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

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

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

revision 1.26 by heimbach, Mon Mar 6 18:25:49 2006 UTC revision 1.41 by heimbach, Wed Apr 18 19:54:06 2007 UTC
# Line 11  C $Name$ Line 11  C $Name$
11  # ifdef ALLOW_KPP  # ifdef ALLOW_KPP
12  #  include "KPP_OPTIONS.h"  #  include "KPP_OPTIONS.h"
13  # endif  # endif
14    # ifdef ALLOW_SEAICE
15    #  include "SEAICE_OPTIONS.h"
16    # endif
17  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
18    
19  CBOP  CBOP
# Line 19  C     !INTERFACE: Line 22  C     !INTERFACE:
22        SUBROUTINE DO_OCEANIC_PHYS(myTime, myIter, myThid)        SUBROUTINE DO_OCEANIC_PHYS(myTime, myIter, myThid)
23  C     !DESCRIPTION: \bv  C     !DESCRIPTION: \bv
24  C     *==========================================================*  C     *==========================================================*
25  C     | SUBROUTINE DO_OCEANIC_PHYS                                  C     | SUBROUTINE DO_OCEANIC_PHYS
26  C     | o Controlling routine for oceanic physics and  C     | o Controlling routine for oceanic physics and
27  C     |   parameterization  C     |   parameterization
28  C     *==========================================================*  C     *==========================================================*
29  C     | o originally, part of S/R thermodynamics  C     | o originally, part of S/R thermodynamics
# Line 59  C     == Global variables === Line 62  C     == Global variables ===
62  # ifdef EXACT_CONSERV  # ifdef EXACT_CONSERV
63  #  include "SURFACE.h"  #  include "SURFACE.h"
64  # endif  # endif
65    # ifdef ALLOW_EXF
66    #  include "ctrl.h"
67    #  include "EXF_FIELDS.h"
68    #  ifdef ALLOW_BULKFORMULAE
69    #   include "EXF_CONSTANTS.h"
70    #  endif
71    # endif
72    # ifdef ALLOW_SEAICE
73    #  include "SEAICE.h"
74    # endif
75  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
76    
77  C     !INPUT/OUTPUT PARAMETERS:  C     !INPUT/OUTPUT PARAMETERS:
# Line 77  C     iMin, iMax    :: Ranges and sub-bl Line 90  C     iMin, iMax    :: Ranges and sub-bl
90  C     jMin, jMax       are applied.  C     jMin, jMax       are applied.
91  C     bi, bj        :: tile indices  C     bi, bj        :: tile indices
92  C     i,j,k         :: loop indices  C     i,j,k         :: loop indices
93          _RL rhokp1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96        _RL sigmaX  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL sigmaX  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
# Line 100  C--   dummy statement to end declaration Line 114  C--   dummy statement to end declaration
114  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
115    
116  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
117        IF ( debugLevel .GE. debLevB )        IF ( debugLevel .GE. debLevB )
118       &    CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid)       &     CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid)
119  #endif  #endif
120    
121        doDiagsRho = 0        doDiagsRho = 0
122  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
123        IF ( useDiagnostics .AND. fluidIsWater ) THEN        IF ( useDiagnostics .AND. fluidIsWater ) THEN
# Line 116  C--   dummy statement to end declaration Line 130  C--   dummy statement to end declaration
130        ENDIF        ENDIF
131  #endif /* ALLOW_DIAGNOSTICS */  #endif /* ALLOW_DIAGNOSTICS */
132    
133  #ifdef ALLOW_THSICE  #ifdef ALLOW_SEAICE
134    C--   Call sea ice model to compute forcing/external data fields.  In
135    C     addition to computing prognostic sea-ice variables and diagnosing the
136    C     forcing/external data fields that drive the ocean model, SEAICE_MODEL
137    C     also sets theta to the freezing point under sea-ice.  The implied
138    C     surface heat flux is then stored in variable surfaceTendencyTice,
139    C     which is needed by KPP package (kpp_calc.F and kpp_transport_t.F)
140    C     to diagnose surface buoyancy fluxes and for the non-local transport
141    C     term.  Because this call precedes model thermodynamics, temperature
142    C     under sea-ice may not be "exactly" at the freezing point by the time
143    C     theta is dumped or time-averaged.
144          IF ( useSEAICE ) THEN
145    #ifdef ALLOW_AUTODIFF_TAMC
146    CADJ STORE atemp,aqh,precip    = comlev1, key = ikey_dynamics
147    CADJ STORE swdown,lwdown       = comlev1, key = ikey_dynamics
148    CADJ STORE theta               = comlev1, key = ikey_dynamics
149    cph# ifdef EXF_READ_EVAP
150    CADJ STORE evap                = comlev1, key = ikey_dynamics
151    cph# endif
152    cph# ifdef SEAICE_ALLOW_DYNAMICS
153    CADJ STORE uvel,vvel           = comlev1, key = ikey_dynamics
154    cph# endif
155    # ifdef SEAICE_CGRID
156    CADJ STORE fu, fv              = comlev1, key = ikey_dynamics
157    CADJ STORE seaicemasku         = comlev1, key = ikey_dynamics
158    CADJ STORE seaicemaskv         = comlev1, key = ikey_dynamics
159    #  ifdef ATMOSPHERIC_LOADING
160    CADJ STORE siceload            = comlev1, key = ikey_dynamics
161    #  endif
162    # endif
163    #endif
164    #ifdef ALLOW_DEBUG
165            IF ( debugLevel .GE. debLevB )
166         &    CALL DEBUG_CALL('SEAICE_MODEL',myThid)
167    #endif
168            CALL TIMER_START('SEAICE_MODEL    [DO_OCEANIC_PHYS]', myThid)
169            CALL SEAICE_MODEL( myTime, myIter, myThid )
170            CALL TIMER_STOP ('SEAICE_MODEL    [DO_OCEANIC_PHYS]', myThid)
171    #ifdef ALLOW_COST_ICE
172            CALL COST_ICE_TEST ( myTime, myIter, myThid )
173    #endif
174          ENDIF
175    #endif /* ALLOW_SEAICE */
176    
177    #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
178        IF ( useThSIce .AND. fluidIsWater ) THEN        IF ( useThSIce .AND. fluidIsWater ) THEN
179  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
180          IF ( debugLevel .GE. debLevB )          IF ( debugLevel .GE. debLevB )
# Line 156  CADJ STORE theta = comlev1, key = ikey_d Line 214  CADJ STORE theta = comlev1, key = ikey_d
214          CALL FREEZE_SURFACE(  myTime, myIter, myThid )          CALL FREEZE_SURFACE(  myTime, myIter, myThid )
215        ENDIF        ENDIF
216    
217  #ifdef COMPONENT_MODULE  #ifdef ALLOW_OCN_COMPON_INTERF
 # ifndef ALLOW_AIM  
