1 |
jmc |
1.17 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/correction_step.F,v 1.16 2001/03/06 17:10:29 jmc Exp $ |
2 |
jmc |
1.15 |
C $Name: $ |
3 |
adcroft |
1.1 |
|
4 |
cnh |
1.10 |
#include "CPP_OPTIONS.h" |
5 |
adcroft |
1.1 |
|
6 |
|
|
C /==========================================================\ |
7 |
|
|
C | S/R CORRECTION_STEP | |
8 |
|
|
C | o Corrects the horizontal flow fields with the surface | |
9 |
cnh |
1.7 |
C | slope. | |
10 |
adcroft |
1.1 |
C \==========================================================/ |
11 |
|
|
SUBROUTINE CORRECTION_STEP( bi, bj, iMin, iMax, jMin, jMax, |
12 |
jmc |
1.17 |
I K, phiSurfX, phiSurfY, |
13 |
adcroft |
1.6 |
I myCurrentTime, myThid ) |
14 |
cnh |
1.7 |
IMPLICIT NONE |
15 |
heimbach |
1.12 |
|
16 |
|
|
C == Global variables == |
17 |
adcroft |
1.1 |
#include "SIZE.h" |
18 |
|
|
#include "DYNVARS.h" |
19 |
cnh |
1.3 |
#include "EEPARAMS.h" |
20 |
adcroft |
1.1 |
#include "PARAMS.h" |
21 |
|
|
#include "GRID.h" |
22 |
adcroft |
1.11 |
#ifdef ALLOW_NONHYDROSTATIC |
23 |
|
|
#include "CG3D.h" |
24 |
|
|
#endif |
25 |
adcroft |
1.1 |
C == Routine Arguments == |
26 |
jmc |
1.17 |
C phiSurfX, phiSurfY - Surface Potential gradient |
27 |
cnh |
1.7 |
C bi,bj,iMin,iMax,jMin,jMax, K - Loop counters |
28 |
cnh |
1.9 |
C myThid - Instance number for |
29 |
|
|
C this call to S/R CORRECTION_STEP |
30 |
|
|
C myCurrentTime - Current simulation time for this instance. |
31 |
jmc |
1.17 |
_RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
32 |
|
|
_RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
33 |
adcroft |
1.1 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
34 |
|
|
INTEGER K |
35 |
|
|
INTEGER myThid |
36 |
adcroft |
1.6 |
_RL myCurrentTime |
37 |
cnh |
1.7 |
|
38 |
adcroft |
1.1 |
C == Local variables == |
39 |
|
|
INTEGER i,j |
40 |
adcroft |
1.11 |
_RL hxFac,hyFac |
41 |
|
|
_RL hx3dFac,hy3dFac |
42 |
adcroft |
1.1 |
|
43 |
|
|
C On/off scaling paramters |
44 |
|
|
hxFac = pfFacMom |
45 |
|
|
hyFac = pfFacMom |
46 |
adcroft |
1.11 |
IF ( nonHydrostatic ) THEN |
47 |
|
|
hx3dFac = pfFacMom |
48 |
|
|
hy3dFac = pfFacMom |
49 |
|
|
ELSE |
50 |
|
|
hx3dFac = 0. |
51 |
|
|
hy3dFac = 0. |
52 |
|
|
ENDIF |
53 |
adcroft |
1.1 |
|
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 |
jmc |
1.17 |
& -deltaTmom*hxFac*implicSurfPress*phiSurfX(i,j) |
59 |
adcroft |
1.11 |
#ifdef ALLOW_NONHYDROSTATIC |
60 |
jmc |
1.17 |
& -deltaTmom*hx3dFac*_recip_dxC(i,j,bi,bj)* |
61 |
adcroft |
1.11 |
& (cg3d_x(i,j,k,bi,bj)-cg3d_x(i-1,j,k,bi,bj)) |
62 |
|
|
#endif |
63 |
adcroft |
1.1 |
& )*_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 |
jmc |
1.17 |
& -deltaTmom*hyFac*implicSurfPress*phiSurfY(i,j) |
73 |
adcroft |
1.11 |
#ifdef ALLOW_NONHYDROSTATIC |
74 |
jmc |
1.17 |
& -deltaTmom*hy3dFac*_recip_dyC(i,j,bi,bj)* |
75 |
adcroft |
1.11 |
& (cg3d_x(i,j,k,bi,bj)-cg3d_x(i,j-1,k,bi,bj)) |
76 |
|
|
#endif |
77 |
adcroft |
1.1 |
& )*_maskS(i,j,k,bi,bj) |
78 |
|
|
gVNm1(i,j,k,bi,bj)=gV(i,j,k,bi,bj) |
79 |
adcroft |
1.4 |
ENDDO |
80 |
|
|
ENDDO |
81 |
adcroft |
1.1 |
|
82 |
|
|
RETURN |
83 |
|
|
END |