FV3 Bundle
type_vbal_blk.F90
Go to the documentation of this file.
1 !----------------------------------------------------------------------
2 ! Module: type_vbal_blk
3 ! Purpose: vertical balance block 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 type_bpar, only: bpar_type
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 ! Vertical balance block derived type
20  integer :: iv ! First variable index
21  integer :: jv ! Second variable index
22  character(len=1024) :: name ! Name
23  real(kind_real),allocatable :: auto(:,:,:) ! Auto-covariance
24  real(kind_real),allocatable :: cross(:,:,:) ! Cross-covariance
25  real(kind_real),allocatable :: auto_inv(:,:,:) ! Inverse auto-covariance
26  real(kind_real),allocatable :: reg(:,:,:) ! Regression
27 contains
28  procedure :: alloc => vbal_blk_alloc
29  procedure :: dealloc => vbal_blk_dealloc
30  procedure :: apply => vbal_blk_apply
31  procedure :: apply_ad => vbal_blk_apply_ad
32 end type vbal_blk_type
33 
34 private
35 public :: vbal_blk_type
36 
37 contains
38 
39 !----------------------------------------------------------------------
40 ! Subroutine: vbal_blk_alloc
41 ! Purpose: vertical balance block data allocation
42 !----------------------------------------------------------------------
43 subroutine vbal_blk_alloc(vbal_blk,nam,geom,nc2b,iv,jv)
44 
45 implicit none
46 
47 ! Passed variables
48 class(vbal_blk_type),intent(inout) :: vbal_blk ! Vertical balance block
49 type(nam_type),intent(in) :: nam ! Namelist
50 type(geom_type),intent(in) :: geom ! Geometry
51 integer,intent(in) :: nc2b ! Subset Sc2 size, halo B
52 integer,intent(in) :: iv ! First variable index
53 integer,intent(in) :: jv ! Second variable index
54 
55 ! Set attributes
56 vbal_blk%iv = iv
57 vbal_blk%jv = jv
58 vbal_blk%name = trim(nam%varname(jv))//'-'//trim(nam%varname(iv))
59 
60 ! Allocation
61 allocate(vbal_blk%auto(nc2b,geom%nl0,geom%nl0))
62 allocate(vbal_blk%cross(nc2b,geom%nl0,geom%nl0))
63 allocate(vbal_blk%auto_inv(nc2b,geom%nl0,geom%nl0))
64 allocate(vbal_blk%reg(nc2b,geom%nl0,geom%nl0))
65 
66 end subroutine vbal_blk_alloc
67 
68 !----------------------------------------------------------------------
69 ! Subroutine: vbal_blk_dealloc
70 ! Purpose: vertical balance block data deallocation
71 !----------------------------------------------------------------------
72 subroutine vbal_blk_dealloc(vbal_blk)
73 
74 implicit none
75 
76 ! Passed variables
77 class(vbal_blk_type),intent(inout) :: vbal_blk ! Vertical balance block
78 
79 ! Allocation
80 if (allocated(vbal_blk%auto)) deallocate(vbal_blk%auto)
81 if (allocated(vbal_blk%cross)) deallocate(vbal_blk%cross)
82 if (allocated(vbal_blk%auto_inv)) deallocate(vbal_blk%auto_inv)
83 if (allocated(vbal_blk%reg)) deallocate(vbal_blk%reg)
84 
85 end subroutine vbal_blk_dealloc
86 
87 !----------------------------------------------------------------------
88 ! Subroutine: vbal_blk_apply
89 ! Purpose: apply vertical balance block
90 !----------------------------------------------------------------------
91 subroutine vbal_blk_apply(vbal_blk,geom,np,h_n_s,h_c2b,h_S,fld)
92 
93 implicit none
94 
95 ! Passed variables
96 class(vbal_blk_type),intent(in) :: vbal_blk ! Vertical balance block
97 type(geom_type),intent(in) :: geom ! Geometry
98 integer,intent(in) :: np ! Maximum number of neighbors
99 integer,intent(in) :: h_n_s(geom%nc0a,geom%nl0i) ! Number of neighbors for the horizontal interpolation
100 integer,intent(in) :: h_c2b(np,geom%nc0a,geom%nl0i) ! Index of neighbors for the horizontal interpolation
101 real(kind_real),intent(in) :: h_S(np,geom%nc0a,geom%nl0i) ! Weight of neighbors for the horizontal interpolation
102 real(kind_real),intent(inout) :: fld(geom%nc0a,geom%nl0) ! Source/destination vector
103 
104 ! Local variables
105 integer :: ic0a,il0,jl0,i_s,ic2b
106 real(kind_real) :: S,fld_tmp(geom%nc0a,geom%nl0)
107 
108 ! Initialization
109 fld_tmp = 0.0
110 
111 do ic0a=1,geom%nc0a
112  do il0=1,geom%nl0
113  do jl0=1,geom%nl0
114  ! Loop over neighbors
115  do i_s=1,h_n_s(ic0a,min(il0,geom%nl0i))
116  ! Find neighbor index and weight
117  ic2b = h_c2b(i_s,ic0a,min(il0,geom%nl0i))
118  s = h_s(i_s,ic0a,min(il0,geom%nl0i))
119 
120  ! Apply regression coefficient weighted by the neighbor weight
121  fld_tmp(ic0a,il0) = fld_tmp(ic0a,il0)+s*vbal_blk%reg(ic2b,il0,jl0)*fld(ic0a,jl0)
122  end do
123  end do
124  end do
125 end do
126 
127 ! Final copy
128 fld = fld_tmp
129 
130 end subroutine vbal_blk_apply
131 
132 !----------------------------------------------------------------------
133 ! Subroutine: vbal_blk_apply_ad
134 ! Purpose: apply adjoint vertical balance block
135 !----------------------------------------------------------------------
136 subroutine vbal_blk_apply_ad(vbal_blk,geom,np,h_n_s,h_c2b,h_S,fld)
138 implicit none
139 
140 ! Passed variables
141 class(vbal_blk_type),intent(in) :: vbal_blk ! Vertical balance block
142 type(geom_type),intent(in) :: geom ! Geometry
143 integer,intent(in) :: np ! Maximum number of neighbors
144 integer,intent(in) :: h_n_s(geom%nc0a,geom%nl0i) ! Number of neighbors for the horizontal interpolation
145 integer,intent(in) :: h_c2b(np,geom%nc0a,geom%nl0i) ! Index of neighbors for the horizontal interpolation
146 real(kind_real),intent(in) :: h_S(np,geom%nc0a,geom%nl0i) ! Weight of neighbors for the horizontal interpolation
147 real(kind_real),intent(inout) :: fld(geom%nc0a,geom%nl0) ! Source/destination vector
148 
149 ! Local variables
150 integer :: ic0a,il0,jl0,i_s,ic2b
151 real(kind_real) :: S,fld_tmp(geom%nc0a,geom%nl0)
152 
153 ! Initialization
154 fld_tmp = 0.0
155 
156 do ic0a=1,geom%nc0a
157  do il0=1,geom%nl0
158  do jl0=1,geom%nl0
159  ! Loop over neighbors
160  do i_s=1,h_n_s(ic0a,min(il0,geom%nl0i))
161  ! Find neighbor index and weight
162  ic2b = h_c2b(i_s,ic0a,min(il0,geom%nl0i))
163  s = h_s(i_s,ic0a,min(il0,geom%nl0i))
164 
165  ! Apply regression coefficient weighted by the neighbor weight
166  fld_tmp(ic0a,jl0) = fld_tmp(ic0a,jl0)+s*vbal_blk%reg(ic2b,il0,jl0)*fld(ic0a,il0)
167  end do
168  end do
169  end do
170 end do
171 
172 ! Final copy
173 fld = fld_tmp
174 
175 end subroutine vbal_blk_apply_ad
176 
177 end module type_vbal_blk
subroutine vbal_blk_alloc(vbal_blk, nam, geom, nc2b, iv, jv)
subroutine vbal_blk_apply_ad(vbal_blk, geom, np, h_n_s, h_c2b, h_S, fld)
subroutine vbal_blk_dealloc(vbal_blk)
double * cross(const double *p1, const double *p2)
integer, parameter, public kind_real
#define min(a, b)
Definition: mosaic_util.h:32
subroutine vbal_blk_apply(vbal_blk, geom, np, h_n_s, h_c2b, h_S, fld)