FV3 Bundle
type_ens.F90
Go to the documentation of this file.
1 !----------------------------------------------------------------------
2 ! Module: type_ens
3 ! Purpose: ensemble 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 !----------------------------------------------------------------------
8 module type_ens
9 
10 use tools_kinds, only: kind_real
11 use tools_missing, only: msi,msr
12 use type_geom, only: geom_type
13 use type_mpl, only: mpl_type
14 use type_nam, only: nam_type
15 
16 implicit none
17 
18 ! Ensemble derived type
20  ! Attributes
21  integer :: ne ! Ensemble size
22  integer :: nsub ! Number of sub-ensembles
23 
24  ! Data
25  real(kind_real),allocatable :: fld(:,:,:,:,:) ! Ensemble perturbation
26  real(kind_real),allocatable :: mean(:,:,:,:,:) ! Ensemble mean
27 contains
28  procedure :: alloc => ens_alloc
29  procedure :: dealloc => ens_dealloc
30  procedure :: copy => ens_copy
31  procedure :: remove_mean => ens_remove_mean
32  procedure :: ens_from
33  procedure :: ens_from_nemovar
34  generic :: from => ens_from,ens_from_nemovar
35 end type ens_type
36 
37 private
38 public :: ens_type
39 
40 contains
41 
42 !----------------------------------------------------------------------
43 ! Subroutine: ens_alloc
44 ! Purpose: ensemble data allocation
45 !----------------------------------------------------------------------
46 subroutine ens_alloc(ens,nam,geom,ne,nsub)
47 
48 implicit none
49 
50 ! Passed variables
51 class(ens_type),intent(inout) :: ens ! Ensemble
52 type(nam_type),intent(in) :: nam ! Namelist
53 type(geom_type),intent(in) :: geom ! Geometry
54 integer,intent(in) :: ne ! Ensemble size
55 integer,intent(in) :: nsub ! Number of sub-ensembles
56 
57 ! Allocate
58 if (ne>0) then
59  allocate(ens%fld(geom%nc0a,geom%nl0,nam%nv,nam%nts,ne))
60  allocate(ens%mean(geom%nc0a,geom%nl0,nam%nv,nam%nts,nsub))
61 end if
62 
63 ! Initialization
64 ens%ne = ne
65 ens%nsub = nsub
66 if (ens%ne>0) then
67  call msr(ens%fld)
68  call msr(ens%mean)
69 end if
70 
71 end subroutine ens_alloc
72 
73 !----------------------------------------------------------------------
74 ! Subroutine: ens_dealloc
75 ! Purpose: ensemble data deallocation
76 !----------------------------------------------------------------------
77 subroutine ens_dealloc(ens)
78 
79 implicit none
80 
81 ! Passed variables
82 class(ens_type),intent(inout) :: ens ! Ensemble
83 
84 ! Release memory
85 if (allocated(ens%fld)) deallocate(ens%fld)
86 if (allocated(ens%mean)) deallocate(ens%mean)
87 
88 end subroutine ens_dealloc
89 
90 !----------------------------------------------------------------------
91 ! Subroutine: ens_copy
92 ! Purpose: ensemble data copy
93 !----------------------------------------------------------------------
94 subroutine ens_copy(ens_out,ens_in)
95 
96 implicit none
97 
98 ! Passed variables
99 class(ens_type),intent(inout) :: ens_out ! Ensemble
100 type(ens_type),intent(in) :: ens_in ! Ensemble
101 
102 ! Allocate
103 if (ens_in%ne>0) then
104  allocate(ens_out%fld(size(ens_in%fld,1),size(ens_in%fld,2),size(ens_in%fld,3),size(ens_in%fld,4),size(ens_in%fld,5)))
105  allocate(ens_out%mean(size(ens_in%mean,1),size(ens_in%mean,2),size(ens_in%mean,3),size(ens_in%mean,4),size(ens_in%mean,5)))
106 end if
107 
108 ! Initialization
109 ens_out%ne = ens_in%ne
110 ens_out%nsub = ens_in%nsub
111 if (ens_in%ne>0) then
112  ens_out%fld = ens_in%fld
113  ens_out%mean = ens_in%mean
114 end if
115 
116 end subroutine ens_copy
117 
118 !----------------------------------------------------------------------
119 ! Subroutine: ens_remove_mean
120 ! Purpose: remove ensemble mean
121 !----------------------------------------------------------------------
122 subroutine ens_remove_mean(ens)
124 implicit none
125 
126 ! Passed variables
127 class(ens_type),intent(inout) :: ens ! Ensemble
128 
129 ! Local variables
130 integer :: isub,ie_sub,ie
131 
132 if (ens%ne>0) then
133  ! Loop over sub-ensembles
134  do isub=1,ens%nsub
135  ! Compute mean
136  ens%mean(:,:,:,:,isub) = 0.0
137  do ie_sub=1,ens%ne/ens%nsub
138  ie = ie_sub+(isub-1)*ens%ne/ens%nsub
139  ens%mean(:,:,:,:,isub) = ens%mean(:,:,:,:,isub)+ens%fld(:,:,:,:,ie)
140  end do
141  ens%mean(:,:,:,:,isub) = ens%mean(:,:,:,:,isub)/(ens%ne/ens%nsub)
142 
143  ! Remove mean
144  do ie_sub=1,ens%ne/ens%nsub
145  ie = ie_sub+(isub-1)*ens%ne/ens%nsub
146  ens%fld(:,:,:,:,ie) = ens%fld(:,:,:,:,ie)-ens%mean(:,:,:,:,isub)
147  end do
148  end do
149 end if
150 
151 end subroutine ens_remove_mean
152 
153 !----------------------------------------------------------------------
154 ! Subroutine: ens_from
155 ! Purpose: copy ensemble array into ensemble data
156 !----------------------------------------------------------------------
157 subroutine ens_from(ens,nam,geom,ne,ens_mga)
159 implicit none
160 
161 ! Passed variables
162 class(ens_type),intent(inout) :: ens ! Ensemble
163 type(nam_type),intent(in) :: nam ! Namelist
164 type(geom_type),intent(in) :: geom ! Geometry
165 integer,intent(in) :: ne ! Ensemble size
166 real(kind_real),intent(in) :: ens_mga(geom%nmga,geom%nl0,nam%nv,nam%nts,ne) ! Ensemble on model grid, halo A
167 
168 ! Local variables
169 integer :: ie,its,iv,il0
170 
171 ! Allocation
172 call ens%alloc(nam,geom,ne,1)
173 
174 if (ens%ne>0) then
175  ! Copy
176  do ie=1,ens%ne
177  do its=1,nam%nts
178  do iv=1,nam%nv
179  do il0=1,geom%nl0
180  ens%fld(:,il0,iv,its,ie) = ens_mga(geom%c0a_to_mga,il0,iv,its,ie)
181  end do
182  end do
183  end do
184  end do
185 
186  ! Remove mean
187  call ens%remove_mean
188 end if
189 
190 end subroutine ens_from
191 
192 !----------------------------------------------------------------------
193 ! Subroutine: ens_from_nemovar
194 ! Purpose: copy 2d NEMOVAR ensemble into ensemble data
195 !----------------------------------------------------------------------
196 subroutine ens_from_nemovar(ens,mpl,nam,geom,nx,ny,nens,ncyc,ens_2d,ens_3d)
198 implicit none
199 
200 ! Passed variables
201 class(ens_type),intent(inout) :: ens ! Ensemble
202 type(mpl_type),intent(in) :: mpl ! MPI data
203 type(nam_type),intent(in) :: nam ! Namelist
204 type(geom_type),intent(in) :: geom ! Geometry
205 integer,intent(in) :: nx ! X-axis size
206 integer,intent(in) :: ny ! Y-axis size
207 integer,intent(in) :: nens ! Ensemble size at each cycle
208 integer,intent(in) :: ncyc ! Number of cycles
209 real(kind_real),intent(in),optional :: ens_2d(nx,ny,nens,ncyc) ! Ensemble on model grid, halo A
210 real(kind_real),intent(in),optional :: ens_3d(nx,ny,geom%nl0,nens,ncyc) ! Ensemble on model grid, halo A
211 
212 ! Local variables
213 integer :: ie,iens,icyc,its,iv,il0
214 real(kind_real) :: tmp(geom%nmga)
215 
216 ! Allocation
217 call ens%alloc(nam,geom,nens*ncyc,ncyc)
218 
219 ! Copy
220 iv = 1
221 its = 1
222 ie = 1
223 do iens=1,nens
224  do icyc=1,ncyc
225  ! Pack and copy
226  if (present(ens_2d)) then
227  il0 = 1
228  tmp = pack(ens_2d(:,:,iens,icyc),.true.)
229  ens%fld(:,il0,iv,its,ie) = tmp(geom%c0a_to_mga)
230  elseif (present(ens_3d)) then
231  do il0=1,geom%nl0
232  ! Pack and copy
233  tmp = pack(ens_3d(:,:,il0,iens,icyc),.true.)
234  ens%fld(:,il0,iv,its,ie) = tmp(geom%c0a_to_mga)
235  end do
236  else
237  call mpl%abort('ens_2d or ens_3d should be provided in ens_from_nemovar')
238  end if
239 
240  ! Update
241  ie = ie+1
242  end do
243 end do
244 
245 ! Remove mean
246 call ens%remove_mean
247 
248 end subroutine ens_from_nemovar
249 
250 end module type_ens
subroutine ens_copy(ens_out, ens_in)
Definition: type_ens.F90:95
subroutine ens_remove_mean(ens)
Definition: type_ens.F90:123
subroutine, public copy(self, rhs)
Definition: qg_fields.F90:225
subroutine ens_alloc(ens, nam, geom, ne, nsub)
Definition: type_ens.F90:47
subroutine ens_from(ens, nam, geom, ne, ens_mga)
Definition: type_ens.F90:158
subroutine ens_dealloc(ens)
Definition: type_ens.F90:78
subroutine ens_from_nemovar(ens, mpl, nam, geom, nx, ny, nens, ncyc, ens_2d, ens_3d)
Definition: type_ens.F90:197
integer, parameter, public kind_real