forked from darwinproject/darwin3
-
Notifications
You must be signed in to change notification settings - Fork 1
/
solve_for_pressure.F
461 lines (431 loc) · 15.8 KB
/
solve_for_pressure.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
#include "PACKAGES_CONFIG.h"
#include "CPP_OPTIONS.h"
CBOP
C !ROUTINE: SOLVE_FOR_PRESSURE
C !INTERFACE:
SUBROUTINE SOLVE_FOR_PRESSURE( myTime, myIter, myThid )
C !DESCRIPTION: \bv
C *==========================================================*
C | SUBROUTINE SOLVE_FOR_PRESSURE
C | o Controls inversion of two and/or three-dimensional
C | elliptic problems for the pressure field.
C *==========================================================*
C \ev
C !USES:
IMPLICIT NONE
C == Global variables
#include "SIZE.h"
#include "EEPARAMS.h"
#include "PARAMS.h"
#include "GRID.h"
#include "SURFACE.h"
#include "FFIELDS.h"
#include "DYNVARS.h"
#include "SOLVE_FOR_PRESSURE.h"
#ifdef ALLOW_NONHYDROSTATIC
#include "SOLVE_FOR_PRESSURE3D.h"
#include "NH_VARS.h"
#endif
#ifdef ALLOW_CD_CODE
#include "CD_CODE_VARS.h"
#endif
#ifdef ALLOW_AUTODIFF
# include "AUTODIFF_PARAMS.h"
#endif
C === Functions ====
LOGICAL DIFFERENT_MULTIPLE
EXTERNAL DIFFERENT_MULTIPLE
#ifdef ALLOW_DIAGNOSTICS
LOGICAL DIAGNOSTICS_IS_ON
EXTERNAL DIAGNOSTICS_IS_ON
#endif /* ALLOW_DIAGNOSTICS */
C !INPUT/OUTPUT PARAMETERS:
C == Routine arguments ==
C myTime :: Current time in simulation
C myIter :: Current iteration number in simulation
C myThid :: Thread number for this instance of SOLVE_FOR_PRESSURE
_RL myTime
INTEGER myIter
INTEGER myThid
C !LOCAL VARIABLES:
C == Local variables ==
INTEGER i,j,k,bi,bj
INTEGER ks
INTEGER numIters, nIterMin
_RL firstResidual, minResidualSq, lastResidual
_RL tmpFac
_RL sumEmP, tileEmP(nSx,nSy)
LOGICAL putPmEinXvector
INTEGER ioUnit
CHARACTER*(MAX_LEN_MBUF) msgBuf
#ifdef ALLOW_NONHYDROSTATIC
LOGICAL zeroPsNH, zeroMeanPnh, oldFreeSurfTerm
#else
_RL cg3d_b(1)
#endif
#ifdef ALLOW_DIAGNOSTICS
CHARACTER*8 diagName
_RL tmpVar (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
#endif /* ALLOW_DIAGNOSTICS */
CEOP
#ifdef ALLOW_NONHYDROSTATIC
zeroPsNH = .FALSE.
c zeroPsNH = use3Dsolver .AND. exactConserv
c & .AND. select_rStar.EQ.0
zeroMeanPnh = .FALSE.
c zeroMeanPnh = use3Dsolver .AND. select_rStar.NE.0
c oldFreeSurfTerm = use3Dsolver .AND. select_rStar.EQ.0
c & .AND. .NOT.zeroPsNH
oldFreeSurfTerm = use3Dsolver .AND. .NOT.exactConserv
#else
cg3d_b(1) = 0.
#endif
C deepAtmosphere & useRealFreshWaterFlux: only valid if deepFac2F(ksurf)=1
C anelastic (always Z-coordinate):
C 1) assume that rhoFacF(1)=1 (and ksurf == 1);
C (this reduces the number of lines of code to modify)
C 2) (a) 2-D continuity eq. compute div. of mass transport (<- add rhoFac)
C (b) gradient of surf.Press in momentum eq. (<- add 1/rhoFac)
C => 2 factors cancel in elliptic eq. for Phi_s ,
C but 1rst factor(a) remains in RHS cg2d_b.
C-- Initialise the Vector solution with etaN + deltaT*Global_mean_PmE
C instead of simply etaN ; This can speed-up the solver convergence in
C the case where |Global_mean_PmE| is large.
putPmEinXvector = .FALSE.
c putPmEinXvector = useRealFreshWaterFlux.AND.fluidIsWater
IF ( myIter.EQ.1+nIter0 .AND. debugLevel .GE. debLevA ) THEN
_BEGIN_MASTER( myThid )
ioUnit = standardMessageUnit
WRITE(msgBuf,'(2A,L5)') 'SOLVE_FOR_PRESSURE:',
& ' putPmEinXvector =', putPmEinXvector
CALL PRINT_MESSAGE( msgBuf, ioUnit, SQUEEZE_RIGHT, myThid )
#ifdef ALLOW_NONHYDROSTATIC
WRITE(msgBuf,'(A,2(A,L5))') 'SOLVE_FOR_PRESSURE:',
& ' zeroPsNH=', zeroPsNH, ' , zeroMeanPnh=', zeroMeanPnh
CALL PRINT_MESSAGE( msgBuf, ioUnit, SQUEEZE_RIGHT, myThid )
WRITE(msgBuf,'(2A,L5)') 'SOLVE_FOR_PRESSURE:',
& ' oldFreeSurfTerm =', oldFreeSurfTerm
CALL PRINT_MESSAGE( msgBuf, ioUnit, SQUEEZE_RIGHT, myThid )
#endif
_END_MASTER( myThid )
ENDIF
C-- Save previous solution & Initialise Vector solution and source term :
sumEmP = 0.
DO bj=myByLo(myThid),myByHi(myThid)
DO bi=myBxLo(myThid),myBxHi(myThid)
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
#ifdef ALLOW_CD_CODE
etaNm1(i,j,bi,bj) = etaN(i,j,bi,bj)
#endif
cg2d_x(i,j,bi,bj) = Bo_surf(i,j,bi,bj)*etaN(i,j,bi,bj)
cg2d_b(i,j,bi,bj) = 0.
ENDDO
ENDDO
IF (useRealFreshWaterFlux.AND.fluidIsWater) THEN
tmpFac = freeSurfFac*mass2rUnit*implicDiv2DFlow
DO j=1,sNy
DO i=1,sNx
cg2d_b(i,j,bi,bj) =
& tmpFac*_rA(i,j,bi,bj)*EmPmR(i,j,bi,bj)/deltaTMom
& *maskInC(i,j,bi,bj)
ENDDO
ENDDO
ENDIF
IF ( putPmEinXvector ) THEN
tileEmP(bi,bj) = 0.
DO j=1,sNy
DO i=1,sNx
tileEmP(bi,bj) = tileEmP(bi,bj)
& + rA(i,j,bi,bj)*EmPmR(i,j,bi,bj)
& *maskInC(i,j,bi,bj)
ENDDO
ENDDO
ENDIF
ENDDO
ENDDO
IF ( putPmEinXvector ) THEN
CALL GLOBAL_SUM_TILE_RL( tileEmP, sumEmP, myThid )
ENDIF
DO bj=myByLo(myThid),myByHi(myThid)
DO bi=myBxLo(myThid),myBxHi(myThid)
IF ( putPmEinXvector ) THEN
tmpFac = 0.
IF (globalArea.GT.0.) tmpFac =
& freeSurfFac*deltaTFreeSurf*mass2rUnit*sumEmP/globalArea
DO j=1,sNy
DO i=1,sNx
cg2d_x(i,j,bi,bj) = cg2d_x(i,j,bi,bj)
& - tmpFac*Bo_surf(i,j,bi,bj)
ENDDO
ENDDO
ENDIF
C- RHS: similar to the divergence of the vertically integrated mass transport:
C del_i { Sum_k [ rhoFac.(dr.hFac).(dy.deepFac).(u*) ] } / deltaT
DO k=Nr,1,-1
CALL CALC_DIV_GHAT(
I bi,bj,k,
U cg2d_b, cg3d_b,
I myThid )
ENDDO
ENDDO
ENDDO
DO bj=myByLo(myThid),myByHi(myThid)
DO bi=myBxLo(myThid),myBxHi(myThid)
#ifdef ALLOW_NONHYDROSTATIC
IF ( oldFreeSurfTerm ) THEN
C-- Add source term arising from w=d/dt (p_s + p_nh)
DO j=1,sNy
DO i=1,sNx
ks = kSurfC(i,j,bi,bj)
IF ( ks.LE.Nr ) THEN
cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj)
& -freeSurfFac*_rA(i,j,bi,bj)*deepFac2F(ks)
& /deltaTMom/deltaTFreeSurf
& *( etaN(i,j,bi,bj)
& +phi_nh(i,j,ks,bi,bj)*recip_Bo(i,j,bi,bj) )
cg3d_b(i,j,ks,bi,bj) = cg3d_b(i,j,ks,bi,bj)
& -freeSurfFac*_rA(i,j,bi,bj)*deepFac2F(ks)
& /deltaTMom/deltaTFreeSurf
& *( etaN(i,j,bi,bj)
& +phi_nh(i,j,ks,bi,bj)*recip_Bo(i,j,bi,bj) )
ENDIF
ENDDO
ENDDO
ELSEIF ( exactConserv ) THEN
#else
C-- Add source term arising from w=d/dt (p_s)
IF ( exactConserv ) THEN
#endif /* ALLOW_NONHYDROSTATIC */
DO j=1,sNy
DO i=1,sNx
ks = kSurfC(i,j,bi,bj)
cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj)
& -freeSurfFac*_rA(i,j,bi,bj)*deepFac2F(ks)
& /deltaTMom/deltaTFreeSurf
& * etaH(i,j,bi,bj)
ENDDO
ENDDO
ELSE
DO j=1,sNy
DO i=1,sNx
ks = kSurfC(i,j,bi,bj)
cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj)
& -freeSurfFac*_rA(i,j,bi,bj)*deepFac2F(ks)
& /deltaTMom/deltaTFreeSurf
& * etaN(i,j,bi,bj)
ENDDO
ENDDO
ENDIF
#ifdef ALLOW_OBCS
C- Note: solver matrix is trivial outside OB region (main diagonal only)
C => no real need to reset RHS (=cg2d_b) & cg2d_x, except that:
C a) normalisation is fct of Max(RHS), which can be large ouside OB region
C (would be different if we were solving for increment of eta/g
C instead of directly for eta/g).
C => need to reset RHS to ensure that interior solution does not depend
C on ouside OB region.
C b) provide directly the trivial solution cg2d_x == 0 for outside OB region
C (=> no residual => no effect on solver convergence and interior solution)
IF (useOBCS) THEN
DO j=1,sNy
DO i=1,sNx
cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj)*maskInC(i,j,bi,bj)
cg2d_x(i,j,bi,bj) = cg2d_x(i,j,bi,bj)*maskInC(i,j,bi,bj)
ENDDO
ENDDO
ENDIF
#endif /* ALLOW_OBCS */
C- end bi,bj loops
ENDDO
ENDDO
#ifdef ALLOW_DEBUG
IF ( debugLevel .GE. debLevD ) THEN
CALL DEBUG_STATS_RL(1,cg2d_b,'cg2d_b (SOLVE_FOR_PRESSURE)',
& myThid)
ENDIF
#endif
IF ( DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock) ) THEN
CALL WRITE_FLD_XY_RL( 'cg2d_b', 'I10', cg2d_b, myIter, myThid )
ENDIF
C-- Find the surface pressure using a two-dimensional conjugate
C gradient solver. See CG2D.h for the interface to this routine.
C In rare cases of a poor solver convergence, better to select the
C solver minimum-residual solution (instead of the last-iter solution)
C by setting cg2dUseMinResSol=1 (<-> nIterMin=0 in input)
numIters = cg2dMaxIters
nIterMin = cg2dUseMinResSol - 1
c CALL TIMER_START('CG2D [SOLVE_FOR_PRESSURE]',myThid)
#ifdef DISCONNECTED_TILES
C-- Call the disconnected-tile (no EXCH) version of cg2d
CALL CG2D_EX0(
U cg2d_b, cg2d_x,
O firstResidual, minResidualSq, lastResidual,
U numIters, nIterMin,
I myThid )
#else /* not DISCONNECTED_TILES = default */
IF ( useSRCGSolver ) THEN
C-- Call the single reduce CG solver
#ifdef ALLOW_SRCG
CALL CG2D_SR(
U cg2d_b, cg2d_x,
O firstResidual, minResidualSq, lastResidual,
U numIters, nIterMin,
I myThid )
#endif /* ALLOW_SRCG */
#ifdef ALLOW_CG2D_NSA
ELSEIF ( useNSACGSolver ) THEN
C-- Call the not-self-adjoint version of cg2d
CALL CG2D_NSA(
U cg2d_b, cg2d_x,
O firstResidual, minResidualSq, lastResidual,
U numIters, nIterMin,
I myThid )
#endif
ELSE
C-- Call the default CG solver
CALL CG2D(
U cg2d_b, cg2d_x,
O firstResidual, minResidualSq, lastResidual,
U numIters, nIterMin,
I myThid )
ENDIF
#endif /* DISCONNECTED_TILES */
_EXCH_XY_RL( cg2d_x, myThid )
#ifdef ALLOW_AUTODIFF
IF ( .NOT. useNSACGSolver .AND. cg2dFullAdjoint )
& CALL CG2D_STORE( cg2d_x, .TRUE., myThid )
#endif
c CALL TIMER_STOP ('CG2D [SOLVE_FOR_PRESSURE]',myThid)
#ifdef ALLOW_DEBUG
IF ( debugLevel .GE. debLevD ) THEN
CALL DEBUG_STATS_RL(1,cg2d_x,'cg2d_x (SOLVE_FOR_PRESSURE)',
& myThid)
ENDIF
#endif
C- dump CG2D output at monitorFreq (to reduce size of STD-OUTPUT files) :
IF ( DIFFERENT_MULTIPLE(monitorFreq,myTime,deltaTClock)
& ) THEN
IF ( debugLevel .GE. debLevA ) THEN
_BEGIN_MASTER( myThid )
WRITE(msgBuf,'(A20,1PE23.14)') 'cg2d_init_res =',firstResidual
CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
WRITE(msgBuf,'(A27,2I8)')
& 'cg2d_iters(min,last) =', nIterMin, numIters
CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
IF ( minResidualSq.GE.0. ) THEN
minResidualSq = SQRT(minResidualSq)
WRITE(msgBuf,'(A20,1PE23.14)') 'cg2d_min_res =',minResidualSq
CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
ENDIF
WRITE(msgBuf,'(A20,1PE23.14)') 'cg2d_last_res =',lastResidual
CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
_END_MASTER( myThid )
ENDIF
ENDIF
#ifdef ALLOW_DIAGNOSTICS
C-- Fill diagnostics
IF ( useDiagnostics .AND. implicSurfPress.NE.oneRL ) THEN
diagName = 'PHI_SURF'
IF ( DIAGNOSTICS_IS_ON(diagName,myThid) ) THEN
DO bj=myByLo(myThid),myByHi(myThid)
DO bi=myBxLo(myThid),myBxHi(myThid)
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
tmpVar(i,j) = implicSurfPress * cg2d_x(i,j,bi,bj)
& + (oneRL - implicSurfPress)* Bo_surf(i,j,bi,bj)
& * etaN(i,j,bi,bj)
ENDDO
ENDDO
CALL DIAGNOSTICS_FILL( tmpVar,diagName,1,1,2,bi,bj,myThid )
ENDDO
ENDDO
ENDIF
ELSEIF ( useDiagnostics ) THEN
CALL DIAGNOSTICS_FILL( cg2d_x,'PHI_SURF', 0,1, 0,1,1, myThid )
ENDIF
#endif /* ALLOW_DIAGNOSTICS */
C-- Transfert the 2D-solution to "etaN" :
DO bj=myByLo(myThid),myByHi(myThid)
DO bi=myBxLo(myThid),myBxHi(myThid)
DO j=1-OLy,sNy+OLy
DO i=1-OLx,sNx+OLx
etaN(i,j,bi,bj) = recip_Bo(i,j,bi,bj)*cg2d_x(i,j,bi,bj)
ENDDO
ENDDO
ENDDO
ENDDO
#ifdef ALLOW_NONHYDROSTATIC
IF ( use3Dsolver ) THEN
IF ( DIFFERENT_MULTIPLE(diagFreq, myTime, deltaTClock) ) THEN
CALL WRITE_FLD_XY_RL( 'cg2d_x','I10', cg2d_x, myIter, myThid )
ENDIF
C-- Solve for a three-dimensional pressure term (NH or IGW or both ).
C see CG3D.h for the interface to this routine.
C-- Finish updating cg3d_b: 1) Add EmPmR contribution to top level cg3d_b:
C 2) Update or Add free-surface contribution
C 3) increment in horiz velocity due to new cg2d_x
C 4) add vertical velocity contribution.
CALL PRE_CG3D(
I oldFreeSurfTerm,
I cg2d_x,
U cg3d_b,
I myTime, myIter, myThid )
#ifdef ALLOW_DEBUG
IF ( debugLevel .GE. debLevD ) THEN
CALL DEBUG_STATS_RL(Nr,cg3d_b,'cg3d_b (SOLVE_FOR_PRESSURE)',
& myThid)
ENDIF
#endif
IF ( DIFFERENT_MULTIPLE( diagFreq, myTime, deltaTClock) ) THEN
CALL WRITE_FLD_XYZ_RL('cg3d_b','I10', cg3d_b, myIter,myThid )
ENDIF
firstResidual=0.
lastResidual=0.
numIters=cg3dMaxIters
CALL TIMER_START('CG3D [SOLVE_FOR_PRESSURE]',myThid)
#ifdef DISCONNECTED_TILES
CALL CG3D_EX0(
U cg3d_b, phi_nh,
O firstResidual, lastResidual,
U numIters,
I myIter, myThid )
#else /* not DISCONNECTED_TILES = default */
CALL CG3D(
U cg3d_b, phi_nh,
O firstResidual, lastResidual,
U numIters,
I myIter, myThid )
#endif /* DISCONNECTED_TILES */
_EXCH_XYZ_RL( phi_nh, myThid )
CALL TIMER_STOP ('CG3D [SOLVE_FOR_PRESSURE]',myThid)
IF ( DIFFERENT_MULTIPLE(monitorFreq,myTime,deltaTClock)
& ) THEN
IF ( debugLevel .GE. debLevA ) THEN
_BEGIN_MASTER( myThid )
WRITE(msgBuf,'(A20,1PE23.14)') 'cg3d_init_res =',firstResidual
CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
WRITE(msgBuf,'(A27,I16)') 'cg3d_iters (last) = ',numIters
CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
WRITE(msgBuf,'(A20,1PE23.14)') 'cg3d_last_res =',lastResidual
CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
_END_MASTER( myThid )
ENDIF
ENDIF
C-- Separate the Hydrostatic Surface Pressure adjusment (=> put it in dPhiNH)
C from the Non-hydrostatic pressure (since cg3d_x contains both contribution)
IF ( nonHydrostatic .AND. exactConserv ) THEN
IF ( DIFFERENT_MULTIPLE( diagFreq, myTime, deltaTClock) ) THEN
CALL WRITE_FLD_XYZ_RL('cg3d_x','I10', phi_nh, myIter,myThid )
ENDIF
CALL POST_CG3D(
I zeroPsNH, zeroMeanPnh,
I myTime, myIter, myThid )
ENDIF
ENDIF
#endif /* ALLOW_NONHYDROSTATIC */
#ifdef ALLOW_SHOWFLOPS
CALL SHOWFLOPS_INSOLVE( myThid)
#endif
RETURN
END