Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
adv_dummy.f90
Go to the documentation of this file.
1! Copyright (c) 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 space, only : space_t
38 use field, only : field_t
39 use coefs, only : coef_t
40 implicit none
41 private
42
44 type, public, extends(advection_t) :: adv_dummy_t
45 contains
47 procedure, pass(this) :: init => init_adv_dummy
49 procedure, pass(this) :: free => free_adv_dummy
52 procedure, pass(this) :: compute => compute_adv_dummy
55 procedure, pass(this) :: compute_scalar => &
59 procedure, pass(this) :: compute_ale => compute_ale_adv_dummy
61 procedure, pass(this) :: recompute_metrics => recompute_metrics_dummy_noop
62 end type adv_dummy_t
63
64contains
65
68 subroutine init_adv_dummy(this, coef)
69 class(adv_dummy_t), intent(inout) :: this
70 type(coef_t), intent(in) :: coef
71
72 end subroutine init_adv_dummy
73
75 subroutine free_adv_dummy(this)
76 class(adv_dummy_t), intent(inout) :: this
77
78 end subroutine free_adv_dummy
79
93 subroutine compute_adv_dummy(this, vx, vy, vz, fx, fy, fz, Xh, &
94 coef, n, dt)
95 class(adv_dummy_t), intent(inout) :: this
96 type(space_t), intent(in) :: Xh
97 type(coef_t), intent(in) :: coef
98 type(field_t), intent(inout) :: vx, vy, vz
99 type(field_t), intent(inout) :: fx, fy, fz
100 integer, intent(in) :: n
101 real(kind=rp), intent(in), optional :: dt
102
103 end subroutine compute_adv_dummy
104
118 subroutine compute_scalar_adv_dummy(this, vx, vy, vz, s, fs, Xh, &
119 coef, n, dt)
120 class(adv_dummy_t), intent(inout) :: this
121 type(field_t), intent(inout) :: vx, vy, vz
122 type(field_t), intent(inout) :: s
123 type(field_t), intent(inout) :: fs
124 type(space_t), intent(in) :: Xh
125 type(coef_t), intent(in) :: coef
126 integer, intent(in) :: n
127 real(kind=rp), intent(in), optional :: dt
128
129 end subroutine compute_scalar_adv_dummy
130
131 subroutine recompute_metrics_dummy_noop(this, coef, moving_boundary)
132 class(adv_dummy_t), intent(inout) :: this
133 type(coef_t), intent(in) :: coef
134 logical, intent(in) :: moving_boundary
135 ! no-op
136 end subroutine recompute_metrics_dummy_noop
137
138 subroutine compute_ale_adv_dummy(this, vx, vy, vz, wm_x, wm_y, wm_z, &
139 fx, fy, fz, Xh, coef, n, dt)
140 class(adv_dummy_t), intent(inout) :: this
141 type(field_t), intent(inout) :: vx, vy, vz
142 type(field_t), intent(inout) :: wm_x, wm_y, wm_z
143 type(field_t), intent(inout) :: fx, fy, fz
144 type(space_t), intent(in) :: Xh
145 type(coef_t), intent(in) :: coef
146 integer, intent(in) :: n
147 real(kind=rp), intent(in), optional :: dt
148 ! no-op
149 end subroutine compute_ale_adv_dummy
150end module adv_dummy
Implements adv_dummy_t
Definition adv_dummy.f90:34
subroutine compute_adv_dummy(this, vx, vy, vz, fx, fy, fz, xh, coef, n, dt)
Add the advection term for the fluid, i.e. to the RHS.
Definition adv_dummy.f90:95
subroutine free_adv_dummy(this)
Destructor.
Definition adv_dummy.f90:76
subroutine init_adv_dummy(this, coef)
Constructor.
Definition adv_dummy.f90:69
subroutine recompute_metrics_dummy_noop(this, coef, moving_boundary)
subroutine compute_ale_adv_dummy(this, vx, vy, vz, wm_x, wm_y, wm_z, fx, fy, fz, xh, coef, n, dt)
subroutine compute_scalar_adv_dummy(this, vx, vy, vz, s, fs, xh, coef, n, dt)
Add the advection term for a scalar, 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
Defines a field.
Definition field.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a function space.
Definition space.f90:34
A zero-valued advection that can be used to kill the advection term.
Definition adv_dummy.f90:44
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