/[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.4 - (show annotations) (download)
Thu Jun 20 23:46:22 2002 UTC (21 years, 11 months ago) by jmc
Branch: MAIN
Changes since 1.3: +8 -1 lines
apply rescaling (NONLIN_FRSURF) also to gNm1 to get a better conservation

1 C $Header: /u/gcmpack/MITgcm/pkg/ptracers/ptracers_integrate.F,v 1.3 2002/06/15 03:37:44 jmc Exp $
2 C $Name: $
3
4 #include "PTRACERS_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: PTRACERS_INTEGERATE
8
9 C !INTERFACE: ==========================================================
10 SUBROUTINE PTRACERS_INTEGERATE(
11 I bi,bj,k,
12 I xA,yA,uTrans,vTrans,rTrans,maskUp,
13 X KappaRtr,
14 I myIter,myTime,myThid )
15
16 C !DESCRIPTION:
17 C Calculates tendancy for passive tracers and integrates forward
18 C in time.
19
20 C !USES: ===============================================================
21 IMPLICIT NONE
22 #include "SIZE.h"
23 #include "EEPARAMS.h"
24 #include "PTRACERS.h"
25 #include "GAD.h"
26
27 C !INPUT PARAMETERS: ===================================================
28 C bi,bj :: tile indices
29 C k :: vertical level number
30 C xA :: face area at U points in level k
31 C yA :: face area at V points in level k
32 C uTrans :: zonal transport in level k
33 C vTrans :: meridional transport in level k
34 C rTrans :: vertical transport across level k
35 C maskUp :: mask for vertical transport
36 C KappaRtr :: vertical diffusion of passive tracers
37 C NOTE! This is infact KappaRS from thermodynamics()
38 C and is being used only temporarily
39 C until we removed the need to store
40 C a "3D" diffusivity
41 C myIter :: time-step number
42 C myTime :: model time
43 C myThid :: thread number
44 INTEGER bi,bj,k
45 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
46 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
47 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
48 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
49 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
50 _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
51 _RL KappaRtr(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
52 INTEGER myIter
53 _RL myTime
54 INTEGER myThid
55
56 C !OUTPUT PARAMETERS: ==================================================
57 C none
58
59 #ifdef ALLOW_PTRACERS
60
61 C !LOCAL VARIABLES: ====================================================
62 C i,j,k,bi,bj,iTracer :: loop indices
63 C iMin,iMax,jMin,jMax :: loop ranges
64 C kUp,kDown :: toggle indices for even/odd level fluxes
65 C km1 :: =min(1,k-1)
66 C rFlx :: vertical flux
67 INTEGER i,j,iTracer
68 INTEGER iMin,iMax,jMin,jMax
69 INTEGER kUp,kDown,km1
70 _RL rFlx(1-OLx:sNx+OLx,1-OLy:sNy+OLy,2,PTRACERS_num)
71 LOGICAL calcAdvection
72 CEOP
73
74 C Loop over tracers
75 DO iTracer=1,PTRACERS_numInUse
76
77 C Initialize vertical flux to zero and set no-flux across k=Nr+1
78 IF (k.EQ.Nr) THEN
79 DO j=1-Oly,sNy+Oly
80 DO i=1-Olx,sNx+Olx
81 rFlx(i,j,1,iTracer)=0.
82 rFlx(i,j,2,iTracer)=0.
83 ENDDO
84 ENDDO
85 ENDIF
86
87 C Loop ranges for daughter routines
88 iMin = 1-OLx+2
89 iMax = sNx+OLx-1
90 jMin = 1-OLy+2
91 jMax = sNy+OLy-1
92
93 km1 = MAX(1,k-1)
94 kUp = 1+MOD(k+1,2)
95 kDown= 1+MOD(k,2)
96
97 C Calculate active tracer tendencies (gPtr) due to internal processes
98 C (advection, [explicit] diffusion, parameterizations,...)
99 calcAdvection = .NOT.multiDimAdvection
100 & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
101 & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
102 & .OR. PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH
103 CALL GAD_CALC_RHS(
104 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
105 I xA,yA,uTrans,vTrans,rTrans,maskUp,
106 I PTRACERS_diffKh(iTracer),
107 I PTRACERS_diffK4(iTracer),
108 I KappaRtr,
109 I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
110 I GAD_TR1,
111 I PTRACERS_advScheme(iTracer),calcAdvection,
112 U rFlx(1-Olx,1-Oly,1,iTracer),
113 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
114 I myThid )
115
116 C External forcing term(s)
117 IF ( forcing_In_AB )
118 & CALL PTRACERS_FORCING(
119 I bi,bj,k,iTracer,
120 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
121 I myIter,myTime,myThid)
122
123 C If using Adams-Bashforth II, then extrapolate tendancies
124 IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
125 & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
126 & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH ) THEN
127 CALL ADAMS_BASHFORTH2(
128 I bi,bj,K,
129 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
130 U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
131 I myIter,myThid )
132 ENDIF
133
134 C External forcing term(s)
135 IF ( .NOT.forcing_In_AB )
136 & CALL PTRACERS_FORCING(
137 I bi,bj,k,iTracer,
138 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
139 I myIter,myTime,myThid)
140
141 #ifdef NONLIN_FRSURF
142 C Account for change in level thickness
143 IF (nonlinFreeSurf.GT.0) THEN
144 CALL FREESURF_RESCALE_G(
145 I bi,bj,K,
146 U gPtr(1-Olx,1-Oly,1,1,1,iTracer),
147 I myThid )
148 IF ( PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_2ND
149 & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_UPWIND_3RD
150 & .OR.PTRACERS_advScheme(iTracer).EQ.ENUM_CENTERED_4TH )
151 & CALL FREESURF_RESCALE_G(
152 I bi,bj,K,
153 U gPtrNm1(1-Olx,1-Oly,1,1,1,iTracer),
154 I myThid )
155 ENDIF
156 #endif /* NONLIN_FRSURF */
157
158 C Integrate forward in time, storing in gPtr: G=T+dt*G
159 CALL TIMESTEP_TRACER(
160 I bi,bj,iMin,iMax,jMin,jMax,k,
161 I PTRACERS_advScheme(iTracer),
162 I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
163 I gPtr(1-Olx,1-Oly,1,1,1,iTracer),
164 I myIter,myThid )
165
166 C end of tracer loop
167 ENDDO
168
169 #endif /* ALLOW_PTRACERS */
170
171 RETURN
172 END

  ViewVC Help
Powered by ViewVC 1.1.22