/[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.39 by heimbach, Mon Apr 16 22:44:29 2007 UTC revision 1.64 by heimbach, Thu May 1 23:52:24 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"
 #  include "exf_clim_fields.h"  
66  #  ifdef ALLOW_BULKFORMULAE  #  ifdef ALLOW_BULKFORMULAE
67  #   include "exf_constants.h"  #   include "EXF_CONSTANTS.h"
68  #  endif  #  endif
69  # endif  # endif
70  # ifdef ALLOW_SEAICE  # ifdef ALLOW_SEAICE
# Line 86  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 122  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  CADJ STORE atemp,aqh,precip    = comlev1, key = ikey_dynamics  CADJ STORE atemp,aqh,precip    = comlev1, key = ikey_dynamics
138  CADJ STORE swdown,lwdown       = comlev1, key = ikey_dynamics  CADJ STORE swdown,lwdown       = comlev1, key = ikey_dynamics
 CADJ STORE theta               = comlev1, key = ikey_dynamics  
139  cph# ifdef EXF_READ_EVAP  cph# ifdef EXF_READ_EVAP
140  CADJ STORE evap                = comlev1, key = ikey_dynamics  CADJ STORE evap                = comlev1, key = ikey_dynamics
141  cph# endif  cph# endif
 cph# ifdef SEAICE_ALLOW_DYNAMICS  
142  CADJ STORE uvel,vvel           = comlev1, key = ikey_dynamics  CADJ STORE uvel,vvel           = comlev1, key = ikey_dynamics
143  cph# endif  #  ifdef SEAICE_ALLOW_DYNAMICS
144  # ifdef SEAICE_CGRID  CADJ STORE uice                = comlev1, key = ikey_dynamics
145  CADJ STORE fu, fv              = comlev1, key = ikey_dynamics  CADJ STORE vice                = comlev1, key = ikey_dynamics
146  CADJ STORE seaicemasku         = comlev1, key = ikey_dynamics  #   ifdef SEAICE_ALLOW_EVP
147  CADJ STORE seaicemaskv         = comlev1, key = ikey_dynamics  CADJ STORE seaice_sigma1       = comlev1, key = ikey_dynamics
148    CADJ STORE seaice_sigma2       = comlev1, key = ikey_dynamics
149    CADJ STORE seaice_sigma12      = comlev1, key = ikey_dynamics
150    #   endif
151    #  endif
152    #  ifdef SEAICE_SALINITY
153    CADJ STORE salt                = comlev1, key = ikey_dynamics
154    #  endif
155  #  ifdef ATMOSPHERIC_LOADING  #  ifdef ATMOSPHERIC_LOADING
156    CADJ STORE pload               = comlev1, key = ikey_dynamics
157  CADJ STORE siceload            = comlev1, key = ikey_dynamics  CADJ STORE siceload            = comlev1, key = ikey_dynamics
158  #  endif  #  endif
159    #  ifdef NONLIN_FRSURF
160    CADJ STORE recip_hfacc         = comlev1, key = ikey_dynamics
161    #  endif
162  # endif  # endif
163  #endif  # ifdef ALLOW_DEBUG
 #ifdef ALLOW_DEBUG  
164          IF ( debugLevel .GE. debLevB )          IF ( debugLevel .GE. debLevB )
165       &    CALL DEBUG_CALL('SEAICE_MODEL',myThid)       &    CALL DEBUG_CALL('SEAICE_MODEL',myThid)
166  #endif  # endif
167          CALL TIMER_START('SEAICE_MODEL    [DO_OCEANIC_PHYS]', myThid)          CALL TIMER_START('SEAICE_MODEL    [DO_OCEANIC_PHYS]', myThid)
168          CALL SEAICE_MODEL( myTime, myIter, myThid )          CALL SEAICE_MODEL( myTime, myIter, myThid )
169          CALL TIMER_STOP ('SEAICE_MODEL    [DO_OCEANIC_PHYS]', myThid)          CALL TIMER_STOP ('SEAICE_MODEL    [DO_OCEANIC_PHYS]', myThid)
170  #ifdef ALLOW_COST_ICE  # ifdef ALLOW_COST
171          CALL COST_ICE_TEST ( myTime, myIter, myThid )          CALL SEAICE_COST_SENSI ( myTime, myIter, myThid )
172  #endif  # endif
173        ENDIF        ENDIF
174  #endif /* ALLOW_SEAICE */  #endif /* ALLOW_SEAICE */
175    
176    #ifdef ALLOW_AUTODIFF_TAMC
177    CADJ STORE sst, sss           = comlev1, key = ikey_dynamics
178    CADJ STORE qsw                = comlev1, key = ikey_dynamics
179    # ifdef ALLOW_SEAICE
180    CADJ STORE area               = comlev1, key = ikey_dynamics
181    # endif
182    #endif
183    
184  #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)  #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
185        IF ( useThSIce .AND. fluidIsWater ) THEN        IF ( useThSIce .AND. fluidIsWater ) THEN
186  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
# Line 195  C       and modify forcing terms includi Line 201  C       and modify forcing terms includi
201          IF ( debugLevel .GE. debLevB )          IF ( debugLevel .GE. debLevB )
202       &    CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid)       &    CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid)
203  #endif  #endif
204  C     compute temperature and (virtual) salt flux at the  C     compute temperature and (virtual) salt flux at the
205  C     shelf-ice ocean interface  C     shelf-ice ocean interface
206         CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',         CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
207       &       myThid)       &       myThid)
# Line 256  CHPF$ INDEPENDENT Line 262  CHPF$ INDEPENDENT
262       &                      + act3*max1*max2       &                      + act3*max1*max2
263       &                      + act4*max1*max2*max3       &                      + act4*max1*max2*max3
264  #else  /* ALLOW_AUTODIFF_TAMC */  #else  /* ALLOW_AUTODIFF_TAMC */
265  C     if fluid is not water, by-pass find_rho, gmredi, surfaceForcing  C     if fluid is not water, by-pass find_rho, gmredi, surfaceForcing
266  C     and all vertical mixing schemes, but keep OBCS_CALC  C     and all vertical mixing schemes, but keep OBCS_CALC
267          IF ( fluidIsWater ) THEN          IF ( fluidIsWater ) THEN
268  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
# Line 269  C     uninitialised but inert locations. Line 275  C     uninitialised but inert locations.
275    
276          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
277           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
278            rhok   (i,j)   = 0. _d 0            rhoK   (i,j)   = 0. _d 0
279            rhoKM1 (i,j)   = 0. _d 0            rhoKm1 (i,j)   = 0. _d 0
280            rhoKP1 (i,j)   = 0. _d 0            rhoKp1 (i,j)   = 0. _d 0
281           ENDDO           ENDDO
282          ENDDO          ENDDO
283    
# Line 306  cph although some of these are re-initia Line 312  cph although some of these are re-initia
312             VisbeckK(i,j,bi,bj)   = 0. _d 0             VisbeckK(i,j,bi,bj)   = 0. _d 0
313  #  endif  #  endif
314  # endif /* ALLOW_GMREDI */  # endif /* ALLOW_GMREDI */
315    # ifdef ALLOW_KPP
316               KPPdiffKzS(i,j,k,bi,bj)  = 0. _d 0
317               KPPdiffKzT(i,j,k,bi,bj)  = 0. _d 0
318    # endif /* ALLOW_KPP */
319  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
320            ENDDO            ENDDO
321           ENDDO           ENDDO
# Line 339  C--     Start of diagnostic loop Line 349  C--     Start of diagnostic loop
349  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?
350  C? Do we still need this?  C? Do we still need this?
351  cph kkey formula corrected.  cph kkey formula corrected.
352  cph Needed for rhok, rhokm1, in the case useGMREDI.  cph Needed for rhoK, rhoKm1, in the case useGMREDI.
353           kkey = (itdkey-1)*Nr + k           kkey = (itdkey-1)*Nr + k
354  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
355    
356  C--       Calculate gradients of potential density for isoneutral  C--       Calculate gradients of potential density for isoneutral
357  C         slope terms (e.g. GM/Redi tensor or IVDC diffusivity)  C         slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
358            IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)            IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
359       &                   .OR. doDiagsRho.GE.1 ) THEN       &         .OR. useSALT_PLUME .OR. doDiagsRho.GE.1 ) THEN
360  #ifdef ALLOW_DEBUG  #ifdef ALLOW_DEBUG
361              IF ( debugLevel .GE. debLevB )              IF ( debugLevel .GE. debLevB )
362       &       CALL DEBUG_CALL('FIND_RHO',myThid)       &       CALL DEBUG_CALL('FIND_RHO',myThid)
# Line 379  CADJ STORE salt (:,:,k-1,bi,bj) = comlev Line 389  CADJ STORE salt (:,:,k-1,bi,bj) = comlev
389  cph Avoid variable aliasing for adjoint !!!  cph Avoid variable aliasing for adjoint !!!
390              DO j=jMin,jMax              DO j=jMin,jMax
391               DO i=iMin,iMax               DO i=iMin,iMax
392                rhoKP1(i,j) = rhoK(i,j)                rhoKp1(i,j) = rhoK(i,j)
393               ENDDO               ENDDO
394              ENDDO              ENDDO
395              CALL GRAD_SIGMA(              CALL GRAD_SIGMA(
# Line 389  cph Avoid variable aliasing for adjoint Line 399  cph Avoid variable aliasing for adjoint
399       I             myThid )       I             myThid )
400            ENDIF            ENDIF
401    
 #ifdef ALLOW_AUTODIFF_TAMC  
 ctest# ifndef GM_EXCLUDE_CLIPPING  
 CADJ STORE rhok   (:,:) = comlev1_bibj_k ,       key=kkey, byte=isbyte  
 ctest# endif  
 CADJ STORE rhokm1 (:,:) = comlev1_bibj_k ,       key=kkey, byte=isbyte  
 #endif /* ALLOW_AUTODIFF_TAMC */  
