Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
interpolation.f90
Go to the documentation of this file.
1! Copyright (c) 2021-2023, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
35 use neko_config
36 use num_types, only : rp
37 use device
38 use fast3d
39 use tensor, only : tnsr3d
40 use tensor_cpu, only : tnsr3d_cpu
41 use space, only : space_t, operator(.eq.), gl, gll
42 use utils, only : neko_error
43 use, intrinsic :: iso_c_binding
44 implicit none
45 private
46
53 type, public :: interpolator_t
55 type(space_t), pointer :: xh
57 type(space_t), pointer :: yh
59 real(kind=rp), allocatable :: xh_to_yh(:,:), xh_to_yht(:,:)
61 real(kind=rp), allocatable :: yh_to_xh(:,:), yh_to_xht(:,:)
63 type(c_ptr) :: xh_yh_d = c_null_ptr
65 type(c_ptr) :: xh_yht_d = c_null_ptr
67 type(c_ptr) :: yh_xh_d = c_null_ptr
69 type(c_ptr) :: yh_xht_d = c_null_ptr
70 contains
71
73 procedure, pass(this) :: init => interpolator_init
75 procedure, pass(this) :: free => interpolator_free
77 procedure, pass(this) :: map => interpolator_map
79 procedure, pass(this) :: map_host => interpolator_map_host
80
81 end type interpolator_t
82
83contains
84
88 subroutine interpolator_init(this, Xh, Yh)
89 class(interpolator_t), intent(inout), target :: this
90 type(space_t), intent(inout), target :: Xh
91 type(space_t), intent(inout), target :: Yh
92 integer :: deg_derivate
93
94 call this%free()
95
96 allocate(this%Xh_to_Yh(yh%lx,xh%lx))
97 allocate(this%Xh_to_YhT(xh%lx,yh%lx))
98 allocate(this%Yh_to_Xh(xh%lx,yh%lx))
99 allocate(this%Yh_to_XhT(yh%lx,xh%lx))
100
101 if (xh%t .eq. gll .and. yh%t .eq. gll) then
102 else if ((xh%t .eq. gl .and. yh%t .eq. gll) .or. &
103 (yh%t .eq. gl .and. xh%t .eq. gll)) then
104 else
105 call neko_error('Unsupported interpolation')
106 end if
107 deg_derivate = 0
108 call setup_intp(this%Xh_to_Yh, this%Xh_to_YhT, &
109 yh%zg, xh%zg, yh%lx, xh%lx, deg_derivate)
110 call setup_intp(this%Yh_to_Xh, this%Yh_to_XhT, &
111 xh%zg, yh%zg, xh%lx, yh%lx, deg_derivate)
112
113 this%Xh => xh
114 this%Yh => yh
115 if (neko_bcknd_device .eq. 1) then
116 call device_map(this%Xh_to_Yh, this%Xh_Yh_d, yh%lx*xh%lx)
117 call device_map(this%Xh_to_YhT, this%Xh_YhT_d, yh%lx*xh%lx)
118 call device_map(this%Yh_to_Xh, this%Yh_Xh_d, yh%lx*xh%lx)
119 call device_map(this%Yh_to_XhT, this%Yh_XhT_d, yh%lx*xh%lx)
120 call device_memcpy(this%Xh_to_Yh, this%Xh_Yh_d, yh%lx*xh%lx, &
121 host_to_device, sync=.false.)
122 call device_memcpy(this%Xh_to_YhT, this%Xh_YhT_d, yh%lx*xh%lx, &
123 host_to_device, sync=.false.)
124 call device_memcpy(this%Yh_to_Xh, this%Yh_Xh_d, yh%lx*xh%lx, &
125 host_to_device, sync=.false.)
126 call device_memcpy(this%Yh_to_XhT, this%Yh_XhT_d, yh%lx*xh%lx, &
127 host_to_device, sync=.false.)
128 end if
129
130 end subroutine interpolator_init
131
132 subroutine interpolator_free(this)
133 class(interpolator_t), intent(inout) :: this
134
135 if (allocated(this%Xh_to_Yh)) then
136 deallocate(this%Xh_to_Yh)
137 end if
138 if (allocated(this%Xh_to_YhT)) then
139 deallocate(this%Xh_to_YhT)
140 end if
141 if (allocated(this%Yh_to_Xh)) then
142 deallocate(this%Yh_to_Xh)
143 end if
144 if (allocated(this%Yh_to_XhT)) then
145 deallocate(this%Yh_to_XhT)
146 end if
147 if (c_associated(this%Yh_Xh_d)) then
148 call device_free(this%Yh_Xh_d)
149 end if
150 if (c_associated(this%Yh_XhT_d)) then
151 call device_free(this%Yh_XhT_d)
152 end if
153 if (c_associated(this%Xh_Yh_d)) then
154 call device_free(this%Xh_Yh_d)
155 end if
156 if (c_associated(this%Xh_YhT_d)) then
157 call device_free(this%Xh_YhT_d)
158 end if
159
160 end subroutine interpolator_free
161
167 subroutine interpolator_map(this, y, x, nel, to_space)
168 class(interpolator_t), intent(inout) :: this
169 integer :: nel
170 type(space_t) :: to_space
171 real(kind=rp), intent(in) :: x(this%Xh%lx, this%Xh%lx, this%Xh%lx, nel)
172 real(kind=rp), intent(inout) :: y(this%Yh%lx, this%Yh%lx, this%Yh%lx, nel)
173 if (to_space .eq. this%Yh) then
174 call tnsr3d(y, this%Yh%lx, x, &
175 this%Xh%lx,this%Yh_to_XhT, &
176 this%Yh_to_Xh, this%Yh_to_Xh, nel)
177 else if (to_space .eq. this%Xh) then
178 call tnsr3d(y, this%Xh%lx, x, &
179 this%Yh%lx,this%Yh_to_Xh, &
180 this%Yh_to_XhT, this%Yh_to_XhT, nel)
181 else
182 call neko_error('Invalid interpolation')
183 end if
184 end subroutine interpolator_map
185
186
193 subroutine interpolator_map_host(this, y, x, nel, to_space)
194 class(interpolator_t), intent(inout) :: this
195 integer :: nel
196 type(space_t) :: to_space
197 real(kind=rp), intent(inout) :: x(this%Xh%lx, this%Xh%lx, this%Xh%lx, nel)
198 real(kind=rp), intent(inout) :: y(this%Yh%lx, this%Yh%lx, this%Yh%lx, nel)
199 if (to_space .eq. this%Yh) then
200 call tnsr3d_cpu(y, this%Yh%lx, x, &
201 this%Xh%lx,this%Yh_to_XhT, &
202 this%Yh_to_Xh, this%Yh_to_Xh, nel)
203 else if (to_space .eq. this%Xh) then
204 call tnsr3d_cpu(y, this%Xh%lx, x, &
205 this%Yh%lx,this%Yh_to_Xh, &
206 this%Yh_to_XhT, this%Yh_to_XhT, nel)
207 else
208 call neko_error('Invalid interpolation')
209 end if
210 end subroutine interpolator_map_host
211
212
213end module interpolation
Map a Fortran array to a device (allocate and associate)
Definition device.F90:71
Copy data between host and device (or device and device)
Definition device.F90:65
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:46
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:200
Fast diagonalization methods from NEKTON.
Definition fast3d.f90:61
subroutine, public setup_intp(jh, jht, z_to, z_from, n_to, n_from, derivative)
Compute interpolation weights for points z_to using values at points z_from.
Definition fast3d.f90:243
Routines to interpolate between different spaces.
subroutine interpolator_map_host(this, y, x, nel, to_space)
Interpolates an array to one of Xh or Yh on host.
subroutine interpolator_map(this, y, x, nel, to_space)
Interpolates an array to one of Xh or Yh.
subroutine interpolator_init(this, xh, yh)
Constructor to initialize with two different spaces.
subroutine interpolator_free(this)
NEKTON map.
Definition map.f90:3
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
integer, parameter, public gll
Definition space.f90:48
integer, parameter, public gl
Definition space.f90:48
subroutine, public tnsr3d_cpu(v, nv, u, nu, a, bt, ct, nelv)
Tensor operations.
Definition tensor.f90:61
subroutine, public tnsr3d(v, nv, u, nu, a, bt, ct, nelv)
Tensor product performed on nelv elements.
Definition tensor.f90:223
Utilities.
Definition utils.f90:35
Interpolation between two space::space_t.
The function space for the SEM solution fields.
Definition space.f90:62