FE-Project
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines | Variables
scale_mesh_base1d Module Reference

module FElib / Mesh / Base 1D More...

Data Types

type  meshbase1d
 
interface  meshbase1d_generate
 

Functions/Subroutines

subroutine, public meshbase1d_init (this, neg, dom_xmin, dom_xmax, refelem, nlocalmeshperprc, nprocs, myrank, fx)
 
subroutine, public meshbase1d_final (this)
 
subroutine, public meshbase1d_setgeometricinfo (lcmesh)
 
subroutine, public meshbase1d_assigndomid (this, tileid_table, panelid_table, pi_table)
 
subroutine, public meshbase1d_setuplocaldom (lcmesh, tileid, panelid, i, nprc, dom_xmin, dom_xmax, ne, fx)
 

Variables

integer, public meshbase1d_dimtype_num = 2
 
integer, public meshbase1d_dimtypeid_x = 1
 
integer, public meshbase1d_dimtypeid_xt = 2
 

Detailed Description

module FElib / Mesh / Base 1D

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

Function/Subroutine Documentation

◆ meshbase1d_init()

subroutine, public scale_mesh_base1d::meshbase1d_init ( class(meshbase1d), intent(inout) this,
integer, intent(in) neg,
real(rp), intent(in) dom_xmin,
real(rp), intent(in) dom_xmax,
class(elementbase1d), intent(in), target refelem,
integer, intent(in) nlocalmeshperprc,
integer, intent(in), optional nprocs,
integer, intent(in), optional myrank,
real(rp), dimension(neg+1), intent(in), optional fx )

Definition at line 79 of file scale_mesh_base1d.F90.

84
85 implicit none
86
87 class(MeshBase1D), intent(inout) :: this
88 integer, intent(in) :: NeG
89 real(RP), intent(in) :: dom_xmin
90 real(RP), intent(in) :: dom_xmax
91 class(ElementBase1D), intent(in), target :: refElem
92 integer, intent(in) :: NLocalMeshPerPrc
93 integer, intent(in), optional :: nprocs
94 integer, intent(in), optional :: myrank
95 real(RP), intent(in), optional :: FX(NeG+1)
96
97 integer :: n
98 integer :: k
99 real(RP) :: dx
100 !-----------------------------------------------------------------------------
101
102 this%NeG = neg
103
104 this%xmin_gl = dom_xmin
105 this%xmax_gl = dom_xmax
106
107 !- Fx
108 allocate( this%FX(neg+1) )
109 if ( present(fx) ) then
110 this%FX(:) = fx(:)
111 else
112 this%FX(1 ) = dom_xmin
113 this%FX(neg+1) = dom_xmax
114 dx = (dom_xmax - dom_xmin) / dble(neg)
115 do k=2, neg
116 this%FX(k) = this%FX(k-1) + dx
117 end do
118 end if
119
120 this%refElem1D => refelem
121 call meshbase_init( this, &
122 meshbase1d_dimtype_num, refelem, &
123 nlocalmeshperprc, 2, &
124 nprocs )
125
126 this%Nprc = this%PRC_NUM
127
128 allocate( this%lcmesh_list(this%LOCAL_MESH_NUM) )
129 do n=1, this%LOCAL_MESH_NUM
130 call localmesh1d_init( this%lcmesh_list(n), n, refelem, myrank )
131 end do
132
133 call this%SetDimInfo( meshbase1d_dimtypeid_x, "x", "m", "X-coordinate" )
134 call this%SetDimInfo( meshbase1d_dimtypeid_xt, "xt", "m", "X-coordinate" )
135
136 return

References scale_localmesh_1d::localmesh1d_init(), meshbase1d_dimtype_num, meshbase1d_dimtypeid_x, meshbase1d_dimtypeid_xt, and scale_mesh_base::meshbase_init().

Referenced by scale_mesh_base1d::meshbase1d_generate::meshbase1d_generate(), and scale_mesh_linedom1d::meshlinedom1d_init().

◆ meshbase1d_final()

subroutine, public scale_mesh_base1d::meshbase1d_final ( class(meshbase1d), intent(inout) this)

