Neko  0.8.1
A portable framework for high-order spectral element flow simulations
non_normal.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 non_normal
35  use symmetry
36  use neko_config
37  use num_types
38  use dirichlet
39  use tuple
40  use device
41  use coefs
42  use math
43  use utils
44  use stack
45  use, intrinsic :: iso_c_binding
46  implicit none
47  private
48 
50  type, public, extends(symmetry_t) :: non_normal_t
51  contains
52  procedure, pass(this) :: init_msk => non_normal_init_msk
53  end type non_normal_t
54 
55 contains
56 
58  subroutine non_normal_init_msk(this)
59  class(non_normal_t), intent(inout) :: this
60  integer :: i, j, k, l
61  type(tuple_i4_t), pointer :: bfp(:)
62  real(kind=rp) :: sx,sy,sz
63  real(kind=rp), parameter :: tol = 1d-3
64  type(tuple_i4_t) :: bc_facet
65  integer :: facet, el
66 
67  call this%bc_x%free()
68  call this%bc_y%free()
69  call this%bc_z%free()
70  call this%bc_x%init(this%coef)
71  call this%bc_y%init(this%coef)
72  call this%bc_z%init(this%coef)
73 
74  associate(c=>this%coef, nx => this%coef%nx, ny => this%coef%ny, &
75  nz => this%coef%nz)
76  bfp => this%marked_facet%array()
77  do i = 1, this%marked_facet%size()
78  bc_facet = bfp(i)
79  facet = bc_facet%x(1)
80  el = bc_facet%x(2)
81  sx = 0d0
82  sy = 0d0
83  sz = 0d0
84  select case (facet)
85  case(1,2)
86  do l = 2, c%Xh%lx - 1
87  do j = 2, c%Xh%lx -1
88  sx = sx + abs(abs(nx(l, j, facet, el)) - 1d0)
89  sy = sy + abs(abs(ny(l, j, facet, el)) - 1d0)
90  sz = sz + abs(abs(nz(l, j, facet, el)) - 1d0)
91  end do
92  end do
93  case(3,4)
94  do l = 2, c%Xh%lx - 1
95  do j = 2, c%Xh%lx - 1
96  sx = sx + abs(abs(nx(l, j, facet, el)) - 1d0)
97  sy = sy + abs(abs(ny(l, j, facet, el)) - 1d0)
98  sz = sz + abs(abs(nz(l, j, facet, el)) - 1d0)
99  end do
100  end do
101  case(5,6)
102  do l = 2, c%Xh%lx - 1
103  do j = 2, c%Xh%lx - 1
104  sx = sx + abs(abs(nx(l, j, facet, el)) - 1d0)
105  sy = sy + abs(abs(ny(l, j, facet, el)) - 1d0)
106  sz = sz + abs(abs(nz(l, j, facet, el)) - 1d0)
107  end do
108  end do
109  end select
110  sx = sx / (c%Xh%lx - 2)**2
111  sy = sy / (c%Xh%lx - 2)**2
112  sz = sz / (c%Xh%lx - 2)**2
113 
114  if (sx .lt. tol) then
115  call this%bc_y%mark_facet(facet, el)
116  call this%bc_z%mark_facet(facet, el)
117  end if
118 
119  if (sy .lt. tol) then
120  call this%bc_x%mark_facet(facet, el)
121  call this%bc_z%mark_facet(facet, el)
122  end if
123 
124  if (sz .lt. tol) then
125  call this%bc_y%mark_facet(facet, el)
126  call this%bc_x%mark_facet(facet, el)
127  end if
128  end do
129  end associate
130  call this%bc_x%finalize()
131  call this%bc_x%set_g(0.0_rp)
132  call this%bc_y%finalize()
133  call this%bc_y%set_g(0.0_rp)
134  call this%bc_z%finalize()
135  call this%bc_z%set_g(0.0_rp)
136  end subroutine non_normal_init_msk
137 
138 
139  subroutine non_normal_free(this)
140  type(non_normal_t), intent(inout) :: this
141 
142  call this%bc_x%free()
143  call this%bc_y%free()
144  call this%bc_z%free()
145 
146  end subroutine non_normal_free
147 end module non_normal
Coefficients.
Definition: coef.f90:34
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
Defines a dirichlet boundary condition.
Definition: dirichlet.f90:34
Definition: math.f90:60
Build configurations.
Definition: neko_config.f90:34
Dirichlet condition on axis aligned plane in the non normal direction.
Definition: non_normal.f90:34
subroutine non_normal_init_msk(this)
Initialize symmetry mask for each axis.
Definition: non_normal.f90:59
subroutine non_normal_free(this)
Definition: non_normal.f90:140
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
Implements a n-tuple.
Definition: tuple.f90:34
Utilities.
Definition: utils.f90:35
Dirichlet condition in non normal direction of a plane.
Definition: non_normal.f90:50
Mixed Dirichlet-Neumann symmetry plane condition.
Definition: symmetry.f90:50
Integer based 2-tuple.
Definition: tuple.f90:51