FV3 Bundle
type_bpar.F90
Go to the documentation of this file.
1 !----------------------------------------------------------------------
2 ! Module: type_bpar
3 ! Purpose: block parameters 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_bpar
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_nam, only: nam_type
14 
15 implicit none
16 
18  ! Block parameters
19  integer :: nb ! Number of blocks
20  integer :: nbe ! Extended number of blocks
21  integer :: nl0rmax ! Maximum effective number of levels
22  integer,allocatable :: nl0r(:) ! Effective number of levels
23  integer,allocatable :: l0rl0b_to_l0(:,:,:) ! Effective level to level
24  integer,allocatable :: il0rz(:,:) ! Effective zero separation level
25  integer,allocatable :: nc3(:) ! Maximum class
26  logical,allocatable :: vbal_block(:,:) ! Vertical balance block
27  logical,allocatable :: diag_block(:) ! HDIAG block
28  logical,allocatable :: avg_block(:) ! Averaging block
29  logical,allocatable :: fit_block(:) ! Fit block
30  logical,allocatable :: b_block(:) ! B-involved block
31  logical,allocatable :: nicas_block(:) ! NICAS block
32  logical,allocatable :: cv_block(:) ! Control variable block
33  character(len=11),allocatable :: blockname(:) ! Block name
34  integer,allocatable :: b_to_v1(:) ! Block to first variable
35  integer,allocatable :: b_to_v2(:) ! Block to second variable
36  integer,allocatable :: b_to_ts1(:) ! Block to first timeslot
37  integer,allocatable :: b_to_ts2(:) ! Block to second timeslot
38 contains
39  procedure :: alloc => bpar_alloc
40  procedure :: dealloc => bpar_dealloc
41 end type bpar_type
42 
43 private
44 public :: bpar_type
45 
46 contains
47 
48 !----------------------------------------------------------------------
49 ! Subroutine: bpar_alloc
50 ! Purpose: allocate general parameters
51 !----------------------------------------------------------------------
52 subroutine bpar_alloc(bpar,nam,geom)
53 
54 implicit none
55 
56 ! Passed variable
57 class(bpar_type),intent(inout) :: bpar ! Block parameters
58 type(nam_type),intent(in) :: nam ! Namelist
59 type(geom_type),intent(in) :: geom ! Geometry
60 
61 ! Local variables
62 integer :: ib,iv,jv,its,jts,il0,jl0r,jl0off
63 
64 ! Number of blocks
65 if (nam%new_lct) then
66  ! Special case for tensors
67  bpar%nb = nam%nv*nam%nts
68  bpar%nbe = bpar%nb
69 else
70  ! General case
71  bpar%nb = nam%nv**2*nam%nts**2
72  if (bpar%nb==1) then
73  bpar%nbe = bpar%nb
74  else
75  bpar%nbe = bpar%nb+1
76  end if
77 end if
78 
79 ! Allocation
80 bpar%nl0rmax = min(nam%nl0r,geom%nl0)
81 allocate(bpar%l0rl0b_to_l0(bpar%nl0rmax,geom%nl0,bpar%nbe))
82 allocate(bpar%il0rz(geom%nl0,bpar%nbe))
83 allocate(bpar%nl0r(bpar%nbe))
84 allocate(bpar%nc3(bpar%nbe))
85 allocate(bpar%vbal_block(nam%nv,nam%nv))
86 allocate(bpar%diag_block(bpar%nbe))
87 allocate(bpar%avg_block(bpar%nbe))
88 allocate(bpar%fit_block(bpar%nbe))
89 allocate(bpar%B_block(bpar%nbe))
90 allocate(bpar%nicas_block(bpar%nbe))
91 allocate(bpar%cv_block(bpar%nbe))
92 allocate(bpar%blockname(bpar%nbe))
93 allocate(bpar%b_to_v1(bpar%nbe))
94 allocate(bpar%b_to_v2(bpar%nbe))
95 allocate(bpar%b_to_ts1(bpar%nbe))
96 allocate(bpar%b_to_ts2(bpar%nbe))
97 
98 ! Initialization
99 call msi(bpar%l0rl0b_to_l0)
100 call msi(bpar%il0rz)
101 
102 if (nam%new_lct) then
103  ! Individual blocks
104  ib = 1
105  its = 1
106  do iv=1,nam%nv
107  ! Classes and levels
108  bpar%nl0r(ib) = bpar%nl0rmax
109  do il0=1,geom%nl0
110  jl0off = il0-(bpar%nl0r(ib)-1)/2-1
111  if (jl0off<1) jl0off = 0
112  if (jl0off+bpar%nl0rmax>geom%nl0) jl0off = geom%nl0-bpar%nl0rmax
113  do jl0r=1,bpar%nl0rmax
114  bpar%l0rl0b_to_l0(jl0r,il0,ib) = jl0off+jl0r
115  if (bpar%l0rl0b_to_l0(jl0r,il0,ib)==il0) bpar%il0rz(il0,ib) = jl0r
116  end do
117  end do
118  bpar%nc3(ib) = nam%nc3
119  bpar%nc3(ib) = nam%nc3
120  do jv=1,nam%nv
121  bpar%vbal_block(iv,jv) = (iv>1).and.(jv<iv).and.nam%vbal_block((iv-1)*(iv-2)/2+jv)
122  end do
123  bpar%diag_block(ib) = .true.
124  bpar%avg_block(ib) = .false.
125  bpar%fit_block(ib) = .false.
126  bpar%B_block(ib) = .true.
127  bpar%nicas_block(ib) = .true.
128  bpar%cv_block(ib) = .true.
129 
130  ! Blocks information
131  write(bpar%blockname(ib),'(i2.2,a,i2.2,a,i2.2,a,i2.2)') iv,'_',iv,'_',its,'_',its
132  bpar%b_to_v1(ib) = iv
133  bpar%b_to_v2(ib) = iv
134  bpar%b_to_ts1(ib) = its
135  bpar%b_to_ts2(ib) = its
136 
137  ! Update block index
138  ib = ib+1
139  end do
140 else
141  ! Individual blocks
142  ib = 1
143  do iv=1,nam%nv
144  do jv=1,nam%nv
145  do its=1,nam%nts
146  do jts=1,nam%nts
147  ! Classes and levels
148  if ((trim(nam%strategy)=='diag_all').or.((iv==jv).and.(its==jts))) then
149  bpar%nl0r(ib) = bpar%nl0rmax
150  do il0=1,geom%nl0
151  jl0off = il0-(bpar%nl0r(ib)-1)/2-1
152  if (jl0off<1) jl0off = 0
153  if (jl0off+bpar%nl0rmax>geom%nl0) jl0off = geom%nl0-bpar%nl0rmax
154  do jl0r=1,bpar%nl0rmax
155  bpar%l0rl0b_to_l0(jl0r,il0,ib) = jl0off+jl0r
156  if (bpar%l0rl0b_to_l0(jl0r,il0,ib)==il0) bpar%il0rz(il0,ib) = jl0r
157  end do
158  end do
159  bpar%nc3(ib) = nam%nc3
160  else
161  do il0=1,geom%nl0
162  bpar%l0rl0b_to_l0(:,il0,ib) = il0
163  end do
164  bpar%il0rz(:,ib) = 1
165  bpar%nl0r(ib) = 1
166  bpar%nc3(ib) = 1
167  end if
168 
169  ! Select blocks
170  bpar%vbal_block(iv,jv) = (iv>1).and.(jv<iv).and.nam%vbal_block((iv-1)*(iv-2)/2+jv)
171  select case (nam%strategy)
172  case ('diag_all')
173  bpar%diag_block(ib) = .true.
174  bpar%avg_block(ib) = .true.
175  bpar%B_block(ib) = .false.
176  bpar%nicas_block(ib) = .false.
177  bpar%cv_block(ib) = .false.
178  case ('common')
179  bpar%diag_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
180  bpar%avg_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
181  bpar%B_block(ib) = (ib==bpar%nbe)
182  bpar%nicas_block(ib) = (ib==bpar%nbe)
183  bpar%cv_block(ib) = (ib==bpar%nbe)
184  case ('common_univariate')
185  bpar%diag_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
186  bpar%avg_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
187  bpar%B_block(ib) = (ib==bpar%nbe)
188  bpar%nicas_block(ib) = (ib==bpar%nbe)
189  bpar%cv_block(ib) = (iv==jv).and.(its==jts)
190  case ('common_weighted')
191  bpar%diag_block(ib) = .true.
192  bpar%avg_block(ib) = (iv==jv).and.(its==jts)
193  bpar%B_block(ib) = .true.
194  bpar%nicas_block(ib) = (bpar%nbe==bpar%nb)
195  bpar%cv_block(ib) = (iv==jv).and.(its==jts)
196  case ('specific_univariate')
197  bpar%diag_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
198  bpar%avg_block(ib) = .false.
199  bpar%B_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
200  bpar%nicas_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
201  bpar%cv_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
202  case ('specific_multivariate')
203  bpar%diag_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
204  bpar%avg_block(ib) = (bpar%nbe==bpar%nb)
205  bpar%B_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
206  bpar%nicas_block(ib) = (iv==jv).and.(its==1).and.(jts==1)
207  bpar%cv_block(ib) = (ib==bpar%nbe)
208  end select
209  bpar%fit_block(ib) = (bpar%nb==1).or.(bpar%diag_block(ib).and.(iv==jv).and.(its==jts) &
210  & .and.(trim(nam%minim_algo)/='none'))
211  if (nam%local_diag) bpar%fit_block(ib) = bpar%fit_block(ib).and.bpar%nicas_block(ib)
212 
213  ! Blocks information
214  write(bpar%blockname(ib),'(i2.2,a,i2.2,a,i2.2,a,i2.2)') iv,'_',jv,'_',its,'_',jts
215  bpar%b_to_v1(ib) = iv
216  bpar%b_to_v2(ib) = jv
217  bpar%b_to_ts1(ib) = its
218  bpar%b_to_ts2(ib) = jts
219 
220  ! Update block index
221  ib = ib+1
222  end do
223  end do
224  end do
225  end do
226 
227  if (bpar%nbe>bpar%nb) then
228  ! Common block
229  ib = bpar%nbe
230 
231  ! Classes and levels
232  bpar%nl0r(ib) = bpar%nl0rmax
233  do il0=1,geom%nl0
234  jl0off = il0-(bpar%nl0r(ib)-1)/2-1
235  if (jl0off<1) jl0off = 0
236  if (jl0off+bpar%nl0rmax>geom%nl0) jl0off = geom%nl0-bpar%nl0rmax
237  do jl0r=1,bpar%nl0rmax
238  bpar%l0rl0b_to_l0(jl0r,il0,ib) = jl0off+jl0r
239  if (bpar%l0rl0b_to_l0(jl0r,il0,ib)==il0) bpar%il0rz(il0,ib) = jl0r
240  end do
241  end do
242  bpar%nc3(ib) = nam%nc3
243 
244  ! Select blocks
245  select case (nam%strategy)
246  case ('diag_all')
247  bpar%diag_block(ib) = .true.
248  bpar%avg_block(ib) = .false.
249  bpar%B_block(ib) = .false.
250  bpar%nicas_block(ib) = .false.
251  bpar%cv_block(ib) = .false.
252  case ('common')
253  bpar%diag_block(ib) = .true.
254  bpar%avg_block(ib) = .false.
255  bpar%B_block(ib) = .true.
256  bpar%nicas_block(ib) = .true.
257  bpar%cv_block(ib) = .true.
258  case ('common_univariate')
259  bpar%diag_block(ib) = .true.
260  bpar%avg_block(ib) = .false.
261  bpar%B_block(ib) = .true.
262  bpar%nicas_block(ib) = .true.
263  bpar%cv_block(ib) = .false.
264  case ('common_weighted')
265  bpar%diag_block(ib) = .true.
266  bpar%avg_block(ib) = .false.
267  bpar%B_block(ib) = .true.
268  bpar%nicas_block(ib) = .true.
269  bpar%cv_block(ib) = .false.
270  case ('specific_univariate')
271  bpar%diag_block(ib) = .false.
272  bpar%avg_block(ib) = .false.
273  bpar%B_block(ib) = .false.
274  bpar%nicas_block(ib) = .false.
275  bpar%cv_block(ib) = .false.
276  case ('specific_multivariate')
277  bpar%diag_block(ib) = .false.
278  bpar%avg_block(ib) = .false.
279  bpar%B_block(ib) = .false.
280  bpar%nicas_block(ib) = .false.
281  bpar%cv_block(ib) = .true.
282  end select
283  bpar%fit_block(ib) = bpar%diag_block(ib).and.(trim(nam%minim_algo)/='none')
284  if (nam%local_diag) bpar%fit_block(ib) = bpar%fit_block(ib).and.bpar%nicas_block(ib)
285 
286  ! Blocks information
287  bpar%blockname(ib) = 'common'
288  bpar%b_to_v1(ib) = 0
289  bpar%b_to_v2(ib) = 0
290  bpar%b_to_ts1(ib) = 0
291  bpar%b_to_ts2(ib) = 0
292  end if
293 end if
294 
295 end subroutine bpar_alloc
296 
297 !----------------------------------------------------------------------
298 ! Subroutine: bpar_dealloc
299 ! Purpose: deallocate general parameters
300 !----------------------------------------------------------------------
301 subroutine bpar_dealloc(bpar)
303 implicit none
304 
305 ! Passed variable
306 class(bpar_type),intent(inout) :: bpar ! Block parameters
307 
308 ! Release memory
309 if (allocated(bpar%nl0r)) deallocate(bpar%nl0r)
310 if (allocated(bpar%l0rl0b_to_l0)) deallocate(bpar%l0rl0b_to_l0)
311 if (allocated(bpar%il0rz)) deallocate(bpar%il0rz)
312 if (allocated(bpar%nc3)) deallocate(bpar%nc3)
313 if (allocated(bpar%vbal_block)) deallocate(bpar%vbal_block)
314 if (allocated(bpar%diag_block)) deallocate(bpar%diag_block)
315 if (allocated(bpar%avg_block)) deallocate(bpar%avg_block)
316 if (allocated(bpar%fit_block)) deallocate(bpar%fit_block)
317 if (allocated(bpar%B_block)) deallocate(bpar%B_block)
318 if (allocated(bpar%nicas_block)) deallocate(bpar%nicas_block)
319 if (allocated(bpar%cv_block)) deallocate(bpar%cv_block)
320 if (allocated(bpar%blockname)) deallocate(bpar%blockname)
321 if (allocated(bpar%b_to_v1)) deallocate(bpar%b_to_v1)
322 if (allocated(bpar%b_to_v2)) deallocate(bpar%b_to_v2)
323 if (allocated(bpar%b_to_ts1)) deallocate(bpar%b_to_ts1)
324 if (allocated(bpar%b_to_ts2)) deallocate(bpar%b_to_ts2)
325 
326 end subroutine bpar_dealloc
327 
328 end module type_bpar
subroutine bpar_dealloc(bpar)
Definition: type_bpar.F90:302
subroutine bpar_alloc(bpar, nam, geom)
Definition: type_bpar.F90:53
integer, parameter, public kind_real
#define min(a, b)
Definition: mosaic_util.h:32