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

module FElib / Mesh / Cubic 3D domain More...

Data Types

type  meshcubedom3d

Functions/Subroutines

subroutine meshcubedom3d_init (this, negx, negy, negz, dom_xmin, dom_xmax, dom_ymin, dom_ymax, dom_zmin, dom_zmax, isperiodicx, isperiodicy, isperiodicz, refelem, nlocalmeshperprc, nprcx, nprcy, nproc, myrank, fz)
 Initialize an object to manage a cubic 3D domain.
subroutine, public meshcubedom3d_coord_conv (x, y, z, xx, xy, xz, yx, yy, yz, zx, zy, zz, vx, vy, vz, elem)

Detailed Description

module FElib / Mesh / Cubic 3D domain

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

Function/Subroutine Documentation

◆ meshcubedom3d_init()

subroutine scale_mesh_cubedom3d::meshcubedom3d_init ( class(meshcubedom3d), intent(inout) this,
integer, intent(in) negx,
integer, intent(in) negy,
integer, intent(in) negz,
real(rp), intent(in) dom_xmin,
real(rp), intent(in) dom_xmax,
real(rp), intent(in) dom_ymin,
real(rp), intent(in) dom_ymax,
real(rp), intent(in) dom_zmin,
real(rp), intent(in) dom_zmax,
logical, intent(in) isperiodicx,
logical, intent(in) isperiodicy,
logical, intent(in) isperiodicz,
type(hexahedralelement), 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,
real(rp), dimension(negz+1), intent(in), optional fz )

Initialize an object to manage a cubic 3D domain.

Definition at line 96 of file scale_mesh_cubedom3d.F90.

103
104 implicit none
105
106 class(MeshCubeDom3D), intent(inout) :: this
107 integer, intent(in) :: NeGX
108 integer, intent(in) :: NeGY
109 integer, intent(in) :: NeGZ
110 real(RP), intent(in) :: dom_xmin
111 real(RP), intent(in) :: dom_xmax
112 real(RP), intent(in) :: dom_ymin
113 real(RP), intent(in) :: dom_ymax
114 real(RP), intent(in) :: dom_Zmin
115 real(RP), intent(in) :: dom_zmax
116 logical, intent(in) :: isPeriodicX
117 logical, intent(in) :: isPeriodicY
118 logical, intent(in) :: isPeriodicZ
119 type(HexahedralElement), intent(in), target :: refElem
120 integer, intent(in) :: NLocalMeshPerPrc
121 integer, intent(in) :: NprcX
122 integer, intent(in) :: NprcY
123 integer, intent(in), optional :: nproc
124 integer, intent(in), optional :: myrank
125 real(RP), intent(in), optional :: FZ(NeGZ+1)
126
127 integer :: k
128 real(RP) :: dz
129 !-----------------------------------------------------------------------------
130
131 this%NeGX = negx
132 this%NeGY = negy
133 this%NeGZ = negz
134
135 this%xmin_gl = dom_xmin
136 this%xmax_gl = dom_xmax
137 this%ymin_gl = dom_ymin
138 this%ymax_gl = dom_ymax
139 this%zmin_gl = dom_zmin
140 this%zmax_gl = dom_zmax
141 this%dom_vol = (this%xmax_gl - this%xmin_gl) * (this%ymax_gl - this%ymin_gl) * (this%zmax_gl - this%zmin_gl)
142
143 this%isPeriodicX = isperiodicx
144 this%isPeriodicY = isperiodicy
145 this%isPeriodicZ = isperiodicz
146
147 this%NprcX = nprcx
148 this%NprcY = nprcy
149 this%NprcZ = 1
150
151 !- Fz
152 allocate( this%FZ(this%NeGZ+1) )
153 if ( present(fz) ) then
154 this%FZ(:) = fz(:)
155 else
156 this%FZ(1 ) = dom_zmin
157 this%FZ(this%NeGZ+1) = dom_zmax
158 dz = (dom_zmax - dom_zmin) / dble(this%NeGZ)
159 do k=2, this%NeGZ
160 this%FZ(k) = this%FZ(k-1) + dz
161 end do
162 end if
163
164 !--
165 call meshbase3d_init( this, refelem, nlocalmeshperprc, 6, &
166 nproc, myrank )
167
168 !--- 2D mesh
169
170 call this%refElem2D%Init( this%refElem3D%PolyOrder_h, refelem%IsLumpedMatrix() )
171
172 call this%mesh2D%Init( this%NeGX, this%NeGY, dom_xmin, dom_xmax, dom_ymin, dom_ymax, &
173 isperiodicx, isperiodicy, this%refElem2D, nlocalmeshperprc, &
174 nprcx, nprcy, nproc, myrank )
175
176 return

