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

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

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


Revision 1.28 - (hide annotations) (download)
Mon Feb 27 20:45:34 2006 UTC (18 years, 3 months ago) by heimbach
Branch: MAIN
Changes since 1.27: +2 -2 lines
Bring call gad_calc_rhs parameter list up-to-date w.r.t. calc_gt, calc_gs

1 heimbach 1.28 C $Header: /u/gcmpack/MITgcm/pkg/ptracers/ptracers_integrate.F,v 1.27 2005/11/06 22:52:49 jmc Exp $
2 adcroft 1.2 C $Name: $
3 adcroft 1.1
4     #include "PTRACERS_OPTIONS.h"
5    
6     CBOP
7 heimbach 1.7 C !ROUTINE: PTRACERS_INTEGRATE
8 adcroft 1.1
9     C !INTERFACE: ==========================================================
10 heimbach 1.7 SUBROUTINE PTRACERS_INTEGRATE(
11 adcroft 1.1 I bi,bj,k,
12 jmc 1.13 I xA,yA,uTrans,vTrans,rTrans,
13     I rTransKp1,maskUp,
14 heimbach 1.11 X rFlx,KappaRtr,
15 adcroft 1.1 I myIter,myTime,myThid )
16    
17     C !DESCRIPTION:
18 jmc 1.26 C Calculates tendency for passive tracers and integrates forward
19 adcroft 1.1 C in time.
20    
21     C !USES: ===============================================================
22     IMPLICIT NONE
23     #include "SIZE.h"
24     #include "EEPARAMS.h"
25 adcroft 1.5 #include "PARAMS.h"
26 jmc 1.13 #include "DYNVARS.h"
27 jmc 1.15 #include "PTRACERS_SIZE.h"
28 adcroft 1.1 #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 jmc 1.13 C rTrans :: vertical volume transport at interface k
39     C rTransKp1 :: vertical volume transport at interface k+1
40 adcroft 1.1 C maskUp :: mask for vertical transport
41 jmc 1.12 C rFlx :: vertical flux
42 jmc 1.18 C KappaRtr :: vertical diffusion of passive tracers, interf k
43 adcroft 1.1 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 jmc 1.13 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
53 adcroft 1.1 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
54 jmc 1.12 _RL rFlx(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num)
55 jmc 1.19 _RL KappaRtr(1-Olx:sNx+Olx,1-Oly:sNy+Oly,PTRACERS_num)
56 adcroft 1.1 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 jmc 1.27 C iTracer :: tracer index
67 adcroft 1.1 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 dimitri 1.8 C GAD_TR :: passive tracer id (GAD_TR1+iTracer-1)
71 jmc 1.27 INTEGER iTracer
72 adcroft 1.1 INTEGER iMin,iMax,jMin,jMax
73     INTEGER kUp,kDown,km1
74 dimitri 1.8 INTEGER GAD_TR
75 jmc 1.3 LOGICAL calcAdvection
76 jmc 1.23 INTEGER iterNb
77 adcroft 1.1 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 jmc 1.3 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 dimitri 1.8 GAD_TR = GAD_TR1 + iTracer - 1
99 jmc 1.13 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 jmc 1.19 I KappaRtr(1-Olx,1-Oly,iTracer),
106 heimbach 1.28 I gPtrNM1(1-Olx,1-Oly,1,1,1,iTracer),
107 jmc 1.26 I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
108 jmc 1.13 I GAD_TR,
109     I PTRACERS_advScheme(iTracer),
110 jmc 1.14 I PTRACERS_advScheme(iTracer),
111 jmc 1.22 I calcAdvection, PTRACERS_ImplVertAdv(iTracer),
112 jmc 1.13 U rFlx(1-Olx,1-Oly,1,iTracer),
113     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
114 jmc 1.17 I myTime, myIter, myThid )
115 adcroft 1.1
116     C External forcing term(s)
117 jmc 1.3 IF ( forcing_In_AB )
118     & CALL PTRACERS_FORCING(
119 mlosch 1.21 I bi,bj,iMin,iMax,jMin,jMax,k,iTracer,
120 heimbach 1.11 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
121 jmc 1.16 I surfaceForcingPtr(1-Olx,1-Oly,1,1,iTracer),
122 heimbach 1.11 I myIter,myTime,myThid)
123 adcroft 1.1
124 jmc 1.27 C If using Adams-Bashforth II, then extrapolate tendencies
125 spk 1.24 C gPtr is now the tracer tendency for explicit advection/diffusion
126 adcroft 1.1 IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
127     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
128     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH ) THEN
129 spk 1.24 #ifdef ALLOW_MATRIX
130     C If matrix is being computed, block call to S/R ADAMS_BASHFORTH2 to
131     C prevent gPtr from being replaced by the average of gPtr and gPtrNm1.
132     IF (.NOT.useMATRIX) THEN
133     #endif
134 jmc 1.23 iterNb = myIter
135     IF (staggerTimeStep) iterNb = myIter - 1
136 adcroft 1.1 CALL ADAMS_BASHFORTH2(
137     I bi,bj,K,
138     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
139     U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
140 jmc 1.23 I iterNb, myThid )
141 spk 1.24 #ifdef ALLOW_MATRIX
142     ENDIF
143     #endif
144 adcroft 1.1 ENDIF
145 jmc 1.3
146     C External forcing term(s)
147     IF ( .NOT.forcing_In_AB )
148     & CALL PTRACERS_FORCING(
149 mlosch 1.21 I bi,bj,iMin,iMax,jMin,jMax,k,iTracer,
150 heimbach 1.11 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
151 jmc 1.16 I surfaceForcingPtr(1-Olx,1-Oly,1,1,iTracer),
152 heimbach 1.11 I myIter,myTime,myThid)
153 adcroft 1.1
154     #ifdef NONLIN_FRSURF
155     C Account for change in level thickness
156     IF (nonlinFreeSurf.GT.0) THEN
157     CALL FREESURF_RESCALE_G(
158     I bi,bj,K,
159     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
160 jmc 1.4 I myThid )
161     IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
162     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
163     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH )
164     & CALL FREESURF_RESCALE_G(
165     I bi,bj,K,
166     U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
167 adcroft 1.1 I myThid )
168     ENDIF
169     #endif /* NONLIN_FRSURF */
170    
171     C Integrate forward in time, storing in gPtr: G=T+dt*G
172     CALL TIMESTEP_TRACER(
173     I bi,bj,iMin,iMax,jMin,jMax,k,
174     I PTRACERS_advScheme(iTracer),
175     I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
176     I gPtr(1-Olx,1-Oly,1,1,1,iTracer),
177     I myIter,myThid )
178    
179 mlosch 1.25 #ifdef ALLOW_OBCS
180     C Apply open boundary conditions
181     IF (useOBCS) THEN
182     CALL OBCS_APPLY_PTRACER(
183     I bi, bj, k, iTracer,
184     U gPtr(1-Olx,1-Oly,k,bi,bj,iTracer),
185     I myThid )
186     END IF
187     #endif /* ALLOW_OBCS */
188 adcroft 1.1 C end of tracer loop
189     ENDDO
190    
191     #endif /* ALLOW_PTRACERS */
192    
193     RETURN
194     END

  ViewVC Help
Powered by ViewVC 1.1.22