/[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.41 by heimbach, Wed Apr 18 19:54:06 2007 UTC revision 1.66 by gforget, Fri May 30 21:10:31 2008 UTC
# Line 49  C     == Global variables === Line 49  C     == Global variables ===
49  # include "tamc.h"  # include "tamc.h"
50  # include "tamc_keys.h"  # include "tamc_keys.h"
51  # include "FFIELDS.h"  # include "FFIELDS.h"
52    # include "SURFACE.h"
53  # include "EOS.h"  # include "EOS.h"
54  # ifdef ALLOW_KPP  # ifdef ALLOW_KPP
55  #  include "KPP.h"  #  include "KPP.h"
# Line 59  C     == Global variables === Line 60  C     == Global variables ===
60  # ifdef ALLOW_EBM  # ifdef ALLOW_EBM
61  #  include "EBM.h"  #  include "EBM.h"
62  # endif  # endif
 # ifdef EXACT_CONSERV  
 #  include "SURFACE.h"  
 # endif  
63  # ifdef ALLOW_EXF  # ifdef ALLOW_EXF
64  #  include "ctrl.h"  #  include "ctrl.h"
65  #  include "EXF_FIELDS.h"  #  include "EXF_FIELDS.h"
# Line 85  C     myThid :: Thread number for this i Line 83  C     myThid :: Thread number for this i
83    
84  C     !LOCAL VARIABLES:  C     !LOCAL VARIABLES:
85  C     == Local variables  C     == Local variables
86  C     rhoK, rhoKM1  :: Density at current level, and level above  C     rhoK, rhoKm1  :: Density at current level, and level above
87  C     iMin, iMax    :: Ranges and sub-block indices on which calculations  C     iMin, iMax    :: Ranges and sub-block indices on which calculations
88  C     jMin, jMax       are applied.  C     jMin, jMax       are applied.
89  C     bi, bj        :: tile indices  C     bi, bj        :: tile indices
90  C     i,j,k         :: loop indices  C     i,j,k         :: loop indices
91        _RL rhokp1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhoKp1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92        _RL rhokm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhoKm1  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93        _RL rhok    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhoK    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94        _RL sigmaX  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL sigmaX  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
95        _RL sigmaY  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL sigmaY  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
96        _RL sigmaR  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL sigmaR  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
# Line 121  C--   dummy statement to end declaration Line 119  C--   dummy statement to end declaration
119        doDiagsRho = 0        doDiagsRho = 0
120  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
121        IF ( useDiagnostics .AND. fluidIsWater ) THEN        IF ( useDiagnostics .AND. fluidIsWater ) THEN
         IF ( DIAGNOSTICS_IS_ON('DRHODR  ',myThid) ) doDiagsRho = 1  
122          IF ( DIAGNOSTICS_IS_ON('RHOANOSQ',myThid) .OR.          IF ( DIAGNOSTICS_IS_ON('RHOANOSQ',myThid) .OR.
123       &       DIAGNOSTICS_IS_ON('URHOMASS',myThid) .OR.       &       DIAGNOSTICS_IS_ON('URHOMASS',myThid) .OR.
124       &       DIAGNOSTICS_IS_ON('VRHOMASS',myThid) .OR.       &       DIAGNOSTICS_IS_ON('VRHOMASS',myThid) .OR.
125       &       DIAGNOSTICS_IS_ON('WRHOMASS',myThid) .OR.       &       DIAGNOSTICS_IS_ON('WRHOMASS',myThid) .OR.
126       &       DIAGNOSTICS_IS_ON('WRHOMASS',myThid) ) doDiagsRho = 2       &       DIAGNOSTICS_IS_ON('WRHOMASS',myThid) ) doDiagsRho = 2
127            IF ( doDiagsRho.EQ.0 .AND.
128         &       DIAGNOSTICS_IS_ON('MXLDEPTH',myThid) ) doDiagsRho = 1
129            IF ( doDiagsRho.EQ.0 .AND.
130         &       DIAGNOSTICS_IS_ON('DRHODR  ',myThid) ) doDiagsRho = 1
131        ENDIF        ENDIF
132  #endif /* ALLOW_DIAGNOSTICS */  #endif /* ALLOW_DIAGNOSTICS */
133    
134  #ifdef ALLOW_SEAICE  #ifdef ALLOW_SEAICE
 C--   Call sea ice model to compute forcing/external data fields.  In  
 C     addition to computing prognostic sea-ice variables and diagnosing the  
 C     forcing/external data fields that drive the ocean model, SEAICE_MODEL  
 C     also sets theta to the freezing point under sea-ice.  The implied  
 C     surface heat flux is then stored in variable surfaceTendencyTice,  
 C     which is needed by KPP package (kpp_calc.F and kpp_transport_t.F)  
 C     to diagnose surface buoyancy fluxes and for the non-local transport  
 C     term.  Because this call precedes model thermodynamics, temperature  
 C     under sea-ice may not be "exactly" at the freezing point by the time  
 C     theta is dumped or time-averaged.  
135        IF ( useSEAICE ) THEN        IF ( useSEAICE ) THEN
136  #ifdef ALLOW_AUTODIFF_TAMC  # ifdef ALLOW_AUTODIFF_TAMC
137    cph-adj-test(
138    CADJ STORE area,empmr,qsw,theta   = comlev1, key = ikey_dynamics
139    cph-adj-test)
140  CADJ STORE atemp,aqh,precip    = comlev1, key = ikey_dynamics  CADJ STORE atemp,aqh,precip    = comlev1, key = ikey_dynamics
141  CADJ STORE swdown,lwdown       = comlev1, key = ikey_dynamics  CADJ STORE swdown,lwdown       = comlev1, key = ikey_dynamics
 CADJ STORE theta               = comlev1, key = ikey_dynamics  
142  cph# ifdef EXF_READ_EVAP  cph# ifdef EXF_READ_EVAP
143  CADJ STORE evap                = comlev1, key = ikey_dynamics  CADJ STORE evap                = comlev1, key = ikey_dynamics
144  cph# endif  cph# endif
 cph# ifdef SEAICE_ALLOW_DYNAMICS  
145  CADJ STORE uvel,vvel           = comlev1, key = ikey_dynamics  CADJ STORE uvel,vvel           = comlev1, key = ikey_dynamics
146  cph# endif  #  ifdef SEAICE_ALLOW_DYNAMICS
147  # ifdef SEAICE_CGRID  CADJ STORE uice                = comlev1, key = ikey_dynamics
148  CADJ STORE fu, fv              = comlev1, key = ikey_dynamics  CADJ STORE vice                = comlev1, key = ikey_dynamics
149  CADJ STORE seaicemasku         = comlev1, key = ikey_dynamics  #   ifdef SEAICE_ALLOW_EVP
150  CADJ STORE seaicemaskv         = comlev1, key = ikey_dynamics  CADJ STORE seaice_sigma1       = comlev1, key = ikey_dynamics
151    CADJ STORE seaice_sigma2       = comlev1, key = ikey_dynamics
152    CADJ STORE seaice_sigma12      = comlev1, key = ikey_dynamics
153    #   endif
154    #  endif
155    #  ifdef SEAICE_SALINITY
156    CADJ STORE salt                = comlev1, key = ikey_dynamics
157    #  endif
158  #  ifdef ATMOSPHERIC_LOADING  #  ifdef ATMOSPHERIC_LOADING
159    CADJ STORE pload               = comlev1, key = ikey_dynamics
160  CADJ STORE siceload            = comlev1, key = ikey_dynamics  CADJ STORE siceload            = comlev1, key = ikey_dynamics
161  #  endif  #  endif
162    #  ifdef NONLIN_FRSURF
163    CADJ STORE recip_hfacc         = comlev1, key = ikey_dynamics
164    #  endif
165  # endif  # endif
166  #endif  # ifdef ALLOW_DEBUG
 #ifdef ALLOW_DEBUG  
167          IF ( debugLevel .GE. debLevB )          IF ( debugLevel .GE. debLevB )
168       &    CALL DEBUG_CALL('SEAICE_MODEL',myThid)       &    CALL DEBUG_CALL('SEAICE_MODEL',myThid)
169  #endif  # endif
170          CALL TIMER_START('SEAICE_MODEL    [DO_OCEANIC_PHYS]', myThid)          CALL TIMER_START('SEAICE_MODEL    [DO_OCEANIC_PHYS]', myThid)
171          CALL SEAICE_MODEL( myTime, myIter, myThid )          CALL SEAICE_MODEL( myTime, myIter, myThid )
172          CALL TIMER_STOP ('SEAICE_MODEL    [DO_OCEANIC_PHYS]', myThid)          CALL TIMER_STOP ('SEAICE_MODEL    [DO_OCEANIC_PHYS]', myThid)
173  #ifdef ALLOW_COST_ICE  # ifdef ALLOW_COST
174          CALL COST_ICE_TEST ( myTime, myIter, myThid )          CALL SEAICE_COST_SENSI ( myTime, myIter, myThid )
175  #endif  # endif
176        ENDIF        ENDIF
177  #endif /* ALLOW_SEAICE */  #endif /* ALLOW_SEAICE */
178    
179    #ifdef ALLOW_AUTODIFF_TAMC
180    CADJ STORE sst, sss           = comlev1, key = ikey_dynamics
181    CADJ STORE qsw                = comlev1, key = ikey_dynamics
182    # ifdef ALLOW_SEAICE
183    CADJ STORE area               = comlev1, key = ikey_dynamics
184    # endif
185    #endif
186    
187  #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)  #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
188        IF ( useThSIce .AND. fluidIsWater ) THEN        IF ( useThSIce .AND. fluidIsWater ) THEN
189  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
# Line 194  C       and modify forcing terms includi Line 204  C       and modify forcing terms includi
204          IF ( debugLevel .GE. debLevB )          IF ( debugLevel .GE. debLevB )
205       &    CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid)       &    CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid)
206  #endif  #endif
207  C     compute temperature and (virtual) salt flux at the  C     compute temperature and (virtual) salt flux at the
208  C     shelf-ice ocean interface  C     shelf-ice ocean interface
209         CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',         CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
210       &       myThid)       &       myThid)
# Line 255  CHPF$ INDEPENDENT Line 265  CHPF$ INDEPENDENT
265       &                      + act3*max1*max2       &                      + act3*max1*max2
266       &                      + act4*max1*max2*max3       &                      + act4*max1*max2*max3
267  #else  /* ALLOW_AUTODIFF_TAMC */  #else  /* ALLOW_AUTODIFF_TAMC */
268  C     if fluid is not water, by-pass find_rho, gmredi, surfaceForcing  C     if fluid is not water, by-pass find_rho, gmredi, surfaceForcing
269  C     and all vertical mixing schemes, but keep OBCS_CALC  C     and all vertical mixing schemes, but keep OBCS_CALC
270          IF ( fluidIsWater ) THEN          IF ( fluidIsWater ) THEN
271  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
# Line 268  C     uninitialised but inert locations. Line 278  C     uninitialised but inert locations.
278    
279          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
280           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
281            rhok   (i,j)   = 0. _d 0            rhoK   (i,j)   = 0. _d 0
282            rhoKM1 (i,j)   = 0. _d 0            rhoKm1 (i,j)   = 0. _d 0
283            rhoKP1 (i,j)   = 0. _d 0            rhoKp1 (i,j)   = 0. _d 0
284           ENDDO           ENDDO
285          ENDDO          ENDDO
286    
# Line 305  cph although some of these are re-initia Line 315  cph although some of these are re-initia
315             VisbeckK(i,j,bi,bj)   = 0. _d 0             VisbeckK(i,j,bi,bj)   = 0. _d 0
316  #  endif  #  endif
317  # endif /* ALLOW_GMREDI */  # endif /* ALLOW_GMREDI */
318    # ifdef ALLOW_KPP
319               KPPdiffKzS(i,j,k,bi,bj)  = 0. _d 0
320               KPPdiffKzT(i,j,k,bi,bj)  = 0. _d 0
321    # endif /* ALLOW_KPP */
322  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
323            ENDDO            ENDDO
324           ENDDO           ENDDO
# Line 338  C--     Start of diagnostic loop Line 352  C--     Start of diagnostic loop
352  C? Patrick, is this formula correct now that we change the loop range?  C? Patrick, is this formula correct now that we change the loop range?
353  C? Do we still need this?  C? Do we still need this?
354  cph kkey formula corrected.  cph kkey formula corrected.
355  cph Needed for rhok, rhokm1, in the case useGMREDI.  cph Needed for rhoK, rhoKm1, in the case useGMREDI.
356           kkey = (itdkey-1)*Nr + k           kkey = (itdkey-1)*Nr + k
357  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
358    
359  C--       Calculate gradients of potential density for isoneutral  C--       Calculate gradients of potential density for isoneutral
360  C         slope terms (e.g. GM/Redi tensor or IVDC diffusivity)  C         slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
361            IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)            IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
362       &                   .OR. doDiagsRho.GE.1 ) THEN       &         .OR. useSALT_PLUME .OR. doDiagsRho.GE.1 ) THEN
363  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
364              IF ( debugLevel .GE. debLevB )              IF ( debugLevel .GE. debLevB )
365       &       CALL DEBUG_CALL('FIND_RHO',myThid)       &       CALL DEBUG_CALL('FIND_RHO',myThid)
# Line 378  CADJ STORE salt (:,:,k-1,bi,bj) = comlev Line 392  CADJ STORE salt (:,:,k-1,bi,bj) = comlev
392  cph Avoid variable aliasing for adjoint !!!  cph Avoid variable aliasing for adjoint !!!
393              DO j=jMin,jMax              DO j=jMin,jMax
394               DO i=iMin,iMax               DO i=iMin,iMax
395                rhoKP1(i,j) = rhoK(i,j)                rhoKp1(i,j) = rhoK(i,j)
396               ENDDO               ENDDO
397              ENDDO              ENDDO
398              CALL GRAD_SIGMA(              CALL GRAD_SIGMA(
# Line 386  cph Avoid variable aliasing for adjoint Line 400  cph Avoid variable aliasing for adjoint
400       I             rhoK, rhoKm1, rhoKp1,       I             rhoK, rhoKm1, rhoKp1,
401       O             sigmaX, sigmaY, sigmaR,       O             sigmaX, sigmaY, sigmaR,
402       I             myThid )       I             myThid )
           ENDIF  
   
