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

Contents of /MITgcm/pkg/ptracers/ptracers_advection.F

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


Revision 1.12 - (show annotations) (download)
Thu May 8 19:50:08 2008 UTC (16 years ago) by jahn
Branch: MAIN
Changes since 1.11: +22 -1 lines
add second-order moment advection schemes (80 and 81);
this uses a dynamically allocated internal state data structure
(#define PTRACERS_ALLOW_DYN_STATE in PTRACERS_OPTIONS.h)
and requires a fortran 90 compiler

1 C $Header: /u/gcmpack/MITgcm/verification/tutorial_advection_in_gyre/code/ptracers_advection.F,v 1.1 2008/01/28 20:09:08 jahn Exp $
2 C $Name: $
3
4 #include "PTRACERS_OPTIONS.h"
5 #include "GAD_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: PTRACERS_ADVECTION
9
10 C !INTERFACE: ==========================================================
11 SUBROUTINE PTRACERS_ADVECTION( bi,bj,myIter,myTime,myThid )
12
13 C !DESCRIPTION:
14 C Calculates tendency for passive tracers and integrates forward
15 C in time.
16
17 C !USES: ===============================================================
18 #include "PTRACERS_MOD.h"
19 IMPLICIT NONE
20 #include "SIZE.h"
21 #include "EEPARAMS.h"
22 #include "PARAMS.h"
23 #include "DYNVARS.h"
24 #include "PTRACERS_SIZE.h"
25 #include "PTRACERS_PARAMS.h"
26 #include "PTRACERS_FIELDS.h"
27 #include "GAD.h"
28 #ifdef ALLOW_AUTODIFF_TAMC
29 # include "tamc.h"
30 # include "tamc_keys.h"
31 #endif
32
33 C !INPUT PARAMETERS: ===================================================
34 C bi,bj :: tile indices
35 C myIter :: time-step number
36 C myTime :: model time
37 C myThid :: thread number
38 INTEGER bi,bj
39 INTEGER myIter
40 _RL myTime
41 INTEGER myThid
42
43 C !OUTPUT PARAMETERS: ==================================================
44 C none
45
46 #ifdef ALLOW_PTRACERS
47
48 C !LOCAL VARIABLES: ====================================================
49 C i,j,k,bi,bj,iTracer :: loop indices
50 C iMin,iMax,jMin,jMax :: loop ranges
51 C kUp,kDown :: toggle indices for even/odd level fluxes
52 C km1 :: =min(1,k-1)
53 C rFlx :: vertical flux
54 INTEGER iTracer
55 CEOP
56
57 C Loop over tracers
58 DO iTracer=1,PTRACERS_numInUse
59
60 #ifdef ALLOW_AUTODIFF_TAMC
61 act0 = iTracer - 1
62 max0 = PTRACERS_num
63 act1 = bi - myBxLo(myThid)
64 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
65 act2 = bj - myByLo(myThid)
66 max2 = myByHi(myThid) - myByLo(myThid) + 1
67 act3 = myThid - 1
68 max3 = nTx*nTy
69 act4 = ikey_dynamics - 1
70 iptrkey = (act0 + 1)
71 & + act1*max0
72 & + act2*max0*max1
73 & + act3*max0*max1*max2
74 & + act4*max0*max1*max2*max3
75 #endif /* ALLOW_AUTODIFF_TAMC */
76
77 #ifdef ALLOW_AUTODIFF_TAMC
78 CADJ STORE pTracer(:,:,:,bi,bj,iTracer)
79 CADJ & = comlev1_bibj_ptracers, key=iptrkey, byte=isbyte
80 #endif /* ALLOW_AUTODIFF_TAMC */
81
82 #if defined(GAD_ALLOW_SOM_ADVECT) && defined(PTRACERS_ALLOW_DYN_STATE)
83 IF ( PTRACERS_SOM_Advection(iTracer) ) THEN
84 #ifdef ALLOW_DEBUG
85 IF ( debugLevel .GE. debLevB )
86 & CALL DEBUG_CALL('GAD_SOM_ADVECT',myThid)
87 #endif
88 CALL GAD_SOM_ADVECT(
89 I PTRACERS_ImplVertAdv(iTracer),
90 I PTRACERS_advScheme(iTracer),
91 I PTRACERS_advScheme(iTracer),
92 I GAD_TR1+iTracer-1,
93 I uVel, vVel, wVel,
94 I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
95 U _Ptracers_som(1-Olx,1-Oly,1,1,1,1,iTracer),
96 O gPtr(1-Olx,1-Oly,1,1,1,iTracer),
97 I bi,bj,myTime,myIter,myThid)
98 ELSEIF ( PTRACERS_MultiDimAdv(iTracer) ) THEN
99 #else /* GAD_ALLOW_SOM_ADVECT && PTRACERS_ALLOW_DYN_STATE */
100 IF ( PTRACERS_MultiDimAdv(iTracer) ) THEN
101 #endif /* GAD_ALLOW_SOM_ADVECT && PTRACERS_ALLOW_DYN_STATE */
102 CALL GAD_ADVECTION(
103 I PTRACERS_ImplVertAdv(iTracer),
104 I PTRACERS_advScheme(iTracer),
105 I PTRACERS_advScheme(iTracer),
106 I GAD_TR1+iTracer-1,
107 I uVel, vVel, wVel,
108 I pTracer(1-Olx,1-Oly,1,1,1,iTracer),
109 O gPtr(1-Olx,1-Oly,1,1,1,iTracer),
110 I bi,bj,myTime,myIter,myThid)
111 ENDIF
112
113
114 C end of tracer loop
115 ENDDO
116
117 #endif /* ALLOW_PTRACERS */
118
119 RETURN
120 END

  ViewVC Help
Powered by ViewVC 1.1.22