References scale_localmesh_base::bctype_interior, scale_mesh_base3d::meshbase3d_final(), scale_mesh_base3d::meshbase3d_init(), scale_mesh_base3d::meshbase3d_setgeometricinfo(), meshcubedom3d_coord_conv(), scale_mesh_rectdom2d::meshrectdom2d_setuplocaldom(), scale_meshutil_3d::meshutil3d_buildglobalmap(), scale_meshutil_3d::meshutil3d_buildinteriormap(), scale_meshutil_3d::meshutil3d_genconnectivity(), scale_meshutil_3d::meshutil3d_gencubedomain(), and scale_meshutil_3d::meshutil3d_genpatchboundarymap().

◆ meshcubedom3d_coord_conv()

subroutine, public scale_mesh_cubedom3d::meshcubedom3d_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) z,
real(rp), dimension(elem%np), intent(out) xx,
real(rp), dimension(elem%np), intent(out) xy,
real(rp), dimension(elem%np), intent(out) xz,
real(rp), dimension(elem%np), intent(out) yx,
real(rp), dimension(elem%np), intent(out) yy,
real(rp), dimension(elem%np), intent(out) yz,
real(rp), dimension(elem%np), intent(out) zx,
real(rp), dimension(elem%np), intent(out) zy,
real(rp), dimension(elem%np), intent(out) zz,
real(rp), dimension(elem%nv), intent(in) vx,
real(rp), dimension(elem%nv), intent(in) vy,
real(rp), dimension(elem%nv), intent(in) vz,
type(elementbase3d), intent(in) elem )

Definition at line 514 of file scale_mesh_cubedom3d.F90.

516
517 implicit none
518
519 type(ElementBase3D), intent(in) :: elem
520 real(RP), intent(out) :: x(elem%Np), y(elem%Np), z(elem%Np)
521 real(RP), intent(out) :: xX(elem%Np), xY(elem%Np), xZ(elem%Np)
522 real(RP), intent(out) :: yX(elem%Np), yY(elem%Np), yZ(elem%Np)
523 real(RP), intent(out) :: zX(elem%Np), zY(elem%Np), zZ(elem%Np)
524 real(RP), intent(in) :: vx(elem%Nv), vy(elem%Nv), vz(elem%Nv)
525
526 !-------------------------------------------------
527
528 x(:) = vx(1) + 0.5_rp*(elem%x1(:) + 1.0_rp)*(vx(2) - vx(1))
529 y(:) = vy(1) + 0.5_rp*(elem%x2(:) + 1.0_rp)*(vy(3) - vy(1))
530 z(:) = vz(1) + 0.5_rp*(elem%x3(:) + 1.0_rp)*(vz(5) - vz(1))
531
532 xx(:) = 0.5_rp*(vx(2) - vx(1)) !matmul(refElem%Dx1,mesh%x1(:,n))
533 xy(:) = 0.0_rp !matmul(refElem%Dx2,mesh%x1(:,n))
534 xz(:) = 0.0_rp !matmul(refElem%Dx3,mesh%x1(:,n))
535 yx(:) = 0.0_rp !matmul(refElem%Dx1,mesh%x2(:,n))
536 yy(:) = 0.5_rp*(vy(3) - vy(1)) !matmul(refElem%Dx2,mesh%x2(:,n))
537 yz(:) = 0.0_rp !matmul(refElem%Dx3,mesh%x2(:,n))
538 zx(:) = 0.0_rp !matmul(refElem%Dx1,mesh%x3(:,n))
539 zy(:) = 0.0_rp !matmul(refElem%Dx2,mesh%x3(:,n))
540 zz(:) = 0.5_rp*(vz(5) - vz(1)) !matmul(refElem%Dx3,mesh%x3(:,n))
541
542 return

Referenced by meshcubedom3d_init(), and scale_mesh_cubedom3d::meshcubedom3d::set_geometric_with_vcoord().