403  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
404  ctest# ifndef GM_EXCLUDE_CLIPPING  #ifdef GMREDI_WITH_STABLE_ADJOINT
405  CADJ STORE rhok   (:,:) = comlev1_bibj_k ,       key=kkey, byte=isbyte  cgf zero out adjoint fields to stabilize pkg/gmredi adjoint
406  ctest# endif  cgf -> cuts adjoint dependency from slope to state
407  CADJ STORE rhokm1 (:,:) = comlev1_bibj_k ,       key=kkey, byte=isbyte        CALL ZERO_ADJ_3D( bi, bj, sigmaX, myThid)
408          CALL ZERO_ADJ_3D( bi, bj, sigmaY, myThid)
409          CALL ZERO_ADJ_3D( bi, bj, sigmaR, myThid)
410    #endif
411  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
412              ENDIF
413    
414  C--       Implicit Vertical Diffusion for Convection  C--       Implicit Vertical Diffusion for Convection
415  c ==> should use sigmaR !!!  c ==> should use sigmaR !!!
416            IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN            IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
# Line 418  c ==> should use sigmaR !!! Line 435  c ==> should use sigmaR !!!
435  C--     end of diagnostic k loop (Nr:1)  C--     end of diagnostic k loop (Nr:1)
436          ENDDO          ENDDO
437    
438    #ifdef ALLOW_AUTODIFF_TAMC
439    CADJ STORE IVDConvCount(:,:,:,bi,bj)
440    CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
441    #endif
442    
443    C--     Diagnose Mixed Layer Depth:
444            IF ( useGMRedi .OR. doDiagsRho.GE.1 ) THEN
445              CALL CALC_OCE_MXLAYER( rhoK, sigmaR,
446         &              bi, bj, myTime, myIter, myThid )
447            ENDIF
448    
449    #ifdef ALLOW_SALT_PLUME
450            IF ( useSALT_PLUME ) THEN
451              CALL SALT_PLUME_CALC_DEPTH( rhoK, sigmaR,
452         &              bi, bj, myTime, myIter, myThid )
453            ENDIF
454    #endif /* ALLOW_SALT_PLUME */
455    
456  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
457          IF ( doDiagsRho.GE.1 ) THEN          IF ( doDiagsRho.GE.1 ) THEN
458            CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR  ', 0, Nr,            CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR  ', 0, Nr,
459       &         2, bi, bj, myThid)       &         2, bi, bj, myThid)
460          ENDIF          ENDIF
461  #endif  #endif /* ALLOW_DIAGNOSTICS */
462    
463  C--     Determines forcing terms based on external fields  C--     Determines forcing terms based on external fields
464  C       relaxation terms, etc.  C       relaxation terms, etc.
# Line 470  CADJ STORE surfaceForcingTice(:,:,bi,bj) Line 505  CADJ STORE surfaceForcingTice(:,:,bi,bj)
505  CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte  CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
506  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
507    
 #ifdef  ALLOW_GMREDI  
   
 #ifdef ALLOW_AUTODIFF_TAMC  
 # ifndef GM_EXCLUDE_CLIPPING  
 cph storing here is needed only for one GMREDI_OPTIONS:  
 cph define GM_BOLUS_ADVEC  
 cph keep it although TAF says you dont need to.  
 cph but I've avoided the #ifdef for now, in case more things change  
 CADJ STORE sigmaX(:,:,:)        = comlev1_bibj, key=itdkey, byte=isbyte  
 CADJ STORE sigmaY(:,:,:)        = comlev1_bibj, key=itdkey, byte=isbyte  
 CADJ STORE sigmaR(:,:,:)        = comlev1_bibj, key=itdkey, byte=isbyte  
 # endif  
 #endif /* ALLOW_AUTODIFF_TAMC */  
   
 C--     Calculate iso-neutral slopes for the GM/Redi parameterisation  
         IF (useGMRedi) THEN  
 #ifdef ALLOW_DEBUG  
           IF ( debugLevel .GE. debLevB )  
      &     CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)  
 #endif  
           CALL GMREDI_CALC_TENSOR(  
      I             bi, bj, iMin, iMax, jMin, jMax,  
      I             sigmaX, sigmaY, sigmaR,  
      I             myThid )  
 #ifdef ALLOW_AUTODIFF_TAMC  
         ELSE  
           CALL GMREDI_CALC_TENSOR_DUMMY(  
      I             bi, bj, iMin, iMax, jMin, jMax,  
      I             sigmaX, sigmaY, sigmaR,  
      I             myThid )  
 #endif /* ALLOW_AUTODIFF_TAMC */  
         ENDIF  
   
 #endif  /* ALLOW_GMREDI */  
   
