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

Contents of /MITgcm/model/src/calc_gs.F

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


Revision 1.52 - (show annotations) (download)
Wed Sep 15 03:41:59 2010 UTC (13 years, 8 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62k, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint63, checkpoint63d, checkpoint63e, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.51: +21 -1 lines
Adjoint compatible with combined AB3 and NLFS.

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_gs.F,v 1.51 2009/06/26 23:10:09 jahn Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: CALC_GS
9 C !INTERFACE:
10 SUBROUTINE CALC_GS(
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 KappaRS,
15 U fVerS,
16 I myTime,myIter,myThid )
17 C !DESCRIPTION: \bv
18 C *==========================================================*
19 C | SUBROUTINE CALC_GS
20 C | o Calculate the salt tendency terms.
21 C *==========================================================*
22 C | A procedure called EXTERNAL_FORCING_S is called from
23 C | here. These procedures can be used to add per problem
24 C | E-P 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, jMin,jMax :: Range of points for which calculation
65 C results will be set.
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 KappaRS :: Vertical diffusion for Salinity
80 C fVerS :: Flux of salt (S) 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
86 INTEGER bi,bj,iMin,iMax,jMin,jMax
87 INTEGER k,kUp,kDown,kM1
88 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92 _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93 _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 _RL KappaRS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
100 _RL myTime
101 INTEGER myIter
102 INTEGER myThid
103 CEOP
104
105 #ifdef ALLOW_GENERIC_ADVDIFF
106 C === Local variables ===
107 LOGICAL calcAdvection
108 INTEGER iterNb
109 #ifdef ALLOW_ADAMSBASHFORTH_3
110 INTEGER m1, m2
111 #endif
112
113 #ifdef ALLOW_AUTODIFF_TAMC
114 act1 = bi - myBxLo(myThid)
115 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
116 act2 = bj - myByLo(myThid)
117 max2 = myByHi(myThid) - myByLo(myThid) + 1
118 act3 = myThid - 1
119 max3 = nTx*nTy
120 act4 = ikey_dynamics - 1
121 itdkey = (act1 + 1) + act2*max1
122 & + act3*max1*max2
123 & + act4*max1*max2*max3
124 kkey = (itdkey-1)*Nr + k
125 #endif /* ALLOW_AUTODIFF_TAMC */
126
127 #ifdef ALLOW_AUTODIFF_TAMC
128 C-- only the kUp part of fverS is set in this subroutine
129 C-- the kDown is still required
130 fVerS(1,1,kDown) = fVerS(1,1,kDown)
131 # ifdef NONLIN_FRSURF
132 CADJ STORE fVerS(:,:,:) =
133 CADJ & comlev1_bibj_k, key=kkey, byte=isbyte,
134 CADJ & kind = isbyte
135 # ifndef ALLOW_ADAMSBASHFORTH_3
136 CADJ STORE gsNm1(:,:,k,bi,bj) =
137 CADJ & comlev1_bibj_k, key=kkey, byte=isbyte,
138 CADJ & kind = isbyte
139 # else
140 CADJ STORE gs(:,:,k,bi,bj) =
141 CADJ & comlev1_bibj_k, key=kkey, byte=isbyte,
142 CADJ & kind = isbyte
143 CADJ STORE gsNm(:,:,k,bi,bj,1) =
144 CADJ & comlev1_bibj_k, key=kkey, byte=isbyte,
145 CADJ & kind = isbyte
146 CADJ STORE gsNm(:,:,k,bi,bj,2) =
147 CADJ & comlev1_bibj_k, key=kkey, byte=isbyte,
148 CADJ & kind = isbyte
149 # endif
150 # endif
151 #endif
152
153 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
154
155 calcAdvection = saltAdvection .AND. .NOT.saltMultiDimAdvec
156 iterNb = myIter
157 IF (staggerTimeStep) iterNb = myIter - 1
158
159 #ifdef ALLOW_ADAMSBASHFORTH_3
160 m1 = 1 + MOD(iterNb+1,2)
161 m2 = 1 + MOD( iterNb ,2)
162 CALL GAD_CALC_RHS(
163 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
164 I xA, yA, maskUp, uFld, vFld, wFld,
165 I uTrans, vTrans, rTrans, rTransKp1,
166 I diffKhS, diffK4S, KappaRS,
167 I gsNm(1-Olx,1-Oly,1,1,1,m2), salt, dTtracerLev,
168 I GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
169 I calcAdvection, saltImplVertAdv, AdamsBashforth_S,
170 I useGMRedi, useKPP,
171 U fVerS, gS,
172 I myTime, myIter, myThid )
173 #else /* ALLOW_ADAMSBASHFORTH_3 */
174 CALL GAD_CALC_RHS(
175 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
176 I xA, yA, maskUp, uFld, vFld, wFld,
177 I uTrans, vTrans, rTrans, rTransKp1,
178 I diffKhS, diffK4S, KappaRS, gsNm1, salt, dTtracerLev,
179 I GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
180 I calcAdvection, saltImplVertAdv, AdamsBashforth_S,
181 I useGMRedi, useKPP,
182 U fVerS, gS,
183 I myTime, myIter, myThid )
184 #endif /* ALLOW_ADAMSBASHFORTH_3 */
185
186 C-- External salinity forcing term(s) inside Adams-Bashforth:
187 IF ( saltForcing .AND. tracForcingOutAB.NE.1 )
188 & CALL EXTERNAL_FORCING_S(
189 I iMin,iMax,jMin,jMax,bi,bj,k,
190 I myTime,myThid)
191
192 IF ( AdamsBashforthGs ) THEN
193 #ifdef ALLOW_ADAMSBASHFORTH_3
194 CALL ADAMS_BASHFORTH3(
195 I bi, bj, k,
196 U gS, gsNm,
197 I saltStartAB, iterNb, myThid )
198 #else
199 CALL ADAMS_BASHFORTH2(
200 I bi, bj, k,
201 U gS, gsNm1,
202 I saltStartAB, iterNb, myThid )
203 #endif
204 ENDIF
205
206 C-- External salinity forcing term(s) outside Adams-Bashforth:
207 IF ( saltForcing .AND. tracForcingOutAB.EQ.1 )
208 & CALL EXTERNAL_FORCING_S(
209 I iMin,iMax,jMin,jMax,bi,bj,k,
210 I myTime,myThid)
211
212 #ifdef NONLIN_FRSURF
213 IF (nonlinFreeSurf.GT.0) THEN
214 CALL FREESURF_RESCALE_G(
215 I bi, bj, k,
216 U gS,
217 I myThid )
218 IF ( AdamsBashforthGs ) THEN
219 #ifdef ALLOW_ADAMSBASHFORTH_3
220 # ifdef ALLOW_AUTODIFF_TAMC
221 CADJ STORE gsNm(:,:,k,bi,bj,1) =
222 CADJ & comlev1_bibj_k, key=kkey, byte=isbyte,
223 CADJ & kind = isbyte
224 CADJ STORE gsNm(:,:,k,bi,bj,2) =
225 CADJ & comlev1_bibj_k, key=kkey, byte=isbyte,
226 CADJ & kind = isbyte
227 # endif
228 CALL FREESURF_RESCALE_G(
229 I bi, bj, k,
230 U gsNm(1-OLx,1-OLy,1,1,1,1),
231 I myThid )
232 CALL FREESURF_RESCALE_G(
233 I bi, bj, k,
234 U gsNm(1-OLx,1-OLy,1,1,1,2),
235 I myThid )
236 #else
237 CALL FREESURF_RESCALE_G(
238 I bi, bj, k,
239 U gsNm1,
240 I myThid )
241 #endif
242 ENDIF
243 ENDIF
244 #endif /* NONLIN_FRSURF */
245
246 #endif /* ALLOW_GENERIC_ADVDIFF */
247
248 RETURN
249 END

  ViewVC Help
Powered by ViewVC 1.1.22