forked from KuangLab-Harvard/SAM_SRCv6.11
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathstepout.f90
328 lines (266 loc) · 9.78 KB
/
stepout.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
subroutine stepout(nstatsteps)
use vars
use rad, only: qrad
use sgs, only: tk
use tracers
use microphysics, only: micro_print, micro_statistics
use sgs, only: sgs_print, sgs_statistics
use stat_moments, only: compmoments
use movies, only: mvmovies
use params
use hbuffer
use instrument_diagnostics, only : isccp_write, modis_write, misr_write, zero_instr_diag
use slm_vars, only: slm_stepout
implicit none
integer i,j,k,ic,jc,nstatsteps
real div, divmax(1), divmin(1), divmax1(1), divmin1(1)
real rdx, rdy, rdz, coef
integer im,jm,km
real wmax, qnmax(1), qnmax1(1)
real(8) buffer(5), buffer1(5)
if(nstatsteps.lt.0) goto 999
if(mod(nstep,nstatis).eq.0) then
call statistics()
call micro_statistics()
call sgs_statistics()
end if
if(mod(nstep,nmovie).eq.0.and.nstep.ge.nmoviestart &
.and.nstep.le.nmovieend) then
call mvmovies()
endif
if(mod(nstep,nstatmom).eq.0.and.nstep.ge.nstatmomstart &
.and.nstep.le.nstatmomend) then
call compmoments()
endif
if(mod(nstep,nstat).eq.0) then
if(masterproc) print *,'Writting statistics:nstatsteps=',nstatsteps
call t_startf ('stat_out')
! write out instrument simulator data and compute domain means for statistics
!bloss(2016-02-06): This needs to be called before hbuf_avg for the simulator
! cloud fraction to be recorded at the first statistics output.
call isccp_write()
call modis_write()
call misr_write()
call hbuf_average(nstatsteps)
call hbuf_write(nstatsteps)
call hbuf_flush()
nstatsteps = 0
call zero_instr_diag()
call t_stopf ('stat_out')
endif
if(mod(nstep,nstat*(1+nrestart_skip)).eq.0.or.nstep.eq.nstop.or.nelapse.eq.0) then
! Kuang Ensemble run: turn on mpi before writing restart (Song Qiyu, 2022)
if(dompiensemble) dompi = .true.
call write_all() ! save restart file
! Kuang Ensemble run: turn off mpi after writing restart (Song Qiyu, 2022)
if(dompiensemble) dompi = .false.
end if
call t_startf ('2D_out')
if(mod(nstep,nsave2D).eq.0.and.nstep.ge.nsave2Dstart &
.and.nstep.le.nsave2Dend) then
! Kuang Ensemble run: turn on mpi before writing 2D output (Song Qiyu, 2022)
if(dompiensemble) dompi = .true.
call write_fields2D()
! Kuang Ensemble run: turn off mpi after writing 2D output (Song Qiyu, 2022)
if(dompiensemble) dompi = .false.
endif
if(.not.save2Davg.or.nstep.eq.nsave2Dstart-nsave2D) call stat_2Dinit(0) ! argument of 0 means storage terms for statistics are preserved
call t_stopf ('2D_out')
call t_startf ('3D_out')
if(mod(nstep,nsave3D).eq.0.and.nstep.ge.nsave3Dstart.and.nstep.le.nsave3Dend ) then
! Kuang Ensemble run: turn on mpi before writing 3D output (Song Qiyu, 2022)
if(dompiensemble) dompi = .true.
! determine if the maximum cloud water exceeds the threshold
! value to save 3D fields:
qnmax(1)=0.
do k=1,nzm
do j=1,ny
do i=1,nx
qnmax(1) = max(qnmax(1),qcl(i,j,k))
qnmax(1) = max(qnmax(1),qci(i,j,k))
end do
enddo
enddo
if(dompi) then
call task_max_real(qnmax(1),qnmax1(1),1)
qnmax(1) = qnmax1(1)
end if
if(qnmax(1).ge.qnsave3D) then
call write_fields3D()
end if
! Kuang Ensemble run: turn off mpi after writing 3D output (Song Qiyu, 2022)
if(dompiensemble) dompi = .false.
endif
call t_stopf ('3D_out')
!------------------------------------------------------------------------------
!------------------------------------------------------------------------------
! Print stuff out:
999 continue
call t_startf ('print_out')
if(masterproc) write(*,'(a,i7,a,i2,a,f6.3,a,f6.3,a,f6.3,a,f6.3)') 'NSTEP = ',nstep,' NCYCLE=',ncycle, &
' CFL=',cfl_adv,' CFLH=',cfl_advh,' CFLZ=',cfl_advz,' CFL_sgs=',cfl_sgs
if(mod(nstep,nprint).eq.0) then
divmin(1)=1.e20
divmax(1)=-1.e20
rdx = 1./dx
rdy = 1./dy
wmax=0.
do k=1,nzm
coef = rho(k)*adz(k)*dz
rdz = 1./coef
if(ny.ne.1) then
do j=1,ny-1*YES3D
jc = j+1*YES3D
do i=1,nx-1
ic = i+1
div = (u(ic,j,k)-u(i,j,k))*rdx + (v(i,jc,k)-v(i,j,k))*rdy + &
(w(i,j,k+1)*rhow(k+1)-w(i,j,k)*rhow(k))*rdz
divmax(1) = max(divmax(1),div)
divmin(1) = min(divmin(1),div)
if(w(i,j,k).gt.wmax) then
wmax=w(i,j,k)
im=i
jm=j
km=k
endif
end do
end do
else
j = 1
do i=1,nx-1
ic = i+1
div = (u(ic,j,k)-u(i,j,k))*rdx +(w(i,j,k+1)*rhow(k+1)-w(i,j,k)*rhow(k))*rdz
divmax(1) = max(divmax(1),div)
divmin(1) = min(divmin(1),div)
if(w(i,j,k).gt.wmax) then
wmax=w(i,j,k)
im=i
jm=j
km=k
endif
end do
endif
end do
if(dompi) then
call task_max_real(divmax(1),divmax1(1),1)
call task_min_real(divmin(1),divmin1(1),1)
divmax(1) = divmax1(1)
divmin(1) = divmin1(1)
end if
! MPI Ensemble run: If MPI ensemble then need to sum up across all domains (Nat, 2022)
if(dompiensemble) dompi = .true.
if(dompi) then
buffer(1) = total_water_before
buffer(2) = total_water_after
buffer(3) = total_water_evap
buffer(4) = total_water_prec
buffer(5) = total_water_ls
call task_sum_real8(buffer, buffer1,5)
total_water_before = buffer1(1)
total_water_after = buffer1(2)
total_water_evap = buffer1(3)
total_water_prec = buffer1(4)
total_water_ls = buffer1(5)
end if
! MPI Ensemble run: If MPI ensemble then need to sum up across all domains (Nat, 2022)
if(dompiensemble) dompi = .false.
!print*,rank,minval(u(1:nx,1:ny,:)),maxval(u(1:nx,1:ny,:))
!print*,rank,'min:',minloc(u(1:nx,1:ny,:))
!print*,rank,'max:',maxloc(u(1:nx,1:ny,:))
!if(masterproc) then
!print*,'--->',tk(27,1,1)
!print*,'tk->:'
!write(6,'(16f7.2)')((tk(i,1,k),i=1,16),k=nzm,1,-1)
!print*,'p->:'
!write(6,'(16f7.2)')((p(i,1,k),i=1,16),k=nzm,1,-1)
!print*,'u->:'
!write(6,'(16f7.2)')((u(i,1,k),i=1,16),k=nzm,1,-1)
!print*,'v->:'
!write(6,'(16f7.2)')((v(i,1,k),i=1,16),k=nzm,1,-1)
!print*,'w->:'
!write(6,'(16f7.2)')((w(i,1,k),i=1,16),k=nzm,1,-1)
!print*,'qcl:'
!write(6,'(16f7.2)')((qcl(i,1,k)*1000.,i=16,31),k=nzm,1,-1)
!print*,'qpl->:'
!write(6,'(16f7.2)')((qpl(i,1,k)*1000.,i=16,31),k=nzm,1,-1)
!print*,'qcl:->'
!write(6,'(16f7.2)')((qci(i,1,k)*1000.,i=16,31),k=nzm,1,-1)
!print*,'qpl->:'
!write(6,'(16f7.2)')((qpi(i,1,k)*1000.,i=16,31),k=nzm,1,-1)
!print*,'qrad->:'
!write(6,'(16f7.2)')((qrad(i,1,k)*3600.,i=16,31),k=nzm,1,-1)
!print*,'qv->:'
!write(6,'(16f7.2)')((qv(i,1,k)*1000.,i=16,31),k=nzm,1,-1)
!print*,'t->:'
!write(6,'(16f7.2)')((t(i,1,k),i=1,16),k=nzm,1,-1)
!print*,'tabs->:'
!write(6,'(16f7.2)')((tabs(i,1,k),i=16,31),k=nzm,1,-1)
!end if
!--------------------------------------------------------
if(masterproc) then
print*,'DAY = ',day
write(6,*) 'NSTEP=',nstep
write(6,*) 'div:',divmax,divmin
if(.not.dodynamicocean) write(6,*) 'SST=',tabs_s
write(6,*) 'surface pressure=',pres0
endif
call fminmax_print('u:',u+ug,dimx1_u,dimx2_u,dimy1_u,dimy2_u,nzm)
call fminmax_print('v:',v+vg,dimx1_v,dimx2_v,dimy1_v,dimy2_v,nzm)
call fminmax_print('w:',w,dimx1_w,dimx2_w,dimy1_w,dimy2_w,nz)
call fminmax_print('p:',p,0,nx,1-YES3D,ny,nzm)
call fminmax_print('t:',t,dimx1_s,dimx2_s,dimy1_s,dimy2_s,nzm)
call fminmax_print('tabs:',tabs,1,nx,1,ny,nzm)
call fminmax_print('qv:',qv,1,nx,1,ny,nzm)
if(dosgs) call sgs_print()
if(docloud) then
call fminmax_print('qcl:',qcl,1,nx,1,ny,nzm)
call fminmax_print('qci:',qci,1,nx,1,ny,nzm)
call micro_print()
end if
if(doprecip) then
call fminmax_print('qpl:',qpl,1,nx,1,ny,nzm)
call fminmax_print('qpi:',qpi,1,nx,1,ny,nzm)
end if
if(dolongwave.or.doshortwave) call fminmax_print('qrad(K/day):',qrad*86400.,1,nx,1,ny,nzm)
if(dotracers) then
do k=1,ntracers
call fminmax_print(trim(tracername(k))//':',tracer(:,:,:,k),dimx1_s,dimx2_s,dimy1_s,dimy2_s,nzm)
end do
end if
call fminmax_print('shf:',fluxbt*cp*rhow(1),1,nx,1,ny,1)
call fminmax_print('lhf:',fluxbq*lcond*rhow(1),1,nx,1,ny,1)
call fminmax_print('uw:',fluxbu,1,nx,1,ny,1)
call fminmax_print('vw:',fluxbv,1,nx,1,ny,1)
call fminmax_print('sst:',sstxy+t00,1,nx,1,ny,1)
call fminmax_print('precinst (mm/h):',precinst*3600.,1,nx,1,ny,1)
if(masterproc.and..not.dosmoke.and.nstatsteps.ne.-1) then
total_water_before = total_water_before/float(nx_gl*ny_gl)
total_water_after = total_water_after/float(nx_gl*ny_gl)
total_water_evap = total_water_evap/float(nx_gl*ny_gl)
total_water_prec = total_water_prec/float(nx_gl*ny_gl)
total_water_ls = total_water_ls/float(nx_gl*ny_gl)
print*,'total water budget:'
write(*,991) total_water_before !'before (mm): ',total_water_before
write(*,992) total_water_after !'after (mm) : ',total_water_after
write(*,993) total_water_evap !'evap (mm/day): ',total_water_evap
write(*,994) total_water_prec !'prec (mm/day): ',total_water_prec
write(*,995) total_water_ls !'ls (mm/day): ',total_water_ls
write(*,996) (total_water_after-(total_water_before+total_water_evap+total_water_ls-total_water_prec))
991 format(' before (mm): ',F16.11)
992 format(' after (mm): ',F16.11)
993 format(' evaporation (mm): ',F16.11)
994 format(' precipitation (mm):',F16.11)
995 format(' large-scale (mm): ',F16.11)
996 format(' Imbalance (mm) ',F16.11)
print*,' imbalance (rel error):', &
(total_water_after-(total_water_before+total_water_evap+total_water_ls-total_water_prec))/total_water_after
print*,'evap (mm/day):',total_water_evap/dt*86400.
print*,'prec (mm/day):',total_water_prec/dt*86400.
print*,'ls (mm/day):',total_water_ls/dt*86400.
print*,'imbalance (mm/day)', &
(total_water_after-total_water_before-total_water_evap-total_water_ls+total_water_prec)/dt*86400.
end if
if(SLM) call slm_stepout()
end if ! (mod(nstep,nprint).eq.0)
call t_stopf ('print_out')
end