/[MITgcm]/MITgcm/model/src/ini_cg2d.F
ViewVC logotype

Annotation of /MITgcm/model/src/ini_cg2d.F

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


Revision 1.1 - (hide annotations) (download)
Wed Apr 22 19:15:30 1998 UTC (26 years, 1 month ago) by cnh
Branch: MAIN
Branch point for: cnh
Initial revision

1 cnh 1.1 C $Id$
2    
3     #include "CPP_EEOPTIONS.h"
4    
5     CStartOfInterface
6     SUBROUTINE INI_CG2D( myThid )
7     C /==========================================================\
8     C | SUBROUTINE INI_CG2D |
9     C | o Initialise 2d conjugate gradient solver operators. |
10     C |==========================================================|
11     C | These arrays are purely a function of the basin geom. |
12     C | We set then here once and them use then repeatedly. |
13     C \==========================================================/
14    
15     C === Global variables ===
16     #include "SIZE.h"
17     #include "EEPARAMS.h"
18     #include "PARAMS.h"
19     #include "GRID.h"
20     #include "CG2D.h"
21    
22     C === Routine arguments ===
23     C myThid - Thread no. that called this routine.
24     INTEGER myThid
25     CEndOfInterface
26    
27     C === Local variables ===
28     C xG, yG - Global coordinate location.
29     C zG
30     C iG, jG - Global coordinate index
31     C bi,bj - Loop counters
32     C faceArea - Temporary used to hold cell face areas.
33     C I,J,K
34     C myNorm - Work variable used in clculating normalisation factor
35     CHARACTER*(MAX_LEN_MBUF) msgBuf
36     INTEGER bi, bj
37     INTEGER I, J, K
38     real faceArea
39     _RL myNorm
40    
41     C-- Initialise laplace operator
42     C aW2d: integral in Z Ax/dX
43     C aS2d: integral in Z Ay/dY
44     myNorm = 0. _d 0
45     DO bj=myByLo(myThid),myByHi(myThid)
46     DO bi=myBxLo(myThid),myBxHi(myThid)
47     DO J=1,sNy
48     DO I=1,sNx
49     aW2d(I,J,bi,bj) = 0. _d 0
50     aS2d(I,J,bi,bj) = 0. _d 0
51     ENDDO
52     ENDDO
53     DO K=1,Nz
54     DO J=1,sNy
55     DO I=1,sNx
56     faceArea = dyG(I,J,bi,bj)*dzF(K)*HFacW(I,J,K,bi,bj)
57     aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) +
58     & gravity*faceArea*rDxC(I,J,bi,bj)
59     faceArea = dxG(I,J,bi,bj)*dzF(K)*HFacS(I,J,K,bi,bj)
60     aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) +
61     & gravity*faceArea*rDyC(I,J,bi,bj)
62     ENDDO
63     ENDDO
64     ENDDO
65     DO J=1,sNy
66     DO I=1,sNx
67     myNorm = MAX(ABS(aW2d(I,J,bi,bj)),myNorm)
68     myNorm = MAX(ABS(aS2d(I,J,bi,bj)),myNorm)
69     ENDDO
70     ENDDO
71     ENDDO
72     ENDDO
73     cg2dNbuf(1,myThid) = myNorm
74     _GLOBAL_MAX_R8( cg2dNbuf, myNorm, myThid )
75     IF ( cg2dNbuf(1,1) .NE. 0. _d 0 ) THEN
76     myNorm = 1. _d 0/cg2dNbuf(1,1)
77     ELSE
78     myNorm = 1. _d 0
79     ENDIF
80     _BEGIN_MASTER( myThid )
81     cg2dNorm = myNorm
82     CcnhDebugStarts
83     WRITE(msgBuf,*) ' CG2D normalisation factor = ', cg2dNorm
84     CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
85     CcnhDebugEnds
86     _END_MASTER( myThid )
87     DO bj=myByLo(myThid),myByHi(myThid)
88     DO bi=myBxLo(myThid),myBxHi(myThid)
89     DO J=1,sNy
90     DO I=1,sNx
91     aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*myNorm
92     aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*myNorm
93     ENDDO
94     ENDDO
95     ENDDO
96     ENDDO
97    
98     C-- Update overlap regions
99     CcnhDebugStarts
100     C CALL PLOT_FIELD_XYR8( aW2d, 'AW2D INI_CG2D.1' , 1, myThid )
101     C CALL PLOT_FIELD_XYR8( aS2d, 'AS2D INI_CG2D.1' , 1, myThid )
102     CcnhDebugEnds
103     _EXCH_XY_R4(aW2d, myThid)
104     _EXCH_XY_R4(aS2d, myThid)
105     CcnhDebugStarts
106     C CALL PLOT_FIELD_XYR8( aW2d, 'AW2D INI_CG2D.2' , 1, myThid )
107     C CALL PLOT_FIELD_XYR8( aS2d, 'AS2D INI_CG2D.2' , 1, myThid )
108     CcnhDebugEnds
109    
110     C-- Initialise preconditioner
111     DO bj=myByLo(myThid),myByHi(myThid)
112     DO bi=myBxLo(myThid),myBxHi(myThid)
113     DO J=1,sNy
114     DO I=1,sNx
115     pC(I,J,bi,bj) = 1. _d 0
116     IF (
117     & aW2d(I,J,bi,bj) + aW2d(I+1,J,bi,bj)
118     & +aS2d(I,J,bi,bj) + aS2D(I,J+1,bi,bj)
119     & .EQ. 0.
120     & ) pC(I,J,bi,bj) = 0. _d 0
121     pW(I,J,bi,bj) = 0.
122     pS(I,J,bi,bj) = 0.
123     ENDDO
124     ENDDO
125     ENDDO
126     ENDDO
127     C-- Update overlap regions
128     _EXCH_XY_R4(pC, myThid)
129     _EXCH_XY_R4(pW, myThid)
130     _EXCH_XY_R4(pS, myThid)
131    
132     C-- Set default values for initial guess
133     DO bj=myByLo(myThid),myByHi(myThid)
134     DO bi=myBxLo(myThid),myBxHi(myThid)
135     DO J=1,sNy
136     DO I=1,sNx
137     cg2d_x(I,J,bi,bj) = 0. _d 0
138     ENDDO
139     ENDDO
140     ENDDO
141     ENDDO
142     C-- Update overlap regions
143     _EXCH_XY_R8(cg2d_x, myThid)
144    
145     RETURN
146     END

  ViewVC Help
Powered by ViewVC 1.1.22