/[MITgcm]/MITgcm/model/src/calc_gs.F
ViewVC logotype

Contents of /MITgcm/model/src/calc_gs.F

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


Revision 1.39 - (show annotations) (download)
Fri Dec 3 15:39:11 2004 UTC (19 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57g_pre, checkpoint57b_post, checkpoint57d_post, checkpoint57, checkpoint57f_post, checkpoint57a_post, checkpoint57c_post, checkpoint57c_pre, checkpoint57e_post, eckpoint57e_pre, checkpoint57f_pre, checkpoint57a_pre
Changes since 1.38: +8 -1 lines
allow to compile without generic_advdiff pkg.

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_gs.F,v 1.38 2004/10/21 21:25:33 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: CALC_GS
9 C !INTERFACE:
10 SUBROUTINE CALC_GS(
11 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
12 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
13 I KappaRS,
14 U fVerS,
15 I myTime,myIter,myThid )
16 C !DESCRIPTION: \bv
17 C *==========================================================*
18 C | SUBROUTINE CALC_GS
19 C | o Calculate the salt tendency terms.
20 C *==========================================================*
21 C | A procedure called EXTERNAL_FORCING_S is called from
22 C | here. These procedures can be used to add per problem
23 C | E-P flux source terms.
24 C | Note: Although it is slightly counter-intuitive the
25 C | EXTERNAL_FORCING routine is not the place to put
26 C | file I/O. Instead files that are required to
27 C | calculate the external source terms are generally
28 C | read during the model main loop. This makes the
29 C | logisitics of multi-processing simpler and also
30 C | makes the adjoint generation simpler. It also
31 C | allows for I/O to overlap computation where that
32 C | is supported by hardware.
33 C | Aside from the problem specific term the code here
34 C | forms the tendency terms due to advection and mixing
35 C | The baseline implementation here uses a centered
36 C | difference form for the advection term and a tensorial
37 C | divergence of a flux form for the diffusive term. The
38 C | diffusive term is formulated so that isopycnal mixing and
39 C | GM-style subgrid-scale terms can be incorporated b simply
40 C | setting the diffusion tensor terms appropriately.
41 C *==========================================================*
42 C \ev
43
44 C !USES:
45 IMPLICIT NONE
46 C == GLobal variables ==
47 #include "SIZE.h"
48 #include "DYNVARS.h"
49 #include "EEPARAMS.h"
50 #include "PARAMS.h"
51 #include "PACKAGES_CONFIG.h"
52 #ifdef ALLOW_GENERIC_ADVDIFF
53 #include "GAD.h"
54 #endif
55
56 C !INPUT/OUTPUT PARAMETERS:
57 C == Routine arguments ==
58 C fVerS :: Flux of salt (S) in the vertical
59 C direction at the upper(U) and lower(D) faces of a cell.
60 C maskUp :: Land mask used to denote base of the domain.
61 C xA :: Tracer cell face area normal to X
62 C yA :: Tracer cell face area normal to X
63 C uTrans :: Zonal volume transport through cell face
64 C vTrans :: Meridional volume transport through cell face
65 C rTrans :: Vertical volume transport at interface k
66 C rTransKp1 :: Vertical volume transport at inteface k+1
67 C bi, bj, iMin, iMax, jMin, jMax :: Range of points for which calculation
68 C results will be set.
69 C myThid :: Instance number for this innvocation of CALC_GT
70 _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
71 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
72 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78 _RL KappaRS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
79 INTEGER k,kUp,kDown,kM1
80 INTEGER bi,bj,iMin,iMax,jMin,jMax
81 _RL myTime
82 INTEGER myIter
83 INTEGER myThid
84
85 CEOP
86
87 #ifdef ALLOW_GENERIC_ADVDIFF
88 C === Local variables ===
89 LOGICAL calcAdvection
90
91 #ifdef ALLOW_AUTODIFF_TAMC
92 C-- only the kUp part of fverS is set in this subroutine
93 C-- the kDown is still required
94 fVerS(1,1,kDown) = fVerS(1,1,kDown)
95 #endif
96
97 calcAdvection = saltAdvection .AND. .NOT.saltMultiDimAdvec
98 CALL GAD_CALC_RHS(
99 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
100 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
101 I uVel, vVel, wVel,
102 I diffKhS, diffK4S, KappaRS, Salt,
103 I GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
104 I calcAdvection, saltImplVertAdv,
105 U fVerS, gS,
106 I myTime, myIter, myThid )
107
108 C-- External salinity forcing term(s) inside Adams-Bashforth:
109 IF ( saltForcing .AND. forcing_In_AB )
110 & CALL EXTERNAL_FORCING_S(
111 I iMin,iMax,jMin,jMax,bi,bj,k,
112 I myTime,myThid)
113
114 IF ( saltAdamsBashforth ) THEN
115 CALL ADAMS_BASHFORTH2(
116 I bi, bj, K,
117 U gS, gSnm1,
118 I myIter, myThid )
119 ENDIF
120
121 C-- External salinity forcing term(s) outside Adams-Bashforth:
122 IF ( saltForcing .AND. .NOT.forcing_In_AB )
123 & CALL EXTERNAL_FORCING_S(
124 I iMin,iMax,jMin,jMax,bi,bj,k,
125 I myTime,myThid)
126
127 #ifdef NONLIN_FRSURF
128 IF (nonlinFreeSurf.GT.0) THEN
129 CALL FREESURF_RESCALE_G(
130 I bi, bj, K,
131 U gS,
132 I myThid )
133 IF ( saltAdamsBashforth )
134 & CALL FREESURF_RESCALE_G(
135 I bi, bj, K,
136 U gSnm1,
137 I myThid )
138 ENDIF
139 #endif /* NONLIN_FRSURF */
140
141 #endif /* ALLOW_GENERIC_ADVDIFF */
142
143 RETURN
144 END

  ViewVC Help
Powered by ViewVC 1.1.22