Definition at line 140 of file scale_mesh_base1d.F90.

141 implicit none
142 class(MeshBase1D), intent(inout) :: this
143
144 integer :: n
145 !-----------------------------------------------------------------------------
146
147 if ( allocated ( this%lcmesh_list ) ) then
148 do n=1, this%LOCAL_MESH_NUM
149 call localmesh1d_final( this%lcmesh_list(n), this%isGenerated )
150 end do
151
152 deallocate( this%lcmesh_list )
153 end if
154
155 call meshbase_final(this)
156
157 return

References scale_localmesh_1d::localmesh1d_final(), and scale_mesh_base::meshbase_final().

Referenced by scale_mesh_base1d::meshbase1d_generate::meshbase1d_generate(), and scale_mesh_linedom1d::meshlinedom1d_final().

◆ meshbase1d_setgeometricinfo()

subroutine, public scale_mesh_base1d::meshbase1d_setgeometricinfo ( type(localmesh1d), intent(inout) lcmesh)

Definition at line 175 of file scale_mesh_base1d.F90.

176 implicit none
177
178 type(LocalMesh1D), intent(inout) :: lcmesh
179
180 class(ElementBase1D), pointer :: refElem
181 integer :: ke
182 integer :: f
183 integer :: i
184
185 integer :: node_ids(lcmesh%refElem%Nv)
186 real(RP) :: vx(lcmesh%refElem%Nv)
187 real(RP) :: xr(lcmesh%refElem%Np)
188 real(DP) :: Escale(1,1,lcmesh%refElem%Np)
189 integer :: fmask(lcmesh%refElem%NfpTot)
190 integer :: fid(lcmesh%refElem1D%Nfp,lcmesh%refElem1D%Nfaces)
191
192 !-----------------------------------------------------------------------------
193
194 call meshbase_setgeometricinfo(lcmesh, 1)
195 refelem => lcmesh%refElem1D
196
197 fmask(:) = reshape(refelem%Fmask, shape(fmask))
198 do f=1, refelem%Nfaces
199 do i=1, refelem%Nfp
200 fid(i,f) = i + (f-1)*refelem%Nfp
201 end do
202 end do
203
204 do ke=1, lcmesh%Ne
205 node_ids(:) = lcmesh%EToV(ke,:)
206 vx(:) = lcmesh%pos_ev(node_ids(:),1)
207 lcmesh%pos_en(:,ke,1) = vx(1) + 0.5_rp*(refelem%x1(:) + 1.0_rp)*(vx(2) - vx(1))
208
209 xr(:) = 0.5_rp*(vx(2) - vx(1)) !matmul(refElem%Dr,mesh%x(:,n))
210
211 lcmesh%J(:,ke) = xr
212 lcmesh%Escale(:,ke,1,1) = 1.0_rp/lcmesh%J(:,ke)
213
214 !* Face
215
216 !
217 !mesh%fx(:,n) = mesh%x(fmask(:),n)
218 !mesh%fy(:,n) = mesh%y(fmask(:),n)
219
220 ! Calculate normal vectors
221 lcmesh%normal_fn(fid(:,1),ke,1) = - 1.0_rp
222 lcmesh%normal_fn(fid(:,2),ke,1) = + 1.0_rp
223 lcmesh%sJ(:,ke) = 1.0_rp
224 lcmesh%Fscale(:,ke) = lcmesh%sJ(:,ke)/lcmesh%J(fmask(:),ke)
225 lcmesh%Gsqrt(:,ke) = 1.0_rp
226 end do
227
228 return

References scale_mesh_base::meshbase_setgeometricinfo().

Referenced by scale_mesh_base1d::meshbase1d_generate::meshbase1d_generate(), and meshbase1d_setuplocaldom().

◆ meshbase1d_assigndomid()

subroutine, public scale_mesh_base1d::meshbase1d_assigndomid ( class(meshbase1d), intent(inout) this,
integer, dimension(this%local_mesh_num, this%prc_num), intent(out) tileid_table,
integer, dimension(this%local_mesh_num*this%prc_num), intent(out) panelid_table,
integer, dimension(this%local_mesh_num*this%prc_num), intent(out) pi_table )

