-
Notifications
You must be signed in to change notification settings - Fork 151
/
Copy pathcu_gf_deep.F90
5796 lines (5436 loc) · 209 KB
/
cu_gf_deep.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 cu_gf_deep.F90
!! This file is the Grell-Freitas deep convection scheme.
module cu_gf_deep
use machine , only : kind_phys
use physcons, only : qamin
real(kind=kind_phys), parameter::g=9.81
real(kind=kind_phys), parameter:: cp=1004.
real(kind=kind_phys), parameter:: xlv=2.5e6
real(kind=kind_phys), parameter::r_v=461.
real(kind=kind_phys), parameter :: tcrit=258.
!> tuning constant for cloudwater/ice detrainment
real(kind=kind_phys), parameter:: c1= 0.003 !.002 ! .0005
!> parameter to turn on or off evaporation of rainwater as done in sas
integer, parameter :: irainevap=1
!> max allowed fractional coverage (frh_thresh)
real(kind=kind_phys), parameter::frh_thresh = .9
!> rh threshold. if fractional coverage ~ frh_thres, do not use cupa any further
real(kind=kind_phys), parameter::rh_thresh = .97
!> tuning constant for j. brown closure (ichoice = 4,5,6)
real(kind=kind_phys), parameter::betajb=1.2
!> tuning for shallow and mid convection. ec uses 1.5
integer, parameter:: use_excess=0
real(kind=kind_phys), parameter :: fluxtune=1.5
!> flag to turn off or modify mom transport by downdrafts
real(kind=kind_phys), parameter :: pgcd = 0.1
!
!> aerosol awareness, do not use yet!
integer, parameter :: autoconv=1 !2
integer, parameter :: aeroevap=1 !3
real(kind=kind_phys), parameter :: scav_factor = 0.5
real(kind=kind_phys), parameter :: dx_thresh = 6500.
!> still 16 ensembles for clousres
integer, parameter:: maxens3=16
!---meltglac-------------------------------------------------
logical, parameter :: melt_glac = .true. !<- turn on/off ice phase/melting
real(kind=kind_phys), parameter :: &
t_0 = 273.16, & !< k
t_ice = 250.16, & !< k
xlf = 0.333e6 !< latent heat of freezing (j k-1 kg-1)
!---meltglac-------------------------------------------------
!-----srf-08aug2017-----begin
real(kind=kind_phys), parameter :: qrc_crit= 2.e-4
!-----srf-08aug2017-----end
contains
!>\defgroup cu_gf_deep_group Grell-Freitas Deep Convection Module
!>\ingroup cu_gf_group
!! This is Grell-Freitas deep convection scheme module
!> @{
integer function my_maxloc1d(A,N)
!$acc routine vector
implicit none
real(kind_phys), intent(in) :: A(:)
integer, intent(in) :: N
real(kind_phys) :: imaxval
integer :: i
imaxval = MAXVAL(A)
my_maxloc1d = 1
!$acc loop
do i = 1, N
if ( A(i) == imaxval ) then
my_maxloc1d = i
return
endif
end do
return
end function my_maxloc1d
!>Driver for the deep or congestus GF routine.
!! \section general_gf_deep Grell-Freitas Deep Convection General Algorithm
subroutine cu_gf_deep_run( &
itf,ktf,its,ite, kts,kte &
,dicycle & ! diurnal cycle flag
,ichoice & ! choice of closure, use "0" for ensemble average
,ipr & ! this flag can be used for debugging prints
,ccn & ! not well tested yet
,ccnclean &
,dtime & ! dt over which forcing is applied
,imid & ! flag to turn on mid level convection
,kpbl & ! level of boundary layer height
,dhdt & ! boundary layer forcing (one closure for shallow)
,xland & ! land mask
,zo & ! heights above surface
,forcing & ! only diagnostic
,t & ! t before forcing
,q & ! q before forcing
,z1 & ! terrain
,tn & ! t including forcing
,qo & ! q including forcing
,po & ! pressure (mb)
,psur & ! surface pressure (mb)
,us & ! u on mass points
,vs & ! v on mass points
,rho & ! density
,hfx & ! w/m2, positive upward
,qfx & ! w/m2, positive upward
,dx & ! dx is grid point dependent here
,mconv & ! integrated vertical advection of moisture
,omeg & ! omega (pa/s)
,csum & ! used to implement memory, set to zero if not avail
,cnvwt & ! gfs needs this
,zuo & ! nomalized updraft mass flux
,zdo & ! nomalized downdraft mass flux
,zdm & ! nomalized downdraft mass flux from mid scheme
,edto & !
,edtm & !
,xmb_out & ! the xmb's may be needed for dicycle
,xmbm_in & !
,xmbs_in & !
,pre & !
,outu & ! momentum tendencies at mass points
,outv & !
,outt & ! temperature tendencies
,outq & ! q tendencies
,outqc & ! ql/qice tendencies
,kbcon & ! lfc of parcel from k22
,ktop & ! cloud top
,cupclw & ! used for direct coupling to radiation, but with tuning factors
,frh_out & ! fractional coverage
,ierr & ! ierr flags are error flags, used for debugging
,ierrc & ! the following should be set to zero if not available
,nchem &
,fscav &
,chem3d &
,wetdpc_deep &
,do_smoke_transport &
,rand_mom & ! for stochastics mom, if temporal and spatial patterns exist
,rand_vmas & ! for stochastics vertmass, if temporal and spatial patterns exist
,rand_clos & ! for stochastics closures, if temporal and spatial patterns exist
,nranflag & ! flag to what you want perturbed
!! 1 = momentum transport
!! 2 = normalized vertical mass flux profile
!! 3 = closures
!! more is possible, talk to developer or
!! implement yourself. pattern is expected to be
!! betwee -1 and +1
,do_capsuppress,cap_suppress_j & !
,k22 & !
,jmin,kdt,tropics) !
implicit none
integer &
,intent (in ) :: &
nranflag,itf,ktf,its,ite, kts,kte,ipr,imid,kdt
integer, intent (in ) :: &
ichoice,nchem
real(kind=kind_phys), dimension (its:ite,4) &
,intent (in ) :: rand_clos
real(kind=kind_phys), dimension (its:ite) &
,intent (in ) :: rand_mom,rand_vmas
!$acc declare copyin(rand_clos,rand_mom,rand_vmas)
integer, intent(in) :: do_capsuppress
real(kind=kind_phys), intent(in), dimension(:) :: cap_suppress_j
!$acc declare create(cap_suppress_j)
!
!
!
real(kind=kind_phys), dimension (its:ite,1:maxens3) :: xf_ens,pr_ens
!$acc declare create(xf_ens,pr_ens)
! outtem = output temp tendency (per s)
! outq = output q tendency (per s)
! outqc = output qc tendency (per s)
! pre = output precip
real(kind=kind_phys), dimension (its:ite,kts:kte) &
,intent (inout ) :: &
cnvwt,outu,outv,outt,outq,outqc,cupclw
real(kind=kind_phys), dimension (its:ite) &
,intent (out ) :: &
frh_out
real(kind=kind_phys), dimension (its:ite) &
,intent (inout ) :: &
pre,xmb_out
!$acc declare copy(cnvwt,outu,outv,outt,outq,outqc,cupclw,frh_out,pre,xmb_out)
real(kind=kind_phys), dimension (its:ite) &
,intent (in ) :: &
hfx,qfx,xmbm_in,xmbs_in
!$acc declare copyin(hfx,qfx,xmbm_in,xmbs_in)
integer, dimension (its:ite) &
,intent (inout ) :: &
kbcon,ktop
!$acc declare copy(kbcon,ktop)
integer, dimension (its:ite) &
,intent (in ) :: &
kpbl,tropics
!$acc declare copyin(kpbl,tropics)
!
! basic environmental input includes moisture convergence (mconv)
! omega (omeg), windspeed (us,vs), and a flag (ierr) to turn off
! convection for this call only and at that particular gridpoint
!
real(kind=kind_phys), dimension (its:ite,kts:kte) &
,intent (in ) :: &
dhdt,rho,t,po,us,vs,tn
!$acc declare copyin(dhdt,rho,t,po,us,vs,tn)
real(kind=kind_phys), dimension (its:ite,kts:kte) &
,intent (inout ) :: &
omeg
!$acc declare copy(omeg)
real(kind=kind_phys), dimension (its:ite,kts:kte) &
,intent (inout) :: &
q,qo,zuo,zdo,zdm
!$acc declare copy(q,qo,zuo,zdo,zdm)
real(kind=kind_phys), dimension (its:ite) &
,intent (in ) :: &
dx,z1,psur,xland
!$acc declare copyin(dx,z1,psur,xland)
real(kind=kind_phys), dimension (its:ite) &
,intent (inout ) :: &
mconv,ccn
!$acc declare copy(mconv,ccn)
real(kind=kind_phys), dimension (:,:,:) &
,intent (inout) :: &
chem3d
logical, intent (in) :: do_smoke_transport
real(kind=kind_phys), dimension (:,:) &
, intent (out) :: wetdpc_deep
real(kind=kind_phys), intent (in) :: fscav(:)
!$acc declare copy(chem3d) copyout(wetdpc_deep) copyin(fscav)
real(kind=kind_phys) &
,intent (in ) :: &
dtime,ccnclean
!
! local ensemble dependent variables in this routine
!
real(kind=kind_phys), dimension (its:ite,1) :: &
xaa0_ens
real(kind=kind_phys), dimension (its:ite,1) :: &
edtc
real(kind=kind_phys), dimension (its:ite,kts:kte,1) :: &
dellat_ens,dellaqc_ens,dellaq_ens,pwo_ens
!$acc declare create(xaa0_ens,edtc,dellat_ens,dellaqc_ens,dellaq_ens,pwo_ens)
!
!
!
!***************** the following are your basic environmental
! variables. they carry a "_cup" if they are
! on model cloud levels (staggered). they carry
! an "o"-ending (z becomes zo), if they are the forced
! variables. they are preceded by x (z becomes xz)
! to indicate modification by some typ of cloud
!
! z = heights of model levels
! q = environmental mixing ratio
! qes = environmental saturation mixing ratio
! t = environmental temp
! p = environmental pressure
! he = environmental moist static energy
! hes = environmental saturation moist static energy
! z_cup = heights of model cloud levels
! q_cup = environmental q on model cloud levels
! qes_cup = saturation q on model cloud levels
! t_cup = temperature (kelvin) on model cloud levels
! p_cup = environmental pressure
! he_cup = moist static energy on model cloud levels
! hes_cup = saturation moist static energy on model cloud levels
! gamma_cup = gamma on model cloud levels
!
!
! hcd = moist static energy in downdraft
! zd normalized downdraft mass flux
! dby = buoancy term
! entr = entrainment rate
! zd = downdraft normalized mass flux
! entr= entrainment rate
! hcd = h in model cloud
! bu = buoancy term
! zd = normalized downdraft mass flux
! gamma_cup = gamma on model cloud levels
! qcd = cloud q (including liquid water) after entrainment
! qrch = saturation q in cloud
! pwd = evaporate at that level
! pwev = total normalized integrated evaoprate (i2)
! entr= entrainment rate
! z1 = terrain elevation
! entr = downdraft entrainment rate
! jmin = downdraft originating level
! kdet = level above ground where downdraft start detraining
! psur = surface pressure
! z1 = terrain elevation
! pr_ens = precipitation ensemble
! xf_ens = mass flux ensembles
! omeg = omega from large scale model
! mconv = moisture convergence from large scale model
! zd = downdraft normalized mass flux
! zu = updraft normalized mass flux
! dir = "storm motion"
! mbdt = arbitrary numerical parameter
! dtime = dt over which forcing is applied
! kbcon = lfc of parcel from k22
! k22 = updraft originating level
! ichoice = flag if only want one closure (usually set to zero!)
! dby = buoancy term
! ktop = cloud top (output)
! xmb = total base mass flux
! hc = cloud moist static energy
! hkb = moist static energy at originating level
real(kind=kind_phys), dimension (its:ite,kts:kte,nchem) :: &
chem
real(kind=kind_phys), dimension (its:ite,kts:kte,nchem) :: &
chem_cup,chem_up,chem_down,dellac,dellac2,chem_c,chem_pw,chem_pwd
real(kind=kind_phys), dimension (its:ite,nchem) :: &
chem_pwav,chem_psum
real(kind=kind_phys):: dtime_max,sum1,sum2
real(kind=kind_phys), dimension (kts:kte) :: trac,trcflx_in,trcflx_out,trc,trco
real(kind=kind_phys), dimension (its:ite,kts:kte) :: pwdper, massflx
integer :: nv
!$acc declare create(chem,chem_cup,chem_up,chem_down,dellac,dellac2,chem_c,chem_pw,chem_pwd, &
!$acc chem_pwav,chem_psum,pwdper,massflux)
real(kind=kind_phys), dimension (its:ite,kts:kte) :: &
entr_rate_2d,mentrd_rate_2d,he,hes,qes,z, heo,heso,qeso,zo, &
xhe,xhes,xqes,xz,xt,xq,qes_cup,q_cup,he_cup,hes_cup,z_cup, &
p_cup,gamma_cup,t_cup, qeso_cup,qo_cup,heo_cup,heso_cup, &
zo_cup,po_cup,gammao_cup,tn_cup, &
xqes_cup,xq_cup,xhe_cup,xhes_cup,xz_cup, &
xt_cup, dby,hc,zu,clw_all, &
dbyo,qco,qrcdo,pwdo,pwo,hcdo,qcdo,dbydo,hco,qrco, &
dbyt,xdby,xhc,xzu, &
! cd = detrainment function for updraft
! cdd = detrainment function for downdraft
! dellat = change of temperature per unit mass flux of cloud ensemble
! dellaq = change of q per unit mass flux of cloud ensemble
! dellaqc = change of qc per unit mass flux of cloud ensemble
cd,cdd,dellah,dellaq,dellat,dellaqc, &
u_cup,v_cup,uc,vc,ucd,vcd,dellu,dellv
!$acc declare create( &
!$acc entr_rate_2d,mentrd_rate_2d,he,hes,qes,z, heo,heso,qeso,zo, &
!$acc xhe,xhes,xqes,xz,xt,xq,qes_cup,q_cup,he_cup,hes_cup,z_cup, &
!$acc p_cup,gamma_cup,t_cup, qeso_cup,qo_cup,heo_cup,heso_cup, &
!$acc zo_cup,po_cup,gammao_cup,tn_cup, &
!$acc xqes_cup,xq_cup,xhe_cup,xhes_cup,xz_cup, &
!$acc xt_cup, dby,hc,zu,clw_all, &
!$acc dbyo,qco,qrcdo,pwdo,pwo,hcdo,qcdo,dbydo,hco,qrco, &
!$acc dbyt,xdby,xhc,xzu, &
!$acc cd,cdd,dellah,dellaq,dellat,dellaqc, &
!$acc u_cup,v_cup,uc,vc,ucd,vcd,dellu,dellv)
! aa0 cloud work function for downdraft
! edt = epsilon
! aa0 = cloud work function without forcing effects
! aa1 = cloud work function with forcing effects
! xaa0 = cloud work function with cloud effects (ensemble dependent)
! edt = epsilon
real(kind=kind_phys), dimension (its:ite) :: &
edt,edto,edtm,aa1,aa0,xaa0,hkb, &
hkbo,xhkb, &
xmb,pwavo,ccnloss, &
pwevo,bu,bud,cap_max, &
cap_max_increment,closure_n,psum,psumh,sig,sigd
real(kind=kind_phys), dimension (its:ite) :: &
axx,edtmax,edtmin,entr_rate
integer, dimension (its:ite) :: &
kzdown,kdet,k22,jmin,kstabi,kstabm,k22x,xland1, &
ktopdby,kbconx,ierr2,ierr3,kbmax
!$acc declare create(edt,edto,edtm,aa1,aa0,xaa0,hkb, &
!$acc hkbo,xhkb, &
!$acc xmb,pwavo,ccnloss, &
!$acc pwevo,bu,bud,cap_max, &
!$acc cap_max_increment,closure_n,psum,psumh,sig,sigd, &
!$acc axx,edtmax,edtmin,entr_rate, &
!$acc kzdown,kdet,k22,jmin,kstabi,kstabm,k22x,xland1, &
!$acc ktopdby,kbconx,ierr2,ierr3,kbmax)
integer, dimension (its:ite), intent(inout) :: ierr
integer, dimension (its:ite), intent(in) :: csum
!$acc declare copy(ierr) copyin(csum)
integer :: &
iloop,nens3,ki,kk,i,k
real(kind=kind_phys) :: &
dz,dzo,mbdt,radius, &
zcutdown,depth_min,zkbmax,z_detr,zktop, &
dh,cap_maxs,trash,trash2,frh,sig_thresh
real(kind=kind_phys), dimension (its:ite) :: pefc
real(kind=kind_phys) entdo,dp,subin,detdo,entup, &
detup,subdown,entdoj,entupk,detupk,totmas
real(kind=kind_phys), dimension (its:ite) :: lambau,flux_tun,zws,ztexec,zqexec
!$acc declare create(lambau,flux_tun,zws,ztexec,zqexec)
integer :: jprnt,jmini,start_k22
logical :: keep_going,flg(its:ite)
!$acc declare create(flg)
character*50 :: ierrc(its:ite)
character*4 :: cumulus
real(kind=kind_phys), dimension (its:ite,kts:kte) :: &
up_massentr,up_massdetr,c1d &
,up_massentro,up_massdetro,dd_massentro,dd_massdetro
real(kind=kind_phys), dimension (its:ite,kts:kte) :: &
up_massentru,up_massdetru,dd_massentru,dd_massdetru
!$acc declare create(up_massentr,up_massdetr,c1d,up_massentro,up_massdetro,dd_massentro,dd_massdetro, &
!$acc up_massentru,up_massdetru,dd_massentru,dd_massdetru)
real(kind=kind_phys) c1_max,buo_flux,pgcon,pgc,blqe
real(kind=kind_phys) :: xff_mid(its:ite,2)
!$acc declare create(xff_mid)
integer :: iversion=1
real(kind=kind_phys) :: denom,h_entr,umean,t_star,dq
integer, intent(in) :: dicycle
real(kind=kind_phys), dimension (its:ite) :: aa1_bl,hkbo_bl,tau_bl,tau_ecmwf,wmean
real(kind=kind_phys), dimension (its:ite,kts:kte) :: tn_bl, qo_bl, qeso_bl, heo_bl, heso_bl &
,qeso_cup_bl,qo_cup_bl, heo_cup_bl,heso_cup_bl &
,gammao_cup_bl,tn_cup_bl,hco_bl,dbyo_bl
real(kind=kind_phys), dimension(its:ite) :: xf_dicycle
!$acc declare create(aa1_bl,hkbo_bl,tau_bl,tau_ecmwf,wmean, &
!$acc tn_bl, qo_bl, qeso_bl, heo_bl, heso_bl, &
!$acc qeso_cup_bl,qo_cup_bl, heo_cup_bl,heso_cup_bl, &
!$acc gammao_cup_bl,tn_cup_bl,hco_bl,dbyo_bl,xf_dicycle)
real(kind=kind_phys), intent(inout), dimension(its:ite,10) :: forcing
!$acc declare copy(forcing)
integer :: turn,pmin_lev(its:ite),start_level(its:ite),ktopkeep(its:ite)
real(kind=kind_phys), dimension (its:ite,kts:kte) :: dtempdz
integer, dimension (its:ite,kts:kte) :: k_inv_layers
real(kind=kind_phys), dimension (its:ite) :: c0 ! HCB
real(kind=kind_phys), dimension (its:ite,kts:kte) :: c0t3d ! hli for smoke/dust wet scavenging
!$acc declare create(pmin_lev,start_level,ktopkeep,dtempdz,k_inv_layers,c0,c0t3d)
! rainevap from sas
real(kind=kind_phys) zuh2(40)
real(kind=kind_phys), dimension (its:ite) :: rntot,delqev,delq2,qevap,rn,qcond
!$acc declare create(zuh2,rntot,delqev,delq2,qevap,rn,qcond)
real(kind=kind_phys) :: rain,t1,q1,elocp,evef,el2orc,evfact,evfactl,g_rain,e_dn,c_up
real(kind=kind_phys) :: pgeoh,dts,fp,fpi,pmin,x_add,beta,beta_u
real(kind=kind_phys) :: cbeg,cmid,cend,const_a,const_b,const_c
!---meltglac-------------------------------------------------
real(kind=kind_phys), dimension (its:ite,kts:kte) :: p_liq_ice,melting_layer,melting
!$acc declare create(p_liq_ice,melting_layer,melting)
integer :: itemp
!---meltglac-------------------------------------------------
!$acc kernels
melting_layer(:,:)=0.
melting(:,:)=0.
flux_tun(:)=fluxtune
!$acc end kernels
! if(imid.eq.1)flux_tun(:)=fluxtune+.5
cumulus='deep'
if(imid.eq.1)cumulus='mid'
pmin=150.
if(imid.eq.1)pmin=75.
!$acc kernels
ktopdby(:)=0
!$acc end kernels
c1_max=c1
elocp=xlv/cp
el2orc=xlv*xlv/(r_v*cp)
evfact=0.25 ! .4
evfactl=0.25 ! .2
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
!proportionality constant to estimate pressure gradient of updraft (zhang and wu, 2003, jas
!
! ecmwf
pgcon=0.
!$acc kernels
lambau(:)=2.0
if(imid.eq.1)lambau(:)=2.0
! here random must be between -1 and 1
if(nranflag == 1)then
lambau(:)=1.5+rand_mom(:)
endif
!$acc end kernels
! sas
! lambau=0.
! pgcon=-.55
!
!---------------------------------------------------- ! HCB
! Set cloud water to rain water conversion rate (c0)
!$acc kernels
c0(:)=0.004
do i=its,itf
xland1(i)=int(xland(i)+.0001) ! 1.
if(xland(i).gt.1.5 .or. xland(i).lt.0.5)then
xland1(i)=0
endif
if(xland1(i).eq.1)c0(i)=0.002
if(imid.eq.1)then
c0(i)=0.002
endif
enddo
!$acc end kernels
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!$acc kernels
ztexec(:) = 0.
zqexec(:) = 0.
zws(:) = 0.
do i=its,itf
!- buoyancy flux (h+le)
buo_flux= (hfx(i)/cp+0.608*t(i,1)*qfx(i)/xlv)/rho(i,1)
pgeoh = zo(i,2)*g
!-convective-scale velocity w*
zws(i) = max(0.,flux_tun(i)*0.41*buo_flux*zo(i,2)*g/t(i,1))
if(zws(i) > tiny(pgeoh)) then
!-convective-scale velocity w*
zws(i) = 1.2*zws(i)**.3333
!- temperature excess
ztexec(i) = max(flux_tun(i)*hfx(i)/(rho(i,1)*zws(i)*cp),0.0)
!- moisture excess
zqexec(i) = max(flux_tun(i)*qfx(i)/xlv/(rho(i,1)*zws(i)),0.)
endif
!- zws for shallow convection closure (grant 2001)
!- height of the pbl
zws(i) = max(0.,.001-flux_tun(i)*0.41*buo_flux*zo(i,kpbl(i))*g/t(i,kpbl(i)))
zws(i) = 1.2*zws(i)**.3333
zws(i) = zws(i)*rho(i,kpbl(i)) !check if zrho is correct
enddo
!$acc end kernels
cap_maxs=75. ! 150.
!$acc kernels
do i=its,itf
edto(i)=0.
closure_n(i)=16.
xmb_out(i)=0.
cap_max(i)=cap_maxs
cap_max_increment(i)=20.
!
! for water or ice
!
if (xland1(i)==0) then
cap_max_increment(i)=20.
else
if(ztexec(i).gt.0.)cap_max(i)=cap_max(i)+25.
if(ztexec(i).lt.0.)cap_max(i)=cap_max(i)-25.
endif
#ifndef _OPENACC
ierrc(i)=" "
#endif
enddo
!$acc end kernels
if(use_excess == 0 )then
!$acc kernels
ztexec(:)=0
zqexec(:)=0
!$acc end kernels
endif
if(do_capsuppress == 1) then
!$acc kernels
do i=its,itf
cap_max(i)=cap_maxs
if (abs(cap_suppress_j(i) - 1.0 ) < 0.1 ) then
cap_max(i)=cap_maxs+75.
elseif (abs(cap_suppress_j(i) - 0.0 ) < 0.1 ) then
cap_max(i)=10.0
endif
enddo
!$acc end kernels
endif
!
!--- initial entrainment rate (these may be changed later on in the
!--- program
!
!$acc kernels
start_level(:)=kte
!$acc end kernels
!$acc kernels
!$acc loop private(radius,frh)
do i=its,ite
c1d(i,:)= 0. !c1 ! 0. ! c1 ! max(.003,c1+float(csum(i))*.0001)
entr_rate(i)=7.e-5 - min(20.,float(csum(i))) * 3.e-6
if(xland1(i) == 0)entr_rate(i)=7.e-5
if(dx(i)<dx_thresh) entr_rate(i)=2.e-4
if(imid.eq.1)entr_rate(i)=3.e-4
radius=.2/entr_rate(i)
frh=min(1.,3.14*radius*radius/dx(i)/dx(i))
if(frh > frh_thresh)then
frh=frh_thresh
radius=sqrt(frh*dx(i)*dx(i)/3.14)
entr_rate(i)=.2/radius
endif
sig(i)=(1.-frh)**2
!frh_out(i) = frh
if(forcing(i,7).eq.0.)sig(i)=1.
if(kdt.le.(3600./dtime))sig(i)=1.
frh_out(i) = frh*sig(i)
enddo
!$acc end kernels
sig_thresh = (1.-frh_thresh)**2
!
!--- entrainment of mass
!
!
!--- initial detrainmentrates
!
!$acc kernels
do k=kts,ktf
do i=its,itf
cnvwt(i,k)=0.
zuo(i,k)=0.
zdo(i,k)=0.
z(i,k)=zo(i,k)
xz(i,k)=zo(i,k)
cupclw(i,k)=0.
cd(i,k)=.1*entr_rate(i) !1.e-9 ! 1.*entr_rate
if(imid.eq.1)cd(i,k)=.5*entr_rate(i)
cdd(i,k)=1.e-9
hcdo(i,k)=0.
qrcdo(i,k)=0.
dellaqc(i,k)=0.
enddo
enddo
!$acc end kernels
!
!--- max/min allowed value for epsilon (ratio downdraft base mass flux/updraft
! base mass flux
!
!$acc kernels
edtmax(:)=1.
! if(imid.eq.1)edtmax(:)=.15
edtmin(:)=.1
! if(imid.eq.1)edtmin(:)=.05
!$acc end kernels
!
!--- minimum depth (m), clouds must have
!
depth_min=3000.
!--- for RRFS allow only very deep convection
if(dx(its)<dx_thresh)depth_min=5000.
if(imid.eq.1)depth_min=2500.
!
!--- maximum depth (mb) of capping
!--- inversion (larger cap = no convection)
!
!$acc kernels
do i=its,itf
kbmax(i)=1
aa0(i)=0.
aa1(i)=0.
edt(i)=0.
kstabm(i)=ktf-1
ierr2(i)=0
ierr3(i)=0
enddo
!$acc end kernels
x_add=0.
!
!--- max height(m) above ground where updraft air can originate
!
zkbmax=4000.
if(imid.eq.1)zkbmax=2000.
!
!--- height(m) above which no downdrafts are allowed to originate
!
zcutdown=4000.
!
!--- depth(m) over which downdraft detrains all its mass
!
z_detr=500.
!
!
!--- environmental conditions, first heights
!
!$acc kernels
do i=its,itf
do k=1,maxens3
xf_ens(i,k)=0.
pr_ens(i,k)=0.
enddo
enddo
!$acc end kernels
!
!> - Call cup_env() to calculate moist static energy, heights, qes
!
call cup_env(z,qes,he,hes,t,q,po,z1, &
psur,ierr,tcrit,-1, &
itf,ktf, &
its,ite, kts,kte)
call cup_env(zo,qeso,heo,heso,tn,qo,po,z1, &
psur,ierr,tcrit,-1, &
itf,ktf, &
its,ite, kts,kte)
!
!> - Call cup_env_clev() to calculate environmental values on cloud levels
!
call cup_env_clev(t,qes,q,he,hes,z,po,qes_cup,q_cup,he_cup, &
hes_cup,z_cup,p_cup,gamma_cup,t_cup,psur, &
ierr,z1, &
itf,ktf, &
its,ite, kts,kte)
call cup_env_clev(tn,qeso,qo,heo,heso,zo,po,qeso_cup,qo_cup, &
heo_cup,heso_cup,zo_cup,po_cup,gammao_cup,tn_cup,psur, &
ierr,z1, &
itf,ktf, &
its,ite, kts,kte)
!---meltglac-------------------------------------------------
!> - Call get_partition_liq_ice() to calculate partition between liq/ice cloud contents
call get_partition_liq_ice(ierr,tn,po_cup,p_liq_ice,melting_layer,&
itf,ktf,its,ite,kts,kte,cumulus)
!---meltglac-------------------------------------------------
!$acc kernels
do i=its,itf
if(ierr(i).eq.0)then
if(kpbl(i).gt.5 .and. imid.eq.1)cap_max(i)=po_cup(i,kpbl(i))
u_cup(i,kts)=us(i,kts)
v_cup(i,kts)=vs(i,kts)
do k=kts+1,ktf
u_cup(i,k)=.5*(us(i,k-1)+us(i,k))
v_cup(i,k)=.5*(vs(i,k-1)+vs(i,k))
enddo
endif
enddo
do i=its,itf
if(ierr(i).eq.0)then
do k=kts,ktf
if(zo_cup(i,k).gt.zkbmax+z1(i))then
kbmax(i)=k
go to 25
endif
enddo
25 continue
!
!> - Compute the level where detrainment for downdraft starts (\p kdet)
!
do k=kts,ktf
if(zo_cup(i,k).gt.z_detr+z1(i))then
kdet(i)=k
go to 26
endif
enddo
26 continue
!
endif
enddo
!$acc end kernels
!
!
!
!> - Determine level with highest moist static energy content (\p k22)
!
start_k22=2
!$acc parallel loop
do 36 i=its,itf
if(ierr(i).eq.0)then
k22(i)=maxloc(heo_cup(i,start_k22:kbmax(i)+2),1)+start_k22-1
if(k22(i).ge.kbmax(i))then
ierr(i)=2
#ifndef _OPENACC
ierrc(i)="could not find k22"
#endif
ktop(i)=0
k22(i)=0
kbcon(i)=0
endif
endif
36 continue
!$acc end parallel
!
!> - call get_cloud_bc() and cup_kbcon() to determine the
!! level of convective cloud base (\p kbcon)
!
!$acc parallel loop private(x_add)
do i=its,itf
if(ierr(i).eq.0)then
x_add = xlv*zqexec(i)+cp*ztexec(i)
call get_cloud_bc(kte,he_cup (i,1:kte),hkb (i),k22(i),x_add)
call get_cloud_bc(kte,heo_cup (i,1:kte),hkbo (i),k22(i),x_add)
endif ! ierr
enddo
!$acc end parallel
jprnt=0
iloop=1
if(imid.eq.1)iloop=5
call cup_kbcon(ierrc,cap_max_increment,iloop,k22,kbcon,heo_cup,heso_cup, &
hkbo,ierr,kbmax,po_cup,cap_max, &
ztexec,zqexec, &
jprnt,itf,ktf, &
its,ite, kts,kte, &
z_cup,entr_rate,heo,imid)
!
!> - Call cup_minimi() to increase detrainment in stable layers
!
call cup_minimi(heso_cup,kbcon,kstabm,kstabi,ierr, &
itf,ktf, &
its,ite, kts,kte)
!$acc parallel loop private(frh,x_add)
do i=its,itf
if(ierr(i) == 0)then
frh = min(qo_cup(i,kbcon(i))/qeso_cup(i,kbcon(i)),1.)
if(frh >= rh_thresh .and. sig(i) <= sig_thresh )then
ierr(i)=231
cycle
endif
!
! never go too low...
!
x_add=0.
!$acc loop seq
do k=kbcon(i)+1,ktf
if(po(i,kbcon(i))-po(i,k) > pmin+x_add)then
pmin_lev(i)=k
exit
endif
enddo
!
!> - Call get_cloud_bc() to initial conditions for updraft
!
start_level(i)=k22(i)
x_add = xlv*zqexec(i)+cp*ztexec(i)
call get_cloud_bc(kte,he_cup (i,1:kte),hkb (i),k22(i),x_add)
endif
enddo
!$acc end parallel
!
!> - Call get_inversion_layer() to get inversion layers for mid level cloud tops
!
if(imid.eq.1)then
call get_inversion_layers(ierr,p_cup,t_cup,z_cup,q_cup,qes_cup,k_inv_layers, &
kbcon,kstabi,dtempdz,itf,ktf,its,ite, kts,kte)
endif
!$acc kernels
do i=its,itf
if(kstabi(i).lt.kbcon(i))then
kbcon(i)=1
ierr(i)=42
endif
do k=kts,ktf
entr_rate_2d(i,k)=entr_rate(i)
enddo
if(ierr(i).eq.0)then
kbcon(i)=max(2,kbcon(i))
do k=kts+1,ktf
frh = min(qo_cup(i,k)/qeso_cup(i,k),1.)
entr_rate_2d(i,k)=entr_rate(i) *(1.3-frh)
enddo
if(imid.eq.1)then
if(k_inv_layers(i,2).gt.0 .and. &
(po_cup(i,k22(i))-po_cup(i,k_inv_layers(i,2))).lt.500.)then
ktop(i)=min(kstabi(i),k_inv_layers(i,2))
ktopdby(i)=ktop(i)
else
!$acc loop seq
do k=kbcon(i)+1,ktf
if((po_cup(i,k22(i))-po_cup(i,k)).gt.500.)then
ktop(i)=k
ktopdby(i)=ktop(i)
exit
endif
enddo
endif ! k_inv_lay
endif
endif
enddo
!$acc end kernels
!
!> - Call rates_up_pdf() to get normalized mass flux, entrainment and detrainmentrates for updraft
!
i=0
!- for mid level clouds we do not allow clouds taller than where stability
!- changes
if(imid.eq.1)then
call rates_up_pdf(rand_vmas,ipr,'mid',ktop,ierr,po_cup,entr_rate_2d,hkbo,heo,heso_cup,zo_cup, &
xland1,kstabi,k22,kbcon,its,ite,itf,kts,kte,ktf,zuo,kpbl,ktopdby,csum,pmin_lev)
else
call rates_up_pdf(rand_vmas,ipr,'deep',ktop,ierr,po_cup,entr_rate_2d,hkbo,heo,heso_cup,zo_cup, &
xland1,kstabi,k22,kbcon,its,ite,itf,kts,kte,ktf,zuo,kbcon,ktopdby,csum,pmin_lev)
endif
!
!
!
!$acc kernels
do i=its,itf
if(ierr(i).eq.0)then
if(k22(i).gt.1)then
!$acc loop independent
do k=1,k22(i) -1
zuo(i,k)=0.
zu (i,k)=0.
xzu(i,k)=0.
enddo
endif
!$acc loop independent
do k=k22(i),ktop(i)
xzu(i,k)= zuo(i,k)
zu (i,k)= zuo(i,k)
enddo
!$acc loop independent
do k=ktop(i)+1,kte
zuo(i,k)=0.
zu (i,k)=0.
xzu(i,k)=0.
enddo
endif
enddo
!$acc end kernels
!
!> - Call get_lateral_massflux() to calculate mass entrainment and detrainment
!
if(imid.eq.1)then
call get_lateral_massflux(itf,ktf, its,ite, kts,kte &
,ierr,ktop,zo_cup,zuo,cd,entr_rate_2d &
,up_massentro, up_massdetro ,up_massentr, up_massdetr &
,3,kbcon,k22,up_massentru,up_massdetru,lambau)
else
call get_lateral_massflux(itf,ktf, its,ite, kts,kte &
,ierr,ktop,zo_cup,zuo,cd,entr_rate_2d &
,up_massentro, up_massdetro ,up_massentr, up_massdetr &
,1,kbcon,k22,up_massentru,up_massdetru,lambau)
endif
!
! note: ktop here already includes overshooting, ktopdby is without
! overshooting
!
!$acc kernels
do k=kts,ktf
do i=its,itf
uc (i,k)=0.
vc (i,k)=0.
hc (i,k)=0.
dby (i,k)=0.
hco (i,k)=0.
dbyo(i,k)=0.
enddo
enddo
do i=its,itf
if(ierr(i).eq.0)then
do k=1,start_level(i)
uc(i,k)=u_cup(i,k)
vc(i,k)=v_cup(i,k)
enddo
do k=1,start_level(i)-1
hc (i,k)=he_cup(i,k)
hco(i,k)=heo_cup(i,k)
enddo
k=start_level(i)
hc (i,k)=hkb(i)
hco(i,k)=hkbo(i)
endif
enddo
!$acc end kernels
!
!---meltglac-------------------------------------------------
!
!--- 1st guess for moist static energy and dbyo (not including ice phase)
!
!$acc parallel loop private(denom,kk,ki)
do i=its,itf
ktopkeep(i)=0
dbyt(i,:)=0.
if(ierr(i) /= 0) cycle
ktopkeep(i)=ktop(i)
!$acc loop seq
do k=start_level(i) +1,ktop(i) !mass cons option
denom=zuo(i,k-1)-.5*up_massdetro(i,k-1)+up_massentro(i,k-1)
if(denom.lt.1.e-8)then
ierr(i)=51
exit
endif
hco(i,k)=(hco(i,k-1)*zuo(i,k-1)-.5*up_massdetro(i,k-1)*hco(i,k-1)+ &
up_massentro(i,k-1)*heo(i,k-1)) / &
(zuo(i,k-1)-.5*up_massdetro(i,k-1)+up_massentro(i,k-1))
dbyo(i,k)=hco(i,k)-heso_cup(i,k)
enddo
! for now no overshooting (only very little)
!$acc loop seq
do k=ktop(i)-1,kbcon(i),-1
if(dbyo(i,k).gt.0.)then
ktopkeep(i)=k+1
exit
endif
enddo
enddo
!$acc end parallel
!$acc kernels
do 37 i=its,itf