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

module FElib / Mesh / Rectangle 2D domain More...

Data Types

type  meshrectdom2d
 

Functions/Subroutines

subroutine meshrectdom2d_init (this, negx, negy, dom_xmin, dom_xmax, dom_ymin, dom_ymax, isperiodicx, isperiodicy, refelem, nlocalmeshperprc, nprcx, nprcy, nproc, myrank)
 
subroutine, public meshrectdom2d_setuplocaldom (lcmesh, tileid, panelid, i, j, nprcx, nprcy, dom_xmin, dom_xmax, dom_ymin, dom_ymax, nex, ney)
 
subroutine, public meshrectdom2d_coord_conv (x, y, xr, xs, yr, ys, vx, vy, elem)
 

Detailed Description

module FElib / Mesh / Rectangle 2D domain

Description
Manage mesh data of rectangle 2D domain for element-based methods
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ meshrectdom2d_init()

subroutine scale_mesh_rectdom2d::meshrectdom2d_init ( class(meshrectdom2d), intent(inout) this,
integer, intent(in) negx,
integer, intent(in) negy,
real(rp), intent(in) dom_xmin,
real(rp), intent(in) dom_xmax,
real(rp), intent(in) dom_ymin,
real(rp), intent(in) dom_ymax,
logical, intent(in) isperiodicx,
logical, intent(in) isperiodicy,
type(quadrilateralelement), intent(in), target refelem,
integer, intent(in) nlocalmeshperprc,
integer, intent(in) nprcx,
integer, intent(in) nprcy,
integer, intent(in), optional nproc,
integer, intent(in), optional myrank )

Definition at line 75 of file scale_mesh_rectdom2d.F90.

82
83 implicit none
84
85 class(MeshRectDom2D), intent(inout) :: this
86 integer, intent(in) :: NeGX
87 integer, intent(in) :: NeGY
88 real(RP), intent(in) :: dom_xmin
89 real(RP), intent(in) :: dom_xmax
90 real(RP), intent(in) :: dom_ymin
91 real(RP), intent(in) :: dom_ymax
92 logical, intent(in) :: isPeriodicX
93 logical, intent(in) :: isPeriodicY
94 type(QuadrilateralElement), intent(in), target :: refElem
95 integer, intent(in) :: NLocalMeshPerPrc
96 integer, intent(in) :: NprcX
97 integer, intent(in) :: NprcY
98 integer, intent(in), optional :: nproc
99 integer, intent(in), optional :: myrank
100
101 !-----------------------------------------------------------------------------
102
103 this%NeGX = negx
104 this%NeGY = negy
105
106 this%xmin_gl = dom_xmin
107 this%xmax_gl = dom_xmax
108 this%ymin_gl = dom_ymin
109 this%ymax_gl = dom_ymax
110 this%dom_vol = (this%xmax_gl - this%xmin_gl) * (this%ymax_gl - this%ymin_gl)
111
112 this%isPeriodicX = isperiodicx
113 this%isPeriodicY = isperiodicy
114
115 this%NprcX = nprcx
116 this%NprcY = nprcy
117
118 call meshbase2d_init( this, refelem, nlocalmeshperprc, &
119 nproc, myrank )
120
121 return

References scale_mesh_base2d::meshbase2d_final(), scale_mesh_base2d::meshbase2d_init(), and meshrectdom2d_setuplocaldom().

◆ meshrectdom2d_setuplocaldom()

subroutine, public scale_mesh_rectdom2d::meshrectdom2d_setuplocaldom ( type(localmesh2d), intent(inout) lcmesh,
integer, intent(in) tileid,
integer, intent(in) panelid,
integer, intent(in) i,
integer, intent(in) j,
integer, intent(in) nprcx,
integer, intent(in) nprcy,
real(rp), intent(in) dom_xmin,
real(rp), intent(in) dom_xmax,
real(rp), intent(in) dom_ymin,
real(rp), intent(in) dom_ymax,
integer, intent(in) nex,
integer, intent(in) ney )

