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

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

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


Revision 1.11 - (show annotations) (download)
Mon Nov 28 23:05:05 2016 UTC (7 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, HEAD
Changes since 1.10: +41 -15 lines
implement fully implicit bottom friction coupled with implicit surface
  pressure (hydrostatic only)

1 C $Header: /u/gcmpack/MITgcm/model/src/update_cg2d.F,v 1.10 2012/06/17 13:13:08 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: UPDATE_CG2D
9 C !INTERFACE:
10 SUBROUTINE UPDATE_CG2D( myTime, myIter, myThid )
11 C !DESCRIPTION: \bv
12 C *==========================================================*
13 C | SUBROUTINE UPDATE_CG2D
14 C | o Update 2d conjugate gradient solver operators
15 C | account for Free-Surf effect on total column thickness
16 C *==========================================================*
17 C | This routine is based on INI_CG2D, and simplified. It is
18 C | used when the non-linear free surface mode is activated
19 C | or when bottom depth is part of the control vector.
20 C *==========================================================*
21 C \ev
22
23 C !USES:
24 IMPLICIT NONE
25 C === Global variables ===
26 #include "SIZE.h"
27 #include "EEPARAMS.h"
28 #include "PARAMS.h"
29 #include "GRID.h"
30 #include "SURFACE.h"
31 #include "CG2D.h"
32 #ifdef ALLOW_SOLVE4_PS_AND_DRAG
33 # include "DYNVARS.h"
34 #endif /* ALLOW_SOLVE4_PS_AND_DRAG */
35
36 C !INPUT/OUTPUT PARAMETERS:
37 C === Routine arguments ===
38 C myTime :: Current time in simulation
39 C myIter :: Current iteration number in simulation
40 C myThid :: Thread number for this instance of the routine.
41 _RL myTime
42 INTEGER myIter
43 INTEGER myThid
44
45 C !LOCAL VARIABLES:
46 C-- Note : compared to "INI_CG2D", no needs to compute again
47 C the solver normalisation factor or the solver tolerance
48 C === Local variables ===
49 C bi,bj :: tile indices
50 C i,j,k :: Loop counters
51 C faceArea :: Temporary used to hold cell face areas.
52 INTEGER bi, bj
53 INTEGER i, j, k, ks
54 _RL faceArea
55 _RL pW_tmp, pS_tmp
56 LOGICAL updatePreCond
57 CEOP
58
59 C-- Decide when to update cg2d Preconditioner :
60 IF ( cg2dPreCondFreq.EQ.0 ) THEN
61 updatePreCond = .FALSE.
62 ELSE
63 updatePreCond = ( myIter.EQ.nIter0 )
64 IF ( MOD(myIter,cg2dPreCondFreq).EQ.0 ) updatePreCond=.TRUE.
65 ENDIF
66
67 C-- Initialise laplace operator
68 C aW2d: integral in Z Ax/dX
69 C aS2d: integral in Z Ay/dY
70 DO bj=myByLo(myThid),myByHi(myThid)
71 DO bi=myBxLo(myThid),myBxHi(myThid)
72 DO j=1-OLy,sNy+OLy
73 DO i=1-OLx,sNx+OLx
74 aW2d(i,j,bi,bj) = 0. _d 0
75 aS2d(i,j,bi,bj) = 0. _d 0
76 #ifdef ALLOW_AUTODIFF
77 aC2d(i,j,bi,bj) = 0. _d 0
78 #endif /* ALLOW_AUTODIFF */
79 ENDDO
80 ENDDO
81 #ifdef ALLOW_SOLVE4_PS_AND_DRAG
82 IF ( selectImplicitDrag.EQ.2 ) THEN
83 DO k=1,Nr
84 DO j=1,sNy+1
85 DO i=1,sNx+1
86 faceArea = _dyG(i,j,bi,bj)*drF(k)*deepFacC(k)*rhoFacC(k)
87 & *_hFacW(i,j,k,bi,bj)
88 aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)
89 & + faceArea*dU_psFacX(i,j,k,bi,bj)
90 & *recip_dxC(i,j,bi,bj)
91 faceArea = _dxG(i,j,bi,bj)*drF(k)*deepFacC(k)*rhoFacC(k)
92 & *_hFacS(i,j,k,bi,bj)
93 aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)
94 & + faceArea*dV_psFacY(i,j,k,bi,bj)
95 & *recip_dyC(i,j,bi,bj)
96 ENDDO
97 ENDDO
98 ENDDO
99 ELSE
100 #endif /* ALLOW_SOLVE4_PS_AND_DRAG */
101 DO k=1,Nr
102 DO j=1,sNy+1
103 DO i=1,sNx+1
104 C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect
105 faceArea = _dyG(i,j,bi,bj)*drF(k)
106 & *_hFacW(i,j,k,bi,bj)
107 aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)
108 & + faceArea*recip_dxC(i,j,bi,bj)
109 faceArea = _dxG(i,j,bi,bj)*drF(k)
110 & *_hFacS(i,j,k,bi,bj)
111 aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)
112 & + faceArea*recip_dyC(i,j,bi,bj)
113 ENDDO
114 ENDDO
115 ENDDO
116 #ifdef ALLOW_SOLVE4_PS_AND_DRAG
117 ENDIF
118 #endif /* ALLOW_SOLVE4_PS_AND_DRAG */
119 DO j=1,sNy+1
120 DO i=1,sNx+1
121 aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)*cg2dNorm
122 & *implicSurfPress*implicDiv2DFlow
123 #ifdef ALLOW_OBCS
124 & *maskInC(i,j,bi,bj)*maskInC(i-1,j,bi,bj)
125 #endif
126 aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)*cg2dNorm
127 & *implicSurfPress*implicDiv2DFlow
128 #ifdef ALLOW_OBCS
129 & *maskInC(i,j,bi,bj)*maskInC(i,j-1,bi,bj)
130 #endif
131 ENDDO
132 ENDDO
133 C-- compute matrix main diagonal :
134 IF ( deepAtmosphere ) THEN
135 DO j=1,sNy
136 DO i=1,sNx
137 ks = kSurfC(i,j,bi,bj)
138 aC2d(i,j,bi,bj) = -(
139 & aW2d(i,j,bi,bj) + aW2d(i+1, j ,bi,bj)
140 & +aS2d(i,j,bi,bj) + aS2d( i ,j+1,bi,bj)
141 & +freeSurfFac*cg2dNorm*recip_Bo(i,j,bi,bj)*deepFac2F(ks)
142 & *rA(i,j,bi,bj)/deltaTMom/deltaTFreeSurf
143 & )
144 ENDDO
145 ENDDO
146 ELSE
147 DO j=1,sNy
148 DO i=1,sNx
149 aC2d(i,j,bi,bj) = -(
150 & aW2d(i,j,bi,bj) + aW2d(i+1, j ,bi,bj)
151 & +aS2d(i,j,bi,bj) + aS2d( i ,j+1,bi,bj)
152 & +freeSurfFac*cg2dNorm*recip_Bo(i,j,bi,bj)
153 & *rA(i,j,bi,bj)/deltaTMom/deltaTFreeSurf
154 & )
155 ENDDO
156 ENDDO
157 ENDIF
158 C- end bi,bj loops
159 ENDDO
160 ENDDO
161
162 IF ( updatePreCond ) THEN
163 C-- Update overlap regions
164 CALL EXCH_XY_RS(aC2d, myThid)
165
166 C-- Initialise preconditioner
167 DO bj=myByLo(myThid),myByHi(myThid)
168 DO bi=myBxLo(myThid),myBxHi(myThid)
169 DO j=1,sNy+1
170 DO i=1,sNx+1
171 IF ( aC2d(i,j,bi,bj) .EQ. 0. ) THEN
172 pC(i,j,bi,bj) = 1. _d 0
173 ELSE
174 pC(i,j,bi,bj) = 1. _d 0 / aC2d(i,j,bi,bj)
175 ENDIF
176 pW_tmp = aC2d(i,j,bi,bj)+aC2d(i-1,j,bi,bj)
177 IF ( pW_tmp .EQ. 0. ) THEN
178 pW(i,j,bi,bj) = 0.
179 ELSE
180 pW(i,j,bi,bj) =
181 & -aW2d(i,j,bi,bj)/((cg2dpcOffDFac *pW_tmp)**2 )
182 ENDIF
183 pS_tmp = aC2d(i,j,bi,bj)+aC2d(i,j-1,bi,bj)
184 IF ( pS_tmp .EQ. 0. ) THEN
185 pS(i,j,bi,bj) = 0.
186 ELSE
187 pS(i,j,bi,bj) =
188 & -aS2d(i,j,bi,bj)/((cg2dpcOffDFac *pS_tmp)**2 )
189 ENDIF
190 ENDDO
191 ENDDO
192 ENDDO
193 ENDDO
194 C- if update Preconditioner : end
195 ENDIF
196
197 RETURN
198 END

  ViewVC Help
Powered by ViewVC 1.1.22