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