-
Notifications
You must be signed in to change notification settings - Fork 153
/
Copy pathGFS_suite_interstitial.F90
972 lines (790 loc) · 38.3 KB
/
GFS_suite_interstitial.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
!> \file GFS_suite_interstitial.f90
!! Contains code related to more than one scheme in the GFS physics suite.
module GFS_suite_interstitial_rad_reset
contains
subroutine GFS_suite_interstitial_rad_reset_init ()
end subroutine GFS_suite_interstitial_rad_reset_init
subroutine GFS_suite_interstitial_rad_reset_finalize()
end subroutine GFS_suite_interstitial_rad_reset_finalize
!> \section arg_table_GFS_suite_interstitial_rad_reset_run Argument Table
!! \htmlinclude GFS_suite_interstitial_rad_reset_run.html
!!
subroutine GFS_suite_interstitial_rad_reset_run (Interstitial, Model, errmsg, errflg)
use machine, only: kind_phys
use GFS_typedefs, only: GFS_control_type, GFS_interstitial_type
implicit none
! interface variables
type(GFS_interstitial_type), intent(inout) :: Interstitial
type(GFS_control_type), intent(in) :: Model
character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
errmsg = ''
errflg = 0
call Interstitial%rad_reset(Model)
end subroutine GFS_suite_interstitial_rad_reset_run
end module GFS_suite_interstitial_rad_reset
module GFS_suite_interstitial_phys_reset
contains
subroutine GFS_suite_interstitial_phys_reset_init ()
end subroutine GFS_suite_interstitial_phys_reset_init
subroutine GFS_suite_interstitial_phys_reset_finalize()
end subroutine GFS_suite_interstitial_phys_reset_finalize
!> \section arg_table_GFS_suite_interstitial_phys_reset_run Argument Table
!! \htmlinclude GFS_suite_interstitial_phys_reset_run.html
!!
subroutine GFS_suite_interstitial_phys_reset_run (Interstitial, Model, errmsg, errflg)
use machine, only: kind_phys
use GFS_typedefs, only: GFS_control_type, GFS_interstitial_type
implicit none
! interface variables
type(GFS_interstitial_type), intent(inout) :: Interstitial
type(GFS_control_type), intent(in ) :: Model
character(len=*), intent( out) :: errmsg
integer, intent( out) :: errflg
errmsg = ''
errflg = 0
call Interstitial%phys_reset(Model)
end subroutine GFS_suite_interstitial_phys_reset_run
end module GFS_suite_interstitial_phys_reset
module GFS_suite_interstitial_1
contains
subroutine GFS_suite_interstitial_1_init ()
end subroutine GFS_suite_interstitial_1_init
subroutine GFS_suite_interstitial_1_finalize()
end subroutine GFS_suite_interstitial_1_finalize
!> \section arg_table_GFS_suite_interstitial_1_run Argument Table
!! \htmlinclude GFS_suite_interstitial_1_run.html
!!
subroutine GFS_suite_interstitial_1_run (im, levs, ntrac, dtf, dtp, slmsk, area, dxmin, dxinv, pgr, &
islmsk, work1, work2, psurf, dudt, dvdt, dtdt, dqdt, errmsg, errflg)
use machine, only: kind_phys
implicit none
! interface variables
integer, intent(in ) :: im, levs, ntrac
real(kind=kind_phys), intent(in ) :: dtf, dtp, dxmin, dxinv
real(kind=kind_phys), intent(in ), dimension(:) :: slmsk, area, pgr
integer, intent(out), dimension(:) :: islmsk
real(kind=kind_phys), intent(out), dimension(:) :: work1, work2, psurf
real(kind=kind_phys), intent(out), dimension(:,:) :: dudt, dvdt, dtdt
real(kind=kind_phys), intent(out), dimension(:,:,:) :: dqdt
character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
! local variables
real(kind=kind_phys), parameter :: zero = 0.0_kind_phys, one = 1.0_kind_phys
integer :: i, k, n
! Initialize CCPP error handling variables
errmsg = ''
errflg = 0
do i = 1, im
islmsk(i) = nint(slmsk(i))
work1(i) = (log(area(i)) - dxmin) * dxinv
work1(i) = max(zero, min(one, work1(i)))
work2(i) = one - work1(i)
psurf(i) = pgr(i)
end do
do k=1,levs
do i=1,im
dudt(i,k) = zero
dvdt(i,k) = zero
dtdt(i,k) = zero
enddo
enddo
do n=1,ntrac
do k=1,levs
do i=1,im
dqdt(i,k,n) = zero
enddo
enddo
enddo
end subroutine GFS_suite_interstitial_1_run
end module GFS_suite_interstitial_1
module GFS_suite_interstitial_2
use machine, only: kind_phys
real(kind=kind_phys), parameter :: one = 1.0_kind_phys
logical :: linit_mod = .false.
contains
subroutine GFS_suite_interstitial_2_init ()
end subroutine GFS_suite_interstitial_2_init
subroutine GFS_suite_interstitial_2_finalize()
end subroutine GFS_suite_interstitial_2_finalize
!> \section arg_table_GFS_suite_interstitial_2_run Argument Table
!! \htmlinclude GFS_suite_interstitial_2_run.html
!!
subroutine GFS_suite_interstitial_2_run (im, levs, lssav, ldiag3d, lsidea, flag_cice, shal_cnv, old_monin, mstrat, &
do_shoc, frac_grid, imfshalcnv, dtf, xcosz, adjsfcdsw, adjsfcdlw, cice, pgr, ulwsfc_cice, lwhd, htrsw, htrlw, xmu, ctei_rm, &
work1, work2, prsi, tgrs, prsl, qgrs_water_vapor, qgrs_cloud_water, cp, hvap, prslk, suntim, adjsfculw, adjsfculw_lnd, &
adjsfculw_ice, adjsfculw_wat, dlwsfc, ulwsfc, psmean, dtend, dtidx, index_of_process_longwave, index_of_process_shortwave, &
index_of_process_pbl, index_of_process_dcnv, index_of_process_scnv, index_of_process_mp, index_of_temperature, &
ctei_rml, ctei_r, kinver, dry, icy, wet, frland, huge, use_LW_jacobian, htrlwu, errmsg, errflg)
implicit none
! interface variables
integer, intent(in ) :: im, levs, imfshalcnv
logical, intent(in ) :: lssav, ldiag3d, lsidea, shal_cnv
logical, intent(in ) :: old_monin, mstrat, do_shoc, frac_grid, use_LW_jacobian
real(kind=kind_phys), intent(in ) :: dtf, cp, hvap
logical, intent(in ), dimension(:) :: flag_cice
real(kind=kind_phys), intent(in ), dimension(:) :: ctei_rm
real(kind=kind_phys), intent(in ), dimension(:) :: xcosz, adjsfcdsw, adjsfcdlw, pgr, xmu, work1, work2
real(kind=kind_phys), intent(in ), dimension(:) :: ulwsfc_cice
real(kind=kind_phys), intent(in ), dimension(:) :: cice
real(kind=kind_phys), intent(in ), dimension(:,:) :: htrsw, htrlw, htrlwu, tgrs, prsl, qgrs_water_vapor, qgrs_cloud_water, prslk
real(kind=kind_phys), intent(in ), dimension(:,:) :: prsi
real(kind=kind_phys), intent(in ), dimension(:,:,:) :: lwhd
integer, intent(inout), dimension(:) :: kinver
real(kind=kind_phys), intent(inout), dimension(:) :: suntim, dlwsfc, ulwsfc, psmean, ctei_rml, ctei_r
real(kind=kind_phys), intent(in ), dimension(:) :: adjsfculw_lnd, adjsfculw_ice, adjsfculw_wat
real(kind=kind_phys), intent(inout), dimension(:) :: adjsfculw
! dtend is only allocated if ldiag3d is .true.
real(kind=kind_phys), optional, intent(inout), dimension(:,:,:) :: dtend
integer, intent(in), dimension(:,:) :: dtidx
integer, intent(in) :: index_of_process_longwave, index_of_process_shortwave, &
index_of_process_pbl, index_of_process_dcnv, index_of_process_scnv, &
index_of_process_mp, index_of_temperature
logical, intent(in ), dimension(:) :: dry, icy, wet
real(kind=kind_phys), intent(in ), dimension(:) :: frland
real(kind=kind_phys), intent(in ) :: huge
character(len=*), intent( out) :: errmsg
integer, intent( out) :: errflg
! local variables
real(kind=kind_phys), parameter :: czmin = 0.0001_kind_phys ! cos(89.994)
integer :: i, k, idtend
real(kind=kind_phys) :: tem1, tem2, tem, hocp
logical, dimension(im) :: invrsn
real(kind=kind_phys), dimension(im) :: tx1, tx2
real(kind=kind_phys), parameter :: zero = 0.0_kind_phys, one = 1.0_kind_phys
real(kind=kind_phys), parameter :: qmin = 1.0e-10_kind_phys, epsln=1.0e-10_kind_phys
! Initialize CCPP error handling variables
errmsg = ''
errflg = 0
hocp = hvap/cp
if (lssav) then ! --- ... accumulate/save output variables
! --- ... sunshine duration time is defined as the length of time (in mdl output
! interval) that solar radiation falling on a plane perpendicular to the
! direction of the sun >= 120 w/m2
do i = 1, im
if ( xcosz(i) >= czmin ) then ! zenth angle > 89.994 deg
tem1 = adjsfcdsw(i) / xcosz(i)
if ( tem1 >= 120.0_kind_phys ) then
suntim(i) = suntim(i) + dtf
endif
endif
enddo
! --- ... sfc lw fluxes used by atmospheric model are saved for output
if (frac_grid) then
do i=1,im
tem = (one - frland(i)) * cice(i) ! tem = ice fraction wrt whole cell
if (flag_cice(i)) then
adjsfculw(i) = adjsfculw_lnd(i) * frland(i) &
+ ulwsfc_cice(i) * tem &
+ adjsfculw_wat(i) * (one - frland(i) - tem)
else
adjsfculw(i) = adjsfculw_lnd(i) * frland(i) &
+ adjsfculw_ice(i) * tem &
+ adjsfculw_wat(i) * (one - frland(i) - tem)
endif
enddo
else
do i=1,im
if (dry(i)) then ! all land
adjsfculw(i) = adjsfculw_lnd(i)
elseif (icy(i)) then ! ice (and water)
tem = one - cice(i)
if (flag_cice(i)) then
if (wet(i) .and. adjsfculw_wat(i) /= huge) then
adjsfculw(i) = ulwsfc_cice(i)*cice(i) + adjsfculw_wat(i)*tem
else
adjsfculw(i) = ulwsfc_cice(i)
endif
else
if (wet(i) .and. adjsfculw_wat(i) /= huge) then
adjsfculw(i) = adjsfculw_ice(i)*cice(i) + adjsfculw_wat(i)*tem
else
adjsfculw(i) = adjsfculw_ice(i)
endif
endif
else ! all water
adjsfculw(i) = adjsfculw_wat(i)
endif
enddo
endif
do i=1,im
dlwsfc(i) = dlwsfc(i) + adjsfcdlw(i)*dtf
ulwsfc(i) = ulwsfc(i) + adjsfculw(i)*dtf
psmean(i) = psmean(i) + pgr(i)*dtf ! mean surface pressure
enddo
if (ldiag3d) then
if (lsidea) then
idtend = dtidx(index_of_temperature,index_of_process_longwave)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + lwhd(:,:,1)*dtf
endif
idtend = dtidx(index_of_temperature,index_of_process_shortwave)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + lwhd(:,:,2)*dtf
endif
idtend = dtidx(index_of_temperature,index_of_process_pbl)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + lwhd(:,:,3)*dtf
endif
idtend = dtidx(index_of_temperature,index_of_process_dcnv)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + lwhd(:,:,4)*dtf
endif
idtend = dtidx(index_of_temperature,index_of_process_scnv)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + lwhd(:,:,5)*dtf
endif
idtend = dtidx(index_of_temperature,index_of_process_mp)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + lwhd(:,:,6)*dtf
endif
else
idtend = dtidx(index_of_temperature,index_of_process_longwave)
if(idtend>=1) then
if (use_LW_jacobian) then
dtend(:,:,idtend) = dtend(:,:,idtend) + htrlwu(:,:)*dtf
else
dtend(:,:,idtend) = dtend(:,:,idtend) + htrlw(:,:)*dtf
endif
endif
idtend = dtidx(index_of_temperature,index_of_process_shortwave)
if(idtend>=1) then
do k=1,levs
do i=1,im
dtend(i,k,idtend) = dtend(i,k,idtend) + htrsw(i,k)*dtf*xmu(i)
enddo
enddo
endif
endif
endif
endif ! end if_lssav_block
do i=1, im
invrsn(i) = .false.
tx1(i) = zero
tx2(i) = 10.0_kind_phys
ctei_r(i) = 10.0_kind_phys
enddo
if ((((imfshalcnv == 0 .and. shal_cnv) .or. old_monin) .and. mstrat) &
.or. do_shoc) then
ctei_rml(:) = ctei_rm(1)*work1(:) + ctei_rm(2)*work2(:)
do k=1,levs/2
do i=1,im
if (prsi(i,1)-prsi(i,k+1) < 0.35_kind_phys*prsi(i,1) &
.and. (.not. invrsn(i))) then
tem = (tgrs(i,k+1) - tgrs(i,k)) &
/ (prsl(i,k) - prsl(i,k+1))
if (((tem > 0.0001_kind_phys) .and. (tx1(i) < zero)) .or. &
((tem-abs(tx1(i)) > zero) .and. (tx2(i) < zero))) then
invrsn(i) = .true.
if (qgrs_water_vapor(i,k) > qgrs_water_vapor(i,k+1)) then
tem1 = tgrs(i,k+1) + hocp*max(qgrs_water_vapor(i,k+1),qmin)
tem2 = tgrs(i,k) + hocp*max(qgrs_water_vapor(i,k),qmin)
tem1 = tem1 / prslk(i,k+1) - tem2 / prslk(i,k)
! --- ... (cp/l)(deltathetae)/(deltatwater) > ctei_rm -> conditon for CTEI
ctei_r(i) = (one/hocp)*tem1/(qgrs_water_vapor(i,k+1)-qgrs_water_vapor(i,k) &
+ qgrs_cloud_water(i,k+1)-qgrs_cloud_water(i,k))
else
ctei_r(i) = 10.0_kind_phys
endif
if ( ctei_rml(i) > ctei_r(i) ) then
kinver(i) = k
else
kinver(i) = levs
endif
endif
tx2(i) = tx1(i)
tx1(i) = tem
endif
enddo
enddo
endif
end subroutine GFS_suite_interstitial_2_run
end module GFS_suite_interstitial_2
module GFS_suite_stateout_reset
contains
subroutine GFS_suite_stateout_reset_init ()
end subroutine GFS_suite_stateout_reset_init
subroutine GFS_suite_stateout_reset_finalize()
end subroutine GFS_suite_stateout_reset_finalize
!> \section arg_table_GFS_suite_stateout_reset_run Argument Table
!! \htmlinclude GFS_suite_stateout_reset_run.html
!!
subroutine GFS_suite_stateout_reset_run (im, levs, ntrac, &
tgrs, ugrs, vgrs, qgrs, &
gt0 , gu0 , gv0 , gq0 , &
errmsg, errflg)
use machine, only: kind_phys
implicit none
! interface variables
integer, intent(in ) :: im
integer, intent(in ) :: levs
integer, intent(in ) :: ntrac
real(kind=kind_phys), intent(in ), dimension(:,:) :: tgrs, ugrs, vgrs
real(kind=kind_phys), intent(in ), dimension(:,:,:) :: qgrs
real(kind=kind_phys), intent(out), dimension(:,:) :: gt0, gu0, gv0
real(kind=kind_phys), intent(out), dimension(:,:,:) :: gq0
character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
! Initialize CCPP error handling variables
errmsg = ''
errflg = 0
gt0(:,:) = tgrs(:,:)
gu0(:,:) = ugrs(:,:)
gv0(:,:) = vgrs(:,:)
gq0(:,:,:) = qgrs(:,:,:)
end subroutine GFS_suite_stateout_reset_run
end module GFS_suite_stateout_reset
module GFS_suite_stateout_update
contains
subroutine GFS_suite_stateout_update_init ()
end subroutine GFS_suite_stateout_update_init
subroutine GFS_suite_stateout_update_finalize()
end subroutine GFS_suite_stateout_update_finalize
!> \section arg_table_GFS_suite_stateout_update_run Argument Table
!! \htmlinclude GFS_suite_stateout_update_run.html
!!
subroutine GFS_suite_stateout_update_run (im, levs, ntrac, dtp, &
tgrs, ugrs, vgrs, qgrs, dudt, dvdt, dtdt, dqdt, &
gt0, gu0, gv0, gq0, ntiw, nqrimef, imp_physics, &
imp_physics_fer_hires, epsq, errmsg, errflg)
use machine, only: kind_phys
implicit none
! Interface variables
integer, intent(in ) :: im
integer, intent(in ) :: levs
integer, intent(in ) :: ntrac
integer, intent(in ) :: imp_physics,imp_physics_fer_hires
integer, intent(in ) :: ntiw, nqrimef
real(kind=kind_phys), intent(in ) :: dtp, epsq
real(kind=kind_phys), intent(in ), dimension(:,:) :: tgrs, ugrs, vgrs
real(kind=kind_phys), intent(in ), dimension(:,:,:) :: qgrs
real(kind=kind_phys), intent(in ), dimension(:,:) :: dudt, dvdt, dtdt
real(kind=kind_phys), intent(in ), dimension(:,:,:) :: dqdt
real(kind=kind_phys), intent(out), dimension(:,:) :: gt0, gu0, gv0
real(kind=kind_phys), intent(out), dimension(:,:,:) :: gq0
character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
integer :: i, k
! Initialize CCPP error handling variables
errmsg = ''
errflg = 0
gt0(:,:) = tgrs(:,:) + dtdt(:,:) * dtp
gu0(:,:) = ugrs(:,:) + dudt(:,:) * dtp
gv0(:,:) = vgrs(:,:) + dvdt(:,:) * dtp
gq0(:,:,:) = qgrs(:,:,:) + dqdt(:,:,:) * dtp
if (imp_physics == imp_physics_fer_hires) then
do k=1,levs
do i=1,im
if(gq0(i,k,ntiw) > epsq) then
gq0(i,k,nqrimef) = max(1., gq0(i,k,nqrimef)/gq0(i,k,ntiw))
else
gq0(i,k,nqrimef) = 1.
end if
end do
end do
end if
end subroutine GFS_suite_stateout_update_run
end module GFS_suite_stateout_update
module GFS_suite_interstitial_3
contains
subroutine GFS_suite_interstitial_3_init ()
end subroutine GFS_suite_interstitial_3_init
subroutine GFS_suite_interstitial_3_finalize()
end subroutine GFS_suite_interstitial_3_finalize
!> \section arg_table_GFS_suite_interstitial_3_run Argument Table
!! \htmlinclude GFS_suite_interstitial_3_run.html
!!
subroutine GFS_suite_interstitial_3_run (im, levs, nn, cscnv, &
satmedmf, trans_trac, do_shoc, ltaerosol, ntrac, ntcw, &
ntiw, ntclamt, ntrw, ntsw, ntrnc, ntsnc, ntgl, ntgnc, &
xlon, xlat, gt0, gq0, imp_physics, imp_physics_mg, &
imp_physics_zhao_carr, imp_physics_zhao_carr_pdf, &
imp_physics_gfdl, imp_physics_thompson, dtidx, ntlnc, &
imp_physics_wsm6, imp_physics_fer_hires, prsi, ntinc, &
prsl, prslk, rhcbot,rhcpbl, rhctop, rhcmax, islmsk, &
work1, work2, kpbl, kinver, ras, me, save_lnc, save_inc, &
ldiag3d, qdiag3d, index_of_process_conv_trans, &
clw, rhc, save_qc, save_qi, save_tcp, errmsg, errflg)
use machine, only: kind_phys
implicit none
! interface variables
integer, intent(in ) :: im, levs, nn, ntrac, ntcw, ntiw, ntclamt, ntrw, ntsw,&
ntrnc, ntsnc, ntgl, ntgnc, imp_physics, imp_physics_mg, imp_physics_zhao_carr, imp_physics_zhao_carr_pdf, &
imp_physics_gfdl, imp_physics_thompson, imp_physics_wsm6,imp_physics_fer_hires, me, index_of_process_conv_trans
integer, intent(in ), dimension(:) :: islmsk, kpbl, kinver
logical, intent(in ) :: cscnv, satmedmf, trans_trac, do_shoc, ltaerosol, ras
integer, intent(in) :: ntinc, ntlnc
logical, intent(in) :: ldiag3d, qdiag3d
integer, dimension(:,:), intent(in) :: dtidx
real, dimension(:,:), intent(out) :: save_lnc, save_inc
real(kind=kind_phys), intent(in ) :: rhcbot, rhcmax, rhcpbl, rhctop
real(kind=kind_phys), intent(in ), dimension(:) :: work1, work2
real(kind=kind_phys), intent(in ), dimension(:,:) :: prsl, prslk
real(kind=kind_phys), intent(in ), dimension(:,:) :: prsi
real(kind=kind_phys), intent(in ), dimension(:) :: xlon, xlat
real(kind=kind_phys), intent(in ), dimension(:,:) :: gt0
real(kind=kind_phys), intent(in ), dimension(:,:,:) :: gq0
real(kind=kind_phys), intent(inout), dimension(:,:) :: rhc, save_qc
! save_qi is not allocated for Zhao-Carr MP
real(kind=kind_phys), intent(inout), dimension(:,:) :: save_qi
real(kind=kind_phys), intent(inout), dimension(:,:) :: save_tcp
real(kind=kind_phys), intent(inout), dimension(:,:,:) :: clw
character(len=*), intent( out) :: errmsg
integer, intent( out) :: errflg
! local variables
integer :: i,k,n,tracers,kk
real(kind=kind_phys) :: tem, tem1, tem2
real(kind=kind_phys), dimension(im) :: tx1, tx2, tx3, tx4
!real(kind=kind_phys),parameter :: slope_mg = 0.02, slope_upmg = 0.04, &
! turnrhcrit = 0.900, turnrhcrit_upper = 0.150
! in the following inverse of slope_mg and slope_upmg are specified
real(kind=kind_phys), parameter :: zero = 0.0_kind_phys, one = 1.0_kind_phys
real(kind=kind_phys), parameter :: slope_mg = 50.0_kind_phys, &
slope_upmg = 25.0_kind_phys
! Initialize CCPP error handling variables
errmsg = ''
errflg = 0
if (cscnv .or. satmedmf .or. trans_trac .or. ras) then
tracers = 2
do n=2,ntrac
if ( n /= ntcw .and. n /= ntiw .and. n /= ntclamt .and. &
n /= ntrw .and. n /= ntsw .and. n /= ntrnc .and. &
n /= ntsnc .and. n /= ntgl .and. n /= ntgnc) then
tracers = tracers + 1
do k=1,levs
do i=1,im
clw(i,k,tracers) = gq0(i,k,n)
enddo
enddo
endif
enddo
endif ! end if_ras or cfscnv or samf
if (ntcw > 0) then
if (imp_physics == imp_physics_mg .and. rhcpbl < 0.5_kind_phys) then ! compute rhc for GMAO macro physics cloud pdf
do i=1,im
tx1(i) = one / prsi(i,1)
tx2(i) = one - rhcmax*work1(i)-rhcbot*work2(i)
kk = min(kinver(i), max(2,kpbl(i)))
tx3(i) = prsi(i,kk)*tx1(i)
tx4(i) = rhcpbl - rhctop*abs(cos(xlat(i)))
enddo
do k = 1, levs
do i = 1, im
tem = prsl(i,k) * tx1(i)
tem1 = min(max((tem-tx3(i))*slope_mg, -20.0_kind_phys), 20.0_kind_phys)
! Using rhcpbl and rhctop from the namelist instead of 0.3 and 0.2
! and rhcbot represents pbl top critical relative humidity
tem2 = min(max((tx4(i)-tem)*slope_upmg, -20.0_kind_phys), 20.0_kind_phys) ! Anning
if (islmsk(i) > 0) then
tem1 = one / (one+exp(tem1+tem1))
else
tem1 = 2.0_kind_phys / (one+exp(tem1+tem1))
endif
tem2 = one / (one+exp(tem2))
rhc(i,k) = min(rhcmax, max(0.7_kind_phys, one-tx2(i)*tem1*tem2))
enddo
enddo
else
do k=1,levs
do i=1,im
kk = max(10,kpbl(i))
if (k < kk) then
tem = rhcbot - (rhcbot-rhcpbl) * (one-prslk(i,k)) / (one-prslk(i,kk))
else
tem = rhcpbl - (rhcpbl-rhctop) * (prslk(i,kk)-prslk(i,k)) / prslk(i,kk)
endif
tem = rhcmax * work1(i) + tem * work2(i)
rhc(i,k) = max(zero, min(one,tem))
enddo
enddo
endif
else
rhc(:,:) = 1.0
endif
if (imp_physics == imp_physics_zhao_carr .or. imp_physics == imp_physics_zhao_carr_pdf) then ! zhao-carr microphysics
!GF* move to GFS_MP_generic_pre (from gscond/precpd)
! do i=1,im
! psautco_l(i) = Model%psautco(1)*work1(i) + Model%psautco(2)*work2(i)
! prautco_l(i) = Model%prautco(1)*work1(i) + Model%prautco(2)*work2(i)
! enddo
!*GF
do k=1,levs
do i=1,im
clw(i,k,1) = gq0(i,k,ntcw)
enddo
enddo
elseif (imp_physics == imp_physics_gfdl) then
clw(1:im,:,1) = gq0(1:im,:,ntcw)
elseif (imp_physics == imp_physics_thompson) then
do k=1,levs
do i=1,im
clw(i,k,1) = gq0(i,k,ntiw) ! ice
clw(i,k,2) = gq0(i,k,ntcw) ! water
save_tcp(i,k) = gt0(i,k)
enddo
enddo
if(ltaerosol) then
save_qi(:,:) = clw(:,:,1)
save_qc(:,:) = clw(:,:,2)
else
save_qi(:,:) = clw(:,:,1)
endif
elseif (imp_physics == imp_physics_wsm6 .or. imp_physics == imp_physics_mg .or. imp_physics == imp_physics_fer_hires) then
do k=1,levs
do i=1,im
clw(i,k,1) = gq0(i,k,ntiw) ! ice
clw(i,k,2) = gq0(i,k,ntcw) ! water
enddo
enddo
endif
if(imp_physics == imp_physics_thompson .and. ldiag3d .and. qdiag3d) then
if(dtidx(100+ntlnc,index_of_process_conv_trans)>0) then
save_lnc = gq0(:,:,ntlnc)
endif
if(dtidx(100+ntinc,index_of_process_conv_trans)>0) then
save_inc = gq0(:,:,ntinc)
endif
endif
end subroutine GFS_suite_interstitial_3_run
end module GFS_suite_interstitial_3
module GFS_suite_interstitial_4
contains
subroutine GFS_suite_interstitial_4_init ()
end subroutine GFS_suite_interstitial_4_init
subroutine GFS_suite_interstitial_4_finalize()
end subroutine GFS_suite_interstitial_4_finalize
!> \section arg_table_GFS_suite_interstitial_4_run Argument Table
!! \htmlinclude GFS_suite_interstitial_4_run.html
!!
subroutine GFS_suite_interstitial_4_run (im, levs, ltaerosol, tracers_total, ntrac, ntcw, ntiw, ntclamt, &
ntrw, ntsw, ntrnc, ntsnc, ntgl, ntgnc, ntlnc, ntinc, nn, imp_physics, imp_physics_gfdl, imp_physics_thompson, &
imp_physics_zhao_carr, imp_physics_zhao_carr_pdf, convert_dry_rho, dtf, save_qc, save_qi, con_pi, dtidx, dtend,&
index_of_process_conv_trans, gq0, clw, prsl, save_tcp, con_rd, con_eps, nwfa, spechum, ldiag3d, &
qdiag3d, save_lnc, save_inc, ntk, ntke, errmsg, errflg)
use machine, only: kind_phys
use module_mp_thompson_make_number_concentrations, only: make_IceNumber, make_DropletNumber
implicit none
! interface variables
integer, intent(in ) :: im, levs, tracers_total, ntrac, ntcw, ntiw, ntclamt, ntrw, &
ntsw, ntrnc, ntsnc, ntgl, ntgnc, ntlnc, ntinc, nn, imp_physics, imp_physics_gfdl, imp_physics_thompson, &
imp_physics_zhao_carr, imp_physics_zhao_carr_pdf
logical, intent(in) :: ltaerosol, convert_dry_rho
real(kind=kind_phys), intent(in ) :: con_pi, dtf
real(kind=kind_phys), intent(in ), dimension(:,:) :: save_qc
! save_qi is not allocated for Zhao-Carr MP
real(kind=kind_phys), intent(in ), dimension(:,:) :: save_qi, save_lnc, save_inc
! dtend and dtidx are only allocated if ldiag3d
logical, intent(in) :: ldiag3d, qdiag3d
real(kind=kind_phys), dimension(:,:,:), intent(inout) :: dtend
integer, dimension(:,:), intent(in) :: dtidx
integer, intent(in) :: index_of_process_conv_trans,ntk,ntke
real(kind=kind_phys), dimension(:,:,:), intent(inout) :: gq0
real(kind=kind_phys), dimension(:,:,:), intent(inout) :: clw
real(kind=kind_phys), dimension(:,:), intent(in) :: prsl
real(kind=kind_phys), intent(in) :: con_rd, con_eps
real(kind=kind_phys), dimension(:,:), intent(in) :: nwfa, save_tcp
real(kind=kind_phys), dimension(:,:), intent(in) :: spechum
character(len=*), intent( out) :: errmsg
integer, intent( out) :: errflg
! local variables
real(kind=kind_phys), parameter :: zero = 0.0_kind_phys, one = 1.0_kind_phys
integer :: i,k,n,tracers,idtend
real(kind=kind_phys) :: rho, orho
real(kind=kind_phys), dimension(im,levs) :: qv_mp !< kg kg-1 (dry mixing ratio)
real(kind=kind_phys), dimension(im,levs) :: qc_mp !< kg kg-1 (dry mixing ratio)
real(kind=kind_phys), dimension(im,levs) :: qi_mp !< kg kg-1 (dry mixing ratio)
real(kind=kind_phys), dimension(im,levs) :: nc_mp !< kg-1 (dry mixing ratio)
real(kind=kind_phys), dimension(im,levs) :: ni_mp !< kg-1 (dry mixing ratio)
! Initialize CCPP error handling variables
errmsg = ''
errflg = 0
! This code was previously in GFS_SCNV_generic_post, but it really belongs
! here, because it fixes the convective transportable_tracers mess for Zhao-Carr
! and GFDL MP from GFS_suite_interstitial_3. This whole code around clw(:,:,2)
! being set to -999 for Zhao-Carr MP (which doesn't have cloud ice) and GFDL-MP
! (which does have cloud ice, but for some reason it was decided to code it up
! in the same way as for Zhao-Carr, nowadays unnecessary and confusing) needs
! to be cleaned up. The convection schemes doing something different internally
! based on clw(i,k,2) being -999.0 or not is not a good idea.
do k=1,levs
do i=1,im
if (clw(i,k,2) <= -999.0) clw(i,k,2) = 0.0
enddo
enddo
if(ldiag3d) then
if(ntk>0 .and. ntk<=size(clw,3)) then
idtend=dtidx(100+ntke,index_of_process_conv_trans)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + clw(:,:,ntk)-gq0(:,:,ntk)
endif
endif
if(ntcw>0) then
if (imp_physics == imp_physics_zhao_carr .or. &
imp_physics == imp_physics_zhao_carr_pdf .or. &
imp_physics == imp_physics_gfdl) then
idtend=dtidx(100+ntcw,index_of_process_conv_trans)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + clw(:,:,1)+clw(:,:,2) - gq0(:,:,ntcw)
endif
else if(ntiw>0) then
idtend=dtidx(100+ntiw,index_of_process_conv_trans)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + clw(:,:,1)-gq0(:,:,ntiw)
endif
idtend=dtidx(100+ntcw,index_of_process_conv_trans)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + clw(:,:,2)-gq0(:,:,ntcw)
endif
else
idtend=dtidx(100+ntcw,index_of_process_conv_trans)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + clw(:,:,1)+clw(:,:,2) - gq0(:,:,ntcw)
endif
endif
endif
endif
! --- update the tracers due to deep & shallow cumulus convective transport
! (except for suspended water and ice)
if (tracers_total > 0) then
tracers = 2
do n=2,ntrac
! if ( n /= ntcw .and. n /= ntiw .and. n /= ntclamt) then
if ( n /= ntcw .and. n /= ntiw .and. n /= ntclamt .and. &
n /= ntrw .and. n /= ntsw .and. n /= ntrnc .and. &
n /= ntsnc .and. n /= ntgl .and. n /= ntgnc ) then
tracers = tracers + 1
if(n/=ntk .and. n/=ntlnc .and. n/=ntinc .and. n /= ntcw .and. n /= ntiw) then
idtend=dtidx(100+n,index_of_process_conv_trans)
if(idtend>=1) then
dtend(:,:,idtend) = dtend(:,:,idtend) + clw(:,:,tracers)-gq0(:,:,n)
endif
endif
do k=1,levs
do i=1,im
gq0(i,k,n) = clw(i,k,tracers)
enddo
enddo
endif
enddo
endif
if (ntcw > 0) then
! for microphysics
if (imp_physics == imp_physics_zhao_carr .or. &
imp_physics == imp_physics_zhao_carr_pdf .or. &
imp_physics == imp_physics_gfdl) then
gq0(1:im,:,ntcw) = clw(1:im,:,1) + clw(1:im,:,2)
elseif (ntiw > 0) then
do k=1,levs
do i=1,im
gq0(i,k,ntiw) = clw(i,k,1) ! ice
gq0(i,k,ntcw) = clw(i,k,2) ! water
enddo
enddo
if (imp_physics == imp_physics_thompson .and. (ntlnc>0 .or. ntinc>0)) then
if_convert_dry_rho: if (convert_dry_rho) then
do k=1,levs
do i=1,im
!> - Convert specific humidity to dry mixing ratio
qv_mp(i,k) = spechum(i,k) / (one-spechum(i,k))
!> - Density of air in kg m-3 and inverse density
rho = con_eps*prsl(i,k) / (con_rd*save_tcp(i,k)*(qv_mp(i,k)+con_eps))
orho = one/rho
if (ntlnc>0) then
!> - Convert moist mixing ratio to dry mixing ratio
qc_mp(i,k) = (clw(i,k,2)-save_qc(i,k)) / (one-spechum(i,k))
!> - Convert number concentration from moist to dry
nc_mp(i,k) = gq0(i,k,ntlnc) / (one-spechum(i,k))
nc_mp(i,k) = max(zero, nc_mp(i,k) + make_DropletNumber(qc_mp(i,k) * rho, nwfa(i,k)*rho) * orho)
!> - Convert number concentrations from dry to moist
gq0(i,k,ntlnc) = nc_mp(i,k) / (one+qv_mp(i,k))
endif
if (ntinc>0) then
!> - Convert moist mixing ratio to dry mixing ratio
qi_mp(i,k) = (clw(i,k,1)-save_qi(i,k)) / (one-spechum(i,k))
!> - Convert number concentration from moist to dry
ni_mp(i,k) = gq0(i,k,ntinc) / (one-spechum(i,k))
ni_mp(i,k) = max(zero, ni_mp(i,k) + make_IceNumber(qi_mp(i,k) * rho, save_tcp(i,k)) * orho)
!> - Convert number concentrations from dry to moist
gq0(i,k,ntinc) = ni_mp(i,k) / (one+qv_mp(i,k))
endif
enddo
enddo
else
do k=1,levs
do i=1,im
!> - Density of air in kg m-3 and inverse density
rho = con_eps*prsl(i,k) / (con_rd*save_tcp(i,k)*(spechum(i,k)+con_eps))
orho = one/rho
if (ntlnc>0) then
!> - Update cloud water mixing ratio
qc_mp(i,k) = (clw(i,k,2)-save_qc(i,k))
!> - Update cloud water number concentration
gq0(i,k,ntlnc) = max(zero, gq0(i,k,ntlnc) + make_DropletNumber(qc_mp(i,k) * rho, nwfa(i,k)*rho) * orho)
endif
if (ntinc>0) then
!> - Update cloud ice mixing ratio
qi_mp(i,k) = (clw(i,k,1)-save_qi(i,k))
!> - Update cloud ice number concentration
gq0(i,k,ntinc) = max(zero, gq0(i,k,ntinc) + make_IceNumber(qi_mp(i,k) * rho, save_tcp(i,k)) * orho)
endif
enddo
enddo
end if if_convert_dry_rho
if(ldiag3d .and. qdiag3d) then
idtend = dtidx(100+ntlnc,index_of_process_conv_trans)
if(idtend>0) then
dtend(:,:,idtend) = dtend(:,:,idtend) + gq0(:,:,ntlnc) - save_lnc
endif
idtend = dtidx(100+ntinc,index_of_process_conv_trans)
if(idtend>0) then
dtend(:,:,idtend) = dtend(:,:,idtend) + gq0(:,:,ntinc) - save_inc
endif
endif
endif
else
do k=1,levs
do i=1,im
gq0(i,k,ntcw) = clw(i,k,1) + clw(i,k,2)
enddo
enddo
endif ! end if_ntiw
else
do k=1,levs
do i=1,im
clw(i,k,1) = clw(i,k,1) + clw(i,k,2)
enddo
enddo
endif ! end if_ntcw
end subroutine GFS_suite_interstitial_4_run
end module GFS_suite_interstitial_4
module GFS_suite_interstitial_5
contains
subroutine GFS_suite_interstitial_5_init ()
end subroutine GFS_suite_interstitial_5_init
subroutine GFS_suite_interstitial_5_finalize()
end subroutine GFS_suite_interstitial_5_finalize
!> \section arg_table_GFS_suite_interstitial_5_run Argument Table
!! \htmlinclude GFS_suite_interstitial_5_run.html
!!
subroutine GFS_suite_interstitial_5_run (im, levs, ntrac, ntcw, ntiw, nn, gq0, clw, errmsg, errflg)
use machine, only: kind_phys
implicit none
! interface variables
integer, intent(in ) :: im, levs, ntrac, ntcw, ntiw, nn
real(kind=kind_phys), intent(in ), dimension(:,:,:) :: gq0
real(kind=kind_phys), intent(out), dimension(:,:,:) :: clw
character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
! local variables
integer :: i,k
! Initialize CCPP error handling variables
errmsg = ''
errflg = 0
do k=1,levs
do i=1,im
clw(i,k,1) = gq0(i,k,ntiw) ! ice
clw(i,k,2) = gq0(i,k,ntcw) ! water
enddo
enddo
end subroutine GFS_suite_interstitial_5_run
end module GFS_suite_interstitial_5