/[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.41 - (show annotations) (download)
Sun Nov 6 22:19:08 2005 UTC (18 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57x_post, checkpoint57y_pre
Changes since 1.40: +33 -9 lines
Allow to apply AB on T,S rather than on AB(gT,gS):
 - implemented within #ifdef ALLOW_ADAMSBASHFORTH_3
 - use the same arrays (gtNm,gsNm) to hold tracer field at previous
   time-steps (if AB(T,S)) and tendencies (if AB(gT,gS)).

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_gs.F,v 1.40 2005/04/15 14:18:50 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,uTrans,vTrans,rTrans,rTransKp1,maskUp,
13 I KappaRS,
14 U fVerS,
15 I myTime,myIter,myThid )
16 C !DESCRIPTION: \bv
17 C *==========================================================*
18 C | SUBROUTINE CALC_GS
19 C | o Calculate the salt tendency terms.
20 C *==========================================================*
21 C | A procedure called EXTERNAL_FORCING_S is called from
22 C | here. These procedures can be used to add per problem
23 C | E-P flux source terms.
24 C | Note: Although it is slightly counter-intuitive the
25 C | EXTERNAL_FORCING routine is not the place to put
26 C | file I/O. Instead files that are required to
27 C | calculate the external source terms are generally
28 C | read during the model main loop. This makes the
29 C | logisitics of multi-processing simpler and also
30 C | makes the adjoint generation simpler. It also
31 C | allows for I/O to overlap computation where that
32 C | is supported by hardware.
33 C | Aside from the problem specific term the code here
34 C | forms the tendency terms due to advection and mixing
35 C | The baseline implementation here uses a centered
36 C | difference form for the advection term and a tensorial
37 C | divergence of a flux form for the diffusive term. The
38 C | diffusive term is formulated so that isopycnal mixing and
39 C | GM-style subgrid-scale terms can be incorporated b simply
40 C | setting the diffusion tensor terms appropriately.
41 C *==========================================================*
42 C \ev
43
44 C !USES:
45 IMPLICIT NONE
46 C == GLobal variables ==
47 #include "SIZE.h"
48 #include "DYNVARS.h"
49 #include "EEPARAMS.h"
50 #include "PARAMS.h"
51 #ifdef ALLOW_GENERIC_ADVDIFF
52 #include "GAD.h"
53 #endif
54
55 C !INPUT/OUTPUT PARAMETERS:
56 C == Routine arguments ==
57 C fVerS :: Flux of salt (S) in the vertical
58 C direction at the upper(U) and lower(D) faces of a cell.
59 C maskUp :: Land mask used to denote base of the domain.
60 C xA :: Tracer cell face area normal to X
61 C yA :: Tracer cell face area normal to X
62 C uTrans :: Zonal volume transport through cell face
63 C vTrans :: Meridional volume transport through cell face
64 C rTrans :: Vertical volume transport at interface k
65 C rTransKp1 :: Vertical volume transport at inteface k+1
66 C bi, bj, iMin, iMax, jMin, jMax :: Range of points for which calculation
67 C results will be set.
68 C myThid :: Instance number for this innvocation of CALC_GT
69 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
70 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
71 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 _RL KappaRS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 INTEGER k,kUp,kDown,kM1
79 INTEGER bi,bj,iMin,iMax,jMin,jMax
80 _RL myTime
81 INTEGER myIter
82 INTEGER myThid
83
84 CEOP
85
86 #ifdef ALLOW_GENERIC_ADVDIFF
87 C === Local variables ===
88 LOGICAL calcAdvection
89 INTEGER iterNb
90 #ifdef ALLOW_ADAMSBASHFORTH_3
91 INTEGER m1, m2
92 #endif
93
94 #ifdef ALLOW_AUTODIFF_TAMC
95 C-- only the kUp part of fverS is set in this subroutine
96 C-- the kDown is still required
97 fVerS(1,1,kDown) = fVerS(1,1,kDown)
98 #endif
99
100 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
101
102 calcAdvection = saltAdvection .AND. .NOT.saltMultiDimAdvec
103 iterNb = myIter
104 IF (staggerTimeStep) iterNb = myIter - 1
105
106 #ifdef ALLOW_ADAMSBASHFORTH_3
107 IF ( AdamsBashforth_S ) THEN
108 m1 = 1 + MOD(iterNb+1,2)
109 m2 = 1 + MOD( iterNb ,2)
110 CALL GAD_CALC_RHS(
111 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
112 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
113 I uVel, vVel, wVel,
114 I diffKhS, diffK4S, KappaRS,
115 I gsNm(1-Olx,1-Oly,1,1,1,m2), salt,
116 I GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
117 I calcAdvection, saltImplVertAdv,
118 U fVerS, gS,
119 I myTime, myIter, myThid )
120 ELSE
121 #endif /* ALLOW_ADAMSBASHFORTH_3 */
122 CALL GAD_CALC_RHS(
123 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
124 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
125 I uVel, vVel, wVel,
126 I diffKhS, diffK4S, KappaRS, salt, salt,
127 I GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
128 I calcAdvection, saltImplVertAdv,
129 U fVerS, gS,
130 I myTime, myIter, myThid )
131 #ifdef ALLOW_ADAMSBASHFORTH_3
132 ENDIF
133 #endif
134
135 C-- External salinity forcing term(s) inside Adams-Bashforth:
136 IF ( saltForcing .AND. forcing_In_AB )
137 & CALL EXTERNAL_FORCING_S(
138 I iMin,iMax,jMin,jMax,bi,bj,k,
139 I myTime,myThid)
140
141 IF ( AdamsBashforthGs ) THEN
142 #ifdef ALLOW_ADAMSBASHFORTH_3
143 CALL ADAMS_BASHFORTH3(
144 I bi, bj, k,
145 U gS, gsNm,
146 I saltStartAB, iterNb, myThid )
147 #else
148 CALL ADAMS_BASHFORTH2(
149 I bi, bj, k,
150 U gS, gsNm1,
151 I iterNb, myThid )
152 #endif
153 ENDIF
154
155 C-- External salinity forcing term(s) outside Adams-Bashforth:
156 IF ( saltForcing .AND. .NOT.forcing_In_AB )
157 & CALL EXTERNAL_FORCING_S(
158 I iMin,iMax,jMin,jMax,bi,bj,k,
159 I myTime,myThid)
160
161 #ifdef NONLIN_FRSURF
162 IF (nonlinFreeSurf.GT.0) THEN
163 CALL FREESURF_RESCALE_G(
164 I bi, bj, k,
165 U gS,
166 I myThid )
167 IF ( AdamsBashforthGs ) THEN
168 #ifdef ALLOW_ADAMSBASHFORTH_3
169 CALL FREESURF_RESCALE_G(
170 I bi, bj, k,
171 U gsNm(1-OLx,1-OLy,1,1,1,1),
172 I myThid )
173 CALL FREESURF_RESCALE_G(
174 I bi, bj, k,
175 U gsNm(1-OLx,1-OLy,1,1,1,2),
176 I myThid )
177 #else
178 CALL FREESURF_RESCALE_G(
179 I bi, bj, k,
180 U gsNm1,
181 I myThid )
182 #endif
183 ENDIF
184 ENDIF
185 #endif /* NONLIN_FRSURF */
186
187 #endif /* ALLOW_GENERIC_ADVDIFF */
188
189 RETURN
190 END

  ViewVC Help
Powered by ViewVC 1.1.22