FV3 Bundle
fv_arrays_tlmadm.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU General Public License *
3 !* This file is a part of fvGFS. *
4 !* *
5 !* fvGFS is free software; you can redistribute it and/or modify it *
6 !* and are expected to follow the terms of the GNU General Public *
7 !* License as published by the Free Software Foundation; either *
8 !* version 2 of the License, or (at your option) any later version. *
9 !* *
10 !* fvGFS is distributed in the hope that it will be useful, but *
11 !* WITHOUT ANY WARRANTY; without even the implied warranty of *
12 !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
13 !* General Public License for more details. *
14 !* *
15 !* For the full text of the GNU General Public License, *
16 !* write to: Free Software Foundation, Inc., *
17 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
18 !* or see: http://www.gnu.org/licenses/gpl.html *
19 !***********************************************************************
21 
22 #include <fms_platform.h>
23 
24 implicit none
25 public
26 
27 logical :: fv_timing_onoff
28 
29 !Compiler definitions that are regular if statements for the tlm/adm, should still match how the code is compiled
31  logical :: fpp_mapl_mode = .true. !GMAO Mode
32  logical :: fpp_overload_r4 = .false.
33 end type fpp_type
34 
35 type(fpp_type) :: fpp
36 
38 
39 ! Traj = pert, true is faster
40  logical :: split_hord = .true. !Force traj to use pert options if false
41 
42 ! TLM/ADM options for horizontal transport
43  integer :: hord_mt_pert = 2 ! Linear options:
44  integer :: hord_vt_pert = 2 ! 1 (all) first order scheme
45  integer :: hord_tm_pert = 2 ! 6 (all) quasi fifth order linear scheme
46  integer :: hord_dp_pert = 2 ! -5 (vt,tm,dp,tr) linear scheme based on PPM 4th order interpolation
47  integer :: hord_tr_pert = 2 ! 333 (vt,tm,dp,tr) linear third order cheme
48 
49 ! Number of sponge layers for the perturbations
50  integer :: n_sponge_pert = 9
51 
52 !Perturbation damping in the sponge
53  logical :: hord_ks_pert = .true. !Use the below for the trajectory in the sponge layer
54  integer :: hord_mt_ks_pert = 1
55  integer :: hord_vt_ks_pert = 1
56  integer :: hord_tm_ks_pert = 1
57  integer :: hord_dp_ks_pert = 1
58  integer :: hord_tr_ks_pert = 1
59 
60 !Trajectory damping in the sponge
61  logical :: hord_ks_traj = .true. !Use the below for the trajectory in the sponge layer
62  integer :: hord_mt_ks_traj = 1
63  integer :: hord_vt_ks_traj = 1
64  integer :: hord_tm_ks_traj = 1
65  integer :: hord_dp_ks_traj = 1
66  integer :: hord_tr_ks_traj = 1
67 
68 ! TLM/ADM options for vertical remapping
69  logical :: split_kord = .true. !Force traj to use pert options if false
70  integer :: kord_mt_pert = 17 ! |kord|>16 is linear remapping, not shape preserving or
71  integer :: kord_wz_pert = 17 ! positive definite, anything else not recommended
72  integer :: kord_tm_pert = 17 ! for perturbation quantities
73  integer :: kord_tr_pert = 17
74 
75 ! TLM/ADM options for damping
76  logical :: split_damp = .true. !Force traj to use pert options if false
77  integer :: nord_pert = 1
78  real :: dddmp_pert = 0.2
79  real :: d2_bg_pert = 0.015
80  real :: d4_bg_pert = 0.150
81  logical :: do_vort_damp_pert = .true.
82  real :: vtdm4_pert = 0.0005
83  real :: d2_bg_k1_pert = 4. ! factor for d2_bg (k=1)
84  real :: d2_bg_k2_pert = 2. ! factor for d2_bg (k=2)
85  real :: d2_bg_ks_pert = 2. ! factor for d2_bg (k=1)
86 
87 ! TLM/ADM options for tracer damping
88  logical :: split_damp_tr = .true. !Force traj to use pert options if false
89  integer :: nord_tr_pert=0 ! Damping order
90  real :: trdm2_pert = 0.0 ! coefficient for damping
91 
92 end type fv_flags_pert_type
93 
95 
96  type(fv_flags_pert_type) :: flagstruct
97 
98  real, allocatable :: up(:,:,:) ! D grid zonal wind (m/s)
99  real, allocatable :: vp(:,:,:) ! D grid meridional wind (m/s)
100  real, allocatable :: ptp(:,:,:) ! temperature (K)
101  real, allocatable :: delpp(:,:,:) ! pressure thickness (pascal)
102  real, allocatable :: qp(:,:,:,:) ! specific humidity and constituents
103 
104  real, allocatable :: wp(:,:,:) ! cell center vertical wind (m/s)
105  real, allocatable :: delzp(:,:,:) ! layer thickness (meters)
106  real, allocatable :: ze0p(:,:,:) ! height at layer edges for remapping
107  real, allocatable :: q_conp(:,:,:) ! total condensates
108 
109  real, allocatable :: psp(:,:) ! Surface pressure (pascal)
110  real, allocatable :: pep(:,:,:) ! edge pressure (pascal)
111  real, allocatable :: pkp(:,:,:) ! pe**cappa
112  real, allocatable :: pelnp(:,:,:) ! ln(pe)
113  real, allocatable :: pkzp(:,:,:) ! finite-volume mean pk
114 
115  real, allocatable :: omgap(:,:,:) ! Vertical pressure velocity (pa/s)
116 
117  real, allocatable :: uap(:,:,:) ! (ua, va) are mostly used as the A grid winds
118  real, allocatable :: vap(:,:,:)
119  real, allocatable :: ucp(:,:,:) ! (uc, vc) are mostly used as the C grid winds
120  real, allocatable :: vcp(:,:,:)
121 
122  real, allocatable :: mfxp(:,:,:) ! Mass fluxes
123  real, allocatable :: mfyp(:,:,:)
124 
125  real, allocatable :: cxp(:,:,:)
126  real, allocatable :: cyp(:,:,:)
127 
128 end type fv_atmos_pert_type
129 
130 contains
131 
132 subroutine allocate_fv_atmos_pert_type(AtmP, isd, ied, jsd, jed, is, ie, js, je, npz, ncnst)
134  implicit none
135  type(fv_atmos_pert_type), intent(INOUT), target :: AtmP
136  integer, intent(IN) :: isd, ied, jsd, jed, is, ie, js, je, npz, ncnst
137 
138  ! Allocate the perturbation variables
139  allocate ( atmp%up(isd:ied ,jsd:jed+1,npz) )
140  allocate ( atmp%vp(isd:ied+1,jsd:jed ,npz) )
141  allocate ( atmp%ptp(isd:ied ,jsd:jed ,npz) )
142  allocate ( atmp%delpp(isd:ied ,jsd:jed ,npz) )
143  allocate ( atmp%qp(isd:ied ,jsd:jed ,npz, ncnst) )
144  allocate ( atmp%wp(isd:ied, jsd:jed ,npz ) )
145  allocate ( atmp%delzp(isd:ied, jsd:jed ,npz) )
146  allocate ( atmp%ze0p(is:ie , js:je ,npz+1) )
147  allocate (atmp%q_conp(isd:ied,jsd:jed,1:npz) )
148  allocate ( atmp%psp(isd:ied ,jsd:jed) )
149  allocate ( atmp%pep(is-1:ie+1, npz+1,js-1:je+1) )
150  allocate ( atmp%pkp(is:ie ,js:je , npz+1) )
151  allocate ( atmp%pelnp(is:ie,npz+1,js:je) )
152  allocate ( atmp%pkzp(is:ie,js:je,npz) )
153  allocate ( atmp%omgap(isd:ied ,jsd:jed ,npz) )
154  allocate ( atmp%uap(isd:ied ,jsd:jed ,npz) )
155  allocate ( atmp%vap(isd:ied ,jsd:jed ,npz) )
156  allocate ( atmp%ucp(isd:ied+1,jsd:jed ,npz) )
157  allocate ( atmp%vcp(isd:ied ,jsd:jed+1,npz) )
158  allocate ( atmp%mfxp(is:ie+1, js:je, npz) )
159  allocate ( atmp%mfyp(is:ie , js:je+1,npz) )
160  allocate ( atmp%cxp(is:ie+1, jsd:jed, npz) )
161  allocate ( atmp%cyp(isd:ied ,js:je+1, npz) )
162 
163  atmp%up = 0.0
164  atmp%vp = 0.0
165  atmp%ptp = 0.0
166  atmp%delpp = 0.0
167  atmp%qp = 0.0
168  atmp%wp = 0.0
169  atmp%delzp = 0.0
170  atmp%ze0p = 0.0
171  atmp%q_conp = 0.0
172  atmp%psp = 0.0
173  atmp%pep = 0.0
174  atmp%pkp = 0.0
175  atmp%pelnp = 0.0
176  atmp%pkzp = 0.0
177  atmp%omgap = 0.0
178  atmp%uap = 0.0
179  atmp%vap = 0.0
180  atmp%ucp = 0.0
181  atmp%vcp = 0.0
182  atmp%mfxp = 0.0
183  atmp%mfyp = 0.0
184  atmp%cxp = 0.0
185  atmp%cyp = 0.0
186 
187 end subroutine allocate_fv_atmos_pert_type
188 
189 subroutine deallocate_fv_atmos_pert_type(AtmP)
191  implicit none
192  type(fv_atmos_pert_type), intent(INOUT) :: AtmP
193 
194  deallocate ( atmp%up )
195  deallocate ( atmp%vp )
196  deallocate ( atmp%ptp )
197  deallocate ( atmp%delpp )
198  deallocate ( atmp%qp )
199  deallocate ( atmp%wp )
200  deallocate ( atmp%delzp )
201  deallocate ( atmp%ze0p )
202  deallocate ( atmp%psp )
203  deallocate ( atmp%pep )
204  deallocate ( atmp%pkp )
205  deallocate ( atmp%pelnp )
206  deallocate ( atmp%pkzp )
207  deallocate ( atmp%omgap )
208  deallocate ( atmp%uap )
209  deallocate ( atmp%vap )
210  deallocate ( atmp%ucp )
211  deallocate ( atmp%vcp )
212  deallocate ( atmp%mfxp )
213  deallocate ( atmp%mfyp )
214  deallocate ( atmp%cxp )
215  deallocate ( atmp%cyp )
216 
217 end subroutine deallocate_fv_atmos_pert_type
218 
219 end module fv_arrays_tlmadm_mod
subroutine allocate_fv_atmos_pert_type(AtmP, isd, ied, jsd, jed, is, ie, js, je, npz, ncnst)
subroutine deallocate_fv_atmos_pert_type(AtmP)
integer, parameter, public up