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

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

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

revision 1.26 by adcroft, Mon May 24 15:42:23 1999 UTC revision 1.38 by jmc, Sun Feb 10 20:04:11 2002 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    C $Name$
3    
4  #include "CPP_OPTIONS.h"  #include "CPP_OPTIONS.h"
5    
6  CStartOfInterface  CBOP
7    C     !ROUTINE: INI_CG2D
8    C     !INTERFACE:
9        SUBROUTINE INI_CG2D( myThid )        SUBROUTINE INI_CG2D( myThid )
 C     /==========================================================\  
 C     | SUBROUTINE INI_CG2D                                      |  
 C     | o Initialise 2d conjugate gradient solver operators.     |  
 C     |==========================================================|  
 C     | These arrays are purely a function of the basin geom.    |  
 C     | We set then here once and them use then repeatedly.      |  
 C     \==========================================================/  
       IMPLICIT NONE  
10    
11    C     !DESCRIPTION: \bv
12    C     *==========================================================*
13    C     | SUBROUTINE INI_CG2D                                      
14    C     | o Initialise 2d 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 ===  C     === Global variables ===
24  #include "SIZE.h"  #include "SIZE.h"
25  #include "EEPARAMS.h"  #include "EEPARAMS.h"
26  #include "PARAMS.h"  #include "PARAMS.h"
27  #include "GRID.h"  #include "GRID.h"
28  #include "DYNVARS.h"  c#include "DYNVARS.h"
29    #include "SURFACE.h"
30  #include "CG2D.h"  #include "CG2D.h"
31  #ifdef ALLOW_OBCS  #ifdef ALLOW_OBCS
32  #include "OBCS.h"  #include "OBCS.h"
33  #endif  #endif
34    
35    C     !INPUT/OUTPUT PARAMETERS:
36  C     === Routine arguments ===  C     === Routine arguments ===
37  C     myThid - Thread no. that called this routine.  C     myThid - Thread no. that called this routine.
38        INTEGER myThid        INTEGER myThid
 CEndOfInterface  