Definition at line 205 of file scale_mesh_rectdom2d.F90.

210
211 use scale_meshutil_2d, only: &
216
218 implicit none
219
220 type(LocalMesh2D), intent(inout) :: lcmesh
221 integer, intent(in) :: tileID
222 integer, intent(in) :: panelID
223 integer, intent(in) :: i, j
224 integer, intent(in) :: NprcX, NprcY
225 real(RP), intent(in) :: dom_xmin, dom_xmax
226 real(RP), intent(in) :: dom_ymin, dom_ymax
227 integer, intent(in) ::NeX, NeY
228
229 class(ElementBase2D), pointer :: elem
230 real(RP) :: delx, dely
231
232 !-----------------------------------------------------------------------------
233
234 elem => lcmesh%refElem2D
235 lcmesh%tileID = tileid
236 lcmesh%panelID = panelid
237
238 !--
239
240 lcmesh%Ne = nex * ney
241 lcmesh%Nv = (nex + 1)*(ney + 1)
242 lcmesh%NeS = 1
243 lcmesh%NeE = lcmesh%Ne
244 lcmesh%NeA = lcmesh%Ne + 2*(nex + ney)
245
246 lcmesh%NeX = nex
247 lcmesh%NeY = ney
248
249 delx = (dom_xmax - dom_xmin)/dble(nprcx)
250 dely = (dom_ymax - dom_ymin)/dble(nprcy)
251
252 lcmesh%xmin = dom_xmin + (i-1)*delx
253 lcmesh%xmax = dom_xmin + i *delx
254 lcmesh%ymin = dom_ymin + (j-1)*dely
255 lcmesh%ymax = dom_ymin + j *dely
256
257 allocate( lcmesh%pos_ev(lcmesh%Nv,2) )
258 allocate( lcmesh%EToV(lcmesh%Ne,elem%Nv) )
259 allocate( lcmesh%EToE(lcmesh%Ne,elem%Nfaces) )
260 allocate( lcmesh%EToF(lcmesh%Ne,elem%Nfaces) )
261 allocate( lcmesh%BCType(lcmesh%refElem%Nfaces,lcmesh%Ne) )
262 allocate( lcmesh%VMapM(elem%NfpTot, lcmesh%Ne) )
263 allocate( lcmesh%VMapP(elem%NfpTot, lcmesh%Ne) )
264 allocate( lcmesh%MapM(elem%NfpTot, lcmesh%Ne) )
265 allocate( lcmesh%MapP(elem%NfpTot, lcmesh%Ne) )
266
267 lcmesh%BCType(:,:) = bctype_interior
268
269 !----
270
271 call meshutil2d_genrectdomain( lcmesh%pos_ev, lcmesh%EToV, & ! (out)
272 lcmesh%NeX, lcmesh%xmin, lcmesh%xmax, & ! (in)
273 lcmesh%NeY, lcmesh%ymin, lcmesh%ymax ) ! (in)
274
275 !---
276 call meshbase2d_setgeometricinfo(lcmesh, meshrectdom2d_coord_conv, meshrectdom2d_calc_normal )
277
278 !---
279
280 call meshutil2d_genconnectivity( lcmesh%EToE, lcmesh%EToF, & ! (out)
281 lcmesh%EToV, lcmesh%Ne, elem%Nfaces ) ! (in)
282
283 !---
284 call meshutil2d_buildinteriormap( lcmesh%VmapM, lcmesh%VMapP, lcmesh%MapM, lcmesh%MapP, &
285 lcmesh%pos_en, lcmesh%pos_ev, lcmesh%EToE, lcmesh%EtoF, lcmesh%EtoV, &
286 elem%Fmask, lcmesh%Ne, elem%Np, elem%Nfp, elem%Nfaces, lcmesh%Nv )
287
288 call meshutil2d_genpatchboundarymap( lcmesh%VMapB, lcmesh%MapB, lcmesh%VMapP, &
289 lcmesh%pos_en, lcmesh%xmin, lcmesh%xmax, lcmesh%ymin, lcmesh%ymax, &
290 elem%Fmask, lcmesh%Ne, elem%Np, elem%Nfp, elem%Nfaces, lcmesh%Nv)
291
292 return
module FElib / Mesh / Local, Base
integer, parameter, public bctype_interior
module FElib / Mesh / utility for 2D mesh
subroutine, public meshutil2d_buildinteriormap(vmapm, vmapp, mapm, mapp, pos_en, pos_ev, etoe, etof, etov, fmask, ne, np, nfp, nfaces, nv)
subroutine, public meshutil2d_genconnectivity(etoe, etof, etov, ne, nfaces)
subroutine, public meshutil2d_genrectdomain(pos_v, etov, ke_x, xmin, xmax, ke_y, ymin, ymax)
subroutine, public meshutil2d_genpatchboundarymap(vmapb, mapb, vmapp, pos_en, xmin, xmax, ymin, ymax, fmask, ne, np, nfp, nfaces, nv)

