| 1 |
C $Header: /u/gcmpack/MITgcm/model/src/calc_gt.F,v 1.54 2009/02/13 21:56:48 heimbach Exp $ |
| 2 |
C $Name: $ |
| 3 |
C added mods to 1.55 7/10/09 |
| 4 |
#include "PACKAGES_CONFIG.h" |
| 5 |
#include "CPP_OPTIONS.h" |
| 6 |
|
| 7 |
CBOP |
| 8 |
C !ROUTINE: CALC_GT |
| 9 |
C !INTERFACE: |
| 10 |
SUBROUTINE CALC_GT( |
| 11 |
I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown, |
| 12 |
I xA, yA, maskUp, uFld, vFld, wFld, |
| 13 |
I uTrans, vTrans, rTrans, rTransKp1, |
| 14 |
I KappaRT, diffKh3d_x, diffKh3d_y, |
| 15 |
U fVerT, |
| 16 |
I myTime,myIter,myThid ) |
| 17 |
C !DESCRIPTION: \bv |
| 18 |
C *==========================================================* |
| 19 |
C | SUBROUTINE CALC_GT |
| 20 |
C | o Calculate the temperature tendency terms. |
| 21 |
C *==========================================================* |
| 22 |
C | A procedure called EXTERNAL_FORCING_T is called from |
| 23 |
C | here. These procedures can be used to add per problem |
| 24 |
C | heat flux source terms. |
| 25 |
C | Note: Although it is slightly counter-intuitive the |
| 26 |
C | EXTERNAL_FORCING routine is not the place to put |
| 27 |
C | file I/O. Instead files that are required to |
| 28 |
C | calculate the external source terms are generally |
| 29 |
C | read during the model main loop. This makes the |
| 30 |
C | logisitics of multi-processing simpler and also |
| 31 |
C | makes the adjoint generation simpler. It also |
| 32 |
C | allows for I/O to overlap computation where that |
| 33 |
C | is supported by hardware. |
| 34 |
C | Aside from the problem specific term the code here |
| 35 |
C | forms the tendency terms due to advection and mixing |
| 36 |
C | The baseline implementation here uses a centered |
| 37 |
C | difference form for the advection term and a tensorial |
| 38 |
C | divergence of a flux form for the diffusive term. The |
| 39 |
C | diffusive term is formulated so that isopycnal mixing and |
| 40 |
C | GM-style subgrid-scale terms can be incorporated b simply |
| 41 |
C | setting the diffusion tensor terms appropriately. |
| 42 |
C *==========================================================* |
| 43 |
C \ev |
| 44 |
|
| 45 |
C !USES: |
| 46 |
IMPLICIT NONE |
| 47 |
C == GLobal variables == |
| 48 |
#include "SIZE.h" |
| 49 |
#include "DYNVARS.h" |
| 50 |
#include "EEPARAMS.h" |
| 51 |
#include "PARAMS.h" |
| 52 |
#include "RESTART.h" |
| 53 |
#ifdef ALLOW_GENERIC_ADVDIFF |
| 54 |
#include "GAD.h" |
| 55 |
#endif |
| 56 |
#ifdef ALLOW_AUTODIFF_TAMC |
| 57 |
# include "tamc.h" |
| 58 |
# include "tamc_keys.h" |
| 59 |
#endif |
| 60 |
|
| 61 |
C !INPUT/OUTPUT PARAMETERS: |
| 62 |
C == Routine arguments == |
| 63 |
C bi, bj, :: tile indices |
| 64 |
C iMin,iMax :: loop range for called routines |
| 65 |
C jMin,jMax :: loop range for called routines |
| 66 |
C k :: vertical index |
| 67 |
C kM1 :: =k-1 for k>1, =1 for k=1 |
| 68 |
C kUp :: index into 2 1/2D array, toggles between 1|2 |
| 69 |
C kDown :: index into 2 1/2D array, toggles between 2|1 |
| 70 |
C xA :: Tracer cell face area normal to X |
| 71 |
C yA :: Tracer cell face area normal to X |
| 72 |
C maskUp :: Land mask used to denote base of the domain. |
| 73 |
C uFld,vFld :: Local copy of horizontal velocity field |
| 74 |
C wFld :: Local copy of vertical velocity field |
| 75 |
C uTrans :: Zonal volume transport through cell face |
| 76 |
C vTrans :: Meridional volume transport through cell face |
| 77 |
C rTrans :: Vertical volume transport at interface k |
| 78 |
C rTransKp1 :: Vertical volume transport at inteface k+1 |
| 79 |
C KappaRT :: Vertical diffusion for Tempertature |
| 80 |
C fVerT :: Flux of temperature (T) in the vertical direction |
| 81 |
C at the upper(U) and lower(D) faces of a cell. |
| 82 |
C myTime :: current time |
| 83 |
C myIter :: current iteration number |
| 84 |
C myThid :: my Thread Id. number |
| 85 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
| 86 |
INTEGER k,kUp,kDown,kM1 |
| 87 |
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 88 |
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 89 |
_RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 90 |
_RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 91 |
_RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 92 |
_RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 93 |
_RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 94 |
_RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 95 |
_RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 96 |
_RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 97 |
_RL KappaRT(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 98 |
_RL diffKh3d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
| 99 |
_RL diffKh3d_y(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
| 100 |
_RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
| 101 |
_RL myTime |
| 102 |
INTEGER myIter |
| 103 |
INTEGER myThid |
| 104 |
CEOP |
| 105 |
|
| 106 |
#ifdef ALLOW_GENERIC_ADVDIFF |
| 107 |
C === Local variables === |
| 108 |
LOGICAL calcAdvection |
| 109 |
INTEGER iterNb |
| 110 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
| 111 |
INTEGER m1, m2 |
| 112 |
#endif |
| 113 |
|
| 114 |
#ifdef ALLOW_AUTODIFF_TAMC |
| 115 |
act1 = bi - myBxLo(myThid) |
| 116 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
| 117 |
act2 = bj - myByLo(myThid) |
| 118 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
| 119 |
act3 = myThid - 1 |
| 120 |
max3 = nTx*nTy |
| 121 |
act4 = ikey_dynamics - 1 |
| 122 |
itdkey = (act1 + 1) + act2*max1 |
| 123 |
& + act3*max1*max2 |
| 124 |
& + act4*max1*max2*max3 |
| 125 |
kkey = (itdkey-1)*Nr + k |
| 126 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
| 127 |
|
| 128 |
#ifdef ALLOW_AUTODIFF_TAMC |
| 129 |
C-- only the kUp part of fverT is set in this subroutine |
| 130 |
C-- the kDown is still required |
| 131 |
fVerT(1,1,kDown) = fVerT(1,1,kDown) |
| 132 |
# ifdef NONLIN_FRSURF |
| 133 |
CADJ STORE fVerT(:,:,:) = |
| 134 |
CADJ & comlev1_bibj_k, key=kkey, byte=isbyte, |
| 135 |
CADJ & kind = isbyte |
| 136 |
CADJ STORE gtNm1(:,:,k,bi,bj) = |
| 137 |
CADJ & comlev1_bibj_k, key=kkey, byte=isbyte, |
| 138 |
CADJ & kind = isbyte |
| 139 |
# endif |
| 140 |
#endif |
| 141 |
|
| 142 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
| 143 |
|
| 144 |
calcAdvection = tempAdvection .AND. .NOT.tempMultiDimAdvec |
| 145 |
iterNb = myIter |
| 146 |
IF (staggerTimeStep) iterNb = myIter -1 |
| 147 |
|
| 148 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
| 149 |
m1 = 1 + MOD(iterNb+1,2) |
| 150 |
m2 = 1 + MOD( iterNb ,2) |
| 151 |
CALL GAD_CALC_RHS( |
| 152 |
I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown, |
| 153 |
I xA, yA, maskUp, uFld, vFld, wFld, |
| 154 |
I uTrans, vTrans, rTrans, rTransKp1, |
| 155 |
I diffKhT, diffK4T, KappaRT, |
| 156 |
I gtNm(1-Olx,1-Oly,1,1,1,m2), theta, dTtracerLev, |
| 157 |
I GAD_TEMPERATURE, tempAdvScheme, tempVertAdvScheme, |
| 158 |
I calcAdvection, tempImplVertAdv, AdamsBashforth_T, |
| 159 |
I useGMRedi, useKPP, |
| 160 |
U fVerT, gT, |
| 161 |
I myTime, myIter, myThid ) |
| 162 |
#else /* ALLOW_ADAMSBASHFORTH_3 */ |
| 163 |
CALL GAD_CALC_RHS_RAF( |
| 164 |
I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown, |
| 165 |
I xA, yA, maskUp, uFld, vFld, wFld, |
| 166 |
I uTrans, vTrans, rTrans, rTransKp1, |
| 167 |
I diffKh3d_x, diffKh3d_y, |
| 168 |
I diffK4T, KappaRT, gtNm1, theta, dTtracerLev, |
| 169 |
I GAD_TEMPERATURE, tempAdvScheme, tempVertAdvScheme, |
| 170 |
I calcAdvection, tempImplVertAdv, AdamsBashforth_T, |
| 171 |
I useGMRedi, useKPP, |
| 172 |
U fVerT, gT, |
| 173 |
I myTime, myIter, myThid ) |
| 174 |
#endif |
| 175 |
|
| 176 |
C-- External thermal forcing term(s) inside Adams-Bashforth: |
| 177 |
IF ( tempForcing .AND. tracForcingOutAB.NE.1 ) |
| 178 |
& CALL EXTERNAL_FORCING_T( |
| 179 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
| 180 |
I myTime,myThid) |
| 181 |
|
| 182 |
IF ( AdamsBashforthGt ) THEN |
| 183 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
| 184 |
CALL ADAMS_BASHFORTH3( |
| 185 |
I bi, bj, k, |
| 186 |
U gT, gtNm, |
| 187 |
I tempStartAB, iterNb, myThid ) |
| 188 |
#else |
| 189 |
CALL ADAMS_BASHFORTH2( |
| 190 |
I bi, bj, k, |
| 191 |
U gT, gtNm1, |
| 192 |
I tempStartAB, iterNb, myThid ) |
| 193 |
#endif |
| 194 |
ENDIF |
| 195 |
|
| 196 |
C-- External thermal forcing term(s) outside Adams-Bashforth: |
| 197 |
IF ( tempForcing .AND. tracForcingOutAB.EQ.1 ) |
| 198 |
& CALL EXTERNAL_FORCING_T( |
| 199 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
| 200 |
I myTime,myThid) |
| 201 |
|
| 202 |
#ifdef NONLIN_FRSURF |
| 203 |
IF (nonlinFreeSurf.GT.0) THEN |
| 204 |
CALL FREESURF_RESCALE_G( |
| 205 |
I bi, bj, k, |
| 206 |
U gT, |
| 207 |
I myThid ) |
| 208 |
IF ( AdamsBashforthGt ) THEN |
| 209 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
| 210 |
CALL FREESURF_RESCALE_G( |
| 211 |
I bi, bj, k, |
| 212 |
U gtNm(1-Olx,1-Oly,1,1,1,1), |
| 213 |
I myThid ) |
| 214 |
CALL FREESURF_RESCALE_G( |
| 215 |
I bi, bj, k, |
| 216 |
U gtNm(1-Olx,1-Oly,1,1,1,2), |
| 217 |
I myThid ) |
| 218 |
#else |
| 219 |
CALL FREESURF_RESCALE_G( |
| 220 |
I bi, bj, k, |
| 221 |
U gtNm1, |
| 222 |
I myThid ) |
| 223 |
#endif |
| 224 |
ENDIF |
| 225 |
ENDIF |
| 226 |
#endif /* NONLIN_FRSURF */ |
| 227 |
|
| 228 |
#endif /* ALLOW_GENERIC_ADVDIFF */ |
| 229 |
|
| 230 |
RETURN |
| 231 |
END |