17 use,
intrinsic :: iso_c_binding
21 use mpi_f08,
only : mpi_allreduce, mpi_integer, mpi_sum
27 integer :: nelv_2d = 0
28 integer :: glb_nelv_2d = 0
29 integer :: offset_el_2d = 0
32 integer,
allocatable :: idx_2d(:)
33 integer,
allocatable :: el_idx_2d(:)
37 type(
coef_t),
pointer :: coef => null()
43 real(kind=
rp) :: domain_height
47 generic :: init => init_int, init_char
50 generic :: average => average_list, average_file
56 class(
map_2d_t),
intent(inout) :: this
57 type(
coef_t),
intent(inout),
target :: coef
58 integer,
intent(in) :: dir
59 real(kind=
rp),
intent(in) :: tol
60 real(kind=
rp) :: el_dim(3,3), el_h
61 integer :: i, e, j, ierr, k, lx, lxy, n
62 call this%map_1d%init(coef,dir,tol)
67 call this%u%init(this%dof)
68 call this%old_u%init(this%dof)
69 call this%avg_u%init(this%dof)
70 call this%el_heights%init(this%dof)
77 do i = 1, this%msh%nelv
78 if (this%map_1d%el_lvl(i) .eq. 1)
then
79 this%nelv_2d = this%nelv_2d + 1
83 call mpi_allreduce(this%nelv_2d, this%glb_nelv_2d, 1, &
87 call mpi_exscan(this%nelv_2d, this%offset_el_2d, 1, &
89 allocate(this%el_idx_2d(this%nelv_2d))
90 do i = 1, this%nelv_2d
91 this%el_idx_2d(i) = this%offset_el_2d + i
93 this%n_2d = this%nelv_2d*lxy
94 allocate(this%idx_2d(this%n_2d))
97 do e = 1, this%msh%nelv
98 if (this%map_1d%el_lvl(e) .eq. 1)
then
99 if (this%map_1d%dir_el(e) .eq. 1)
then
104 if (this%map_1d%dir_el(e) .eq. 2)
then
107 this%idx_2d(j*lxy+k+lx*(i-1)) =
linear_index(k,1,i,e,lx,lx,lx)
111 if (this%map_1d%dir_el(e) .eq. 3)
then
119 do i = 1, this%msh%nelv
124 el_dim(1,:) = abs(this%msh%elements(i)%e%pts(1)%p%x-&
125 this%msh%elements(i)%e%pts(2)%p%x)
126 el_dim(2,:) = abs(this%msh%elements(i)%e%pts(1)%p%x-&
127 this%msh%elements(i)%e%pts(3)%p%x)
128 el_dim(3,:) = abs(this%msh%elements(i)%e%pts(1)%p%x-&
129 this%msh%elements(i)%e%pts(5)%p%x)
131 el_h = el_dim(this%map_1d%dir_el(i),dir)
132 this%el_heights%x(:,:,:,i) = el_h
137 call copy(this%u%x,this%el_heights%x,n)
138 call copy(this%old_u%x,this%el_heights%x,n)
139 call copy(this%avg_u%x,this%el_heights%x,n)
141 this%map_1d%n_el_lvls, &
142 this%map_1d%dir_el,this%coef%gs_h, this%coef%mult, this%msh%nelv, lx)
143 this%domain_height = this%u%x(1,1,1,1)
149 type(
coef_t),
intent(inout),
target :: coef
150 character(len=*),
intent(in) :: dir
151 real(kind=
rp),
intent(in) :: tol
154 if (trim(dir) .eq.
'x')
then
156 else if (trim(dir) .eq.
'y')
then
158 else if (trim(dir) .eq.
'z')
then
161 call neko_error(
'Direction not supported, map_2d')
164 call this%init(coef,idir,tol)
174 class(
map_2d_t),
intent(inout) :: this
177 real(kind=
rp),
pointer,
dimension(:,:,:,:) :: x_ptr, y_ptr
181 integer :: i, j, lx, lxy
183 call fld_data2d%init(this%nelv_2d,this%offset_el_2d)
185 fld_data2d%lx = this%dof%Xh%lx
186 fld_data2d%ly = this%dof%Xh%ly
189 fld_data2d%glb_nelv = this%glb_nelv_2d
190 lxy = fld_data2d%lx*fld_data2d%ly
191 n_2d = lxy*this%nelv_2d
193 call fld_data2d%x%init(n_2d)
194 call fld_data2d%y%init(n_2d)
195 allocate(fld_data2d%idx(this%nelv_2d))
197 if (this%dir .eq. 1)
then
201 if (this%dir .eq. 2)
then
205 if (this%dir .eq. 3)
then
209 do j = 1, this%nelv_2d
210 fld_data2d%idx(j) = this%el_idx_2d(j)
213 fld_data2d%x%x(j) = x_ptr(this%idx_2d(j),1,1,1)
214 fld_data2d%y%x(j) = y_ptr(this%idx_2d(j),1,1,1)
216 allocate(fields2d(fld_data3d%size()))
218 call fld_data2d%init_n_fields(fld_data3d%size(),n_2d)
222 do i = 1, fld_data3d%size()
223 call copy(this%old_u%x,fld_data3d%items(i)%ptr%x,n)
225 this%el_heights, this%domain_height, &
226 this%map_1d%dir_el, this%coef, this%msh%nelv, lx)
227 call copy(this%old_u%x,this%u%x,n)
228 call copy(this%avg_u%x,this%u%x,n)
230 this%old_u, this%map_1d%n_el_lvls, &
231 this%map_1d%dir_el,this%coef%gs_h, this%coef%mult, &
233 call copy(fld_data3d%items(i)%ptr%x,this%u%x,n)
235 call fld_data2d%get_list(fields2d,fld_data2d%size())
236 do i = 1, fld_data3d%size()
238 fields2d(i)%ptr%x(j) = fld_data3d%items(i)%ptr%x(this%idx_2d(j),1,1,1)
250 class(
map_2d_t),
intent(inout) :: this
253 real(kind=
rp),
pointer,
dimension(:,:,:,:) :: x_ptr, y_ptr
255 type(
vector_ptr_t),
allocatable :: fields3d(:), fields2d(:)
257 integer :: i, j, lx, lxy
259 call fld_data2d%init(this%nelv_2d,this%offset_el_2d)
261 fld_data2d%lx = this%dof%Xh%lx
262 fld_data2d%ly = this%dof%Xh%ly
265 fld_data2d%glb_nelv = this%glb_nelv_2d
266 lxy = fld_data2d%lx*fld_data2d%ly
267 n_2d = lxy*this%nelv_2d
269 call fld_data2d%x%init(n_2d)
270 call fld_data2d%y%init(n_2d)
271 allocate(fld_data2d%idx(n_2d))
273 if (this%dir .eq. 1)
then
277 if (this%dir .eq. 2)
then
281 if (this%dir .eq. 3)
then
286 fld_data2d%idx(j) = this%idx_2d(j)
287 fld_data2d%x%x(j) = x_ptr(this%idx_2d(j),1,1,1)
288 fld_data2d%y%x(j) = y_ptr(this%idx_2d(j),1,1,1)
290 allocate(fields3d(fld_data3d%size()))
291 allocate(fields2d(fld_data3d%size()))
293 call fld_data2d%init_n_fields(fld_data3d%size(),n_2d)
294 call fld_data3d%get_list(fields3d,fld_data3d%size())
299 do i = 1, fld_data3d%size()
300 call copy(this%old_u%x,fields3d(i)%ptr%x,n)
302 this%el_heights, this%domain_height, &
303 this%map_1d%dir_el, this%coef, this%msh%nelv, lx)
304 call copy(this%old_u%x,this%u%x,n)
305 call copy(this%avg_u%x,this%u%x,n)
307 this%old_u, this%map_1d%n_el_lvls, &
308 this%map_1d%dir_el,this%coef%gs_h,&
309 this%coef%mult, this%msh%nelv, lx)
310 call copy(fields3d(i)%ptr%x,this%u%x,n)
312 call fld_data2d%get_list(fields2d,fld_data2d%size())
313 do i = 1, fld_data3d%size()
315 fields2d(i)%ptr%x(j) = fields3d(i)%ptr%x(this%idx_2d(j))
321 hom_dir_el, gs_h, mult, nelv, lx)
322 type(
field_t),
intent(inout) :: u, avg_u, old_u
323 type(
gs_t),
intent(inout) :: gs_h
324 integer,
intent(in) :: n_levels, nelv, lx
325 integer,
intent(in) :: hom_dir_el(nelv)
326 real(kind=
rp),
intent(in) :: mult(nelv*lx**3)
327 real(kind=
rp) :: temp_el(lx,lx,lx)
328 integer :: n, i, j, e
337 call gs_h%op(u,gs_op_add)
340 call col2(u%x,mult,n)
343 temp_el = 2.0*u%x(:,:,:,e)-old_u%x(:,:,:,e)
344 if (hom_dir_el(e) .eq. 1)
then
345 u%x(1,:,:,e) = temp_el(lx,:,:)
346 avg_u%x(1,:,:,e) = avg_u%x(1,:,:,e)+temp_el(1,:,:)
347 u%x(lx,:,:,e) = temp_el(1,:,:)
348 else if (hom_dir_el(e) .eq. 2)
then
349 u%x(:,1,:,e) = temp_el(:,lx,:)
350 avg_u%x(:,1,:,e) = avg_u%x(:,1,:,e)+temp_el(:,1,:)
351 u%x(:,lx,:,e) = temp_el(:,1,:)
352 else if (hom_dir_el(e) .eq. 3)
then
353 u%x(:,:,1,e) = temp_el(:,:,lx)
354 avg_u%x(:,:,1,e) = avg_u%x(:,:,1,e)+temp_el(:,:,1)
355 u%x(:,:,lx,e) = temp_el(:,:,1)
357 old_u%x(:,:,:,e) = u%x(:,:,:,e)
363 if (hom_dir_el(e) .eq. 1)
then
364 u%x(:,i,j,e) = avg_u%x(1,i,j,e)
365 else if (hom_dir_el(e) .eq. 2)
then
366 u%x(i,:,j,e) = avg_u%x(i,1,j,e)
367 else if (hom_dir_el(e) .eq. 3)
then
368 u%x(i,j,:,e) = avg_u%x(i,j,1,e)
376 hom_dir_el, coef, nelv, lx)
377 type(
field_t),
intent(inout) :: u, u_out, el_heights
378 type(
coef_t),
intent(inout) :: coef
379 integer,
intent(in) :: nelv, lx
380 integer,
intent(in) :: hom_dir_el(nelv)
381 real(kind=
rp),
intent(in) :: domain_height
383 integer :: n, i, j, e, k
387 call col2(u%x,el_heights%x,n)
388 call cmult(u%x, 1.0_rp/(2.0*domain_height),n)
389 call rzero(u_out%x,n)
397 if (hom_dir_el(e) .eq. 1)
then
398 u_out%x(1,i,j,e) = u_out%x(1,i,j,e)+wt*u%x(k,i,j,e)
399 else if (hom_dir_el(e) .eq. 2)
then
400 u_out%x(i,1,j,e) = u_out%x(i,1,j,e)+wt*u%x(i,k,j,e)
401 else if (hom_dir_el(e) .eq. 3)
then
402 u_out%x(i,j,1,e) = u_out%x(i,j,1,e)+wt*u%x(i,j,k,e)
412 if (hom_dir_el(e) .eq. 1)
then
413 u_out%x(:,i,j,e) = u_out%x(1,i,j,e)
414 else if (hom_dir_el(e) .eq. 2)
then
415 u_out%x(i,:,j,e) = u_out%x(i,1,j,e)
416 else if (hom_dir_el(e) .eq. 3)
then
417 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...