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

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

Data Types

type  meshbase2d
 
interface  meshbase2d_generate
 

Functions/Subroutines

subroutine, public meshbase2d_init (this, refelem, nlocalmeshperprc, nprocs, myrank)
 
subroutine, public meshbase2d_final (this)
 
subroutine, public meshbase2d_setgeometricinfo (lcmesh, coord_conv, calc_normal)
 

Variables

integer, public meshbase2d_dimtype_num = 4
 
integer, public meshbase2d_dimtypeid_x = 1
 
integer, public meshbase2d_dimtypeid_y = 2
 
integer, public meshbase2d_dimtypeid_xy = 3
 
integer, public meshbase2d_dimtypeid_xyt = 4
 

Detailed Description

module FElib / Mesh / Base 2D

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

Function/Subroutine Documentation

◆ meshbase2d_init()

subroutine, public scale_mesh_base2d::meshbase2d_init ( class(meshbase2d), intent(inout) this,
class(elementbase2d), intent(in), target refelem,
integer, intent(in) nlocalmeshperprc,
integer, intent(in), optional nprocs,
integer, intent(in), optional myrank )

Definition at line 74 of file scale_mesh_base2d.F90.

77
78 implicit none
79
80 class(MeshBase2D), intent(inout) :: this
81 class(ElementBase2D), intent(in), target :: refElem
82 integer, intent(in) :: NLocalMeshPerPrc
83 integer, intent(in), optional :: nprocs
84 integer, intent(in), optional :: myrank
85
86 integer :: n
87 !-----------------------------------------------------------------------------
88
89 this%refElem2D => refelem
90 call meshbase_init( this, &
91 meshbase2d_dimtype_num, refelem, &
92 nlocalmeshperprc, 4, &
93 nprocs )
94
95 allocate( this%lcmesh_list(this%LOCAL_MESH_NUM) )
96 do n=1, this%LOCAL_MESH_NUM
97 call localmesh2d_init( this%lcmesh_list(n), n, refelem, myrank )
98 end do
99
100 call this%SetDimInfo( meshbase2d_dimtypeid_x, "x", "m", "X-coordinate" )
101 call this%SetDimInfo( meshbase2d_dimtypeid_y, "y", "m", "Y-coordinate" )
102 call this%SetDimInfo( meshbase2d_dimtypeid_xy, "xy", "m", "XY-coordinate" )
103 call this%SetDimInfo( meshbase2d_dimtypeid_xyt, "xyt", "m", "XY-coordinate" )
104
105 return

References scale_localmesh_2d::localmesh2d_init(), meshbase2d_dimtype_num, meshbase2d_dimtypeid_x, meshbase2d_dimtypeid_xy, meshbase2d_dimtypeid_xyt, meshbase2d_dimtypeid_y, and scale_mesh_base::meshbase_init().

Referenced by scale_mesh_base2d::meshbase2d_generate::meshbase2d_generate(), scale_mesh_cubedom3d::meshcubedom3d_init(), scale_mesh_cubedspheredom2d::meshcubedspheredom2d_init(), and scale_mesh_rectdom2d::meshrectdom2d_init().

◆ meshbase2d_final()

subroutine, public scale_mesh_base2d::meshbase2d_final ( class(meshbase2d), intent(inout) this)

Definition at line 109 of file scale_mesh_base2d.F90.

110 implicit none
111
112 class(MeshBase2D), intent(inout) :: this
113 integer :: n
114 !-----------------------------------------------------------------------------
115
116 if ( allocated ( this%lcmesh_list ) ) then
117 do n=1, this%LOCAL_MESH_NUM
118 call localmesh2d_final( this%lcmesh_list(n), this%isGenerated )
119 end do
120
121 deallocate( this%lcmesh_list )
122 end if
123
124 call meshbase_final(this)
125
126 return

References scale_localmesh_2d::localmesh2d_final(), and scale_mesh_base::meshbase_final().

Referenced by scale_mesh_base2d::meshbase2d_generate::meshbase2d_generate(), scale_mesh_cubedspheredom2d::meshcubedspheredom2d_init(), and scale_mesh_rectdom2d::meshrectdom2d_init().

