/[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.10 - (show annotations) (download)
Sun Jun 17 13:13:08 2012 UTC (11 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint64, checkpoint65, checkpoint66a, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63o, checkpoint65o, checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f
Changes since 1.9: +6 -3 lines
for TAF: (re-)initialise aC2d to zero

1 C $Header: /u/gcmpack/MITgcm/model/src/update_cg2d.F,v 1.9 2012/06/15 21:54:51 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
33 C !INPUT/OUTPUT PARAMETERS:
34 C === Routine arguments ===
35 C myTime :: Current time in simulation
36 C myIter :: Current iteration number in simulation
37 C myThid :: Thread number for this instance of the routine.
38 _RL myTime
39 INTEGER myIter
40 INTEGER myThid
41
42 C !LOCAL VARIABLES:
43 C-- Note : compared to "INI_CG2D", no needs to compute again
44 C the solver normalisation factor or the solver tolerance
45 C === Local variables ===
46 C bi,bj :: tile indices
47 C i,j,k :: Loop counters
48 C faceArea :: Temporary used to hold cell face areas.
49 INTEGER bi, bj
50 INTEGER i, j, k, ks
51 _RL faceArea
52 _RL pW_tmp, pS_tmp
53 LOGICAL updatePreCond
54 CEOP
55
56 C-- Decide when to update cg2d Preconditioner :
57 IF ( cg2dPreCondFreq.EQ.0 ) THEN
58 updatePreCond = .FALSE.
59 ELSE
60 updatePreCond = ( myIter.EQ.nIter0 )
61 IF ( MOD(myIter,cg2dPreCondFreq).EQ.0 ) updatePreCond=.TRUE.
62 ENDIF
63
64 C-- Initialise laplace operator
65 C aW2d: integral in Z Ax/dX
66 C aS2d: integral in Z Ay/dY
67 DO bj=myByLo(myThid),myByHi(myThid)
68 DO bi=myBxLo(myThid),myBxHi(myThid)
69 DO j=1-OLy,sNy+OLy
70 DO i=1-OLx,sNx+OLx
71 aW2d(i,j,bi,bj) = 0. _d 0
72 aS2d(i,j,bi,bj) = 0. _d 0
73 #ifdef ALLOW_AUTODIFF
74 aC2d(i,j,bi,bj) = 0. _d 0
75 #endif /* ALLOW_AUTODIFF */
76 ENDDO
77 ENDDO
78 DO k=1,Nr
79 DO j=1,sNy+1
80 DO i=1,sNx+1
81 C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect
82 faceArea = _dyG(i,j,bi,bj)*drF(k)
83 & *_hFacW(i,j,k,bi,bj)
84 aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)
85 & + faceArea*recip_dxC(i,j,bi,bj)
86 faceArea = _dxG(i,j,bi,bj)*drF(k)
87 & *_hFacS(i,j,k,bi,bj)
88 aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)
89 & + faceArea*recip_dyC(i,j,bi,bj)
90 ENDDO
91 ENDDO
92 ENDDO
93 DO j=1,sNy+1
94 DO i=1,sNx+1
95 aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)*cg2dNorm
96 & *implicSurfPress*implicDiv2DFlow
97 #ifdef ALLOW_OBCS
98 & *maskInC(i,j,bi,bj)*maskInC(i-1,j,bi,bj)
99 #endif
100 aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)*cg2dNorm
101 & *implicSurfPress*implicDiv2DFlow
102 #ifdef ALLOW_OBCS
103 & *maskInC(i,j,bi,bj)*maskInC(i,j-1,bi,bj)
104 #endif
105 ENDDO
106 ENDDO
107 C-- compute matrix main diagonal :
108 IF ( deepAtmosphere ) THEN
109 DO j=1,sNy
110 DO i=1,sNx
111 ks = kSurfC(i,j,bi,bj)
112 aC2d(i,j,bi,bj) = -(
113 & aW2d(i,j,bi,bj) + aW2d(i+1, j ,bi,bj)
114 & +aS2d(i,j,bi,bj) + aS2d( i ,j+1,bi,bj)
115 & +freeSurfFac*cg2dNorm*recip_Bo(i,j,bi,bj)*deepFac2F(ks)
116 & *rA(i,j,bi,bj)/deltaTMom/deltaTfreesurf
117 & )
118 ENDDO
119 ENDDO
120 ELSE
121 DO j=1,sNy
122 DO i=1,sNx
123 aC2d(i,j,bi,bj) = -(
124 & aW2d(i,j,bi,bj) + aW2d(i+1, j ,bi,bj)
125 & +aS2d(i,j,bi,bj) + aS2d( i ,j+1,bi,bj)
126 & +freeSurfFac*cg2dNorm*recip_Bo(i,j,bi,bj)
127 & *rA(i,j,bi,bj)/deltaTMom/deltaTfreesurf
128 & )
129 ENDDO
130 ENDDO
131 ENDIF
132 C- end bi,bj loops
133 ENDDO
134 ENDDO
135
136 IF ( updatePreCond ) THEN
137 C-- Update overlap regions
138 CALL EXCH_XY_RS(aC2d, myThid)
139
140 C-- Initialise preconditioner
141 DO bj=myByLo(myThid),myByHi(myThid)
142 DO bi=myBxLo(myThid),myBxHi(myThid)
143 DO j=1,sNy+1
144 DO i=1,sNx+1
145 IF ( aC2d(i,j,bi,bj) .EQ. 0. ) THEN
146 pC(i,j,bi,bj) = 1. _d 0
147 ELSE
148 pC(i,j,bi,bj) = 1. _d 0 / aC2d(i,j,bi,bj)
149 ENDIF
150 pW_tmp = aC2d(i,j,bi,bj)+aC2d(i-1,j,bi,bj)
151 IF ( pW_tmp .EQ. 0. ) THEN
152 pW(i,j,bi,bj) = 0.
153 ELSE
154 pW(i,j,bi,bj) =
155 & -aW2d(i,j,bi,bj)/((cg2dpcOffDFac *pW_tmp)**2 )
156 ENDIF
157 pS_tmp = aC2d(i,j,bi,bj)+aC2d(i,j-1,bi,bj)
158 IF ( pS_tmp .EQ. 0. ) THEN
159 pS(i,j,bi,bj) = 0.
160 ELSE
161 pS(i,j,bi,bj) =
162 & -aS2d(i,j,bi,bj)/((cg2dpcOffDFac *pS_tmp)**2 )
163 ENDIF
164 ENDDO
165 ENDDO
166 ENDDO
167 ENDDO
168 C- if update Preconditioner : end
169 ENDIF
170
171 RETURN
172 END

  ViewVC Help
Powered by ViewVC 1.1.22