Neko  0.8.99
A portable framework for high-order spectral element flow simulations
tri_mesh.f90
Go to the documentation of this file.
1 ! Copyright (c) 2022, 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 module tri_mesh
36  use tri, only : tri_t, neko_tri_npts
37  use point, only : point_t
38  implicit none
39  private
40 
41  type, public :: tri_mesh_t
42  type(tri_t), allocatable :: el(:)
43  type(point_t), allocatable :: points(:)
44  integer :: nelv
45  integer :: mpts
46  contains
47  procedure, pass(this) :: init => tri_mesh_init
48  procedure, pass(this) :: free => tri_mesh_free
49  procedure, pass(this) :: add_element => tri_mesh_add_element
50  end type tri_mesh_t
51 
52 contains
53 
55  subroutine tri_mesh_init(this, nelv)
56  class(tri_mesh_t), intent(inout) :: this
57  integer, intent(in) :: nelv
58 
59  call this%free()
60 
61  this%nelv = nelv
62  allocate(this%el(this%nelv))
63  allocate(this%points(nelv * neko_tri_npts))
64 
65  this%mpts = 0
66  this%nelv = 0
67 
68  end subroutine tri_mesh_init
69 
71  subroutine tri_mesh_free(this)
72  class(tri_mesh_t), intent(inout) :: this
73 
74  if (allocated(this%el)) then
75  deallocate(this%el)
76  end if
77 
78  if (allocated(this%points)) then
79  deallocate(this%points)
80  end if
81 
82  end subroutine tri_mesh_free
83 
85  subroutine tri_mesh_add_element(this, p1, p2, p3)
86  class(tri_mesh_t), intent(inout) :: this
87  type(point_t), target, intent(inout) :: p1, p2, p3
88 
89  this%points(this%mpts + 1) = p1
90  this%points(this%mpts + 2) = p2
91  this%points(this%mpts + 3) = p3
92 
93  this%nelv = this%nelv + 1
94  call this%el(this%nelv)%init(this%nelv, &
95  this%points(this%mpts + 1), &
96  this%points(this%mpts + 2), &
97  this%points(this%mpts + 3))
98  this%mpts = this%mpts + 3
99 
100  end subroutine tri_mesh_add_element
101 
102 end module tri_mesh
Implements a point.
Definition: point.f90:35
Defines a triangular surface mesh.
Definition: tri_mesh.f90:35
subroutine tri_mesh_init(this, nelv)
Initialise a triangular surface mesh.
Definition: tri_mesh.f90:56
subroutine tri_mesh_add_element(this, p1, p2, p3)
Add an element to a mesh.
Definition: tri_mesh.f90:86
subroutine tri_mesh_free(this)
Deallocate a triangular surface mesh.
Definition: tri_mesh.f90:72
Defines a triangular element.
Definition: tri.f90:34
integer, parameter, public neko_tri_npts
Number of points.
Definition: tri.f90:42
A point in with coordinates .
Definition: point.f90:43
Triangular element.
Definition: tri.f90:60