402  C--       Implicit Vertical Diffusion for Convection  C--       Implicit Vertical Diffusion for Convection
403  c ==> should use sigmaR !!!  c ==> should use sigmaR !!!
404            IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN            IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
# Line 419  c ==> should use sigmaR !!! Line 423  c ==> should use sigmaR !!!
423  C--     end of diagnostic k loop (Nr:1)  C--     end of diagnostic k loop (Nr:1)
424          ENDDO          ENDDO
425    
426    #ifdef ALLOW_AUTODIFF_TAMC
427    CADJ STORE IVDConvCount(:,:,:,bi,bj)
428    CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
429    #endif
430    
431    C--     Diagnose Mixed Layer Depth:
432            IF ( useGMRedi .OR. doDiagsRho.GE.1 ) THEN
433              CALL CALC_OCE_MXLAYER( rhoK, sigmaR,
434         &              bi, bj, myTime, myIter, myThid )
435            ENDIF
436    
437    #ifdef ALLOW_SALT_PLUME
438            IF ( useSALT_PLUME ) THEN
439              CALL SALT_PLUME_CALC_DEPTH( rhoK, sigmaR,
440         &              bi, bj, myTime, myIter, myThid )
441            ENDIF
442    #endif /* ALLOW_SALT_PLUME */
443    
444  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
445          IF ( doDiagsRho.GE.1 ) THEN          IF ( doDiagsRho.GE.1 ) THEN
446            CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR  ', 0, Nr,            CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR  ', 0, Nr,
447       &         2, bi, bj, myThid)       &         2, bi, bj, myThid)
448          ENDIF          ENDIF
449  #endif  #endif /* ALLOW_DIAGNOSTICS */
450    
451  C--     Determines forcing terms based on external fields  C--     Determines forcing terms based on external fields
452  C       relaxation terms, etc.  C       relaxation terms, etc.
# Line 471  CADJ STORE surfaceForcingTice(:,:,bi,bj) Line 493  CADJ STORE surfaceForcingTice(:,:,bi,bj)
493  CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte  CADJ &     = comlev1_bibj, key=itdkey, byte=isbyte
494  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
495    
 #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 */  
   
