1 |
#include "CPP_OPTIONS.h" |
2 |
|
3 |
subroutine ecco_check_exp( |
4 |
& mythid, mycurrentiter, mycurrenttime, yprefix ) |
5 |
|
6 |
c ================================================================= |
7 |
c SUBROUTINE ecco_check_exp |
8 |
c ================================================================= |
9 |
c |
10 |
c o Check details of the model run |
11 |
c |
12 |
c This routine dumps a collection of model fields for diagnostic |
13 |
c or testimg purposes, respectively. |
14 |
c |
15 |
c Variables for experiment 06: |
16 |
c |
17 |
c Dynamical core: |
18 |
c Potential temperature theta [C] |
19 |
c Salinity salt [psu] |
20 |
c Zonal velocity uvel [m/s] |
21 |
c Meridional velocity vvel [m/s] |
22 |
c Vertical velocity ( --> check_fld) rvel [m/s] |
23 |
c Surface pressure cg2d_x [m] |
24 |
c Surface heat flux qnet [K/s] |
25 |
c Qnet contrib. from external forcing tflux [K/s] |
26 |
c Qnet contrib. from relaxation to Levitus qlev [K/s] |
27 |
c Qnet contrib. from relaxation to Reynolds qrey [K/s] |
28 |
c Surface virtual salt flux empmr [psu/s] |
29 |
c Surface zonal wind stress fu [m/s^2] |
30 |
c Surface meridional wind stress fv [m/s^2] |
31 |
c |
32 |
c Control vector contributions: |
33 |
c Heat flux correction xx_hflux [W/m^2] |
34 |
c Virtual salt flux correction xx_sflux [psu/s/m^2] |
35 |
c Zonal wind stress correction xx_tauu [N/m^2] |
36 |
c Meridional wind stress correction xx_tauv [N/m^2] |
37 |
c |
38 |
c Bulk formulae: |
39 |
c Atmospheric zonal wind uwind [m/s] |
40 |
c Atmospheric meridional wind vwind [m/s] |
41 |
c Air temperature atemp [K] |
42 |
c Specific humidity aqh [kg/kg] |
43 |
c Precipitation precip [kg/s/m^2] |
44 |
c Short wave radiative flux swflux/qsw [W/m^2] |
45 |
c Long wave radiative flux lwflux/qlw [W/m^2] |
46 |
c |
47 |
c Non-local K-Profile Parameterization (KPP): |
48 |
c Short wave radiative flux swflux/qsw [W/m^2] |
49 |
c Boundary layer depth kpphbl [m] |
50 |
c |
51 |
c |
52 |
c Beta Version: Christian Eckert (MIT) 15-Nov-1999 |
53 |
c |
54 |
c ================================================================= |
55 |
c SUBROUTINE check_exp |
56 |
c ================================================================= |
57 |
|
58 |
implicit none |
59 |
|
60 |
c-- == global variables == |
61 |
|
62 |
cph#ifdef ALLOW_SNAPSHOTS |
63 |
|
64 |
#include "EEPARAMS.h" |
65 |
#include "SIZE.h" |
66 |
#include "PARAMS.h" |
67 |
cph#include "CG2D_EXTERNAL.h" |
68 |
#include "DYNVARS.h" |
69 |
#include "FFIELDS.h" |
70 |
#include "GRID.h" |
71 |
#include "ctrl.h" |
72 |
cph#include "cal.h" |
73 |
cph#include "exf_clim_param.h" |
74 |
cph#include "exf_fields.h" |
75 |
|
76 |
#ifdef ALLOW_KPP |
77 |
# include "KPP_OPTIONS.h" |
78 |
# include "KPP_PARAMS.h" |
79 |
# include "KPP.h" |
80 |
#endif |
81 |
|
82 |
cph#endif |
83 |
|
84 |
c == routine arguments == |
85 |
c mythid - thread number for this instance of the routine. |
86 |
|
87 |
integer mythid |
88 |
integer mycurrentiter |
89 |
_RL mycurrenttime |
90 |
character yprefix*3 |
91 |
|
92 |
cph#ifdef ALLOW_SNAPSHOTS |
93 |
|
94 |
c-- == local variables == |
95 |
|
96 |
INTEGER bi,bj,i,j |
97 |
integer irec |
98 |
integer mydate(4) |
99 |
character yfname*128 |
100 |
|
101 |
_RS tmpflux (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
102 |
|
103 |
c == end of interface == |
104 |
|
105 |
irec = 0 |
106 |
|
107 |
if ( mod((mycurrentiter),1) .eq. 0 ) then |
108 |
irec = (mycurrentiter)/1 + 1 |
109 |
|
110 |
cph( |
111 |
cph call cal_GetDate( |
112 |
cph I mycurrentiter, |
113 |
cph I mycurrenttime, |
114 |
cph O mydate, |
115 |
cph I mythid |
116 |
cph & ) |
117 |
|
118 |
print *, 'pathei: in check_exp: iter/time/rec/yprefix ', |
119 |
& mycurrentiter, mycurrenttime, irec, ' ', yprefix |
120 |
print *, 'pathei: in check_exp: date ', mycurrentiter |
121 |
print *, 'pathei: in check_exp: theta ', theta(10,10,1,1,1) |
122 |
print *, 'pathei: in check_exp: salt ', salt(10,10,1,1,1) |
123 |
print *, 'pathei: in check_exp: uvel ', uvel(10,10,1,1,1) |
124 |
print *, 'pathei: in check_exp: vvel ', vvel(10,10,1,1,1) |
125 |
print *, 'pathei: in check_exp: qnet ', qnet(10,10,1,1) |
126 |
print *, 'pathei: in check_exp: empmr ', empmr(10,10,1,1) |
127 |
print *, 'pathei: in check_exp: fu ', fu(10,10,1,1) |
128 |
print *, 'pathei: in check_exp: fv ', fv(10,10,1,1) |
129 |
cph) |
130 |
|
131 |
c-- Potential temperature: |
132 |
write(yfname,'(128a)') ' ' |
133 |
write(yfname,'(2a)') yprefix, 'snapshot_theta' |
134 |
call mdswritefield( yfname, 32, .false., |
135 |
& 'RL', nr, theta, irec, |
136 |
& mycurrentiter, mythid ) |
137 |
|
138 |
c-- Salinity: |
139 |
write(yfname,'(128a)') ' ' |
140 |
write(yfname,'(2a)') yprefix, 'snapshot_salt' |
141 |
call mdswritefield( yfname, 32, .false., |
142 |
& 'RL', nr, salt, irec, |
143 |
& mycurrentiter, mythid ) |
144 |
|
145 |
c-- Zonal velocity: |
146 |
write(yfname,'(128a)') ' ' |
147 |
write(yfname,'(2a)') yprefix, 'snapshot_uvel' |
148 |
call mdswritefield( yfname, 32, .false., |
149 |
& 'RL', nr, uvel, irec, |
150 |
& mycurrentiter, mythid ) |
151 |
|
152 |
c-- Meridional velocity: |
153 |
write(yfname,'(128a)') ' ' |
154 |
write(yfname,'(2a)') yprefix, 'snapshot_vvel' |
155 |
call mdswritefield( yfname, 32, .false., |
156 |
& 'RL', nr, vvel, irec, |
157 |
& mycurrentiter, mythid ) |
158 |
|
159 |
c-- Surface pressure: |
160 |
cph write(yfname,'(128a)') ' ' |
161 |
cph write(yfname,'(2a)') yprefix, 'snapshot_cg2d_x' |
162 |
cph call mdswritefield( yfname, 32, .false., |
163 |
cph & 'RL', 1, cg2d_x, irec, |
164 |
cph & mycurrentiter, mythid ) |
165 |
|
166 |
c-- Surface heat flux: |
167 |
DO bj = myByLo(myThid), myByHi(myThid) |
168 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
169 |
DO j=1-oLy,sNy+oLy |
170 |
DO i=1-oLx,sNx+oLx |
171 |
tmpflux(i,j,bi,bj) = |
172 |
& - Qnet(i,j,bi,bj)*HeatCapacity_Cp*rhoNil*dRf(1) |
173 |
ENDDO |
174 |
ENDDO |
175 |
ENDDO |
176 |
ENDDO |
177 |
write(yfname,'(128a)') ' ' |
178 |
write(yfname,'(2a)') yprefix, 'snapshot_qnet' |
179 |
call mdswritefield( yfname, 32, .false., |
180 |
& 'RS', 1, tmpflux, irec, |
181 |
& mycurrentiter, mythid ) |
182 |
|
183 |
c-- Surface virtual salt flux: |
184 |
DO bj = myByLo(myThid), myByHi(myThid) |
185 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
186 |
DO j=1-oLy,sNy+oLy |
187 |
DO i=1-oLx,sNx+oLx |
188 |
tmpflux(i,j,bi,bj) = |
189 |
& EmPmR(i,j,bi,bj)*dRf(1)/35. |
190 |
ENDDO |
191 |
ENDDO |
192 |
ENDDO |
193 |
ENDDO |
194 |
write(yfname,'(128a)') ' ' |
195 |
write(yfname,'(2a)') yprefix, 'snapshot_empmr' |
196 |
call mdswritefield( yfname, 32, .false., |
197 |
& 'RS', 1, tmpflux, irec, |
198 |
& mycurrentiter, mythid ) |
199 |
|
200 |
c-- Surface zonal wind stress: |
201 |
DO bj = myByLo(myThid), myByHi(myThid) |
202 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
203 |
DO j=1-oLy,sNy+oLy |
204 |
DO i=1-oLx,sNx+oLx |
205 |
tmpflux(i,j,bi,bj) = |
206 |
& -fu(i,j,bi,bj)*rhoNil*dRf(1)/horiVertRatio |
207 |
ENDDO |
208 |
ENDDO |
209 |
ENDDO |
210 |
ENDDO |
211 |
write(yfname,'(128a)') ' ' |
212 |
write(yfname,'(2a)') yprefix, 'snapshot_fu' |
213 |
call mdswritefield( yfname, 32, .false., |
214 |
& 'RS', 1, tmpflux, irec, |
215 |
& mycurrentiter, mythid ) |
216 |
|
217 |
c-- Surface meridional wind stress: |
218 |
DO bj = myByLo(myThid), myByHi(myThid) |
219 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
220 |
DO j=1-oLy,sNy+oLy |
221 |
DO i=1-oLx,sNx+oLx |
222 |
tmpflux(i,j,bi,bj) = |
223 |
& -fv(i,j,bi,bj)*rhoNil*dRf(1)/horiVertRatio |
224 |
ENDDO |
225 |
ENDDO |
226 |
ENDDO |
227 |
ENDDO |
228 |
write(yfname,'(128a)') ' ' |
229 |
write(yfname,'(2a)') yprefix, 'snapshot_fv' |
230 |
call mdswritefield( yfname, 32, .false., |
231 |
& 'RS', 1, tmpflux, irec, |
232 |
& mycurrentiter, mythid ) |
233 |
|
234 |
c-- Control vector contributions: |
235 |
|
236 |
c-- Heat flux (control): |
237 |
cph call mdswritefield( yprefix//'snapshot_xx_hfl', 32, .false., |
238 |
cph & 'RS', 1, xx_hfl, irec, |
239 |
cph & mycurrentiter, mythid ) |
240 |
|
241 |
c-- Virtual salt flux (control): |
242 |
cph call mdswritefield( yprefix//'snapshot_xx_sfl', 32, .false., |
243 |
cph & 'RS', 1, xx_sfl, irec, |
244 |
cph & mycurrentiter, mythid ) |
245 |
|
246 |
c-- Zonal wind stress (control): |
247 |
cph call mdswritefield( yprefix//'snapshot_xx_tauu', 32, .false., |
248 |
cph & 'RS', 1, xx_tauu, irec, |
249 |
cph & mycurrentiter, mythid ) |
250 |
|
251 |
c-- Meridional wind stress (control): |
252 |
cph call mdswritefield( yprefix//'snapshot_xx_tauv', 32, .false., |
253 |
cph & 'RS', 1, xx_tauv, irec, |
254 |
cph & mycurrentiter, mythid ) |
255 |
|
256 |
#ifdef ALLOW_BULKFORMULAE |
257 |
c-- Atmospheric zonal wind: |
258 |
write(yfname,'(128a)') ' ' |
259 |
write(yfname,'(2a)') yprefix, 'snapshot_uwind' |
260 |
call mdswritefield( yfname, 32, .false., |
261 |
& 'RS', 1, uwind, irec, |
262 |
& mycurrentiter, mythid ) |
263 |
|
264 |
c-- Atmospheric meridional wind: |
265 |
write(yfname,'(128a)') ' ' |
266 |
write(yfname,'(2a)') yprefix, 'snapshot_vwind' |
267 |
call mdswritefield( yfname, 32, .false., |
268 |
& 'RS', 1,vwind, irec, |
269 |
& mycurrentiter, mythid ) |
270 |
|
271 |
c-- Air temperature: |
272 |
write(yfname,'(128a)') ' ' |
273 |
write(yfname,'(2a)') yprefix, 'snapshot_atemp' |
274 |
call mdswritefield( yfname, 32, .false., |
275 |
& 'RS', 1, atemp, irec, |
276 |
& mycurrentiter, mythid ) |
277 |
|
278 |
c-- Relative humidity: |
279 |
write(yfname,'(128a)') ' ' |
280 |
write(yfname,'(2a)') yprefix, 'snapshot_aqh' |
281 |
call mdswritefield( yfname, 32, .false., |
282 |
& 'RS', 1, aqh, irec, |
283 |
& mycurrentiter, mythid ) |
284 |
|
285 |
c-- Precipitation: |
286 |
write(yfname,'(128a)') ' ' |
287 |
write(yfname,'(2a)') yprefix, 'snapshot_precip' |
288 |
call mdswritefield( yfname, 32, .false., |
289 |
& 'RS', 1, precip, irec, |
290 |
& mycurrentiter, mythid ) |
291 |
|
292 |
#ifndef ALLOW_KPP |
293 |
c-- Short wave radiative flux: |
294 |
write(yfname,'(128a)') ' ' |
295 |
write(yfname,'(2a)') yprefix, 'snapshot_swflux' |
296 |
call mdswritefield( yfname, 32, .false., |
297 |
& 'RS', 1, swflux, irec, |
298 |
& mycurrentiter, mythid ) |
299 |
#endif |
300 |
|
301 |
c-- Long wave radiative flux: |
302 |
write(yfname,'(128a)') ' ' |
303 |
write(yfname,'(2a)') yprefix, 'snapshot_lwflux' |
304 |
call mdswritefield( yfname, 32, .false., |
305 |
& 'RS', 1, lwflux, irec, |
306 |
& mycurrentiter, mythid ) |
307 |
#endif / * ALLOW_BULKFORMULAE * / |
308 |
|
309 |
#ifdef ALLOW_KPP |
310 |
c-- Short wave radiative flux: |
311 |
DO bj = myByLo(myThid), myByHi(myThid) |
312 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
313 |
DO j=1-oLy,sNy+oLy |
314 |
DO i=1-oLx,sNx+oLx |
315 |
tmpflux(i,j,bi,bj) = |
316 |
& -Qsw(i,j,bi,bj)*HeatCapacity_Cp*rhoNil*dRf(1) |
317 |
ENDDO |
318 |
ENDDO |
319 |
ENDDO |
320 |
ENDDO |
321 |
write(yfname,'(128a)') ' ' |
322 |
write(yfname,'(2a)') yprefix, 'snapshot_swflux' |
323 |
call mdswritefield( yfname, 32, .false., |
324 |
& 'RS', 1, tmpflux, irec, |
325 |
& mycurrentiter, mythid ) |
326 |
|
327 |
c-- Boundary layer depth: |
328 |
write(yfname,'(128a)') ' ' |
329 |
write(yfname,'(2a)') yprefix, 'snapshot_kpphbl' |
330 |
call mdswritefield( yfname, 32, .false., |
331 |
& 'RL', 1, kpphbl, irec, |
332 |
& mycurrentiter, mythid ) |
333 |
#endif / * ALLOW_KPP * / |
334 |
|
335 |
#ifdef ALLOW_CLIMSST_RELAXATION |
336 |
c-- SST climatology: |
337 |
write(yfname,'(128a)') ' ' |
338 |
write(yfname,'(2a)') yprefix, 'snapshot_sst' |
339 |
call mdswritefield( yfname, 32, .false., |
340 |
& 'RS', 1, sst, irec, |
341 |
& mycurrentiter, mythid ) |
342 |
#endif / * ALLOW_CLIMSST_RELAXATION * / |
343 |
|
344 |
#ifdef ALLOW_CLIMSSS_RELAXATION |
345 |
c-- SSS climatology: |
346 |
write(yfname,'(128a)') ' ' |
347 |
write(yfname,'(2a)') yprefix, 'snapshot_sss' |
348 |
call mdswritefield( yfname, 32, .false., |
349 |
& 'RS', 1, sss, irec, |
350 |
& mycurrentiter, mythid ) |
351 |
#endif / * ALLOW_CLIMSSS_RELAXATION * / |
352 |
|
353 |
endif |
354 |
|
355 |
cph#endif / * ALLOW_SNAPSHOTS * / |
356 |
|
357 |
return |
358 |
end |
359 |
|
360 |
|