218  C--    Apply imported data (from coupled interface) to forcing fields  C--    Apply imported data (from coupled interface) to forcing fields
219  C jmc: do not know precisely where to put this call (bf or af thSIce ?)  C jmc: do not know precisely where to put this call (bf or af thSIce ?)
220         IF ( useCoupler ) THEN        IF ( useCoupler ) THEN
221           CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )           CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )
222         ENDIF        ENDIF
223  # endif  #endif /* ALLOW_OCN_COMPON_INTERF */
 #endif /* COMPONENT_MODULE */  
224    
225  #ifdef ALLOW_BALANCE_FLUXES  #ifdef ALLOW_BALANCE_FLUXES
226  C     balance fluxes  C     balance fluxes
227         IF ( balanceEmPmR )        IF ( balanceEmPmR )
228       &        CALL REMOVE_MEAN_RS( 1, EmPmR, maskH, maskH, rA, drF,       &        CALL REMOVE_MEAN_RS( 1, EmPmR, maskH, maskH, rA, drF,
229       &        'EmPmR', myTime, myThid )       &        'EmPmR', myTime, myThid )
230         IF ( balanceQnet )        IF ( balanceQnet )
231       &        CALL REMOVE_MEAN_RS( 1, Qnet,  maskH, maskH, rA, drF,       &        CALL REMOVE_MEAN_RS( 1, Qnet,  maskH, maskH, rA, drF,
232       &        'Qnet ', myTime, myThid )       &        'Qnet ', myTime, myThid )
233  #endif /* ALLOW_BALANCE_FLUXES */  #endif /* ALLOW_BALANCE_FLUXES */
# Line 198  CHPF$ INDEPENDENT Line 254  CHPF$ INDEPENDENT
254            itdkey = (act1 + 1) + act2*max1            itdkey = (act1 + 1) + act2*max1
255       &                      + act3*max1*max2       &                      + act3*max1*max2
256       &                      + act4*max1*max2*max3       &                      + act4*max1*max2*max3
257    #else  /* ALLOW_AUTODIFF_TAMC */
258    C     if fluid is not water, by-pass find_rho, gmredi, surfaceForcing
259    C     and all vertical mixing schemes, but keep OBCS_CALC
260            IF ( fluidIsWater ) THEN
261  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
262    
263  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
# Line 210  C     uninitialised but inert locations. Line 270  C     uninitialised but inert locations.
270           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
271            rhok   (i,j)   = 0. _d 0            rhok   (i,j)   = 0. _d 0
272            rhoKM1 (i,j)   = 0. _d 0            rhoKM1 (i,j)   = 0. _d 0
273              rhoKP1 (i,j)   = 0. _d 0
274           ENDDO           ENDDO
275          ENDDO          ENDDO
276    
# Line 266  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_ Line 327  CADJ STORE vvel (:,:,:,bi,bj) = comlev1_
327  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
328    
329  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
330          IF ( debugLevel .GE. debLevB )          IF ( debugLevel .GE. debLevB )
331       &    CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)       &    CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
332  #endif  #endif
333    
# Line 283  cph Needed for rhok, rhokm1, in the case Line 344  cph Needed for rhok, rhokm1, in the case
344    
345  C--       Calculate gradients of potential density for isoneutral  C--       Calculate gradients of potential density for isoneutral
346  C         slope terms (e.g. GM/Redi tensor or IVDC diffusivity)  C         slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
 c         IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN  
