/[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.11 - (show annotations) (download)
Mon Mar 22 15:54:03 1999 UTC (25 years, 2 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint28, checkpoint29, checkpoint20, checkpoint21, checkpoint22, checkpoint23, checkpoint24, checkpoint25, checkpoint27, checkpoint26, checkpoint30
Changes since 1.10: +21 -2 lines
Modifications for non-hydrostatic ability + updates for open-boundaries.

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

  ViewVC Help
Powered by ViewVC 1.1.22