forked from szy21/pycles_GCM
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Microphysics.pyx
937 lines (754 loc) · 46.6 KB
/
Microphysics.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
#!python
#cython: boundscheck=False
#cython: wraparound=False
#cython: initializedcheck=False
#cython: cdivision=True
cimport numpy as np
import numpy as np
cimport Grid
cimport Lookup
cimport PrognosticVariables
cimport DiagnosticVariables
cimport ReferenceState
cimport ParallelMPI
cimport TimeStepping
from NetCDFIO cimport NetCDFIO_Stats
from Thermodynamics cimport LatentHeat, ClausiusClapeyron
from Microphysics_Arctic_1M cimport Microphysics_Arctic_1M
from libc.math cimport fmax, fmin, fabs, pow
from thermodynamic_functions cimport cpm_c, pv_c, pd_c
include 'parameters.pxi'
cdef extern from "microphysics.h":
void microphysics_stokes_sedimentation_velocity(Grid.DimStruct *dims, double* density, double ccn, double* ql, double* qt_velocity)
cdef extern from "scalar_advection.h":
void compute_advective_fluxes_a(Grid.DimStruct *dims, double *rho0, double *rho0_half, double *velocity, double *scalar, double* flux, int d, int scheme) nogil
cdef extern from "microphysics_sb.h":
void sb_sedimentation_velocity_liquid(Grid.DimStruct *dims, double* density, double ccn, double* ql, double* qt_velocity)nogil
cdef class No_Microphysics_Dry:
def __init__(self, ParallelMPI.ParallelMPI Par, LatentHeat LH, namelist):
LH.Lambda_fp = lambda_constant
LH.L_fp = latent_heat_constant
self.thermodynamics_type = 'dry'
return
cpdef initialize(self, Grid.Grid Gr, PrognosticVariables.PrognosticVariables PV,DiagnosticVariables.DiagnosticVariables DV, NetCDFIO_Stats NS, ParallelMPI.ParallelMPI Pa):
return
cpdef update(self, Grid.Grid Gr, ReferenceState.ReferenceState Ref, Th, PrognosticVariables.PrognosticVariables PV, DiagnosticVariables.DiagnosticVariables DV, TimeStepping.TimeStepping TS,ParallelMPI.ParallelMPI Pa):
return
cpdef stats_io(self, Grid.Grid Gr, ReferenceState.ReferenceState Ref, Th, PrognosticVariables.PrognosticVariables PV, DiagnosticVariables.DiagnosticVariables DV, NetCDFIO_Stats NS, ParallelMPI.ParallelMPI Pa):
return
cdef class No_Microphysics_SA:
def __init__(self, ParallelMPI.ParallelMPI Par, LatentHeat LH, namelist):
LH.Lambda_fp = lambda_constant
LH.L_fp = latent_heat_variable
self.thermodynamics_type = 'SA'
#also set local versions
self.Lambda_fp = lambda_constant
self.L_fp = latent_heat_variable
# Extract case-specific parameter values from the namelist
# Get number concentration of cloud condensation nuclei (1/m^3)
try:
self.ccn = namelist['microphysics']['ccn']
except:
self.ccn = 100.0e6
try:
self.order = namelist['scalar_transport']['order_sedimentation']
except:
self.order = namelist['scalar_transport']['order']
try:
self.cloud_sedimentation = namelist['microphysics']['cloud_sedimentation']
except:
self.cloud_sedimentation = False
if namelist['meta']['casename'] == 'DYCOMS_RF02':
self.stokes_sedimentation = True
else:
self.stokes_sedimentation = False
return
cpdef initialize(self, Grid.Grid Gr, PrognosticVariables.PrognosticVariables PV,DiagnosticVariables.DiagnosticVariables DV, NetCDFIO_Stats NS, ParallelMPI.ParallelMPI Pa):
if self.cloud_sedimentation:
DV.add_variables('w_qt', 'm/s', r'w_ql', 'cloud liquid water sedimentation velocity', 'sym', Pa)
NS.add_profile('qt_sedimentation_flux', Gr, Pa)
NS.add_profile('s_qt_sedimentation_source',Gr,Pa)
return
cpdef update(self, Grid.Grid Gr, ReferenceState.ReferenceState Ref, Th, PrognosticVariables.PrognosticVariables PV, DiagnosticVariables.DiagnosticVariables DV, TimeStepping.TimeStepping TS,ParallelMPI.ParallelMPI Pa):
cdef:
Py_ssize_t wqt_shift
Py_ssize_t ql_shift = DV.get_varshift(Gr,'ql')
if self.cloud_sedimentation:
wqt_shift = DV.get_varshift(Gr, 'w_qt')
if self.stokes_sedimentation:
microphysics_stokes_sedimentation_velocity(&Gr.dims, &Ref.rho0_half[0], self.ccn, &DV.values[ql_shift], &DV.values[wqt_shift])
else:
sb_sedimentation_velocity_liquid(&Gr.dims, &Ref.rho0_half[0], self.ccn, &DV.values[ql_shift], &DV.values[wqt_shift])
return
cpdef stats_io(self, Grid.Grid Gr, ReferenceState.ReferenceState Ref, Th, PrognosticVariables.PrognosticVariables PV, DiagnosticVariables.DiagnosticVariables DV, NetCDFIO_Stats NS, ParallelMPI.ParallelMPI Pa):
cdef:
Py_ssize_t gw = Gr.dims.gw
double[:] dummy = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
Py_ssize_t ql_shift = DV.get_varshift(Gr, 'ql')
Py_ssize_t qv_shift = DV.get_varshift(Gr, 'qv')
Py_ssize_t t_shift = DV.get_varshift(Gr, 'temperature')
Py_ssize_t qt_shift = PV.get_varshift(Gr, 'qt')
Py_ssize_t s_shift
Py_ssize_t wqt_shift
double[:] s_src = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
double[:] tmp
if self.cloud_sedimentation:
s_shift = PV.get_varshift(Gr, 's')
wqt_shift = DV.get_varshift(Gr,'w_qt')
compute_advective_fluxes_a(&Gr.dims, &Ref.rho0[0], &Ref.rho0_half[0], &DV.values[wqt_shift], &DV.values[ql_shift], &dummy[0], 2, self.order)
tmp = Pa.HorizontalMean(Gr, &dummy[0])
NS.write_profile('qt_sedimentation_flux', tmp[gw:-gw], Pa)
compute_qt_sedimentation_s_source(&Gr.dims, &Ref.p0_half[0], &Ref.rho0_half[0], &dummy[0],
&PV.values[qt_shift], &DV.values[qv_shift],&DV.values[t_shift], &s_src[0], self.Lambda_fp,
self.L_fp, Gr.dims.dx[2], 2)
tmp = Pa.HorizontalMean(Gr, &s_src[0])
NS.write_profile('s_qt_sedimentation_source', tmp[gw:-gw], Pa)
return
cdef extern from "microphysics_sb.h":
double sb_rain_shape_parameter_0(double density, double qr, double Dm) nogil
double sb_rain_shape_parameter_1(double density, double qr, double Dm) nogil
double sb_rain_shape_parameter_2(double density, double qr, double Dm) nogil
double sb_rain_shape_parameter_4(double density, double qr, double Dm) nogil
double sb_droplet_nu_0(double density, double ql) nogil
double sb_droplet_nu_1(double density, double ql) nogil
double sb_droplet_nu_2(double density, double ql) nogil
void sb_sedimentation_velocity_rain(Grid.DimStruct *dims, double (*rain_mu)(double,double,double),
double* density, double* nr, double* qr, double* nr_velocity, double* qr_velocity) nogil
# void sb_sedimentation_velocity_liquid(Grid.DimStruct *dims, double* density, double ccn, double* ql, double* qt_velocity)nogil
void sb_microphysics_sources(Grid.DimStruct *dims, Lookup.LookupStruct *LT, double (*lam_fp)(double), double (*L_fp)(double, double),
double (*rain_mu)(double,double,double), double (*droplet_nu)(double,double),
double* density, double* p0, double* temperature, double* qt, double ccn,
double* ql, double* nr, double* qr, double dt, double* nr_tendency_micro, double* qr_tendency_micro,
double* nr_tendency, double* qr_tendency) nogil
void sb_qt_source_formation(Grid.DimStruct *dims,double* qr_tendency, double* qt_tendency )nogil
void sb_entropy_source_formation(Grid.DimStruct *dims, Lookup.LookupStruct *LT, double (*lam_fp)(double),
double (*L_fp)(double, double), double* p0, double* T, double* Twet, double* qt,
double* qv, double* qr_tendency, double* entropy_tendency)nogil
void sb_entropy_source_heating(Grid.DimStruct *dims, double* T, double* Twet, double* qr, double* w_qr, double* w,
double* entropy_tendency)nogil
void sb_entropy_source_drag(Grid.DimStruct *dims, double* T, double* qr, double* w_qr, double* entropy_tendency)nogil
void sb_autoconversion_rain_wrapper(Grid.DimStruct *dims, double (*droplet_nu)(double,double), double* density,
double ccn, double* ql, double* qr, double* nr_tendency, double* qr_tendency) nogil
void sb_accretion_rain_wrapper(Grid.DimStruct *dims, double* density, double* ql, double* qr, double* qr_tendency)nogil
void sb_selfcollection_breakup_rain_wrapper(Grid.DimStruct *dims, double (*rain_mu)(double,double,double),
double* density, double* nr, double* qr, double* nr_tendency)nogil
void sb_evaporation_rain_wrapper(Grid.DimStruct *dims, Lookup.LookupStruct *LT, double (*lam_fp)(double), double (*L_fp)(double, double),
double (*rain_mu)(double,double,double), double* density, double* p0, double* temperature, double* qt,
double* ql, double* nr, double* qr, double* nr_tendency, double* qr_tendency)nogil
cdef extern from "scalar_advection.h":
void compute_qt_sedimentation_s_source(Grid.DimStruct *dims, double *p0_half, double* rho0_half, double *flux,
double* qt, double* qv, double* T, double* tendency, double (*lam_fp)(double),
double (*L_fp)(double, double), double dx, ssize_t d)nogil
cdef extern from "microphysics.h":
void microphysics_wetbulb_temperature(Grid.DimStruct *dims, Lookup.LookupStruct *LT, double* p0, double* s,
double* qt, double* T, double* Twet )nogil
cdef class Microphysics_SB_Liquid:
def __init__(self, ParallelMPI.ParallelMPI Par, LatentHeat LH, namelist):
# Create the appropriate linkages to the bulk thermodynamics
LH.Lambda_fp = lambda_constant
LH.L_fp = latent_heat_variable
self.thermodynamics_type = 'SA'
#also set local versions
self.Lambda_fp = lambda_constant
self.L_fp = latent_heat_variable
self.CC = ClausiusClapeyron()
self.CC.initialize(namelist, LH, Par)
# Extract case-specific parameter values from the namelist
# Set the number concentration of cloud condensation nuclei (1/m^3)
# First set a default value, then set a case specific value, which can then be overwritten using namelist options
self.ccn = 100.0e6
if namelist['meta']['casename'] == 'DYCOMS_RF02':
self.ccn = 55.0e6
elif namelist['meta']['casename'] == 'Rico':
self.ccn = 70.0e6
try:
self.ccn = namelist['microphysics']['ccn']
except:
pass
# Set option for calculation of mu (distribution shape parameter)
try:
mu_opt = namelist['microphysics']['SB_Liquid']['mu_rain']
if mu_opt == 1:
self.compute_rain_shape_parameter = sb_rain_shape_parameter_1
elif mu_opt == 2:
self.compute_rain_shape_parameter = sb_rain_shape_parameter_2
elif mu_opt == 4:
self.compute_rain_shape_parameter = sb_rain_shape_parameter_4
elif mu_opt == 0:
self.compute_rain_shape_parameter = sb_rain_shape_parameter_0
else:
Par.root_print("SB_Liquid mu_rain option not recognized, defaulting to option 1")
self.compute_rain_shape_parameter = sb_rain_shape_parameter_1
except:
Par.root_print("SB_Liquid mu_rain option not selected, defaulting to option 1")
self.compute_rain_shape_parameter = sb_rain_shape_parameter_1
# Set option for calculation of nu parameter of droplet distribution
try:
nu_opt = namelist['microphysics']['SB_Liquid']['nu_droplet']
if nu_opt == 0:
self.compute_droplet_nu = sb_droplet_nu_0
elif nu_opt == 1:
self.compute_droplet_nu = sb_droplet_nu_1
elif nu_opt ==2:
self.compute_droplet_nu = sb_droplet_nu_2
else:
Par.root_print("SB_Liquid nu_droplet_option not recognized, defaulting to option 0")
self.compute_droplet_nu = sb_droplet_nu_0
except:
Par.root_print("SB_Liquid nu_droplet_option not selected, defaulting to option 0")
self.compute_droplet_nu = sb_droplet_nu_0
try:
self.order = namelist['scalar_transport']['order_sedimentation']
except:
self.order = namelist['scalar_transport']['order']
try:
self.cloud_sedimentation = namelist['microphysics']['cloud_sedimentation']
except:
self.cloud_sedimentation = False
if namelist['meta']['casename'] == 'DYCOMS_RF02':
self.stokes_sedimentation = True
else:
self.stokes_sedimentation = False
return
cpdef initialize(self, Grid.Grid Gr, PrognosticVariables.PrognosticVariables PV, DiagnosticVariables.DiagnosticVariables DV, NetCDFIO_Stats NS, ParallelMPI.ParallelMPI Pa):
# add prognostic variables for mass and number of rain
PV.add_variable('nr', '1/kg', r'n_r', 'rain droplet number concentration','sym','scalar',Pa)
PV.add_variable('qr', 'kg/kg', r'q_r', 'rain water specific humidity','sym','scalar',Pa)
# add sedimentation velocities as diagnostic variables
DV.add_variables('w_qr', 'm/s', r'w_{qr}', 'rain mass sedimentation veloctiy', 'sym', Pa)
DV.add_variables('w_nr', 'm/s', r'w_{nr}', 'rain number sedimentation velocity', 'sym', Pa)
if self.cloud_sedimentation:
DV.add_variables('w_qt', 'm/s', r'w_ql', 'cloud liquid water sedimentation velocity', 'sym', Pa)
NS.add_profile('qt_sedimentation_flux', Gr, Pa)
NS.add_profile('s_qt_sedimentation_source',Gr,Pa)
# add wet bulb temperature
DV.add_variables('temperature_wb', 'K', r'T_{wb}','wet bulb temperature','sym', Pa)
# add statistical output for the class
NS.add_profile('qr_sedimentation_flux', Gr, Pa)
NS.add_profile('nr_sedimentation_flux', Gr, Pa)
NS.add_profile('qr_autoconversion', Gr, Pa)
NS.add_profile('nr_autoconversion', Gr, Pa)
NS.add_profile('s_autoconversion', Gr, Pa)
NS.add_profile('nr_selfcollection', Gr, Pa)
NS.add_profile('qr_accretion', Gr, Pa)
NS.add_profile('s_accretion', Gr, Pa)
NS.add_profile('nr_evaporation', Gr, Pa)
NS.add_profile('qr_evaporation', Gr,Pa)
NS.add_profile('s_evaporation', Gr,Pa)
NS.add_profile('s_precip_heating', Gr, Pa)
NS.add_profile('s_precip_drag', Gr, Pa)
return
cpdef update(self, Grid.Grid Gr, ReferenceState.ReferenceState Ref, Th, PrognosticVariables.PrognosticVariables PV, DiagnosticVariables.DiagnosticVariables DV, TimeStepping.TimeStepping TS, ParallelMPI.ParallelMPI Pa):
cdef:
Py_ssize_t t_shift = DV.get_varshift(Gr, 'temperature')
Py_ssize_t ql_shift = DV.get_varshift(Gr,'ql')
Py_ssize_t qv_shift = DV.get_varshift(Gr,'qv')
Py_ssize_t nr_shift = PV.get_varshift(Gr, 'nr')
Py_ssize_t qr_shift = PV.get_varshift(Gr, 'qr')
Py_ssize_t qt_shift = PV.get_varshift(Gr, 'qt')
Py_ssize_t w_shift = PV.get_varshift(Gr, 'w')
double dt = TS.dt
Py_ssize_t wqr_shift = DV.get_varshift(Gr, 'w_qr')
Py_ssize_t wnr_shift = DV.get_varshift(Gr, 'w_nr')
Py_ssize_t wqt_shift
double[:] qr_tend_micro = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
double[:] nr_tend_micro = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
sb_microphysics_sources(&Gr.dims, &self.CC.LT.LookupStructC, self.Lambda_fp, self.L_fp, self.compute_rain_shape_parameter,
self.compute_droplet_nu, &Ref.rho0_half[0], &Ref.p0_half[0], &DV.values[t_shift],
&PV.values[qt_shift], self.ccn, &DV.values[ql_shift], &PV.values[nr_shift],
&PV.values[qr_shift], dt, &nr_tend_micro[0], &qr_tend_micro[0], &PV.tendencies[nr_shift], &PV.tendencies[qr_shift] )
sb_sedimentation_velocity_rain(&Gr.dims,self.compute_rain_shape_parameter,
&Ref.rho0_half[0],&PV.values[nr_shift], &PV.values[qr_shift],
&DV.values[wnr_shift], &DV.values[wqr_shift])
if self.cloud_sedimentation:
wqt_shift = DV.get_varshift(Gr, 'w_qt')
if self.stokes_sedimentation:
microphysics_stokes_sedimentation_velocity(&Gr.dims, &Ref.rho0_half[0], self.ccn, &DV.values[ql_shift], &DV.values[wqt_shift])
else:
sb_sedimentation_velocity_liquid(&Gr.dims, &Ref.rho0_half[0], self.ccn, &DV.values[ql_shift], &DV.values[wqt_shift])
# update the Boundary conditions and ghost cells of the sedimentation velocities
# wnr_nv = DV.name_index['w_nr']
# wqr_nv = DV.name_index['w_qr']
# DV.communicate_variable(Gr,Pa,wnr_nv)
# DV.communicate_variable(Gr,Pa,wqr_nv )
sb_qt_source_formation(&Gr.dims, &qr_tend_micro[0], &PV.tendencies[qt_shift])
cdef:
Py_ssize_t tw_shift = DV.get_varshift(Gr, 'temperature_wb')
Py_ssize_t s_shift
Py_ssize_t thli_shift
if 's' in PV.name_index:
s_shift = PV.get_varshift(Gr, 's')
microphysics_wetbulb_temperature(&Gr.dims, &self.CC.LT.LookupStructC, &Ref.p0_half[0], &PV.values[s_shift],
&PV.values[qt_shift], &DV.values[t_shift], &DV.values[tw_shift])
sb_entropy_source_formation(&Gr.dims, &self.CC.LT.LookupStructC, self.Lambda_fp, self.L_fp, &Ref.p0_half[0],
&DV.values[t_shift], &DV.values[tw_shift], &PV.values[qt_shift], &DV.values[qv_shift],
&qr_tend_micro[0], &PV.tendencies[s_shift])
sb_entropy_source_heating(&Gr.dims, &DV.values[t_shift], &DV.values[tw_shift], &PV.values[qr_shift],
&DV.values[wqr_shift], &PV.values[w_shift], &PV.tendencies[s_shift])
sb_entropy_source_drag(&Gr.dims, &DV.values[t_shift], &PV.values[qr_shift], &DV.values[wqr_shift], &PV.tendencies[s_shift])
else:
thli_shfit = PV.get_varshift(Gr, 'thli')
s_shift = DV.get_varshift(Gr, 's')
microphysics_wetbulb_temperature(&Gr.dims, &self.CC.LT.LookupStructC, &Ref.p0_half[0], &DV.values[s_shift],
&PV.values[qt_shift], &DV.values[t_shift], &DV.values[tw_shift])
return
#
cpdef stats_io(self, Grid.Grid Gr, ReferenceState.ReferenceState Ref, Th, PrognosticVariables.PrognosticVariables PV, DiagnosticVariables.DiagnosticVariables DV, NetCDFIO_Stats NS, ParallelMPI.ParallelMPI Pa):
cdef:
Py_ssize_t i, j, k, ijk
Py_ssize_t gw = Gr.dims.gw
Py_ssize_t imax = Gr.dims.nlg[0]
Py_ssize_t jmax = Gr.dims.nlg[1]
Py_ssize_t kmax = Gr.dims.nlg[2]
Py_ssize_t istride = Gr.dims.nlg[1] * Gr.dims.nlg[2]
Py_ssize_t jstride = Gr.dims.nlg[2]
Py_ssize_t ishift, jshift
Py_ssize_t t_shift = DV.get_varshift(Gr, 'temperature')
Py_ssize_t tw_shift = DV.get_varshift(Gr, 'temperature_wb')
Py_ssize_t qv_shift = DV.get_varshift(Gr, 'qv')
Py_ssize_t ql_shift = DV.get_varshift(Gr,'ql')
Py_ssize_t nr_shift = PV.get_varshift(Gr, 'nr')
Py_ssize_t qr_shift = PV.get_varshift(Gr, 'qr')
Py_ssize_t qt_shift = PV.get_varshift(Gr, 'qt')
Py_ssize_t w_shift = PV.get_varshift(Gr, 'w')
double[:] qr_tendency = np.empty((Gr.dims.npg,), dtype=np.double, order='c')
double[:] nr_tendency = np.empty((Gr.dims.npg,), dtype=np.double, order='c')
double[:] tmp
double[:] dummy = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
Py_ssize_t wqr_shift = DV.get_varshift(Gr, 'w_qr')
Py_ssize_t wnr_shift = DV.get_varshift(Gr, 'w_nr')
Py_ssize_t wqt_shift
cdef double[:] s_src = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
if self.cloud_sedimentation:
wqt_shift = DV.get_varshift(Gr,'w_qt')
compute_advective_fluxes_a(&Gr.dims, &Ref.rho0[0], &Ref.rho0_half[0], &DV.values[wqt_shift], &DV.values[ql_shift], &dummy[0], 2, self.order)
tmp = Pa.HorizontalMean(Gr, &dummy[0])
NS.write_profile('qt_sedimentation_flux', tmp[gw:-gw], Pa)
compute_qt_sedimentation_s_source(&Gr.dims, &Ref.p0_half[0], &Ref.rho0_half[0], &dummy[0],
&PV.values[qt_shift], &DV.values[qv_shift],&DV.values[t_shift], &s_src[0], self.Lambda_fp,
self.L_fp, Gr.dims.dx[2], 2)
tmp = Pa.HorizontalMean(Gr, &s_src[0])
NS.write_profile('s_qt_sedimentation_source', tmp[gw:-gw], Pa)
#compute sedimentation flux only of nr
compute_advective_fluxes_a(&Gr.dims, &Ref.rho0[0], &Ref.rho0_half[0], &DV.values[wnr_shift], &PV.values[nr_shift], &dummy[0], 2, self.order)
tmp = Pa.HorizontalMean(Gr, &dummy[0])
NS.write_profile('nr_sedimentation_flux', tmp[gw:-gw], Pa)
#compute sedimentation flux only of qr
compute_advective_fluxes_a(&Gr.dims, &Ref.rho0[0], &Ref.rho0_half[0], &DV.values[wqr_shift], &PV.values[qr_shift], &dummy[0], 2, self.order)
tmp = Pa.HorizontalMean(Gr, &dummy[0])
NS.write_profile('qr_sedimentation_flux', tmp[gw:-gw], Pa)
#note we can re-use nr_tendency and qr_tendency because they are overwritten in each function
#must have a zero array to pass as entropy tendency and need to send a dummy variable for qt tendency
# Autoconversion tendencies of qr, nr, s
sb_autoconversion_rain_wrapper(&Gr.dims, self.compute_droplet_nu, &Ref.rho0_half[0], self.ccn,
&DV.values[ql_shift], &PV.values[qr_shift], &nr_tendency[0], &qr_tendency[0])
tmp = Pa.HorizontalMean(Gr, &nr_tendency[0])
NS.write_profile('nr_autoconversion', tmp[gw:-gw], Pa)
tmp = Pa.HorizontalMean(Gr, &qr_tendency[0])
NS.write_profile('qr_autoconversion', tmp[gw:-gw], Pa)
cdef double[:] s_auto = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
sb_entropy_source_formation(&Gr.dims, &self.CC.LT.LookupStructC, self.Lambda_fp, self.L_fp, &Ref.p0_half[0],
&DV.values[t_shift], &DV.values[tw_shift],&PV.values[qt_shift], &DV.values[qv_shift],
&qr_tendency[0], &s_auto[0])
tmp = Pa.HorizontalMean(Gr, &s_auto[0])
NS.write_profile('s_autoconversion', tmp[gw:-gw], Pa)
# Accretion tendencies of qr, s
sb_accretion_rain_wrapper(&Gr.dims, &Ref.rho0_half[0], &DV.values[ql_shift], &PV.values[qr_shift], &qr_tendency[0])
tmp = Pa.HorizontalMean(Gr, &qr_tendency[0])
NS.write_profile('qr_accretion', tmp[gw:-gw], Pa)
cdef double[:] s_accr = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
sb_entropy_source_formation(&Gr.dims, &self.CC.LT.LookupStructC, self.Lambda_fp, self.L_fp, &Ref.p0_half[0],
&DV.values[t_shift], &DV.values[tw_shift],&PV.values[qt_shift], &DV.values[qv_shift],
&qr_tendency[0], &s_accr[0])
tmp = Pa.HorizontalMean(Gr, &s_accr[0])
NS.write_profile('s_accretion', tmp[gw:-gw], Pa)
# Self-collection and breakup tendencies (lumped) of nr
sb_selfcollection_breakup_rain_wrapper(&Gr.dims, self.compute_rain_shape_parameter, &Ref.rho0_half[0],
&PV.values[nr_shift], &PV.values[qr_shift], &nr_tendency[0])
tmp = Pa.HorizontalMean(Gr, &nr_tendency[0])
NS.write_profile('nr_selfcollection', tmp[gw:-gw], Pa)
# Evaporation tendencies of qr, nr, s
sb_evaporation_rain_wrapper(&Gr.dims, &self.CC.LT.LookupStructC, self.Lambda_fp, self.L_fp,
self.compute_rain_shape_parameter, &Ref.rho0_half[0], &Ref.p0_half[0],
&DV.values[t_shift], &PV.values[qt_shift], &DV.values[ql_shift],
&PV.values[nr_shift], &PV.values[qr_shift], &nr_tendency[0], &qr_tendency[0])
tmp = Pa.HorizontalMean(Gr, &nr_tendency[0])
NS.write_profile('nr_evaporation', tmp[gw:-gw], Pa)
tmp = Pa.HorizontalMean(Gr, &qr_tendency[0])
NS.write_profile('qr_evaporation', tmp[gw:-gw], Pa)
cdef double[:] s_evp = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
sb_entropy_source_formation(&Gr.dims, &self.CC.LT.LookupStructC, self.Lambda_fp, self.L_fp, &Ref.p0_half[0],
&DV.values[t_shift], &DV.values[tw_shift],&PV.values[qt_shift], &DV.values[qv_shift],
&qr_tendency[0], &s_evp[0])
tmp = Pa.HorizontalMean(Gr, &s_evp[0])
NS.write_profile('s_evaporation', tmp[gw:-gw], Pa)
cdef double[:] s_heat = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
sb_entropy_source_heating(&Gr.dims, &DV.values[t_shift], &DV.values[tw_shift], &PV.values[qr_shift],
&DV.values[wqr_shift], &PV.values[w_shift], &s_heat[0])
tmp = Pa.HorizontalMean(Gr, &s_heat[0])
NS.write_profile('s_precip_heating', tmp[gw:-gw], Pa)
cdef double[:] s_drag = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
sb_entropy_source_drag(&Gr.dims, &DV.values[t_shift], &PV.values[qr_shift], &DV.values[wqr_shift], &s_drag[0])
tmp = Pa.HorizontalMean(Gr, &s_drag[0])
NS.write_profile('s_precip_drag', tmp[gw:-gw], Pa)
return
cdef extern from "microphysics_CLIMA.h":
void CLIMA_sedimentation_velocity_rain(Grid.DimStruct *dims,\
double* density, double* qr,\
double* qr_velocity) nogil
void CLIMA_microphysics_sources(Grid.DimStruct *dims, Lookup.LookupStruct *LT,\
double (*lam_fp)(double), double (*L_fp)(double, double),\
double* density, double* p0, double* temperature,\
double* qt, double* ql, double* qr, double dt,\
double* qr_tendency_micro, double* qr_tendency) nogil
void CLIMA_qt_source_formation(Grid.DimStruct *dims, double* qr_tendency,\
double* qt_tendency) nogil
void CLIMA_entropy_source_formation(Grid.DimStruct *dims, Lookup.LookupStruct *LT,\
double (*lam_fp)(double), double (*L_fp)(double, double),\
double* p0, double* T, double* Twet, double* qt,\
double* qv, double* qr_tendency,\
double* entropy_tendency) nogil
void CLIMA_entropy_source_heating(Grid.DimStruct *dims, double* T, double* Twet,\
double* qr, double* w_qr, double* w,\
double* entropy_tendency) nogil
void CLIMA_entropy_source_drag(Grid.DimStruct *dims, double* T,\
double* qr, double* w_qr,\
double* entropy_tendency) nogil
void CLIMA_autoconversion_rain_wrapper(Grid.DimStruct *dims,\
double* ql,\
double* qr_tendency) nogil
void CLIMA_accretion_rain_wrapper(Grid.DimStruct *dims, double* density,\
double* ql, double* qr,\
double* qr_tendency) nogil
void CLIMA_evaporation_rain_wrapper(Grid.DimStruct *dims, Lookup.LookupStruct *LT,\
double (*lam_fp)(double),\
double (*L_fp)(double, double),\
double* density, double* p0,\
double* temperature, double* qt,\
double* ql, double* qr,\
double* qr_tendency) nogil
cdef class Microphysics_CLIMA_Liquid_1M:
def __init__(self, ParallelMPI.ParallelMPI Par, LatentHeat LH, namelist):
# Create the appropriate linkages to the bulk thermodynamics
LH.Lambda_fp = lambda_constant
LH.L_fp = latent_heat_variable
self.thermodynamics_type = 'SA'
#also set local versions
self.Lambda_fp = lambda_constant
self.L_fp = latent_heat_variable
self.CC = ClausiusClapeyron()
self.CC.initialize(namelist, LH, Par)
try:
self.order = namelist['scalar_transport']['order_sedimentation']
except:
self.order = namelist['scalar_transport']['order']
return
cpdef initialize(self, Grid.Grid Gr, PrognosticVariables.PrognosticVariables PV,\
DiagnosticVariables.DiagnosticVariables DV, NetCDFIO_Stats NS,\
ParallelMPI.ParallelMPI Pa):
PV.add_variable('qr', 'kg/kg', r'q_r', 'rain water specific humidity','sym','scalar',Pa)
DV.add_variables('w_qr', 'm/s', r'w_{qr}', 'rain mass sedimentation veloctiy', 'sym', Pa)
DV.add_variables('temperature_wb', 'K', r'T_{wb}','wet bulb temperature','sym', Pa)
# add statistical output for the class
NS.add_profile('qr_sedimentation_flux', Gr, Pa)
NS.add_profile('qr_autoconversion', Gr, Pa)
NS.add_profile('qr_accretion', Gr, Pa)
NS.add_profile('qr_evaporation', Gr,Pa)
NS.add_profile('s_autoconversion', Gr, Pa)
NS.add_profile('s_accretion', Gr, Pa)
NS.add_profile('s_evaporation', Gr,Pa)
NS.add_profile('s_precip_heating', Gr, Pa)
NS.add_profile('s_precip_drag', Gr, Pa)
return
cpdef update(self, Grid.Grid Gr, ReferenceState.ReferenceState Ref, Th,\
PrognosticVariables.PrognosticVariables PV,\
DiagnosticVariables.DiagnosticVariables DV,\
TimeStepping.TimeStepping TS, ParallelMPI.ParallelMPI Pa):
cdef:
Py_ssize_t qt_shift = PV.get_varshift(Gr, 'qt')
Py_ssize_t qr_shift = PV.get_varshift(Gr, 'qr')
Py_ssize_t w_shift = PV.get_varshift(Gr, 'w')
Py_ssize_t t_shift = DV.get_varshift(Gr, 'temperature')
Py_ssize_t ql_shift = DV.get_varshift(Gr,'ql')
Py_ssize_t qv_shift = DV.get_varshift(Gr,'qv')
double dt = TS.dt
Py_ssize_t wqr_shift = DV.get_varshift(Gr, 'w_qr')
double[:] qr_tend_micro = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
CLIMA_microphysics_sources(&Gr.dims, &self.CC.LT.LookupStructC,\
self.Lambda_fp, self.L_fp,\
&Ref.rho0_half[0], &Ref.p0_half[0],\
&DV.values[t_shift], &PV.values[qt_shift],\
&DV.values[ql_shift], &PV.values[qr_shift],\
dt, &qr_tend_micro[0], &PV.tendencies[qr_shift])
CLIMA_sedimentation_velocity_rain(&Gr.dims, &Ref.rho0_half[0],\
&PV.values[qr_shift],\
&DV.values[wqr_shift])
CLIMA_qt_source_formation(&Gr.dims, &qr_tend_micro[0],\
&PV.tendencies[qt_shift])
cdef:
Py_ssize_t tw_shift = DV.get_varshift(Gr, 'temperature_wb')
Py_ssize_t s_shift = PV.get_varshift(Gr, 's')
microphysics_wetbulb_temperature(&Gr.dims, &self.CC.LT.LookupStructC,\
&Ref.p0_half[0], &PV.values[s_shift],\
&PV.values[qt_shift],\
&DV.values[t_shift], &DV.values[tw_shift])
CLIMA_entropy_source_formation(&Gr.dims, &self.CC.LT.LookupStructC,\
self.Lambda_fp, self.L_fp, &Ref.p0_half[0],\
&DV.values[t_shift], &DV.values[tw_shift],\
&PV.values[qt_shift], &DV.values[qv_shift],\
&qr_tend_micro[0], &PV.tendencies[s_shift])
CLIMA_entropy_source_heating(&Gr.dims, &DV.values[t_shift],\
&DV.values[tw_shift], &PV.values[qr_shift],\
&DV.values[wqr_shift], &PV.values[w_shift],\
&PV.tendencies[s_shift])
CLIMA_entropy_source_drag(&Gr.dims, &DV.values[t_shift],\
&PV.values[qr_shift], &DV.values[wqr_shift],\
&PV.tendencies[s_shift])
return
cpdef stats_io(self, Grid.Grid Gr, ReferenceState.ReferenceState Ref, Th,\
PrognosticVariables.PrognosticVariables PV,\
DiagnosticVariables.DiagnosticVariables DV,\
NetCDFIO_Stats NS, ParallelMPI.ParallelMPI Pa):
cdef:
Py_ssize_t i, j, k, ijk
Py_ssize_t gw = Gr.dims.gw
Py_ssize_t imax = Gr.dims.nlg[0]
Py_ssize_t jmax = Gr.dims.nlg[1]
Py_ssize_t kmax = Gr.dims.nlg[2]
Py_ssize_t istride = Gr.dims.nlg[1] * Gr.dims.nlg[2]
Py_ssize_t jstride = Gr.dims.nlg[2]
Py_ssize_t ishift, jshift
Py_ssize_t qr_shift = PV.get_varshift(Gr, 'qr')
Py_ssize_t qt_shift = PV.get_varshift(Gr, 'qt')
Py_ssize_t w_shift = PV.get_varshift(Gr, 'w')
Py_ssize_t t_shift = DV.get_varshift(Gr, 'temperature')
Py_ssize_t tw_shift = DV.get_varshift(Gr, 'temperature_wb')
Py_ssize_t qv_shift = DV.get_varshift(Gr, 'qv')
Py_ssize_t ql_shift = DV.get_varshift(Gr, 'ql')
Py_ssize_t wqr_shift = DV.get_varshift(Gr, 'w_qr')
double[:] qr_tendency = np.empty((Gr.dims.npg,), dtype=np.double, order='c')
double[:] tmp
double[:] dummy = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
double[:] s_src = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
# copied from SB microphysics
# sedimentation flux of qr
compute_advective_fluxes_a(&Gr.dims, &Ref.rho0[0], &Ref.rho0_half[0],\
&DV.values[wqr_shift], &PV.values[qr_shift],\
&dummy[0], 2, self.order)
tmp = Pa.HorizontalMean(Gr, &dummy[0])
NS.write_profile('qr_sedimentation_flux', tmp[gw:-gw], Pa)
# autoconversion tendencies of qr
CLIMA_autoconversion_rain_wrapper(&Gr.dims, &DV.values[ql_shift], &qr_tendency[0])
tmp = Pa.HorizontalMean(Gr, &qr_tendency[0])
NS.write_profile('qr_autoconversion', tmp[gw: -gw], Pa)
# autoconversion tendencies of s
cdef double[:] s_auto = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
CLIMA_entropy_source_formation(&Gr.dims, &self.CC.LT.LookupStructC,\
self.Lambda_fp, self.L_fp, &Ref.p0_half[0],\
&DV.values[t_shift], &DV.values[tw_shift],\
&PV.values[qt_shift], &DV.values[qv_shift],\
&qr_tendency[0], &s_auto[0])
tmp = Pa.HorizontalMean(Gr, &s_auto[0])
NS.write_profile('s_autoconversion', tmp[gw: -gw], Pa)
# accretion tendencies of qr
CLIMA_accretion_rain_wrapper(&Gr.dims, &Ref.rho0_half[0],\
&DV.values[ql_shift], &PV.values[qr_shift],\
&qr_tendency[0])
tmp = Pa.HorizontalMean(Gr, &qr_tendency[0])
NS.write_profile('qr_accretion', tmp[gw: -gw], Pa)
# accretion tendencies of s
cdef double[:] s_accr = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
CLIMA_entropy_source_formation(&Gr.dims, &self.CC.LT.LookupStructC,\
self.Lambda_fp, self.L_fp, &Ref.p0_half[0],\
&DV.values[t_shift], &DV.values[tw_shift],\
&PV.values[qt_shift], &DV.values[qv_shift],\
&qr_tendency[0], &s_accr[0])
tmp = Pa.HorizontalMean(Gr, &s_accr[0])
NS.write_profile('s_accretion', tmp[gw: -gw], Pa)
# evaporation tendencies of qr
CLIMA_evaporation_rain_wrapper(&Gr.dims, &self.CC.LT.LookupStructC,\
self.Lambda_fp, self.L_fp,\
&Ref.rho0_half[0], &Ref.p0_half[0],\
&DV.values[t_shift], &PV.values[qt_shift],\
&DV.values[ql_shift], &PV.values[qr_shift],\
&qr_tendency[0])
tmp = Pa.HorizontalMean(Gr, &qr_tendency[0])
NS.write_profile('qr_evaporation', tmp[gw: -gw], Pa)
# evaporation tendencies of s
cdef double[:] s_evp = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
CLIMA_entropy_source_formation(&Gr.dims, &self.CC.LT.LookupStructC,\
self.Lambda_fp, self.L_fp, &Ref.p0_half[0],\
&DV.values[t_shift], &DV.values[tw_shift],\
&PV.values[qt_shift], &DV.values[qv_shift],\
&qr_tendency[0], &s_evp[0])
tmp = Pa.HorizontalMean(Gr, &s_evp[0])
NS.write_profile('s_evaporation', tmp[gw: -gw], Pa)
# precip heating tendencies of s
cdef double[:] s_heat = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
CLIMA_entropy_source_heating(&Gr.dims, &DV.values[t_shift],\
&DV.values[tw_shift], &PV.values[qr_shift],\
&DV.values[wqr_shift], &PV.values[w_shift],\
&s_heat[0])
tmp = Pa.HorizontalMean(Gr, &s_heat[0])
NS.write_profile('s_precip_heating', tmp[gw: -gw], Pa)
# precip drag tendencies of s
cdef double[:] s_drag = np.zeros((Gr.dims.npg,), dtype=np.double, order='c')
CLIMA_entropy_source_drag(&Gr.dims, &DV.values[t_shift],\
&PV.values[qr_shift], &DV.values[wqr_shift],\
&s_drag[0])
tmp = Pa.HorizontalMean(Gr, &s_drag[0])
NS.write_profile('s_precip_drag', tmp[gw: -gw], Pa)
return
cdef extern from "entropies.h":
inline double sd_c(double pd, double T) nogil
inline double sv_c(double pv, double T) nogil
cdef extern from "thermodynamic_functions.h":
inline double qv_star_c(const double p0, const double qt, const double pv)nogil
cdef cython_wetbulb(Grid.DimStruct *dims, Lookup.LookupStruct *LT, double *p0, double *s, double *qt, double *T, double *Twet):
cdef:
Py_ssize_t imin = 0
Py_ssize_t jmin = 0
Py_ssize_t kmin = 0
Py_ssize_t imax = dims.nlg[0]
Py_ssize_t jmax = dims.nlg[1]
Py_ssize_t kmax = dims.nlg[2]
Py_ssize_t istride = dims.nlg[1] * dims.nlg[2]
Py_ssize_t jstride = dims.nlg[2]
Py_ssize_t ishift, jshift, ijk, i,j,k, iter = 0
cdef:
double T_1, T_2, T_n, pv_star_1, pv_star_2, qv_star_1, qv_star_2
double pd_1, pd_2, s_1, s_2, f_1, f_2, delta_T
cdef Py
with nogil:
for i in xrange(imin,imax):
ishift = i*istride
for j in xrange(jmin,jmax):
jshift = j*jstride
for k in xrange(kmin,kmax):
ijk = ishift + jshift + k
T_1 = T[ijk]
pv_star_1 = Lookup.lookup(LT, T_1)
qv_star_1 = qv_star_c(p0[k], qt[ijk], pv_star_1)
if qt[ijk] >= qv_star_1:
Twet[ijk] = T_1
else:
T_2 = T_1 + 1.0
delta_T = fabs(T_2 - T_1)
qv_star_1 = pv_star_1/(eps_vi * (p0[k] - pv_star_1) + pv_star_1)
pd_1 = p0[k] - pv_star_1
s_1 = sd_c(pd_1,T_1) * (1.0 - qv_star_1) + sv_c(pv_star_1,T_1) * qv_star_1
f_1 = s[ijk] - s_1
iter = 0
while delta_T > 1.0e-3:
pv_star_2 = Lookup.lookup(LT, T_2)
qv_star_2 = pv_star_2/(eps_vi * (p0[k] - pv_star_2) + pv_star_2)
pd_2 = p0[k] - pv_star_2
s_2 = sd_c(pd_2,T_2) * (1.0 - qv_star_2) + sv_c(pv_star_2,T_2) * qv_star_2
f_2 = s[ijk] - s_2
T_n = T_2 - f_2*(T_2 - T_1)/(f_2 - f_1)
T_1 = T_2
T_2 = T_n
f_1 = f_2
delta_T = fabs(T_2 - T_1)
iter += 1
Twet[ijk] = T_2
with gil:
print(T[ijk]-Twet[ijk], iter)
print('leaving wetbulb')
return
cdef class Microphysics_T_Liquid:
def __init__(self, ParallelMPI.ParallelMPI Par, LatentHeat LH, namelist):
LH.Lambda_fp = lambda_T
LH.L_fp = latent_heat_T
self.thermodynamics_type = 'SA'
#also set local versions
self.Lambda_fp = lambda_T
self.L_fp = latent_heat_T
self.CC = ClausiusClapeyron()
self.CC.initialize(namelist, LH, Par)
return
cpdef initialize(self, Grid.Grid Gr, PrognosticVariables.PrognosticVariables PV,
DiagnosticVariables.DiagnosticVariables DV, NetCDFIO_Stats NS, ParallelMPI.ParallelMPI Pa):
DV.add_variables('dqtdt_precip', 'kgkg^-1s^-1',
r'\left(\frac{dq_t}{d_t}\right)_{auto}', 'qt autoconversion tendency', 'sym', Pa)
DV.add_variables('dsdt_precip', 'J kg^-1 K^-1s^-1',
r'\left(\frac{ds}{d_t}\right)_{auto}', 's autoconversion tendency', 'sym', Pa)
return
cpdef update(self, Grid.Grid Gr, ReferenceState.ReferenceState Ref, Th, PrognosticVariables.PrognosticVariables PV,
DiagnosticVariables.DiagnosticVariables DV, TimeStepping.TimeStepping TS, ParallelMPI.ParallelMPI Pa):
cdef:
Py_ssize_t i, j, k, ijk, ishift, jshift
Py_ssize_t istride = Gr.dims.nlg[1] * Gr.dims.nlg[2]
Py_ssize_t jstride = Gr.dims.nlg[2]
Py_ssize_t imin = 0
Py_ssize_t jmin = 0
Py_ssize_t kmin = 0
Py_ssize_t imax = Gr.dims.nlg[0]
Py_ssize_t jmax = Gr.dims.nlg[1]
Py_ssize_t kmax = Gr.dims.nlg[2]
Py_ssize_t count
Py_ssize_t s_shift = PV.get_varshift(Gr, 's')
Py_ssize_t qt_shift = PV.get_varshift(Gr, 'qt')
Py_ssize_t qv_shift = DV.get_varshift(Gr,'qv')
Py_ssize_t ql_shift = DV.get_varshift(Gr,'ql')
Py_ssize_t qi_shift = DV.get_varshift(Gr,'qi')
Py_ssize_t t_shift = DV.get_varshift(Gr,'temperature')
Py_ssize_t dqtdt_shift = DV.get_varshift(Gr,'dqtdt_precip')
Py_ssize_t dsdt_shift = DV.get_varshift(Gr,'dsdt_precip')
double lam, t, p0, rho0, qt, qv, pd, pv
double lv
double qco = 1e-3
double alpha_kessler = 0.001
double dqtdt_rain_auto
double dqtdt_snow_auto
double tau_a
# Ouput profiles of relative humidity
with nogil:
count = 0
for i in range(imin, imax):
ishift = i * istride
for j in range(jmin, jmax):
jshift = j * jstride
for k in range(kmin, kmax):
ijk = ishift + jshift + k
dqtdt_rain_auto = 0.0
dqtdt_snow_auto = 0.0
#Zero the diagnotic tendencies
DV.values[dsdt_shift + ijk] = 0.0
DV.values[dqtdt_shift + ijk] = 0.0
if DV.values[ql_shift + ijk] > 0.0:
p0 = Ref.p0_half[k]
rho0 = Ref.rho0_half[k]
qt = PV.values[qt_shift + ijk]
qv = qt - DV.values[ql_shift + ijk]
pd = pd_c(p0,qt,qv)
pv = pv_c(p0,qt,qv)
t = DV.values[t_shift + ijk]
lam = self.Lambda_fp(t)
lv = self.L_fp(t, lam)
tau_a = 32.0/9.0*pow((t-258.15),2.0) + 200.0;
if tau_a > 1000.0:
tau_a = 1000.0
dqtdt_rain_auto = -fmax(0.0, alpha_kessler * (DV.values[ql_shift + ijk] - qco))
dqtdt_snow_auto = - rho0 * DV.values[qi_shift + ijk]/ tau_a
DV.values[dqtdt_shift + ijk] = dqtdt_rain_auto + dqtdt_snow_auto #-fmax(0.0, DV.values[ql_shift + ijk] - 0.02 * DV.values[qv_shift+ijk])/TS.dt
DV.values[dsdt_shift + ijk] = -(sd_c(pd,t) - sv_c(pv,t) + lv/t ) * DV.values[dqtdt_shift + ijk]
PV.tendencies[qt_shift + ijk] += DV.values[dqtdt_shift + ijk]
PV.tendencies[s_shift + ijk] += DV.values[dsdt_shift + ijk]
return
cpdef stats_io(self, Grid.Grid Gr, ReferenceState.ReferenceState Ref, Th, PrognosticVariables.PrognosticVariables PV,
DiagnosticVariables.DiagnosticVariables DV, NetCDFIO_Stats NS, ParallelMPI.ParallelMPI Pa):
return
def MicrophysicsFactory(namelist, LatentHeat LH, ParallelMPI.ParallelMPI Par):
if(namelist['microphysics']['scheme'] == 'None_Dry'):
return No_Microphysics_Dry(Par, LH, namelist)
elif(namelist['microphysics']['scheme'] == 'None_SA'):
return No_Microphysics_SA(Par, LH, namelist)
elif(namelist['microphysics']['scheme'] == 'SB_Liquid'):
return Microphysics_SB_Liquid(Par, LH, namelist)
elif(namelist['microphysics']['scheme'] == 'T_Liquid'):
return Microphysics_T_Liquid(Par, LH, namelist)
elif(namelist['microphysics']['scheme'] == 'Arctic_1M'):
return Microphysics_Arctic_1M(Par, LH, namelist)
elif(namelist['microphysics']['scheme'] == 'CLIMA_liquid_1M'):
return Microphysics_CLIMA_Liquid_1M(Par, LH, namelist)