Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
adv_no_dealias.f90
Go to the documentation of this file.
1! Copyright (c) 2021-2024, 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 advection, only : advection_t
36 use num_types, only : rp
37 use math, only : subcol3, rzero, addcol3
38 use field_math, only : field_col3
39 use space, only : space_t
40 use field, only : field_t
41 use coefs, only : coef_t
44 use operators, only : conv1, div
47 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr, c_associated
48 implicit none
49 private
50
52 type, public, extends(advection_t) :: adv_no_dealias_t
53 real(kind=rp), allocatable :: temp(:)
54 type(c_ptr) :: temp_d = c_null_ptr
55 contains
57 procedure, pass(this) :: init => init_no_dealias
59 procedure, pass(this) :: free => free_no_dealias
62 procedure, pass(this) :: compute => compute_advection_no_dealias
65 procedure, pass(this) :: compute_scalar => &
69 procedure, pass(this) :: compute_ale => compute_ale_advection_no_dealias
71 procedure, pass(this) :: recompute_metrics => recompute_metrics_no_dealias
72 end type adv_no_dealias_t
73
74contains
75
78 subroutine init_no_dealias(this, coef)
79 class(adv_no_dealias_t), intent(inout) :: this
80 type(coef_t), intent(in) :: coef
81
82 allocate(this%temp(coef%dof%size()))
83
84 if (neko_bcknd_device .eq. 1) then
85 call device_map(this%temp, this%temp_d, coef%dof%size())
86 end if
87
88 end subroutine init_no_dealias
89
91 subroutine free_no_dealias(this)
92 class(adv_no_dealias_t), intent(inout) :: this
93
94 if (allocated(this%temp)) then
95 deallocate(this%temp)
96 end if
97 if (c_associated(this%temp_d)) then
98 call device_free(this%temp_d)
99 end if
100 end subroutine free_no_dealias
101
114 subroutine compute_advection_no_dealias(this, vx, vy, vz, fx, fy, fz, Xh, &
115 coef, n, dt)
116 class(adv_no_dealias_t), intent(inout) :: this
117 type(space_t), intent(in) :: Xh
118 type(coef_t), intent(in) :: coef
119 type(field_t), intent(inout) :: vx, vy, vz
120 type(field_t), intent(inout) :: fx, fy, fz
121 integer, intent(in) :: n
122 real(kind=rp), intent(in), optional :: dt
123
124 if (neko_bcknd_device .eq. 1) then
125
126 call conv1(this%temp, vx%x, vx%x, vy%x, vz%x, xh, coef)
127 call device_subcol3 (fx%x_d, coef%B_d, this%temp_d, n)
128 call conv1(this%temp, vy%x, vx%x, vy%x, vz%x, xh, coef)
129 call device_subcol3 (fy%x_d, coef%B_d, this%temp_d, n)
130 if (coef%Xh%lz .eq. 1) then
131 call device_rzero (this%temp_d, n)
132 else
133 call conv1(this%temp, vz%x, vx%x, vy%x, vz%x, xh, coef)
134 call device_subcol3(fz%x_d, coef%B_d, this%temp_d, n)
135 end if
136 else
137 call conv1(this%temp, vx%x, vx%x, vy%x, vz%x, xh, coef)
138 call subcol3 (fx%x, coef%B, this%temp, n)
139 call conv1(this%temp, vy%x, vx%x, vy%x, vz%x, xh, coef)
140 call subcol3 (fy%x, coef%B, this%temp, n)
141 if (coef%Xh%lz .eq. 1) then
142 call rzero (this%temp, n)
143 else
144 call conv1(this%temp, vz%x, vx%x, vy%x, vz%x, xh, coef)
145 call subcol3(fz%x, coef%B, this%temp, n)
146 end if
147 end if
148
149 end subroutine compute_advection_no_dealias
150
163 subroutine compute_scalar_advection_no_dealias(this, vx, vy, vz, s, fs, Xh, &
164 coef, n, dt)
165 class(adv_no_dealias_t), intent(inout) :: this
166 type(field_t), intent(inout) :: vx, vy, vz
167 type(field_t), intent(inout) :: s
168 type(field_t), intent(inout) :: fs
169 type(space_t), intent(in) :: Xh
170 type(coef_t), intent(in) :: coef
171 integer, intent(in) :: n
172 real(kind=rp), intent(in), optional :: dt
173
174 if (neko_bcknd_device .eq. 1) then
175
176 call conv1(this%temp, s%x, vx%x, vy%x, vz%x, xh, coef)
177 call device_subcol3 (fs%x_d, coef%B_d, this%temp_d, n)
178 if (coef%Xh%lz .eq. 1) then
179 call device_rzero (this%temp_d, n)
180 end if
181 else
182 ! temp will hold vx*ds/dx + vy*ds/dy + vz*ds/sz
183 call conv1(this%temp, s%x, vx%x, vy%x, vz%x, xh, coef)
184
185 ! fs = fs - B*temp
186 call subcol3 (fs%x, coef%B, this%temp, n)
187 if (coef%Xh%lz .eq. 1) then
188 call rzero (this%temp, n)
189 end if
190 end if
191
193
194 subroutine recompute_metrics_no_dealias(this, coef, moving_boundary)
195 class(adv_no_dealias_t), intent(inout) :: this
196 type(coef_t), intent(in) :: coef
197 logical, intent(in) :: moving_boundary
198 ! no-op
199 end subroutine recompute_metrics_no_dealias
200
201 subroutine compute_ale_advection_no_dealias(this, vx, vy, vz, &
202 wm_x, wm_y, wm_z, fx, fy, fz, Xh, coef, n, dt)
203 class(adv_no_dealias_t), intent(inout) :: this
204 type(field_t), intent(inout) :: vx, vy, vz
205 type(field_t), intent(inout) :: wm_x, wm_y, wm_z
206 type(field_t), intent(inout) :: fx, fy, fz
207 type(space_t), intent(in) :: Xh
208 type(coef_t), intent(in) :: coef
209 type(field_t), pointer :: work_x, work_y, work_z
210 integer :: id_x, id_y, id_z
211 integer, intent(in) :: n
212 real(kind=rp), intent(in), optional :: dt
213
214 call neko_scratch_registry%request_field(work_x, id_x, .false.)
215 call neko_scratch_registry%request_field(work_y, id_y, .false.)
216 call neko_scratch_registry%request_field(work_z, id_z, .false.)
217
218 ! u.wm_*
219 call field_col3(work_x, vx, wm_x)
220 call field_col3(work_y, vx, wm_y)
221 call field_col3(work_z, vx, wm_z)
222 call div(this%temp, work_x%x, work_y%x, work_z%x, coef)
223
224 if (neko_bcknd_device .eq. 1) then
225 call device_addcol3(fx%x_d, coef%B_d, this%temp_d, n)
226 else
227 call addcol3(fx%x, coef%B, this%temp, n)
228 end if
229
230 ! v.wm_*
231 call field_col3(work_x, vy, wm_x)
232 call field_col3(work_y, vy, wm_y)
233 call field_col3(work_z, vy, wm_z)
234 call div(this%temp, work_x%x, work_y%x, work_z%x, coef)
235
236 if (neko_bcknd_device .eq. 1) then
237 call device_addcol3(fy%x_d, coef%B_d, this%temp_d, n)
238 else
239 call addcol3(fy%x, coef%B, this%temp, n)
240 end if
241
242 ! w.wm_*
243 call field_col3(work_x, vz, wm_x)
244 call field_col3(work_y, vz, wm_y)
245 call field_col3(work_z, vz, wm_z)
246 call div(this%temp, work_x%x, work_y%x, work_z%x, coef)
247
248 if (neko_bcknd_device .eq. 1) then
249 call device_addcol3(fz%x_d, coef%B_d, this%temp_d, n)
250 else
251 call addcol3(fz%x, coef%B, this%temp, n)
252 end if
253
254 call neko_scratch_registry%relinquish_field(id_x)
255 call neko_scratch_registry%relinquish_field(id_y)
256 call neko_scratch_registry%relinquish_field(id_z)
257
259end module adv_no_dealias
Return the device pointer for an associated Fortran array.
Definition device.F90:107
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Subroutines to add advection terms to the RHS of a transport equation.
subroutine init_no_dealias(this, coef)
Constructor.
subroutine free_no_dealias(this)
Destructor.
subroutine compute_ale_advection_no_dealias(this, vx, vy, vz, wm_x, wm_y, wm_z, fx, fy, fz, xh, coef, n, dt)
subroutine compute_scalar_advection_no_dealias(this, vx, vy, vz, s, fs, xh, coef, n, dt)
Add the advection term for a scalar, i.e. , to the RHS.
subroutine recompute_metrics_no_dealias(this, coef, moving_boundary)
subroutine compute_advection_no_dealias(this, vx, vy, vz, fx, fy, fz, xh, coef, n, dt)
Add the advection term for the fluid, i.e. to the RHS.
Subroutines to add advection terms to the RHS of a transport equation.
Definition advection.f90:34
Coefficients.
Definition coef.f90:34
subroutine, public device_addcol3(a_d, b_d, c_d, n, strm)
Returns .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_subcol3(a_d, b_d, c_d, n, strm)
Returns .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:225
subroutine, public field_col3(a, b, c, n)
Vector multiplication with 3 vectors .
Defines a field.
Definition field.f90:34
Definition math.f90:60
subroutine, public subcol3(a, b, c, n)
Returns .
Definition math.f90:882
subroutine, public addcol3(a, b, c, n)
Returns .
Definition math.f90:959
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:206
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public div(res, ux, uy, uz, coef)
Compute the divergence of a vector field.
subroutine, public conv1(du, u, vx, vy, vz, xh, coef, es, ee)
Compute the advection term.
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Defines a function space.
Definition space.f90:34
Type encapsulating advection routines with no dealiasing applied.
Base abstract type for computing the advection operator.
Definition advection.f90:46
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:62
The function space for the SEM solution fields.
Definition space.f90:63