-
Notifications
You must be signed in to change notification settings - Fork 55
/
Copy pathmodgenstat.f90
1823 lines (1563 loc) · 61.6 KB
/
modgenstat.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
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!> \file modgenstat.f90
!! Genstat calculates slab averages of several variables
!>
!! Genstat calculates slab averages of several variables
!>
!! Written to fields.expnr, moments.expnr and flux1.expnr and flux2.expnr
!! If netcdf is true, this module leads the profiles.expnr.nc output
!! \author Hans Cuijpers, K.N.M.I.
!! \author Pier Siebesma, K.N.M.I.
!! \author Stephan de Roode,TU Delft
!! \author Chiel van Heerwaarden, Wageningen U.R.
!! \author Thijs Heus,MPI-M
!! \par Revision list
!! \todo Documentation
! This file is part of DALES.
!
! DALES is free software; you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation; either version 3 of the License, or
! (at your option) any later version.
!
! DALES is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with this program. If not, see <http://www.gnu.org/licenses/>.
!
! Copyright 1993-2009 Delft University of Technology, Wageningen University, Utrecht University, KNMI
!
module modgenstat
!-----------------------------------------------------------------|
! |
!*** *stattend* calculates generic slabaveraged statistics |
! |
! Pier Siebesma K.N.M.I. 12/05/1995 |
! Hans Cuijpers I.M.A.U. 21/06/1995 |
! |
! purpose. |
! -------- |
! |
! stattend.f calculates: |
! |
! * time averaged fieldsof theta_l, theta, theta_v, qt, qv |
! ql, u and v |
! * time averaged tendencies of theta_l, theta, qt, qv, ql, |
! u and v |
! * time averaged turbulent fluxes of theta_l, theta_v, |
! theta, qt, qv, ql, u and v |
! * variances of qt, w, u, theta_l and theta_v |
! * skewness of qt and w |
!*** *genstat* calculates timeseries of several variables |
! |
!____________________SETTINGS_AND_SWITCHES________________________|
! IN &NAMTIMESTAT |
! |
! dtav SAMPLING INTERVAL |
! |
! timeav INTERVAL OF WRITING |
! |
! lstat SWITCH TO ENABLE TIMESERIES |
!-----------------------------------------------------------------|
use modprecision
implicit none
! private
PUBLIC :: initgenstat, genstat, exitgenstat
save
!NetCDF variables
integer :: nvar = 47
integer :: ncid,nrec = 0
character(80) :: fname = 'profiles.xxx.nc'
character(80),allocatable, dimension(:,:) :: ncname
character(80),dimension(1,4) :: tncname
real :: dtav, timeav
integer(kind=longint) :: idtav,itimeav,tnext,tnextwrite
logical :: lstat= .false. ! switch for conditional sampling cloud (on/off)
integer :: nsamples
! ---- total fields ---
real, allocatable :: umn (:) ,vmn (:)
real, allocatable :: thlmn (:) ,thvmn (:)
real, allocatable :: qtmn (:) ,qlmn (:), qlhmn(:),cfracmn(:),hurmn(:),tamn(:)
real, allocatable :: clwmn(:), climn(:), plwmn(:), plimn(:)
! real, allocatable :: --- fluxes (resolved, subgrid and total) ---
real, allocatable :: wthlsmn (:),wthlrmn (:),wthltmn(:)
real, allocatable :: wthvsmn (:),wthvrmn (:),wthvtmn(:)
real, allocatable :: wqlsmn (:),wqlrmn (:),wqltmn(:)
real, allocatable :: wqtsmn (:),wqtrmn (:),wqttmn(:)
real, allocatable :: wsvsmn (:,:),wsvrmn(:,:),wsvtmn(:,:)
real, allocatable :: uwtmn (:),vwtmn (:) !total uw, vw
real, allocatable :: uwsmn (:),vwsmn (:) !resolved uw, vw
real, allocatable :: uwrmn (:),vwrmn (:) !subgrid uw, vw
! real, allocatable :: --- various moments ---
! real, allocatable :: rmn (:), r2mn (:), r3mn (:), rhmn (:)
real, allocatable :: w2mn (:), skewmn (:)
real, allocatable :: w2submn (:)
real, allocatable :: u2mn (:), v2mn (:), qt2mn(:)
real, allocatable :: thl2mn (:), thv2mn(:), th2mn(:), ql2mn(:)
! real, allocatable :: qs2mn (:), qsmn (:)
real, allocatable :: svmmn(:,:),svptmn(:,:),svplsmn(:,:),svpmn(:,:)
real, allocatable :: sv2mn(:,:)
real(field_r), allocatable :: umav (:) ! slab averaged ql_0 at full level
real(field_r), allocatable :: vmav (:) ! slab averaged ql_0 at full level
real(field_r), allocatable :: thlmav (:) ! slab averaged ql_0 at full level
real(field_r), allocatable :: thmav (:) ! slab averaged ql_0 at full level
real(field_r), allocatable :: qtmav (:) ! slab averaged ql_0 at full level
real(field_r), allocatable :: qlmav (:) ! slab averaged ql_0 at full level
real, allocatable :: cfracav (:) ! slab averaged ql_0 at full level
real, allocatable :: hurav (:)
real, allocatable :: clwav(:), cliav(:), plwav(:), pliav(:)
real(field_r), allocatable :: svmav (:,:) ! slab averaged ql_0 at full level
real(field_r), allocatable :: taav (:)
real, allocatable :: svpav(:,:) ! slab average total tendency of sv(n)
real, allocatable :: svptav(:,:) ! slab average tendency of sv(n) due to turb.
real, allocatable :: uptav(:) ! slab averaged tendency for u
real, allocatable :: vptav(:) ! slab averaged tendency for v
real, allocatable :: thptav(:) ! slab averaged turbulence tendency of theta
real, allocatable :: qlptav(:) ! slab averaged turbulence tendency of q_liq
real, allocatable :: uwtot (:) ! slab averaged tot w-u flux at half levels
real, allocatable :: vwtot (:) ! slab averaged tot w-v flux at half levels
real, allocatable :: uwsub (:) ! slab averaged sgs w-u flux at half levels
real, allocatable :: vwsub (:) ! slab averaged sgs w-v flux at half levels
real, allocatable :: uwres (:) ! slab averaged res w-u flux at half levels
real, allocatable :: vwres (:) ! slab averaged res w-v flux at half levels
real, allocatable :: thvhav(:)
real, allocatable :: th0av(:)
real, allocatable :: wthlsub(:) ! slab averaged sub w-theta_l flux at half levels
real, allocatable :: wthlres(:) ! slab averaged res w-theta_l flux at half levels
real, allocatable :: wthltot(:) ! slab averaged tot w-theta_l flux at half levels
real, allocatable :: wqtsub(:) ! slab averaged sub w-qtot flux at half levels
real, allocatable :: wqtres(:) ! slab averaged res w-qtot flux at half levels
real, allocatable :: wqttot(:) ! slab averaged tot w-qtot flux at half levels
real, allocatable :: wqlsub(:) ! slab averaged sub w-ql flux at half levels
real, allocatable :: wqlres(:) ! slab averaged tot w-ql flux at half levels
real, allocatable :: wqltot(:) ! slab averaged tot w-ql flux at half levels
real, allocatable :: wthvtot(:) ! slab averaged total wthv-flux
real, allocatable :: wthvres(:) ! slab averaged res wthv-flux
real, allocatable :: wthvsub(:) ! slab averaged sub wthv-flux
real, allocatable :: wsvsub(:,:)! slab averaged sub w-sv(n) flux
real, allocatable :: wsvres(:,:)! slab averaged res w-sv(n) flux
real, allocatable :: wsvtot(:,:)! slab averaged tot w-sv(n) flux
real, allocatable :: cszmn(:) ! Smagorinsky constant
real, allocatable :: cszav(:) ! Smagorinsky constant
real, allocatable :: qlmnlast(:)
real, allocatable :: wthvtmnlast(:)
contains
subroutine initgenstat
use modmpi, only : myid,mpierr, comm3d, D_MPI_BCAST
use modglobal, only : kmax,k1, nsv,ifnamopt,fname_options, ifoutput,&
cexpnr,dtav_glob,timeav_glob,dt_lim,btime,tres,lwarmstart,checknamelisterror
use modstat_nc, only : lnetcdf, open_nc,define_nc,ncinfo,nctiminfo,writestat_dims_nc
use modsurfdata, only : isurf,ksoilmax
implicit none
integer n, ierr
character(40) :: name
character(3) :: csvname
namelist/NAMGENSTAT/ &
dtav,timeav,lstat
dtav=dtav_glob;timeav=timeav_glob
if(myid==0)then
open(ifnamopt,file=fname_options,status='old',iostat=ierr)
read (ifnamopt,NAMGENSTAT,iostat=ierr)
call checknamelisterror(ierr, ifnamopt, 'NAMGENSTAT')
write(6 ,NAMGENSTAT)
close(ifnamopt)
end if
call D_MPI_BCAST(timeav ,1,0,comm3d,mpierr)
call D_MPI_BCAST(dtav ,1,0,comm3d,mpierr)
call D_MPI_BCAST(lstat ,1,0,comm3d,mpierr)
idtav = int(dtav/tres,kind=longint)
itimeav = int(timeav/tres,kind=longint)
tnext = idtav +btime
tnextwrite = itimeav +btime
nsamples = int(itimeav/idtav)
if(.not.(lstat)) return
dt_lim = min(dt_lim,tnext)
if (abs(timeav/dtav-nsamples)>1e-4) then
stop 'timeav must be a integer multiple of dtav'
end if
allocate(umn(k1) ,vmn (k1))
allocate(thlmn (k1) ,thvmn (k1))
allocate(qtmn (k1) ,qlmn (k1), qlhmn(k1),cfracmn(k1), hurmn(k1), tamn(k1))
allocate(clwmn(k1), climn(k1), plwmn(k1), plimn(k1))
allocate(wthlsmn (k1),wthlrmn (k1),wthltmn(k1))
allocate(wthvsmn (k1),wthvrmn (k1),wthvtmn(k1))
allocate(wqlsmn (k1),wqlrmn (k1),wqltmn(k1))
allocate(wqtsmn (k1),wqtrmn (k1),wqttmn(k1))
allocate(wsvsmn (k1,nsv),wsvrmn(k1,nsv),wsvtmn(k1,nsv))
allocate(uwtmn (k1),vwtmn (k1),uwrmn (k1), vwrmn(k1),uwsmn(k1),vwsmn(k1))
! allocate(rmn(k1), r2mn (k1), r3mn (k1), rhmn (k1))
allocate(w2mn (k1), skewmn (k1))
allocate(w2submn (k1))
allocate(u2mn (k1), v2mn (k1), qt2mn(k1))
allocate(thl2mn (k1), thv2mn(k1), th2mn(k1), ql2mn(k1))
! allocate(qs2mn (k1), qsmn (k1))
allocate(svmmn(k1,nsv),svptmn(k1,nsv),svplsmn(k1,nsv),svpmn(k1,nsv))
allocate(sv2mn(k1,nsv))
allocate(umav (k1))
allocate(vmav (k1))
allocate(thlmav (k1))
allocate(thmav (k1))
allocate(qtmav (k1))
allocate(qlmav (k1))
allocate(cfracav(k1))
allocate(hurav(k1))
allocate(clwav(k1), cliav(k1), plwav(k1), pliav(k1))
allocate(taav(k1))
allocate(svmav (k1,nsv))
allocate(uptav(k1))
allocate(vptav(k1))
allocate(thptav(k1))
allocate(qlptav(k1))
allocate(uwtot (k1))
allocate(vwtot (k1))
allocate(uwsub (k1))
allocate(vwsub (k1))
allocate(uwres (k1))
allocate(vwres (k1))
allocate(wthlsub(k1))
allocate(wthlres(k1))
allocate(wthltot(k1))
allocate(wqtsub(k1))
allocate(wqtres(k1))
allocate(wqttot(k1))
allocate(wqlsub(k1))
allocate(wqlres(k1))
allocate(wqltot(k1))
allocate(wthvtot(k1))
allocate(wthvres(k1))
allocate(wthvsub(k1))
allocate(wsvsub(k1,nsv))
allocate(wsvres(k1,nsv))
allocate(wsvtot(k1,nsv))
allocate(thvhav(k1))
allocate(th0av(k1))
allocate(svptav(k1,nsv))
allocate(svpav(k1,nsv))
allocate(cszmn(k1), cszav(k1))
allocate(qlmnlast(k1))
allocate(wthvtmnlast(k1))
umn = 0.
vmn = 0.
thlmn = 0.
thvmn = 0.
qtmn = 0.
qlmn = 0.
qlhmn = 0.
cfracmn = 0.
hurmn = 0.
clwmn = 0.
climn = 0.
plwmn = 0.
plimn = 0.
tamn = 0.
wthlsmn = 0.
wthlrmn = 0.
wthltmn = 0.
wthvsmn = 0.
wthvrmn = 0.
wthvtmn = 0.
wqtsmn = 0.
wqtrmn = 0.
wqttmn = 0.
wqlsmn = 0.
wqlrmn = 0.
wqltmn = 0.
uwtmn = 0.
vwtmn = 0.
uwrmn = 0.
vwrmn = 0.
uwsmn = 0.
vwsmn = 0.
u2mn = 0.
v2mn = 0.
w2mn = 0.
w2submn = 0.
skewmn = 0.
qt2mn = 0.
thl2mn = 0.
thv2mn = 0.
th2mn = 0.
ql2mn = 0.
! qs2mn = 0.
! qsmn = 0.
! rhmn = 0.
! rmn = 0.
! r2mn = 0.
! r3mn = 0.
svmmn = 0.
svpmn = 0.
svpav = 0.
! svplsmn = 0.
svptmn = 0.
svptav = 0.
sv2mn = 0.
wsvsmn = 0.
wsvrmn = 0.
wsvtmn = 0.
cszav = 0.
cszmn = 0.
qlmnlast = 0.
wthvtmnlast = 0.
if(myid==0)then
if (.not. lwarmstart) then
open (ifoutput,file='field.'//cexpnr,status='replace')
close (ifoutput)
open (ifoutput,file='flux1.'//cexpnr,status='replace')
close (ifoutput)
open (ifoutput,file='flux2.'//cexpnr,status='replace')
close (ifoutput)
open (ifoutput,file='moments.'//cexpnr,status='replace')
close (ifoutput)
do n=1,nsv
name = 'svnnnfld.'//cexpnr
write (name(3:5),'(i3.3)') n
open (ifoutput,file=name,status='replace')
close (ifoutput)
end do
do n=1,nsv
name = 'svnnnflx.'//cexpnr
write (name(3:5),'(i3.3)') n
open (ifoutput,file=name,status='replace')
close (ifoutput)
end do
end if
if (lnetcdf) then
fname(10:12) = cexpnr
nvar = nvar + 7*nsv
allocate(ncname(nvar,4))
call nctiminfo(tncname(1,:))
call ncinfo(ncname( 1,:),'rhof','Full level slab averaged density','kg/m^3','tt')
call ncinfo(ncname( 2,:),'rhobf','Full level base-state density','kg/m^3','tt')
call ncinfo(ncname( 3,:),'rhobh','Half level base-state density','kg/m^3','mt')
call ncinfo(ncname( 4,:),'presh','Pressure at cell center','Pa','tt')
call ncinfo(ncname( 5,:),'u','West-East velocity','m/s','tt')
call ncinfo(ncname( 6,:),'v','South-North velocity','m/s','tt')
call ncinfo(ncname( 7,:),'thl','Liquid water potential temperature','K','tt')
call ncinfo(ncname( 8,:),'thv','Virtual potential temperature','K','tt')
call ncinfo(ncname( 9,:),'qt','Total water specific humidity','kg/kg','tt')
call ncinfo(ncname(10,:),'ql','Liquid water specific humidity','kg/kg','tt')
call ncinfo(ncname(11,:),'wthls','SFS-Theta_l flux','Km/s','mt')
call ncinfo(ncname(12,:),'wthlr','Resolved Theta_l flux','Km/s','mt')
call ncinfo(ncname(13,:),'wthlt','Total Theta_l flux','Km/s','mt')
call ncinfo(ncname(14,:),'wthvs','SFS-buoyancy flux','Km/s','mt')
call ncinfo(ncname(15,:),'wthvr','Resolved buoyancy flux','Km/s','mt')
call ncinfo(ncname(16,:),'wthvt','Total buoyancy flux','Km/s','mt')
call ncinfo(ncname(17,:),'wqts','SFS-moisture flux','kg/kg m/s','mt')
call ncinfo(ncname(18,:),'wqtr','Resolved moisture flux','kg/kg m/s','mt')
call ncinfo(ncname(19,:),'wqtt','Total moisture flux','kg/kg m/s','mt')
call ncinfo(ncname(20,:),'wqls','SFS-liquid water flux','kg/kg m/s','mt')
call ncinfo(ncname(21,:),'wqlr','Resolved liquid water flux','kg/kg m/s','mt')
call ncinfo(ncname(22,:),'wqlt','Total liquid water flux','kg/kg m/s','mt')
call ncinfo(ncname(23,:),'uws','SFS-momentum flux (uw)','m^2/s^2','mt')
call ncinfo(ncname(24,:),'uwr','Resolved momentum flux (uw)','m^2/s^2','mt')
call ncinfo(ncname(25,:),'uwt','Total momentum flux (uw)','m^2/s^2','mt')
call ncinfo(ncname(26,:),'vws','SFS-momentum flux (vw)','m^2/s^2','mt')
call ncinfo(ncname(27,:),'vwr','Resolved momentum flux (vw)','m^2/s^2','mt')
call ncinfo(ncname(28,:),'vwt','Total momentum flux (vw)','m^2/s^2','mt')
call ncinfo(ncname(29,:),'w2s','SFS-TKE','m^2/s^2','mt')
call ncinfo(ncname(30,:),'w2r','Resolved vertical velocity variance','m^2/s^2','mt')
!call ncinfo(ncname(31,:),'w2t','Total vertical velocity variance','m^2/s^2','mt')
call ncinfo(ncname(31,:),'skew','vertical velocity skewness','-','mt')
call ncinfo(ncname(32,:),'u2r','Resolved horizontal velocity variance (u)','m^2/s^2','tt')
call ncinfo(ncname(33,:),'v2r','Resolved horizontal velocity variance (v)','m^2/s^2','tt')
call ncinfo(ncname(34,:),'thl2r','Resolved theta_l variance','K^2','tt')
call ncinfo(ncname(35,:),'thv2r','Resolved buoyancy variance','K^2','tt')
call ncinfo(ncname(36,:),'th2r','Resolved theta variance','K^2','tt')
call ncinfo(ncname(37,:),'qt2r','Resolved total water variance','(kg/kg)^2','tt')
call ncinfo(ncname(38,:),'ql2r','Resolved liquid water variance','(kg/kg)^2','tt')
call ncinfo(ncname(39,:),'cs','Smagorinsky constant','-','tt')
call ncinfo(ncname(40,:),'cfrac','Cloud fraction','-','tt')
call ncinfo(ncname(41,:),'hur','Relative humidity','%','tt')
call ncinfo(ncname(42,:),'hus','Specific humidity','kg/kg','tt')
call ncinfo(ncname(43,:),'ta', 'Temperature','K','tt')
call ncinfo(ncname(44,:),'clw', 'Specific cloud liquid water content','kg/kg','tt')
call ncinfo(ncname(45,:),'cli', 'Specific cloud ice content','kg/kg','tt')
call ncinfo(ncname(46,:),'plw', 'Specific precipitation liquid water content','kg/kg','tt')
call ncinfo(ncname(47,:),'pli', 'Specific precipitation ice content','kg/kg','tt')
do n=1,nsv
write (csvname(1:3),'(i3.3)') n
call ncinfo(ncname(47+7*(n-1)+1,:),'sv'//csvname,'Scalar '//csvname//' specific mixing ratio','(kg/kg)','tt')
call ncinfo(ncname(47+7*(n-1)+2,:),'svp'//csvname,'Scalar '//csvname//' tendency','(kg/kg/s)','tt')
call ncinfo(ncname(47+7*(n-1)+3,:),'svpt'//csvname,'Scalar '//csvname//' turbulence tendency','(kg/kg/s)','tt')
call ncinfo(ncname(47+7*(n-1)+4,:),'sv'//csvname//'2r','Resolved scalar '//csvname//' variance','(kg/kg)^2','tt')
call ncinfo(ncname(47+7*(n-1)+5,:),'wsv'//csvname//'s','SFS scalar '//csvname//' flux','kg/kg m/s','mt')
call ncinfo(ncname(47+7*(n-1)+6,:),'wsv'//csvname//'r','Resolved scalar '//csvname//' flux','kg/kg m/s','mt')
call ncinfo(ncname(47+7*(n-1)+7,:),'wsv'//csvname//'t','Total scalar '//csvname//' flux','kg/kg m/s','mt')
end do
if (isurf==1) then
call open_nc(fname, ncid,nrec,n3=kmax,ns=ksoilmax)
else
call open_nc(fname, ncid,nrec,n3=kmax)
endif
if (nrec == 0) then
call define_nc( ncid, 1, tncname)
call writestat_dims_nc(ncid)
end if
call define_nc( ncid, NVar, ncname)
end if
end if
end subroutine initgenstat
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine genstat
use modglobal, only : rk3step,timee,dt_lim
implicit none
if (.not. lstat) return
if (rk3step/=3) return
if(timee<tnext .and. timee<tnextwrite) then
dt_lim = minval((/dt_lim,tnext-timee,tnextwrite-timee/))
return
end if
if (timee>=tnext) then
tnext = tnext+idtav
call do_genstat
end if
if (timee>=tnextwrite) then
tnextwrite = tnextwrite+itimeav
call writestat
end if
dt_lim = minval((/dt_lim,tnext-timee,tnextwrite-timee/))
end subroutine genstat
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine do_genstat
use modfields, only : u0,v0,w0,thl0,qt0,qt0h,e120, &
ql0,ql0h,thl0h,thv0h,sv0,exnf,exnh,tmp0,presf
use modsurfdata,only: thls,qts,svs,ustar,thlflux,qtflux,svflux
use modsubgriddata,only : ekm, ekh, csz
use modglobal, only : i1,ih,j1,jh,k1,kmax,nsv,dzf,dzh,rlv,rv,rd,cp, &
ijtot,cu,cv,iadv_sv,iadv_kappa,eps1,dxi,dyi,tup,tdn
use modmpi, only : comm3d,mpi_sum,mpierr,slabsum,D_MPI_ALLREDUCE
use advec_kappa, only : halflev_kappa
use modmicrodata, only: tuprsg, tdnrsg, imicro, imicro_sice, imicro_sice2
use modthermodynamics, only: qsat_tab
implicit none
real cthl,cqt,den
real,allocatable, dimension(:) :: &
qlhavl , & ! slab averaged ql_0 at half level &
u2avl , &
v2avl , &
w2avl , &
w3avl , &
w2subavl , &
qt2avl , &
thl2avl , &
thv2avl , &
ql2avl , &
th2avl
real,allocatable, dimension(:,:) :: &
wsvsubl,& ! slab averaged sub w-sv(n) flux &
wsvresl,& ! slab averaged res w-sv(n) flux &
sv2avl,&
sv2av
real,allocatable, dimension(:) :: wqlsubl
real,allocatable, dimension(:) :: wqlresl
real,allocatable, dimension(:):: wthlsubl
real,allocatable, dimension(:):: wthlresl
real,allocatable, dimension(:):: wqtsubl
real,allocatable, dimension(:):: wqtresl
real,allocatable, dimension(:):: wthvsubl
real,allocatable, dimension(:) ::wthvresl
real,allocatable, dimension(:):: cfracavl ! cloudfraction at full level
real,allocatable, dimension(:):: huravl
real,allocatable, dimension(:):: qlptavl ! slab averaged turbulence tendency of q_liq
real,allocatable, dimension(:):: uwsubl
real,allocatable, dimension(:):: vwsubl
real,allocatable, dimension(:):: uwresl
real,allocatable, dimension(:):: vwresl
real,allocatable, dimension(:):: qlhav
real,allocatable, dimension(:):: u2av , &
v2av , &
w2av , &
w3av , &
w2subav , &
qt2av , &
thl2av , &
thv2av , &
th2av , &
ql2av
real(field_r),allocatable, dimension(:,:,:):: thv0
real(field_r),allocatable, dimension(:):: thvmav
real(field_r),allocatable, dimension(:,:,:):: sv0h
integer i, j, k, n, km
real tsurf, qsat_, c1, c2
real qs0h, t0h, ekhalf, euhalf, evhalf
real wthls, wthlr, wqts, wqtr, wqls, wqlr, wthvs, wthvr
real uws,vws,uwr,vwr
real upcu, vpcv
real qls
real(field_r) :: ilratio
allocate( &
qlhavl (k1), & ! slab averaged ql_0 at half level &
wsvsubl(k1,nsv),& ! slab averaged sub w-sv(n) flux &
wsvresl(k1,nsv),& ! slab averaged res w-sv(n) flux &
u2avl (k1), &
v2avl (k1), &
w2avl (k1), &
w3avl (k1), &
w2subavl (k1), &
qt2avl (k1), &
thl2avl (k1), &
thv2avl (k1), &
th2avl (k1))
allocate( &
ql2avl (k1), &
sv2avl (k1,nsv))
allocate( wqlsubl (k1))
allocate( wqlresl (k1))
allocate( wthlsubl (k1))
allocate( wthlresl (k1))
allocate( wqtsubl (k1))
allocate( wqtresl (k1))
allocate( wthvsubl (k1))
allocate( wthvresl (k1))
allocate( cfracavl(k1)) ! slab averaged cloud fraction
allocate( huravl(k1))
allocate( qlptavl(k1)) ! slab averaged turbulence tendency of q_liq
allocate( uwsubl(k1))
allocate( vwsubl(k1))
allocate( uwresl(k1))
allocate( vwresl(k1))
allocate( qlhav(k1))
allocate( u2av (k1), &
v2av (k1), &
w2av (k1), &
w3av (k1), &
w2subav (k1), &
qt2av (k1), &
thl2av (k1), &
thv2av (k1), &
th2av (k1), &
ql2av (k1), &
sv2av (k1,nsv))
allocate(thv0(2-ih:i1+ih,2-jh:j1+jh,k1))
allocate(thvmav(k1))
allocate(sv0h(2-ih:i1+ih,2-jh:j1+jh,k1))
!-----------------------------------------------------------------------
! 1. INITIALISE LOCAL CONSTANTS
! -- --------------------------
! --------------------------------------------------------
! 3.0 RESET ARRAYS FOR SLAB AVERAGES
! --- ------------------------------
! --------------------------------------------------------
qlhavl = 0.0
cfracavl = 0.0
huravl = 0.0
qlptavl = 0.0
wqlsubl = 0.0
wqlresl = 0.0
wqltot = 0.0
wthlsubl = 0.0
wthlresl = 0.0
wthltot = 0.0
wqtsubl = 0.0
wqtresl = 0.0
wqttot = 0.0
wthvsubl = 0.0
wthvresl = 0.0
wthvtot = 0.0
wsvsubl = 0.
wsvresl = 0.
sv2avl = 0.
uwresl = 0.
vwresl = 0.
uwtot = 0.
uwsubl = 0.
vwsubl = 0.
vwtot = 0.
u2avl = 0.0
v2avl = 0.0
w2avl = 0.0
w3avl = 0.0
w2subavl = 0.0
qt2avl = 0.0
thl2avl = 0.0
thv2avl = 0.0
th2avl = 0.0
ql2avl = 0.0
thvmav = 0.0
sv2av = 0.0
umav = 0.0
vmav = 0.0
thlmav = 0.0
thmav = 0.0
qtmav = 0.0
qlmav = 0.0
cfracav= 0.0
hurav = 0.0
clwav = 0.0
cliav = 0.0
plwav = 0.0
pliav = 0.0
taav = 0.0
svmav = 0.
cszav = 0.
do k=1,k1
do j=2,j1
do i=2,i1
thv0(i,j,k) = (thl0(i,j,k)+rlv*ql0(i,j,k)/(cp*exnf(k))) &
*(1+(rv/rd-1)*qt0(i,j,k)-rv/rd*ql0(i,j,k))
huravl(k) = huravl(k) + 100 * (qt0(i,j,k) - ql0(i,j,k)) / qsat_tab(tmp0(i,j,k), presf(k))
enddo
enddo
enddo
do k=1,k1
do j=2,j1
do i=2,i1
ilratio = max(0._field_r,min(1._field_r,(tmp0(i,j,k)-tdn)/(tup-tdn)))
clwav(k) = clwav(k) + ql0(i,j,k) * ilratio
cliav(k) = cliav(k) + ql0(i,j,k) * (1-ilratio)
end do
end do
end do
if (imicro == imicro_sice .or. imicro == imicro_sice2) then
do k=1,k1
do j=2,j1
do i=2,i1
ilratio = max(0._field_r,min(1._field_r,(tmp0(i,j,k)-tdnrsg)/(tuprsg-tdnrsg)))
plwav(k) = plwav(k) + ql0(i,j,k) * ilratio
pliav(k) = pliav(k) + ql0(i,j,k) * (1-ilratio)
end do
end do
end do
end if
do k=1,k1
cfracavl(k) = cfracavl(k)+count(ql0(2:i1,2:j1,k)>0)
end do
call D_MPI_ALLREDUCE(clwav,k1,MPI_SUM,comm3d,mpierr)
call D_MPI_ALLREDUCE(cliav,k1,MPI_SUM,comm3d,mpierr)
call D_MPI_ALLREDUCE(plwav,k1,MPI_SUM,comm3d,mpierr)
call D_MPI_ALLREDUCE(pliav,k1,MPI_SUM,comm3d,mpierr)
clwav = clwav / ijtot
cliav = cliav / ijtot
plwav = plwav / ijtot
pliav = pliav / ijtot
call D_MPI_ALLREDUCE(cfracavl,cfracav,k1,MPI_SUM,comm3d,mpierr)
call D_MPI_ALLREDUCE(huravl,hurav,k1,MPI_SUM,comm3d,mpierr)
call slabsum(umav ,1,k1,u0 ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(vmav ,1,k1,v0 ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(thlmav,1,k1,thl0,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(qtmav ,1,k1,qt0 ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(qlmav ,1,k1,ql0 ,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(thvmav,1,k1,thv0,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
call slabsum(taav ,1,k1,tmp0,2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
umav = umav /ijtot + cu
vmav = vmav /ijtot + cv
thlmav = thlmav/ijtot
qtmav = qtmav /ijtot
qlmav = qlmav /ijtot
cfracav = cfracav / ijtot
hurav = hurav / ijtot
taav = taav / ijtot
thmav = thlmav + (rlv/cp)*qlmav/exnf
thvmav = thvmav/ijtot
cszav = csz
!
do n=1,nsv
call slabsum(svmav(1:1,n),1,k1,sv0(:,:,:,n),2-ih,i1+ih,2-jh,j1+jh,1,k1,2,i1,2,j1,1,k1)
enddo
svmav = svmav/ijtot
!------------------------------------------------------------------
! 4 CALCULATE SLAB AVERAGED OF FLUXES AND SEVERAL MOMENTS
! -------------------------------------------------------------
! 4.1 special treatment for lowest level
! -------------------------------------------------
qls = 0.0 ! hj: no liquid water at the surface
tsurf = thls*exnh(1)+(rlv/cp)*qls
qsat_ = qts - qls
if (qls< eps1) then ! TH: Should always be true
c1 = 1.+(rv/rd-1)*qts
c2 = (rv/rd-1)
else
c1 = (1.-qts+rv/rd*qsat_*(1.+rlv/(rv*tsurf))) &
/(1.+rlv/(rv*tsurf)*rlv/(cp*tsurf)*qsat_)
c2 = c1*rlv/(tsurf*cp)-1.
end if
den = 1. + (rlv**2)*qsat_/(rv*cp*(tsurf**2))
cthl = (exnh(1)*cp/rlv)*((1-den)/den)
cqt = 1./den
do j=2,j1
do i=2,i1
qlhavl(1) = qlhavl(1) + ql0h(i,j,1)
! thv(1) = thlm(i,j,1) * (1.+(rv/rd-1)*qtm(i,j,1))
wthlsubl(1) = wthlsubl(1) + thlflux(i,j)
wqtsubl(1) = wqtsubl(1) + qtflux (i,j)
wqlsubl(1) = 0
wthvsubl(1) = wthvsubl(1) + ( c1*thlflux(i,j)+c2*thls*qtflux(i,j) ) !hj: thv0 replaced by thls
!Momentum flux
if (abs(u0(i,j,1)+cu)<eps1) then
upcu = sign(eps1,u0(i,j,1)+cu)
else
upcu = u0(i,j,1)+cu
end if
uwsubl(1) = uwsubl(1) - ( 0.5*( ustar(i,j)+ustar(i-1,j) ) )**2 * &
upcu/sqrt(upcu**2 + &
((v0(i,j,1)+v0(i-1,j,1)+v0(i,j+1,1)+v0(i-1,j+1,1))/4.+cv)**2)
if (abs(v0(i,j,1)+cv)<eps1) then
vpcv = sign(eps1,v0(i,j,1)+cv)
else
vpcv = v0(i,j,1)+cv
end if
vwsubl(1) = vwsubl(1) - ( 0.5*( ustar(i,j)+ustar(i,j-1) ) )**2 * &
vpcv/sqrt(vpcv**2 + &
((u0(i,j,1)+u0(i+1,j,1)+u0(i,j-1,1)+u0(i+1,j-1,1))/4.+cu)**2)
!Higher order moments
u2avl (1) = u2avl (1) + (u0 (i,j,1)+cu - umav(1))**2
v2avl (1) = v2avl (1) + (v0 (i,j,1)+cv - vmav(1))**2
w2avl (1) = w2avl (1) + (w0 (i,j,1)**2)
w3avl (1) = w3avl (1) + (w0 (i,j,1)**3)
w2subavl (1) = w2subavl (1) + (e120(i,j,1)**2)
qt2avl (1) = qt2avl (1) + (qt0 (i,j,1) - qtmav (1))**2
thl2avl (1) = thl2avl (1) + (thl0(i,j,1) - thlmav(1))**2
thv2avl (1) = thv2avl (1) + (thv0(i,j,1) - thvmav(1))**2
th2avl (1) = th2avl (1) + (thl0(i,j,1) - thmav (1))**2
ql2avl (1) = ql2avl (1) + (ql0(i,j,1) - qlmav (1))**2
! qs2avl (1) = qs2avl (1) + qs0**2
! qsavl (1) = qsavl (1) + qs0
! rhavl (1) = rhavl (1) + qtm (i,j,1)/qs0
! ravl (1) = ravl (1) + ((qt0(i,j,1) - qs0))
! r2avl (1) = r2avl (1) + ((qt0(i,j,1) - qs0)**2)
! r3avl (1) = r3avl (1) + ((qt0(i,j,1) - qs0)**3)
do n=1,nsv
wsvsubl(1,n) = wsvsubl(1,n) + svflux(i,j,n)
sv2avl(1,n) = sv2avl(1,n) + (sv0(i,j,1,n)-svmav(1,n))**2
end do
end do
end do
! --------------------------
! 4.2 higher levels
! --------------------------
do j=2,j1
do i=2,i1
! --------------------------------------------------------
! Calculate half level fields for thl and qt consistent
! with the used advection scheme ( kappa or cent. diff.)
! ---------------
! --------------------------------------------------------
do k=2,kmax
km = k-1
! ------------------------------------------------------
! calculate ql and thv at time t0 at full and half level
! ----------------------------------------------------
qlhavl(k) = qlhavl(k) + ql0h(i,j,k)
! -----------------------------------------------------------
! calculate prefactors for subgrid wthv and wql fluxes
! at half levels
! -----------------------------------------------------------
qs0h = (qt0h(i,j,k) - ql0h(i,j,k))
t0h = exnh(k)*thl0h(i,j,k) + (rlv/cp)*ql0h(i,j,k)
den = 1. + (rlv**2)*qs0h/(rv*cp*(t0h**2))
cthl = (exnh(k)*cp/rlv)*((1-den)/den)
cqt = (1./den)
if (ql0h(i,j,k)>0) then
c1 = (1.-qt0h(i,j,k)+rv/rd*qs0h &
* (1.+rd/rv*rlv/(rd*t0h)))/den
c2 = c1*rlv/(t0h*cp)-1.
else
c1 = 1. + (rv/rd-1)*qt0h(i,j,k)
c2 = (rv/rd-1)
end if
! -----------------------------------------------------------
! calculate resolved and subgrid fluxes at half levels
! -----------------------------------------------------------
ekhalf = (ekh(i,j,k)*dzf(km)+ekh(i,j,km)*dzf(k))/(2*dzh(k))
euhalf = ( dzf(km) * ( ekm(i,j,k) + ekm(i-1,j,k) ) + &
dzf(k) * ( ekm(i,j,km) + ekm(i-1,j,km) ) ) / &
( 4. * dzh(k) )
evhalf = ( dzf(km) * ( ekm(i,j,k) + ekm(i,j-1,k) ) + &
dzf(k) * ( ekm(i,j,km) + ekm(i,j-1,km) ) ) / &
( 4. * dzh(k) )
wthls = -ekhalf*(thl0(i,j,k)-thl0(i,j,km))/dzh(k)
wthlr = w0(i,j,k)*thl0h(i,j,k)
wqts = -ekhalf*(qt0(i,j,k)-qt0(i,j,km))/dzh(k)
wqtr = w0(i,j,k)*qt0h(i,j,k)
wqls = cthl*wthls+ cqt*wqts
wqlr = w0(i,j,k)*ql0h(i,j,k)
wthvs = c1*wthls + c2*thl0h(i,j,k)*wqts
wthvr = w0(i,j,k)*thv0h(i,j,k)
uwr = (w0(i,j,k)+w0(i-1,j,k)) &
*((u0(i,j,k-1)+cu)*dzf(k)+(u0(i,j,k)+cu)*dzf(k-1))/(4*dzh(k))
vwr = (w0(i,j,k)+w0(i,j-1,k)) &
*((v0(i,j,k-1)+cv)*dzf(k)+(v0(i,j,k)+cv)*dzf(k-1))/(4*dzh(k))
uws = -euhalf &
*((u0(i,j,k)-u0(i,j,k-1))/dzh(k)+(w0(i,j,k)-w0(i-1,j,k))*dxi)
vws = -evhalf &
*((v0(i,j,k)-v0(i,j,k-1))/dzh(k)+(w0(i,j,k)-w0(i,j-1,k))*dyi)
if (ql0h(i,j,k)>0) then
wqlsubl(k) = wqlsubl(k) + wqls
end if
wqlresl(k) = wqlresl(k) + wqlr
wthlsubl(k) = wthlsubl(k) + wthls
wthlresl(k) = wthlresl(k) + wthlr
wthvsubl(k) = wthvsubl(k) + wthvs
wthvresl(k) = wthvresl(k) + wthvr
wqtsubl(k) = wqtsubl(k) + wqts
wqtresl(k) = wqtresl(k) + wqtr
uwresl(k) = uwresl(k) + uwr
vwresl(k) = vwresl(k) + vwr
uwsubl(k) = uwsubl(k) + uws
vwsubl(k) = vwsubl(k) + vws
! -----------------------------------------------------------
! calculate various moments
! -----------------------------------------------------------
u2avl (k) = u2avl (k) + (u0 (i,j,k)+cu - umav(k))**2
v2avl (k) = v2avl (k) + (v0 (i,j,k)+cv - vmav(k))**2
w2avl (k) = w2avl (k) + (w0 (i,j,k)**2)
w3avl (k) = w3avl (k) + (w0 (i,j,k)**3)
w2subavl (k) = w2subavl (k) + (e120(i,j,k)**2)
qt2avl (k) = qt2avl (k) + (qt0 (i,j,k) - qtmav (k))**2
thl2avl (k) = thl2avl (k) + (thl0(i,j,k) - thlmav(k))**2
thv2avl (k) = thv2avl (k) + (thv0(i,j,k) - thvmav(k))**2
th2avl (k) = th2avl (k) + (thl0(i,j,k) - thmav (k))**2 !thlm, no thm !?!
ql2avl (k) = ql2avl (k) + (ql0(i,j,k) - qlmav (k))**2
! qs2avl (k) = qs2avl (k) + qs0**2
! qsavl (k) = qsavl (k) + qs0
! rhavl (k) = rhavl (k) + qtm (i,j,k)/qs0
! ravl (k) = ravl (k) + (qt0(i,j,k) - qs0)
! r2avl (k) = r2avl (k) + (qt0(i,j,k) - qs0)**2
! r3avl (k) = r3avl (k) + (qt0(i,j,k) - qs0)**3
end do
end do
end do
! -------------------
do n=1,nsv
do k=2,kmax
do j=2,j1
do i=2,i1
sv2avl(k,n) = sv2avl(k,n) + (sv0(i,j,k,n)-svmav(k,n))**2
end do
end do
end do
end do
do n=1,nsv
if (iadv_sv(n)==iadv_kappa) then
call halflev_kappa(sv0(2-ih:i1+ih,2-jh:j1+jh,1:k1,n),sv0h)
else
do k=2,k1
do j=2,j1
do i=2,i1
sv0h(i,j,k) = (sv0(i,j,k,n)*dzf(k-1)+sv0(i,j,k-1,n)*dzf(k))/(2*dzh(k))
enddo
enddo
enddo
sv0h(2:i1,2:j1,1) = svs(n)
end if
do k=2,kmax
do j=2,j1
do i=2,i1
wsvresl(k,n) = wsvresl(k,n) + w0(i,j,k)*sv0h(i,j,k)
end do
end do
end do
do k=2,kmax
km = k-1
do j=2,j1
do i=2,i1
ekhalf = (ekh(i,j,k)*dzf(km)+ekh(i,j,km)*dzf(k))/(2*dzh(k))
wsvsubl(k,n)= wsvsubl(k,n)-ekhalf*(sv0(i,j,k,n)-sv0(i,j,km,n)) &
/dzh(k)
end do