508  #ifdef  ALLOW_KPP  #ifdef  ALLOW_KPP
509  C--     Compute KPP mixing coefficients  C--     Compute KPP mixing coefficients
510          IF (useKPP) THEN          IF (useKPP) THEN
# Line 513  C--     Compute KPP mixing coefficients Line 513  C--     Compute KPP mixing coefficients
513       &     CALL DEBUG_CALL('KPP_CALC',myThid)       &     CALL DEBUG_CALL('KPP_CALC',myThid)
514  #endif  #endif
515            CALL KPP_CALC(            CALL KPP_CALC(
516       I                  bi, bj, myTime, myThid )       I                  bi, bj, myTime, myIter, myThid )
517  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
518          ELSE          ELSE
519            CALL KPP_CALC_DUMMY(            CALL KPP_CALC_DUMMY(
520       I                  bi, bj, myTime, myThid )       I                  bi, bj, myTime, myIter, myThid )
521  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
522          ENDIF          ENDIF
523    
# Line 569  C--     Compute GGL90 mixing coefficient Line 569  C--     Compute GGL90 mixing coefficient
569          ENDIF          ENDIF
570  #endif /* ALLOW_TIMEAVE */  #endif /* ALLOW_TIMEAVE */
571    
572    #ifdef  ALLOW_GMREDI
573    #ifdef ALLOW_AUTODIFF_TAMC
574    # ifndef GM_EXCLUDE_CLIPPING
575    cph storing here is needed only for one GMREDI_OPTIONS:
576    cph define GM_BOLUS_ADVEC
577    cph keep it although TAF says you dont need to.
578    cph but I've avoided the #ifdef for now, in case more things change
579    CADJ STORE sigmaX(:,:,:)        = comlev1_bibj, key=itdkey, byte=isbyte
580    CADJ STORE sigmaY(:,:,:)        = comlev1_bibj, key=itdkey, byte=isbyte
581    CADJ STORE sigmaR(:,:,:)        = comlev1_bibj, key=itdkey, byte=isbyte
582    # endif
583    #endif /* ALLOW_AUTODIFF_TAMC */
584    
585    C--     Calculate iso-neutral slopes for the GM/Redi parameterisation
586            IF (useGMRedi) THEN
587    #ifdef ALLOW_DEBUG
588              IF ( debugLevel .GE. debLevB )
589         &     CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
590    #endif
591              CALL GMREDI_CALC_TENSOR(
592    c    I             bi, bj, iMin, iMax, jMin, jMax,
593    c    I             sigmaX, sigmaY, sigmaR,
594    c    I             myThid )
595         I             iMin, iMax, jMin, jMax,
596         I             sigmaX, sigmaY, sigmaR,
597         I             bi, bj, myTime, myIter, myThid )
598    #ifdef ALLOW_AUTODIFF_TAMC
599            ELSE
600              CALL GMREDI_CALC_TENSOR_DUMMY(
601    c    I             bi, bj, iMin, iMax, jMin, jMax,
602    c    I             sigmaX, sigmaY, sigmaR,
603    c    I             myThid )
604         I             iMin, iMax, jMin, jMax,
605         I             sigmaX, sigmaY, sigmaR,
606         I             bi, bj, myTime, myIter, myThid )
607    #endif /* ALLOW_AUTODIFF_TAMC */
608            ENDIF
609    #endif  /* ALLOW_GMREDI */
610    
611  #ifndef ALLOW_AUTODIFF_TAMC  #ifndef ALLOW_AUTODIFF_TAMC
612  C---  if fluid Is Water: end  C---  if fluid Is Water: end
613          ENDIF          ENDIF
# Line 591  C--   end bi,bj loops. Line 630  C--   end bi,bj loops.
630         ENDDO         ENDDO
631        ENDDO        ENDDO
632    
633    #ifdef  ALLOW_KPP
634          IF (useKPP) THEN
635            CALL KPP_DO_EXCH( myThid )
636          ENDIF
637    #endif  /* ALLOW_KPP */
638    
639  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
640        IF ( fluidIsWater .AND. useDiagnostics ) THEN        IF ( fluidIsWater .AND. useDiagnostics ) THEN
641          CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )          CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )

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

  ViewVC Help
Powered by ViewVC 1.1.22