1 |
C $Header: /u/gcmpack/MITgcm/model/src/diags_rho.F,v 1.3 2006/06/07 01:55:12 heimbach 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 k, bi, bj, |
18 |
I rhoK, rhoKm1, wFld, |
19 |
I myTime, myIter, myThid ) |
20 |
C !DESCRIPTION: \bv |
21 |
C *==========================================================* |
22 |
C | S/R DIAGS_RHO_L |
23 |
C | o Vertical buoyancy Flux 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 k :: level index |
40 |
C bi, bj :: tile indices |
41 |
C rhoK :: in-situ density anomaly |
42 |
C rhoKm1 :: density of water @ level above (k-1), evaluated at pressure level k |
43 |
C wFld :: vertical velocity |
44 |
C myTime :: Current time |
45 |
C myIter :: Iteration number |
46 |
C myThid :: my Thread Id number |
47 |
INTEGER k, bi, bj |
48 |
_RL rhoK (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
49 |
_RL rhoKm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
50 |
_RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
51 |
_RL myTime |
52 |
INTEGER myIter, myThid |
53 |
CEOP |
54 |
|
55 |
#ifdef ALLOW_DIAGNOSTICS |
56 |
|
57 |
C !LOCAL VARIABLES: |
58 |
C == Local variables == |
59 |
C i,j :: Loop counters |
60 |
INTEGER i,j |
61 |
_RL tmpFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
62 |
|
63 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
64 |
|
65 |
IF ( k.EQ.1 ) THEN |
66 |
DO j=1,sNy |
67 |
DO i=1,sNx |
68 |
tmpFld(i,j) = wFld(i,j,k,bi,bj)*rhoK(i,j) |
69 |
ENDDO |
70 |
ENDDO |
71 |
ELSE |
72 |
DO j=1,sNy |
73 |
DO i=1,sNx |
74 |
tmpFld(i,j) = wFld(i,j,k,bi,bj) |
75 |
& *(rhoKm1(i,j)+rhoK(i,j))*0.5 _d 0 |
76 |
ENDDO |
77 |
ENDDO |
78 |
ENDIF |
79 |
CALL DIAGNOSTICS_FILL(tmpFld,'WRHOMASS',k,1,2,bi,bj,myThid) |
80 |
|
81 |
#endif /* ALLOW_DIAGNOSTICS */ |
82 |
|
83 |
RETURN |
84 |
END |
85 |
|
86 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
87 |
CBOP |
88 |
C !ROUTINE: DIAGS_RHO_G |
89 |
C !INTERFACE: |
90 |
SUBROUTINE DIAGS_RHO_G( |
91 |
I rho3d, uFld, vFld, |
92 |
I myTime, myIter, myThid ) |
93 |
C !DESCRIPTION: \bv |
94 |
C *==========================================================* |
95 |
C | S/R DIAGS_RHO_G |
96 |
C | o Buoyancy & horizontal buoyancy Flux diagnostics |
97 |
C *==========================================================* |
98 |
C | works with global arrays; k,bi,bj loops are done here |
99 |
C *==========================================================* |
100 |
C \ev |
101 |
|
102 |
C !USES: |
103 |
IMPLICIT NONE |
104 |
C == Global variables == |
105 |
#include "SIZE.h" |
106 |
#include "EEPARAMS.h" |
107 |
#include "PARAMS.h" |
108 |
#include "GRID.h" |
109 |
|
110 |
C !INPUT/OUTPUT PARAMETERS: |
111 |
C == Routine Arguments == |
112 |
C rho3d :: in-situ density anomaly |
113 |
C uFld :: zonal velocity |
114 |
C vFld :: meridional velocity |
115 |
C myTime :: Current time |
116 |
C myIter :: Iteration number |
117 |
C myThid :: my Thread Id number |
118 |
_RL rho3d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
119 |
_RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
120 |
_RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
121 |
_RL myTime |
122 |
INTEGER myIter, myThid |
123 |
CEOP |
124 |
|
125 |
#ifdef ALLOW_DIAGNOSTICS |
126 |
|
127 |
C !LOCAL VARIABLES: |
128 |
C == Local variables == |
129 |
C i,j :: Loop counters |
130 |
C k, bi,bj :: level & tile indices |
131 |
INTEGER i,j |
132 |
INTEGER k, bi,bj |
133 |
_RL tmpFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
134 |
_RL tmpFac |
135 |
LOGICAL DIAGNOSTICS_IS_ON |
136 |
EXTERNAL DIAGNOSTICS_IS_ON |
137 |
|
138 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
139 |
|
140 |
CALL DIAGNOSTICS_FILL( rho3d, 'RHOAnoma', |
141 |
& 0, Nr, 0, 1, 1, myThid ) |
142 |
tmpFac = 1. _d 0 |
143 |
CALL DIAGNOSTICS_SCALE_FILL( rho3d, tmpFac, 2, |
144 |
& 'RHOANOSQ', 0, Nr, 0, 1, 1, myThid ) |
145 |
|
146 |
IF ( DIAGNOSTICS_IS_ON('URHOMASS',myThid) ) THEN |
147 |
DO bj=myByLo(myThid),myByHi(myThid) |
148 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
149 |
DO k=1,Nr |
150 |
DO j=1,sNy |
151 |
DO i=1,sNx+1 |
152 |
tmpFld(i,j) = uFld(i,j,k,bi,bj)*_hFacW(i,j,k,bi,bj) |
153 |
& *(rho3d(i-1,j,k,bi,bj)+rho3d(i,j,k,bi,bj)) |
154 |
& *0.5 _d 0 |
155 |
ENDDO |
156 |
ENDDO |
157 |
CALL DIAGNOSTICS_FILL(tmpFld,'URHOMASS',k,1,2,bi,bj,myThid) |
158 |
ENDDO |
159 |
ENDDO |
160 |
ENDDO |
161 |
ENDIF |
162 |
|
163 |
IF ( DIAGNOSTICS_IS_ON('VRHOMASS',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+1 |
168 |
DO i=1,sNx |
169 |
tmpFld(i,j) = vFld(i,j,k,bi,bj)*_hFacS(i,j,k,bi,bj) |
170 |
& *(rho3d(i,j-1,k,bi,bj)+rho3d(i,j,k,bi,bj)) |
171 |
& *0.5 _d 0 |
172 |
ENDDO |
173 |
ENDDO |
174 |
CALL DIAGNOSTICS_FILL(tmpFld,'VRHOMASS',k,1,2,bi,bj,myThid) |
175 |
ENDDO |
176 |
ENDDO |
177 |
ENDDO |
178 |
ENDIF |
179 |
|
180 |
#endif /* ALLOW_DIAGNOSTICS */ |
181 |
|
182 |
RETURN |
183 |
END |