CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
cu_gf_driver_pre.F90
1
3
6
7 implicit none
8
9 private
10
11 public :: cu_gf_driver_pre_run
12
13 contains
14
19 subroutine cu_gf_driver_pre_run (flag_init, flag_restart, gf_coldstart, kdt, fhour, dtp, t, q, prevst, prevsq, &
20 forcet, forceq, cactiv, cactiv_m, conv_act, conv_act_m, &
21 ntsmoke, ntdust, ntcoarsepm, chem3d, gq0, errmsg, errflg)
22
23 use machine, only: kind_phys
24
25 implicit none
26
27 logical, intent(in) :: flag_init
28 logical, intent(in) :: flag_restart
29 logical, intent(in) :: gf_coldstart
30 integer, intent(in) :: kdt
31 real(kind_phys), intent(in) :: fhour
32 real(kind_phys), intent(in) :: dtp
33 real(kind_phys), intent(in) :: t(:,:)
34 real(kind_phys), intent(in) :: q(:,:)
35 real(kind_phys), intent(in) :: prevst(:,:)
36 real(kind_phys), intent(in) :: prevsq(:,:)
37!$acc declare copyin(t,q,prevst,prevsq)
38 real(kind_phys), intent(out) :: forcet(:,:)
39 real(kind_phys), intent(out) :: forceq(:,:)
40 integer, intent(out) :: cactiv(:)
41 integer, intent(out) :: cactiv_m(:)
42 integer, intent(in) :: ntsmoke, ntdust, ntcoarsepm
43!$acc declare copyout(forcet,forceq,cactiv,cactiv_m)
44 real(kind_phys), intent(in) :: conv_act(:)
45 real(kind_phys), intent(in) :: conv_act_m(:)
46 real(kind_phys), intent(inout), optional :: chem3d(:,:,:)
47 real(kind_phys), intent(inout) :: gq0(:,:,:)
48!$acc declare copyin(conv_act,conv_act_m) copy(chem3d,gq0)
49 character(len=*), intent(out) :: errmsg
50 integer, intent(out) :: errflg
51
52 ! local variables
53 real(kind=kind_phys) :: dtdyn
54
55 ! Initialize CCPP error handling variables
56 errmsg = ''
57 errflg = 0
58
59 ! For restart runs, can assume that prevst and prevsq
60 ! are read from the restart files beforehand, same
61 ! for conv_act.
62 if((flag_init .and. .not.flag_restart) .or. gf_coldstart) then
63!$acc kernels
64 forcet(:,:)=0.0
65 forceq(:,:)=0.0
66!$acc end kernels
67 else
68 dtdyn=3600.0*(fhour)/kdt
69 if(dtp > dtdyn) then
70!$acc kernels
71 forcet(:,:)=(t(:,:) - prevst(:,:))/dtp
72 forceq(:,:)=(q(:,:) - prevsq(:,:))/dtp
73!$acc end kernels
74 else
75!$acc kernels
76 forcet(:,:)=(t(:,:) - prevst(:,:))/dtdyn
77 forceq(:,:)=(q(:,:) - prevsq(:,:))/dtdyn
78!$acc end kernels
79 endif
80 endif
81
82!$acc kernels
83 cactiv(:)=nint(conv_act(:))
84 cactiv_m(:)=nint(conv_act_m(:))
85
86 if (present(chem3d)) then
87 chem3d(:,:,1) = gq0(:,:,ntsmoke)
88 chem3d(:,:,2) = gq0(:,:,ntdust)
89 chem3d(:,:,3) = gq0(:,:,ntcoarsepm)
90 endif
91!$acc end kernels
92
93 end subroutine cu_gf_driver_pre_run
94
95end module cu_gf_driver_pre
subroutine, public cu_gf_driver_pre_run(flag_init, flag_restart, gf_coldstart, kdt, fhour, dtp, t, q, prevst, prevsq, forcet, forceq, cactiv, cactiv_m, conv_act, conv_act_m, ntsmoke, ntdust, ntcoarsepm, chem3d, gq0, errmsg, errflg)
This module contains code related to GF convective schemes to be used within the GFS physics suite.