/[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.11 - (hide annotations) (download)
Sun Oct 26 01:10:34 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint52d_pre, checkpoint52, checkpoint51t_post, checkpoint51s_post, checkpoint52b_pre, checkpoint51q_post, checkpoint52b_post, checkpoint52c_post, checkpoint51r_post, checkpoint52a_pre, branch-netcdf, checkpoint51o_post, checkpoint52a_post, ecco_c52_e35, checkpoint51p_post, checkpoint51u_post
Branch point for: branch-nonh
Changes since 1.10: +10 -18 lines
o initialisation of rFlx extended to full array (required by TAF)
  and shifted to thermodynamics
o removed PTRACERS.h in ptracers routine
o added surfacetendencyPtr to S/R parameter list pracers_forcing

1 heimbach 1.11 C $Header: /u/gcmpack/MITgcm/pkg/ptracers/ptracers_integrate.F,v 1.10 2003/10/06 20:25:54 stephd Exp $
2 adcroft 1.2 C $Name: $
3 adcroft 1.1
4     #include "PTRACERS_OPTIONS.h"
5 stephd 1.6 cswdptr -- add ---
6 heimbach 1.7 #ifdef ALLOW_GCHEM
7     # include "GCHEM_OPTIONS.h"
8     #endif
9 stephd 1.6 cswdptr -- end add ---
10 adcroft 1.1
11     CBOP
12 heimbach 1.7 C !ROUTINE: PTRACERS_INTEGRATE
13 adcroft 1.1
14     C !INTERFACE: ==========================================================
15 heimbach 1.7 SUBROUTINE PTRACERS_INTEGRATE(
16 adcroft 1.1 I bi,bj,k,
17     I xA,yA,uTrans,vTrans,rTrans,maskUp,
18 heimbach 1.11 X rFlx,KappaRtr,
19 adcroft 1.1 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 adcroft 1.5 #include "PARAMS.h"
30 adcroft 1.1 #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 dimitri 1.8 C GAD_TR :: passive tracer id (GAD_TR1+iTracer-1)
74 adcroft 1.1 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 dimitri 1.8 INTEGER GAD_TR
79 jmc 1.3 LOGICAL calcAdvection
80 adcroft 1.1 CEOP
81    
82     C Loop over tracers
83     DO iTracer=1,PTRACERS_numInUse
84    
85     C Loop ranges for daughter routines
86     iMin = 1-OLx+2
87     iMax = sNx+OLx-1
88     jMin = 1-OLy+2
89     jMax = sNy+OLy-1
90    
91     km1 = MAX(1,k-1)
92     kUp = 1+MOD(k+1,2)
93     kDown= 1+MOD(k,2)
94    
95     C Calculate active tracer tendencies (gPtr) due to internal processes
96     C (advection, [explicit] diffusion, parameterizations,...)
97 jmc 1.3 calcAdvection = .NOT.multiDimAdvection
98     & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
99     & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
100     & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH
101 dimitri 1.8 GAD_TR = GAD_TR1 + iTracer - 1
102 adcroft 1.1 CALL GAD_CALC_RHS(
103     I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
104     I xA,yA,uTrans,vTrans,rTrans,maskUp,
105     I PTRACERS_diffKh(iTracer),
106     I PTRACERS_diffK4(iTracer),
107     I KappaRtr,
108     I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
109 dimitri 1.8 I GAD_TR,
110 jmc 1.3 I PTRACERS_advScheme(iTracer),calcAdvection,
111 adcroft 1.1 U rFlx(1-Olx,1-Oly,1,iTracer),
112     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
113     I myThid )
114    
115     C External forcing term(s)
116 stephd 1.10 cswdptr --add --
117     #ifdef ALLOW_GCHEM
118     #ifndef PTRACERS_SEPARATE_FORCING
119     IF ( forcing_In_AB )
120     & CALL GCHEM_FORCING_INT(
121     I bi,bj,iMin,iMax,jMin,jMax,k,
122     I iTracer,
123     I myTime,myIter, myThid)
124     #endif
125     #else
126 stephd 1.6 cswdptr - end add ---
127 jmc 1.3 IF ( forcing_In_AB )
128     & CALL PTRACERS_FORCING(
129 heimbach 1.11 I bi,bj,iMin,iMax,jMin,jMax,k,
130     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
131     I surfaceTendencyPtr(1-Olx,1-Oly,1,1,iTracer),
132     I myIter,myTime,myThid)
133 stephd 1.6 cswdptr --add---
134     #endif
135     cswdptr -- end add ---
136 adcroft 1.1
137     C If using Adams-Bashforth II, then extrapolate tendancies
138     IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
139     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
140     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH ) THEN
141     CALL ADAMS_BASHFORTH2(
142     I bi,bj,K,
143     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
144     U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
145     I myIter,myThid )
146     ENDIF
147 jmc 1.3
148     C External forcing term(s)
149 stephd 1.6 cswdptr - add--
150 stephd 1.10 #ifdef ALLOW_GCHEM
151     #ifndef PTRACERS_SEPARATE_FORCING
152     IF ( .NOT.forcing_In_AB )
153     & CALL GCHEM_FORCING_INT(
154     I bi,bj,iMin,iMax,jMin,jMax,k,
155     I iTracer,
156     I myTime,myIter, myThid)
157     #endif
158     #else
159 stephd 1.6 cswdptr - end add ---
160 jmc 1.3 IF ( .NOT.forcing_In_AB )
161     & CALL PTRACERS_FORCING(
162 heimbach 1.11 I bi,bj,iMin,iMax,jMin,jMax,k,
163     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
164     I surfaceTendencyPtr(1-Olx,1-Oly,1,1,iTracer),
165     I myIter,myTime,myThid)
166 stephd 1.6 cswdptr - add--
167     #endif
168     cswdptr -- end add ---
169 adcroft 1.1
170     #ifdef NONLIN_FRSURF
171     C Account for change in level thickness
172     IF (nonlinFreeSurf.GT.0) THEN
173     CALL FREESURF_RESCALE_G(
174     I bi,bj,K,
175     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
176 jmc 1.4 I myThid )
177     IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
178     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
179     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH )
180     & CALL FREESURF_RESCALE_G(
181     I bi,bj,K,
182     U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
183 adcroft 1.1 I myThid )
184     ENDIF
185     #endif /* NONLIN_FRSURF */
186    
187     C Integrate forward in time, storing in gPtr: G=T+dt*G
188     CALL TIMESTEP_TRACER(
189     I bi,bj,iMin,iMax,jMin,jMax,k,
190     I PTRACERS_advScheme(iTracer),
191     I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
192     I gPtr(1-Olx,1-Oly,1,1,1,iTracer),
193     I myIter,myThid )
194    
195     C end of tracer loop
196     ENDDO
197    
198     #endif /* ALLOW_PTRACERS */
199    
200     RETURN
201     END

  ViewVC Help
Powered by ViewVC 1.1.22