/[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.25 - (hide annotations) (download)
Mon Oct 10 05:53:48 2005 UTC (18 years, 7 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint57v_post, checkpint57u_post, checkpoint57w_post
Changes since 1.24: +10 -1 lines
o OBCS and PTRACERS: add open boundary support for passive tracers
  - either use homogenous (pseudo) v.Neumann conditions or prescribe
    OB-values from file; this is not different from the way theta and salinity
    are treated
  - however, Orlanski-radiation conditions are not supported, and the model
    will stop if you use pTracers and Orlanski at the same time.
  - beefed up the rountine obcs_external_fields_load: now only those open
    boundary values are overwritten with values from files for which there
    are really files, otherwise the OB-fields remain untouched. This makes
    it possible to use different OBs at different ends of the domain (as
    with EXF)
  - TODO: add support for OB?w and OB?eta, which can currently not be read
    from a file.

1 mlosch 1.25 C $Header: /u/gcmpack/MITgcm/pkg/ptracers/ptracers_integrate.F,v 1.24 2005/04/20 15:54:57 spk 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     C Calculates tendancy 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 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     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 dimitri 1.8 C GAD_TR :: passive tracer id (GAD_TR1+iTracer-1)
71 adcroft 1.1 INTEGER i,j,iTracer
72     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 jmc 1.13 I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
107     I GAD_TR,
108     I PTRACERS_advScheme(iTracer),
109 jmc 1.14 I PTRACERS_advScheme(iTracer),
110 jmc 1.22 I calcAdvection, PTRACERS_ImplVertAdv(iTracer),
111 jmc 1.13 U rFlx(1-Olx,1-Oly,1,iTracer),
112     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
113 jmc 1.17 I myTime, myIter, myThid )
114 adcroft 1.1
115     C External forcing term(s)
116 jmc 1.3 IF ( forcing_In_AB )
117     & CALL PTRACERS_FORCING(
118 mlosch 1.21 I bi,bj,iMin,iMax,jMin,jMax,k,iTracer,
119 heimbach 1.11 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
120 jmc 1.16 I surfaceForcingPtr(1-Olx,1-Oly,1,1,iTracer),
121 heimbach 1.11 I myIter,myTime,myThid)
122 adcroft 1.1
123     C If using Adams-Bashforth II, then extrapolate tendancies
124 spk 1.24 C gPtr is now the tracer tendency for explicit advection/diffusion
125 adcroft 1.1 IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
126     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
127     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH ) THEN
128 spk 1.24 #ifdef ALLOW_MATRIX
129     C If matrix is being computed, block call to S/R ADAMS_BASHFORTH2 to
130     C prevent gPtr from being replaced by the average of gPtr and gPtrNm1.
131     IF (.NOT.useMATRIX) THEN
132     #endif
133 jmc 1.23 iterNb = myIter
134     IF (staggerTimeStep) iterNb = myIter - 1
135 adcroft 1.1 CALL ADAMS_BASHFORTH2(
136     I bi,bj,K,
137     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
138     U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
139 jmc 1.23 I iterNb, myThid )
140 spk 1.24 #ifdef ALLOW_MATRIX
141     ENDIF
142     #endif
143 adcroft 1.1 ENDIF
144 jmc 1.3
145     C External forcing term(s)
146     IF ( .NOT.forcing_In_AB )
147     & CALL PTRACERS_FORCING(
148 mlosch 1.21 I bi,bj,iMin,iMax,jMin,jMax,k,iTracer,
149 heimbach 1.11 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
150 jmc 1.16 I surfaceForcingPtr(1-Olx,1-Oly,1,1,iTracer),
151 heimbach 1.11 I myIter,myTime,myThid)
152 adcroft 1.1
153     #ifdef NONLIN_FRSURF
154     C Account for change in level thickness
155     IF (nonlinFreeSurf.GT.0) THEN
156     CALL FREESURF_RESCALE_G(
157     I bi,bj,K,
158     U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
159 jmc 1.4 I myThid )
160     IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
161     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
162     & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH )
163     & CALL FREESURF_RESCALE_G(
164     I bi,bj,K,
165     U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
166 adcroft 1.1 I myThid )
167     ENDIF
168     #endif /* NONLIN_FRSURF */
169    
170     C Integrate forward in time, storing in gPtr: G=T+dt*G
171     CALL TIMESTEP_TRACER(
172     I bi,bj,iMin,iMax,jMin,jMax,k,
173     I PTRACERS_advScheme(iTracer),
174     I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
175     I gPtr(1-Olx,1-Oly,1,1,1,iTracer),
176     I myIter,myThid )
177    
178 mlosch 1.25 #ifdef ALLOW_OBCS
179     C Apply open boundary conditions
180     IF (useOBCS) THEN
181     CALL OBCS_APPLY_PTRACER(
182     I bi, bj, k, iTracer,
183     U gPtr(1-Olx,1-Oly,k,bi,bj,iTracer),
184     I myThid )
185     END IF
186     #endif /* ALLOW_OBCS */
187 adcroft 1.1 C end of tracer loop
188     ENDDO
189    
190     #endif /* ALLOW_PTRACERS */
191    
192     RETURN
193     END

  ViewVC Help
Powered by ViewVC 1.1.22