CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
GFS_SCNV_generic_pre.F90
1
3
5
6 contains
7
11 subroutine gfs_scnv_generic_pre_run (im, levs, ldiag3d, qdiag3d, gu0, gv0, gt0, gq0, &
12 save_u, save_v, save_t, save_q, ntqv, nsamftrac, flag_for_scnv_generic_tend, &
13 dtidx, index_of_process_scnv, ntcw,ntiw,ntclamt,ntrw,ntsw,ntrnc,ntsnc,ntgl,ntgnc, &
14 ntsigma, cscnv, satmedmf, trans_trac, ras, ntrac, clw, errmsg, errflg)
15
16 use machine, only: kind_phys
17
18 implicit none
19
20 integer, intent(in) :: im, levs, ntqv, nsamftrac, index_of_process_scnv, dtidx(:,:)
21 integer, intent(in) :: ntcw,ntiw,ntclamt,ntrw,ntsw,ntrnc,ntsnc,ntgl,ntgnc, ntsigma,ntrac
22 logical, intent(in) :: ldiag3d, qdiag3d, flag_for_scnv_generic_tend
23 real(kind=kind_phys), dimension(:,:), intent(in) :: gu0, gv0, gt0
24 real(kind=kind_phys), dimension(:,:,:), intent(in) :: gq0
25 real(kind=kind_phys), dimension(:,:,:), intent(inout) :: save_q
26 real(kind=kind_phys), dimension(:,:), intent(inout) :: save_u, save_v, save_t
27 character(len=*), intent(out) :: errmsg
28 integer, intent(out) :: errflg
29 logical, intent(in) :: cscnv, satmedmf, trans_trac, ras
30 real(kind=kind_phys), dimension(:,:,:), intent(in) :: clw
31
32 integer :: i, k, n, tracers
33
34 ! Initialize CCPP error handling variables
35 errmsg = ''
36 errflg = 0
37
38 if (ldiag3d .and. flag_for_scnv_generic_tend) then
39 do k=1,levs
40 do i=1,im
41 save_u(i,k) = gu0(i,k)
42 save_v(i,k) = gv0(i,k)
43 save_t(i,k) = gt0(i,k)
44 enddo
45 enddo
46 if (qdiag3d) then
47 if (cscnv .or. satmedmf .or. trans_trac .or. ras) then
48 tracers = 2
49 do n=2,ntrac
50 if ( n /= ntcw .and. n /= ntiw .and. n /= ntclamt .and. &
51 n /= ntrw .and. n /= ntsw .and. n /= ntrnc .and. &
52 n /= ntsnc .and. n /= ntgl .and. n /= ntgnc .and. n /= ntsigma) then
53 tracers = tracers + 1
54 if(dtidx(100+n,index_of_process_scnv)>0) then
55 save_q(:,:,n) = clw(:,:,tracers)
56 endif
57 endif
58 enddo
59 else
60 do n=2,ntrac
61 if(dtidx(100+n,index_of_process_scnv)>0) then
62 save_q(:,:,n) = gq0(:,:,n)
63 endif
64 enddo
65 endif ! end if_ras or cfscnv or samf
66 save_q(:,:,ntqv) = gq0(:,:,ntqv)
67 endif
68 endif
69
70 end subroutine gfs_scnv_generic_pre_run
71
72
73 end module gfs_scnv_generic_pre