496  #ifdef  ALLOW_KPP  #ifdef  ALLOW_KPP
497  C--     Compute KPP mixing coefficients  C--     Compute KPP mixing coefficients
498          IF (useKPP) THEN          IF (useKPP) THEN
# Line 514  C--     Compute KPP mixing coefficients Line 501  C--     Compute KPP mixing coefficients
501       &     CALL DEBUG_CALL('KPP_CALC',myThid)       &     CALL DEBUG_CALL('KPP_CALC',myThid)
502  #endif  #endif
503            CALL KPP_CALC(            CALL KPP_CALC(
504       I                  bi, bj, myTime, myThid )       I                  bi, bj, myTime, myIter, myThid )
505  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
506          ELSE          ELSE
507            CALL KPP_CALC_DUMMY(            CALL KPP_CALC_DUMMY(
508       I                  bi, bj, myTime, myThid )       I                  bi, bj, myTime, myIter, myThid )
509  #endif /* ALLOW_AUTODIFF_TAMC */  #endif /* ALLOW_AUTODIFF_TAMC */
510          ENDIF          ENDIF
511    
# Line 570  C--     Compute GGL90 mixing coefficient Line 557  C--     Compute GGL90 mixing coefficient
557          ENDIF          ENDIF
558  #endif /* ALLOW_TIMEAVE */  #endif /* ALLOW_TIMEAVE */
559    
560    #ifdef  ALLOW_GMREDI
561    #ifdef ALLOW_AUTODIFF_TAMC
562    # ifndef GM_EXCLUDE_CLIPPING
563    cph storing here is needed only for one GMREDI_OPTIONS:
564    cph define GM_BOLUS_ADVEC
565    cph keep it although TAF says you dont need to.
566    cph but I've avoided the #ifdef for now, in case more things change
567    CADJ STORE sigmaX(:,:,:)        = comlev1_bibj, key=itdkey, byte=isbyte
568    CADJ STORE sigmaY(:,:,:)        = comlev1_bibj, key=itdkey, byte=isbyte
569    CADJ STORE sigmaR(:,:,:)        = comlev1_bibj, key=itdkey, byte=isbyte
570    # endif
571    #endif /* ALLOW_AUTODIFF_TAMC */
572    
573    C--     Calculate iso-neutral slopes for the GM/Redi parameterisation
574            IF (useGMRedi) THEN
575    #ifdef ALLOW_DEBUG
576              IF ( debugLevel .GE. debLevB )
577         &     CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
578    #endif
579              CALL GMREDI_CALC_TENSOR(
580    c    I             bi, bj, iMin, iMax, jMin, jMax,
581    c    I             sigmaX, sigmaY, sigmaR,
582    c    I             myThid )
583         I             iMin, iMax, jMin, jMax,
584         I             sigmaX, sigmaY, sigmaR,
585         I             bi, bj, myTime, myIter, myThid )
586    #ifdef ALLOW_AUTODIFF_TAMC
587            ELSE
588              CALL GMREDI_CALC_TENSOR_DUMMY(
589    c    I             bi, bj, iMin, iMax, jMin, jMax,
590    c    I             sigmaX, sigmaY, sigmaR,
591    c    I             myThid )
592         I             iMin, iMax, jMin, jMax,
593         I             sigmaX, sigmaY, sigmaR,
594         I             bi, bj, myTime, myIter, myThid )
595    #endif /* ALLOW_AUTODIFF_TAMC */
596            ENDIF
597    #endif  /* ALLOW_GMREDI */
598    
599  #ifndef ALLOW_AUTODIFF_TAMC  #ifndef ALLOW_AUTODIFF_TAMC
600  C---  if fluid Is Water: end  C---  if fluid Is Water: end
601          ENDIF          ENDIF
# Line 592  C--   end bi,bj loops. Line 618  C--   end bi,bj loops.
618         ENDDO         ENDDO
619        ENDDO        ENDDO
620    
621    #ifdef  ALLOW_KPP
622          IF (useKPP) THEN
623            CALL KPP_DO_EXCH( myThid )
624          ENDIF
625    #endif  /* ALLOW_KPP */
626    
627  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
628        IF ( fluidIsWater .AND. useDiagnostics ) THEN        IF ( fluidIsWater .AND. useDiagnostics ) THEN
629          CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )          CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )

Legend:
Removed from v.1.39  
changed lines
  Added in v.1.64

  ViewVC Help
Powered by ViewVC 1.1.22