/[MITgcm]/MITgcm_contrib/jscott/code_rafmod/calc_gs.F
ViewVC logotype

Contents of /MITgcm_contrib/jscott/code_rafmod/calc_gs.F

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


Revision 1.1 - (show annotations) (download)
Tue Aug 21 16:34:17 2007 UTC (17 years, 11 months ago) by jscott
Branch: MAIN
code directory for crude ML horizontal mixing scheme

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_gs.F,v 1.47 2006/06/18 23:22:43 jmc 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, diffKh3d_x, diffKh3d_y,
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 #ifdef ALLOW_GENERIC_ADVDIFF
53 #include "GAD.h"
54 #endif
55 #ifdef ALLOW_AUTODIFF_TAMC
56 # include "tamc.h"
57 # include "tamc_keys.h"
58 #endif
59
60 C !INPUT/OUTPUT PARAMETERS:
61 C == Routine arguments ==
62 C bi, bj, :: tile indices
63 C iMin,iMax, jMin,jMax :: Range of points for which calculation
64 C results will be set.
65 C k :: vertical index
66 C kM1 :: =k-1 for k>1, =1 for k=1
67 C kUp :: index into 2 1/2D array, toggles between 1|2
68 C kDown :: index into 2 1/2D array, toggles between 2|1
69 C xA :: Tracer cell face area normal to X
70 C yA :: Tracer cell face area normal to X
71 C maskUp :: Land mask used to denote base of the domain.
72 C uFld,vFld :: Local copy of horizontal velocity field
73 C wFld :: Local copy of vertical velocity field
74 C uTrans :: Zonal volume transport through cell face
75 C vTrans :: Meridional volume transport through cell face
76 C rTrans :: Vertical volume transport at interface k
77 C rTransKp1 :: Vertical volume transport at inteface k+1
78 C KappaRS :: Vertical diffusion for Salinity
79 C fVerS :: Flux of salt (S) in the vertical direction
80 C at the upper(U) and lower(D) faces of a cell.
81 C myTime :: current time
82 C myIter :: current iteration number
83 C myThid :: my Thread Id. number
84
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 KappaRS(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 fVerS (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 fverS is set in this subroutine
130 C-- the kDown is still required
131 fVerS(1,1,kDown) = fVerS(1,1,kDown)
132 # ifdef NONLIN_FRSURF
133 CADJ STORE fVerS(:,:,:) =
134 CADJ & comlev1_bibj_k, key=kkey, byte=isbyte
135 CADJ STORE gsNm1(:,:,k,bi,bj) =
136 CADJ & comlev1_bibj_k, key=kkey, byte=isbyte
137 # endif
138 #endif
139
140 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
141
142 calcAdvection = saltAdvection .AND. .NOT.saltMultiDimAdvec
143 iterNb = myIter
144 IF (staggerTimeStep) iterNb = myIter - 1
145
146 #ifdef ALLOW_ADAMSBASHFORTH_3
147 m1 = 1 + MOD(iterNb+1,2)
148 m2 = 1 + MOD( iterNb ,2)
149 CALL GAD_CALC_RHS(
150 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
151 I xA, yA, maskUp, uFld, vFld, wFld,
152 I uTrans, vTrans, rTrans, rTransKp1,
153 I diffKhS, diffK4S, KappaRS,
154 I gsNm(1-Olx,1-Oly,1,1,1,m2), salt,
155 I GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
156 I calcAdvection, saltImplVertAdv, AdamsBashforth_S,
157 U fVerS, gS,
158 I myTime, myIter, myThid )
159 #else /* ALLOW_ADAMSBASHFORTH_3 */
160 CALL GAD_CALC_RHS_RAF(
161 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
162 I xA, yA, maskUp, uFld, vFld, wFld,
163 I uTrans, vTrans, rTrans, rTransKp1,
164 I diffKh3d_x, diffKh3d_y,
165 I diffK4S, KappaRS, gsNm1, salt,
166 I GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
167 I calcAdvection, saltImplVertAdv, AdamsBashforth_S,
168 U fVerS, gS,
169 I myTime, myIter, myThid )
170 #endif /* ALLOW_ADAMSBASHFORTH_3 */
171
172 C-- External salinity forcing term(s) inside Adams-Bashforth:
173 IF ( saltForcing .AND. tracForcingOutAB.NE.1 )
174 & CALL EXTERNAL_FORCING_S(
175 I iMin,iMax,jMin,jMax,bi,bj,k,
176 I myTime,myThid)
177
178 IF ( AdamsBashforthGs ) THEN
179 #ifdef ALLOW_ADAMSBASHFORTH_3
180 CALL ADAMS_BASHFORTH3(
181 I bi, bj, k,
182 U gS, gsNm,
183 I saltStartAB, iterNb, myThid )
184 #else
185 CALL ADAMS_BASHFORTH2(
186 I bi, bj, k,
187 U gS, gsNm1,
188 I iterNb, myThid )
189 #endif
190 ENDIF
191
192 C-- External salinity forcing term(s) outside Adams-Bashforth:
193 IF ( saltForcing .AND. tracForcingOutAB.EQ.1 )
194 & CALL EXTERNAL_FORCING_S(
195 I iMin,iMax,jMin,jMax,bi,bj,k,
196 I myTime,myThid)
197
198 #ifdef NONLIN_FRSURF
199 IF (nonlinFreeSurf.GT.0) THEN
200 CALL FREESURF_RESCALE_G(
201 I bi, bj, k,
202 U gS,
203 I myThid )
204 IF ( AdamsBashforthGs ) THEN
205 #ifdef ALLOW_ADAMSBASHFORTH_3
206 CALL FREESURF_RESCALE_G(
207 I bi, bj, k,
208 U gsNm(1-OLx,1-OLy,1,1,1,1),
209 I myThid )
210 CALL FREESURF_RESCALE_G(
211 I bi, bj, k,
212 U gsNm(1-OLx,1-OLy,1,1,1,2),
213 I myThid )
214 #else
215 CALL FREESURF_RESCALE_G(
216 I bi, bj, k,
217 U gsNm1,
218 I myThid )
219 #endif
220 ENDIF
221 ENDIF
222 #endif /* NONLIN_FRSURF */
223
224 #endif /* ALLOW_GENERIC_ADVDIFF */
225
226 RETURN
227 END

  ViewVC Help
Powered by ViewVC 1.1.22