/[MITgcm]/MITgcm/model/src/diags_phi_rlow.F
ViewVC logotype

Annotation of /MITgcm/model/src/diags_phi_rlow.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.4 - (hide annotations) (download)
Mon Oct 26 21:48:43 2009 UTC (14 years, 7 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint62, checkpoint61z, checkpoint61y
Changes since 1.3: +1 -7 lines
initialize phihydlow to avoid TAF recomputation
(in the context of bottom pressure (GRACE) cost term)
done correctly, this time...

1 gforget 1.4 C $Header: /u/gcmpack/MITgcm/model/src/diags_phi_rlow.F,v 1.3 2009/10/26 00:43:38 gforget Exp $
2 jmc 1.1 C $Name: $
3    
4     #include "CPP_OPTIONS.h"
5    
6     CBOP
7     C !ROUTINE: DIAGS_PHI_RLOW
8     C !INTERFACE:
9     SUBROUTINE DIAGS_PHI_RLOW(
10     I k, bi, bj, iMin,iMax, jMin,jMax,
11 jmc 1.2 I phiHydF, phiHydC, alphRho, tFld, sFld,
12 jmc 1.1 I myTime, myIter, myThid)
13     C !DESCRIPTION: \bv
14     C *==========================================================*
15     C | S/R DIAGS_PHI_RLOW
16     C | o Diagnose Phi-Hydrostatic at r-lower boundary
17     C | = bottom pressure (ocean in z-coord) ;
18     C | = sea surface elevation (ocean in p-coord) ;
19     C | = height at the top of atmosphere (in p-coord) ;
20     C *==========================================================*
21     C \ev
22    
23     C !USES:
24     IMPLICIT NONE
25     C == Global variables ==
26     #include "SIZE.h"
27     #include "EEPARAMS.h"
28     #include "PARAMS.h"
29     #include "GRID.h"
30     #include "SURFACE.h"
31     #include "DYNVARS.h"
32    
33     C !INPUT/OUTPUT PARAMETERS:
34     C == Routine Arguments ==
35     C bi,bj :: tile index
36     C iMin,iMax,jMin,jMax :: Loop counters
37 jmc 1.2 C phiHydF :: hydrostatic potential anomaly at middle between
38     C 2 centers k & k+1 (interface k+1)
39     C phiHydC :: hydrostatic potential anomaly at cell center
40 jmc 1.1 C (atmos: =Geopotential ; ocean-z: =Pressure/rho)
41     C alphRho :: Density (z-coord) or specific volume (p-coord)
42     C tFld :: Potential temp.
43     C sFld :: Salinity
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 jmc 1.2 _RL phiHydF(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
49     _RL phiHydC(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
50 jmc 1.1 _RL alphRho(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
51     _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
52     _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
53     _RL myTime
54     INTEGER myIter, myThid
55    
56     #ifdef INCLUDE_PHIHYD_CALCULATION_CODE
57    
58     C !LOCAL VARIABLES:
59     C == Local variables ==
60     C i,j :: Loop counters
61     INTEGER i,j
62     _RL zero, one, half
63 jmc 1.2 _RL ddRloc, ratioRm, ratioRp
64 jmc 1.1 PARAMETER ( zero= 0. _d 0 , one= 1. _d 0 , half= .5 _d 0 )
65     CEOP
66    
67 jmc 1.2 IF ( buoyancyRelation .EQ. 'OCEANIC' ) THEN
68 jmc 1.1
69 jmc 1.2 C----- Compute bottom pressure deviation from gravity*rho0*H
70     C Start from phiHyd at the (bottom) tracer point and add Del_h*g*rho'
71     C with Del_h = distance from the bottom up to tracer point
72 jmc 1.1
73     IF (integr_GeoPot.EQ.1) THEN
74     C -- Finite Volume Form
75    
76     DO j=jMin,jMax
77     DO i=iMin,iMax
78     IF ( k .EQ. kLowC(i,j,bi,bj) ) THEN
79 jmc 1.2 ddRloc = rC(k)-R_low(i,j,bi,bj)
80     phiHydLow(i,j,bi,bj) = phiHydC(i,j)
81     & + ddRloc*gravity*alphRho(i,j)*recip_rhoConst
82 jmc 1.1 ENDIF
83     ENDDO
84     ENDDO
85    
86     ELSE
87     C -- Finite Difference Form
88    
89 jmc 1.2 ratioRm = one
90     ratioRp = one
91     IF (k.GT.1 ) ratioRm = half*drC(k)/(rF(k)-rC(k))
92     IF (k.LT.Nr) ratioRp = half*drC(k+1)/(rC(k)-rF(k+1))
93 jmc 1.1
94     DO j=jMin,jMax
95     DO i=iMin,iMax
96 jmc 1.2 IF ( k .EQ. kLowC(i,j,bi,bj) ) THEN
97     ddRloc = rC(k)-R_low(i,j,bi,bj)
98     phiHydLow(i,j,bi,bj) = phiHydC(i,j)
99     & +( MIN(zero,ddRloc)*ratioRm
100     & +MAX(zero,ddRloc)*ratioRp
101     & )*gravity*alphRho(i,j)*recip_rhoConst
102 jmc 1.1 ENDIF
103     ENDDO
104     ENDDO
105    
106     C -- end if integr_GeoPot = ...
107     ENDIF
108    
109 jmc 1.2 C -- end buoyancyR = Oceanic (z)
110     ENDIF
111    
112     IF (k.EQ.Nr) THEN
113 jmc 1.1 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
114 jmc 1.2 C -- last level (bottom): rescale (r*) and add surface contribution
115 jmc 1.1
116 jmc 1.2 IF ( buoyancyRelation .EQ. 'OCEANICP' .OR.
117     & buoyancyRelation .EQ. 'ATMOSPHERIC' ) THEN
118     C -- P coordinate : Phi(R_low) is simply at the top :
119     DO j=jMin,jMax
120     DO i=iMin,iMax
121     phiHydLow(i,j,bi,bj) = phiHydF(i,j)
122 jmc 1.1 ENDDO
123 jmc 1.2 ENDDO
124     ENDIF
125 jmc 1.1
126 jmc 1.2 DO j=jMin,jMax
127     DO i=iMin,iMax
128     phiHydLow(i,j,bi,bj) = phiHydLow(i,j,bi,bj)
129     & + Bo_surf(i,j,bi,bj)*etaN(i,j,bi,bj)
130     & + phi0surf(i,j,bi,bj)
131 jmc 1.1 ENDDO
132 jmc 1.2 ENDDO
133 jmc 1.1
134     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
135 jmc 1.2 C -- end if k=Nr.
136 jmc 1.1 ENDIF
137    
138     #endif /* INCLUDE_PHIHYD_CALCULATION_CODE */
139    
140     RETURN
141     END

  ViewVC Help
Powered by ViewVC 1.1.22