/[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.29 - (show annotations) (download)
Wed Mar 1 03:07:18 2006 UTC (18 years, 2 months ago) by jmc
Branch: MAIN
Changes since 1.28: +3 -2 lines
added 1 argument to S/R GAD_CALC_RHS

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

  ViewVC Help
Powered by ViewVC 1.1.22