/[MITgcm]/MITgcm/pkg/ptracers/ptracers_integrate.F
ViewVC logotype

Contents of /MITgcm/pkg/ptracers/ptracers_integrate.F

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


Revision 1.14 - (show annotations) (download)
Sat Jun 26 02:39:31 2004 UTC (19 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint54, checkpoint54b_post, checkpoint54a_pre, checkpoint54a_post, checkpoint53g_post
Changes since 1.13: +2 -1 lines
T & S: separate Vert.Advec.Scheme from horizontal Advec.Scheme

1 C $Header: /u/gcmpack/MITgcm/pkg/ptracers/ptracers_integrate.F,v 1.13 2004/01/07 21:21:26 jmc Exp $
2 C $Name: $
3
4 #include "PTRACERS_OPTIONS.h"
5 cswdptr -- add ---
6 #ifdef ALLOW_GCHEM
7 # include "GCHEM_OPTIONS.h"
8 #endif
9 cswdptr -- end add ---
10
11 CBOP
12 C !ROUTINE: PTRACERS_INTEGRATE
13
14 C !INTERFACE: ==========================================================
15 SUBROUTINE PTRACERS_INTEGRATE(
16 I bi,bj,k,
17 I xA,yA,uTrans,vTrans,rTrans,
18 I rTransKp1,maskUp,
19 X rFlx,KappaRtr,
20 I myIter,myTime,myThid )
21
22 C !DESCRIPTION:
23 C Calculates tendancy for passive tracers and integrates forward
24 C in time.
25
26 C !USES: ===============================================================
27 IMPLICIT NONE
28 #include "SIZE.h"
29 #include "EEPARAMS.h"
30 #include "PARAMS.h"
31 #include "DYNVARS.h"
32 #include "PTRACERS.h"
33 #include "GAD.h"
34
35 C !INPUT PARAMETERS: ===================================================
36 C bi,bj :: tile indices
37 C k :: vertical level number
38 C xA :: face area at U points in level k
39 C yA :: face area at V points in level k
40 C uTrans :: zonal transport in level k
41 C vTrans :: meridional transport in level k
42 C rTrans :: vertical volume transport at interface k
43 C rTransKp1 :: vertical volume transport at interface k+1
44 C maskUp :: mask for vertical transport
45 C rFlx :: vertical flux
46 C KappaRtr :: vertical diffusion of passive tracers
47 C NOTE! This is infact KappaRS from thermodynamics()
48 C and is being used only temporarily
49 C until we removed the need to store
50 C a "3D" diffusivity
51 C myIter :: time-step number
52 C myTime :: model time
53 C myThid :: thread number
54 INTEGER bi,bj,k
55 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
56 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
57 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
58 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
59 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
60 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
61 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
62 _RL rFlx(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num)
63 _RL KappaRtr(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
64 INTEGER myIter
65 _RL myTime
66 INTEGER myThid
67
68 C !OUTPUT PARAMETERS: ==================================================
69 C none
70
71 #ifdef ALLOW_PTRACERS
72
73 C !LOCAL VARIABLES: ====================================================
74 C i,j,k,bi,bj,iTracer :: loop indices
75 C iMin,iMax,jMin,jMax :: loop ranges
76 C kUp,kDown :: toggle indices for even/odd level fluxes
77 C km1 :: =min(1,k-1)
78 C GAD_TR :: passive tracer id (GAD_TR1+iTracer-1)
79 INTEGER i,j,iTracer
80 INTEGER iMin,iMax,jMin,jMax
81 INTEGER kUp,kDown,km1
82 INTEGER GAD_TR
83 LOGICAL calcAdvection
84 CEOP
85
86 C Loop over tracers
87 DO iTracer=1,PTRACERS_numInUse
88
89 C Loop ranges for daughter routines
90 iMin = 1-OLx+2
91 iMax = sNx+OLx-1
92 jMin = 1-OLy+2
93 jMax = sNy+OLy-1
94
95 km1 = MAX(1,k-1)
96 kUp = 1+MOD(k+1,2)
97 kDown= 1+MOD(k,2)
98
99 C Calculate active tracer tendencies (gPtr) due to internal processes
100 C (advection, [explicit] diffusion, parameterizations,...)
101 calcAdvection = .NOT.multiDimAdvection
102 & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
103 & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
104 & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH
105 GAD_TR = GAD_TR1 + iTracer - 1
106 CALL GAD_CALC_RHS(
107 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
108 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
109 I uVel, vVel, wVel,
110 I PTRACERS_diffKh(iTracer),
111 I PTRACERS_diffK4(iTracer),
112 I KappaRtr,
113 I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
114 I GAD_TR,
115 I PTRACERS_advScheme(iTracer),
116 I PTRACERS_advScheme(iTracer),
117 I calcAdvection, .FALSE.,
118 U rFlx(1-Olx,1-Oly,1,iTracer),
119 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
120 I myThid )
121
122 C External forcing term(s)
123 cswdptr --add --
124 #ifdef ALLOW_GCHEM
125 #ifndef PTRACERS_SEPARATE_FORCING
126 IF ( forcing_In_AB )
127 & CALL GCHEM_FORCING_INT(
128 I bi,bj,iMin,iMax,jMin,jMax,k,
129 I iTracer,
130 I myTime,myIter, myThid)
131 #endif
132 #else
133 cswdptr - end add ---
134 IF ( forcing_In_AB )
135 & CALL PTRACERS_FORCING(
136 I bi,bj,iMin,iMax,jMin,jMax,k,
137 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
138 I surfaceTendencyPtr(1-Olx,1-Oly,1,1,iTracer),
139 I myIter,myTime,myThid)
140 cswdptr --add---
141 #endif
142 cswdptr -- end add ---
143
144 C If using Adams-Bashforth II, then extrapolate tendancies
145 IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
146 & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
147 & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH ) THEN
148 CALL ADAMS_BASHFORTH2(
149 I bi,bj,K,
150 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
151 U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
152 I myIter,myThid )
153 ENDIF
154
155 C External forcing term(s)
156 cswdptr - add--
157 #ifdef ALLOW_GCHEM
158 #ifndef PTRACERS_SEPARATE_FORCING
159 IF ( .NOT.forcing_In_AB )
160 & CALL GCHEM_FORCING_INT(
161 I bi,bj,iMin,iMax,jMin,jMax,k,
162 I iTracer,
163 I myTime,myIter, myThid)
164 #endif
165 #else
166 cswdptr - end add ---
167 IF ( .NOT.forcing_In_AB )
168 & CALL PTRACERS_FORCING(
169 I bi,bj,iMin,iMax,jMin,jMax,k,
170 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
171 I surfaceTendencyPtr(1-Olx,1-Oly,1,1,iTracer),
172 I myIter,myTime,myThid)
173 cswdptr - add--
174 #endif
175 cswdptr -- end add ---
176
177 #ifdef NONLIN_FRSURF
178 C Account for change in level thickness
179 IF (nonlinFreeSurf.GT.0) THEN
180 CALL FREESURF_RESCALE_G(
181 I bi,bj,K,
182 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
183 I myThid )
184 IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
185 & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
186 & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH )
187 & CALL FREESURF_RESCALE_G(
188 I bi,bj,K,
189 U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
190 I myThid )
191 ENDIF
192 #endif /* NONLIN_FRSURF */
193
194 C Integrate forward in time, storing in gPtr: G=T+dt*G
195 CALL TIMESTEP_TRACER(
196 I bi,bj,iMin,iMax,jMin,jMax,k,
197 I PTRACERS_advScheme(iTracer),
198 I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
199 I gPtr(1-Olx,1-Oly,1,1,1,iTracer),
200 I myIter,myThid )
201
202 C end of tracer loop
203 ENDDO
204
205 #endif /* ALLOW_PTRACERS */
206
207 RETURN
208 END

  ViewVC Help
Powered by ViewVC 1.1.22