FV3 Bundle
type_cmat_blk.F90
Go to the documentation of this file.
1 !----------------------------------------------------------------------
2 ! Module: type_cmat_blk
3 ! Purpose: correlation matrix derived type
4 ! Author: Benjamin Menetrier
5 ! Licensing: this code is distributed under the CeCILL-C license
6 ! Copyright © 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
7 !----------------------------------------------------------------------
9 
10 use tools_kinds, only: kind_real
11 use tools_missing, only: msr
12 use type_bpar, only: bpar_type
13 use type_geom, only: geom_type
14 use type_nam, only: nam_type
15 
16 implicit none
17 
18 ! C matrix block data derived type
20  ! Block index and name
21  integer :: ib ! Block index
22  character(len=1024) :: name ! Name
23  logical :: double_fit ! Double fit
24  logical :: anisotropic ! Anisoptropic tensor
25 
26  ! Read data
27  real(kind_real),allocatable :: oops_coef_ens(:,:) ! OOPS ensemble coefficient
28  real(kind_real),allocatable :: oops_coef_sta(:,:) ! OOPS static coefficient
29  real(kind_real),allocatable :: oops_rh(:,:) ! OOPS horizontal fit support radius
30  real(kind_real),allocatable :: oops_rv(:,:) ! OOPS vertical fit support radius
31  real(kind_real),allocatable :: oops_rv_rfac(:,:) ! OOPS vertical fit support radius factor
32  real(kind_real),allocatable :: oops_rv_coef(:,:) ! OOPS vertical fit coefficient
33 
34  ! Data
35  real(kind_real),allocatable :: coef_ens(:,:) ! Ensemble coefficient
36  real(kind_real),allocatable :: coef_sta(:,:) ! Static coefficient
37  real(kind_real),allocatable :: rh(:,:) ! Horizontal fit support radius
38  real(kind_real),allocatable :: rv(:,:) ! Vertical fit support radius
39  real(kind_real),allocatable :: rv_rfac(:,:) ! Vertical fit support radius factor
40  real(kind_real),allocatable :: rv_coef(:,:) ! Vertical fit coefficient
41  real(kind_real),allocatable :: rhs(:,:) ! Fit support radius for sampling
42  real(kind_real),allocatable :: rvs(:,:) ! Fit support radius, for sampling
43  real(kind_real) :: wgt ! Block weight
44  real(kind_real),allocatable :: h11(:,:) ! LCT component 11
45  real(kind_real),allocatable :: h22(:,:) ! LCT component 22
46  real(kind_real),allocatable :: h33(:,:) ! LCT component 33
47  real(kind_real),allocatable :: h12(:,:) ! LCT component 12
48  real(kind_real),allocatable :: hcoef(:,:) ! LCT scales coefficients
49  real(kind_real),allocatable :: displ_lon(:,:,:) ! Displaced longitude
50  real(kind_real),allocatable :: displ_lat(:,:,:) ! Displaced latitude
51 contains
52  procedure :: alloc => cmat_blk_alloc
53  procedure :: dealloc => cmat_blk_dealloc
54 end type cmat_blk_type
55 
56 private
57 public :: cmat_blk_type
58 
59 contains
60 
61 !----------------------------------------------------------------------
62 ! Subroutine: cmat_blk_alloc
63 ! Purpose: C matrix block data allocation
64 !----------------------------------------------------------------------
65 subroutine cmat_blk_alloc(cmat_blk,nam,geom,bpar)
66 
67 implicit none
68 
69 ! Passed variables
70 class(cmat_blk_type),intent(inout) :: cmat_blk ! C matrix data block
71 type(nam_type),target,intent(in) :: nam ! Namelist
72 type(geom_type),target,intent(in) :: geom ! Geometry
73 type(bpar_type),intent(in) :: bpar ! Block parameters
74 
75 ! Associate
76 associate(ib=>cmat_blk%ib)
77 
78 if (bpar%diag_block(ib)) then
79  ! Allocation
80  allocate(cmat_blk%coef_ens(geom%nc0a,geom%nl0))
81  allocate(cmat_blk%coef_sta(geom%nc0a,geom%nl0))
82  allocate(cmat_blk%rh(geom%nc0a,geom%nl0))
83  allocate(cmat_blk%rv(geom%nc0a,geom%nl0))
84  if (cmat_blk%double_fit) then
85  allocate(cmat_blk%rv_rfac(geom%nc0a,geom%nl0))
86  allocate(cmat_blk%rv_coef(geom%nc0a,geom%nl0))
87  end if
88  allocate(cmat_blk%rhs(geom%nc0a,geom%nl0))
89  allocate(cmat_blk%rvs(geom%nc0a,geom%nl0))
90  if (cmat_blk%anisotropic) then
91  allocate(cmat_blk%H11(geom%nc0a,geom%nl0))
92  allocate(cmat_blk%H22(geom%nc0a,geom%nl0))
93  allocate(cmat_blk%H33(geom%nc0a,geom%nl0))
94  allocate(cmat_blk%H12(geom%nc0a,geom%nl0))
95  allocate(cmat_blk%Hcoef(geom%nc0a,geom%nl0))
96  end if
97 
98  ! Initialization
99  call msr(cmat_blk%coef_ens)
100  call msr(cmat_blk%coef_sta)
101  call msr(cmat_blk%rh)
102  call msr(cmat_blk%rv)
103  if (cmat_blk%double_fit) then
104  call msr(cmat_blk%rv_rfac)
105  call msr(cmat_blk%rv_coef)
106  end if
107  call msr(cmat_blk%rhs)
108  call msr(cmat_blk%rvs)
109  call msr(cmat_blk%wgt)
110  if (cmat_blk%anisotropic) then
111  call msr(cmat_blk%H11)
112  call msr(cmat_blk%H22)
113  call msr(cmat_blk%H33)
114  call msr(cmat_blk%H12)
115  call msr(cmat_blk%Hcoef)
116  end if
117 end if
118 
119 if ((ib==bpar%nbe).and.nam%displ_diag) then
120  ! Allocation
121  allocate(cmat_blk%displ_lon(geom%nc0a,geom%nl0,2:nam%nts))
122  allocate(cmat_blk%displ_lat(geom%nc0a,geom%nl0,2:nam%nts))
123 
124  ! Initialization
125  if (nam%displ_diag) then
126  call msr(cmat_blk%displ_lon)
127  call msr(cmat_blk%displ_lat)
128  end if
129 end if
130 
131 ! End associate
132 end associate
133 
134 end subroutine cmat_blk_alloc
135 
136 !----------------------------------------------------------------------
137 ! Subroutine: cmat_blk_dealloc
138 ! Purpose: C matrix block data deallocation
139 !----------------------------------------------------------------------
140 subroutine cmat_blk_dealloc(cmat_blk)
142 implicit none
143 
144 ! Passed variables
145 class(cmat_blk_type),intent(inout) :: cmat_blk ! C matrix data block
146 
147 ! Release memory
148 if (allocated(cmat_blk%oops_coef_ens)) deallocate(cmat_blk%oops_coef_ens)
149 if (allocated(cmat_blk%oops_coef_sta)) deallocate(cmat_blk%oops_coef_sta)
150 if (allocated(cmat_blk%oops_rh)) deallocate(cmat_blk%oops_rh)
151 if (allocated(cmat_blk%oops_rv)) deallocate(cmat_blk%oops_rv)
152 if (allocated(cmat_blk%oops_rv_rfac)) deallocate(cmat_blk%oops_rv_rfac)
153 if (allocated(cmat_blk%oops_rv_coef)) deallocate(cmat_blk%oops_rv_coef)
154 if (allocated(cmat_blk%coef_ens)) deallocate(cmat_blk%coef_ens)
155 if (allocated(cmat_blk%coef_sta)) deallocate(cmat_blk%coef_sta)
156 if (allocated(cmat_blk%rh)) deallocate(cmat_blk%rh)
157 if (allocated(cmat_blk%rv)) deallocate(cmat_blk%rv)
158 if (allocated(cmat_blk%rv_rfac)) deallocate(cmat_blk%rv_rfac)
159 if (allocated(cmat_blk%rv_coef)) deallocate(cmat_blk%rv_coef)
160 if (allocated(cmat_blk%rhs)) deallocate(cmat_blk%rhs)
161 if (allocated(cmat_blk%rvs)) deallocate(cmat_blk%rvs)
162 if (allocated(cmat_blk%H11)) deallocate(cmat_blk%H11)
163 if (allocated(cmat_blk%H22)) deallocate(cmat_blk%H22)
164 if (allocated(cmat_blk%H33)) deallocate(cmat_blk%H33)
165 if (allocated(cmat_blk%H12)) deallocate(cmat_blk%H12)
166 if (allocated(cmat_blk%Hcoef)) deallocate(cmat_blk%Hcoef)
167 if (allocated(cmat_blk%displ_lon)) deallocate(cmat_blk%displ_lon)
168 if (allocated(cmat_blk%displ_lat)) deallocate(cmat_blk%displ_lat)
169 
170 end subroutine cmat_blk_dealloc
171 
172 end module type_cmat_blk
subroutine cmat_blk_alloc(cmat_blk, nam, geom, bpar)
subroutine cmat_blk_dealloc(cmat_blk)
integer, parameter, public kind_real