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

Contents of /MITgcm/model/src/correction_step.F

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


Revision 1.15 - (show annotations) (download)
Tue Feb 20 15:06:21 2001 UTC (23 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint36
Changes since 1.14: +4 -4 lines
implement a Crank-Nickelson barotropic time-stepping

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/correction_step.F,v 1.14 2001/02/04 14:38:46 cnh Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 C /==========================================================\
7 C | S/R CORRECTION_STEP |
8 C | o Corrects the horizontal flow fields with the surface |
9 C | slope. |
10 C \==========================================================/
11 SUBROUTINE CORRECTION_STEP( bi, bj, iMin, iMax, jMin, jMax,
12 I K, etaSurfX, etaSurfY,
13 I myCurrentTime, myThid )
14 IMPLICIT NONE
15
16 C == Global variables ==
17 #include "SIZE.h"
18 #include "DYNVARS.h"
19 #include "EEPARAMS.h"
20 #include "PARAMS.h"
21 #include "GRID.h"
22 #include "CG2D.h"
23 #ifdef ALLOW_NONHYDROSTATIC
24 #include "CG3D.h"
25 #endif
26 C == Routine Arguments ==
27 C etaSurfX, etaSurfY - Surface slope
28 C bi,bj,iMin,iMax,jMin,jMax, K - Loop counters
29 C myThid - Instance number for
30 C this call to S/R CORRECTION_STEP
31 C myCurrentTime - Current simulation time for this instance.
32 _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
33 _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
34 INTEGER bi,bj,iMin,iMax,jMin,jMax
35 INTEGER K
36 INTEGER myThid
37 _RL myCurrentTime
38
39 C == Local variables ==
40 INTEGER i,j
41 _RL hxFac,hyFac
42 _RL hx3dFac,hy3dFac
43
44 C On/off scaling paramters
45 hxFac = pfFacMom
46 hyFac = pfFacMom
47 IF ( nonHydrostatic ) THEN
48 hx3dFac = pfFacMom
49 hy3dFac = pfFacMom
50 ELSE
51 hx3dFac = 0.
52 hy3dFac = 0.
53 ENDIF
54
55 C Step forward zonal velocity
56 DO j=jMin,jMax
57 DO i=iMin,iMax
58 uVel(i,j,k,bi,bj)=( gUNm1(i,j,k,bi,bj)
59 & -deltaTmom*hxFac*gBaro*implicSurfPress*etaSurfX(i,j)
60 #ifdef ALLOW_NONHYDROSTATIC
61 & -deltaTmom*hx3dFac*gravity*_recip_dxC(i,j,bi,bj)*
62 & (cg3d_x(i,j,k,bi,bj)-cg3d_x(i-1,j,k,bi,bj))
63 #endif
64 & )*_maskW(i,j,k,bi,bj)
65 gUNm1(i,j,k,bi,bj)=gU(i,j,k,bi,bj)
66 ENDDO
67 ENDDO
68
69 C Step forward meridional velocity
70 DO j=jMin,jMax
71 DO i=iMin,iMax
72 vVel(i,j,k,bi,bj)=( gVNm1(i,j,k,bi,bj)
73 & -deltaTmom*hyFac*gBaro*implicSurfPress*etaSurfY(i,j)
74 #ifdef ALLOW_NONHYDROSTATIC
75 & -deltaTmom*hy3dFac*gravity*_recip_dyC(i,j,bi,bj)*
76 & (cg3d_x(i,j,k,bi,bj)-cg3d_x(i,j-1,k,bi,bj))
77 #endif
78 & )*_maskS(i,j,k,bi,bj)
79 gVNm1(i,j,k,bi,bj)=gV(i,j,k,bi,bj)
80 ENDDO
81 ENDDO
82
83 RETURN
84 END

  ViewVC Help
Powered by ViewVC 1.1.22