forked from cmkaul/SCAMPy
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathTurbulence_PrognosticTKE.pyx
executable file
·2360 lines (2036 loc) · 128 KB
/
Turbulence_PrognosticTKE.pyx
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
#!python
#cython: boundscheck=False
#cython: wraparound=False
#cython: initializedcheck=True
#cython: cdivision=False
import numpy as np
include "parameters.pxi"
import cython
import sys
from Grid cimport Grid
cimport EDMF_Updrafts
cimport EDMF_Environment
cimport EDMF_Rain
from Variables cimport VariablePrognostic, VariableDiagnostic, GridMeanVariables
from Surface cimport SurfaceBase
from Cases cimport CasesBase
from ReferenceState cimport ReferenceState
from TimeStepping cimport TimeStepping
from NetCDFIO cimport NetCDFIO_Stats
from thermodynamic_functions cimport *
from microphysics_functions cimport *
from turbulence_functions cimport *
from utility_functions cimport *
from libc.math cimport fmax, sqrt, exp, pow, cbrt, fmin, fabs
from cpython.mem cimport PyMem_Malloc, PyMem_Realloc, PyMem_Free
cdef class EDMF_PrognosticTKE(ParameterizationBase):
# Initialize the class
def __init__(self, namelist, paramlist, Grid Gr, ReferenceState Ref):
# Initialize the base parameterization class
ParameterizationBase.__init__(self, paramlist, Gr, Ref)
# Set the number of updrafts (1)
try:
self.n_updrafts = namelist['turbulence']['EDMF_PrognosticTKE']['updraft_number']
except:
self.n_updrafts = 1
print('Turbulence--EDMF_PrognosticTKE: defaulting to single updraft')
# TODO - steady updrafts have not been tested!
try:
self.use_steady_updrafts = namelist['turbulence']['EDMF_PrognosticTKE']['use_steady_updrafts']
except:
self.use_steady_updrafts = False
try:
self.calc_tke = namelist['turbulence']['EDMF_PrognosticTKE']['calculate_tke']
except:
self.calc_tke = True
try:
self.use_const_plume_spacing = namelist['turbulence']['EDMF_PrognosticTKE']['use_constant_plume_spacing']
except:
self.use_const_plume_spacing = False
try:
self.calc_scalar_var = namelist['turbulence']['EDMF_PrognosticTKE']['calc_scalar_var']
except:
self.calc_scalar_var = False
if (self.calc_scalar_var==True and self.calc_tke==False):
sys.exit('Turbulence--EDMF_PrognosticTKE: >>calculate_tke<< must be set to True when >>calc_scalar_var<< is True (to calculate the mixing length for the variance and covariance calculations')
try:
if str(namelist['turbulence']['EDMF_PrognosticTKE']['entrainment']) == 'inverse_z':
self.entr_detr_fp = entr_detr_inverse_z
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['entrainment']) == 'dry':
self.entr_detr_fp = entr_detr_dry
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['entrainment']) == 'inverse_w':
self.entr_detr_fp = entr_detr_inverse_w
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['entrainment']) == 'b_w2':
self.entr_detr_fp = entr_detr_b_w2
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['entrainment']) == 'entr_detr_tke':
self.entr_detr_fp = entr_detr_tke
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['entrainment']) == 'suselj':
self.entr_detr_fp = entr_detr_suselj
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['entrainment']) == 'buoyancy_sorting':
self.entr_detr_fp = entr_detr_buoyancy_sorting
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['entrainment']) == 'moisture_deficit':
self.entr_detr_fp = entr_detr_env_moisture_deficit
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['entrainment']) == 'none':
self.entr_detr_fp = entr_detr_none
else:
print('Turbulence--EDMF_PrognosticTKE: Entrainment rate namelist option is not recognized')
except:
self.entr_detr_fp = entr_detr_b_w2
print('Turbulence--EDMF_PrognosticTKE: defaulting to cloudy entrainment formulation')
if(self.calc_tke == False and 'tke' in str(namelist['turbulence']['EDMF_PrognosticTKE']['entrainment'])):
sys.exit('Turbulence--EDMF_PrognosticTKE: >>calc_tke<< must be set to True when entrainment is using tke')
try:
if str(namelist['turbulence']['EDMF_PrognosticTKE']['pressure_closure_buoy']) == 'tan18':
self.pressure_func_buoy = pressure_tan18_buoy
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['pressure_closure_buoy']) == 'normalmode':
self.pressure_func_buoy = pressure_normalmode_buoy
else:
print('Turbulence--EDMF_PrognosticTKE: pressure closure in namelist option is not recognized')
except:
self.pressure_func_buoy = pressure_tan18_buoy
print('Turbulence--EDMF_PrognosticTKE: defaulting to pressure closure Tan2018')
self.drag_sign = False
try:
if str(namelist['turbulence']['EDMF_PrognosticTKE']['pressure_closure_drag']) == 'tan18':
self.pressure_func_drag = pressure_tan18_drag
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['pressure_closure_drag']) == 'normalmode':
self.pressure_func_drag = pressure_normalmode_drag
elif str(namelist['turbulence']['EDMF_PrognosticTKE']['pressure_closure_drag']) == 'normalmode_signdf':
self.pressure_func_drag = pressure_normalmode_drag
self.drag_sign = True
else:
print('Turbulence--EDMF_PrognosticTKE: pressure closure in namelist option is not recognized')
except:
self.pressure_func_drag = pressure_tan18_drag
print('Turbulence--EDMF_PrognosticTKE: defaulting to pressure closure Tan2018')
try:
self.asp_label = str(namelist['turbulence']['EDMF_PrognosticTKE']['pressure_closure_asp_label'])
except:
self.asp_label = 'const'
print('Turbulence--EDMF_PrognosticTKE: H/2R defaulting to constant')
try:
self.similarity_diffusivity = namelist['turbulence']['EDMF_PrognosticTKE']['use_similarity_diffusivity']
except:
self.similarity_diffusivity = False
print('Turbulence--EDMF_PrognosticTKE: defaulting to TKE-based eddy diffusivity')
if(self.similarity_diffusivity == False and self.calc_tke ==False):
sys.exit('Turbulence--EDMF_PrognosticTKE: either >>use_similarity_diffusivity<< or >>calc_tke<< flag is needed to get the eddy diffusivities')
if(self.similarity_diffusivity == True and self.calc_tke == True):
print("TKE will be calculated but not used for eddy diffusivity calculation")
try:
self.extrapolate_buoyancy = namelist['turbulence']['EDMF_PrognosticTKE']['extrapolate_buoyancy']
except:
self.extrapolate_buoyancy = True
print('Turbulence--EDMF_PrognosticTKE: defaulting to extrapolation of updraft buoyancy along a pseudoadiabat')
try:
self.mixing_scheme = str(namelist['turbulence']['EDMF_PrognosticTKE']['mixing_length'])
except:
self.mixing_scheme = 'default'
print 'Using (Tan et al, 2018) default'
# Get values from paramlist
# set defaults at some point?
self.surface_area = paramlist['turbulence']['EDMF_PrognosticTKE']['surface_area']
self.max_area = paramlist['turbulence']['EDMF_PrognosticTKE']['max_area']
self.entrainment_factor = paramlist['turbulence']['EDMF_PrognosticTKE']['entrainment_factor']
self.updraft_mixing_frac = paramlist['turbulence']['EDMF_PrognosticTKE']['updraft_mixing_frac']
self.entrainment_sigma = paramlist['turbulence']['EDMF_PrognosticTKE']['entrainment_sigma']
self.entrainment_smin_tke_coeff = paramlist['turbulence']['EDMF_PrognosticTKE']['entrainment_smin_tke_coeff']
self.entrainment_ed_mf_sigma = paramlist['turbulence']['EDMF_PrognosticTKE']['entrainment_smin_tke_coeff']
self.entrainment_scale = paramlist['turbulence']['EDMF_PrognosticTKE']['entrainment_scale']
self.constant_plume_spacing = paramlist['turbulence']['EDMF_PrognosticTKE']['constant_plume_spacing']
self.detrainment_factor = paramlist['turbulence']['EDMF_PrognosticTKE']['detrainment_factor']
self.sorting_power = paramlist['turbulence']['EDMF_PrognosticTKE']['sorting_power']
self.turbulent_entrainment_factor = paramlist['turbulence']['EDMF_PrognosticTKE']['turbulent_entrainment_factor']
self.pressure_buoy_coeff = paramlist['turbulence']['EDMF_PrognosticTKE']['pressure_buoy_coeff']
self.aspect_ratio = paramlist['turbulence']['EDMF_PrognosticTKE']['aspect_ratio']
if str(namelist['turbulence']['EDMF_PrognosticTKE']['pressure_closure_buoy']) == 'normalmode':
try:
self.pressure_normalmode_buoy_coeff1 = paramlist['turbulence']['EDMF_PrognosticTKE']['pressure_normalmode_buoy_coeff1']
self.pressure_normalmode_buoy_coeff2 = paramlist['turbulence']['EDMF_PrognosticTKE']['pressure_normalmode_buoy_coeff2']
except:
self.pressure_normalmode_buoy_coeff1 = self.pressure_buoy_coeff
self.pressure_normalmode_buoy_coeff2 = 0.0
print 'Using (Tan et al, 2018) parameters as default for Normal Mode pressure formula buoyancy term'
if str(namelist['turbulence']['EDMF_PrognosticTKE']['pressure_closure_drag']) == 'normalmode':
try:
self.pressure_normalmode_adv_coeff = paramlist['turbulence']['EDMF_PrognosticTKE']['pressure_normalmode_adv_coeff']
self.pressure_normalmode_drag_coeff = paramlist['turbulence']['EDMF_PrognosticTKE']['pressure_normalmode_drag_coeff']
except:
self.pressure_normalmode_adv_coeff = 0.0
self.pressure_normalmode_drag_coeff = 1.0
print 'Using (Tan et al, 2018) parameters as default for Normal Mode pressure formula drag term'
# "Legacy" coefficients used by the steady updraft routine
self.vel_buoy_coeff = 1.0-self.pressure_buoy_coeff
if self.calc_tke == True:
self.tke_ed_coeff = paramlist['turbulence']['EDMF_PrognosticTKE']['tke_ed_coeff']
self.tke_diss_coeff = paramlist['turbulence']['EDMF_PrognosticTKE']['tke_diss_coeff']
self.static_stab_coeff = paramlist['turbulence']['EDMF_PrognosticTKE']['static_stab_coeff']
# Latent heat stability effect
self.lambda_stab = paramlist['turbulence']['EDMF_PrognosticTKE']['lambda_stab']
# Need to code up as paramlist option?
self.minimum_area = 1e-5
# Create the class for rain
self.Rain = EDMF_Rain.RainVariables(namelist, Gr)
if self.use_steady_updrafts == True and self.Rain.rain_model != "None":
sys.exit('PrognosticTKE: rain model is available for prognostic updrafts only')
self.RainPhysics = EDMF_Rain.RainPhysics(Gr, Ref)
# Create the updraft variable class (major diagnostic and prognostic variables)
self.UpdVar = EDMF_Updrafts.UpdraftVariables(self.n_updrafts, namelist,paramlist, Gr)
# Create the class for updraft thermodynamics
self.UpdThermo = EDMF_Updrafts.UpdraftThermodynamics(self.n_updrafts, Gr, Ref, self.UpdVar, self.Rain)
# Create the environment variable class (major diagnostic and prognostic variables)
self.EnvVar = EDMF_Environment.EnvironmentVariables(namelist,Gr)
# Create the class for environment thermodynamics
self.EnvThermo = EDMF_Environment.EnvironmentThermodynamics(namelist, Gr, Ref, self.EnvVar, self.Rain)
# Entrainment rates
self.entr_sc = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
#self.press = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
# Detrainment rates
self.detr_sc = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.sorting_function = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.b_mix = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
# turbulent entrainment
self.frac_turb_entr = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.frac_turb_entr_full = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.turb_entr_W = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.turb_entr_H = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.turb_entr_QT = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
# Pressure term in updraft vertical momentum equation
self.nh_pressure = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.nh_pressure_b = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.nh_pressure_adv = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.nh_pressure_drag = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.asp_ratio = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.b_coeff = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
# Mass flux
self.m = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double, order='c')
# mixing length
self.mixing_length = np.zeros((Gr.nzg,),dtype=np.double, order='c')
self.horizontal_KM = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
self.horizontal_KH = np.zeros((self.n_updrafts, Gr.nzg),dtype=np.double,order='c')
# diagnosed tke budget terms
self.tke_transport = np.zeros((Gr.nzg,),dtype=np.double, order='c')
self.tke_advection = np.zeros((Gr.nzg,),dtype=np.double, order='c')
# Near-surface BC of updraft area fraction
self.area_surface_bc= np.zeros((self.n_updrafts,),dtype=np.double, order='c')
self.w_surface_bc= np.zeros((self.n_updrafts,),dtype=np.double, order='c')
self.h_surface_bc= np.zeros((self.n_updrafts,),dtype=np.double, order='c')
self.qt_surface_bc= np.zeros((self.n_updrafts,),dtype=np.double, order='c')
self.pressure_plume_spacing = np.zeros((self.n_updrafts,),dtype=np.double,order='c')
# Mass flux tendencies of mean scalars (for output)
self.massflux_tendency_h = np.zeros((Gr.nzg,),dtype=np.double,order='c')
self.massflux_tendency_qt = np.zeros((Gr.nzg,),dtype=np.double,order='c')
# (Eddy) diffusive tendencies of mean scalars (for output)
self.diffusive_tendency_h = np.zeros((Gr.nzg,),dtype=np.double,order='c')
self.diffusive_tendency_qt = np.zeros((Gr.nzg,),dtype=np.double,order='c')
# Vertical fluxes for output
self.massflux_h = np.zeros((Gr.nzg,),dtype=np.double,order='c')
self.massflux_qt = np.zeros((Gr.nzg,),dtype=np.double,order='c')
self.diffusive_flux_h = np.zeros((Gr.nzg,),dtype=np.double,order='c')
self.diffusive_flux_qt = np.zeros((Gr.nzg,),dtype=np.double,order='c')
self.diffusive_flux_u = np.zeros((Gr.nzg,),dtype=np.double,order='c')
self.diffusive_flux_v = np.zeros((Gr.nzg,),dtype=np.double,order='c')
if self.calc_tke:
self.massflux_tke = np.zeros((Gr.nzg,),dtype=np.double,order='c')
# Added by Ignacio : Length scheme in use (mls), and smooth min effect (ml_ratio)
# Variable Prandtl number initialized as neutral value.
self.prandtl_nvec = np.multiply( self.prandtl_number, np.ones((Gr.nzg,),dtype=np.double, order='c'))
self.mls = np.zeros((Gr.nzg,),dtype=np.double, order='c')
self.ml_ratio = np.zeros((Gr.nzg,),dtype=np.double, order='c')
self.l_entdet = np.zeros((Gr.nzg,),dtype=np.double, order='c')
self.b = np.zeros((Gr.nzg,),dtype=np.double, order='c')
return
cpdef initialize(self, CasesBase Case, GridMeanVariables GMV, ReferenceState Ref):
if Case.casename == 'DryBubble':
print 'updraft initialized for Dry Bubble'
self.UpdVar.initialize_DryBubble(GMV, Ref)
else:
self.UpdVar.initialize(GMV)
return
# Initialize the IO pertaining to this class
cpdef initialize_io(self, NetCDFIO_Stats Stats):
self.UpdVar.initialize_io(Stats)
self.EnvVar.initialize_io(Stats)
self.Rain.initialize_io(Stats)
Stats.add_profile('eddy_viscosity')
Stats.add_profile('eddy_diffusivity')
Stats.add_profile('entrainment_sc')
Stats.add_profile('detrainment_sc')
Stats.add_profile('nh_pressure')
Stats.add_profile('nh_pressure_adv')
Stats.add_profile('nh_pressure_drag')
Stats.add_profile('nh_pressure_b')
Stats.add_profile('asp_ratio')
Stats.add_profile('b_coeff')
Stats.add_profile('horizontal_KM')
Stats.add_profile('horizontal_KH')
Stats.add_profile('sorting_function')
Stats.add_profile('b_mix')
Stats.add_ts('rd')
Stats.add_profile('turbulent_entrainment')
Stats.add_profile('turbulent_entrainment_full')
Stats.add_profile('turbulent_entrainment_W')
Stats.add_profile('turbulent_entrainment_H')
Stats.add_profile('turbulent_entrainment_QT')
Stats.add_profile('massflux')
Stats.add_profile('massflux_h')
Stats.add_profile('massflux_qt')
Stats.add_profile('massflux_tendency_h')
Stats.add_profile('massflux_tendency_qt')
Stats.add_profile('diffusive_flux_h')
Stats.add_profile('diffusive_flux_u')
Stats.add_profile('diffusive_flux_v')
Stats.add_profile('diffusive_flux_qt')
Stats.add_profile('diffusive_tendency_h')
Stats.add_profile('diffusive_tendency_qt')
Stats.add_profile('total_flux_h')
Stats.add_profile('total_flux_qt')
Stats.add_profile('mixing_length')
Stats.add_profile('updraft_qt_precip')
Stats.add_profile('updraft_thetal_precip')
# Diff mixing lengths: Ignacio
Stats.add_profile('ed_length_scheme')
Stats.add_profile('mixing_length_ratio')
Stats.add_profile('entdet_balance_length')
Stats.add_profile('interdomain_tke_t')
if self.calc_tke:
Stats.add_profile('tke_buoy')
Stats.add_profile('tke_dissipation')
Stats.add_profile('tke_entr_gain')
Stats.add_profile('tke_detr_loss')
Stats.add_profile('tke_shear')
Stats.add_profile('tke_pressure')
Stats.add_profile('tke_interdomain')
Stats.add_profile('tke_transport')
Stats.add_profile('tke_advection')
if self.calc_scalar_var:
Stats.add_profile('Hvar_dissipation')
Stats.add_profile('QTvar_dissipation')
Stats.add_profile('HQTcov_dissipation')
Stats.add_profile('Hvar_entr_gain')
Stats.add_profile('QTvar_entr_gain')
Stats.add_profile('Hvar_detr_loss')
Stats.add_profile('QTvar_detr_loss')
Stats.add_profile('HQTcov_detr_loss')
Stats.add_profile('HQTcov_entr_gain')
Stats.add_profile('Hvar_shear')
Stats.add_profile('QTvar_shear')
Stats.add_profile('HQTcov_shear')
Stats.add_profile('Hvar_rain')
Stats.add_profile('QTvar_rain')
Stats.add_profile('HQTcov_rain')
Stats.add_profile('Hvar_interdomain')
Stats.add_profile('QTvar_interdomain')
Stats.add_profile('HQTcov_interdomain')
return
cpdef io(self, NetCDFIO_Stats Stats, TimeStepping TS):
cdef:
Py_ssize_t k, i
Py_ssize_t kmin = self.Gr.gw
Py_ssize_t kmax = self.Gr.nzg-self.Gr.gw
double [:] mean_entr_sc = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_nh_pressure = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_nh_pressure_adv = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_nh_pressure_drag = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_nh_pressure_b = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_asp_ratio = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_b_coeff = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_detr_sc = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] massflux = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mf_h = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mf_qt = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_frac_turb_entr = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_frac_turb_entr_full = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_turb_entr_W = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_turb_entr_H = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_turb_entr_QT = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_horizontal_KM = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_horizontal_KH = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_sorting_function = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
double [:] mean_b_mix = np.zeros((self.Gr.nzg,), dtype=np.double, order='c')
self.UpdVar.io(Stats, self.Ref)
self.EnvVar.io(Stats, self.Ref)
self.Rain.io(Stats, self.Ref, self.UpdThermo, self.EnvThermo, TS)
Stats.write_profile('eddy_viscosity', self.KM.values[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('eddy_diffusivity', self.KH.values[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_ts('rd', np.mean(self.pressure_plume_spacing))
with nogil:
for k in xrange(self.Gr.gw, self.Gr.nzg-self.Gr.gw):
mf_h[k] = interp2pt(self.massflux_h[k], self.massflux_h[k-1])
mf_qt[k] = interp2pt(self.massflux_qt[k], self.massflux_qt[k-1])
if self.UpdVar.Area.bulkvalues[k] > 0.0:
for i in xrange(self.n_updrafts):
massflux[k] += interp2pt(self.m[i,k], self.m[i,k-1])
mean_entr_sc[k] += self.UpdVar.Area.values[i,k] * self.entr_sc[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_detr_sc[k] += self.UpdVar.Area.values[i,k] * self.detr_sc[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_nh_pressure[k] += self.UpdVar.Area.values[i,k] * self.nh_pressure[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_nh_pressure_b[k] += self.UpdVar.Area.values[i,k] * self.nh_pressure_b[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_nh_pressure_adv[k] += self.UpdVar.Area.values[i,k] * self.nh_pressure_adv[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_nh_pressure_drag[k] += self.UpdVar.Area.values[i,k] * self.nh_pressure_drag[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_asp_ratio[k] += self.UpdVar.Area.values[i,k] * self.asp_ratio[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_b_coeff[k] += self.UpdVar.Area.values[i,k] * self.b_coeff[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_frac_turb_entr_full[k] += self.UpdVar.Area.values[i,k] * self.frac_turb_entr_full[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_frac_turb_entr[k] += self.UpdVar.Area.values[i,k] * self.frac_turb_entr[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_turb_entr_W[k] += self.UpdVar.Area.values[i,k] * self.turb_entr_W[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_turb_entr_H[k] += self.UpdVar.Area.values[i,k] * self.turb_entr_H[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_turb_entr_QT[k] += self.UpdVar.Area.values[i,k] * self.turb_entr_QT[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_horizontal_KM[k] += self.UpdVar.Area.values[i,k] * self.horizontal_KM[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_horizontal_KH[k] += self.UpdVar.Area.values[i,k] * self.horizontal_KH[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_sorting_function[k] += self.UpdVar.Area.values[i,k] * self.sorting_function[i,k]/self.UpdVar.Area.bulkvalues[k]
mean_b_mix[k] += self.UpdVar.Area.values[i,k] * self.b_mix[i,k]/self.UpdVar.Area.bulkvalues[k]
Stats.write_profile('turbulent_entrainment', mean_frac_turb_entr[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('turbulent_entrainment_full', mean_frac_turb_entr_full[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('turbulent_entrainment_W', mean_turb_entr_W[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('turbulent_entrainment_H', mean_turb_entr_H[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('turbulent_entrainment_QT', mean_turb_entr_QT[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('horizontal_KM', mean_horizontal_KM[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('horizontal_KH', mean_horizontal_KH[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('entrainment_sc', mean_entr_sc[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('detrainment_sc', mean_detr_sc[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('sorting_function', mean_sorting_function[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('b_mix', mean_b_mix[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('nh_pressure', mean_nh_pressure[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('nh_pressure_adv', mean_nh_pressure_adv[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('nh_pressure_drag', mean_nh_pressure_drag[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('nh_pressure_b', mean_nh_pressure_b[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('asp_ratio', mean_asp_ratio[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('b_coeff', mean_b_coeff[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('massflux', massflux[self.Gr.gw:self.Gr.nzg-self.Gr.gw ])
Stats.write_profile('massflux_h', mf_h[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('massflux_qt', mf_qt[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('massflux_tendency_h', self.massflux_tendency_h[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('massflux_tendency_qt', self.massflux_tendency_qt[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('diffusive_flux_h', self.diffusive_flux_h[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('diffusive_flux_qt', self.diffusive_flux_qt[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('diffusive_flux_u', self.diffusive_flux_u[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('diffusive_flux_v', self.diffusive_flux_v[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('diffusive_tendency_h', self.diffusive_tendency_h[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('diffusive_tendency_qt', self.diffusive_tendency_qt[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('total_flux_h', np.add(mf_h[self.Gr.gw:self.Gr.nzg-self.Gr.gw],
self.diffusive_flux_h[self.Gr.gw:self.Gr.nzg-self.Gr.gw]))
Stats.write_profile('total_flux_qt', np.add(mf_qt[self.Gr.gw:self.Gr.nzg-self.Gr.gw],
self.diffusive_flux_qt[self.Gr.gw:self.Gr.nzg-self.Gr.gw]))
Stats.write_profile('mixing_length', self.mixing_length[kmin:kmax])
Stats.write_profile('updraft_qt_precip', self.UpdThermo.prec_source_qt_tot[kmin:kmax])
Stats.write_profile('updraft_thetal_precip', self.UpdThermo.prec_source_h_tot[kmin:kmax])
#Different mixing lengths : Ignacio
Stats.write_profile('ed_length_scheme', self.mls[kmin:kmax])
Stats.write_profile('mixing_length_ratio', self.ml_ratio[kmin:kmax])
Stats.write_profile('entdet_balance_length', self.l_entdet[kmin:kmax])
Stats.write_profile('interdomain_tke_t', self.b[kmin:kmax])
if self.calc_tke:
self.compute_covariance_dissipation(self.EnvVar.TKE)
Stats.write_profile('tke_dissipation', self.EnvVar.TKE.dissipation[kmin:kmax])
Stats.write_profile('tke_entr_gain', self.EnvVar.TKE.entr_gain[kmin:kmax])
self.compute_covariance_detr(self.EnvVar.TKE)
Stats.write_profile('tke_detr_loss', self.EnvVar.TKE.detr_loss[kmin:kmax])
Stats.write_profile('tke_shear', self.EnvVar.TKE.shear[kmin:kmax])
Stats.write_profile('tke_buoy', self.EnvVar.TKE.buoy[kmin:kmax])
Stats.write_profile('tke_pressure', self.EnvVar.TKE.press[kmin:kmax])
Stats.write_profile('tke_interdomain', self.EnvVar.TKE.interdomain[kmin:kmax])
self.compute_tke_transport()
Stats.write_profile('tke_transport', self.tke_transport[kmin:kmax])
self.compute_tke_advection()
Stats.write_profile('tke_advection', self.tke_advection[kmin:kmax])
if self.calc_scalar_var:
self.compute_covariance_dissipation(self.EnvVar.Hvar)
Stats.write_profile('Hvar_dissipation', self.EnvVar.Hvar.dissipation[kmin:kmax])
self.compute_covariance_dissipation(self.EnvVar.QTvar)
Stats.write_profile('QTvar_dissipation', self.EnvVar.QTvar.dissipation[kmin:kmax])
self.compute_covariance_dissipation(self.EnvVar.HQTcov)
Stats.write_profile('HQTcov_dissipation', self.EnvVar.HQTcov.dissipation[kmin:kmax])
Stats.write_profile('Hvar_entr_gain', self.EnvVar.Hvar.entr_gain[kmin:kmax])
Stats.write_profile('QTvar_entr_gain', self.EnvVar.QTvar.entr_gain[kmin:kmax])
Stats.write_profile('HQTcov_entr_gain', self.EnvVar.HQTcov.entr_gain[kmin:kmax])
self.compute_covariance_detr(self.EnvVar.Hvar)
self.compute_covariance_detr(self.EnvVar.QTvar)
self.compute_covariance_detr(self.EnvVar.HQTcov)
Stats.write_profile('Hvar_detr_loss', self.EnvVar.Hvar.detr_loss[kmin:kmax])
Stats.write_profile('QTvar_detr_loss', self.EnvVar.QTvar.detr_loss[kmin:kmax])
Stats.write_profile('HQTcov_detr_loss', self.EnvVar.HQTcov.detr_loss[kmin:kmax])
Stats.write_profile('Hvar_shear', self.EnvVar.Hvar.shear[kmin:kmax])
Stats.write_profile('QTvar_shear', self.EnvVar.QTvar.shear[kmin:kmax])
Stats.write_profile('HQTcov_shear', self.EnvVar.HQTcov.shear[kmin:kmax])
Stats.write_profile('Hvar_rain', self.EnvVar.Hvar.rain_src[kmin:kmax])
Stats.write_profile('QTvar_rain', self.EnvVar.QTvar.rain_src[kmin:kmax])
Stats.write_profile('HQTcov_rain', self.EnvVar.HQTcov.rain_src[kmin:kmax])
Stats.write_profile('Hvar_interdomain', self.EnvVar.Hvar.interdomain[kmin:kmax])
Stats.write_profile('QTvar_interdomain', self.EnvVar.QTvar.interdomain[kmin:kmax])
Stats.write_profile('HQTcov_interdomain', self.EnvVar.HQTcov.interdomain[kmin:kmax])
return
# Perform the update of the scheme
cpdef update(self, GridMeanVariables GMV, CasesBase Case, TimeStepping TS):
cdef:
Py_ssize_t k
Py_ssize_t kmin = self.Gr.gw
Py_ssize_t kmax = self.Gr.nzg - self.Gr.gw
self.update_inversion(GMV, Case.inversion_option)
self.compute_pressure_plume_spacing(GMV, Case)
self.wstar = get_wstar(Case.Sur.bflux, self.zi)
if TS.nstep == 0:
self.decompose_environment(GMV, 'values')
if Case.casename == 'DryBubble':
self.EnvThermo.saturation_adjustment(self.EnvVar)
self.UpdThermo.buoyancy(self.UpdVar, self.EnvVar, GMV, self.extrapolate_buoyancy)
self.EnvThermo.microphysics(self.EnvVar, self.Rain, TS.dt)
self.initialize_covariance(GMV, Case)
with nogil:
for k in xrange(self.Gr.nzg):
if self.calc_tke:
self.EnvVar.TKE.values[k] = GMV.TKE.values[k]
if self.calc_scalar_var:
self.EnvVar.Hvar.values[k] = GMV.Hvar.values[k]
self.EnvVar.QTvar.values[k] = GMV.QTvar.values[k]
self.EnvVar.HQTcov.values[k] = GMV.HQTcov.values[k]
self.decompose_environment(GMV, 'values')
if self.use_steady_updrafts:
self.compute_diagnostic_updrafts(GMV, Case)
else:
self.compute_prognostic_updrafts(GMV, Case, TS)
# TODO -maybe not needed? - both diagnostic and prognostic updrafts end with decompose_environment
# But in general ok here without thermodynamics because MF doesnt depend directly on buoyancy
self.decompose_environment(GMV, 'values')
self.update_GMV_MF(GMV, TS)
# (###)
# decompose_environment + EnvThermo.saturation_adjustment + UpdThermo.buoyancy should always be used together
# This ensures that:
# - the buoyancy of updrafts and environment is up to date with the most recent decomposition,
# - the buoyancy of updrafts and environment is updated such that
# the mean buoyancy with repect to reference state alpha_0 is zero.
self.decompose_environment(GMV, 'mf_update')
self.EnvThermo.microphysics(self.EnvVar, self.Rain, TS.dt) # saturation adjustment + rain creation
# Sink of environmental QT and H due to rain creation is applied in tridiagonal solver
self.UpdThermo.buoyancy(self.UpdVar, self.EnvVar, GMV, self.extrapolate_buoyancy)
self.compute_eddy_diffusivities_tke(GMV, Case)
self.update_GMV_ED(GMV, Case, TS)
self.compute_covariance(GMV, Case, TS)
if self.Rain.rain_model == "clima_1m":
# sum updraft and environment rain into bulk rain
self.Rain.sum_subdomains_rain(self.UpdThermo, self.EnvThermo)
# rain fall (all three categories are assumed to be falling though "grid-mean" conditions
self.RainPhysics.solve_rain_fall(GMV, TS, self.Rain.QR, self.Rain.RainArea)
self.RainPhysics.solve_rain_fall(GMV, TS, self.Rain.Upd_QR, self.Rain.Upd_RainArea)
self.RainPhysics.solve_rain_fall(GMV, TS, self.Rain.Env_QR, self.Rain.Env_RainArea)
# rain evaporation (all three categories are assumed to be evaporating in "grid-mean" conditions
self.RainPhysics.solve_rain_evap(GMV, TS, self.Rain.QR, self.Rain.RainArea)
self.RainPhysics.solve_rain_evap(GMV, TS, self.Rain.Upd_QR, self.Rain.Upd_RainArea)
self.RainPhysics.solve_rain_evap(GMV, TS, self.Rain.Env_QR, self.Rain.Env_RainArea)
# update grid-mean cloud fraction and cloud cover
for k in xrange(self.Gr.nzg):
self.EnvVar.Area.values[k] = 1.0 - self.UpdVar.Area.bulkvalues[k]
GMV.cloud_fraction.values[k] = \
self.EnvVar.Area.values[k] * self.EnvVar.cloud_fraction.values[k] +\
self.UpdVar.Area.bulkvalues[k] * self.UpdVar.cloud_fraction[k]
GMV.cloud_cover = min(self.EnvVar.cloud_cover + np.sum(self.UpdVar.cloud_cover), 1)
# Back out the tendencies of the grid mean variables for the whole timestep
# by differencing GMV.new and GMV.values
ParameterizationBase.update(self, GMV, Case, TS)
return
cpdef compute_prognostic_updrafts(self, GridMeanVariables GMV, CasesBase Case, TimeStepping TS):
cdef:
Py_ssize_t iter_
double time_elapsed = 0.0
self.set_subdomain_bcs()
self.UpdVar.set_new_with_values()
self.UpdVar.set_old_with_values()
self.set_updraft_surface_bc(GMV, Case)
self.dt_upd = np.minimum(TS.dt, 0.5 * self.Gr.dz/fmax(np.max(self.UpdVar.W.values),1e-10))
self.UpdThermo.clear_precip_sources()
while time_elapsed < TS.dt:
self.compute_entrainment_detrainment(GMV, Case)
if self.turbulent_entrainment_factor > 1.0e-6:
self.compute_horizontal_eddy_diffusivities(GMV)
self.compute_turbulent_entrainment(GMV,Case)
self.compute_nh_pressure()
self.solve_updraft_velocity_area()
self.solve_updraft_scalars(GMV)
self.UpdThermo.microphysics(self.UpdVar, self.Rain, TS.dt) # causes division error in dry bubble first time step
self.UpdVar.set_values_with_new()
self.zero_area_fraction_cleanup(GMV)
time_elapsed += self.dt_upd
self.dt_upd = np.minimum(TS.dt-time_elapsed, 0.5 * self.Gr.dz/fmax(np.max(self.UpdVar.W.values),1e-10))
# (####)
# TODO - see comment (###)
# It would be better to have a simple linear rule for updating environment here
# instead of calling EnvThermo saturation adjustment scheme for every updraft.
self.decompose_environment(GMV, 'values')
self.EnvThermo.saturation_adjustment(self.EnvVar)
self.UpdThermo.buoyancy(self.UpdVar, self.EnvVar, GMV, self.extrapolate_buoyancy)
self.set_subdomain_bcs()
self.UpdThermo.update_total_precip_sources()
return
cpdef compute_diagnostic_updrafts(self, GridMeanVariables GMV, CasesBase Case):
cdef:
Py_ssize_t i, k
Py_ssize_t gw = self.Gr.gw
double dz = self.Gr.dz
double dzi = self.Gr.dzi
eos_struct sa
mph_struct mph
entr_struct ret
entr_in_struct input
double a,b,c, w, w_km, w_mid, w_low, denom, arg
double entr_w, detr_w, B_k, area_k, w2
self.set_updraft_surface_bc(GMV, Case)
self.compute_entrainment_detrainment(GMV, Case)
with nogil:
for i in xrange(self.n_updrafts):
self.UpdVar.H.values[i,gw] = self.h_surface_bc[i]
self.UpdVar.QT.values[i,gw] = self.qt_surface_bc[i]
# do saturation adjustment
sa = eos(
self.UpdThermo.t_to_prog_fp,self.UpdThermo.prog_to_t_fp,
self.Ref.p0_half[gw], self.UpdVar.QT.values[i,gw],
self.UpdVar.H.values[i,gw]
)
self.UpdVar.QL.values[i,gw] = sa.ql
self.UpdVar.T.values[i, gw] = sa.T
for k in xrange(gw+1, self.Gr.nzg-gw):
denom = 1.0 + self.entr_sc[i,k] * dz
self.UpdVar.H.values[i,k] = (self.UpdVar.H.values[i, k-1] + self.entr_sc[i,k] * dz * GMV.H.values[k])/denom
self.UpdVar.QT.values[i,k] = (self.UpdVar.QT.values[i,k-1] + self.entr_sc[i,k] * dz * GMV.QT.values[k])/denom
# do saturation adjustment
sa = eos(
self.UpdThermo.t_to_prog_fp, self.UpdThermo.prog_to_t_fp,
self.Ref.p0_half[k], self.UpdVar.QT.values[i,k],
self.UpdVar.H.values[i,k]
)
self.UpdVar.QL.values[i,k] = sa.ql
self.UpdVar.T.values[i, k] = sa.T
self.UpdVar.QT.set_bcs(self.Gr)
self.UpdVar.H.set_bcs(self.Gr)
# TODO - see comment (####)
self.decompose_environment(GMV, 'values')
self.EnvThermo.saturation_adjustment(self.EnvVar)
self.UpdThermo.buoyancy(self.UpdVar, self.EnvVar, GMV, self.extrapolate_buoyancy)
# Solve updraft velocity equation
with nogil:
for i in xrange(self.n_updrafts):
self.UpdVar.W.values[i, self.Gr.gw-1] = self.w_surface_bc[i]
self.entr_sc[i,gw] = 2.0 /dz # 0.0 ?
self.detr_sc[i,gw] = 0.0
for k in range(self.Gr.gw, self.Gr.nzg-self.Gr.gw):
area_k = interp2pt(self.UpdVar.Area.values[i,k], self.UpdVar.Area.values[i,k+1])
if area_k >= self.minimum_area:
w_km = self.UpdVar.W.values[i,k-1]
entr_w = interp2pt(self.entr_sc[i,k], self.entr_sc[i,k+1])
detr_w = interp2pt(self.detr_sc[i,k], self.detr_sc[i,k+1])
B_k = interp2pt(self.UpdVar.B.values[i,k], self.UpdVar.B.values[i,k+1])
w2 = ((self.vel_buoy_coeff * B_k + 0.5 * w_km * w_km * dzi)
/(0.5 * dzi +entr_w + (1.0/self.pressure_plume_spacing[i])))
if w2 > 0.0:
self.UpdVar.W.values[i,k] = sqrt(w2)
else:
self.UpdVar.W.values[i,k:] = 0
break
else:
self.UpdVar.W.values[i,k:] = 0
self.UpdVar.W.set_bcs(self.Gr)
cdef double au_lim
with nogil:
for i in xrange(self.n_updrafts):
au_lim = self.max_area
self.UpdVar.Area.values[i,gw] = self.area_surface_bc[i]
w_mid = 0.5* (self.UpdVar.W.values[i,gw])
for k in xrange(gw+1, self.Gr.nzg):
w_low = w_mid
w_mid = interp2pt(self.UpdVar.W.values[i,k],self.UpdVar.W.values[i,k-1])
if w_mid > 0.0:
if self.entr_sc[i,k]>(0.9/dz):
self.entr_sc[i,k] = 0.9/dz
self.UpdVar.Area.values[i,k] = (self.Ref.rho0_half[k-1]*self.UpdVar.Area.values[i,k-1]*w_low/
(1.0-(self.entr_sc[i,k]-self.detr_sc[i,k])*dz)/w_mid/self.Ref.rho0_half[k])
# # Limit the increase in updraft area when the updraft decelerates
if self.UpdVar.Area.values[i,k] > au_lim:
self.UpdVar.Area.values[i,k] = au_lim
self.detr_sc[i,k] =(self.Ref.rho0_half[k-1] * self.UpdVar.Area.values[i,k-1]
* w_low / au_lim / w_mid / self.Ref.rho0_half[k] + self.entr_sc[i,k] * dz -1.0)/dz
else:
# the updraft has terminated so set its area fraction to zero at this height and all heights above
self.UpdVar.Area.values[i,k] = 0.0
self.UpdVar.H.values[i,k] = GMV.H.values[k]
self.UpdVar.QT.values[i,k] = GMV.QT.values[k]
#TODO wouldnt it be more consistent to have here?
#self.UpdVar.QL.values[i,k] = GMV.QL.values[k]
#self.UpdVar.T.values[i,k] = GMV.T.values[k]
sa = eos(self.UpdThermo.t_to_prog_fp,self.UpdThermo.prog_to_t_fp, self.Ref.p0_half[k],
self.UpdVar.QT.values[i,k], self.UpdVar.H.values[i,k])
self.UpdVar.QL.values[i,k] = sa.ql
self.UpdVar.T.values[i,k] = sa.T
# TODO - see comment (####)
self.decompose_environment(GMV, 'values')
self.EnvThermo.saturation_adjustment(self.EnvVar)
self.UpdThermo.buoyancy(self.UpdVar, self.EnvVar, GMV, self.extrapolate_buoyancy)
self.UpdVar.Area.set_bcs(self.Gr)
return
cpdef update_inversion(self,GridMeanVariables GMV, option):
ParameterizationBase.update_inversion(self, GMV,option)
return
cpdef compute_mixing_length(self, double obukhov_length, double ustar, GridMeanVariables GMV):
cdef:
Py_ssize_t k
Py_ssize_t gw = self.Gr.gw
double tau = get_mixing_tau(self.zi, self.wstar)
double l1, l2, l3, z_, N
double l[3]
double ri_grad, shear2
double qt_dry, th_dry, t_cloudy, qv_cloudy, qt_cloudy, th_cloudy
double lh, cpm, prefactor, d_buoy_thetal_dry, d_buoy_qt_dry
double d_buoy_thetal_cloudy, d_buoy_qt_cloudy, d_buoy_thetal_total, d_buoy_qt_total
double grad_thl_plus=0.0, grad_qt_plus=0.0, grad_thv_plus=0.0
double thv, grad_qt, grad_qt_low, grad_thv_low, grad_thv
double grad_b_thl, grad_b_qt
double m_eps = 1.0e-9 # Epsilon to avoid zero
double a, c_neg, wc_upd_nn, wc_env, frac_turb_entr_half
if (self.mixing_scheme == 'sbl'):
for k in xrange(gw, self.Gr.nzg-gw):
z_ = self.Gr.z_half[k]
# kz scale (surface layer)
if obukhov_length < 0.0: #unstable
l2 = vkb * z_ /(sqrt(self.EnvVar.TKE.values[self.Gr.gw]/ustar/ustar)*self.tke_ed_coeff) * fmin(
(1.0 - 100.0 * z_/obukhov_length)**0.2, 1.0/vkb )
else: # neutral or stable
l2 = vkb * z_ /(sqrt(self.EnvVar.TKE.values[self.Gr.gw]/ustar/ustar)*self.tke_ed_coeff)
# Shear-dissipation TKE equilibrium scale (Stable)
shear2 = pow((GMV.U.values[k+1] - GMV.U.values[k-1]) * 0.5 * self.Gr.dzi, 2) + \
pow((GMV.V.values[k+1] - GMV.V.values[k-1]) * 0.5 * self.Gr.dzi, 2) + \
pow((self.EnvVar.W.values[k] - self.EnvVar.W.values[k-1]) * self.Gr.dzi, 2)
qt_dry = self.EnvThermo.qt_dry[k]
th_dry = self.EnvThermo.th_dry[k]
t_cloudy = self.EnvThermo.t_cloudy[k]
qv_cloudy = self.EnvThermo.qv_cloudy[k]
qt_cloudy = self.EnvThermo.qt_cloudy[k]
th_cloudy = self.EnvThermo.th_cloudy[k]
lh = latent_heat(t_cloudy)
cpm = cpm_c(qt_cloudy)
grad_thl_low = grad_thl_plus
grad_qt_low = grad_qt_plus
grad_thl_plus = (self.EnvVar.THL.values[k+1] - self.EnvVar.THL.values[k]) * self.Gr.dzi
grad_qt_plus = (self.EnvVar.QT.values[k+1] - self.EnvVar.QT.values[k]) * self.Gr.dzi
grad_thl = interp2pt(grad_thl_low, grad_thl_plus)
grad_qt = interp2pt(grad_qt_low, grad_qt_plus)
# g/theta_ref
prefactor = g * ( Rd / self.Ref.alpha0_half[k] /self.Ref.p0_half[k]) * exner_c(self.Ref.p0_half[k])
d_buoy_thetal_dry = prefactor * (1.0 + (eps_vi-1.0) * qt_dry)
d_buoy_qt_dry = prefactor * th_dry * (eps_vi-1.0)
if self.EnvVar.cloud_fraction.values[k] > 0.0:
d_buoy_thetal_cloudy = (prefactor * (1.0 + eps_vi * (1.0 + lh / Rv / t_cloudy) * qv_cloudy - qt_cloudy )
/ (1.0 + lh * lh / cpm / Rv / t_cloudy / t_cloudy * qv_cloudy))
d_buoy_qt_cloudy = (lh / cpm / t_cloudy * d_buoy_thetal_cloudy - prefactor) * th_cloudy
else:
d_buoy_thetal_cloudy = 0.0
d_buoy_qt_cloudy = 0.0
d_buoy_thetal_total = (self.EnvVar.cloud_fraction.values[k] * d_buoy_thetal_cloudy
+ (1.0-self.EnvVar.cloud_fraction.values[k]) * d_buoy_thetal_dry)
d_buoy_qt_total = (self.EnvVar.cloud_fraction.values[k] * d_buoy_qt_cloudy
+ (1.0-self.EnvVar.cloud_fraction.values[k]) * d_buoy_qt_dry)
# Partial buoyancy gradients
grad_b_thl = grad_thl * d_buoy_thetal_total
grad_b_qt = grad_qt * d_buoy_qt_total
ri_grad = fmin( grad_b_thl/fmax(shear2, m_eps) + grad_b_qt/fmax(shear2, m_eps) , 0.25)
# Turbulent Prandtl number:
if obukhov_length > 0.0 and ri_grad>0.0: #stable
# CSB (Dan Li, 2019), with Pr_neutral=0.74 and w1=40.0/13.0
self.prandtl_nvec[k] = self.prandtl_number*( 2.0*ri_grad/
(1.0+(53.0/13.0)*ri_grad -sqrt( (1.0+(53.0/13.0)*ri_grad)**2.0 - 4.0*ri_grad ) ) )
else:
self.prandtl_nvec[k] = self.prandtl_number
l3 = sqrt(self.tke_diss_coeff/fmax(self.tke_ed_coeff, m_eps)) * sqrt(self.EnvVar.TKE.values[k])
l3 /= sqrt(fmax(shear2 - grad_b_thl/self.prandtl_nvec[k] - grad_b_qt/self.prandtl_nvec[k], m_eps))
if ( shear2 - grad_b_thl/self.prandtl_nvec[k] - grad_b_qt/self.prandtl_nvec[k] < m_eps):
l3 = 1.0e6
# Limiting stratification scale (Deardorff, 1976)
thv = theta_virt_c(self.Ref.p0_half[k], self.EnvVar.T.values[k], self.EnvVar.QT.values[k],
self.EnvVar.QL.values[k])
grad_thv_low = grad_thv_plus
grad_thv_plus = ( theta_virt_c(self.Ref.p0_half[k+1], self.EnvVar.T.values[k+1], self.EnvVar.QT.values[k+1],
self.EnvVar.QL.values[k+1]) - thv) * self.Gr.dzi
grad_thv = interp2pt(grad_thv_low, grad_thv_plus)
N = sqrt(fmax(g/thv*grad_thv, 0.0))
if N > 0.0:
l1 = fmin(sqrt(fmax(self.static_stab_coeff*self.EnvVar.TKE.values[k],0.0))/N, 1.0e6)
else:
l1 = 1.0e6
l[0]=l1; l[1]=l3; l[2]=l2;
j = 0
while(j<len(l)):
if l[j]<m_eps or l[j]>1.0e6:
l[j] = 1.0e6
j += 1
self.mls[k] = np.argmin(l)
self.mixing_length[k] = auto_smooth_minimum(l, 0.1)
self.ml_ratio[k] = self.mixing_length[k]/l[int(self.mls[k])]
elif (self.mixing_scheme == 'sbtd_eq'):
for k in xrange(gw, self.Gr.nzg-gw):
z_ = self.Gr.z_half[k]
# kz scale (surface layer)
if obukhov_length < 0.0: #unstable
l2 = vkb * z_ /(sqrt(self.EnvVar.TKE.values[self.Gr.gw]/ustar/ustar)*self.tke_ed_coeff) * fmin(
(1.0 - 100.0 * z_/obukhov_length)**0.2, 1.0/vkb )
else: # neutral or stable
l2 = vkb * z_ /(sqrt(self.EnvVar.TKE.values[self.Gr.gw]/ustar/ustar)*self.tke_ed_coeff)
# Buoyancy-shear-subdomain exchange-dissipation TKE equilibrium scale
shear2 = pow((GMV.U.values[k+1] - GMV.U.values[k-1]) * 0.5 * self.Gr.dzi, 2) + \
pow((GMV.V.values[k+1] - GMV.V.values[k-1]) * 0.5 * self.Gr.dzi, 2) + \
pow((self.EnvVar.W.values[k] - self.EnvVar.W.values[k-1]) * self.Gr.dzi, 2)
qt_dry = self.EnvThermo.qt_dry[k]
th_dry = self.EnvThermo.th_dry[k]
t_cloudy = self.EnvThermo.t_cloudy[k]
qv_cloudy = self.EnvThermo.qv_cloudy[k]
qt_cloudy = self.EnvThermo.qt_cloudy[k]
th_cloudy = self.EnvThermo.th_cloudy[k]
lh = latent_heat(t_cloudy)
cpm = cpm_c(qt_cloudy)
grad_thl_low = grad_thl_plus
grad_qt_low = grad_qt_plus
grad_thl_plus = (self.EnvVar.THL.values[k+1] - self.EnvVar.THL.values[k]) * self.Gr.dzi
grad_qt_plus = (self.EnvVar.QT.values[k+1] - self.EnvVar.QT.values[k]) * self.Gr.dzi
grad_thl = interp2pt(grad_thl_low, grad_thl_plus)
grad_qt = interp2pt(grad_qt_low, grad_qt_plus)
# g/theta_ref
prefactor = g * ( Rd / self.Ref.alpha0_half[k] /self.Ref.p0_half[k]) * exner_c(self.Ref.p0_half[k])
d_buoy_thetal_dry = prefactor * (1.0 + (eps_vi-1.0) * qt_dry)
d_buoy_qt_dry = prefactor * th_dry * (eps_vi-1.0)
if self.EnvVar.cloud_fraction.values[k] > 0.0:
d_buoy_thetal_cloudy = (prefactor * (1.0 + eps_vi * (1.0 + lh / Rv / t_cloudy) * qv_cloudy - qt_cloudy )
/ (1.0 + lh * lh / cpm / Rv / t_cloudy / t_cloudy * qv_cloudy))
d_buoy_qt_cloudy = (lh / cpm / t_cloudy * d_buoy_thetal_cloudy - prefactor) * th_cloudy
else:
d_buoy_thetal_cloudy = 0.0
d_buoy_qt_cloudy = 0.0
d_buoy_thetal_total = (self.EnvVar.cloud_fraction.values[k] * d_buoy_thetal_cloudy
+ (1.0-self.EnvVar.cloud_fraction.values[k]) * d_buoy_thetal_dry)
d_buoy_qt_total = (self.EnvVar.cloud_fraction.values[k] * d_buoy_qt_cloudy
+ (1.0-self.EnvVar.cloud_fraction.values[k]) * d_buoy_qt_dry)
# Partial buoyancy gradients
grad_b_thl = grad_thl * d_buoy_thetal_total
grad_b_qt = grad_qt * d_buoy_qt_total
ri_grad = fmin( grad_b_thl/fmax(shear2, m_eps) + grad_b_qt/fmax(shear2, m_eps) , 0.25)
# Turbulent Prandtl number:
if obukhov_length > 0.0 and ri_grad>0.0: #stable
# CSB (Dan Li, 2019), with Pr_neutral=0.74 and w1=40.0/13.0
self.prandtl_nvec[k] = self.prandtl_number*( 2.0*ri_grad/
(1.0+(53.0/13.0)*ri_grad -sqrt( (1.0+(53.0/13.0)*ri_grad)**2.0 - 4.0*ri_grad ) ) )
else:
self.prandtl_nvec[k] = self.prandtl_number
# Production/destruction terms
a = self.tke_ed_coeff*(shear2 - grad_b_thl/self.prandtl_nvec[k] - grad_b_qt/self.prandtl_nvec[k])* sqrt(self.EnvVar.TKE.values[k])
# Dissipation term
c_neg = self.tke_diss_coeff*self.EnvVar.TKE.values[k]*sqrt(self.EnvVar.TKE.values[k])
# Subdomain exchange term
self.b[k] = 0.0
for nn in xrange(self.n_updrafts):
wc_upd_nn = (self.UpdVar.W.values[nn,k] + self.UpdVar.W.values[nn,k-1])/2.0
wc_env = (self.EnvVar.W.values[k] + self.EnvVar.W.values[k-1])/2.0
frac_turb_entr_half = interp2pt(self.frac_turb_entr_full[nn,k],self.frac_turb_entr_full[nn,k-1])
self.b[k] += self.UpdVar.Area.values[nn,k]*wc_upd_nn*self.detr_sc[nn,k]/(1.0-self.UpdVar.Area.bulkvalues[k])*(
(wc_upd_nn-wc_env)*(wc_upd_nn-wc_env)/2.0-self.EnvVar.TKE.values[k]) - self.UpdVar.Area.values[nn,k]*wc_upd_nn*(
wc_upd_nn-wc_env)*frac_turb_entr_half*wc_env/(1.0-self.UpdVar.Area.bulkvalues[k])
if abs(a) > m_eps and 4.0*a*c_neg > - self.b[k]*self.b[k]:
self.l_entdet[k] = fmax( -self.b[k]/2.0/a + sqrt( self.b[k]*self.b[k] + 4.0*a*c_neg )/2.0/a, 0.0)
elif abs(a) < m_eps and abs(self.b[k]) > m_eps:
self.l_entdet[k] = c_neg/self.b[k]
l3 = self.l_entdet[k]
# Limiting stratification scale (Deardorff, 1976)
thv = theta_virt_c(self.Ref.p0_half[k], self.EnvVar.T.values[k], self.EnvVar.QT.values[k],
self.EnvVar.QL.values[k])
grad_thv_low = grad_thv_plus
grad_thv_plus = ( theta_virt_c(self.Ref.p0_half[k+1], self.EnvVar.T.values[k+1], self.EnvVar.QT.values[k+1],
self.EnvVar.QL.values[k+1]) - thv) * self.Gr.dzi
grad_thv = interp2pt(grad_thv_low, grad_thv_plus)
# Effective static stability using environmental mean.
# Set lambda for now to environmental cloud_fraction (TBD: Rain)
grad_th_eff = (1.0-self.EnvVar.cloud_fraction.values[k])*grad_thv + self.EnvVar.cloud_fraction.values[k]*(
1.0/exp( - latent_heat(self.EnvVar.T.values[k]) * self.EnvVar.QL.values[k]
/ cpm_c(self.EnvVar.QT.values[k]) / self.EnvVar.T.values[k] )*(
(1.0+ (eps_vi-1.0)*self.EnvVar.QT.values[k])*grad_thl+(eps_vi-1.0)*self.EnvVar.THL.values[k]*grad_qt))
N = sqrt(fmax(g/thv*grad_th_eff, 0.0))
if N > 0.0:
l1 = fmin(sqrt(fmax(self.static_stab_coeff*self.EnvVar.TKE.values[k],0.0))/N, 1.0e6)
else:
l1 = 1.0e6
l[0]=l1; l[1]=l3; l[2]=l2;
j = 0
while(j<len(l)):
if l[j]<m_eps or l[j]>1.0e6:
l[j] = 1.0e6
j += 1
self.mls[k] = np.argmin(l)
self.mixing_length[k] = lamb_smooth_minimum(l, 0.1, 1.5)
self.ml_ratio[k] = self.mixing_length[k]/l[int(self.mls[k])]
else:
# default mixing scheme , see Tan et al. (2018)
with nogil:
for k in xrange(gw, self.Gr.nzg-gw):
l1 = tau * sqrt(fmax(self.EnvVar.TKE.values[k],0.0))
z_ = self.Gr.z_half[k]
if obukhov_length < 0.0: #unstable
l2 = vkb * z_ * ( (1.0 - 100.0 * z_/obukhov_length)**0.2 )
elif obukhov_length > 0.0: #stable
l2 = vkb * z_ / (1. + 2.7 *z_/obukhov_length)
l1 = 1.0/m_eps
else:
l2 = vkb * z_
self.mixing_length[k] = fmax( 1.0/(1.0/fmax(l1,m_eps) + 1.0/l2), 1e-3)
self.prandtl_nvec[k] = 1.0