1 |
jmc |
1.1 |
C $Header: /u/gcmpack/MITgcm/model/src/calc_grad_phi_hyd.F,v 1.14 2011/12/14 01:12:59 jmc Exp $ |
2 |
|
|
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "CPP_OPTIONS.h" |
5 |
|
|
|
6 |
|
|
CBOP |
7 |
|
|
C !ROUTINE: CALC_GRAD_PHI_FV |
8 |
|
|
C !INTERFACE: |
9 |
|
|
SUBROUTINE CALC_GRAD_PHI_FV( |
10 |
|
|
I k, bi, bj, iMin,iMax, jMin,jMax, |
11 |
|
|
I phiHydF, phiHydU, pKappaF, pKappaU, |
12 |
|
|
O dPhiHydX, dPhiHydY, |
13 |
|
|
I myTime, myIter, myThid) |
14 |
|
|
C !DESCRIPTION: \bv |
15 |
|
|
C *==========================================================* |
16 |
|
|
C | S/R CALC_GRAD_PHI_FV |
17 |
|
|
C | o Calculate the gradient of Hydrostatic pressure anomaly |
18 |
|
|
C | using Finite-Volume method from S.-J. Lin (QJRMS, 1997) |
19 |
|
|
C *==========================================================* |
20 |
|
|
C | o used with sigma-coords - might be useful also with r* |
21 |
|
|
C *==========================================================* |
22 |
|
|
C \ev |
23 |
|
|
|
24 |
|
|
C !USES: |
25 |
|
|
IMPLICIT NONE |
26 |
|
|
C == Global variables == |
27 |
|
|
#include "SIZE.h" |
28 |
|
|
#include "EEPARAMS.h" |
29 |
|
|
#include "PARAMS.h" |
30 |
|
|
#include "GRID.h" |
31 |
|
|
#include "SURFACE.h" |
32 |
|
|
c#include "DYNVARS.h" |
33 |
|
|
|
34 |
|
|
C !INPUT/OUTPUT PARAMETERS: |
35 |
|
|
C == Routine Arguments == |
36 |
|
|
C bi,bj :: tile index |
37 |
|
|
C iMin,iMax,jMin,jMax :: Loop counters |
38 |
|
|
C phiHydF :: hydrostatic potential anomaly at interface k |
39 |
|
|
C (atmos: =Geopotential ; ocean-z: =Pressure/rho) |
40 |
|
|
C phiHydU :: hydrostatic potential anomaly at interface k+1 (Upper k) |
41 |
|
|
C pKappaF :: (p/Po)^kappa at interface k |
42 |
|
|
C pKappaU :: (p/Po)^kappa at interface k+1 (Upper k) |
43 |
|
|
C dPhiHydX,Y :: Gradient (X & Y directions) of Hyd. Potential |
44 |
|
|
C myTime :: Current time |
45 |
|
|
C myIter :: Current iteration number |
46 |
|
|
C myThid :: Instance number for this call of the routine. |
47 |
|
|
INTEGER k, bi,bj, iMin,iMax, jMin,jMax |
48 |
|
|
_RL phiHydF(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
49 |
|
|
_RL phiHydU(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
50 |
|
|
_RL pKappaF(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
51 |
|
|
_RL pKappaU(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
52 |
|
|
_RL dPhiHydX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
53 |
|
|
_RL dPhiHydY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
54 |
|
|
_RL myTime |
55 |
|
|
INTEGER myIter, myThid |
56 |
|
|
|
57 |
|
|
#ifdef INCLUDE_PHIHYD_CALCULATION_CODE |
58 |
|
|
|
59 |
|
|
C !LOCAL VARIABLES: |
60 |
|
|
C == Local variables == |
61 |
|
|
C i,j :: Loop counters |
62 |
|
|
INTEGER i,j |
63 |
|
|
_RL dpk_dip, dpk_dim |
64 |
|
|
_RL dpk_djp, dpk_djm |
65 |
|
|
c CHARACTER*(MAX_LEN_MBUF) msgBuf |
66 |
|
|
CEOP |
67 |
|
|
|
68 |
|
|
C-- Zonal & Meridional gradient of potential anomaly |
69 |
|
|
DO j=1-OLy,sNy+OLy |
70 |
|
|
DO i=1-OLx,sNx+OLx |
71 |
|
|
dPhiHydX(i,j) = 0. _d 0 |
72 |
|
|
dPhiHydY(i,j) = 0. _d 0 |
73 |
|
|
ENDDO |
74 |
|
|
ENDDO |
75 |
|
|
DO j=jMin,jMax |
76 |
|
|
DO i=iMin+1,iMax |
77 |
|
|
dpk_dip = pKappaF( i ,j) - pKappaU(i-1,j) |
78 |
|
|
dpk_dim = pKappaF(i-1,j) - pKappaU( i ,j) |
79 |
|
|
dPhiHydX(i,j) = ( phi0surf(i,j,bi,bj) - phi0surf(i-1,j,bi,bj) ) |
80 |
|
|
& + ( dpk_dip*( phiHydU(i,j) - phiHydF(i-1,j) ) |
81 |
|
|
& + dpk_dim*( phiHydF(i,j) - phiHydU(i-1,j) ) |
82 |
|
|
& )/( dpk_dip + dpk_dim ) |
83 |
|
|
dPhiHydX(i,j) = _recip_dxC(i,j,bi,bj)*dPhiHydX(i,j) |
84 |
|
|
c & * recip_deepFacC(k)*recip_rhoFacC(k) |
85 |
|
|
ENDDO |
86 |
|
|
ENDDO |
87 |
|
|
DO j=jMin+1,jMax |
88 |
|
|
DO i=iMin,iMax |
89 |
|
|
dpk_djp = pKappaF(i, j ) - pKappaU(i,j-1) |
90 |
|
|
dpk_djm = pKappaF(i,j-1) - pKappaU(i, j ) |
91 |
|
|
dPhiHydY(i,j) = ( phi0surf(i,j,bi,bj) - phi0surf(i,j-1,bi,bj) ) |
92 |
|
|
& + ( dpk_djp*( phiHydU(i,j) - phiHydF(i,j-1) ) |
93 |
|
|
& + dpk_djm*( phiHydF(i,j) - phiHydU(i,j-1) ) |
94 |
|
|
& )/( dpk_djp + dpk_djm ) |
95 |
|
|
dPhiHydY(i,j) = _recip_dyC(i,j,bi,bj)*dPhiHydY(i,j) |
96 |
|
|
c & * recip_deepFacC(k)*recip_rhoFacC(k) |
97 |
|
|
ENDDO |
98 |
|
|
ENDDO |
99 |
|
|
|
100 |
|
|
C-- Apply mask: |
101 |
|
|
DO j=1-OLy,sNy+OLy |
102 |
|
|
DO i=1-OLx,sNx+OLx |
103 |
|
|
dPhiHydX(i,j) = dPhiHydX(i,j)*_maskW(i,j,k,bi,bj) |
104 |
|
|
dPhiHydY(i,j) = dPhiHydY(i,j)*_maskS(i,j,k,bi,bj) |
105 |
|
|
ENDDO |
106 |
|
|
ENDDO |
107 |
|
|
|
108 |
|
|
#endif /* INCLUDE_PHIHYD_CALCULATION_CODE */ |
109 |
|
|
|
110 |
|
|
RETURN |
111 |
|
|
END |