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