Definition at line 232 of file scale_mesh_base1d.F90.

235
236 use scale_meshutil_1d, only: &
238 implicit none
239
240 class(MeshBase1D), intent(inout) :: this
241 integer, intent(out) :: tileID_table(this%LOCAL_MESH_NUM, this%PRC_NUM)
242 integer, intent(out) :: panelID_table(this%LOCAL_MESH_NUM*this%PRC_NUM)
243 integer, intent(out) :: pi_table(this%LOCAL_MESH_NUM*this%PRC_NUM)
244
245 integer :: n
246 integer :: p
247 integer :: tileID
248 !-----------------------------------------------------------------------------
249
251 panelid_table, pi_table, & ! (out)
252 this%tileID_globalMap, this%tileFaceID_globalMap, this%tilePanelID_globalMap, & ! (out)
253 this%LOCAL_MESH_NUM_global ) ! (in)
254
255 do p=1, this%PRC_NUM
256 do n=1, this%LOCAL_MESH_NUM
257 tileid = n + (p-1)*this%LOCAL_MESH_NUM
258
259 tileid_table(n,p) = tileid
260 this%tileID_global2localMap(tileid) = n
261 this%PRCRank_globalMap(tileid) = p - 1
262 end do
263 end do
264
265 return
module FElib / Mesh / utility for 1D mesh
subroutine, public meshutil1d_buildglobalmap(panelid_table, pi_table, tileid_map, tilefaceid_map, tilepanelid_map, ntile)

References scale_meshutil_1d::meshutil1d_buildglobalmap().

Referenced by scale_mesh_base1d::meshbase1d_generate::meshbase1d_generate(), and scale_mesh_linedom1d::meshlinedom1d_generate().

◆ meshbase1d_setuplocaldom()

subroutine, public scale_mesh_base1d::meshbase1d_setuplocaldom ( type(localmesh1d), intent(inout) lcmesh,
integer, intent(in) tileid,
integer, intent(in) panelid,
integer, intent(in) i,
integer, intent(in) nprc,
real(rp) dom_xmin,
real(rp) dom_xmax,
integer, intent(in) ne,
real(rp), dimension(ne*nprc+1), intent(in) fx )

Definition at line 269 of file scale_mesh_base1d.F90.

