/[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.36 - (show annotations) (download)
Sat Jun 26 02:38:09 2004 UTC (19 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint54d_post, checkpoint54e_post, checkpoint55, checkpoint54, checkpoint54f_post, checkpoint55c_post, checkpoint54b_post, checkpoint54a_pre, checkpoint54a_post, checkpoint55b_post, checkpoint53g_post, checkpoint55a_post, checkpoint54c_post
Changes since 1.35: +2 -2 lines
T & S: separate Vert.Advec.Scheme from horizontal Advec.Scheme

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

  ViewVC Help
Powered by ViewVC 1.1.22