Changeset 5238


Ignore:
Timestamp:
Apr 24, 2008, 10:41:14 AM (16 years ago)
Author:
ole
Message:

Cleanup of C-code and coding style.

File:
1 edited

Legend:

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

    r5224 r5238  
    77// or use python compile.py
    88//
    9 // See the module shallow_water.py
     9// See the module shallow_water_domain.py for more documentation on
     10// how to use this module
    1011//
    1112//
     
    6768      else
    6869        *qmax=dq0;
    69       if ((*qmin=dq0+dq2)<0)
    70         ;// qmin is already set to correct value
    71       else
    72         *qmin=0.0;
     70       
     71      *qmin=dq0+dq2;
     72      if (*qmin>=0.0) *qmin = 0.0;
    7373    }
    7474    else{// dq1<dq2
     
    7777      else
    7878        *qmax=dq0;
    79       if ((*qmin=dq0+dq1)<0)
    80         ;// qmin is the correct value
    81       else
    82         *qmin=0.0;
     79       
     80      *qmin=dq0+dq1;   
     81      if (*qmin>=0.0) *qmin=0.0;
    8382    }
    8483  }
     
    8988      else
    9089        *qmin=dq0;
    91       if ((*qmax=dq0+dq2)>0.0)
    92         ;// qmax is already set to the correct value
    93       else
    94         *qmax=0.0;
     90       
     91      *qmax=dq0+dq2;   
     92      if (*qmax<=0.0) *qmax=0.0;
    9593    }
    9694    else{// dq1>dq2
     
    9997      else
    10098        *qmin=dq0;
    101       if ((*qmax=dq0+dq1)>0.0)
    102         ;// qmax is already set to the correct value
    103       else
    104         *qmax=0.0;
     99       
     100      *qmax=dq0+dq1;
     101      if (*qmax<=0.0) *qmax=0.0;
    105102    }
    106103  }
     
    18491846      if (tri_full_flag[k] == 1) {
    18501847        if (max_speed > epsilon) {
     1848       
     1849          // Original CFL calculation
    18511850          timestep = min(timestep, radii[k]/max_speed);
    18521851          if (n>=0)
    18531852            timestep = min(timestep, radii[n]/max_speed);
    18541853         
    1855           // Ted Rigby's suggestion
     1854          // Ted Rigby's suggested less conservative version
    18561855          //if (n>=0) {             
    1857           //  timestep = min(timestep, 0.8*(radii[k]+radii[n])/max_speed);
     1856          //  timestep = min(timestep, (radii[k]+radii[n])/max_speed);
    18581857          //} else {
    18591858          //  timestep = min(timestep, radii[k]/max_speed);
    18601859          // }
    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          
    18711860        }
    18721861      }
Note: See TracChangeset for help on using the changeset viewer.