1 |
#include "AUTODIFF_OPTIONS.h" |
2 |
|
3 |
subroutine adecco_check_exp( |
4 |
& mythid, mycurrentiter, mycurrenttime, yprefix ) |
5 |
|
6 |
c ================================================================= |
7 |
c SUBROUTINE adecco_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 adecco_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 |
#ifdef ALLOW_AUTODIFF_MONITOR |
73 |
# include "adcommon.h" |
74 |
#endif |
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 |
#ifdef ALLOW_AUTODIFF_MONITOR |
95 |
|
96 |
c-- == local variables == |
97 |
|
98 |
INTEGER bi,bj,i,j |
99 |
integer irec |
100 |
integer mydate(4) |
101 |
character yfname*128 |
102 |
|
103 |
_RS tmpflux (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
104 |
|
105 |
c == end of interface == |
106 |
|
107 |
irec = 0 |
108 |
|
109 |
if ( mod((mycurrentiter),1) .eq. 0 ) then |
110 |
irec = (mycurrentiter)/1 + 1 |
111 |
|
112 |
cph( |
113 |
cph call cal_GetDate( |
114 |
cph I mycurrentiter, |
115 |
cph I mycurrenttime, |
116 |
cph O mydate, |
117 |
cph I mythid |
118 |
cph & ) |
119 |
|
120 |
print *, 'pathei: in check_exp: iter/time/rec/yprefix ', |
121 |
& mycurrentiter, mycurrenttime, irec, ' ', yprefix |
122 |
print *, 'pathei: in check_exp: date ', mycurrentiter |
123 |
print *, 'pathei: in check_exp: adtheta ', adtheta(10,10,1,1,1) |
124 |
print *, 'pathei: in check_exp: adsalt ', adsalt(10,10,1,1,1) |
125 |
print *, 'pathei: in check_exp: aduvel ', aduvel(10,10,1,1,1) |
126 |
print *, 'pathei: in check_exp: advvel ', advvel(10,10,1,1,1) |
127 |
print *, 'pathei: in check_exp: adqnet ', adqnet(10,10,1,1) |
128 |
print *, 'pathei: in check_exp: adempmr ', adempmr(10,10,1,1) |
129 |
print *, 'pathei: in check_exp: adfu ', adfu(10,10,1,1) |
130 |
print *, 'pathei: in check_exp: adfv ', adfv(10,10,1,1) |
131 |
print *, 'pathei: in check_exp: adetan ', adetan(10,10,1,1) |
132 |
#ifdef ALLOW_PASSIVE_TRACER |
133 |
print *, 'pathei: in check_exp: adtr1 ', adtr1(10,10,1,1,1) |
134 |
#endif |
135 |
#ifdef ALLOW_DIFFKR_CONTROL |
136 |
print *, 'pathei: in check_exp: addiffkr', addiffkr(10,10,2,1,1) |
137 |
#endif |
138 |
#ifdef ALLOW_KAPGM_CONTROL |
139 |
print *, 'pathei: in check_exp: adkapgm ', adkapgm(10,10,1,1,1) |
140 |
#endif |
141 |
cph) |
142 |
|
143 |
c-- Potential temperature: |
144 |
write(yfname,'(128a)') ' ' |
145 |
write(yfname,'(2a)') yprefix, 'snapshot_adtheta' |
146 |
call mdswritefield( yfname, 32, .false., |
147 |
& 'RL', nr, adtheta, irec, |
148 |
& mycurrentiter, mythid ) |
149 |
|
150 |
c-- Salinity: |
151 |
write(yfname,'(128a)') ' ' |
152 |
write(yfname,'(2a)') yprefix, 'snapshot_adsalt' |
153 |
call mdswritefield( yfname, 32, .false., |
154 |
& 'RL', nr, adsalt, irec, |
155 |
& mycurrentiter, mythid ) |
156 |
|
157 |
c-- Zonal velocity: |
158 |
write(yfname,'(128a)') ' ' |
159 |
write(yfname,'(2a)') yprefix, 'snapshot_aduvel' |
160 |
call mdswritefield( yfname, 32, .false., |
161 |
& 'RL', nr, aduvel, irec, |
162 |
& mycurrentiter, mythid ) |
163 |
|
164 |
c-- Meridional velocity: |
165 |
write(yfname,'(128a)') ' ' |
166 |
write(yfname,'(2a)') yprefix, 'snapshot_advvel' |
167 |
call mdswritefield( yfname, 32, .false., |
168 |
& 'RL', nr, advvel, irec, |
169 |
& mycurrentiter, mythid ) |
170 |
|
171 |
c-- Surface pressure: |
172 |
write(yfname,'(128a)') ' ' |
173 |
write(yfname,'(2a)') yprefix, 'snapshot_adetan' |
174 |
call mdswritefield( yfname, 32, .false., |
175 |
& 'RL', 1, adetan, irec, |
176 |
& mycurrentiter, mythid ) |
177 |
|
178 |
#ifdef ALLOW_PASSIVE_TRACER |
179 |
c-- Passive tracer: |
180 |
write(yfname,'(128a)') ' ' |
181 |
write(yfname,'(2a)') yprefix, 'snapshot_adtr1' |
182 |
call mdswritefield( yfname, 32, .false., |
183 |
& 'RL', nr, adtr1, irec, |
184 |
& mycurrentiter, mythid ) |
185 |
#endif |
186 |
|
187 |
#ifdef ALLOW_DIFFKR_CONTROL |
188 |
c-- diapycnal diffusion: |
189 |
write(yfname,'(128a)') ' ' |
190 |
write(yfname,'(2a)') yprefix, 'snapshot_addiffkr' |
191 |
call mdswritefield( yfname, 32, .false., |
192 |
& 'RL', nr, addiffkr, irec, |
193 |
& mycurrentiter, mythid ) |
194 |
#endif |
195 |
|
196 |
#ifdef ALLOW_KAPGM_CONTROL |
197 |
c-- isopycnal diffusion: |
198 |
write(yfname,'(128a)') ' ' |
199 |
write(yfname,'(2a)') yprefix, 'snapshot_adkapgm' |
200 |
call mdswritefield( yfname, 32, .false., |
201 |
& 'RL', nr, adkapgm, irec, |
202 |
& mycurrentiter, mythid ) |
203 |
#endif |
204 |
|
205 |
c-- Surface heat flux: |
206 |
DO bj = myByLo(myThid), myByHi(myThid) |
207 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
208 |
DO j=1-oLy,sNy+oLy |
209 |
DO i=1-oLx,sNx+oLx |
210 |
tmpflux(i,j,bi,bj) = |
211 |
& - adQnet(i,j,bi,bj)*HeatCapacity_Cp*rhoConst*dRf(1) |
212 |
ENDDO |
213 |
ENDDO |
214 |
ENDDO |
215 |
ENDDO |
216 |
write(yfname,'(128a)') ' ' |
217 |
write(yfname,'(2a)') yprefix, 'snapshot_adqnet' |
218 |
call mdswritefield( yfname, 32, .false., |
219 |
& 'RS', 1, tmpflux, irec, |
220 |
& mycurrentiter, mythid ) |
221 |
|
222 |
c-- Surface virtual salt flux: |
223 |
DO bj = myByLo(myThid), myByHi(myThid) |
224 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
225 |
DO j=1-oLy,sNy+oLy |
226 |
DO i=1-oLx,sNx+oLx |
227 |
tmpflux(i,j,bi,bj) = |
228 |
& adEmPmR(i,j,bi,bj)*dRf(1)/35. |
229 |
ENDDO |
230 |
ENDDO |
231 |
ENDDO |
232 |
ENDDO |
233 |
write(yfname,'(128a)') ' ' |
234 |
write(yfname,'(2a)') yprefix, 'snapshot_adempmr' |
235 |
call mdswritefield( yfname, 32, .false., |
236 |
& 'RS', 1, tmpflux, irec, |
237 |
& mycurrentiter, mythid ) |
238 |
|
239 |
c-- Surface zonal wind stress: |
240 |
DO bj = myByLo(myThid), myByHi(myThid) |
241 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
242 |
DO j=1-oLy,sNy+oLy |
243 |
DO i=1-oLx,sNx+oLx |
244 |
tmpflux(i,j,bi,bj) = |
245 |
& -adfu(i,j,bi,bj)*rhoConst*dRf(1)/horiVertRatio |
246 |
ENDDO |
247 |
ENDDO |
248 |
ENDDO |
249 |
ENDDO |
250 |
write(yfname,'(128a)') ' ' |
251 |
write(yfname,'(2a)') yprefix, 'snapshot_adfu' |
252 |
call mdswritefield( yfname, 32, .false., |
253 |
& 'RS', 1, tmpflux, irec, |
254 |
& mycurrentiter, mythid ) |
255 |
|
256 |
c-- Surface meridional wind stress: |
257 |
DO bj = myByLo(myThid), myByHi(myThid) |
258 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
259 |
DO j=1-oLy,sNy+oLy |
260 |
DO i=1-oLx,sNx+oLx |
261 |
tmpflux(i,j,bi,bj) = |
262 |
& -adfv(i,j,bi,bj)*rhoConst*dRf(1)/horiVertRatio |
263 |
ENDDO |
264 |
ENDDO |
265 |
ENDDO |
266 |
ENDDO |
267 |
write(yfname,'(128a)') ' ' |
268 |
write(yfname,'(2a)') yprefix, 'snapshot_adfv' |
269 |
call mdswritefield( yfname, 32, .false., |
270 |
& 'RS', 1, tmpflux, irec, |
271 |
& mycurrentiter, mythid ) |
272 |
|
273 |
c-- Control vector contributions: |
274 |
|
275 |
c-- Heat flux (control): |
276 |
cph call mdswritefield( yprefix//'snapshot_xx_hfl', 32, .false., |
277 |
cph & 'RS', 1, xx_hfl, irec, |
278 |
cph & mycurrentiter, mythid ) |
279 |
|
280 |
c-- Virtual salt flux (control): |
281 |
cph call mdswritefield( yprefix//'snapshot_xx_sfl', 32, .false., |
282 |
cph & 'RS', 1, xx_sfl, irec, |
283 |
cph & mycurrentiter, mythid ) |
284 |
|
285 |
c-- Zonal wind stress (control): |
286 |
cph call mdswritefield( yprefix//'snapshot_xx_tauu', 32, .false., |
287 |
cph & 'RS', 1, xx_tauu, irec, |
288 |
cph & mycurrentiter, mythid ) |
289 |
|
290 |
c-- Meridional wind stress (control): |
291 |
cph call mdswritefield( yprefix//'snapshot_xx_tauv', 32, .false., |
292 |
cph & 'RS', 1, xx_tauv, irec, |
293 |
cph & mycurrentiter, mythid ) |
294 |
|
295 |
#ifdef I_HAVE_FIXED_THIS |
296 |
|
297 |
#ifdef ALLOW_BULKFORMULAE |
298 |
c-- Atmospheric zonal wind: |
299 |
write(yfname,'(128a)') ' ' |
300 |
write(yfname,'(2a)') yprefix, 'snapshot_uwind' |
301 |
call mdswritefield( yfname, 32, .false., |
302 |
& 'RS', 1, uwind, irec, |
303 |
& mycurrentiter, mythid ) |
304 |
|
305 |
c-- Atmospheric meridional wind: |
306 |
write(yfname,'(128a)') ' ' |
307 |
write(yfname,'(2a)') yprefix, 'snapshot_vwind' |
308 |
call mdswritefield( yfname, 32, .false., |
309 |
& 'RS', 1,vwind, irec, |
310 |
& mycurrentiter, mythid ) |
311 |
|
312 |
c-- Air temperature: |
313 |
write(yfname,'(128a)') ' ' |
314 |
write(yfname,'(2a)') yprefix, 'snapshot_atemp' |
315 |
call mdswritefield( yfname, 32, .false., |
316 |
& 'RS', 1, atemp, irec, |
317 |
& mycurrentiter, mythid ) |
318 |
|
319 |
c-- Relative humidity: |
320 |
write(yfname,'(128a)') ' ' |
321 |
write(yfname,'(2a)') yprefix, 'snapshot_aqh' |
322 |
call mdswritefield( yfname, 32, .false., |
323 |
& 'RS', 1, aqh, irec, |
324 |
& mycurrentiter, mythid ) |
325 |
|
326 |
c-- Precipitation: |
327 |
write(yfname,'(128a)') ' ' |
328 |
write(yfname,'(2a)') yprefix, 'snapshot_precip' |
329 |
call mdswritefield( yfname, 32, .false., |
330 |
& 'RS', 1, precip, irec, |
331 |
& mycurrentiter, mythid ) |
332 |
|
333 |
#ifndef ALLOW_KPP |
334 |
c-- Short wave radiative flux: |
335 |
write(yfname,'(128a)') ' ' |
336 |
write(yfname,'(2a)') yprefix, 'snapshot_swflux' |
337 |
call mdswritefield( yfname, 32, .false., |
338 |
& 'RS', 1, swflux, irec, |
339 |
& mycurrentiter, mythid ) |
340 |
#endif |
341 |
|
342 |
c-- Long wave radiative flux: |
343 |
write(yfname,'(128a)') ' ' |
344 |
write(yfname,'(2a)') yprefix, 'snapshot_lwflux' |
345 |
call mdswritefield( yfname, 32, .false., |
346 |
& 'RS', 1, lwflux, irec, |
347 |
& mycurrentiter, mythid ) |
348 |
#endif / * ALLOW_BULKFORMULAE * / |
349 |
|
350 |
#endif /* I_HAVE_FIXED_THIS */ |
351 |
|
352 |
endif |
353 |
|
354 |
#endif /* ALLOW_AUTODIFF_MONITOR */ |
355 |
|
356 |
cph#endif / * ALLOW_SNAPSHOTS * / |
357 |
|
358 |
return |
359 |
end |
360 |
|
361 |
|