Neko  0.8.1
A portable framework for high-order spectral element flow simulations
symmetry.f90
Go to the documentation of this file.
1 ! Copyright (c) 2020-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 !
34 module symmetry
35  use device_symmetry
36  use neko_config
37  use num_types
38  use dirichlet
39  use bc
40  use math
41  use utils
42  use stack
43  use tuple
44  use coefs, only : coef_t
45  use, intrinsic :: iso_c_binding, only : c_ptr
46  implicit none
47  private
48 
50  type, public, extends(bc_t) :: symmetry_t
51  type(dirichlet_t) :: bc_x
52  type(dirichlet_t) :: bc_y
53  type(dirichlet_t) :: bc_z
54  contains
55  procedure, pass(this) :: init_msk => symmetry_init_msk
56  procedure, pass(this) :: apply_scalar => symmetry_apply_scalar
57  procedure, pass(this) :: apply_vector => symmetry_apply_vector
58  procedure, pass(this) :: apply_scalar_dev => symmetry_apply_scalar_dev
59  procedure, pass(this) :: apply_vector_dev => symmetry_apply_vector_dev
60  end type symmetry_t
61 
62 contains
63 
65  subroutine symmetry_init_msk(this)
66  class(symmetry_t), intent(inout) :: this
67  integer :: i, m, j, l
68  type(tuple_i4_t), pointer :: bfp(:)
69  real(kind=rp) :: sx,sy,sz
70  real(kind=rp), parameter :: tol = 1d-3
71  type(tuple_i4_t) :: bc_facet
72  integer :: facet, el
73 
74  call this%bc_x%free()
75  call this%bc_y%free()
76  call this%bc_z%free()
77  call this%bc_x%init(this%coef)
78  call this%bc_y%init(this%coef)
79  call this%bc_z%init(this%coef)
80 
81  associate(c=>this%coef, nx => this%coef%nx, ny => this%coef%ny, &
82  nz => this%coef%nz)
83  bfp => this%marked_facet%array()
84  do i = 1, this%marked_facet%size()
85  bc_facet = bfp(i)
86  facet = bc_facet%x(1)
87  el = bc_facet%x(2)
88  sx = 0d0
89  sy = 0d0
90  sz = 0d0
91  select case (facet)
92  case(1,2)
93  do l = 2, c%Xh%lx - 1
94  do j = 2, c%Xh%lx -1
95  sx = sx + abs(abs(nx(l, j, facet, el)) - 1d0)
96  sy = sy + abs(abs(ny(l, j, facet, el)) - 1d0)
97  sz = sz + abs(abs(nz(l, j, facet, el)) - 1d0)
98  end do
99  end do
100  case(3,4)
101  do l = 2, c%Xh%lx - 1
102  do j = 2, c%Xh%lx - 1
103  sx = sx + abs(abs(nx(l, j, facet, el)) - 1d0)
104  sy = sy + abs(abs(ny(l, j, facet, el)) - 1d0)
105  sz = sz + abs(abs(nz(l, j, facet, el)) - 1d0)
106  end do
107  end do
108  case(5,6)
109  do l = 2, c%Xh%lx - 1
110  do j = 2, c%Xh%lx - 1
111  sx = sx + abs(abs(nx(l, j, facet, el)) - 1d0)
112  sy = sy + abs(abs(ny(l, j, facet, el)) - 1d0)
113  sz = sz + abs(abs(nz(l, j, facet, el)) - 1d0)
114  end do
115  end do
116  end select
117  sx = sx / (c%Xh%lx - 2)**2
118  sy = sy / (c%Xh%lx - 2)**2
119  sz = sz / (c%Xh%lx - 2)**2
120 
121  if (sx .lt. tol) then
122  call this%bc_x%mark_facet(facet, el)
123  end if
124 
125  if (sy .lt. tol) then
126  call this%bc_y%mark_facet(facet, el)
127  end if
128 
129  if (sz .lt. tol) then
130  call this%bc_z%mark_facet(facet, el)
131  end if
132  end do
133  end associate
134  call this%bc_x%finalize()
135  call this%bc_x%set_g(0.0_rp)
136  call this%bc_y%finalize()
137  call this%bc_y%set_g(0.0_rp)
138  call this%bc_z%finalize()
139  call this%bc_z%set_g(0.0_rp)
140 
141  end subroutine symmetry_init_msk
142 
143  subroutine symmetry_free(this)
144  type(symmetry_t), intent(inout) :: this
145 
146  call this%bc_x%free()
147  call this%bc_y%free()
148  call this%bc_z%free()
149 
150  end subroutine symmetry_free
151 
153  subroutine symmetry_apply_scalar(this, x, n, t, tstep)
154  class(symmetry_t), intent(inout) :: this
155  integer, intent(in) :: n
156  real(kind=rp), intent(inout), dimension(n) :: x
157  real(kind=rp), intent(in), optional :: t
158  integer, intent(in), optional :: tstep
159  end subroutine symmetry_apply_scalar
160 
162  subroutine symmetry_apply_vector(this, x, y, z, n, t, tstep)
163  class(symmetry_t), intent(inout) :: this
164  integer, intent(in) :: n
165  real(kind=rp), intent(inout), dimension(n) :: x
166  real(kind=rp), intent(inout), dimension(n) :: y
167  real(kind=rp), intent(inout), dimension(n) :: z
168  real(kind=rp), intent(in), optional :: t
169  integer, intent(in), optional :: tstep
170  integer :: i, m, k
171 
172  call this%bc_x%apply_scalar(x,n)
173  call this%bc_y%apply_scalar(y,n)
174  call this%bc_z%apply_scalar(z,n)
175 
176  end subroutine symmetry_apply_vector
177 
179  subroutine symmetry_apply_scalar_dev(this, x_d, t, tstep)
180  class(symmetry_t), intent(inout), target :: this
181  type(c_ptr) :: x_d
182  real(kind=rp), intent(in), optional :: t
183  integer, intent(in), optional :: tstep
184  end subroutine symmetry_apply_scalar_dev
185 
187  subroutine symmetry_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
188  class(symmetry_t), intent(inout), target :: this
189  type(c_ptr) :: x_d
190  type(c_ptr) :: y_d
191  type(c_ptr) :: z_d
192  real(kind=rp), intent(in), optional :: t
193  integer, intent(in), optional :: tstep
194 
195  call device_symmetry_apply_vector(this%bc_x%msk_d, this%bc_y%msk_d, &
196  this%bc_z%msk_d, x_d, y_d, z_d, &
197  this%bc_x%msk(0), &
198  this%bc_y%msk(0), &
199  this%bc_z%msk(0))
200 
201  end subroutine symmetry_apply_vector_dev
202 
203 end module symmetry
Defines a boundary condition.
Definition: bc.f90:34
Coefficients.
Definition: coef.f90:34
Defines a dirichlet boundary condition.
Definition: dirichlet.f90:34
Definition: math.f90:60
Build configurations.
Definition: neko_config.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Implements a dynamic stack ADT.
Definition: stack.f90:35
Mixed Dirichlet-Neumann axis aligned symmetry plane.
Definition: symmetry.f90:34
subroutine symmetry_init_msk(this)
Initialize symmetry mask for each axis.
Definition: symmetry.f90:66
subroutine symmetry_apply_vector(this, x, y, z, n, t, tstep)
Apply symmetry conditions (axis aligned)
Definition: symmetry.f90:163
subroutine symmetry_apply_scalar(this, x, n, t, tstep)
No-op scalar apply.
Definition: symmetry.f90:154
subroutine symmetry_apply_scalar_dev(this, x_d, t, tstep)
No-op scalar apply (device version)
Definition: symmetry.f90:180
subroutine symmetry_free(this)
Definition: symmetry.f90:144
subroutine symmetry_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
Apply symmetry conditions (axis aligned) (device version)
Definition: symmetry.f90:188
Implements a n-tuple.
Definition: tuple.f90:34
Utilities.
Definition: utils.f90:35
Base type for a boundary condition.
Definition: bc.f90:51
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:54
Generic Dirichlet boundary condition on .
Definition: dirichlet.f90:44
Mixed Dirichlet-Neumann symmetry plane condition.
Definition: symmetry.f90:50
Integer based 2-tuple.
Definition: tuple.f90:51