39    
40    C     !LOCAL VARIABLES:
41  C     === Local variables ===  C     === Local variables ===
42  C     xG, yG - Global coordinate location.  C     xG, yG - Global coordinate location.
43  C     zG  C     zG
# Line 36  C     iG, jG - Global coordinate index Line 45  C     iG, jG - Global coordinate index
45  C     bi,bj  - Loop counters  C     bi,bj  - Loop counters
46  C     faceArea - Temporary used to hold cell face areas.  C     faceArea - Temporary used to hold cell face areas.
47  C     I,J,K  C     I,J,K
48  C     myNorm - Work variable used in clculating normalisation factor  C     myNorm - Work variable used in calculating normalisation factor
49    C     sumArea - Work variable used to compute the total Domain Area
50        CHARACTER*(MAX_LEN_MBUF) msgBuf        CHARACTER*(MAX_LEN_MBUF) msgBuf
51        INTEGER bi, bj        INTEGER bi, bj
52        INTEGER I,  J, K        INTEGER I,  J, K
53        _RL     faceArea        _RL     faceArea
54        _RS     myNorm        _RS     myNorm
55          _RL     sumArea
56        _RL     aC, aCw, aCs        _RL     aC, aCw, aCs
57    CEOP
58    
59    C--   Initialize arrays in common blocs (CG2D.h) ; not really necessary
60    C     but safer when EXCH do not fill all the overlap regions.
61          DO bj=myByLo(myThid),myByHi(myThid)
62           DO bi=myBxLo(myThid),myBxHi(myThid)
63            DO J=1-OLy,sNy+OLy
64             DO I=1-OLx,sNx+OLx
65              aW2d(I,J,bi,bj) = 0. _d 0
66              aS2d(I,J,bi,bj) = 0. _d 0
67              pW(I,J,bi,bj) = 0. _d 0
68              pS(I,J,bi,bj) = 0. _d 0
69              pC(I,J,bi,bj) = 0. _d 0
70              cg2d_q(I,J,bi,bj) = 0. _d 0
71             ENDDO
72            ENDDO
73            DO J=1-1,sNy+1
74             DO I=1-1,sNx+1
75              cg2d_r(I,J,bi,bj) = 0. _d 0
76              cg2d_s(I,J,bi,bj) = 0. _d 0
77             ENDDO
78            ENDDO
79           ENDDO
80          ENDDO
81    
82  C--   Initialise laplace operator  C--   Initialise laplace operator
83  C     aW2d: integral in Z Ax/dX  C     aW2d: integral in Z Ax/dX
# Line 62  C     aS2d: integral in Z Ay/dY Line 97  C     aS2d: integral in Z Ay/dY
97             faceArea = _dyG(I,J,bi,bj)*drF(K)             faceArea = _dyG(I,J,bi,bj)*drF(K)
98       &               *_hFacW(I,J,K,bi,bj)       &               *_hFacW(I,J,K,bi,bj)
99             aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) +             aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) +
100       &      gBaro*faceArea*recip_dxC(I,J,bi,bj)       &      implicSurfPress*implicDiv2DFlow
101         &           *faceArea*recip_dxC(I,J,bi,bj)
102             faceArea = _dxG(I,J,bi,bj)*drF(K)             faceArea = _dxG(I,J,bi,bj)*drF(K)
103       &               *_hFacS(I,J,K,bi,bj)       &               *_hFacS(I,J,K,bi,bj)
104             aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) +             aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) +
105       &      gBaro*faceArea*recip_dyC(I,J,bi,bj)       &      implicSurfPress*implicDiv2DFlow
106         &           *faceArea*recip_dyC(I,J,bi,bj)
107            ENDDO            ENDDO
108           ENDDO           ENDDO
109          ENDDO          ENDDO
110  #ifdef ALLOW_OBCS  #ifdef ALLOW_OBCS
111          IF (openBoundaries) THEN          IF (useOBCS) THEN
112           DO I=1,sNx           DO I=1,sNx
113            IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0.            IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0.
114            IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0.            IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0.
# Line 126  CcnhDebugStarts Line 163  CcnhDebugStarts
163  C     CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid )  C     CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid )
164  C     CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid )  C     CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid )
165  CcnhDebugEnds  CcnhDebugEnds
166        _EXCH_XY_R4(aW2d, myThid)  c     _EXCH_XY_R4(aW2d, myThid)
167        _EXCH_XY_R4(aS2d, myThid)  c     _EXCH_XY_R4(aS2d, myThid)
168          CALL EXCH_UV_XY_RS(aW2d,aS2d,.FALSE.,myThid)
169  CcnhDebugStarts  CcnhDebugStarts
170  C     CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid )  C     CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid )
171  C     CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid )  C     CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid )
172  CcnhDebugEnds  CcnhDebugEnds
173    
174    C-- Define the solver tolerance in the appropriate Unit :
175          cg2dNormaliseRHS = cg2dTargetResWunit.LE.0
176          IF (cg2dNormaliseRHS) THEN
177    C-  when using a normalisation of RHS, tolerance has no unit => no conversion
178            cg2dTolerance = cg2dTargetResidual
179          ELSE
180    C-  compute the total Area of the domain :
181            sumArea = 0.
182            DO bj=myByLo(myThid),myByHi(myThid)
183             DO bi=myBxLo(myThid),myBxHi(myThid)
184              DO j=1,sNy
185               DO i=1,sNx
186                 IF (Ro_surf(i,j,bi,bj).GT.R_low(i,j,bi,bj)) THEN
187                   sumArea = sumArea + rA(i,j,bi,bj)
188                 ENDIF
189               ENDDO
190              ENDDO
191             ENDDO
192            ENDDO
193    c       WRITE(*,*) ' mythid, sumArea = ', mythid, sumArea
194            _GLOBAL_SUM_R8( sumArea, myThid )
195    C-  convert Target-Residual (in W unit) to cg2d-solver residual unit [m^2/s^2]
196            cg2dTolerance = cg2dNorm * cg2dTargetResWunit
197         &                           * sumArea / deltaTMom
198            WRITE(*,'(2A,1P2E22.14)') ' ini_cg2d: ',
199         &          'sumArea,cg2dTolerance =', sumArea,cg2dTolerance
200          ENDIF
201        
202  C--   Initialise preconditioner  C--   Initialise preconditioner
203  C     Note. 20th May 1998  C     Note. 20th May 1998
204  C           I made a weird discovery! In the model paper we argue  C           I made a weird discovery! In the model paper we argue
# Line 155  C           defaults to 0.51 but can be Line 221  C           defaults to 0.51 but can be
221            aC = -(            aC = -(
222       &     aW2d(I,J,bi,bj) + aW2d(I+1,J  ,bi,bj)       &     aW2d(I,J,bi,bj) + aW2d(I+1,J  ,bi,bj)
223       &    +aS2d(I,J,bi,bj) + aS2d(I  ,J+1,bi,bj)       &    +aS2d(I,J,bi,bj) + aS2d(I  ,J+1,bi,bj)
224       &    +freeSurfFac*myNorm* horiVertRatio*       &    +freeSurfFac*myNorm*recip_Bo(I,J,bi,bj)*
225       &     rA(I,J,bi,bj)/deltaTMom/deltaTMom       &     rA(I,J,bi,bj)/deltaTMom/deltaTMom
226       &    )       &    )
227            aCs = -(            aCs = -(
228       &     aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj)       &     aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj)
229       &    +aS2d(I,J-1,bi,bj) + aS2d(I  ,J  ,bi,bj)       &    +aS2d(I,J-1,bi,bj) + aS2d(I  ,J  ,bi,bj)
230       &    +freeSurfFac*myNorm* horiVertRatio*       &    +freeSurfFac*myNorm*recip_Bo(I,J-1,bi,bj)*
231       &     rA(I,J-1,bi,bj)/deltaTMom/deltaTMom       &     rA(I,J-1,bi,bj)/deltaTMom/deltaTMom
232       &    )       &    )
233            aCw = -(            aCw = -(
234       &     aW2d(I-1,J,bi,bj) + aW2d(I  ,J  ,bi,bj)       &     aW2d(I-1,J,bi,bj) + aW2d(I  ,J  ,bi,bj)
235       &    +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj)       &    +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj)
236       &    +freeSurfFac*myNorm* horiVertRatio*       &    +freeSurfFac*myNorm*recip_Bo(I-1,J,bi,bj)*
237       &     rA(I-1,J,bi,bj)/deltaTMom/deltaTMom       &     rA(I-1,J,bi,bj)/deltaTMom/deltaTMom
238       &    )       &    )
239            IF ( aC .EQ. 0. ) THEN            IF ( aC .EQ. 0. ) THEN
240              pC(I,J,bi,bj) = 0. _d 0              pC(I,J,bi,bj) = 1. _d 0
241            ELSE            ELSE
242             pC(I,J,bi,bj) =  1. _d 0 / aC             pC(I,J,bi,bj) =  1. _d 0 / aC
243            ENDIF            ENDIF
# Line 196  C         pS(I,J,bi,bj) = 0. Line 262  C         pS(I,J,bi,bj) = 0.
262        ENDDO        ENDDO
263  C--   Update overlap regions  C--   Update overlap regions
264        _EXCH_XY_R4(pC, myThid)        _EXCH_XY_R4(pC, myThid)
265        _EXCH_XY_R4(pW, myThid)  c     _EXCH_XY_R4(pW, myThid)
266        _EXCH_XY_R4(pS, myThid)  c     _EXCH_XY_R4(pS, myThid)
267          CALL EXCH_UV_XY_RS(pW,pS,.FALSE.,myThid)
268  CcnhDebugStarts  CcnhDebugStarts
269  C     CALL PLOT_FIELD_XYRS( pC, 'pC   INI_CG2D.2' , 1, myThid )  C     CALL PLOT_FIELD_XYRS( pC, 'pC   INI_CG2D.2' , 1, myThid )
270  C     CALL PLOT_FIELD_XYRS( pW, 'pW   INI_CG2D.2' , 1, myThid )  C     CALL PLOT_FIELD_XYRS( pW, 'pW   INI_CG2D.2' , 1, myThid )
271  C     CALL PLOT_FIELD_XYRS( pS, 'pS   INI_CG2D.2' , 1, myThid )  C     CALL PLOT_FIELD_XYRS( pS, 'pS   INI_CG2D.2' , 1, myThid )
272  CcnhDebugEnds  CcnhDebugEnds
273    
 C--   Set default values for initial guess and RHS  
       IF ( startTime .EQ. 0 ) THEN  
        DO bj=myByLo(myThid),myByHi(myThid)  
         DO bi=myBxLo(myThid),myBxHi(myThid)  
          DO J=1,sNy  
           DO I=1,sNx  
            cg2d_x(I,J,bi,bj) = 0. _d 0  
            cg2d_b(I,J,bi,bj) = 0. _d 0  
 #ifdef INCLUDE_CD_CODE  
            cg2d_xNM1(I,J,bi,bj) = 0. _d 0  
 #endif  
           ENDDO  
          ENDDO  
         ENDDO  
        ENDDO  
 C--    Update overlap regions  
        _EXCH_XY_R8(cg2d_x, myThid)  
        _EXCH_XY_R8(cg2d_b, myThid)  
 #ifdef INCLUDE_CD_CODE  
        _EXCH_XY_R8(cg2d_xNM1, myThid)  
 #endif  
       ENDIF  
   
274        RETURN        RETURN
275        END        END

Legend:
Removed from v.1.26  
changed lines
  Added in v.1.38

  ViewVC Help
Powered by ViewVC 1.1.22