1 |
C $Header: /u/gcmpack/MITgcm/model/src/diags_rho.F,v 1.4 2008/09/22 17:55:16 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
|
7 |
C-- File diags_rho.F: density & density advection diagnostics |
8 |
C-- Contents |
9 |
C-- o DIAGS_RHO_L |
10 |
C-- o DIAGS_RHO_G |
11 |
|
12 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
13 |
CBOP |
14 |
C !ROUTINE: DIAGS_RHO_L |
15 |
C !INTERFACE: |
16 |
SUBROUTINE DIAGS_RHO_L( |
17 |
I diagRho, k, bi, bj, |
18 |
I rho3d, rhoKm1, wFld, |
19 |
I myTime, myIter, myThid ) |
20 |
C !DESCRIPTION: \bv |
21 |
C *==========================================================* |
22 |
C | S/R DIAGS_RHO_L |
23 |
C | o Density vertical advective term diagnostics |
24 |
C *==========================================================* |
25 |
C | works with local arrays, and called inside k,bi,bj loops |
26 |
C *==========================================================* |
27 |
C \ev |
28 |
|
29 |
C !USES: |
30 |
IMPLICIT NONE |
31 |
C == Global variables == |
32 |
#include "SIZE.h" |
33 |
#include "EEPARAMS.h" |
34 |
#include "PARAMS.h" |
35 |
#include "GRID.h" |
36 |
|
37 |
C !INPUT/OUTPUT PARAMETERS: |
38 |
C == Routine Arguments == |
39 |
C diagRho :: select which diags to fill |
40 |
C k :: level index |
41 |
C bi, bj :: tile indices |
42 |
C rho3d :: in-situ density anomaly |
43 |
C rhoKm1 :: density of water @ level above (k-1), evaluated at pressure level k |
44 |
C wFld :: vertical velocity |
45 |
C myTime :: Current time |
46 |
C myIter :: Iteration number |
47 |
C myThid :: my Thread Id number |
48 |
INTEGER diagRho |
49 |
INTEGER k, bi, bj |
50 |
_RL rho3d (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
51 |
_RL rhoKm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
52 |
_RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
53 |
_RL myTime |
54 |
INTEGER myIter, myThid |
55 |
CEOP |
56 |
|
57 |
#ifdef ALLOW_DIAGNOSTICS |
58 |
|
59 |
C !LOCAL VARIABLES: |
60 |
C == Local variables == |
61 |
C i,j :: Loop counters |
62 |
C tmpFld :: temporary working array |
63 |
INTEGER i,j |
64 |
_RL tmpFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
65 |
|
66 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
67 |
|
68 |
IF ( k.GE.2 .AND. MOD(diagRho,8).GE.4 ) THEN |
69 |
C-- Diagnose Vertical velocity times vertical difference |
70 |
C of potential density reference at level below (i.e. level k) |
71 |
DO j=1,sNy |
72 |
DO i=1,sNx |
73 |
tmpFld(i,j) = wFld(i,j,k,bi,bj) |
74 |
& *( rho3d(i,j,k) - rhoKm1(i,j) )*rkSign |
75 |
ENDDO |
76 |
ENDDO |
77 |
CALL DIAGNOSTICS_FILL(tmpFld,'WdRHO_P ',k,1,2,bi,bj,myThid) |
78 |
IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT('WdRHO_P ',bi,bj,myThid) |
79 |
ENDIF |
80 |
|
81 |
IF ( k.GE.2 .AND. diagRho.GE.8 ) THEN |
82 |
C-- Diagnose Vertical velocity times vertical difference |
83 |
C of density at fixed Temp & Salt (from level above, i.e. level k-1) |
84 |
DO j=1,sNy |
85 |
DO i=1,sNx |
86 |
tmpFld(i,j) = wFld(i,j,k,bi,bj) |
87 |
& *( rhoKm1(i,j) - rho3d(i,j,k-1) )*rkSign |
88 |
ENDDO |
89 |
ENDDO |
90 |
CALL DIAGNOSTICS_FILL(tmpFld,'WdRHOdP ',k,1,2,bi,bj,myThid) |
91 |
IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT('WdRHOdP ',bi,bj,myThid) |
92 |
ENDIF |
93 |
|
94 |
#endif /* ALLOW_DIAGNOSTICS */ |
95 |
|
96 |
RETURN |
97 |
END |
98 |
|
99 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
100 |
CBOP |
101 |
C !ROUTINE: DIAGS_RHO_G |
102 |
C !INTERFACE: |
103 |
SUBROUTINE DIAGS_RHO_G( |
104 |
I rho3d, uFld, vFld, wFld, |
105 |
I myTime, myIter, myThid ) |
106 |
C !DESCRIPTION: \bv |
107 |
C *==========================================================* |
108 |
C | S/R DIAGS_RHO_G |
109 |
C | o Density & Density advective Flux diagnostics |
110 |
C *==========================================================* |
111 |
C | works with global arrays; k,bi,bj loops are done here |
112 |
C *==========================================================* |
113 |
C \ev |
114 |
|
115 |
C !USES: |
116 |
IMPLICIT NONE |
117 |
C == Global variables == |
118 |
#include "SIZE.h" |
119 |
#include "EEPARAMS.h" |
120 |
#include "PARAMS.h" |
121 |
#include "GRID.h" |
122 |
|
123 |
C !INPUT/OUTPUT PARAMETERS: |
124 |
C == Routine Arguments == |
125 |
C rho3d :: in-situ density anomaly |
126 |
C uFld :: zonal velocity |
127 |
C vFld :: meridional velocity |
128 |
C wFld :: vertical velocity |
129 |
C myTime :: Current time |
130 |
C myIter :: Iteration number |
131 |
C myThid :: my Thread Id number |
132 |
_RL rho3d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
133 |
_RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
134 |
_RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
135 |
_RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
136 |
_RL myTime |
137 |
INTEGER myIter, myThid |
138 |
CEOP |
139 |
|
140 |
#ifdef ALLOW_DIAGNOSTICS |
141 |
C !FUNCTIONS: |
142 |
LOGICAL DIAGNOSTICS_IS_ON |
143 |
EXTERNAL DIAGNOSTICS_IS_ON |
144 |
|
145 |
C !LOCAL VARIABLES: |
146 |
C == Local variables == |
147 |
C i,j :: Loop counters |
148 |
C k, bi,bj :: level & tile indices |
149 |
C tmpFld :: temporary working array |
150 |
INTEGER i,j |
151 |
INTEGER k, bi,bj |
152 |
_RL tmpFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
153 |
_RL tmpFac |
154 |
|
155 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
156 |
|
157 |
CALL DIAGNOSTICS_FILL( rho3d, 'RHOAnoma', |
158 |
& 0, Nr, 0, 1, 1, myThid ) |
159 |
tmpFac = 1. _d 0 |
160 |
CALL DIAGNOSTICS_SCALE_FILL( rho3d, tmpFac, 2, |
161 |
& 'RHOANOSQ', 0, Nr, 0, 1, 1, myThid ) |
162 |
|
163 |
IF ( DIAGNOSTICS_IS_ON('URHOMASS',myThid) ) THEN |
164 |
DO bj=myByLo(myThid),myByHi(myThid) |
165 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
166 |
DO k=1,Nr |
167 |
DO j=1,sNy |
168 |
DO i=1,sNx+1 |
169 |
tmpFld(i,j) = uFld(i,j,k,bi,bj)*_hFacW(i,j,k,bi,bj) |
170 |
& *(rho3d(i-1,j,k,bi,bj)+rho3d(i,j,k,bi,bj)) |
171 |
& *0.5 _d 0 |
172 |
ENDDO |
173 |
ENDDO |
174 |
CALL DIAGNOSTICS_FILL(tmpFld,'URHOMASS',k,1,2,bi,bj,myThid) |
175 |
ENDDO |
176 |
ENDDO |
177 |
ENDDO |
178 |
ENDIF |
179 |
|
180 |
IF ( DIAGNOSTICS_IS_ON('VRHOMASS',myThid) ) THEN |
181 |
DO bj=myByLo(myThid),myByHi(myThid) |
182 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
183 |
DO k=1,Nr |
184 |
DO j=1,sNy+1 |
185 |
DO i=1,sNx |
186 |
tmpFld(i,j) = vFld(i,j,k,bi,bj)*_hFacS(i,j,k,bi,bj) |
187 |
& *(rho3d(i,j-1,k,bi,bj)+rho3d(i,j,k,bi,bj)) |
188 |
& *0.5 _d 0 |
189 |
ENDDO |
190 |
ENDDO |
191 |
CALL DIAGNOSTICS_FILL(tmpFld,'VRHOMASS',k,1,2,bi,bj,myThid) |
192 |
ENDDO |
193 |
ENDDO |
194 |
ENDDO |
195 |
ENDIF |
196 |
|
197 |
IF ( DIAGNOSTICS_IS_ON('WRHOMASS',myThid) ) THEN |
198 |
DO bj=myByLo(myThid),myByHi(myThid) |
199 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
200 |
DO k=1,Nr |
201 |
IF ( k.EQ.1 ) THEN |
202 |
DO j=1,sNy |
203 |
DO i=1,sNx |
204 |
tmpFld(i,j) = wFld(i,j,k,bi,bj)*rho3d(i,j,k,bi,bj) |
205 |
ENDDO |
206 |
ENDDO |
207 |
ELSE |
208 |
DO j=1,sNy |
209 |
DO i=1,sNx |
210 |
tmpFld(i,j) = wFld(i,j,k,bi,bj) |
211 |
& *(rho3d(i,j,k-1,bi,bj)+rho3d(i,j,k,bi,bj)) |
212 |
& *0.5 _d 0 |
213 |
ENDDO |
214 |
ENDDO |
215 |
ENDIF |
216 |
CALL DIAGNOSTICS_FILL(tmpFld,'WRHOMASS',k,1,2,bi,bj,myThid) |
217 |
ENDDO |
218 |
ENDDO |
219 |
ENDDO |
220 |
ENDIF |
221 |
|
222 |
#endif /* ALLOW_DIAGNOSTICS */ |
223 |
|
224 |
RETURN |
225 |
END |