/[MITgcm]/MITgcm/pkg/ecco/ecco_check_exp.F
ViewVC logotype

Annotation of /MITgcm/pkg/ecco/ecco_check_exp.F

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


Revision 1.3 - (hide annotations) (download)
Wed Sep 25 19:36:50 2002 UTC (21 years, 8 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint46n_post, checkpoint51k_post, checkpoint47j_post, checkpoint48d_pre, checkpoint51j_post, branch-exfmods-tag, checkpoint47e_post, checkpoint47f_post, checkpoint46j_post, checkpoint47c_post, checkpoint50e_post, checkpoint50c_post, checkpoint46i_post, checkpoint51n_pre, checkpoint47d_post, checkpoint47a_post, checkpoint47h_post, checkpoint48e_post, checkpoint48h_post, checkpoint48d_post, checkpoint50c_pre, branchpoint-genmake2, checkpoint46k_post, checkpoint50d_pre, checkpoint48f_post, checkpoint51r_post, checkpoint47i_post, checkpoint51o_pre, checkpoint46j_pre, checkpoint51i_post, checkpoint48c_post, checkpoint46l_post, checkpoint51e_post, checkpoint51b_post, checkpoint51l_pre, checkpoint51c_post, checkpoint47d_pre, checkpoint47, checkpoint48, checkpoint49, checkpoint47b_post, checkpoint51l_post, checkpoint48i_post, checkpoint51o_post, checkpoint46l_pre, checkpoint51f_pre, checkpoint51q_post, checkpoint51, checkpoint50, checkpoint50d_post, checkpoint51b_pre, checkpoint47g_post, checkpoint46m_post, checkpoint48a_post, checkpoint51h_pre, checkpoint50g_post, checkpoint50b_pre, checkpoint51g_post, checkpoint51f_post, checkpoint48b_post, checkpoint50b_post, checkpoint50f_post, checkpoint50a_post, checkpoint50f_pre, checkpoint51d_post, checkpoint48c_pre, checkpoint51m_post, checkpoint51t_post, checkpoint50h_post, checkpoint51a_post, checkpoint46h_post, checkpoint50e_pre, checkpoint50i_post, checkpoint51p_post, checkpoint51n_post, checkpoint48g_post, checkpoint51i_pre, checkpoint51s_post
Branch point for: branch-genmake2, branch-exfmods-curt, branch-nonh, tg2-branch, checkpoint51n_branch
Changes since 1.2: +4 -4 lines
o cleaned up the use of rhoNil and rhoConst.
  - rhoNil should only appear in the LINEAR equation of state, everywhere
    else rhoNil is replaced by rhoConst, e.g. find_rho computes rho-rhoConst
    and the dynamical equations are all divided by rhoConst
o introduced new parameter rhoConstFresh, a reference density of fresh
  water, to remove the fresh water flux's dependence on rhoNil. The default
  value is 999.8 kg/m^3
o cleanup up external_forcing.F and external_forcing_surf.F
  - can now be used by both OCEANIC and OCEANICP

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

  ViewVC Help
Powered by ViewVC 1.1.22