/[MITgcm]/MITgcm_contrib/dgoldberg/code_cg3d_petsc/ini_cg3d.F
ViewVC logotype

Annotation of /MITgcm_contrib/dgoldberg/code_cg3d_petsc/ini_cg3d.F

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


Revision 1.1 - (hide annotations) (download)
Fri Jul 1 20:40:30 2016 UTC (9 years ago) by dgoldberg
Branch: MAIN
test files for using cg3d with petsc

1 dgoldberg 1.1 C $Header: /u/gcmpack/MITgcm/model/src/ini_cg3d.F,v 1.29 2016/05/04 22:09:54 jmc Exp $
2     C $Name: $
3    
4     #include "PACKAGES_CONFIG.h"
5     #include "CPP_OPTIONS.h"
6    
7     CBOP
8     C !ROUTINE: INI_CG3D
9     C !INTERFACE:
10     SUBROUTINE INI_CG3D( myThid )
11     C !DESCRIPTION: \bv
12     C *==========================================================*
13     C | SUBROUTINE INI_CG3D
14     C | o Initialise 3d conjugate gradient solver operators.
15     C *==========================================================*
16     C | These arrays are purely a function of the basin geom.
17     C | We set then here once and them use then repeatedly.
18     C *==========================================================*
19     C \ev
20    
21     C !USES:
22     IMPLICIT NONE
23     C === Global variables ===
24     #include "SIZE.h"
25     #include "EEPARAMS.h"
26     #include "PARAMS.h"
27     #include "GRID.h"
28     #include "SURFACE.h"
29     #include "CG3D.h"
30     #include "SOLVE_FOR_PRESSURE3D.h"
31    
32     C !INPUT/OUTPUT PARAMETERS:
33     C === Routine arguments ===
34     C myThid :: My Thread Id number
35     INTEGER myThid
36    
37     #ifdef ALLOW_NONHYDROSTATIC
38    
39     C !LOCAL VARIABLES:
40     C === Local variables ===
41     C bi,bj :: tile indices
42     C i,j,k :: Loop counters
43     C faceArea :: Temporary used to hold cell face areas.
44     C myNorm :: Work variable used in clculating normalisation factor
45     CHARACTER*(MAX_LEN_MBUF) msgBuf
46     INTEGER bi, bj
47     INTEGER i, j, k, ks
48     _RL faceArea
49     _RS myNorm
50     _RL theRecip_Dr
51     _RL aU, aL, aW, aE, aN, aS
52     _RL tmpFac, nh_Fac, igwFac
53     _RL locGamma
54     CEOP
55    
56     CcnhDebugStarts
57     c _RL phi(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
58     CcnhDebugEnds
59    
60     C-- Initialise to zero over the full range of indices
61     DO bj=myByLo(myThid),myByHi(myThid)
62     DO bi=myBxLo(myThid),myBxHi(myThid)
63     DO k=1,Nr
64     C- From common bloc CG3D_R:
65     DO j=1-OLy,sNy+OLy
66     DO i=1-OLx,sNx+OLx
67     aW3d(i,j,k,bi,bj) = 0.
68     aS3d(i,j,k,bi,bj) = 0.
69     aV3d(i,j,k,bi,bj) = 0.
70     aC3d(i,j,k,bi,bj) = 0.
71     zMC (i,j,k,bi,bj) = 0.
72     zML (i,j,k,bi,bj) = 0.
73     zMU (i,j,k,bi,bj) = 0.
74     ENDDO
75     ENDDO
76     C- From common bloc CG3D_WK_R:
77     DO j=0,sNy+1
78     DO i=0,sNx+1
79     cg3d_q(i,j,k,bi,bj) = 0.
80     cg3d_r(i,j,k,bi,bj) = 0.
81     cg3d_s(i,j,k,bi,bj) = 0.
82     ENDDO
83     ENDDO
84     C- From common bloc SFP3D_COMMON_R:
85     DO j=1-OLy,sNy+OLy
86     DO i=1-OLx,sNx+OLx
87     cg3d_b(i,j,k,bi,bj) = 0.
88     ENDDO
89     ENDDO
90     ENDDO
91     ENDDO
92     ENDDO
93    
94     nh_Fac = 0.
95     igwFac = 0.
96     IF ( nonHydrostatic
97     & .AND. nh_Am2.NE.0. ) nh_Fac = 1. _d 0 / nh_Am2
98     IF ( implicitIntGravWave ) igwFac = 1. _d 0
99    
100     IF ( use3Dsolver ) THEN
101     C-- Initialise laplace operator
102     C aW3d: Ax/dX
103     C aS3d: Ay/dY
104     C aV3d: Ar/dR
105     myNorm = 0. _d 0
106     DO bj=myByLo(myThid),myByHi(myThid)
107     DO bi=myBxLo(myThid),myBxHi(myThid)
108     DO k=1,Nr
109     DO j=1,sNy
110     DO i=1,sNx+1
111     faceArea = _dyG(i,j,bi,bj)*drF(k)
112     & *_hFacW(i,j,k,bi,bj)
113     #ifdef ALLOW_OBCS
114     & *maskInC(i,j,bi,bj)*maskInC(i-1,j,bi,bj)
115     #endif
116     aW3d(i,j,k,bi,bj) = faceArea*recip_dxC(i,j,bi,bj)
117     & *implicitNHPress*implicDiv2DFlow
118     myNorm = MAX(ABS(aW3d(i,j,k,bi,bj)),myNorm)
119     ENDDO
120     ENDDO
121     C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect
122     DO j=1,sNy+1
123     DO i=1,sNx
124     faceArea = _dxG(i,j,bi,bj)*drF(K)
125     & *_hFacS(i,j,k,bi,bj)
126     #ifdef ALLOW_OBCS
127     & *maskInC(i,j,bi,bj)*maskInC(i,j-1,bi,bj)
128     #endif
129     aS3d(i,j,k,bi,bj) = faceArea*recip_dyC(i,j,bi,bj)
130     & *implicitNHPress*implicDiv2DFlow
131     myNorm = MAX(ABS(aS3d(i,j,k,bi,bj)),myNorm)
132     ENDDO
133     ENDDO
134     ENDDO
135     DO k=1,1
136     DO j=1,sNy
137     DO i=1,sNx
138     aV3d(i,j,k,bi,bj) = 0.
139     ENDDO
140     ENDDO
141     ENDDO
142     DO k=2,Nr
143     tmpFac = nh_Fac*rVel2wUnit(k)*rVel2wUnit(k)
144     & + igwFac*dBdrRef(k)*deltaTMom*dTtracerLev(k)
145     IF (tmpFac.GT.0. ) tmpFac = 1. _d 0 / tmpFac
146     DO j=1,sNy
147     DO i=1,sNx
148     faceArea = _rA(i,j,bi,bj)*maskC(i,j, k ,bi,bj)
149     & *maskC(i,j,k-1,bi,bj)
150     & *deepFac2F(k)
151     #ifdef ALLOW_OBCS
152     & *maskInC(i,j,bi,bj)
153     #endif
154     theRecip_Dr = recip_drC(k)
155     c theRecip_Dr =
156     caja & drF(k )*_hFacC(i,j,k ,bi,bj)*0.5
157     caja & +drF(k-1)*_hFacC(i,j,k-1,bi,bj)*0.5
158     c IF ( theRecip_Dr .NE. 0. )
159     c & theRecip_Dr = 1. _d 0/theRecip_Dr
160     aV3d(i,j,k,bi,bj) = faceArea*theRecip_Dr*tmpFac
161     & *implicitNHPress*implicDiv2DFlow
162     myNorm = MAX(ABS(aV3d(i,j,k,bi,bj)),myNorm)
163     ENDDO
164     ENDDO
165     ENDDO
166     ENDDO
167     ENDDO
168     _GLOBAL_MAX_RS( myNorm, myThid )
169     IF ( myNorm .NE. 0. _d 0 ) THEN
170     myNorm = 1. _d 0/myNorm
171     ELSE
172     myNorm = 1. _d 0
173     ENDIF
174    
175     _BEGIN_MASTER( myThid )
176     C-- set global parameter in common block:
177     cg3dNorm = myNorm
178     WRITE(msgBuf,'(2A,1PE23.16)') 'INI_CG3D: ',
179     & 'CG3D normalisation factor = ', cg3dNorm
180     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
181     & SQUEEZE_RIGHT, myThid )
182     WRITE(msgBuf,*) ' '
183     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
184     & SQUEEZE_RIGHT, myThid )
185     _END_MASTER( myThid )
186    
187     DO bj=myByLo(myThid),myByHi(myThid)
188     DO bi=myBxLo(myThid),myBxHi(myThid)
189     C- Set solver main diagonal term
190     DO k=1,Nr
191     DO j=1,sNy
192     DO i=1,sNx
193     aW = aW3d( i, j, k, bi,bj)
194     aE = aW3d(i+1,j, k, bi,bj)
195     aN = aS3d( i,j+1,k, bi,bj)
196     aS = aS3d( i, j, k, bi,bj)
197     aU = aV3d( i, j, k, bi,bj)
198     IF ( k .NE. Nr ) THEN
199     aL = aV3d(i, j,k+1,bi,bj)
200     ELSE
201     aL = 0.
202     ENDIF
203     aC3d(i,j,k,bi,bj) = -aW-aE-aN-aS-aU-aL
204     ENDDO
205     ENDDO
206     ENDDO
207     C- Add free-surface source term
208     IF ( selectNHfreeSurf.GE.1 ) THEN
209     DO j=1,sNy
210     DO i=1,sNx
211     locGamma = drC(1)*recip_Bo(i,j,bi,bj)
212     & /( deltaTMom*deltaTFreeSurf
213     & *implicitNHPress*implicDiv2DFlow )
214     ks = 1
215     c ks = kSurfC(i,j,bi,bj)
216     c IF ( ks.LE.Nr ) THEN
217     aC3d(i,j,ks,bi,bj) = aC3d(i,j,ks,bi,bj)
218     & - freeSurfFac*recip_Bo(i,j,bi,bj)
219     & *rA(i,j,bi,bj)*deepFac2F(ks)/deltaTMom/deltaTFreeSurf
220     & / (1. _d 0 + locGamma )
221     c ENDIF
222     ENDDO
223     ENDDO
224     ELSE
225     DO j=1,sNy
226     DO i=1,sNx
227     ks = kSurfC(i,j,bi,bj)
228     IF ( ks.LE.Nr ) THEN
229     aC3d(i,j,ks,bi,bj) = aC3d(i,j,ks,bi,bj)
230     & - freeSurfFac*recip_Bo(i,j,bi,bj)
231     & *rA(i,j,bi,bj)*deepFac2F(ks)/deltaTMom/deltaTFreeSurf
232     ENDIF
233     ENDDO
234     ENDDO
235     ENDIF
236     C- Matrix solver normalisation
237     DO k=1,Nr
238     DO j=1,sNy
239     DO i=1,sNx
240     aW3d(i,j,k,bi,bj) = aW3d(i,j,k,bi,bj)*myNorm
241     aS3d(i,j,k,bi,bj) = aS3d(i,j,k,bi,bj)*myNorm
242     aV3d(i,j,k,bi,bj) = aV3d(i,j,k,bi,bj)*myNorm
243     aC3d(i,j,k,bi,bj) = aC3d(i,j,k,bi,bj)*myNorm
244     ENDDO
245     ENDDO
246     ENDDO
247     ENDDO
248     ENDDO
249    
250     C-- Update overlap regions
251     CALL EXCH_UV_XYZ_RS(aW3d,aS3d,.FALSE.,myThid)
252     _EXCH_XYZ_RS(aV3d, myThid)
253     _EXCH_XYZ_RS(aC3d, myThid)
254     CcnhDebugStarts
255     C CALL PLOT_FIELD_XYZRS( aW3d, 'AW3D INI_CG3D.1' , Nr, 1, myThid )
256     C CALL PLOT_FIELD_XYZRS( aS3d, 'AS3D INI_CG3D.1' , Nr, 1, myThid )
257     CcnhDebugEnds
258    
259     C-- Initialise preconditioner
260     C For now PC is just the identity. Change to
261     C be LU factorization of d2/dz2 later. Note
262     C check for consistency with S/R CG3D before
263     C assuming zML is lower and zMU is upper!
264     DO bj=myByLo(myThid),myByHi(myThid)
265     DO bi=myBxLo(myThid),myBxHi(myThid)
266     DO k=1,Nr
267     DO j=1,sNy
268     DO i=1,sNx
269     IF ( aC3d(i,j,k,bi,bj) .NE. 0. ) THEN
270     zMC(i,j,k,bi,bj) = aC3d(i,j,k,bi,bj)
271     zML(i,j,k,bi,bj) = aV3d(i,j,k,bi,bj)
272     IF ( k.NE.Nr ) THEN
273     zMU(i,j,k,bi,bj)= aV3d(i,j,k+1,bi,bj)
274     ELSE
275     zMU(i,j,k,bi,bj)= 0.
276     ENDIF
277     CcnhDebugStarts
278     C zMC(i,j,k,bi,bj) = 1.
279     C zMU(i,j,k,bi,bj) = 0.
280     C zML(i,j,k,bi,bj) = 0.
281     CcnhDebugEnds
282     ELSE
283     zMC(i,j,k,bi,bj) = 1. _d 0
284     zMU(i,j,k,bi,bj) = 0.
285     zML(i,j,k,bi,bj) = 0.
286     ENDIF
287     ENDDO
288     ENDDO
289     ENDDO
290     k = 1
291     DO j=1,sNy
292     DO i=1,sNx
293     zMC(i,j,k,bi,bj) = 1. _d 0 / zMC(i,j,k,bi,bj)
294     zMU(i,j,k,bi,bj) = zMU(i,j,k,bi,bj)*zMC(i,j,k,bi,bj)
295     ENDDO
296     ENDDO
297     DO k=2,Nr
298     DO j=1,sNy
299     DO i=1,sNx
300     zMC(i,j,k,bi,bj) = 1. _d 0 /
301     & (zMC(i,j,k,bi,bj)-zML(i,j,k,bi,bj)*zMU(i,j,k-1,bi,bj))
302     zMU(i,j,k,bi,bj) = zMU(i,j,k,bi,bj)*zMC(i,j,k,bi,bj)
303     ENDDO
304     ENDDO
305     ENDDO
306     DO k=1,Nr
307     DO j=1,sNy
308     DO i=1,sNx
309     IF ( aC3d(i,j,k,bi,bj) .EQ. 0. ) THEN
310     zMC(i,j,k,bi,bj) = 1.
311     zML(i,j,k,bi,bj) = 0.
312     zMU(i,j,k,bi,bj) = 0.
313     CcnhDebugStarts
314     C ELSE
315     C zMC(i,j,k,bi,bj) = 1.
316     C zML(i,j,k,bi,bj) = 0.
317     C zMU(i,j,k,bi,bj) = 0.
318     CcnhDEbugEnds
319     ENDIF
320     ENDDO
321     ENDDO
322     ENDDO
323     ENDDO
324     ENDDO
325     C-- Update overlap regions
326     _EXCH_XYZ_RS(zMC, myThid)
327     _EXCH_XYZ_RS(zML, myThid)
328     _EXCH_XYZ_RS(zMU, myThid)
329    
330     IF ( debugLevel .GE. debLevC ) THEN
331     CALL WRITE_FLD_XYZ_RS( 'zMC',' ',zMC, 0, myThid )
332     CALL WRITE_FLD_XYZ_RS( 'zML',' ',zML, 0, myThid )
333     CALL WRITE_FLD_XYZ_RS( 'zMU',' ',zMU, 0, myThid )
334     ENDIF
335     CcnhDebugStarts
336     c DO k=1,Nr
337     c DO j=1-OLy,sNy+OLy
338     c DO i=1-OLx,sNx+OLx
339     c phi(i,j,1,1) = zMc(i,j,k,1,1)
340     c ENDDO
341     c ENDDO
342     C CALL PLOT_FIELD_XYRS( phi, 'zMC INI_CG3D.1' , 1, myThid )
343     c ENDDO
344     C CALL PLOT_FIELD_XYRS( zMU, 'zMU INI_CG3D.1' , Nr, 1, myThid )
345     C CALL PLOT_FIELD_XYRS( zML, 'zML INI_CG3D.1' , Nr, 1, myThid )
346     CcnhDebugEnds
347    
348     C-- end if (use3Dsolver)
349     ENDIF
350    
351     #ifdef ALLOW_PETSC
352     IF (use_cg3d_petsc) THEN
353     CALL CG3D_INITIALIZE_PETSC
354     CALL CG3D_PETSC_NUMERATE (mythid)
355     ENDIF
356     #endif
357    
358     #endif /* ALLOW_NONHYDROSTATIC */
359    
360     RETURN
361     END

  ViewVC Help
Powered by ViewVC 1.1.22