Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
symmetry.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2021, 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!
36 use dirichlet, only : dirichlet_t
37 use num_types, only : rp
38 use bc, only : bc_t
39 use tuple, only : tuple_i4_t
40 use coefs, only : coef_t
41 use json_module, only : json_file
43 use, intrinsic :: iso_c_binding, only : c_ptr
44 implicit none
45 private
46
49 type, public, extends(bc_t) :: symmetry_t
50 type(zero_dirichlet_t) :: bc_x
51 type(zero_dirichlet_t) :: bc_y
52 type(zero_dirichlet_t) :: bc_z
53 contains
54 procedure, pass(this) :: apply_scalar => symmetry_apply_scalar
55 procedure, pass(this) :: apply_vector => symmetry_apply_vector
56 procedure, pass(this) :: apply_scalar_dev => symmetry_apply_scalar_dev
57 procedure, pass(this) :: apply_vector_dev => symmetry_apply_vector_dev
59 procedure, pass(this) :: init => symmetry_init
61 procedure, pass(this) :: init_from_components => &
64 procedure, pass(this) :: free => symmetry_free
66 procedure, pass(this) :: get_normal_axis => symmetry_get_normal_axis
68 procedure, pass(this) :: finalize => symmetry_finalize
69 end type symmetry_t
70
71contains
75 subroutine symmetry_init(this, coef, json)
76 class(symmetry_t), intent(inout), target :: this
77 type(coef_t), intent(in) :: coef
78 type(json_file), intent(inout) ::json
79
80 call this%init_from_components(coef)
81
82 end subroutine symmetry_init
83
86 subroutine symmetry_init_from_components(this, coef)
87 class(symmetry_t), intent(inout), target :: this
88 type(coef_t), intent(in) :: coef
89
90 call this%free()
91 this%strong = .false.
92
93 call this%init_base(coef)
94 call this%bc_x%init_from_components(this%coef)
95 call this%bc_y%init_from_components(this%coef)
96 call this%bc_z%init_from_components(this%coef)
97
99
103 subroutine symmetry_finalize(this)
104 class(symmetry_t), target, intent(inout) :: this
105 integer :: i, m, j, l
106 type(tuple_i4_t), pointer :: bfp(:)
107 real(kind=rp) :: sx, sy, sz
108 real(kind=rp), parameter :: tol = 1e-3_rp
109 type(tuple_i4_t) :: bc_facet
110 integer :: facet, el
111
112 call this%finalize_base()
113
114 associate(c => this%coef, nx => this%coef%nx, ny => this%coef%ny, &
115 nz => this%coef%nz)
116 bfp => this%marked_facet%array()
117 do i = 1, this%marked_facet%size()
118 bc_facet = bfp(i)
119 facet = bc_facet%x(1)
120 el = bc_facet%x(2)
121 call this%get_normal_axis(sx, sy, sz, facet, el)
122
123 if (sx .lt. tol) then
124 call this%bc_x%mark_facet(facet, el)
125 end if
126
127 if (sy .lt. tol) then
128 call this%bc_y%mark_facet(facet, el)
129 end if
130
131 if (sz .lt. tol) then
132 call this%bc_z%mark_facet(facet, el)
133 end if
134 end do
135 end associate
136 call this%bc_x%finalize()
137 call this%bc_y%finalize()
138 call this%bc_z%finalize()
139 end subroutine symmetry_finalize
140
144 subroutine symmetry_get_normal_axis(this, sx, sy, sz, facet, el)
145 class(symmetry_t), target, intent(inout) :: this
146 real(kind=rp), intent(out) :: sx, sy, sz
147 integer, intent(in) :: facet
148 integer, intent(in) :: el
149 integer :: i, m, j, l
150 type(tuple_i4_t), pointer :: bfp(:)
151 ! TODO is the tolerance too large?
152 real(kind=rp), parameter :: tol = 1d-3
153 type(tuple_i4_t) :: bc_facet
154
155 associate(c => this%coef, nx => this%coef%nx, ny => this%coef%ny, &
156 nz => this%coef%nz)
157 sx = 0.0_rp
158 sy = 0d0
159 sz = 0d0
160 select case (facet)
161 case (1, 2)
162 do l = 2, c%Xh%lx - 1
163 do j = 2, c%Xh%lx -1
164 sx = sx + abs(abs(nx(l, j, facet, el)) - 1d0)
165 sy = sy + abs(abs(ny(l, j, facet, el)) - 1d0)
166 sz = sz + abs(abs(nz(l, j, facet, el)) - 1d0)
167 end do
168 end do
169 case (3, 4)
170 do l = 2, c%Xh%lx - 1
171 do j = 2, c%Xh%lx - 1
172 sx = sx + abs(abs(nx(l, j, facet, el)) - 1d0)
173 sy = sy + abs(abs(ny(l, j, facet, el)) - 1d0)
174 sz = sz + abs(abs(nz(l, j, facet, el)) - 1d0)
175 end do
176 end do
177 case (5, 6)
178 do l = 2, c%Xh%lx - 1
179 do j = 2, c%Xh%lx - 1
180 sx = sx + abs(abs(nx(l, j, facet, el)) - 1d0)
181 sy = sy + abs(abs(ny(l, j, facet, el)) - 1d0)
182 sz = sz + abs(abs(nz(l, j, facet, el)) - 1d0)
183 end do
184 end do
185 end select
186 sx = sx / (c%Xh%lx - 2)**2
187 sy = sy / (c%Xh%lx - 2)**2
188 sz = sz / (c%Xh%lx - 2)**2
189 end associate
190 end subroutine symmetry_get_normal_axis
191
198 subroutine symmetry_apply_scalar(this, x, n, t, tstep, strong)
199 class(symmetry_t), intent(inout) :: this
200 integer, intent(in) :: n
201 real(kind=rp), intent(inout), dimension(n) :: x
202 real(kind=rp), intent(in), optional :: t
203 integer, intent(in), optional :: tstep
204 logical, intent(in), optional :: strong
205 end subroutine symmetry_apply_scalar
206
215 subroutine symmetry_apply_vector(this, x, y, z, n, t, tstep, strong)
216 class(symmetry_t), intent(inout) :: this
217 integer, intent(in) :: n
218 real(kind=rp), intent(inout), dimension(n) :: x
219 real(kind=rp), intent(inout), dimension(n) :: y
220 real(kind=rp), intent(inout), dimension(n) :: z
221 real(kind=rp), intent(in), optional :: t
222 integer, intent(in), optional :: tstep
223 logical, intent(in), optional :: strong
224 logical :: strong_ = .true.
225
226 if (present(strong)) strong_ = strong
227
228 if (strong_) then
229 call this%bc_x%apply_scalar(x, n)
230 call this%bc_y%apply_scalar(y, n)
231 call this%bc_z%apply_scalar(z, n)
232 end if
233
234 end subroutine symmetry_apply_vector
235
241 subroutine symmetry_apply_scalar_dev(this, x_d, t, tstep, strong)
242 class(symmetry_t), intent(inout), target :: this
243 type(c_ptr) :: x_d
244 real(kind=rp), intent(in), optional :: t
245 integer, intent(in), optional :: tstep
246 logical, intent(in), optional :: strong
247 end subroutine symmetry_apply_scalar_dev
248
256 subroutine symmetry_apply_vector_dev(this, x_d, y_d, z_d, t, tstep, strong)
257 class(symmetry_t), intent(inout), target :: this
258 type(c_ptr) :: x_d
259 type(c_ptr) :: y_d
260 type(c_ptr) :: z_d
261 real(kind=rp), intent(in), optional :: t
262 integer, intent(in), optional :: tstep
263 logical, intent(in), optional :: strong
264 logical :: strong_ = .true.
265
266 if (present(strong)) strong_ = strong
267
268 if (strong_ .and. (this%msk(0) .gt. 0)) then
269 call device_symmetry_apply_vector(this%bc_x%msk_d, this%bc_y%msk_d, &
270 this%bc_z%msk_d, x_d, y_d, z_d, &
271 this%bc_x%msk(0), this%bc_y%msk(0), this%bc_z%msk(0))
272 end if
273 end subroutine symmetry_apply_vector_dev
274
276 subroutine symmetry_free(this)
277 class(symmetry_t), target, intent(inout) :: this
278
279 call this%free_base()
280 call this%bc_x%free()
281 call this%bc_y%free()
282 call this%bc_z%free()
283
284 end subroutine symmetry_free
285end module symmetry
Defines a boundary condition.
Definition bc.f90:34
Coefficients.
Definition coef.f90:34
subroutine device_symmetry_apply_vector(xmsk, ymsk, zmsk, x, y, z, m, n, l)
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Mixed Dirichlet-Neumann axis aligned symmetry plane.
Definition symmetry.f90:34
subroutine symmetry_init_from_components(this, coef)
Constructor from components.
Definition symmetry.f90:87
subroutine symmetry_apply_scalar(this, x, n, t, tstep, strong)
No-op scalar apply.
Definition symmetry.f90:199
subroutine symmetry_init(this, coef, json)
Constructor.
Definition symmetry.f90:76
subroutine symmetry_apply_scalar_dev(this, x_d, t, tstep, strong)
No-op scalar apply (device version)
Definition symmetry.f90:242
subroutine symmetry_finalize(this)
Finalize. Marks the appropriate faces for application of a homogeneous Dirchlet condition based on th...
Definition symmetry.f90:104
subroutine symmetry_apply_vector_dev(this, x_d, y_d, z_d, t, tstep, strong)
Apply symmetry conditions (axis aligned) (device version)
Definition symmetry.f90:257
subroutine symmetry_apply_vector(this, x, y, z, n, t, tstep, strong)
Apply symmetry conditions (axis aligned)
Definition symmetry.f90:216
subroutine symmetry_free(this)
Destructor.
Definition symmetry.f90:277
subroutine symmetry_get_normal_axis(this, sx, sy, sz, facet, el)
Compute the average normal for a facet of an element.
Definition symmetry.f90:145
Implements a n-tuple.
Definition tuple.f90:34
Defines a zero-valued Dirichlet boundary condition.
Base type for a boundary condition.
Definition bc.f90:54
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:47
Mixed Dirichlet-Neumann symmetry plane condition.
Definition symmetry.f90:49
Integer based 2-tuple.
Definition tuple.f90:51
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...