/[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.17 - (show annotations) (download)
Thu Mar 8 20:35:39 2001 UTC (23 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint38, pre38tag1, c37_adj, pre38-close, checkpoint39, checkpoint37
Branch point for: pre38
Changes since 1.16: +9 -9 lines
Use directly gradient of PhiSurf as argument (instead of grad.eta)

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/correction_step.F,v 1.16 2001/03/06 17:10:29 jmc 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, phiSurfX, phiSurfY,
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 #ifdef ALLOW_NONHYDROSTATIC
23 #include "CG3D.h"
24 #endif
25 C == Routine Arguments ==
26 C phiSurfX, phiSurfY - Surface Potential gradient
27 C bi,bj,iMin,iMax,jMin,jMax, K - Loop counters
28 C myThid - Instance number for
29 C this call to S/R CORRECTION_STEP
30 C myCurrentTime - Current simulation time for this instance.
31 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
32 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
33 INTEGER bi,bj,iMin,iMax,jMin,jMax
34 INTEGER K
35 INTEGER myThid
36 _RL myCurrentTime
37
38 C == Local variables ==
39 INTEGER i,j
40 _RL hxFac,hyFac
41 _RL hx3dFac,hy3dFac
42
43 C On/off scaling paramters
44 hxFac = pfFacMom
45 hyFac = pfFacMom
46 IF ( nonHydrostatic ) THEN
47 hx3dFac = pfFacMom
48 hy3dFac = pfFacMom
49 ELSE
50 hx3dFac = 0.
51 hy3dFac = 0.
52 ENDIF
53
54 C Step forward zonal velocity
55 DO j=jMin,jMax
56 DO i=iMin,iMax
57 uVel(i,j,k,bi,bj)=( gUNm1(i,j,k,bi,bj)
58 & -deltaTmom*hxFac*implicSurfPress*phiSurfX(i,j)
59 #ifdef ALLOW_NONHYDROSTATIC
60 & -deltaTmom*hx3dFac*_recip_dxC(i,j,bi,bj)*
61 & (cg3d_x(i,j,k,bi,bj)-cg3d_x(i-1,j,k,bi,bj))
62 #endif
63 & )*_maskW(i,j,k,bi,bj)
64 gUNm1(i,j,k,bi,bj)=gU(i,j,k,bi,bj)
65 ENDDO
66 ENDDO
67
68 C Step forward meridional velocity
69 DO j=jMin,jMax
70 DO i=iMin,iMax
71 vVel(i,j,k,bi,bj)=( gVNm1(i,j,k,bi,bj)
72 & -deltaTmom*hyFac*implicSurfPress*phiSurfY(i,j)
73 #ifdef ALLOW_NONHYDROSTATIC
74 & -deltaTmom*hy3dFac*_recip_dyC(i,j,bi,bj)*
75 & (cg3d_x(i,j,k,bi,bj)-cg3d_x(i,j-1,k,bi,bj))
76 #endif
77 & )*_maskS(i,j,k,bi,bj)
78 gVNm1(i,j,k,bi,bj)=gV(i,j,k,bi,bj)
79 ENDDO
80 ENDDO
81
82 RETURN
83 END

  ViewVC Help
Powered by ViewVC 1.1.22