347            IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)            IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
348       &                   .OR. doDiagsRho.GE.1 ) THEN       &                   .OR. doDiagsRho.GE.1 ) THEN
349  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
350              IF ( debugLevel .GE. debLevB )              IF ( debugLevel .GE. debLevB )
351       &       CALL DEBUG_CALL('FIND_RHO',myThid)       &       CALL DEBUG_CALL('FIND_RHO',myThid)
352  #endif  #endif
353  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
# Line 312  CADJ STORE salt (:,:,k-1,bi,bj) = comlev Line 372  CADJ STORE salt (:,:,k-1,bi,bj) = comlev
372       I        myThid )       I        myThid )
373              ENDIF              ENDIF
374  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
375              IF ( debugLevel .GE. debLevB )              IF ( debugLevel .GE. debLevB )
376       &       CALL DEBUG_CALL('GRAD_SIGMA',myThid)       &       CALL DEBUG_CALL('GRAD_SIGMA',myThid)
377  #endif  #endif
378    cph Avoid variable aliasing for adjoint !!!
379                DO j=jMin,jMax
380                 DO i=iMin,iMax
381                  rhoKP1(i,j) = rhoK(i,j)
382                 ENDDO
383                ENDDO
384              CALL GRAD_SIGMA(              CALL GRAD_SIGMA(
385       I             bi, bj, iMin, iMax, jMin, jMax, k,       I             bi, bj, iMin, iMax, jMin, jMax, k,
386       I             rhoK, rhoKm1, rhoK,       I             rhoK, rhoKm1, rhoKp1,
387       O             sigmaX, sigmaY, sigmaR,       O             sigmaX, sigmaY, sigmaR,
388       I             myThid )       I             myThid )
389            ENDIF            ENDIF
# Line 332  C--       Implicit Vertical Diffusion fo Line 398  C--       Implicit Vertical Diffusion fo
398  c ==> should use sigmaR !!!  c ==> should use sigmaR !!!
399            IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN            IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
400  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
401              IF ( debugLevel .GE. debLevB )              IF ( debugLevel .GE. debLevB )
402       &       CALL DEBUG_CALL('CALC_IVDC',myThid)       &       CALL DEBUG_CALL('CALC_IVDC',myThid)
403  #endif  #endif
404              CALL CALC_IVDC(              CALL CALC_IVDC(
# Line 353  C--     end of diagnostic k loop (Nr:1) Line 419  C--     end of diagnostic k loop (Nr:1)
419          ENDDO          ENDDO
420    
421  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
 c       IF ( useDiagnostics .AND.  
 c    &       (useGMRedi .OR. ivdc_kappa.NE.0.) ) THEN  
422          IF ( doDiagsRho.GE.1 ) THEN          IF ( doDiagsRho.GE.1 ) THEN
423            CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR  ', 0, Nr,            CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR  ', 0, Nr,
424       &         2, bi, bj, myThid)       &         2, bi, bj, myThid)
425          ENDIF          ENDIF
426  #endif  #endif
427    
 #ifdef  ALLOW_OBCS  
 C--     Calculate future values on open boundaries  
         IF (useOBCS) THEN  
 #ifdef ALLOW_DEBUG  
           IF ( debugLevel .GE. debLevB )  
      &     CALL DEBUG_CALL('OBCS_CALC',myThid)  
 #endif  
           CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1,  
      I            uVel, vVel, wVel, theta, salt,  
      I            myThid )  
         ENDIF  
 #endif  /* ALLOW_OBCS */  
   
 #ifndef ALLOW_AUTODIFF_TAMC  
         IF ( fluidIsWater ) THEN  
 #endif  
