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