Neko  0.8.99
A portable framework for high-order spectral element flow simulations
facet_zone.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 !
34 module facet_zone
35  use tuple, only : tuple_i4_t, tuple4_i4_t
36  use stack, only : stack_i4t2_t, stack_i4t4_t
37  use utils, only : neko_error
38  implicit none
39  private
40 
41  type, public :: facet_zone_t
42  type(tuple_i4_t), allocatable :: facet_el(:)
43  integer :: size
44  logical, private :: finalized = .false.
45  type(stack_i4t2_t), private :: scratch
46  contains
47  procedure, pass(z) :: init => facet_zone_init
48  procedure, pass(z) :: free => facet_zone_free
49  procedure, pass(z) :: finalize => facet_zone_finalize
50  procedure, pass(z) :: add_facet => facet_zone_add_facet
51  end type facet_zone_t
52 
53  type, public, extends(facet_zone_t) :: facet_zone_periodic_t
54  type(tuple_i4_t), allocatable :: p_facet_el(:)
55  type(stack_i4t2_t), private :: p_scratch
56  type(tuple4_i4_t), allocatable :: p_ids(:)
57  type(stack_i4t4_t), private :: p_id_scratch
58  type(tuple4_i4_t), allocatable :: org_ids(:)
59  type(stack_i4t4_t), private :: org_id_scratch
60  contains
61  procedure, pass(z) :: init => facet_zone_periodic_init
62  procedure, pass(z) :: free => facet_zone_periodic_free
63  procedure, pass(z) :: finalize => facet_zone_periodic_finalize
64  procedure, pass(z) :: add_periodic_facet => facet_zone_periodic_add_facet
65  end type facet_zone_periodic_t
66 
67 contains
68 
70  subroutine facet_zone_init(z, size)
71  class(facet_zone_t), intent(inout) :: z
72  integer, optional :: size
73 
74  call facet_zone_free(z)
75 
76  if (present(size)) then
77  call z%scratch%init(size)
78  else
79  call z%scratch%init()
80  end if
81 
82  end subroutine facet_zone_init
83 
85  subroutine facet_zone_free(z)
86  class(facet_zone_t), intent(inout) :: z
87  if (allocated(z%facet_el)) then
88  deallocate(z%facet_el)
89  end if
90 
91  z%finalized = .false.
92  z%size = 0
93 
94  call z%scratch%free()
95 
96  end subroutine facet_zone_free
97 
100  subroutine facet_zone_finalize(z)
101  class(facet_zone_t), intent(inout) :: z
102  type(tuple_i4_t), pointer :: tp(:)
103  integer :: i
104 
105  if (.not. z%finalized) then
106 
107  allocate(z%facet_el(z%scratch%size()))
108 
109  tp => z%scratch%array()
110  do i = 1, z%scratch%size()
111  z%facet_el(i) = tp(i)
112  end do
113 
114  z%size = z%scratch%size()
115 
116  call z%scratch%clear()
117 
118  z%finalized = .true.
119 
120  end if
121 
122  end subroutine facet_zone_finalize
123 
125  subroutine facet_zone_add_facet(z, facet, el)
126  class(facet_zone_t), intent(inout) :: z
127  integer, intent(in) :: facet
128  integer, intent(in) :: el
129  type(tuple_i4_t) :: t
130 
131  if (z%finalized) then
132  call neko_error('Facet zone already finalized')
133  end if
134 
135  t%x = (/ facet, el /)
136  call z%scratch%push(t)
137 
138  end subroutine facet_zone_add_facet
139 
141  subroutine facet_zone_periodic_init(z, size)
142  class(facet_zone_periodic_t), intent(inout) :: z
143  integer, optional :: size
144 
145  call z%free()
146 
147  if (present(size)) then
148  call facet_zone_init(z, size)
149  call z%p_scratch%init(size)
150  call z%p_id_scratch%init(size)
151  call z%org_id_scratch%init(size)
152  else
153  call facet_zone_init(z)
154  call z%p_scratch%init()
155  call z%p_id_scratch%init()
156  call z%org_id_scratch%init()
157  end if
158 
159  end subroutine facet_zone_periodic_init
160 
163  class(facet_zone_periodic_t), intent(inout) :: z
164 
165  call facet_zone_free(z)
166 
167  if (allocated(z%p_facet_el)) then
168  deallocate(z%p_facet_el)
169  end if
170 
171  if (allocated(z%p_ids)) then
172  deallocate(z%p_ids)
173  end if
174  if (allocated(z%org_ids)) then
175  deallocate(z%org_ids)
176  end if
177 
178  call z%p_scratch%free()
179  call z%p_id_scratch%free()
180  call z%org_id_scratch%free()
181 
182  end subroutine facet_zone_periodic_free
183 
187  class(facet_zone_periodic_t), intent(inout) :: z
188  type(tuple_i4_t), pointer :: tp(:)
189  type(tuple4_i4_t), pointer :: tp2(:)
190  type(tuple4_i4_t), pointer :: tp3(:)
191  integer :: i
192 
193  if (.not. z%finalized) then
194 
195  call facet_zone_finalize(z)
196 
197  if (z%size .ne. z%p_scratch%size()) then
198  call neko_error('Zone size mismatch')
199  end if
200 
201  allocate(z%p_facet_el(z%size))
202  allocate(z%p_ids(z%size))
203  allocate(z%org_ids(z%size))
204 
205  tp => z%p_scratch%array()
206  do i = 1, z%size
207  z%p_facet_el(i) = tp(i)
208  end do
209  tp2 => z%p_id_scratch%array()
210  do i = 1, z%size
211  z%p_ids(i) = tp2(i)
212  end do
213  tp3 => z%org_id_scratch%array()
214  do i = 1, z%size
215  z%org_ids(i) = tp3(i)
216  end do
217 
218  call z%p_scratch%clear()
219  call z%p_id_scratch%clear()
220  call z%org_id_scratch%clear()
221 
222  end if
223 
224  end subroutine facet_zone_periodic_finalize
225 
227  subroutine facet_zone_periodic_add_facet(z, facet, el, p_facet, p_el, pids, org_ids)
228  class(facet_zone_periodic_t), intent(inout) :: z
229  integer, intent(in) :: facet
230  integer, intent(in) :: el
231  integer, intent(in) :: p_facet
232  integer, intent(in) :: p_el
233  integer, intent(in) :: pids(4)
234  integer, intent(in) :: org_ids(4)
235  type(tuple_i4_t) :: t
236  type(tuple4_i4_t) :: t2
237  type(tuple4_i4_t) :: t3
238 
239  if (z%finalized) then
240  call neko_error('Facet zone already finalized')
241  end if
242 
243  call z%add_facet(facet, el)
244 
245  t%x = (/ p_facet, p_el /)
246  call z%p_scratch%push(t)
247  t2%x = pids
248  call z%p_id_scratch%push(t2)
249  t3%x = org_ids
250  call z%org_id_scratch%push(t3)
251 
252  end subroutine facet_zone_periodic_add_facet
253 
254 end module facet_zone
Defines a zone as a subset of facets in a mesh.
Definition: facet_zone.f90:34
subroutine facet_zone_init(z, size)
Initialize a facet zone.
Definition: facet_zone.f90:71
subroutine facet_zone_periodic_finalize(z)
Finalize a periodic zone list.
Definition: facet_zone.f90:187
subroutine facet_zone_periodic_init(z, size)
Initialize a periodic zone.
Definition: facet_zone.f90:142
subroutine facet_zone_add_facet(z, facet, el)
Add a (facet, el) tuple to an unfinalized zone.
Definition: facet_zone.f90:126
subroutine facet_zone_periodic_add_facet(z, facet, el, p_facet, p_el, pids, org_ids)
Add a (facet, el) tuple to an unfinalized zone.
Definition: facet_zone.f90:228
subroutine facet_zone_free(z)
Deallocate a facet zone.
Definition: facet_zone.f90:86
subroutine facet_zone_periodic_free(z)
Deallocate a zone.
Definition: facet_zone.f90:163
subroutine facet_zone_finalize(z)
Finalize a zone list.
Definition: facet_zone.f90:101
Implements a dynamic stack ADT.
Definition: stack.f90:35
Implements a n-tuple.
Definition: tuple.f90:34
Utilities.
Definition: utils.f90:35
Integer 2-tuple based stack.
Definition: stack.f90:84
Integer 4-tuple based stack.
Definition: stack.f90:91
Integer based 4-tuple.
Definition: tuple.f90:69
Integer based 2-tuple.
Definition: tuple.f90:51