428  C--     Determines forcing terms based on external fields  C--     Determines forcing terms based on external fields
429  C       relaxation terms, etc.  C       relaxation terms, etc.
430  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
431          IF ( debugLevel .GE. debLevB )          IF ( debugLevel .GE. debLevB )
432       &    CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)       &    CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
433  #endif  #endif
434  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
# Line 390  CADJ &     = comlev1_bibj, key=itdkey, b Line 438  CADJ &     = comlev1_bibj, key=itdkey, b
438  CADJ STORE PmEpR(:,:,bi,bj)  CADJ STORE PmEpR(:,:,bi,bj)
439  CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte  CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
440  # endif  # endif
441    # ifdef NONLIN_FRSURF
442    CADJ STORE hFac_surfC(:,:,bi,bj)
443    CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
444    CADJ STORE recip_hFacC(:,:,:,bi,bj)
445    CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
446    # endif
447  #endif  #endif
448           CALL EXTERNAL_FORCING_SURF(          CALL EXTERNAL_FORCING_SURF(
449       I             bi, bj, iMin, iMax, jMin, jMax,       I             bi, bj, iMin, iMax, jMin, jMax,
450       I             myTime, myIter, myThid )       I             myTime, myIter, myThid )
451  #ifndef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
452          ENDIF  # ifdef EXACT_CONSERV
453    cph-test
454    cphCADJ STORE PmEpR(:,:,bi,bj)
455    cphCADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
456    # endif
457  #endif  #endif
458    
459  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
# Line 429  CADJ STORE sigmaR(:,:,:)        = comlev Line 487  CADJ STORE sigmaR(:,:,:)        = comlev
487  C--     Calculate iso-neutral slopes for the GM/Redi parameterisation  C--     Calculate iso-neutral slopes for the GM/Redi parameterisation
488          IF (useGMRedi) THEN          IF (useGMRedi) THEN
489  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
490            IF ( debugLevel .GE. debLevB )            IF ( debugLevel .GE. debLevB )
491       &     CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)       &     CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
492  #endif  #endif
493            CALL GMREDI_CALC_TENSOR(            CALL GMREDI_CALC_TENSOR(
# Line 451  C--     Calculate iso-neutral slopes for Line 509  C--     Calculate iso-neutral slopes for
509  C--     Compute KPP mixing coefficients  C--     Compute KPP mixing coefficients
510          IF (useKPP) THEN          IF (useKPP) THEN
511  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
512            IF ( debugLevel .GE. debLevB )            IF ( debugLevel .GE. debLevB )
513       &     CALL DEBUG_CALL('KPP_CALC',myThid)       &     CALL DEBUG_CALL('KPP_CALC',myThid)
514  #endif  #endif
515            CALL KPP_CALC(            CALL KPP_CALC(
# Line 469  C--     Compute KPP mixing coefficients Line 527  C--     Compute KPP mixing coefficients
527  C--     Compute PP81 mixing coefficients  C--     Compute PP81 mixing coefficients
528          IF (usePP81) THEN          IF (usePP81) THEN
529  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
530            IF ( debugLevel .GE. debLevB )            IF ( debugLevel .GE. debLevB )
531       &     CALL DEBUG_CALL('PP81_CALC',myThid)       &     CALL DEBUG_CALL('PP81_CALC',myThid)
532  #endif  #endif
533            CALL PP81_CALC(            CALL PP81_CALC(
# Line 481  C--     Compute PP81 mixing coefficients Line 539  C--     Compute PP81 mixing coefficients
539  C--     Compute MY82 mixing coefficients  C--     Compute MY82 mixing coefficients
540          IF (useMY82) THEN          IF (useMY82) THEN
541  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
542            IF ( debugLevel .GE. debLevB )            IF ( debugLevel .GE. debLevB )
543       &     CALL DEBUG_CALL('MY82_CALC',myThid)       &     CALL DEBUG_CALL('MY82_CALC',myThid)
544  #endif  #endif
545            CALL MY82_CALC(            CALL MY82_CALC(
# Line 493  C--     Compute MY82 mixing coefficients Line 551  C--     Compute MY82 mixing coefficients
551  C--     Compute GGL90 mixing coefficients  C--     Compute GGL90 mixing coefficients
552          IF (useGGL90) THEN          IF (useGGL90) THEN
553  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
554            IF ( debugLevel .GE. debLevB )            IF ( debugLevel .GE. debLevB )
555       &     CALL DEBUG_CALL('GGL90_CALC',myThid)       &     CALL DEBUG_CALL('GGL90_CALC',myThid)
556  #endif  #endif
557            CALL GGL90_CALC(            CALL GGL90_CALC(
# Line 502  C--     Compute GGL90 mixing coefficient Line 560  C--     Compute GGL90 mixing coefficient
560  #endif /* ALLOW_GGL90 */  #endif /* ALLOW_GGL90 */
561    
562  #ifdef ALLOW_TIMEAVE  #ifdef ALLOW_TIMEAVE
563          IF ( taveFreq.GT. 0. _d 0 .AND. fluidIsWater ) THEN          IF ( taveFreq.GT. 0. _d 0 ) THEN
564            CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)            CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
565          ENDIF          ENDIF
566          IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN          IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
# Line 511  C--     Compute GGL90 mixing coefficient Line 569  C--     Compute GGL90 mixing coefficient
569          ENDIF          ENDIF
570  #endif /* ALLOW_TIMEAVE */  #endif /* ALLOW_TIMEAVE */
571    
572    #ifndef ALLOW_AUTODIFF_TAMC
573    C---  if fluid Is Water: end
574            ENDIF
575    #endif
576    
577    #ifdef  ALLOW_OBCS
578    C--     Calculate future values on open boundaries
579            IF (useOBCS) THEN
580    #ifdef ALLOW_DEBUG
581              IF ( debugLevel .GE. debLevB )
582         &     CALL DEBUG_CALL('OBCS_CALC',myThid)
583    #endif
584              CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1,
585         I            uVel, vVel, wVel, theta, salt,
586         I            myThid )
587            ENDIF
588    #endif  /* ALLOW_OBCS */
589    
590  C--   end bi,bj loops.  C--   end bi,bj loops.
591         ENDDO         ENDDO
592        ENDDO        ENDDO
# Line 526  C--   end bi,bj loops. Line 602  C--   end bi,bj loops.
602  #endif  #endif
603    
604  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
605           IF ( debugLevel .GE. debLevB )        IF ( debugLevel .GE. debLevB )
606       &    CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid)       &     CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid)
607  #endif  #endif
608    
609        RETURN        RETURN

Legend:
Removed from v.1.26  
changed lines
  Added in v.1.41

  ViewVC Help
Powered by ViewVC 1.1.22