/[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.10 - (show annotations) (download)
Mon Oct 6 20:25:54 2003 UTC (20 years, 7 months ago) by stephd
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint51o_pre, checkpoint51l_post, checkpoint51n_post, checkpoint51j_post, checkpoint51n_pre, checkpoint51l_pre, checkpoint51h_pre, checkpoint51i_post, checkpoint51i_pre, checkpoint51m_post
Branch point for: tg2-branch, checkpoint51n_branch
Changes since 1.9: +20 -4 lines
modifications so gchem pkg is more versatile

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

  ViewVC Help
Powered by ViewVC 1.1.22