/[MITgcm]/MITgcm/model/src/dynamics.F
ViewVC logotype

Diff of /MITgcm/model/src/dynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph | View Patch Patch

revision 1.3 by adcroft, Wed Apr 29 21:31:09 1998 UTC revision 1.16 by cnh, Tue Jun 9 16:34:03 1998 UTC
# Line 2  C $Header$ Line 2  C $Header$
2    
3  #include "CPP_EEOPTIONS.h"  #include "CPP_EEOPTIONS.h"
4    
5        SUBROUTINE DYNAMICS(myThid)        SUBROUTINE DYNAMICS(myTime, myIter, myThid)
6  C     /==========================================================\  C     /==========================================================\
7  C     | SUBROUTINE DYNAMICS                                      |  C     | SUBROUTINE DYNAMICS                                      |
8  C     | o Controlling routine for the explicit part of the model |  C     | o Controlling routine for the explicit part of the model |
# Line 25  C     == Global variables === Line 25  C     == Global variables ===
25  #include "SIZE.h"  #include "SIZE.h"
26  #include "EEPARAMS.h"  #include "EEPARAMS.h"
27  #include "CG2D.h"  #include "CG2D.h"
28    #include "PARAMS.h"
29  #include "DYNVARS.h"  #include "DYNVARS.h"
30    
31  C     == Routine arguments ==  C     == Routine arguments ==
32    C     myTime - Current time in simulation
33    C     myIter - Current iteration number in simulation
34  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
35        INTEGER myThid        INTEGER myThid
36          _RL myTime
37          INTEGER myIter
38    
39  C     == Local variables  C     == Local variables
40  C     xA, yA                 - Per block temporaries holding face areas  C     xA, yA                 - Per block temporaries holding face areas
41  C     uTrans, vTrans, wTrans - Per block temporaries holding flow transport  C     uTrans, vTrans, wTrans - Per block temporaries holding flow transport
42  C                              o uTrans: Zonal transport  C     wVel                     o uTrans: Zonal transport
43  C                              o vTrans: Meridional transport  C                              o vTrans: Meridional transport
44  C                              o wTrans: Vertical transport  C                              o wTrans: Vertical transport
45    C                              o wVel:   Vertical velocity at upper and lower
46    C                                        cell faces.
47  C     maskC,maskUp             o maskC: land/water mask for tracer cells  C     maskC,maskUp             o maskC: land/water mask for tracer cells
48  C                              o maskUp: land/water mask for W points  C                              o maskUp: land/water mask for W points
49  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in  C     aTerm, xTerm, cTerm    - Work arrays for holding separate terms in
# Line 63  C                          into fVerTerm Line 70  C                          into fVerTerm
70        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73          _RL wVel  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
74        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
# Line 79  C                          into fVerTerm Line 87  C                          into fVerTerm
87        _RL pH    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL pH    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
88        _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89        _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90          _RL rhotmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91          _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92          _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93          _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
94          _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
95          _RL K33   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
96          _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97          _RL KappaZT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)
98    
99        INTEGER iMin, iMax        INTEGER iMin, iMax
100        INTEGER jMin, jMax        INTEGER jMin, jMax
101        INTEGER bi, bj        INTEGER bi, bj
102        INTEGER i, j        INTEGER i, j
103        INTEGER k, kM1, kUp, kDown        INTEGER k, kM1, kUp, kDown
104    
105    C---    The algorithm...
106    C
107    C       "Correction Step"
108    C       =================
109    C       Here we update the horizontal velocities with the surface
110    C       pressure such that the resulting flow is either consistent
111    C       with the free-surface evolution or the rigid-lid:
112    C         U[n] = U* + dt x d/dx P
113    C         V[n] = V* + dt x d/dy P
114    C
115    C       "Calculation of Gs"
116    C       ===================
117    C       This is where all the accelerations and tendencies (ie.
118    C       physics, parameterizations etc...) are calculated
119    C         w = sum_z ( div. u[n] )
120    C         rho = rho ( theta[n], salt[n] )
121    C         K31 = K31 ( rho )
122    C         Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... )
123    C         Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... )
124    C         Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... )
125    C         Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... )
126    C
127    C       "Time-stepping" or "Prediction"
128    C       ================================
129    C       The models variables are stepped forward with the appropriate
130    C       time-stepping scheme (currently we use Adams-Bashforth II)
131    C       - For momentum, the result is always *only* a "prediction"
132    C       in that the flow may be divergent and will be "corrected"
133    C       later with a surface pressure gradient.
134    C       - Normally for tracers the result is the new field at time
135    C       level [n+1} *BUT* in the case of implicit diffusion the result
136    C       is also *only* a prediction.
137    C       - We denote "predictors" with an asterisk (*).
138    C         U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
139    C         V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
140    C         theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
141    C         salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
142    C       With implicit diffusion:
143    C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
144    C         salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
145    C         (1 + dt * K * d_zz) theta[n] = theta*
146    C         (1 + dt * K * d_zz) salt[n] = salt*
147    C---
148    
149  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
150  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
151  C     just ensure that all memory references are to valid floating  C     just ensure that all memory references are to valid floating
# Line 92  C     point numbers. This prevents spuri Line 153  C     point numbers. This prevents spuri
153  C     uninitialised but inert locations.  C     uninitialised but inert locations.
154        DO j=1-OLy,sNy+OLy        DO j=1-OLy,sNy+OLy
155         DO i=1-OLx,sNx+OLx         DO i=1-OLx,sNx+OLx
156          xA(i,j)      = 0.*1. _d 37          xA(i,j)      = 0. _d 0
157          yA(i,j)      = 0.*1. _d 37          yA(i,j)      = 0. _d 0
158          uTrans(i,j)  = 0.*1. _d 37          uTrans(i,j)  = 0. _d 0
159          vTrans(i,j)  = 0.*1. _d 37          vTrans(i,j)  = 0. _d 0
160          aTerm(i,j)   = 0.*1. _d 37          aTerm(i,j)   = 0. _d 0
161          xTerm(i,j)   = 0.*1. _d 37          xTerm(i,j)   = 0. _d 0
162          cTerm(i,j)   = 0.*1. _d 37          cTerm(i,j)   = 0. _d 0
163          mTerm(i,j)   = 0.*1. _d 37          mTerm(i,j)   = 0. _d 0
164          pTerm(i,j)   = 0.*1. _d 37          pTerm(i,j)   = 0. _d 0
165          fZon(i,j)    = 0.*1. _d 37          fZon(i,j)    = 0. _d 0
166          fMer(i,j)    = 0.*1. _d 37          fMer(i,j)    = 0. _d 0
167          DO K=1,nZ          DO K=1,nZ
168           pH (i,j,k)  = 0.*1. _d 37           pH (i,j,k)  = 0. _d 0
169             K13(i,j,k) = 0. _d 0
170             K23(i,j,k) = 0. _d 0
171             K33(i,j,k) = 0. _d 0
172             KappaZT(i,j,k) = 0. _d 0
173          ENDDO          ENDDO
174          rhokm1(i,j)    = 0. _d 0          rhokm1(i,j)  = 0. _d 0
175          rhokp1(i,j)    = 0. _d 0          rhokp1(i,j)  = 0. _d 0
176         ENDDO          rhotmp(i,j)  = 0. _d 0
177        ENDDO          maskC (i,j)  = 0. _d 0
 C--   Set up work arrays that need valid initial values  
       DO j=1-OLy,sNy+OLy  
        DO i=1-OLx,sNx+OLx  
         wTrans(i,j)  = 0. _d 0  
         fVerT(i,j,1) = 0. _d 0  
         fVerT(i,j,2) = 0. _d 0  
         fVerS(i,j,1) = 0. _d 0  
         fVerS(i,j,2) = 0. _d 0  
         fVerU(i,j,1) = 0. _d 0  
         fVerU(i,j,2) = 0. _d 0  
         fVerV(i,j,1) = 0. _d 0  
         fVerV(i,j,2) = 0. _d 0  
