/[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.6 - (hide annotations) (download)
Wed Jun 25 21:03:28 2003 UTC (20 years, 11 months ago) by stephd
Branch: MAIN
CVS Tags: checkpoint51b_pre, checkpoint51a_post
Changes since 1.5: +16 -1 lines
files needed for biogeochemistry package (gchem)

1 stephd 1.6 C $Header: /u/u0/gcmpack/MITgcm/pkg/ptracers/ptracers_integrate.F,v 1.5 2002/06/24 18:47:10 adcroft Exp $
2 adcroft 1.2 C $Name: $
3 adcroft 1.1
4     #include "PTRACERS_OPTIONS.h"
5 stephd 1.6 cswdptr -- add ---
6     #include "GCHEM_OPTIONS.h"
7     cswdptr -- end add ---
8 adcroft 1.1
9     CBOP
10     C !ROUTINE: PTRACERS_INTEGERATE
11    
12     C !INTERFACE: ==========================================================
13     SUBROUTINE PTRACERS_INTEGERATE(
14     I bi,bj,k,
15     I xA,yA,uTrans,vTrans,rTrans,maskUp,
16     X KappaRtr,
17     I myIter,myTime,myThid )
18    
19     C !DESCRIPTION:
20     C Calculates tendancy for passive tracers and integrates forward
21     C in time.
22    
23     C !USES: ===============================================================
24     IMPLICIT NONE
25     #include "SIZE.h"
26     #include "EEPARAMS.h"
27 adcroft 1.5 #include "PARAMS.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     C rTrans :: vertical transport across level k
39     C maskUp :: mask for vertical transport
40     C KappaRtr :: vertical diffusion of passive tracers
41     C NOTE! This is infact KappaRS from thermodynamics()
42     C and is being used only temporarily
43     C until we removed the need to store
44     C a "3D" diffusivity
45     C myIter :: time-step number
46     C myTime :: model time
47     C myThid :: thread number
48     INTEGER bi,bj,k
49     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
50     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
51     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
52     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
53     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
54     _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
55     _RL KappaRtr(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
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 i,j,k,bi,bj,iTracer :: loop indices
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 rFlx :: vertical flux
71     INTEGER i,j,iTracer
72     INTEGER iMin,iMax,jMin,jMax
73     INTEGER kUp,kDown,km1
74     _RL rFlx(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num)
75 jmc 1.3 LOGICAL calcAdvection
76 adcroft 1.1 CEOP
77    
78     C Loop over tracers
79     DO iTracer=1,PTRACERS_numInUse
80    
81     C Initialize vertical flux to zero and set no-flux across k=Nr+1
82     IF (k.EQ.Nr) THEN
83     DO j=1-Oly,sNy+Oly
84     DO i=1-Olx,sNx+Olx
85     rFlx(i,j,1,iTracer)=0.
86     rFlx(i,j,2,iTracer)=0.
87     ENDDO
88     ENDDO
89     ENDIF
90    
91     C Loop ranges for daughter routines
92     iMin = 1-OLx+2
93     iMax = sNx+OLx-1
94     jMin = 1-OLy+2
95     jMax = sNy+OLy-1
96    
97     km1 = MAX(1,k-1)
98     kUp = 1+MOD(k+1,2)
99     kDown= 1+MOD(k,2)
100    
101     C Calculate active tracer tendencies (gPtr) due to internal processes
102     C (advection, [explicit] diffusion, parameterizations,...)
103 jmc 1.3 calcAdvection = .NOT.multiDimAdvection
104     & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
105     & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
106     & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH
107 adcroft 1.1 CALL GAD_CALC_RHS(
108     I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
109     I xA,yA,uTrans,vTrans,rTrans,maskUp,
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_TR1,
115 jmc 1.3 I PTRACERS_advScheme(iTracer),calcAdvection,
116 adcroft 1.1 U rFlx(1-Olx,1-Oly,1,iTracer),
117     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
118     I myThid )
119    
120     C External forcing term(s)
121 stephd 1.6 cswdptr - add--
122     #ifndef PTRACERS_SEPERATE_FORCING
123     cswdptr - end add ---
124 jmc 1.3 IF ( forcing_In_AB )
125     & CALL PTRACERS_FORCING(
126 adcroft 1.1 I bi,bj,k,iTracer,
127 adcroft 1.2 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
128 adcroft 1.1 I myIter,myTime,myThid)
129 stephd 1.6 cswdptr --add---
130     #endif
131     cswdptr -- end add ---
132 adcroft 1.1
133     C If using Adams-Bashforth II, then extrapolate tendancies
134     IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
135     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
136     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH ) THEN
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 myIter,myThid )
142     ENDIF
143 jmc 1.3
144     C External forcing term(s)
145 stephd 1.6 cswdptr - add--
146     #ifndef PTRACERS_SEPERATE_FORCING
147     cswdptr - end add ---
148 jmc 1.3 IF ( .NOT.forcing_In_AB )
149     & CALL PTRACERS_FORCING(
150     I bi,bj,k,iTracer,
151     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
152     I myIter,myTime,myThid)
153 stephd 1.6 cswdptr - add--
154     #endif
155     cswdptr -- end add ---
156 adcroft 1.1
157     #ifdef NONLIN_FRSURF
158     C Account for change in level thickness
159     IF (nonlinFreeSurf.GT.0) THEN
160     CALL FREESURF_RESCALE_G(
161     I bi,bj,K,
162     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
163 jmc 1.4 I myThid )
164     IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
165     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
166     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH )
167     & CALL FREESURF_RESCALE_G(
168     I bi,bj,K,
169     U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
170 adcroft 1.1 I myThid )
171     ENDIF
172     #endif /* NONLIN_FRSURF */
173    
174     C Integrate forward in time, storing in gPtr: G=T+dt*G
175     CALL TIMESTEP_TRACER(
176     I bi,bj,iMin,iMax,jMin,jMax,k,
177     I PTRACERS_advScheme(iTracer),
178     I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
179     I gPtr(1-Olx,1-Oly,1,1,1,iTracer),
180     I myIter,myThid )
181    
182     C end of tracer loop
183     ENDDO
184    
185     #endif /* ALLOW_PTRACERS */
186    
187     RETURN
188     END

  ViewVC Help
Powered by ViewVC 1.1.22