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

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

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


Revision 1.4 - (hide annotations) (download)
Thu Oct 9 04:19:18 2003 UTC (20 years, 7 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint52d_pre, checkpoint51o_pre, checkpoint51l_post, checkpoint52, checkpoint51t_post, checkpoint51n_post, checkpoint51s_post, checkpoint51n_pre, checkpoint52b_pre, checkpoint51l_pre, checkpoint51q_post, checkpoint52b_post, checkpoint52c_post, checkpoint51r_post, checkpoint51i_post, checkpoint52d_post, checkpoint52a_pre, branch-netcdf, checkpoint51o_post, checkpoint52a_post, ecco_c52_e35, checkpoint51m_post, checkpoint51p_post, checkpoint51u_post
Branch point for: branch-nonh, tg2-branch, netcdf-sm0, checkpoint51n_branch
Changes since 1.3: +3 -1 lines
 o first check-in for the "branch-genmake2" merge
 o verification suite as run on shelley (gcc 3.2.2):

Wed Oct  8 23:42:29 EDT 2003
                T           S           U           V
G D M    c        m  s        m  s        m  s        m  s
E p a R  g  m  m  e  .  m  m  e  .  m  m  e  .  m  m  e  .
N n k u  2  i  a  a  d  i  a  a  d  i  a  a  d  i  a  a  d
2 d e n  d  n  x  n  .  n  x  n  .  n  x  n  .  n  x  n  .

OPTFILE=NONE

Y Y Y Y 13 16 16 16  0 16 16 16 16 16 16 16 16 13 12  0  0 pass  adjustment.128x64x1
Y Y Y Y 16 16 16 16  0 16 16 16 16 16 16  0  0 16 16  0  0 pass  adjustment.cs-32x32x1
Y Y Y Y 16 16 16 16  0 16 16 16 16 16 16 22  0 16 16 22  0 pass  adjust_nlfs.cs-32x32x1
Y Y Y Y -- 13 13 16 16 13 13 13 13 16 16 16 16 16 16 16 16 N/O   advect_cs
Y Y Y Y -- 22 16 16 16 16 16 16 13 16 16 16 16 16 16 16 16 N/O   advect_xy
Y Y Y Y -- 13 16 13 16 16 16 16 16 16 16 22 16 16 16 16 16 N/O   advect_xz
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 pass  aim.5l_cs
Y Y Y Y 14 16 16 16 16 16 16 16 16 13 16 16 16 16 16 13 16 pass  aim.5l_Equatorial_Channel
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 13 16 16 13 13 16 pass  aim.5l_LatLon
Y Y Y Y 13 16 16 16 16 16 16 16 16 16 13 12 13 13 16 13 16 pass  exp0
Y Y Y Y 14 16 16 16 16 16 16 16 22 16 16 16 13 16 16 22 16 pass  exp1
Y Y Y Y 13 13 16 13 16 16 16 16 16 13 13 16 16 13 13 13 13 pass  exp2
Y Y Y Y 16 16 16 16 16 16 16 16 22 16 16 16 16 16 16 16 16 pass  exp4
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 22 16 16 16 22 16 pass  exp5
Y Y Y Y 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 16 pass  front_relax
Y Y Y Y 14 16 16 13 13 16 16 13 13 16 13 13 16 12 13 13 16 pass  global_ocean.90x40x15
Y Y Y Y 10 16 16 13 13 16 13 16 16 13 13 13 13 16 16 13 16 FAIL  global_ocean.cs32x15
Y Y Y Y  6 11 12 13 13 12 13 16 13  9  9  9  9 10  9  9 11 FAIL  global_ocean_pressure
Y Y Y Y 14 16 16 13 16 16 16 13 13 13 13 13 16 12 16 13 16 pass  global_with_exf
Y Y Y Y 14 16 16 16 16 16 16 16 16 11 13 22 13 16 16  9 16 pass  hs94.128x64x5
Y Y Y Y 13 16 16 16 16 16 16 16 16 11 16 16 16 13 16 22 13 pass  hs94.1x64x5
Y Y Y Y 14 16 16 16 16 16 16 16 16 13 16 13 13 16 16 22 13 pass  hs94.cs-32x32x5
Y Y Y Y 10 10 16 13 13 16 16 16 22 16 13 13 13 13 13 22 13 FAIL  ideal_2D_oce
Y Y Y Y  8 16 16 16 16 16 16 16 16 13 13  8 16 16 16 16 16 FAIL  internal_wave
Y Y Y Y 14 16 16 16 16 16 16 16 16 13 13 22 13 13 13 22 16 pass  inverted_barometer
Y Y Y Y 12 16 16 16 16 16 16 16 16 16 13 12 13 13 13 13 13 FAIL  lab_sea
Y Y Y Y 11 16 16 16 16 16 16 16 13 13 13 12 13 16 13 12 13 FAIL  natl_box
Y Y Y Y 16 16 16 16 16 16 16 16 22 16 16 16 16 16 16 16 16 pass  plume_on_slope
Y Y Y Y 13 16 16 16 16 13 16 16 16 16 16 16 16 13 16 16 16 pass  solid-body.cs-32x32x1

1 edhill 1.4 C $Header: /u/u3/gcmpack/MITgcm/model/src/update_cg2d.F,v 1.3.6.1 2003/10/02 18:10:45 edhill Exp $
2 jmc 1.1 C $Name: $
3    
4 edhill 1.4 #include "PACKAGES_CONFIG.h"
5 jmc 1.1 #include "CPP_OPTIONS.h"
6 edhill 1.4
7 cnh 1.2 CBOP
8     C !ROUTINE: UPDATE_CG2D
9     C !INTERFACE:
10     SUBROUTINE UPDATE_CG2D( myTime, myIter, myThid )
11     C !DESCRIPTION: \bv
12     C *==========================================================*
13     C | SUBROUTINE UPDATE_CG2D
14     C | o Update 2d conjugate gradient solver operators
15     C | account for Free-Surf effect on total column thickness
16     C *==========================================================*
17     C | This routine is based on INI_CG2D, and simplified. It is
18     C | only active when the non-linear free surface mode of
19     C | equations is active.
20     C *==========================================================*
21     C \ev
22 jmc 1.1
23 cnh 1.2 C !USES:
24 jmc 1.1 IMPLICIT NONE
25     C === Global variables ===
26     #include "SIZE.h"
27     #include "EEPARAMS.h"
28     #include "PARAMS.h"
29     #include "GRID.h"
30     #include "SURFACE.h"
31     #include "CG2D.h"
32     #ifdef ALLOW_OBCS
33     #include "OBCS.h"
34     #endif
35    
36 cnh 1.2 C !INPUT/OUTPUT PARAMETERS:
37 jmc 1.1 C === Routine arguments ===
38     C myTime - Current time in simulation
39     C myIter - Current iteration number in simulation
40     C myThid - Thread number for this instance of the routine.
41     _RL myTime
42     INTEGER myIter
43     INTEGER myThid
44 cnh 1.2
45     C !LOCAL VARIABLES:
46 jmc 1.1 #ifdef NONLIN_FRSURF
47     C-- Note : compared to "INI_CG2D", no needs to compute again
48     C the solver norn=malisation factor of the solver tolerance
49     C === Local variables ===
50     C bi,bj - Loop counters
51     C I,J,K
52     C faceArea - Temporary used to hold cell face areas.
53     INTEGER bi, bj
54     INTEGER I, J, K
55     _RL faceArea
56 jmc 1.3 _RL pW_tmp, pS_tmp
57 cnh 1.2 CEOP
58 jmc 1.1
59     C-- Initialise laplace operator
60     C aW2d: integral in Z Ax/dX
61     C aS2d: integral in Z Ay/dY
62     DO bj=myByLo(myThid),myByHi(myThid)
63     DO bi=myBxLo(myThid),myBxHi(myThid)
64 jmc 1.3 DO J=1,sNy+1
65     DO I=1,sNx+1
66 jmc 1.1 aW2d(I,J,bi,bj) = 0. _d 0
67     aS2d(I,J,bi,bj) = 0. _d 0
68     ENDDO
69     ENDDO
70     DO K=1,Nr
71 jmc 1.3 DO J=1,sNy+1
72     DO I=1,sNx+1
73 jmc 1.1 faceArea = _dyG(I,J,bi,bj)*drF(K)
74     & *_hFacW(I,J,K,bi,bj)
75     aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) +
76     & faceArea*recip_dxC(I,J,bi,bj)
77     faceArea = _dxG(I,J,bi,bj)*drF(K)
78     & *_hFacS(I,J,K,bi,bj)
79     aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) +
80     & faceArea*recip_dyC(I,J,bi,bj)
81     ENDDO
82     ENDDO
83     ENDDO
84     #ifdef ALLOW_OBCS
85     IF (useOBCS) THEN
86     DO I=1,sNx
87     IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0.
88     IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0.
89     IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj)+1,bi,bj)=0.
90     IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj),bi,bj)=0.
91     ENDDO
92     DO J=1,sNy
93     IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj),J,bi,bj)=0.
94     IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj)+1,J,bi,bj)=0.
95     IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj)+1,J,bi,bj)=0.
96     IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj),J,bi,bj)=0.
97     ENDDO
98     ENDIF
99     #endif
100 jmc 1.3 DO J=1,sNy+1
101     DO I=1,sNx+1
102 jmc 1.1 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*cg2dNorm
103     & *implicSurfPress*implicDiv2DFlow
104     aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*cg2dNorm
105     & *implicSurfPress*implicDiv2DFlow
106     ENDDO
107     ENDDO
108 jmc 1.3 C-- Start to compute preconditioner matrix (use cg2d_q as temporary array)
109     DO J=1,sNy
110     DO I=1,sNx
111     cg2d_q(I,J,bi,bj) = -(
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     ENDDO
118     ENDDO
119 jmc 1.1 ENDDO
120     ENDDO
121    
122     C-- Update overlap regions
123 jmc 1.3 _EXCH_XY_R8(cg2d_q, myThid)
124 jmc 1.1
125     C-- Initialise preconditioner
126     DO bj=myByLo(myThid),myByHi(myThid)
127     DO bi=myBxLo(myThid),myBxHi(myThid)
128 jmc 1.3 DO J=1,sNy+1
129     DO I=1,sNx+1
130     IF ( cg2d_q(I,J,bi,bj) .EQ. 0. ) THEN
131 jmc 1.1 pC(I,J,bi,bj) = 1. _d 0
132     ELSE
133 jmc 1.3 pC(I,J,bi,bj) = 1. _d 0 / cg2d_q(I,J,bi,bj)
134 jmc 1.1 ENDIF
135 jmc 1.3 pW_tmp = cg2d_q(I,J,bi,bj)+cg2d_q(I-1,J,bi,bj)
136     IF ( pW_tmp .EQ. 0. ) THEN
137 jmc 1.1 pW(I,J,bi,bj) = 0.
138     ELSE
139     pW(I,J,bi,bj) =
140 jmc 1.3 & -aW2d(I,J,bi,bj)/((cg2dpcOffDFac *pW_tmp)**2 )
141 jmc 1.1 ENDIF
142 jmc 1.3 pS_tmp = cg2d_q(I,J,bi,bj)+cg2d_q(I,J-1,bi,bj)
143     IF ( pS_tmp .EQ. 0. ) THEN
144 jmc 1.1 pS(I,J,bi,bj) = 0.
145     ELSE
146     pS(I,J,bi,bj) =
147 jmc 1.3 & -aS2d(I,J,bi,bj)/((cg2dpcOffDFac *pS_tmp)**2 )
148 jmc 1.1 ENDIF
149     ENDDO
150     ENDDO
151     ENDDO
152     ENDDO
153    
154     #endif /* NONLIN_FRSURF */
155    
156     RETURN
157     END

  ViewVC Help
Powered by ViewVC 1.1.22