32 character(len=H_SHORT) :: name
33 character(len=H_MID) :: desc
34 character(len=H_SHORT) :: unit
35 logical :: positive_down
39 integer :: local_mesh_num
41 integer :: local_mesh_num_global
43 integer,
allocatable :: tileid_globalmap(:,:)
44 integer,
allocatable :: tilefaceid_globalmap(:,:)
45 integer,
allocatable :: tilepanelid_globalmap(:,:)
46 integer,
allocatable :: tileid_global2localmap(:)
47 integer,
allocatable :: prcrank_globalmap(:)
54 logical :: isgenerated
57 procedure :: setdiminfo => meshbase_setdiminfo
64 class(
meshbase),
target,
intent(in) :: this
65 integer,
intent(in) :: id
91 ndimtype, refElem, NLocalMeshPerPrc, NsideTile, &
94 use scale_prc,
only: &
98 class(
meshbase),
intent(inout) :: this
99 integer,
intent(in) :: ndimtype
101 integer,
intent(in) :: nlocalmeshperprc
102 integer,
intent(in) :: nsidetile
103 integer,
intent(in),
optional :: nprocs
108 if (
present(nprocs))
then
109 this%PRC_NUM = nprocs
111 this%PRC_NUM = prc_nprocs
114 this%LOCAL_MESH_NUM = nlocalmeshperprc
115 this%LOCAL_MESH_NUM_global = this%PRC_NUM * this%LOCAL_MESH_NUM
117 this%refElem => refelem
119 allocate( this%tileID_globalMap(nsidetile, this%LOCAL_MESH_NUM_global) )
120 allocate( this%tileFaceID_globalMap(nsidetile, this%LOCAL_MESH_NUM_global) )
121 allocate( this%tilePanelID_globalMap(nsidetile, this%LOCAL_MESH_NUM_global) )
122 allocate( this%tileID_global2localMap(this%LOCAL_MESH_NUM_global) )
123 allocate( this%PRCRank_globalMap(this%LOCAL_MESH_NUM_global) )
124 allocate( this%dimInfo(ndimtype) )
126 this%isGenerated = .false.
134 class(
meshbase),
intent(inout) :: this
137 if (
allocated(this%tileID_globalMap) )
then
138 deallocate( this%tileID_globalMap )
139 deallocate( this%tileFaceID_globalMap )
140 deallocate( this%tilePanelID_globalMap )
141 deallocate( this%tileID_global2localMap )
142 deallocate( this%PRCRank_globalMap )
143 deallocate( this%dimInfo )
154 integer,
intent(in) :: ndim
159 refelem => mesh%refElem
161 allocate( mesh%pos_en(refelem%Np,mesh%Ne,ndim) )
164 allocate( mesh%normal_fn(refelem%NfpTot,mesh%Ne,ndim) )
165 allocate( mesh%sJ(refelem%NfpTot,mesh%Ne) )
166 allocate( mesh%J(refelem%Np,mesh%Ne) )
167 allocate( mesh%Fscale(refelem%NfpTot,mesh%Ne) )
168 allocate( mesh%Escale(refelem%Np,mesh%Ne,ndim,ndim) )
169 allocate( mesh%Gsqrt(refelem%Np,mesh%Ne) )
170 allocate( mesh%G_ij(refelem%Np,mesh%Ne, ndim,ndim) )
171 allocate( mesh%GIJ (refelem%Np,mesh%Ne, ndim,ndim) )
177 subroutine meshbase_setdiminfo( this, &
178 dimID, name, unit, desc, positive_down )
180 class(
meshbase),
intent(inout) :: this
181 integer,
intent(in) :: dimid
182 character(len=*),
intent(in) :: name
183 character(len=*),
intent(in) :: unit
184 character(len=*),
intent(in) :: desc
185 logical,
intent(in),
optional :: positive_down
189 this%dimInfo(dimid)%name = name
190 this%dimInfo(dimid)%unit = unit
191 this%dimInfo(dimid)%desc = desc
192 if (
present(positive_down) )
then
193 this%dimInfo(dimid)%positive_down = positive_down
195 this%dimInfo(dimid)%positive_down = .false.
199 end subroutine meshbase_setdiminfo
201end module scale_mesh_base
module FElib / Element / Base
module FElib / Mesh / Local, Base
module FElib / Mesh / Base
subroutine, public meshbase_final(this)
subroutine, public meshbase_setgeometricinfo(mesh, ndim)
subroutine, public meshbase_init(this, ndimtype, refelem, nlocalmeshperprc, nsidetile, nprocs)