178         ENDDO         ENDDO
179        ENDDO        ENDDO
180    
181        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
182         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
183    
184  C--   Boundary condition on hydrostatic pressure is pH(z=0)=0  C--     Set up work arrays that need valid initial values
185          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
186           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
187              wTrans(i,j)  = 0. _d 0
188              wVel  (i,j,1) = 0. _d 0
189              wVel  (i,j,2) = 0. _d 0
190              fVerT(i,j,1) = 0. _d 0
191              fVerT(i,j,2) = 0. _d 0
192              fVerS(i,j,1) = 0. _d 0
193              fVerS(i,j,2) = 0. _d 0
194              fVerU(i,j,1) = 0. _d 0
195              fVerU(i,j,2) = 0. _d 0
196              fVerV(i,j,1) = 0. _d 0
197              fVerV(i,j,2) = 0. _d 0
198            pH(i,j,1) = 0. _d 0            pH(i,j,1) = 0. _d 0
199              K13(i,j,1) = 0. _d 0
200              K23(i,j,1) = 0. _d 0
201              K33(i,j,1) = 0. _d 0
202              KapGM(i,j) = 0. _d 0
203           ENDDO           ENDDO
204          ENDDO          ENDDO
205    
# Line 140  C--   Boundary condition on hydrostatic Line 208  C--   Boundary condition on hydrostatic
208          jMin = 1-OLy+1          jMin = 1-OLy+1
209          jMax = sNy+OLy          jMax = sNy+OLy
210    
211  C--     Update fields according to tendency terms  C--     Calculate gradient of surface pressure
212          CALL TIMESTEP(          CALL GRAD_PSURF(
213       I       bi,bj,iMin,iMax,jMin,jMax,myThid)       I       bi,bj,iMin,iMax,jMin,jMax,
214         O       pSurfX,pSurfY,
215         I       myThid)
216    
217    C--     Update fields in top level according to tendency terms
218            CALL CORRECTION_STEP(
219         I       bi,bj,iMin,iMax,jMin,jMax,1,pSurfX,pSurfY,myThid)
220    
221    C--     Density of 1st level (below W(1)) reference to level 1
222            CALL FIND_RHO(
223         I     bi, bj, iMin, iMax, jMin, jMax, 1, 1, eosType,
224         O     rhoKm1,
225         I     myThid )
226    C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0
227            CALL CALC_PH(
228         I      bi,bj,iMin,iMax,jMin,jMax,1,rhoKm1,rhoKm1,
229         U      pH,
230         I      myThid )
231            DO J=jMin,jMax
232             DO I=iMin,iMax
233              rhoKp1(I,J)=rhoKm1(I,J)
234             ENDDO
235            ENDDO
236    
237          DO K=2,Nz          DO K=2,Nz
238  C Density of K-1 level (above W(K)) reference to K level  C--     Update fields in Kth level according to tendency terms
239           CALL FIND_RHO(          CALL CORRECTION_STEP(
240       I      bi, bj, iMin, iMax, jMin, jMax,  K-1, K, 'LINEAR',       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myThid)
241       O      rhoKm1,  C--     Density of K-1 level (above W(K)) reference to K-1 level
242       I      myThid )  copt    CALL FIND_RHO(
243  C Density of K level (below W(K)) reference to K level  copt I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, eosType,
244           CALL FIND_RHO(  copt O     rhoKm1,
245       I      bi, bj, iMin, iMax, jMin, jMax,  K, K, 'LINEAR',  copt I     myThid )
246       O      rhoKp1,  C       rhoKm1=rhoKp1
247       I      myThid )          DO J=jMin,jMax
248             DO I=iMin,iMax
249              rhoKm1(I,J)=rhoKp1(I,J)
250             ENDDO
251            ENDDO
252    C--     Density of K level (below W(K)) reference to K level
253            CALL FIND_RHO(
254         I     bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,
255         O     rhoKp1,
256         I     myThid )
257    C--     Density of K-1 level (above W(K)) reference to K level
258            CALL FIND_RHO(
259         I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K, eosType,
260         O     rhotmp,
261         I     myThid )
262    C--     Calculate iso-neutral slopes for the GM/Redi parameterisation
263            CALL CALC_ISOSLOPES(
264         I            bi, bj, iMin, iMax, jMin, jMax, K,
265         I            rhoKm1, rhoKp1, rhotmp,
266         O            K13, K23, K33, KapGM,
267         I            myThid )
268  C--     Calculate static stability and mix where convectively unstable  C--     Calculate static stability and mix where convectively unstable
269           CALL CONVECT(          CALL CONVECT(
270       I       bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,myThid)       I      bi,bj,iMin,iMax,jMin,jMax,K,rhotmp,rhoKp1,
271  C Density of K-1 level (above W(K)) reference to K-1 level       I      myTime,myIter,myThid)
272           CALL FIND_RHO(  C--     Density of K-1 level (above W(K)) reference to K-1 level
273       I      bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, 'LINEAR',          CALL FIND_RHO(
274       O      rhoKm1,       I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, eosType,
275       I      myThid )       O     rhoKm1,
276         I     myThid )
277    C--     Density of K level (below W(K)) referenced to K level
278            CALL FIND_RHO(
279         I     bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,
280         O     rhoKp1,
281         I     myThid )
282  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0
283           CALL CALC_PH(          CALL CALC_PH(
284       I       bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,       I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,
285       U       pH,       U      pH,
286       I       myThid )       I      myThid )
287    
288          ENDDO ! K          ENDDO ! K
289    
290  C Density of Nz level (bottom level) reference to Nz level  C--     Initial boundary condition on barotropic divergence integral
291           CALL FIND_RHO(          DO j=1-OLy,sNy+OLy
292       I      bi, bj, iMin, iMax, jMin, jMax,  Nz, Nz, 'LINEAR',           DO i=1-OLx,sNx+OLx
293       O      rhoKm1,            cg2d_b(i,j,bi,bj) = 0. _d 0
294       I      myThid )           ENDDO
295  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0          ENDDO
          CALL CALC_PH(  
      I       bi,bj,iMin,iMax,jMin,jMax,Nz+1,rhoKm1,  
      U       pH,  
      I       myThid )  
