Changeset 5831


Ignore:
Timestamp:
Oct 10, 2008, 4:50:04 PM (17 years ago)
Author:
steve
Message:

Added cfl condition to compute fluxes in comp_flux_ext.c

File:
1 edited

Legend:

Unmodified
Added
Removed
  • anuga_work/development/anuga_1d/comp_flux_ext.c

    r5727 r5831  
    131131       
    132132// Computational function for flux computation
    133 double _compute_fluxes_ext(double timestep,
     133double _compute_fluxes_ext(
     134                           double cfl,
     135                           double timestep,
    134136                           double epsilon,
    135137                           double g,
     
    152154                double zl, zr, max_speed, normal;
    153155                int k, i, ki, n, m, nm=0;
     156               
    154157               
    155158                for (k=0; k<number_of_elements; k++) {
     
    187190                        // Original CFL calculation
    188191                                       
    189                                         timestep = min(timestep, 0.5*areas[k]/max_speed); //Here, CFL=1.0 is assumed. ?????????????????????????????????????????????
     192                                        timestep = min(timestep, 0.5*cfl*areas[k]/max_speed); //Here, CFL=1.0 is assumed. ?????????????????????????????????????????????
    190193                                        if (n>=0) {
    191                                                 timestep = min(timestep, 0.5*areas[n]/max_speed); //Here, CFL=1.0 is assumed. ?????????????????????????????????????????????
     194                                                timestep = min(timestep, 0.5*cfl*areas[n]/max_speed); //Here, CFL=1.0 is assumed. ?????????????????????????????????????????????
    192195                                        }
    193196                                }
     
    230233        *max_speed_array;
    231234   
    232   double timestep, epsilon, g, h0;
     235  double timestep, epsilon, g, h0, cfl;
    233236  int number_of_elements;
    234237   
    235238  // Convert Python arguments to C
    236   if (!PyArg_ParseTuple(args, "ddddOOOOOOOOOOOiO",
     239  if (!PyArg_ParseTuple(args, "dddddOOOOOOOOOOOiO",
     240                        &cfl,
    237241                        &timestep,
    238242                        &epsilon,
     
    259263  // Call underlying flux computation routine and update
    260264  // the explicit update arrays
    261   timestep = _compute_fluxes_ext(timestep,
     265  timestep = _compute_fluxes_ext(
     266                 cfl,
     267                                 timestep,
    262268                                 epsilon,
    263269                                 g,
     
    304310        *max_speed_array;
    305311   
    306   double timestep, epsilon, g, h0;
     312  double timestep, epsilon, g, h0, cfl;
    307313  int number_of_elements;
    308314
     
    323329    g                 = get_python_double(domain,"g");
    324330    h0                = get_python_double(domain,"h0");
     331        cfl                = get_python_double(domain,"cfl");
    325332 
    326333   
     
    350357    // Call underlying flux computation routine and update
    351358    // the explicit update arrays
    352     timestep = _compute_fluxes_ext(timestep,
     359    timestep = _compute_fluxes_ext(
     360                       cfl,
     361                                   timestep,
    353362                                   epsilon,
    354363                                   g,
Note: See TracChangeset for help on using the changeset viewer.