/[MITgcm]/MITgcm/pkg/autodiff/adecco_check_exp.F
ViewVC logotype

Contents of /MITgcm/pkg/autodiff/adecco_check_exp.F

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


Revision 1.3 - (show annotations) (download)
Tue Nov 20 23:27:29 2001 UTC (22 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint44b_post, checkpoint43a-release1mods, checkpoint44h_pre, checkpoint44f_pre, release1-branch_tutorials, checkpoint44g_post, checkpoint44h_post, checkpoint44e_post, checkpoint45d_post, checkpoint45b_post, chkpt44a_pre, release1-branch-end, release1_final_v1, checkpoint44, chkpt44d_post, chkpt44a_post, checkpoint44b_pre, checkpoint45c_post, checkpoint44e_pre, checkpoint45a_post, release1-branch_branchpoint, checkpoint45, chkpt44c_pre, checkpoint44f_post, chkpt44c_post
Branch point for: release1_final, release1-branch
Changes since 1.2: +2 -2 lines
o Bugfix in adcommon.h: commen blocks were adjusted to latest
  common block structure in DYNVARS.h
o placed a do_field_blocking_exchanges after dummy_in_stepping
  to ensure that addummy_in_stepping is preceded by exchanges.

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

  ViewVC Help
Powered by ViewVC 1.1.22