Changeset 5224


Ignore:
Timestamp:
Apr 21, 2008, 5:16:38 PM (17 years ago)
Author:
ole
Message:

Work done during Water Down Under 2008.
Tried alternative CFL conditions as suggested by Ted Rigby.
They are kept commented out.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • anuga_core/source/anuga/shallow_water/shallow_water_ext.c

    r5175 r5224  
    231231  _rotate(q_right_rotated, n1, n2);
    232232
    233   z = (z_left+z_right)/2; // Average elevation values
     233  z = (z_left+z_right)/2; // Average elevation values.
     234                          // Even though this will nominally allow for discontinuities
     235                          // in the elevation data, there is currently no numerical
     236                          // support for this so results may be strange near jumps in the bed.
    234237
    235238  // Compute speeds in x-direction
     
    17691772      n = neighbours[ki];
    17701773      if (n < 0) {
    1771         m = -n-1; // Convert negative flag to index
     1774        // Neighbour is a boundary condition
     1775        m = -n-1; // Convert negative flag to boundary index
    17721776       
    17731777        qr[0] = stage_boundary_values[m];
     
    17761780        zr = zl; // Extend bed elevation to boundary
    17771781      } else {
     1782        // Neighbour is a real element
    17781783        m = neighbour_edges[ki];
    17791784        nm = n*3+m; // Linear index (triangle n, edge m)
     
    18471852          if (n>=0)
    18481853            timestep = min(timestep, radii[n]/max_speed);
     1854         
     1855          // Ted Rigby's suggestion
     1856          //if (n>=0) {             
     1857          //  timestep = min(timestep, 0.8*(radii[k]+radii[n])/max_speed);
     1858          //} else {
     1859          //  timestep = min(timestep, radii[k]/max_speed);
     1860          // }
     1861         
     1862         
     1863          // Ole's modification
     1864          //if (n>=0) {             
     1865          //  timestep = min(timestep, max(radii[k],radii[n])/max_speed);
     1866          //} else {
     1867          //  timestep = min(timestep, radii[k]/max_speed);
     1868          //}
     1869         
     1870         
    18491871        }
    18501872      }
    18511873     
    1852     } // End edge i
     1874    } // End edge i (and neighbour n)
    18531875   
    18541876   
     
    20102032    in domain.
    20112033
    2012     Compute total flux for each conserved quantity using "flux_function_central"
     2034    THIS IS AN EXPERIMENTAL FUNCTION - NORMALLY flux_function_central IS USED.
     2035   
     2036    Compute total flux for each conserved quantity using "flux_function_kinetic"
    20132037
    20142038    Fluxes across each edge are scaled by edgelengths and summed up
Note: See TracChangeset for help on using the changeset viewer.