/[MITgcm]/MITgcm_contrib/ksnow/press_release/code/update_cg2d.F
ViewVC logotype

Annotation of /MITgcm_contrib/ksnow/press_release/code/update_cg2d.F

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


Revision 1.1 - (hide annotations) (download)
Fri Dec 16 15:23:18 2016 UTC (8 years, 7 months ago) by ksnow
Branch: MAIN
Adding press_release core code files
C: ----------------------------------------------------------------------

1 ksnow 1.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     #ifdef ALLOW_PRESSURE_RELEASE_CODE
7     # include "SHELFICE_OPTIONS.h"
8     #endif
9    
10    
11     CBOP
12     C !ROUTINE: UPDATE_CG2D
13     C !INTERFACE:
14     SUBROUTINE UPDATE_CG2D( myTime, myIter, myThid )
15     C !DESCRIPTION: \bv
16     C *==========================================================*
17     C | SUBROUTINE UPDATE_CG2D
18     C | o Update 2d conjugate gradient solver operators
19     C | account for Free-Surf effect on total column thickness
20     C *==========================================================*
21     C | This routine is based on INI_CG2D, and simplified. It is
22     C | used when the non-linear free surface mode is activated
23     C | or when bottom depth is part of the control vector.
24     C *==========================================================*
25     C \ev
26    
27     C !USES:
28     IMPLICIT NONE
29     C === Global variables ===
30     #include "SIZE.h"
31     #include "EEPARAMS.h"
32     #include "PARAMS.h"
33     #include "GRID.h"
34     #include "SURFACE.h"
35     #include "CG2D.h"
36     #ifdef ALLOW_PRESSURE_RELEASE_CODE
37     # include "SHELFICE.h"
38     #endif
39     #ifdef IMPLICIT_BOTTOMSIDEDRAG
40     # include "MOM_VISC.h"
41     #endif
42    
43    
44     C !INPUT/OUTPUT PARAMETERS:
45     C === Routine arguments ===
46     C myTime :: Current time in simulation
47     C myIter :: Current iteration number in simulation
48     C myThid :: Thread number for this instance of the routine.
49     _RL myTime
50     INTEGER myIter
51     INTEGER myThid
52    
53     C !LOCAL VARIABLES:
54     C-- Note : compared to "INI_CG2D", no needs to compute again
55     C the solver normalisation factor or the solver tolerance
56     C === Local variables ===
57     C bi,bj :: tile indices
58     C i,j,k :: Loop counters
59     C faceArea :: Temporary used to hold cell face areas.
60     INTEGER bi, bj
61     INTEGER i, j, k, ks
62     _RL faceArea
63     _RL pW_tmp, pS_tmp
64     _RL depth, eff_depth, depth_fac, pinit
65     #ifdef ALLOW_PRESSURE_RELEASE_CODE
66     _RL curPrelVisc, avgGrd
67     #endif
68     LOGICAL updatePreCond
69     CEOP
70    
71     C-- Decide when to update cg2d Preconditioner :
72     IF ( cg2dPreCondFreq.EQ.0 ) THEN
73     updatePreCond = .FALSE.
74     ELSE
75     updatePreCond = ( myIter.EQ.nIter0 )
76     IF ( MOD(myIter,cg2dPreCondFreq).EQ.0 ) updatePreCond=.TRUE.
77     ENDIF
78    
79     C-- Initialise laplace operator
80     C aW2d: integral in Z Ax/dX
81     C aS2d: integral in Z Ay/dY
82     DO bj=myByLo(myThid),myByHi(myThid)
83     DO bi=myBxLo(myThid),myBxHi(myThid)
84     DO j=1-OLy,sNy+OLy
85     DO i=1-OLx,sNx+OLx
86     aW2d(i,j,bi,bj) = 0. _d 0
87     aS2d(i,j,bi,bj) = 0. _d 0
88     #ifdef ALLOW_AUTODIFF
89     aC2d(i,j,bi,bj) = 0. _d 0
90     #endif /* ALLOW_AUTODIFF */
91     ENDDO
92     ENDDO
93     DO k=1,Nr
94     DO j=1,sNy+1
95     DO i=1,sNx+1
96     C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect
97     faceArea = _dyG(i,j,bi,bj)*drF(k)
98     & *_hFacW(i,j,k,bi,bj)
99     aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)
100     #ifndef IMPLICIT_BOTTOMSIDEDRAG
101     & + faceArea*recip_dxC(i,j,bi,bj)
102     #else
103     & + faceArea/
104     & (1 - uDragTermsCommon(i,j,k,bi,bj)*deltaTmom) *
105     & recip_dxC(i,j,bi,bj)
106     #endif
107    
108     faceArea = _dxG(i,j,bi,bj)*drF(k)
109     & *_hFacS(i,j,k,bi,bj)
110     aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)
111     #ifndef IMPLICIT_BOTTOMSIDEDRAG
112     & + faceArea*recip_dyC(i,j,bi,bj)
113     #else
114     & + faceArea/
115     & (1 - vDragTermsCommon(i,j,k,bi,bj)*deltaTmom) *
116     & recip_dyC(i,j,bi,bj)
117     #endif
118     ENDDO
119     ENDDO
120     ENDDO
121     DO j=1,sNy+1
122     DO i=1,sNx+1
123     aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)*cg2dNorm
124     & *implicSurfPress*implicDiv2DFlow
125     #ifdef ALLOW_OBCS
126     & *maskInC(i,j,bi,bj)*maskInC(i-1,j,bi,bj)
127     #endif
128     #ifdef ALLOW_PRESSURE_RELEASE_CODE
129     &
130     #endif
131     aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)*cg2dNorm
132     & *implicSurfPress*implicDiv2DFlow
133     #ifdef ALLOW_OBCS
134     & *maskInC(i,j,bi,bj)*maskInC(i,j-1,bi,bj)
135     #endif
136     ENDDO
137     ENDDO
138     ! call write_fld_xy_rl('as2d0','',as2d,0,mythid)
139     #ifdef ALLOW_PRESSURE_RELEASE_CODE
140     DO j=1,sNy+1
141     DO i=1,sNx+1
142     IF (depthColW(i,j,bi,bj).lt.cg2dminColumnEps) THEN
143     IF (maskInC(i,j,bi,bj)*maskInC(i-1,j,bi,bj).eq.1.) THEN
144    
145     avgGrd = min (grdFactor(i,j,bi,bj),grdFactor(i-1,j,bi,bj))
146     curPrelVisc = max((0.55+avgGrd*0.45),0. _d 0)*pReleaseVisc
147    
148     depth = depthColW(i,j,bi,bj)
149     eff_depth = cg2dminColumnEps *
150     & (1.-2./PI * COS (PI*depth / (2.*cg2dminColumnEps)))
151     aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj) +
152     & curPrelVisc*(eff_depth-depth)
153     & *_dyG(i,j,bi,bj)*
154     & recip_dxC(i,j,bi,bj) * cg2dNorm *
155     & implicSurfPress*implicDiv2DFlow /
156     & deltaTmom
157    
158     ENDIF
159     ENDIF
160     IF (depthColS(i,j,bi,bj).lt.cg2dminColumnEps) THEN
161     IF (maskInC(i,j,bi,bj)*maskInC(i,j-1,bi,bj).eq.1.) THEN
162    
163     avgGrd = min (grdFactor(i,j,bi,bj),grdFactor(i,j-1,bi,bj))
164     curPrelVisc = max((0.55+avgGrd*0.45),0. _d 0)*pReleaseVisc
165    
166     depth = depthColS(i,j,bi,bj)
167     eff_depth = cg2dminColumnEps *
168     & (1.-2./PI * COS (PI*depth / (2.*cg2dminColumnEps)))
169     aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj) +
170     & curPrelVisc*(eff_depth-depth)
171     & *_dxG(i,j,bi,bj)*
172     & recip_dyC(i,j,bi,bj) * cg2dNorm *
173     & implicSurfPress*implicDiv2DFlow /
174     & deltaTmom
175     ENDIF
176     ENDIF
177     ENDDO
178     ENDDO
179     #endif
180     ! call write_fld_xy_rl('as2d0','',as2d,0,mythid)
181     C-- compute matrix main diagonal :
182     IF ( deepAtmosphere ) THEN
183     DO j=1,sNy
184     DO i=1,sNx
185     ks = kSurfC(i,j,bi,bj)
186     aC2d(i,j,bi,bj) = -(
187     & aW2d(i,j,bi,bj) + aW2d(i+1, j ,bi,bj)
188     & +aS2d(i,j,bi,bj) + aS2d( i ,j+1,bi,bj)
189     & +freeSurfFac*cg2dNorm*recip_Bo(i,j,bi,bj)*deepFac2F(ks)
190     & *rA(i,j,bi,bj)/deltaTMom/deltaTfreesurf
191     & )
192     ENDDO
193     ENDDO
194     ELSE
195     DO j=1,sNy
196     DO i=1,sNx
197     aC2d(i,j,bi,bj) = -(
198     & aW2d(i,j,bi,bj) + aW2d(i+1, j ,bi,bj)
199     & +aS2d(i,j,bi,bj) + aS2d( i ,j+1,bi,bj)
200     & +freeSurfFac*cg2dNorm*recip_Bo(i,j,bi,bj)
201     & *rA(i,j,bi,bj)/deltaTMom/deltaTfreesurf
202     & )
203     ENDDO
204     ENDDO
205     ENDIF
206     C- end bi,bj loops
207     ENDDO
208     ENDDO
209    
210     IF ( updatePreCond ) THEN
211     C-- Update overlap regions
212     CALL EXCH_XY_RS(aC2d, myThid)
213    
214     C-- Initialise preconditioner
215     DO bj=myByLo(myThid),myByHi(myThid)
216     DO bi=myBxLo(myThid),myBxHi(myThid)
217     DO j=1,sNy+1
218     DO i=1,sNx+1
219     IF ( aC2d(i,j,bi,bj) .EQ. 0. ) THEN
220     pC(i,j,bi,bj) = 1. _d 0
221     ELSE
222     pC(i,j,bi,bj) = 1. _d 0 / aC2d(i,j,bi,bj)
223     ENDIF
224     pW_tmp = aC2d(i,j,bi,bj)+aC2d(i-1,j,bi,bj)
225     IF ( pW_tmp .EQ. 0. ) THEN
226     pW(i,j,bi,bj) = 0.
227     ELSE
228     pW(i,j,bi,bj) =
229     & -aW2d(i,j,bi,bj)/((cg2dpcOffDFac *pW_tmp)**2 )
230     ENDIF
231     pS_tmp = aC2d(i,j,bi,bj)+aC2d(i,j-1,bi,bj)
232     IF ( pS_tmp .EQ. 0. ) THEN
233     pS(i,j,bi,bj) = 0.
234     ELSE
235     pS(i,j,bi,bj) =
236     & -aS2d(i,j,bi,bj)/((cg2dpcOffDFac *pS_tmp)**2 )
237     ENDIF
238     ENDDO
239     ENDDO
240     ENDDO
241     ENDDO
242     C- if update Preconditioner : end
243     ENDIF
244    
245    
246     RETURN
247     END

  ViewVC Help
Powered by ViewVC 1.1.22