FE-Project
Loading...
Searching...
No Matches
scale_localmesh_base.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2!> module FElib / Mesh / Local, Base
3!!
4!! @par Description
5!! Base module to manage meshes for element-based methods
6!!
7!! @author Yuta Kawai, Team SCALE
8!<
9#include "scaleFElib.h"
11
12 !-----------------------------------------------------------------------------
13 !
14 !++ used modules
15 !
16 use scale_precision
17 use scale_io
18
20
21 !-----------------------------------------------------------------------------
22 implicit none
23 private
24
25 !-----------------------------------------------------------------------------
26 !
27 !++ Public type & procedure
28 !
29
30 !> Derived type to manage a local computational domain (base type)
31 type, public :: localmeshbase
32 integer :: ne !< Number of finite elements in the local computational domain
33 integer :: nes !< Start index of finite element indices in the local computational domain
34 integer :: nee !< End index of finite element indices in the local computational domain
35 integer :: nea !< Total number of finite elements (Ne + halo buffer)
36 integer :: nv !< Number of vertices in the local computational domain
37
38 class(elementbase), pointer :: refelem !< Pointer to an object with a reference element
39
40 real(rp), allocatable :: pos_ev(:,:) !< Position of vertices in finite elements in the local computational domain
41 real(rp), allocatable :: pos_en(:,:,:) !< Position of nodes in finite elements in the local computational domain
42 real(rp), allocatable :: normal_fn(:,:,:)
43
44 real(rp), allocatable :: sj(:,:)
45 real(rp), allocatable :: j(:,:)
46
47 real(rp), allocatable :: escale(:,:,:,:)
48 real(rp), allocatable :: fscale(:,:)
49
50 integer, allocatable :: etov(:,:)
51 integer, allocatable :: etoe(:,:)
52 integer, allocatable :: etof(:,:)
53 integer, allocatable :: vmapm(:,:)
54 integer, allocatable :: vmapp(:,:)
55 integer, allocatable :: mapm(:,:)
56 integer, allocatable :: mapp(:,:)
57
58 integer, allocatable :: bctype(:,:)
59 integer, allocatable :: mapb(:)
60 integer, allocatable :: vmapb(:)
61 ! integer, allocatable :: MapD(:)
62 ! integer, allocatable :: VMapD(:)
63 ! integer, allocatable :: MapN(:)
64 ! integer, allocatable :: VMapN(:)
65
66 integer :: tileid
67 integer :: panelid
68 integer :: prc_myrank
69 integer :: lcdomid !< ID of local computational mesh
70
71 real(rp), allocatable :: g_ij(:,:,:,:) !< The covariant component of metric tensor with horizontal general curvilinear coordinate
72 real(rp), allocatable :: gij(:,:,:,:) !< The contravariant component of metric tensor with horizontal general curvilinear coordinate
73 real(rp), allocatable :: gsqrt(:,:) !< The Jacobian of 3D transformation in the computational coordinate (=GsqrtH * GsqrtV)
74 end type localmeshbase
75
76 public :: localmeshbase_init
77 public :: localmeshbase_final
78
79 !-----------------------------------------------------------------------------
80 !
81 !++ Public parameters & variables
82 !
83 integer, public, parameter :: bctype_interior = 0
84 integer, public, parameter :: bctype_undefbc = 1
85 integer, public, parameter :: bctype_periodic = 2
86 integer, public, parameter :: bctype_dirchlet = 3
87 integer, public, parameter :: bctype_neuman = 4
88 integer, public, parameter :: bctype_shoreline = 5
89
90 !-----------------------------------------------------------------------------
91 !
92 !++ Private procedure
93 !
94
95 !-----------------------------------------------------------------------------
96 !
97 !++ Private parameters & variables
98 !
99
100contains
101!> Setup an object to manage a local computational mesh
102 subroutine localmeshbase_init( this, lcdomID, refElem, ndim, myrank )
103
104 use scale_prc, only: prc_myrank
105 implicit none
106
107 class(localmeshbase), intent(inout) :: this
108 integer, intent(in) :: lcdomid
109 class(elementbase), intent(in), target :: refelem
110 integer, intent(in) :: ndim
111 integer, intent(in), optional :: myrank
112 !-----------------------------------------------------------------------------
113
114 this%lcdomID = lcdomid
115 this%refElem => refelem
116
117 if (present(myrank)) then
118 this%PRC_myrank = myrank
119 else
120 this%PRC_myrank = prc_myrank
121 end if
122
123 return
124 end subroutine localmeshbase_init
125
126!> Finalize an object to manage a local computational mesh
127 subroutine localmeshbase_final( this, is_generated )
128 implicit none
129
130 class(localmeshbase), intent(inout) :: this
131 logical, intent(in) :: is_generated
132 !-----------------------------------------------------------------------------
133
134 if ( is_generated ) then
135 deallocate( this%pos_ev, this%pos_en, this%normal_fn )
136 deallocate( this%sJ, this%J )
137
138 deallocate( this%Escale, this%Fscale )
139
140 deallocate( this%EToV, this%EToE, this%EToF )
141
142 if ( allocated(this%VMapM) ) then
143 deallocate( this%VMapM, this%VMapP, this%MapM, this%MapP )
144 end if
145 if ( allocated(this%VMapB) ) then
146 deallocate( this%BCType )
147 deallocate( this%VMapB, this%MapB )
148 end if
149 if ( allocated(this%G_ij) ) then
150 deallocate( this%G_ij, this%GIJ )
151 deallocate( this%Gsqrt )
152 end if
153 end if
154
155 return
156 end subroutine localmeshbase_final
157end module scale_localmesh_base
module FElib / Element / Base
module FElib / Mesh / Local, Base
integer, parameter, public bctype_interior
integer, parameter, public bctype_shoreline
integer, parameter, public bctype_neuman
subroutine, public localmeshbase_final(this, is_generated)
Finalize an object to manage a local computational mesh.
integer, parameter, public bctype_periodic
integer, parameter, public bctype_dirchlet
subroutine, public localmeshbase_init(this, lcdomid, refelem, ndim, myrank)
Setup an object to manage a local computational mesh.
integer, parameter, public bctype_undefbc
Derived type representing an arbitrary finite element.
Derived type to manage a local computational domain (base type)