16 use,
intrinsic :: iso_c_binding
20 use mpi_f08,
only : mpi_allreduce, mpi_integer, mpi_sum, mpi_exscan
26 integer :: nelv_2d = 0
27 integer :: glb_nelv_2d = 0
28 integer :: offset_el_2d = 0
31 integer,
allocatable :: idx_2d(:)
32 integer,
allocatable :: el_idx_2d(:)
36 type(
coef_t),
pointer :: coef => null()
42 real(kind=
rp) :: domain_height
46 generic :: init => init_int, init_char
49 generic :: average => average_list, average_file
55 class(
map_2d_t),
intent(inout) :: this
56 type(
coef_t),
intent(inout),
target :: coef
57 integer,
intent(in) :: dir
58 real(kind=
rp),
intent(in) :: tol
59 real(kind=
rp) :: el_dim(3,3), el_h
60 integer :: i, e, j, ierr, k, lx, lxy, n
61 call this%map_1d%init(coef,dir,tol)
66 call this%u%init(this%dof)
67 call this%old_u%init(this%dof)
68 call this%avg_u%init(this%dof)
69 call this%el_heights%init(this%dof)
76 do i = 1, this%msh%nelv
77 if (this%map_1d%el_lvl(i) .eq. 1)
then
78 this%nelv_2d = this%nelv_2d + 1
82 call mpi_allreduce(this%nelv_2d, this%glb_nelv_2d, 1, &
86 call mpi_exscan(this%nelv_2d, this%offset_el_2d, 1, &
88 allocate(this%el_idx_2d(this%nelv_2d))
89 do i = 1, this%nelv_2d
90 this%el_idx_2d(i) = this%offset_el_2d + i
92 this%n_2d = this%nelv_2d*lxy
93 allocate(this%idx_2d(this%n_2d))
96 do e = 1, this%msh%nelv
97 if (this%map_1d%el_lvl(e) .eq. 1)
then
98 if (this%map_1d%dir_el(e) .eq. 1)
then
103 if (this%map_1d%dir_el(e) .eq. 2)
then
106 this%idx_2d(j*lxy+k+lx*(i-1)) =
linear_index(k,1,i,e,lx,lx,lx)
110 if (this%map_1d%dir_el(e) .eq. 3)
then
118 do i = 1, this%msh%nelv
123 el_dim(1,:) = abs(this%msh%elements(i)%e%pts(1)%p%x-&
124 this%msh%elements(i)%e%pts(2)%p%x)
125 el_dim(2,:) = abs(this%msh%elements(i)%e%pts(1)%p%x-&
126 this%msh%elements(i)%e%pts(3)%p%x)
127 el_dim(3,:) = abs(this%msh%elements(i)%e%pts(1)%p%x-&
128 this%msh%elements(i)%e%pts(5)%p%x)
130 el_h = el_dim(this%map_1d%dir_el(i),dir)
131 this%el_heights%x(:,:,:,i) = el_h
136 call copy(this%u%x,this%el_heights%x,n)
137 call copy(this%old_u%x,this%el_heights%x,n)
138 call copy(this%avg_u%x,this%el_heights%x,n)
140 this%map_1d%n_el_lvls, &
141 this%map_1d%dir_el,this%coef%gs_h, this%coef%mult, this%msh%nelv, lx)
142 this%domain_height = this%u%x(1,1,1,1)
148 type(
coef_t),
intent(inout),
target :: coef
149 character(len=*),
intent(in) :: dir
150 real(kind=
rp),
intent(in) :: tol
153 if (trim(dir) .eq.
'x')
then
155 else if (trim(dir) .eq.
'y')
then
157 else if (trim(dir) .eq.
'z')
then
160 call neko_error(
'Direction not supported, map_2d')
163 call this%init(coef,idir,tol)
173 class(
map_2d_t),
intent(inout) :: this
176 real(kind=
rp),
pointer,
dimension(:,:,:,:) :: x_ptr, y_ptr
180 integer :: i, j, lx, lxy
182 call fld_data2d%init(this%nelv_2d, this%offset_el_2d)
184 fld_data2d%lx = this%dof%Xh%lx
185 fld_data2d%ly = this%dof%Xh%ly
188 fld_data2d%glb_nelv = this%glb_nelv_2d
189 lxy = fld_data2d%lx*fld_data2d%ly
190 n_2d = lxy*this%nelv_2d
192 call fld_data2d%x%init(n_2d)
193 call fld_data2d%y%init(n_2d)
194 allocate(fld_data2d%idx(this%nelv_2d))
196 if (this%dir .eq. 1)
then
200 if (this%dir .eq. 2)
then
204 if (this%dir .eq. 3)
then
208 do j = 1, this%nelv_2d
209 fld_data2d%idx(j) = this%el_idx_2d(j)
212 fld_data2d%x%x(j) = x_ptr(this%idx_2d(j),1,1,1)
213 fld_data2d%y%x(j) = y_ptr(this%idx_2d(j),1,1,1)
215 allocate(fields2d(fld_data3d%size()))
217 call fld_data2d%init_n_fields(fld_data3d%size(), n_2d)
221 do i = 1, fld_data3d%size()
222 call copy(this%old_u%x, fld_data3d%items(i)%ptr%x, n)
224 this%el_heights, this%domain_height, &
225 this%map_1d%dir_el, this%coef, this%msh%nelv, lx)
226 call copy(this%old_u%x, this%u%x, n)
227 call copy(this%avg_u%x, this%u%x, n)
229 this%old_u, this%map_1d%n_el_lvls, &
230 this%map_1d%dir_el,this%coef%gs_h, this%coef%mult, &
232 call copy(fld_data3d%items(i)%ptr%x, this%u%x, n)
234 call fld_data2d%get_list(fields2d, fld_data2d%size())
235 do i = 1, fld_data3d%size()
237 fields2d(i)%ptr%x(j) = fld_data3d%items(i)%ptr%x(this%idx_2d(j),1,1,1)
249 class(
map_2d_t),
intent(inout) :: this
252 real(kind=
rp),
pointer,
dimension(:,:,:,:) :: x_ptr, y_ptr
254 type(
vector_ptr_t),
allocatable :: fields3d(:), fields2d(:)
256 integer :: i, j, lx, lxy
258 call fld_data2d%init(this%nelv_2d, this%offset_el_2d)
260 fld_data2d%lx = this%dof%Xh%lx
261 fld_data2d%ly = this%dof%Xh%ly
264 fld_data2d%glb_nelv = this%glb_nelv_2d
265 lxy = fld_data2d%lx*fld_data2d%ly
266 n_2d = lxy*this%nelv_2d
268 call fld_data2d%x%init(n_2d)
269 call fld_data2d%y%init(n_2d)
270 allocate(fld_data2d%idx(n_2d))
272 if (this%dir .eq. 1)
then
276 if (this%dir .eq. 2)
then
280 if (this%dir .eq. 3)
then
285 fld_data2d%idx(j) = this%idx_2d(j)
286 fld_data2d%x%x(j) = x_ptr(this%idx_2d(j),1,1,1)
287 fld_data2d%y%x(j) = y_ptr(this%idx_2d(j),1,1,1)
289 allocate(fields3d(fld_data3d%size()))
290 allocate(fields2d(fld_data3d%size()))
292 call fld_data2d%init_n_fields(fld_data3d%size(), n_2d)
293 call fld_data3d%get_list(fields3d,fld_data3d%size())
298 do i = 1, fld_data3d%size()
299 call copy(this%old_u%x, fields3d(i)%ptr%x, n)
301 this%el_heights, this%domain_height, &
302 this%map_1d%dir_el, this%coef, this%msh%nelv, lx)
303 call copy(this%old_u%x, this%u%x,n)
304 call copy(this%avg_u%x, this%u%x,n)
306 this%old_u, this%map_1d%n_el_lvls, &
307 this%map_1d%dir_el, this%coef%gs_h,&
308 this%coef%mult, this%msh%nelv, lx)
309 call copy(fields3d(i)%ptr%x, this%u%x, n)
311 call fld_data2d%get_list(fields2d, fld_data2d%size())
312 do i = 1, fld_data3d%size()
314 fields2d(i)%ptr%x(j) = fields3d(i)%ptr%x(this%idx_2d(j))
320 hom_dir_el, gs_h, mult, nelv, lx)
321 type(
field_t),
intent(inout) :: u, avg_u, old_u
322 type(
gs_t),
intent(inout) :: gs_h
323 integer,
intent(in) :: n_levels, nelv, lx
324 integer,
intent(in) :: hom_dir_el(nelv)
325 real(kind=
rp),
intent(in) :: mult(nelv*lx**3)
326 real(kind=
rp) :: temp_el(lx,lx,lx)
327 integer :: n, i, j, e
336 call gs_h%op(u,gs_op_add)
339 call col2(u%x,mult,n)
342 temp_el = 2.0*u%x(:,:,:,e)-old_u%x(:,:,:,e)
343 if (hom_dir_el(e) .eq. 1)
then
344 u%x(1,:,:,e) = temp_el(lx,:,:)
345 avg_u%x(1,:,:,e) = avg_u%x(1,:,:,e)+temp_el(1,:,:)
346 u%x(lx,:,:,e) = temp_el(1,:,:)
347 else if (hom_dir_el(e) .eq. 2)
then
348 u%x(:,1,:,e) = temp_el(:,lx,:)
349 avg_u%x(:,1,:,e) = avg_u%x(:,1,:,e)+temp_el(:,1,:)
350 u%x(:,lx,:,e) = temp_el(:,1,:)
351 else if (hom_dir_el(e) .eq. 3)
then
352 u%x(:,:,1,e) = temp_el(:,:,lx)
353 avg_u%x(:,:,1,e) = avg_u%x(:,:,1,e)+temp_el(:,:,1)
354 u%x(:,:,lx,e) = temp_el(:,:,1)
356 old_u%x(:,:,:,e) = u%x(:,:,:,e)
362 if (hom_dir_el(e) .eq. 1)
then
363 u%x(:,i,j,e) = avg_u%x(1,i,j,e)
364 else if (hom_dir_el(e) .eq. 2)
then
365 u%x(i,:,j,e) = avg_u%x(i,1,j,e)
366 else if (hom_dir_el(e) .eq. 3)
then
367 u%x(i,j,:,e) = avg_u%x(i,j,1,e)
375 hom_dir_el, coef, nelv, lx)
376 type(
field_t),
intent(inout) :: u, u_out, el_heights
377 type(
coef_t),
intent(inout) :: coef
378 integer,
intent(in) :: nelv, lx
379 integer,
intent(in) :: hom_dir_el(nelv)
380 real(kind=
rp),
intent(in) :: domain_height
382 integer :: n, i, j, e, k
386 call col2(u%x,el_heights%x,n)
387 call cmult(u%x, 1.0_rp/(2.0*domain_height),n)
388 call rzero(u_out%x,n)
396 if (hom_dir_el(e) .eq. 1)
then
397 u_out%x(1,i,j,e) = u_out%x(1,i,j,e)+wt*u%x(k,i,j,e)
398 else if (hom_dir_el(e) .eq. 2)
then
399 u_out%x(i,1,j,e) = u_out%x(i,1,j,e)+wt*u%x(i,k,j,e)
400 else if (hom_dir_el(e) .eq. 3)
then
401 u_out%x(i,j,1,e) = u_out%x(i,j,1,e)+wt*u%x(i,j,k,e)
411 if (hom_dir_el(e) .eq. 1)
then
412 u_out%x(:,i,j,e) = u_out%x(1,i,j,e)
413 else if (hom_dir_el(e) .eq. 2)
then
414 u_out%x(i,:,j,e) = u_out%x(i,1,j,e)
415 else if (hom_dir_el(e) .eq. 3)
then
416 u_out%x(i,j,:,e) = u_out%x(i,j,1,e)
Copy data between host and device (or device and device)
type(mpi_comm), public neko_comm
MPI communicator.
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
integer, parameter, public device_to_host
Defines a mapping of the degrees of freedom.
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
Creates a 1d GLL point map along a specified direction based on the connectivity in the mesh.
Maps a 3D dofmap to a 2D spectral element grid.
subroutine perform_local_summation(u_out, u, el_heights, domain_height, hom_dir_el, coef, nelv, lx)
subroutine perform_global_summation(u, avg_u, old_u, n_levels, hom_dir_el, gs_h, mult, nelv, lx)
subroutine map_2d_init(this, coef, dir, tol)
subroutine map_2d_init_char(this, coef, dir, tol)
subroutine map_2d_average_field_list(this, fld_data2d, fld_data3d)
Computes average if field list in one direction and outputs 2D field with averaged values.
subroutine map_2d_average(this, fld_data2d, fld_data3d)
Computes average if field list in one direction and outputs 2D field with averaged values.
subroutine, public cmult(a, c, n)
Multiplication by constant c .
subroutine, public col2(a, b, n)
Vector multiplication .
subroutine, public copy(a, b, n)
Copy a vector .
subroutine, public rzero(a, n)
Zero a real vector.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
pure integer function, public linear_index(i, j, k, l, lx, ly, lz)
Compute the address of a (i,j,k,l) array with sizes (1:lx, 1:ly, 1:lz, :)
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
field_list_t, To be able to group fields together
Type that encapsulates a mapping from each gll point in the mesh to its corresponding (global) GLL po...