Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
map_1d.f90
Go to the documentation of this file.
1
3module map_1d
5 use num_types, only : rp
6 use space, only : space_t
7 use dofmap, only : dofmap_t
8 use gather_scatter, only : gs_op_add
9 use mesh, only : mesh_t
13 use coefs, only : coef_t
14 use field_list, only : field_list_t
15 use matrix, only : matrix_t
16 use vector, only : vector_ptr_t
17 use utils, only : neko_error, neko_warning
18 use math, only : glmax, glmin, glimax, relcmp, cmult, add2s1, col2
19 use mpi_f08, only : mpi_allreduce, mpi_sum, mpi_barrier, mpi_in_place
20 use, intrinsic :: iso_c_binding
21 implicit none
22 private
27 type, public :: map_1d_t
30 integer, allocatable :: dir_el(:)
32 integer, allocatable :: el_lvl(:)
34 integer, allocatable :: pt_lvl(:, :, :, :)
36 integer :: n_el_lvls
38 integer :: n_gll_lvls
41 type(dofmap_t), pointer :: dof => null()
43 type(coef_t), pointer :: coef => null()
45 type(mesh_t), pointer :: msh => null()
47 integer :: dir
49 real(kind=rp) :: tol = 1e-7
51 real(kind=rp), allocatable :: volume_per_gll_lvl(:)
52 contains
54 procedure, pass(this) :: init_int => map_1d_init
55
56 procedure, pass(this) :: init_char => map_1d_init_char
57 generic :: init => init_int, init_char
59 procedure, pass(this) :: free => map_1d_free
61 procedure, pass(this) :: average_planes_fld_lst => &
63
64 procedure, pass(this) :: average_planes_vec_ptr => &
66 generic :: average_planes => average_planes_fld_lst, average_planes_vec_ptr
67 end type map_1d_t
68
69
70contains
71
72 subroutine map_1d_init(this, coef, dir, tol)
73 class(map_1d_t) :: this
74 type(coef_t), intent(inout), target :: coef
75 integer, intent(in) :: dir
76 real(kind=rp), intent(in) :: tol
77 integer :: nelv, lx, n, i, e, lvl, ierr
78 real(kind=rp), contiguous, pointer :: line(:, :, :, :)
79 real(kind=rp), allocatable :: min_vals(:, :, :, :)
80 real(kind=rp), allocatable :: min_temp(:, :, :, :)
81 type(c_ptr) :: min_vals_d = c_null_ptr
82 real(kind=rp) :: el_dim(3, 3), glb_min, glb_max, el_min
83
84 call this%free()
85
86 if (neko_bcknd_device .eq. 1) then
87 if (pe_rank .eq. 0) then
88 call neko_warning('map_1d does not copy indices to device, ' // &
89 ' but ok if used on cpu and for io')
90 end if
91 end if
92
93 this%dir = dir
94 this%dof => coef%dof
95 this%coef => coef
96 this%msh => coef%msh
97 nelv = this%msh%nelv
98 lx = this%dof%Xh%lx
99 n = this%dof%size()
100
101 if (dir .eq. 1) then
102 line => this%dof%x
103 else if (dir .eq. 2) then
104 line => this%dof%y
105 else if (dir .eq. 3) then
106 line => this%dof%z
107 else
108 call neko_error('Invalid dir for geopmetric comm')
109 end if
110 allocate(this%dir_el(nelv))
111 allocate(this%el_lvl(nelv))
112 allocate(this%pt_lvl(lx, lx, lx, nelv))
113 allocate(min_vals(lx, lx, lx, nelv))
114 allocate(min_temp(lx, lx, lx, nelv))
115 call mpi_barrier(neko_comm)
116 if (neko_bcknd_device .eq. 1) then
117
118 call device_map(min_vals, min_vals_d, n)
119 end if
120 call mpi_barrier(neko_comm)
121
122 do i = 1, nelv
123 !store which direction r,s,t corresponds to speciifed direction, x,y,z
124 !we assume elements are stacked on each other...
125 ! Check which one of the normalized vectors are closest to dir
126 ! If we want to incorporate other directions, we should look here
127 el_dim(1, :) = abs(this%msh%elements(i)%e%pts(1)%p%x - &
128 this%msh%elements(i)%e%pts(2)%p%x)
129 el_dim(1, :) = el_dim(1, :)/norm2(el_dim(1, :))
130 el_dim(2, :) = abs(this%msh%elements(i)%e%pts(1)%p%x - &
131 this%msh%elements(i)%e%pts(3)%p%x)
132 el_dim(2, :) = el_dim(2, :)/norm2(el_dim(2, :))
133 el_dim(3, :) = abs(this%msh%elements(i)%e%pts(1)%p%x - &
134 this%msh%elements(i)%e%pts(5)%p%x)
135 el_dim(3, :) = el_dim(3, :)/norm2(el_dim(3, :))
136 ! Checks which directions in rst the xyz corresponds to
137 ! 1 corresponds to r, 2 to s, 3 to t and are stored in dir_el
138 this%dir_el(i) = maxloc(el_dim(:, this%dir), dim = 1)
139 end do
140 glb_min = glmin(line, n)
141 glb_max = glmax(line, n)
142
143 i = 1
144 this%el_lvl = -1
145 ! Check what the mimum value in each element and put in min_vals
146 do e = 1, nelv
147 el_min = minval(line(:, :, :, e))
148 min_vals(:, :, :, e) = el_min
149 ! Check if this element is on the bottom,
150 ! in this case assign el_lvl = i = 1
151 if (relcmp(el_min, glb_min, this%tol)) then
152 if (this%el_lvl(e) .eq. -1) this%el_lvl(e) = i
153 end if
154 end do
155 ! While loop where at each iteation the global maximum value
156 ! propagates down one level.
157 ! When the minumum value has propagated to the highest level this stops.
158 ! Only works when the bottom plate of the domain is flat.
159 do while (.not. relcmp(glmax(min_vals, n), glb_min, this%tol))
160 i = i + 1
161 do e = 1, nelv
162 !Sets the value at the bottom of each element to glb_max
163 if (this%dir_el(e) .eq. 1) then
164 if (line(1, 1, 1, e) .gt. line(lx, 1, 1, e)) then
165 min_vals(lx, :, :, e) = glb_max
166 else
167 min_vals(1, :, :, e) = glb_max
168 end if
169 end if
170 if (this%dir_el(e) .eq. 2) then
171 if (line(1, 1, 1, e) .gt. line(1, lx, 1, e)) then
172 min_vals(:, lx, :, e) = glb_max
173 else
174 min_vals(:, 1, :, e) = glb_max
175 end if
176 end if
177 if (this%dir_el(e) .eq. 3) then
178 if (line(1, 1, 1, e) .gt. line(1, 1, lx, e)) then
179 min_vals(:, :, lx, e) = glb_max
180 else
181 min_vals(:, :, 1, e) = glb_max
182 end if
183 end if
184 end do
185 !Make sketchy min as GS_OP_MIN is not supported with device mpi
186 min_temp = min_vals
187 if (neko_bcknd_device .eq. 1) &
188 call device_memcpy(min_vals, min_vals_d, n, &
189 host_to_device, sync = .false.)
190 !Propagates the minumum value along the element boundary.
191 call coef%gs_h%op(min_vals, n, gs_op_add)
192 if (neko_bcknd_device .eq. 1) &
193 call device_memcpy(min_vals, min_vals_d, n, &
194 device_to_host, sync = .true.)
195 !Obtain average along boundary
196
197 call col2(min_vals, coef%mult, n)
198 call cmult(min_temp, -1.0_rp, n)
199 call add2s1(min_vals, min_temp, 2.0_rp, n)
200
201
202 !Checks the new minimum value on each element
203 !Assign this value to all points in this element in min_val
204 !If the element has not already been assinged a level,
205 !and it has obtained the minval, set el_lvl = i
206 do e = 1, nelv
207 el_min = minval(min_vals(:, :, :, e))
208 min_vals(:, :, :, e) = el_min
209 if (relcmp(el_min, glb_min, this%tol)) then
210 if (this%el_lvl(e) .eq. -1) this%el_lvl(e) = i
211 end if
212 end do
213 end do
214 this%n_el_lvls = glimax(this%el_lvl, nelv)
215 this%n_gll_lvls = this%n_el_lvls*lx
216 !Numbers the points in each element based on the element level
217 !and its orientation
218 do e = 1, nelv
219 do i = 1, lx
220 lvl = lx * (this%el_lvl(e) - 1) + i
221 if (this%dir_el(e) .eq. 1) then
222 if (line(1, 1, 1, e) .gt. line(lx, 1, 1, e)) then
223 this%pt_lvl(lx-i+1, :, :, e) = lvl
224 else
225 this%pt_lvl(i, :, :, e) = lvl
226 end if
227 end if
228 if (this%dir_el(e) .eq. 2) then
229 if (line(1, 1, 1, e) .gt. line(1, lx, 1, e)) then
230 this%pt_lvl(:, lx-i+1, :, e) = lvl
231 else
232 this%pt_lvl(:, i, :, e) = lvl
233 end if
234 end if
235 if (this%dir_el(e) .eq. 3) then
236 if (line(1, 1, 1, e) .gt. line(1, 1, lx, e)) then
237 this%pt_lvl(:, :, lx-i+1, e) = lvl
238 else
239 this%pt_lvl(:, :, i, e) = lvl
240 end if
241 end if
242 end do
243 end do
244 if (allocated(min_vals)) deallocate(min_vals)
245 if (c_associated(min_vals_d)) call device_free(min_vals_d)
246 if (allocated(min_temp)) deallocate(min_temp)
247 allocate(this%volume_per_gll_lvl(this%n_gll_lvls))
248 this%volume_per_gll_lvl = 0.0_rp
249 do i = 1, n
250 this%volume_per_gll_lvl(this%pt_lvl(i, 1, 1, 1)) = &
251 this%volume_per_gll_lvl(this%pt_lvl(i, 1, 1, 1)) + &
252 coef%B(i, 1, 1, 1)
253 end do
254 call mpi_allreduce(mpi_in_place, this%volume_per_gll_lvl, this%n_gll_lvls, &
255 mpi_real_precision, mpi_sum, neko_comm, ierr)
256
257 end subroutine map_1d_init
258
259 subroutine map_1d_init_char(this, coef, dir, tol)
260 class(map_1d_t) :: this
261 type(coef_t), intent(inout), target :: coef
262 character(len=*), intent(in) :: dir
263 real(kind=rp), intent(in) :: tol
264 integer :: idir
265
266 if (trim(dir) .eq. 'yz' .or. trim(dir) .eq. 'zy') then
267 idir = 1
268 else if (trim(dir) .eq. 'xz' .or. trim(dir) .eq. 'zx') then
269 idir = 2
270 else if (trim(dir) .eq. 'xy' .or. trim(dir) .eq. 'yx') then
271 idir = 3
272 else
273 call neko_error('homogenous direction not supported')
274 end if
275
276 call this%init(coef, idir, tol)
277
278 end subroutine map_1d_init_char
279
280 subroutine map_1d_free(this)
281 class(map_1d_t) :: this
282
283 if (allocated(this%dir_el)) deallocate(this%dir_el)
284 if (allocated(this%el_lvl)) deallocate(this%el_lvl)
285 if (allocated(this%pt_lvl)) deallocate(this%pt_lvl)
286 if (associated(this%dof)) nullify(this%dof)
287 if (associated(this%msh)) nullify(this%msh)
288 if (associated(this%coef)) nullify(this%coef)
289 if (allocated(this%volume_per_gll_lvl)) deallocate(this%volume_per_gll_lvl)
290 this%dir = 0
291 this%n_el_lvls = 0
292 this%n_gll_lvls = 0
293
294 end subroutine map_1d_free
295
296
302 subroutine map_1d_average_field_list(this, avg_planes, field_list)
303 class(map_1d_t), intent(inout) :: this
304 type(field_list_t), intent(inout) :: field_list
305 type(matrix_t), intent(inout) :: avg_planes
306 integer :: n, ierr, j, i
307 real(kind=rp) :: coord
308 call avg_planes%free()
309 call avg_planes%init(this%n_gll_lvls, field_list%size() + 1)
310 avg_planes = 0.0_rp
311 !ugly way of getting coordinates, computes average
312 n = this%dof%size()
313 do i = 1, n
314 if (this%dir .eq. 1) coord = this%dof%x(i,1,1,1)
315 if (this%dir .eq. 2) coord = this%dof%y(i,1,1,1)
316 if (this%dir .eq. 3) coord = this%dof%z(i,1,1,1)
317 avg_planes%x(this%pt_lvl(i,1,1,1), 1) = &
318 avg_planes%x(this%pt_lvl(i,1,1,1), 1) + &
319 coord * this%coef%B(i,1,1,1) / &
320 this%volume_per_gll_lvl(this%pt_lvl(i,1,1,1))
321 end do
322 do j = 2, field_list%size() + 1
323 do i = 1, n
324 avg_planes%x(this%pt_lvl(i,1,1,1), j) = &
325 avg_planes%x(this%pt_lvl(i,1,1,1), j) + &
326 field_list%items(j-1)%ptr%x(i,1,1,1) * this%coef%B(i,1,1,1) &
327 /this%volume_per_gll_lvl(this%pt_lvl(i,1,1,1))
328 end do
329 end do
330 if (pe_size .gt. 1) then
331 call mpi_allreduce(mpi_in_place, avg_planes%x, &
332 (field_list%size() + 1) * this%n_gll_lvls, mpi_real_precision, &
333 mpi_sum, neko_comm, ierr)
334 end if
335
336
337 end subroutine map_1d_average_field_list
338
344 subroutine map_1d_average_vector_ptr(this, avg_planes, vector_ptr)
345 class(map_1d_t), intent(inout) :: this
346 !Observe is an array...
347 type(vector_ptr_t), intent(inout) :: vector_ptr(:)
348 type(matrix_t), intent(inout) :: avg_planes
349 integer :: n, ierr, j, i
350 real(kind=rp) :: coord
351
352 call avg_planes%free()
353 call avg_planes%init(this%n_gll_lvls, size(vector_ptr) + 1)
354 !ugly way of getting coordinates, computes average
355 avg_planes = 0.0_rp
356
357 n = this%dof%size()
358 do i = 1, n
359 if (this%dir .eq. 1) coord = this%dof%x(i,1,1,1)
360 if (this%dir .eq. 2) coord = this%dof%y(i,1,1,1)
361 if (this%dir .eq. 3) coord = this%dof%z(i,1,1,1)
362 avg_planes%x(this%pt_lvl(i,1,1,1), 1) = &
363 avg_planes%x(this%pt_lvl(i,1,1,1), 1) + &
364 coord * this%coef%B(i,1,1,1) / &
365 this%volume_per_gll_lvl(this%pt_lvl(i,1,1,1))
366 end do
367 do j = 2, size(vector_ptr) + 1
368 do i = 1, n
369 avg_planes%x(this%pt_lvl(i,1,1,1), j) = &
370 avg_planes%x(this%pt_lvl(i,1,1,1), j) + &
371 vector_ptr(j-1)%ptr%x(i)*this%coef%B(i,1,1,1) &
372 /this%volume_per_gll_lvl(this%pt_lvl(i,1,1,1))
373 end do
374 end do
375 call mpi_allreduce(mpi_in_place, avg_planes%x, &
376 (size(vector_ptr) + 1) * this%n_gll_lvls, mpi_real_precision, &
377 mpi_sum, neko_comm, ierr)
378
379
380 end subroutine map_1d_average_vector_ptr
381
382end module map_1d
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Copy data between host and device (or device and device)
Definition device.F90:71
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:51
integer, public pe_size
MPI size of communicator.
Definition comm.F90:59
integer, public pe_rank
MPI rank.
Definition comm.F90:56
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:219
integer, parameter, public device_to_host
Definition device.F90:47
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Gather-scatter.
Creates a 1d GLL point map along a specified direction based on the connectivity in the mesh.
Definition map_1d.f90:3
subroutine map_1d_init(this, coef, dir, tol)
Definition map_1d.f90:73
subroutine map_1d_free(this)
Definition map_1d.f90:281
subroutine map_1d_average_field_list(this, avg_planes, field_list)
Computes average if field list in two directions and outputs matrix with averaged values avg_planes c...
Definition map_1d.f90:303
subroutine map_1d_init_char(this, coef, dir, tol)
Definition map_1d.f90:260
subroutine map_1d_average_vector_ptr(this, avg_planes, vector_ptr)
Computes average if vector_pt in two directions and outputs matrix with averaged values avg_planes co...
Definition map_1d.f90:345
Definition math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition math.f90:411
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition math.f90:797
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:516
integer function, public glimax(a, n)
Max of an integer vector of length n.
Definition math.f90:531
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:546
Defines a matrix.
Definition matrix.f90:34
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
Defines a function space.
Definition space.f90:34
Utilities.
Definition utils.f90:35
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:346
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:56
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
The function space for the SEM solution fields.
Definition space.f90:63