◆ meshbase2d_setgeometricinfo()

subroutine, public scale_mesh_base2d::meshbase2d_setgeometricinfo ( type(localmesh2d), intent(inout) lcmesh,
external subroutine(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) coord_conv,
external subroutine(real(rp), dimension(elem%nfptot,2), intent(out) normal_fn, real(rp), dimension(elem%nfptot,2,2), intent(in) escale_f, integer, dimension(elem%nfp,elem%nfaces), intent(in) fid, type(elementbase2d), intent(in) elem) calc_normal )

Definition at line 144 of file scale_mesh_base2d.F90.

145
146 implicit none
147
148 type(LocalMesh2D), intent(inout) :: lcmesh
149 interface
150 subroutine coord_conv( x, y, xr, xs, yr, ys, &
151 vx, vy, elem )
152 import elementbase2d
153 import rp
154 type(ElementBase2D), intent(in) :: elem
155 real(RP), intent(out) :: x(elem%Np), y(elem%Np)
156 real(RP), intent(out) :: xr(elem%Np), xs(elem%Np), yr(elem%Np), ys(elem%Np)
157 real(RP), intent(in) :: vx(elem%Nv), vy(elem%Nv)
158 end subroutine coord_conv
159 subroutine calc_normal( normal_fn, &
160 Escale_f, fid, elem )
161 import elementbase2d
162 import rp
163 type(ElementBase2D), intent(in) :: elem
164 real(RP), intent(out) :: normal_fn(elem%NfpTot,2)
165 integer, intent(in) :: fid(elem%Nfp,elem%Nfaces)
166 real(RP), intent(in) :: Escale_f(elem%NfpTot,2,2)
167 end subroutine calc_normal
168 end interface
169
170 class(ElementBase2D), pointer :: refElem
171 integer :: ke
172 integer :: f
173 integer :: i, j
174 integer :: d
175 integer :: node_ids(lcmesh%refElem%Nv)
176 real(RP) :: vx(lcmesh%refElem%Nv), vy(lcmesh%refElem%Nv)
177 real(RP) :: xr(lcmesh%refElem%Np), xs(lcmesh%refElem%Np)
178 real(RP) :: yr(lcmesh%refElem%Np), ys(lcmesh%refElem%Np)
179 integer :: fmask(lcmesh%refElem%NfpTot)
180 integer :: fid(lcmesh%refElem2D%Nfp,lcmesh%refElem2D%Nfaces)
181 real(RP) :: Escale_f(lcmesh%refElem%NfpTot,2,2)
182
183 !-----------------------------------------------------------------------------
184
185 refelem => lcmesh%refElem2D
186
187 allocate( lcmesh%pos_en(refelem%Np,lcmesh%Ne,2) )
188 !allocate( mesh%fx(refElem%Nfaces*refElem%Nfp,mesh%Ne) )
189 !allocate( mesh%fy(refElem%Nfaces*refElem%Nfp,mesh%Ne) )
190 allocate( lcmesh%normal_fn(refelem%NfpTot,lcmesh%Ne,2) )
191 allocate( lcmesh%sJ(refelem%NfpTot,lcmesh%Ne) )
192 allocate( lcmesh%J(refelem%Np,lcmesh%Ne) )
193 allocate( lcmesh%Fscale(refelem%NfpTot,lcmesh%Ne) )
194 allocate( lcmesh%Escale(refelem%Np,lcmesh%Ne,2,2) )
195 allocate( lcmesh%Gsqrt(refelem%Np,lcmesh%NeA) )
196 allocate( lcmesh%G_ij(refelem%Np,lcmesh%Ne, 2, 2) )
197 allocate( lcmesh%GIJ (refelem%Np,lcmesh%Ne, 2, 2) )
198 allocate( lcmesh%lon(refelem%Np,lcmesh%Ne) )
199 allocate( lcmesh%lat(refelem%Np,lcmesh%Ne) )
200
201 fmask(:) = reshape(refelem%Fmask, shape(fmask))
202 do f=1, refelem%Nfaces
203 do i=1, refelem%Nfp
204 fid(i,f) = i + (f-1)*refelem%Nfp
205 end do
206 end do
207
208 do ke=1, lcmesh%Ne
209 node_ids(:) = lcmesh%EToV(ke,:)
210 vx(:) = lcmesh%pos_ev(node_ids(:),1)
211 vy(:) = lcmesh%pos_ev(node_ids(:),2)
212 call coord_conv( &
213 lcmesh%pos_en(:,ke,1), lcmesh%pos_en(:,ke,2), xr, xs, yr, ys, & ! (out)
214 vx, vy, refelem ) ! (in)
215
216 lcmesh%J(:,ke) = - xs*yr + xr*ys
217
218 lcmesh%Escale(:,ke,1,1) = ys/lcmesh%J(:,ke)
219 lcmesh%Escale(:,ke,1,2) = - xs/lcmesh%J(:,ke)
220 lcmesh%Escale(:,ke,2,1) = - yr/lcmesh%J(:,ke)
221 lcmesh%Escale(:,ke,2,2) = xr/lcmesh%J(:,ke)
222
223 !* Face
224
225 !
226 !mesh%fx(:,n) = mesh%x(fmask(:),n)
227 !mesh%fy(:,n) = mesh%y(fmask(:),n)
228
229 ! Calculate normal vectors
230 do j=1, 2
231 do i=1, 2
232 escale_f(:,i,j) = lcmesh%Escale(fmask(:),ke,i,j)
233 end do
234 end do
235
236 call calc_normal( lcmesh%normal_fn(:,ke,:), & ! (out)
237 escale_f, fid, refelem ) ! (in)
238
239 lcmesh%sJ(:,ke) = sqrt( lcmesh%normal_fn(:,ke,1)**2 + lcmesh%normal_fn(:,ke,2)**2 )
240 do d=1, 2
241 lcmesh%normal_fn(:,ke,d) = lcmesh%normal_fn(:,ke,d)/lcmesh%sJ(:,ke)
242 end do
243 lcmesh%sJ(:,ke) = lcmesh%sJ(:,ke)*lcmesh%J(fmask(:),ke)
244
245 lcmesh%Fscale(:,ke) = lcmesh%sJ(:,ke)/lcmesh%J(fmask(:),ke)
246 end do
247
248 !$omp parallel
249 !$omp workshare
250 lcmesh%Gsqrt (:,:) = 1.0_rp
251 lcmesh%GIJ (:,:,1,1) = 1.0_rp
252 lcmesh%GIJ (:,:,2,1) = 0.0_rp
253 lcmesh%GIJ (:,:,1,2) = 0.0_rp
254 lcmesh%GIJ (:,:,2,2) = 1.0_rp
255 lcmesh%G_ij (:,:,1,1) = 1.0_rp
256 lcmesh%G_ij (:,:,2,1) = 0.0_rp
257 lcmesh%G_ij (:,:,1,2) = 0.0_rp
258 lcmesh%G_ij (:,:,2,2) = 1.0_rp
259 !$omp end workshare
260 !$omp end parallel
261
262 return

Referenced by scale_mesh_base2d::meshbase2d_generate::meshbase2d_generate(), scale_mesh_cubedspheredom2d::meshcubedspheredom2d_setuplocaldom(), and scale_mesh_rectdom2d::meshrectdom2d_setuplocaldom().

Variable Documentation

◆ meshbase2d_dimtype_num

integer, public scale_mesh_base2d::meshbase2d_dimtype_num = 4

Definition at line 56 of file scale_mesh_base2d.F90.

56 integer, public :: MESHBASE2D_DIMTYPE_NUM = 4

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

◆ meshbase2d_dimtypeid_x

integer, public scale_mesh_base2d::meshbase2d_dimtypeid_x = 1

◆ meshbase2d_dimtypeid_y

integer, public scale_mesh_base2d::meshbase2d_dimtypeid_y = 2

◆ meshbase2d_dimtypeid_xy

integer, public scale_mesh_base2d::meshbase2d_dimtypeid_xy = 3

◆ meshbase2d_dimtypeid_xyt

integer, public scale_mesh_base2d::meshbase2d_dimtypeid_xyt = 4