274
275 use scale_meshutil_1d, only: &
280
282 implicit none
283
284 type(LocalMesh1D), intent(inout) :: lcmesh
285 integer, intent(in) :: tileID
286 integer, intent(in) :: panelID
287 integer, intent(in) :: i
288 integer, intent(in) :: Nprc
289 real(RP) :: dom_xmin, dom_xmax
290 integer, intent(in) ::Ne
291 real(RP), intent(in) :: FX(Ne*Nprc+1)
292
293 class(ElementBase1D), pointer :: elem
294 real(RP) :: delx
295 real(RP) :: FX_lc(Ne+1)
296 !-----------------------------------------------------------------------------
297
298 elem => lcmesh%refElem1D
299 lcmesh%tileID = tileid
300 lcmesh%panelID = panelid
301
302 !--
303
304 lcmesh%Ne = ne
305 lcmesh%Nv = ne + 1
306 lcmesh%NeS = 1
307 lcmesh%NeE = lcmesh%Ne
308 lcmesh%NeA = lcmesh%Ne + 2
309
310 !delx = (dom_xmax - dom_xmin)/dble(Nprc)
311 fx_lc(:) = fx((i-1)*ne+1:i*ne)
312 lcmesh%xmin = fx_lc(1)
313 lcmesh%xmax = fx_lc(ne+1)
314
315 allocate(lcmesh%pos_ev(lcmesh%Nv,1))
316 allocate( lcmesh%EToV(lcmesh%Ne,2) )
317 allocate( lcmesh%EToE(lcmesh%Ne,elem%Nfaces) )
318 allocate( lcmesh%EToF(lcmesh%Ne,elem%Nfaces) )
319 allocate( lcmesh%BCType(lcmesh%refElem%Nfaces,lcmesh%Ne) )
320 allocate( lcmesh%VMapM(elem%NfpTot, lcmesh%Ne) )
321 allocate( lcmesh%VMapP(elem%NfpTot, lcmesh%Ne) )
322 allocate( lcmesh%MapM(elem%NfpTot, lcmesh%Ne) )
323 allocate( lcmesh%MapP(elem%NfpTot, lcmesh%Ne) )
324
325 lcmesh%BCType(:,:) = bctype_interior
326
327 !----
328
329 call meshutil1d_genlinedomain( lcmesh%pos_ev, lcmesh%EToV, & ! (out)
330 lcmesh%Ne, lcmesh%xmin, lcmesh%xmax, fx=fx_lc ) ! (in)
331
332
333 !---
334 call meshbase1d_setgeometricinfo( lcmesh )
335
336 !---
337
338 call meshutil1d_genconnectivity( lcmesh%EToE, lcmesh%EToF, & ! (out)
339 & lcmesh%EToV, lcmesh%Ne, elem%Nfaces ) ! (in)
340
341 !---
342 call meshutil1d_buildinteriormap( lcmesh%VmapM, lcmesh%VMapP, lcmesh%MapM, lcmesh%MapP, &
343 & lcmesh%pos_en, lcmesh%pos_ev, lcmesh%EToE, lcmesh%EtoF, lcmesh%EtoV, &
344 & elem%Fmask, lcmesh%Ne, elem%Np, elem%Nfp, elem%Nfaces, lcmesh%Nv )
345
346 call meshutil1d_genpatchboundarymap( lcmesh%VMapB, lcmesh%MapB, lcmesh%VMapP, &
347 & lcmesh%pos_en, lcmesh%xmin, lcmesh%xmax, &
348 & elem%Fmask, lcmesh%Ne, elem%Np, elem%Nfp, elem%Nfaces, lcmesh%Nv)
349
350 return
module FElib / Mesh / Local, Base
integer, parameter, public bctype_interior
subroutine, public meshutil1d_buildinteriormap(vmapm, vmapp, mapm, mapp, pos_en, pos_ev, etoe, etof, etov, fmask, ne, np, nfp, nfaces, nv)
subroutine, public meshutil1d_genconnectivity(etoe, etof, etov, ne, nfaces)
subroutine, public meshutil1d_genpatchboundarymap(vmapb, mapb, vmapp, pos_en, xmin, xmax, fmask, ne, np, nfp, nfaces, nv)
subroutine, public meshutil1d_genlinedomain(pos_v, etov, ke_x, xmin, xmax, fx)

References scale_localmesh_base::bctype_interior, meshbase1d_setgeometricinfo(), scale_meshutil_1d::meshutil1d_buildinteriormap(), scale_meshutil_1d::meshutil1d_genconnectivity(), scale_meshutil_1d::meshutil1d_genlinedomain(), and scale_meshutil_1d::meshutil1d_genpatchboundarymap().

Referenced by scale_mesh_base1d::meshbase1d_generate::meshbase1d_generate(), and scale_mesh_linedom1d::meshlinedom1d_generate().

Variable Documentation

◆ meshbase1d_dimtype_num

integer, public scale_mesh_base1d::meshbase1d_dimtype_num = 2

Definition at line 63 of file scale_mesh_base1d.F90.

63 integer, public :: MESHBASE1D_DIMTYPE_NUM = 2

Referenced by scale_file_history_meshfield::file_history_meshfield_in::file_history_meshfield_in3d(), and meshbase1d_init().

◆ meshbase1d_dimtypeid_x

integer, public scale_mesh_base1d::meshbase1d_dimtypeid_x = 1

◆ meshbase1d_dimtypeid_xt

integer, public scale_mesh_base1d::meshbase1d_dimtypeid_xt = 2

Definition at line 65 of file scale_mesh_base1d.F90.

65 integer, public :: MESHBASE1D_DIMTYPEID_XT = 2

Referenced by scale_file_common_meshfield::file_common_meshfield_get_dims::file_common_meshfield_get_dims1d(), and meshbase1d_init().