/[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.1 - (show annotations) (download)
Mon Aug 27 18:50:10 2001 UTC (22 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint40pre9, checkpoint40
New routines relative to NonLin-FreeSurf

1 C $Header: $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 SUBROUTINE UPDATE_CG2D( myTime, myIter, myThid )
7 C /==========================================================\
8 C | SUBROUTINE UPDATE_CG2D |
9 C | o Update 2d conjugate gradient solver operators |
10 C | account for Free-Surf effect on total column thickness |
11 C |==========================================================|
12 C | This routine is based on INI_CG2D, and simplified |
13 C \==========================================================/
14 IMPLICIT NONE
15
16 C === Global variables ===
17 #include "SIZE.h"
18 #include "EEPARAMS.h"
19 #include "PARAMS.h"
20 #include "GRID.h"
21 #include "SURFACE.h"
22 #include "CG2D.h"
23 #ifdef ALLOW_OBCS
24 #include "OBCS.h"
25 #endif
26
27 C === Routine arguments ===
28 C myTime - Current time in simulation
29 C myIter - Current iteration number in simulation
30 C myThid - Thread number for this instance of the routine.
31 _RL myTime
32 INTEGER myIter
33 INTEGER myThid
34
35 #ifdef NONLIN_FRSURF
36 C-- Note : compared to "INI_CG2D", no needs to compute again
37 C the solver norn=malisation factor of the solver tolerance
38
39 C === Local variables ===
40 C bi,bj - Loop counters
41 C I,J,K
42 C faceArea - Temporary used to hold cell face areas.
43 INTEGER bi, bj
44 INTEGER I, J, K
45 _RL faceArea
46 _RL aC, aCw, aCs
47
48 C-- Initialise laplace operator
49 C aW2d: integral in Z Ax/dX
50 C aS2d: integral in Z Ay/dY
51 DO bj=myByLo(myThid),myByHi(myThid)
52 DO bi=myBxLo(myThid),myBxHi(myThid)
53 DO J=1,sNy
54 DO I=1,sNx
55 aW2d(I,J,bi,bj) = 0. _d 0
56 aS2d(I,J,bi,bj) = 0. _d 0
57 ENDDO
58 ENDDO
59 DO K=1,Nr
60 DO J=1,sNy
61 DO I=1,sNx
62 faceArea = _dyG(I,J,bi,bj)*drF(K)
63 & *_hFacW(I,J,K,bi,bj)
64 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) +
65 & faceArea*recip_dxC(I,J,bi,bj)
66 faceArea = _dxG(I,J,bi,bj)*drF(K)
67 & *_hFacS(I,J,K,bi,bj)
68 aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) +
69 & faceArea*recip_dyC(I,J,bi,bj)
70 ENDDO
71 ENDDO
72 ENDDO
73 #ifdef ALLOW_OBCS
74 IF (useOBCS) THEN
75 DO I=1,sNx
76 IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0.
77 IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0.
78 IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj)+1,bi,bj)=0.
79 IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj),bi,bj)=0.
80 ENDDO
81 DO J=1,sNy
82 IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj),J,bi,bj)=0.
83 IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj)+1,J,bi,bj)=0.
84 IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj)+1,J,bi,bj)=0.
85 IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj),J,bi,bj)=0.
86 ENDDO
87 ENDIF
88 #endif
89 DO J=1,sNy
90 DO I=1,sNx
91 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*cg2dNorm
92 & *implicSurfPress*implicDiv2DFlow
93 aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*cg2dNorm
94 & *implicSurfPress*implicDiv2DFlow
95 ENDDO
96 ENDDO
97 ENDDO
98 ENDDO
99
100 C-- Update overlap regions
101 c _EXCH_XY_R4(aW2d, myThid)
102 c _EXCH_XY_R4(aS2d, myThid)
103 CALL EXCH_UV_XY_RS(aW2d,aS2d,.FALSE.,myThid)
104
105 C-- Initialise preconditioner
106 DO bj=myByLo(myThid),myByHi(myThid)
107 DO bi=myBxLo(myThid),myBxHi(myThid)
108 DO J=1,sNy
109 DO I=1,sNx
110 pC(I,J,bi,bj) = 1. _d 0
111 aC = -(
112 & aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj)
113 & +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj)
114 & +freeSurfFac*cg2dNorm*recip_Bo(I,J,bi,bj)*
115 & rA(I,J,bi,bj)/deltaTMom/deltaTMom
116 & )
117 aCs = -(
118 & aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj)
119 & +aS2d(I,J-1,bi,bj) + aS2d(I ,J ,bi,bj)
120 & +freeSurfFac*cg2dNorm*recip_Bo(I,J-1,bi,bj)*
121 & rA(I,J-1,bi,bj)/deltaTMom/deltaTMom
122 & )
123 aCw = -(
124 & aW2d(I-1,J,bi,bj) + aW2d(I ,J ,bi,bj)
125 & +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj)
126 & +freeSurfFac*cg2dNorm*recip_Bo(I-1,J,bi,bj)*
127 & rA(I-1,J,bi,bj)/deltaTMom/deltaTMom
128 & )
129 IF ( aC .EQ. 0. ) THEN
130 pC(I,J,bi,bj) = 1. _d 0
131 ELSE
132 pC(I,J,bi,bj) = 1. _d 0 / aC
133 ENDIF
134 IF ( aC + aCw .EQ. 0. ) THEN
135 pW(I,J,bi,bj) = 0.
136 ELSE
137 pW(I,J,bi,bj) =
138 & -aW2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 )
139 ENDIF
140 IF ( aC + aCs .EQ. 0. ) THEN
141 pS(I,J,bi,bj) = 0.
142 ELSE
143 pS(I,J,bi,bj) =
144 & -aS2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 )
145 ENDIF
146 ENDDO
147 ENDDO
148 ENDDO
149 ENDDO
150 C-- Update overlap regions
151 _EXCH_XY_R4(pC, myThid)
152 c _EXCH_XY_R4(pW, myThid)
153 c _EXCH_XY_R4(pS, myThid)
154 CALL EXCH_UV_XY_RS(pW,pS,.FALSE.,myThid)
155
156 #endif /* NONLIN_FRSURF */
157
158 RETURN
159 END

  ViewVC Help
Powered by ViewVC 1.1.22