296    
297          DO K = Nz, 1, -1          DO K = Nz, 1, -1
298           kM1  =max(1,k-1)   ! Points to level above k (=k-1)           kM1  =max(1,k-1)   ! Points to level above k (=k-1)
# Line 193  C--     Integrate hydrostatic balance fo Line 306  C--     Integrate hydrostatic balance fo
306  C--      Get temporary terms used by tendency routines  C--      Get temporary terms used by tendency routines
307           CALL CALC_COMMON_FACTORS (           CALL CALC_COMMON_FACTORS (
308       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
309       O        xA,yA,uTrans,vTrans,wTrans,maskC,maskUp,       O        xA,yA,uTrans,vTrans,wTrans,wVel,maskC,maskUp,
310       I        myThid)       I        myThid)
311    
312  C--      Calculate accelerations in the momentum equations  C--      Calculate the total vertical diffusivity
313           CALL CALC_MOM_RHS(           CALL CALC_DIFFUSIVITY(
314       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,K,
315       I        xA,yA,uTrans,vTrans,wTrans,maskC,       I        maskC,maskUp,KapGM,K33,
316       I        pH,       O        KappaZT,
      U        aTerm,xTerm,cTerm,mTerm,pTerm,  
      U        fZon, fMer, fVerU, fVerV,  
317       I        myThid)       I        myThid)
318    
319    C--      Calculate accelerations in the momentum equations
320             IF ( momStepping ) THEN
321              CALL CALC_MOM_RHS(
322         I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
323         I         xA,yA,uTrans,vTrans,wTrans,wVel,maskC,
324         I         pH,
325         U         aTerm,xTerm,cTerm,mTerm,pTerm,
326         U         fZon, fMer, fVerU, fVerV,
327         I         myThid)
328             ENDIF
329    
330  C--      Calculate active tracer tendencies  C--      Calculate active tracer tendencies
331           CALL CALC_GT(           IF ( tempStepping ) THEN
332       I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,            CALL CALC_GT(
333       I        xA,yA,uTrans,vTrans,wTrans,maskUp,       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
334       U        aTerm,xTerm,fZon,fMer,fVerT,       I         xA,yA,uTrans,vTrans,wTrans,maskUp,
335       I        myThid)       I         K13,K23,KappaZT,KapGM,
336         U         aTerm,xTerm,fZon,fMer,fVerT,
337         I         myThid)
338             ENDIF
339  Cdbg     CALL CALC_GS(  Cdbg     CALL CALC_GS(
340  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
341  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,
342    Cdbg I        K13,K23,K33,KapGM,
343  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,
344  Cdbg I        myThid)  Cdbg I        myThid)
345    
346          ENDDO  C--      Prediction step (step forward all model variables)
347             CALL TIMESTEP(
348         I       bi,bj,iMin,iMax,jMin,jMax,K,
349         I       myThid)
350    
351    C--      Diagnose barotropic divergence of predicted fields
352             CALL DIV_G(
353         I       bi,bj,iMin,iMax,jMin,jMax,K,
354         I       xA,yA,
355         I       myThid)
356    
357            ENDDO ! K
358    
359    C--     Implicit diffusion
360            IF (implicitDiffusion) THEN
361             CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,
362         I                  KappaZT,
363         I                  myThid )
364            ENDIF
365    
366         ENDDO         ENDDO
367        ENDDO        ENDDO
368    
369          write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),
370         &                           maxval(cg2d_x(1:sNx,1:sNy,:,:))
371          write(0,*) 'dynamics: U  ',minval(uVel(1:sNx,1:sNy,:,:,:)),
372         &                           maxval(uVel(1:sNx,1:sNy,:,:,:))
373          write(0,*) 'dynamics: V  ',minval(vVel(1:sNx,1:sNy,:,:,:)),
374         &                           maxval(vVel(1:sNx,1:sNy,:,:,:))
375    cblk  write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
376    cblk &                           maxval(K13(1:sNx,1:sNy,:))
377    cblk  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
378    cblk &                           maxval(K23(1:sNx,1:sNy,:))
379    cblk  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
380    cblk &                           maxval(K33(1:sNx,1:sNy,:))
381          write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),
382         &                           maxval(gT(1:sNx,1:sNy,:,:,:))
383          write(0,*) 'dynamics: T  ',minval(Theta(1:sNx,1:sNy,:,:,:)),
384         &                           maxval(Theta(1:sNx,1:sNy,:,:,:))
385    cblk  write(0,*) 'dynamics: pH ',minval(pH/(Gravity*Rhonil)),
386    cblk &                           maxval(pH/(Gravity*Rhonil))
387    
388        RETURN        RETURN
389        END        END

Legend:
Removed from v.1.3  
changed lines
  Added in v.1.16

  ViewVC Help
Powered by ViewVC 1.1.22