References scale_localmesh_base::bctype_interior, scale_mesh_base2d::meshbase2d_setgeometricinfo(), meshrectdom2d_coord_conv(), scale_meshutil_2d::meshutil2d_buildglobalmap(), scale_meshutil_2d::meshutil2d_buildinteriormap(), scale_meshutil_2d::meshutil2d_genconnectivity(), scale_meshutil_2d::meshutil2d_genpatchboundarymap(), and scale_meshutil_2d::meshutil2d_genrectdomain().

Referenced by scale_mesh_rectdom2d::meshrectdom2d::assigndomid(), scale_mesh_cubedom3d::meshcubedom3d_init(), and meshrectdom2d_init().

◆ meshrectdom2d_coord_conv()

subroutine, public scale_mesh_rectdom2d::meshrectdom2d_coord_conv ( real(rp), dimension(elem%np), intent(out) x,
real(rp), dimension(elem%np), intent(out) y,
real(rp), dimension(elem%np), intent(out) xr,
real(rp), dimension(elem%np), intent(out) xs,
real(rp), dimension(elem%np), intent(out) yr,
real(rp), dimension(elem%np), intent(out) ys,
real(rp), dimension(elem%nv), intent(in) vx,
real(rp), dimension(elem%nv), intent(in) vy,
type(elementbase2d), intent(in) elem )

Definition at line 359 of file scale_mesh_rectdom2d.F90.

361
362 implicit none
363
364 type(ElementBase2D), intent(in) :: elem
365 real(RP), intent(out) :: x(elem%Np), y(elem%Np)
366 real(RP), intent(out) :: xr(elem%Np), xs(elem%Np), yr(elem%Np), ys(elem%Np)
367 real(RP), intent(in) :: vx(elem%Nv), vy(elem%Nv)
368
369 !-------------------------------------------------
370
371 x(:) = vx(1) + 0.5_rp*(elem%x1(:) + 1.0_rp)*(vx(2) - vx(1))
372 y(:) = vy(1) + 0.5_rp*(elem%x2(:) + 1.0_rp)*(vy(3) - vy(1))
373
374 xr(:) = 0.5_rp*(vx(2) - vx(1)) !matmul(refElem%Dx1,mesh%x1(:,n))
375 xs(:) = 0.0_rp !matmul(refElem%Dx2,mesh%x1(:,n))
376 yr(:) = 0.0_rp !matmul(refElem%Dx1,mesh%x2(:,n))
377 ys(:) = 0.5_rp*(vy(3) - vy(1)) !matmul(refElem%Dx2,mesh%x2(:,n))
378
379 return

Referenced by scale_mesh_rectdom2d::meshrectdom2d::assigndomid(), and meshrectdom2d_setuplocaldom().