FE-Project
Loading...
Searching...
No Matches
scale_mesh_base Module Reference

module FElib / Mesh / Base More...

Data Types

type  meshbase
interface  meshbase_get_localmesh
type  meshdiminfo

Functions/Subroutines

subroutine, public meshbase_init (this, ndimtype, refelem, nlocalmeshperprc, nsidetile, nprocs)
subroutine, public meshbase_final (this)
subroutine, public meshbase_setgeometricinfo (mesh, ndim)

Detailed Description

module FElib / Mesh / Base

Description
Base module to manage meshes for element-based methods
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ meshbase_init()

subroutine, public scale_mesh_base::meshbase_init ( class(meshbase), intent(inout) this,
integer, intent(in) ndimtype,
class(elementbase), intent(in), target refelem,
integer, intent(in) nlocalmeshperprc,
integer, intent(in) nsidetile,
integer, intent(in), optional nprocs )

Definition at line 91 of file scale_mesh_base.F90.

94
95 use scale_prc, only: &
96 prc_nprocs
97 implicit none
98
99 class(MeshBase), intent(inout) :: this
100 integer, intent(in) :: ndimtype
101 class(ElementBase), intent(in), target :: refElem
102 integer, intent(in) :: NLocalMeshPerPrc
103 integer, intent(in) :: NsideTile
104 integer, intent(in), optional :: nprocs
105
106 integer :: n
107 !-----------------------------------------------------------------------------
108
109 if ( present(nprocs) ) then
110 this%PRC_NUM = nprocs
111 else
112 this%PRC_NUM = prc_nprocs
113 end if
114
115 this%LOCAL_MESH_NUM = nlocalmeshperprc
116 this%LOCAL_MESH_NUM_global = this%PRC_NUM * this%LOCAL_MESH_NUM
117
118 this%refElem => refelem
119
120 allocate( this%tileID_globalMap(nsidetile, this%LOCAL_MESH_NUM_global) )
121 allocate( this%tileFaceID_globalMap(nsidetile, this%LOCAL_MESH_NUM_global) )
122 allocate( this%tilePanelID_globalMap(nsidetile, this%LOCAL_MESH_NUM_global) )
123 allocate( this%tileID_global2localMap(this%LOCAL_MESH_NUM_global) )
124 allocate( this%PRCRank_globalMap(this%LOCAL_MESH_NUM_global) )
125 allocate( this%dimInfo(ndimtype) )
126
127 this%isGenerated = .false.
128
129 return

Referenced by scale_mesh_base1d::meshbase1d_init(), scale_mesh_base2d::meshbase2d_init(), scale_mesh_base3d::meshbase3d_init(), and scale_mesh_base::meshbase_get_localmesh::meshbase_get_localmesh().

◆ meshbase_final()

subroutine, public scale_mesh_base::meshbase_final ( class(meshbase), intent(inout) this)

Definition at line 133 of file scale_mesh_base.F90.

134 implicit none
135 class(MeshBase), intent(inout) :: this
136 !-----------------------------------------------------------------------------
137
138 if ( allocated(this%tileID_globalMap) ) then
139 deallocate( this%tileID_globalMap )
140 deallocate( this%tileFaceID_globalMap )
141 deallocate( this%tilePanelID_globalMap )
142 deallocate( this%tileID_global2localMap )
143 deallocate( this%PRCRank_globalMap )
144 deallocate( this%dimInfo )
145 end if
146
147 return

Referenced by scale_mesh_base1d::meshbase1d_final(), scale_mesh_base2d::meshbase2d_final(), scale_mesh_base3d::meshbase3d_final(), and scale_mesh_base::meshbase_get_localmesh::meshbase_get_localmesh().

◆ meshbase_setgeometricinfo()

subroutine, public scale_mesh_base::meshbase_setgeometricinfo ( class(localmeshbase), intent(inout) mesh,
integer, intent(in) ndim )

Definition at line 151 of file scale_mesh_base.F90.

152 implicit none
153
154 class(LocalMeshBase), intent(inout) :: mesh
155 integer, intent(in) :: ndim
156
157 class(ElementBase), pointer :: refElem
158 !-----------------------------------------------------------------------------
159
160 refelem => mesh%refElem
161
162 allocate( mesh%pos_en(refelem%Np,mesh%Ne,ndim) )
163 !allocate( mesh%fx(refElem%Nfaces*refElem%Nfp,mesh%Ne) )
164 !allocate( mesh%fy(refElem%Nfaces*refElem%Nfp,mesh%Ne) )
165 allocate( mesh%normal_fn(refelem%NfpTot,mesh%Ne,ndim) )
166 allocate( mesh%sJ(refelem%NfpTot,mesh%Ne) )
167 allocate( mesh%J(refelem%Np,mesh%Ne) )
168 allocate( mesh%Fscale(refelem%NfpTot,mesh%Ne) )
169 allocate( mesh%Escale(refelem%Np,mesh%Ne,ndim,ndim) )
170 allocate( mesh%Gsqrt(refelem%Np,mesh%Ne) )
171 allocate( mesh%G_ij(refelem%Np,mesh%Ne, ndim,ndim) )
172 allocate( mesh%GIJ (refelem%Np,mesh%Ne, ndim,ndim) )
173
174 return

Referenced by scale_mesh_base1d::meshbase1d_setgeometricinfo(), and scale_mesh_base::meshbase_get_localmesh::meshbase_get_localmesh().