Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
map_2d.f90
Go to the documentation of this file.
1
2
3module map_2d
5 use num_types, only: rp
6 use dofmap, only: dofmap_t
7 use map_1d, only : map_1d_t
8 use gather_scatter, only : gs_t, gs_op_add
9 use mesh, only: mesh_t
10 use field_list, only : field_list_t
11 use coefs, only: coef_t
12 use field_list, only: field_list_t
13 use vector, only: vector_ptr_t
14 use utils, only: neko_error, linear_index
15 use math, only: cmult, col2, copy, rzero
16 use field, only : field_t
17 use, intrinsic :: iso_c_binding
19 use comm, only : neko_comm
21 use mpi_f08, only : mpi_allreduce, mpi_integer, mpi_sum
23 implicit none
24 private
25
26 type, public :: map_2d_t
27 integer :: nelv_2d = 0
28 integer :: glb_nelv_2d = 0
29 integer :: offset_el_2d = 0
30 integer :: lxy = 0
31 integer :: n_2d = 0
32 integer, allocatable :: idx_2d(:)
33 integer, allocatable :: el_idx_2d(:)
35 type(mesh_t), pointer :: msh
36 type(dofmap_t), pointer :: dof => null()
37 type(coef_t), pointer :: coef => null()
38 type(field_t) :: u
39 type(field_t) :: old_u
40 type(field_t) :: avg_u
41 type(field_t) :: el_heights
42 integer :: dir
43 real(kind=rp) :: domain_height
44 contains
45 procedure, pass(this) :: init_int => map_2d_init
46 procedure, pass(this) :: init_char => map_2d_init_char
47 generic :: init => init_int, init_char
48 procedure, pass(this) :: average_file => map_2d_average
49 procedure, pass(this) :: average_list => map_2d_average_field_list
50 generic :: average => average_list, average_file
51 end type map_2d_t
52
53contains
54
55 subroutine map_2d_init(this,coef, dir, tol)
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)
63 this%msh => coef%msh
64 this%coef => coef
65 this%dof => coef%dof
66
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)
71 this%dir = dir
72
73 n = this%dof%size()
74 lx = this%dof%Xh%lx
75 lxy = this%dof%Xh%lxy
76 this%nelv_2d = 0
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
80 end if
81 end do
82 this%glb_nelv_2d = 0
83 call mpi_allreduce(this%nelv_2d, this%glb_nelv_2d, 1, &
84 mpi_integer, mpi_sum, neko_comm, ierr)
85
86 this%offset_el_2d = 0
87 call mpi_exscan(this%nelv_2d, this%offset_el_2d, 1, &
88 mpi_integer, mpi_sum, neko_comm, ierr)
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
92 end do
93 this%n_2d = this%nelv_2d*lxy
94 allocate(this%idx_2d(this%n_2d))
95
96 j = 0
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
100 do i = 1, lxy
101 this%idx_2d(j*lxy+i) = linear_index(1,i,1,e,lx,lx,lx)
102 end do
103 end if
104 if (this%map_1d%dir_el(e) .eq. 2) then
105 do i = 1, lx
106 do k = 1, lx
107 this%idx_2d(j*lxy+k+lx*(i-1)) = linear_index(k,1,i,e,lx,lx,lx)
108 end do
109 end do
110 end if
111 if (this%map_1d%dir_el(e) .eq. 3) then
112 do i = 1, lxy
113 this%idx_2d(j*lxy+i) = linear_index(i,1,1,e,lx,lx,lx)
114 end do
115 end if
116 j = j +1
117 end if
118 end do
119 do i = 1, this%msh%nelv
120 !find height in hom-dir
121 !direction in local coords (r,s,t) that is hom is stored in map_1d%dir_el
122 !set element to height
123 !we assume elements are stacked on eachother...
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)
130 ! 1 corresponds to x, 2 to y, 3 to z
131 el_h = el_dim(this%map_1d%dir_el(i),dir)
132 this%el_heights%x(:,:,:,i) = el_h
133 end do
134 !Need to compute mapping for 3d to 2d
135 !Does order matter? Think its ok as long as values written in same order
136
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)
140 call perform_global_summation(this%u, this%avg_u, this%old_u, &
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)
144
145 end subroutine map_2d_init
146
147 subroutine map_2d_init_char(this, coef, dir, tol)
148 class(map_2d_t) :: this
149 type(coef_t), intent(inout), target :: coef
150 character(len=*), intent(in) :: dir
151 real(kind=rp), intent(in) :: tol
152 integer :: idir
153
154 if (trim(dir) .eq. 'x') then
155 idir = 1
156 else if (trim(dir) .eq. 'y') then
157 idir = 2
158 else if (trim(dir) .eq. 'z') then
159 idir = 3
160 else
161 call neko_error('Direction not supported, map_2d')
162 end if
163
164 call this%init(coef,idir,tol)
165
166 end subroutine map_2d_init_char
167
168
173 subroutine map_2d_average_field_list(this,fld_data2D,fld_data3D)
174 class(map_2d_t), intent(inout) :: this
175 type(fld_file_data_t), intent(inout) :: fld_data2D
176 type(field_list_t), intent(inout) :: fld_data3D
177 real(kind=rp), pointer, dimension(:,:,:,:) :: x_ptr, y_ptr
178
179 type(vector_ptr_t), allocatable :: fields2d(:)
180 integer :: n_2d, n
181 integer :: i, j, lx, lxy
182
183 call fld_data2d%init(this%nelv_2d,this%offset_el_2d)
184 fld_data2d%gdim = 2
185 fld_data2d%lx = this%dof%Xh%lx
186 fld_data2d%ly = this%dof%Xh%ly
187 fld_data2d%lz = 1
188 lx = this%dof%Xh%lx
189 fld_data2d%glb_nelv = this%glb_nelv_2d
190 lxy = fld_data2d%lx*fld_data2d%ly
191 n_2d = lxy*this%nelv_2d
192 n = this%dof%size()
193 call fld_data2d%x%init(n_2d)
194 call fld_data2d%y%init(n_2d)
195 allocate(fld_data2d%idx(this%nelv_2d))
196
197 if (this%dir .eq. 1) then
198 x_ptr => this%dof%z
199 y_ptr => this%dof%y
200 end if
201 if (this%dir .eq. 2) then
202 x_ptr => this%dof%x
203 y_ptr => this%dof%z
204 end if
205 if (this%dir .eq. 3) then
206 x_ptr => this%dof%x
207 y_ptr => this%dof%y
208 end if
209 do j = 1, this%nelv_2d
210 fld_data2d%idx(j) = this%el_idx_2d(j)
211 end do
212 do j = 1, n_2d
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)
215 end do
216 allocate(fields2d(fld_data3d%size()))
217
218 call fld_data2d%init_n_fields(fld_data3d%size(),n_2d)
219 this%u = 0.0_rp
220 this%old_u = 0.0_rp
221 this%avg_u = 0.0_rp
222 do i = 1, fld_data3d%size()
223 call copy(this%old_u%x,fld_data3d%items(i)%ptr%x,n)
224 call perform_local_summation(this%u,this%old_u, &
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)
229 call perform_global_summation(this%u, this%avg_u, &
230 this%old_u, this%map_1d%n_el_lvls, &
231 this%map_1d%dir_el,this%coef%gs_h, this%coef%mult, &
232 this%msh%nelv, lx)
233 call copy(fld_data3d%items(i)%ptr%x,this%u%x,n)
234 end do
235 call fld_data2d%get_list(fields2d,fld_data2d%size())
236 do i = 1, fld_data3d%size()
237 do j = 1, n_2d
238 fields2d(i)%ptr%x(j) = fld_data3d%items(i)%ptr%x(this%idx_2d(j),1,1,1)
239 end do
240 end do
241
242 end subroutine map_2d_average_field_list
243
244
249 subroutine map_2d_average(this,fld_data2D,fld_data3D)
250 class(map_2d_t), intent(inout) :: this
251 type(fld_file_data_t), intent(inout) :: fld_data2D
252 type(fld_file_data_t), intent(inout) :: fld_data3D
253 real(kind=rp), pointer, dimension(:,:,:,:) :: x_ptr, y_ptr
254
255 type(vector_ptr_t), allocatable :: fields3d(:), fields2d(:)
256 integer :: n_2d, n
257 integer :: i, j, lx, lxy
258
259 call fld_data2d%init(this%nelv_2d,this%offset_el_2d)
260 fld_data2d%gdim = 2
261 fld_data2d%lx = this%dof%Xh%lx
262 fld_data2d%ly = this%dof%Xh%ly
263 fld_data2d%lz = 1
264 lx = this%dof%Xh%lx
265 fld_data2d%glb_nelv = this%glb_nelv_2d
266 lxy = fld_data2d%lx*fld_data2d%ly
267 n_2d = lxy*this%nelv_2d
268 n = this%dof%size()
269 call fld_data2d%x%init(n_2d)
270 call fld_data2d%y%init(n_2d)
271 allocate(fld_data2d%idx(n_2d))
272
273 if (this%dir .eq. 1) then
274 x_ptr => this%dof%z
275 y_ptr => this%dof%y
276 end if
277 if (this%dir .eq. 2) then
278 x_ptr => this%dof%x
279 y_ptr => this%dof%z
280 end if
281 if (this%dir .eq. 3) then
282 x_ptr => this%dof%x
283 y_ptr => this%dof%y
284 end if
285 do j = 1, n_2d
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)
289 end do
290 allocate(fields3d(fld_data3d%size()))
291 allocate(fields2d(fld_data3d%size()))
292
293 call fld_data2d%init_n_fields(fld_data3d%size(),n_2d)
294 call fld_data3d%get_list(fields3d,fld_data3d%size())
295
296 this%u = 0.0_rp
297 this%old_u = 0.0_rp
298 this%avg_u = 0.0_rp
299 do i = 1, fld_data3d%size()
300 call copy(this%old_u%x,fields3d(i)%ptr%x,n)
301 call perform_local_summation(this%u,this%old_u,&
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)
306 call perform_global_summation(this%u, this%avg_u, &
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)
311 end do
312 call fld_data2d%get_list(fields2d,fld_data2d%size())
313 do i = 1, fld_data3d%size()
314 do j = 1, n_2d
315 fields2d(i)%ptr%x(j) = fields3d(i)%ptr%x(this%idx_2d(j))
316 end do
317 end do
318 end subroutine map_2d_average
319
320 subroutine perform_global_summation(u, avg_u, old_u, n_levels, &
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
329
330 n = u%dof%size()
331
332 do i = 1, n_levels-1
333 !compute average
334 if (neko_bcknd_device .eq. 1) &
335 call device_memcpy(u%x, u%x_d, n, &
336 host_to_device, sync=.true.)
337 call gs_h%op(u,gs_op_add)
338 if (neko_bcknd_device .eq. 1) &
339 call device_memcpy(u%x, u%x_d, n, device_to_host, sync=.true.)
340 call col2(u%x,mult,n)
341 !Assumes sugarcube I think
342 do e = 1, nelv
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)
356 end if
357 old_u%x(:,:,:,e) = u%x(:,:,:,e)
358 end do
359 end do
360 do e = 1, nelv
361 do i = 1,lx
362 do j = 1, lx
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)
369 end if
370 end do
371 end do
372 end do
373 end subroutine perform_global_summation
374
375 subroutine perform_local_summation(u_out, u, el_heights, domain_height,&
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
382 real(kind=rp) :: wt
383 integer :: n, i, j, e, k
384
385 n = nelv*lx**3
386
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)
390
391
392 do e = 1, nelv
393 do i = 1,lx
394 do j = 1, lx
395 do k = 1, lx
396 wt = coef%Xh%wx(k)
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)
403 end if
404 end do
405 end do
406 end do
407 end do
408
409 do e = 1, nelv
410 do i = 1,lx
411 do j = 1, lx
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)
418 end if
419 end do
420 end do
421 end do
422 end subroutine perform_local_summation
423
424end module map_2d
Copy data between host and device (or device and device)
Definition device.F90:66
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
integer, parameter, public device_to_host
Definition device.F90:47
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a field.
Definition field.f90:34
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
Gather-scatter.
Creates a 1d GLL point map along a specified direction based on the connectivity in the mesh.
Definition map_1d.f90:3
Maps a 3D dofmap to a 2D spectral element grid.
Definition map_2d.f90:3
subroutine perform_local_summation(u_out, u, el_heights, domain_height, hom_dir_el, coef, nelv, lx)
Definition map_2d.f90:377
subroutine perform_global_summation(u, avg_u, old_u, n_levels, hom_dir_el, gs_h, mult, nelv, lx)
Definition map_2d.f90:322
subroutine map_2d_init(this, coef, dir, tol)
Definition map_2d.f90:56
subroutine map_2d_init_char(this, coef, dir, tol)
Definition map_2d.f90:148
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.
Definition map_2d.f90:174
subroutine map_2d_average(this, fld_data2d, fld_data3d)
Computes average if field list in one direction and outputs 2D field with averaged values.
Definition map_2d.f90:250
Definition math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition math.f90:417
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:860
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:255
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:211
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Utilities.
Definition utils.f90:35
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, :)
Definition utils.f90:175
Defines a vector.
Definition vector.f90:34
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
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...
Definition map_1d.f90:27