| 1 |
ksnow |
1.1 |
C $Header: /u/gcmpack/MITgcm/model/src/calc_div_ghat.F,v 1.27 2012/11/09 22:37:05 jmc Exp $ |
| 2 |
|
|
C $Name: $ |
| 3 |
|
|
|
| 4 |
|
|
#include "CPP_OPTIONS.h" |
| 5 |
|
|
#ifdef ALLOW_PRESSURE_RELEASE_CODE |
| 6 |
|
|
# include "SHELFICE_OPTIONS.h" |
| 7 |
|
|
#endif |
| 8 |
|
|
|
| 9 |
|
|
|
| 10 |
|
|
CBOP |
| 11 |
|
|
C !ROUTINE: CALC_DIV_GHAT |
| 12 |
|
|
C !INTERFACE: |
| 13 |
|
|
SUBROUTINE ADJUST_GHAT_PRESS_RELEASE( |
| 14 |
|
|
U cg2d_b,bi,bj, |
| 15 |
|
|
I myThid) |
| 16 |
|
|
C !DESCRIPTION: \bv |
| 17 |
|
|
C *==========================================================* |
| 18 |
|
|
C | S/R CALC_DIV_GHAT |
| 19 |
|
|
C | o Form the right hand-side of the surface pressure eqn. |
| 20 |
|
|
C *==========================================================* |
| 21 |
|
|
C | Right hand side of pressure equation is divergence |
| 22 |
|
|
C | of veclocity tendency (GHAT) term along with a relaxation |
| 23 |
|
|
C | term equal to the barotropic flow field divergence. |
| 24 |
|
|
C *==========================================================* |
| 25 |
|
|
C \ev |
| 26 |
|
|
|
| 27 |
|
|
C !USES: |
| 28 |
|
|
IMPLICIT NONE |
| 29 |
|
|
C == Global variables == |
| 30 |
|
|
#include "SIZE.h" |
| 31 |
|
|
#include "EEPARAMS.h" |
| 32 |
|
|
#include "PARAMS.h" |
| 33 |
|
|
#include "GRID.h" |
| 34 |
|
|
#include "DYNVARS.h" |
| 35 |
|
|
#include "SURFACE.h" |
| 36 |
|
|
#ifdef ALLOW_ADDFLUID |
| 37 |
|
|
# include "FFIELDS.h" |
| 38 |
|
|
#endif |
| 39 |
|
|
#ifdef ALLOW_PRESSURE_RELEASE_CODE |
| 40 |
|
|
# include "SHELFICE.h" |
| 41 |
|
|
#endif |
| 42 |
|
|
|
| 43 |
|
|
C !INPUT/OUTPUT PARAMETERS: |
| 44 |
|
|
C == Routine arguments == |
| 45 |
|
|
C bi, bj :: tile indices |
| 46 |
|
|
C k :: Index of layer. |
| 47 |
|
|
C cg2d_b :: Conjugate Gradient 2-D solver : Right-hand side vector |
| 48 |
|
|
C cg3d_b :: Conjugate Gradient 3-D solver : Right-hand side vector |
| 49 |
|
|
C myThid :: Instance number for this call of CALC_DIV_GHAT |
| 50 |
|
|
INTEGER bi,bj |
| 51 |
|
|
_RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
| 52 |
|
|
INTEGER myThid |
| 53 |
|
|
|
| 54 |
|
|
C !LOCAL VARIABLES: |
| 55 |
|
|
C == Local variables == |
| 56 |
|
|
C i,j :: Loop counters |
| 57 |
|
|
C xA, yA :: Cell vertical face areas |
| 58 |
|
|
C pf :: Intermediate array for building RHS source term. |
| 59 |
|
|
INTEGER i,j, k, km1 |
| 60 |
|
|
_RL pf (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 61 |
|
|
#ifdef ALLOW_PRESSURE_RELEASE_CODE |
| 62 |
|
|
_RL curPrelVisc, avgGrd |
| 63 |
|
|
_RL depth, eff_depth, pinit, phiLow, phiLowM1 |
| 64 |
|
|
#endif |
| 65 |
|
|
CEOP |
| 66 |
|
|
|
| 67 |
|
|
#ifdef ALLOW_PRESSURE_RELEASE_CODE |
| 68 |
|
|
#ifndef ALLOW_NONHYDROSTATIC |
| 69 |
|
|
|
| 70 |
|
|
IF (implicDiv2Dflow.EQ.1) THEN |
| 71 |
|
|
|
| 72 |
|
|
C Calculate vertical face areas |
| 73 |
|
|
! xA(i,j) = _dyG(i,j,bi,bj)*depthColW(i,j,bi,bj)*rhoFacC(k) |
| 74 |
|
|
! yA(i,j) = _dxG(i,j,bi,bj)*depthColS(i,j,bi,bj)*rhoFacC(k) |
| 75 |
|
|
|
| 76 |
|
|
C-- Pressure equation source term |
| 77 |
|
|
C Term is the vertical integral of the divergence of the |
| 78 |
|
|
C time tendency terms along with a relaxation term that |
| 79 |
|
|
C pulls div(U) + dh/dt back toward zero. |
| 80 |
|
|
|
| 81 |
|
|
DO j=1,sNy+1 |
| 82 |
|
|
DO i=1,sNx+1 |
| 83 |
|
|
pf(i,j)=0.0 |
| 84 |
|
|
ENDDO |
| 85 |
|
|
ENDDO |
| 86 |
|
|
DO j=1,sNy |
| 87 |
|
|
DO i=1,sNx+1 |
| 88 |
|
|
IF (maskInC(i,j,bi,bj)*maskInC(i-1,j,bi,bj).eq.1.) THEN |
| 89 |
|
|
! pf(i,j) = xA(i,j)*gU(i,j,k,bi,bj) / deltaTMom |
| 90 |
|
|
depth = depthColW(i,j,bi,bj) |
| 91 |
|
|
IF (depth.lt.cg2dminColumnEps) THEN |
| 92 |
|
|
|
| 93 |
|
|
phiLow = phiHydLowC (i,j,bi,bj) |
| 94 |
|
|
phiLowM1 = phiHydLowC (i-1,j,bi,bj) |
| 95 |
|
|
|
| 96 |
|
|
avgGrd = min (grdFactor(i,j,bi,bj),grdFactor(i-1,j,bi,bj)) |
| 97 |
|
|
curPrelVisc = max((0.55+avgGrd*0.45),0. _d 0)*pReleaseVisc |
| 98 |
|
|
|
| 99 |
|
|
eff_depth = cg2dminColumnEps - |
| 100 |
|
|
& 2. * cg2dminColumnEps / PI * |
| 101 |
|
|
& COS (PI * depth / (2.*cg2dminColumnEps)) |
| 102 |
|
|
|
| 103 |
|
|
pinit = -1. * curPrelVisc * |
| 104 |
|
|
& _recip_dxC(i,j,bi,bj) * _dyG(i,j,bi,bj) |
| 105 |
|
|
& *(phi0surf(i,j,bi,bj)-phi0surf(i-1,j,bi,bj)+ |
| 106 |
|
|
& phiLow-phiLowM1) |
| 107 |
|
|
& /deltaTmom |
| 108 |
|
|
|
| 109 |
|
|
|
| 110 |
|
|
|
| 111 |
|
|
! pinit = -1. * pReleaseVisc * |
| 112 |
|
|
! & _recip_dxC(i,j,bi,bj) * _dyG(i,j,bi,bj) |
| 113 |
|
|
! & *(phiLow-phiLowM1) |
| 114 |
|
|
! & /deltaTmom |
| 115 |
|
|
|
| 116 |
|
|
pf(i,j) = pinit * (eff_depth - depth) |
| 117 |
|
|
|
| 118 |
|
|
ENDIF |
| 119 |
|
|
ENDIF |
| 120 |
|
|
ENDDO |
| 121 |
|
|
ENDDO |
| 122 |
|
|
DO j=1,sNy |
| 123 |
|
|
DO i=1,sNx |
| 124 |
|
|
cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) + |
| 125 |
|
|
& pf(i+1,j)-pf(i,j) |
| 126 |
|
|
ENDDO |
| 127 |
|
|
ENDDO |
| 128 |
|
|
|
| 129 |
|
|
! call write_fld_xy_rl('cg2db_1a','',cg2d_b,0,mythid) |
| 130 |
|
|
|
| 131 |
|
|
DO j=1,sNy+1 |
| 132 |
|
|
DO i=1,sNx+1 |
| 133 |
|
|
pf(i,j)=0.0 |
| 134 |
|
|
ENDDO |
| 135 |
|
|
ENDDO |
| 136 |
|
|
|
| 137 |
|
|
DO j=1,sNy+1 |
| 138 |
|
|
DO i=1,sNx |
| 139 |
|
|
IF (maskInC(i,j,bi,bj)*maskInC(i,j-1,bi,bj).eq.1.) THEN |
| 140 |
|
|
! pf(i,j) = xA(i,j)*gU(i,j,k,bi,bj) / deltaTMom |
| 141 |
|
|
depth = depthColS(i,j,bi,bj) |
| 142 |
|
|
IF (depth.lt.cg2dminColumnEps) THEN |
| 143 |
|
|
|
| 144 |
|
|
phiLow = phiHydLowC (i,j,bi,bj) |
| 145 |
|
|
phiLowM1 = phiHydLowC (i,j-1,bi,bj) |
| 146 |
|
|
|
| 147 |
|
|
avgGrd = min (grdFactor(i,j,bi,bj),grdFactor(i,j-1,bi,bj)) |
| 148 |
|
|
curPrelVisc = max((0.55+avgGrd*0.45),0. _d 0)*pReleaseVisc |
| 149 |
|
|
|
| 150 |
|
|
eff_depth = cg2dminColumnEps - |
| 151 |
|
|
& 2. * cg2dminColumnEps / PI * |
| 152 |
|
|
& COS (PI * depth / (2.*cg2dminColumnEps)) |
| 153 |
|
|
|
| 154 |
|
|
pinit = -1. * curPrelVisc * |
| 155 |
|
|
& _recip_dyC(i,j,bi,bj) * _dxG(i,j,bi,bj) |
| 156 |
|
|
& *(phi0surf(i,j,bi,bj)-phi0surf(i,j-1,bi,bj)+ |
| 157 |
|
|
& phiLow-phiLowM1) |
| 158 |
|
|
& /deltaTmom |
| 159 |
|
|
|
| 160 |
|
|
! pinit = -1. * pReleaseVisc * |
| 161 |
|
|
! & _recip_dyC(i,j,bi,bj) * _dxG(i,j,bi,bj) |
| 162 |
|
|
! & *(phiLow-phiLowM1) |
| 163 |
|
|
! & /deltaTmom |
| 164 |
|
|
|
| 165 |
|
|
pf(i,j) = pinit * (eff_depth - depth) |
| 166 |
|
|
|
| 167 |
|
|
! print *, "GOT HERE ADJUST", i,j,k,philow, philowm1 |
| 168 |
|
|
|
| 169 |
|
|
ENDIF |
| 170 |
|
|
ENDIF |
| 171 |
|
|
ENDDO |
| 172 |
|
|
ENDDO |
| 173 |
|
|
DO j=1,sNy |
| 174 |
|
|
DO i=1,sNx |
| 175 |
|
|
cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) + |
| 176 |
|
|
& pf(i,j+1)-pf(i,j) |
| 177 |
|
|
ENDDO |
| 178 |
|
|
ENDDO |
| 179 |
|
|
|
| 180 |
|
|
ENDIF |
| 181 |
|
|
|
| 182 |
|
|
#endif |
| 183 |
|
|
#endif |
| 184 |
|
|
|
| 185 |
|
|
RETURN |
| 186 |
|
|
END |