FE-Project
Loading...
Searching...
No Matches
scale_mesh_base3d.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
9#include "scaleFElib.h"
11
12 !-----------------------------------------------------------------------------
13 !
14 !++ used modules
15 !
16 use scale_precision
17
18 use scale_localmesh_3d, only: &
20 use scale_localmesh_2d, only: &
22 use scale_mesh_base, only: &
24 use scale_mesh_base2d, only: &
27
28 !-----------------------------------------------------------------------------
29 implicit none
30 private
31
32 !-----------------------------------------------------------------------------
33 !
34 !++ Public type & procedure
35 !
36 type, abstract, public, extends(meshbase) :: meshbase3d
37 type(localmesh3d), allocatable :: lcmesh_list(:)
38 type(elementbase3d), pointer :: refelem3d
39 contains
40 procedure(meshbase3d_generate), deferred :: generate
41 procedure(meshbase3d_getmesh2d), deferred :: getmesh2d
42 procedure(meshbase3d_set_geometric_with_vcoord), deferred :: set_geometric_with_vcoord
43 procedure :: getlocalmesh => meshbase3d_get_localmesh
44 end type meshbase3d
45
46 interface
47 subroutine meshbase3d_generate(this)
48 import meshbase3d
49 class(meshbase3d), intent(inout), target :: this
50 end subroutine meshbase3d_generate
51
52 subroutine meshbase3d_getmesh2d(this, ptr_mesh2D)
53 import meshbase3d
54 import meshbase2d
55 class(meshbase3d), intent(in), target :: this
56 class(meshbase2d), pointer, intent(out) :: ptr_mesh2D
57 end subroutine meshbase3d_getmesh2d
58
59 subroutine meshbase3d_set_geometric_with_vcoord(this, lcdomID, GsqrtV_lc, zlev_lc, G13_lc, G23_lc)
60 import rp
61 import meshbase3d
62 class(meshbase3d), intent(inout), target :: this
63 integer, intent(in) :: lcdomID
64 real(RP), intent(in) :: GsqrtV_lc(this%refElem3D%Np,this%lcmesh_list(lcdomID)%NeA)
65 real(RP), intent(in) :: zlev_lc(this%refElem3D%Np,this%lcmesh_list(lcdomID)%NeA)
66 real(RP), intent(in) :: G13_lc(this%refElem3D%Np,this%lcmesh_list(lcdomID)%NeA)
67 real(RP), intent(in) :: G23_lc(this%refElem3D%Np,this%lcmesh_list(lcdomID)%NeA)
68 end subroutine meshbase3d_set_geometric_with_vcoord
69 end interface
70
73
74 !-----------------------------------------------------------------------------
75 !
76 !++ Public parameters & variables
77 !
78 integer, public :: meshbase3d_dimtype_num = 8
79 integer, public :: meshbase3d_dimtypeid_x = 1
80 integer, public :: meshbase3d_dimtypeid_y = 2
81 integer, public :: meshbase3d_dimtypeid_z = 3
82 integer, public :: meshbase3d_dimtypeid_zt = 4
83 integer, public :: meshbase3d_dimtypeid_xy = 5
84 integer, public :: meshbase3d_dimtypeid_xyt = 6
85 integer, public :: meshbase3d_dimtypeid_xyz = 7
86 integer, public :: meshbase3d_dimtypeid_xyzt = 8
87
88 !-----------------------------------------------------------------------------
89 !
90 !++ Private procedure
91 !
92
93 !-----------------------------------------------------------------------------
94 !
95 !++ Private parameters & variables
96 !
97
98contains
99!OCL SERIAL
100 subroutine meshbase3d_init(this, &
101 refElem, NLocalMeshPerPrc, NsideTile, &
102 nproc, myrank )
103
104 implicit none
105
106 class(meshbase3d), intent(inout) :: this
107 class(elementbase3d), intent(in), target :: refelem
108 integer, intent(in) :: nlocalmeshperprc
109 integer, intent(in) :: nsidetile
110 integer, intent(in), optional :: nproc
111 integer, intent(in), optional :: myrank
112
113 integer :: n
114 !-----------------------------------------------------------------------------
115
116 this%refElem3D => refelem
117 call meshbase_init( this, &
118 meshbase3d_dimtype_num, refelem, &
119 nlocalmeshperprc, nsidetile, nproc )
120
121 allocate( this%lcmesh_list(this%LOCAL_MESH_NUM) )
122 do n=1, this%LOCAL_MESH_NUM
123 call localmesh3d_init( this%lcmesh_list(n), n, refelem, myrank )
124 end do
125
126 call this%SetDimInfo( meshbase3d_dimtypeid_x, "x", "m", "X-coordinate" )
127 call this%SetDimInfo( meshbase3d_dimtypeid_y, "y", "m", "Y-coordinate" )
128 call this%SetDimInfo( meshbase3d_dimtypeid_z, "z", "m", "Z-coordinate" )
129 call this%SetDimInfo( meshbase3d_dimtypeid_zt, "z", "m", "Z-coordinate" )
130 call this%SetDimInfo( meshbase3d_dimtypeid_xy, "xy", "m", "XY-coordinate" )
131 call this%SetDimInfo( meshbase3d_dimtypeid_xyt, "xyt", "m", "XY-coordinate" )
132 call this%SetDimInfo( meshbase3d_dimtypeid_xyz, "xyz", "m", "XYZ-coordinate" )
133 call this%SetDimInfo( meshbase3d_dimtypeid_xyzt, "xyzt", "m", "XYZ-coordinate" )
134
135 return
136 end subroutine meshbase3d_init
137
138!OCL SERIAL
139 subroutine meshbase3d_final( this )
140 implicit none
141 class(meshbase3d), intent(inout) :: this
142
143 integer :: n
144 !-----------------------------------------------------------------------------
145
146 if ( allocated ( this%lcmesh_list ) ) then
147 do n=1, this%LOCAL_MESH_NUM
148 call localmesh3d_final( this%lcmesh_list(n), this%isGenerated )
149 end do
150
151 deallocate( this%lcmesh_list )
152 end if
153
154 call meshbase_final(this)
155
156 return
157 end subroutine meshbase3d_final
158
159!OCL SERIAL
160 subroutine meshbase3d_get_localmesh( this, id, ptr_lcmesh )
162 implicit none
163
164 class(meshbase3d), target, intent(in) :: this
165 integer, intent(in) :: id
166 class(localmeshbase), pointer, intent(out) :: ptr_lcmesh
167 !-------------------------------------------------------------
168
169 ptr_lcmesh => this%lcmesh_list(id)
170 return
171 end subroutine meshbase3d_get_localmesh
172
173!OCL SERIAL
174 subroutine meshbase3d_setgeometricinfo( lcmesh, coord_conv, calc_normal )
175 implicit none
176
177 type(localmesh3d), intent(inout) :: lcmesh
178 interface
179 subroutine coord_conv( x, y, z, xX, xY, xZ, yX, yY, yZ, zX, zY, zZ, &
180 vx, vy, vz, elem )
181 import elementbase3d
182 import rp
183 type(elementbase3d), intent(in) :: elem
184 real(rp), intent(out) :: x(elem%np), y(elem%np), z(elem%np)
185 real(rp), intent(out) :: xx(elem%np), xy(elem%np), xz(elem%np)
186 real(rp), intent(out) :: yx(elem%np), yy(elem%np), yz(elem%np)
187 real(rp), intent(out) :: zx(elem%np), zy(elem%np), zz(elem%np)
188 real(rp), intent(in) :: vx(elem%nv), vy(elem%nv), vz(elem%nv)
189 end subroutine coord_conv
190 subroutine calc_normal( normal_fn, &
191 Escale_f, fid_h, fid_v, elem )
192 import elementbase3d
193 import rp
194 type(elementbase3d), intent(in) :: elem
195 real(rp), intent(out) :: normal_fn(elem%nfptot,3)
196 integer, intent(in) :: fid_h(elem%nfp_h,elem%nfaces_h)
197 integer, intent(in) :: fid_v(elem%nfp_v,elem%nfaces_v)
198 real(rp), intent(in) :: escale_f(elem%nfptot,3,3)
199 end subroutine calc_normal
200 end interface
201
202 class(elementbase3d), pointer :: refelem
203 integer :: ke, ke2d
204 integer :: f
205 integer :: i, j
206 integer :: d
207 integer :: fmask(lcmesh%refelem%nfptot)
208 integer :: fid_h(lcmesh%refelem3d%nfp_h,lcmesh%refelem3d%nfaces_h)
209 integer :: fid_v(lcmesh%refelem3d%nfp_v,lcmesh%refelem3d%nfaces_v)
210 real(rp) :: escale_f(lcmesh%refelem%nfptot,3,3)
211
212 integer :: node_ids(lcmesh%refelem%nv)
213 real(rp) :: vx(lcmesh%refelem%nv), vy(lcmesh%refelem%nv), vz(lcmesh%refelem%nv)
214 real(rp) :: xx(lcmesh%refelem%np), xy(lcmesh%refelem%np), xz(lcmesh%refelem%np)
215 real(rp) :: yx(lcmesh%refelem%np), yy(lcmesh%refelem%np), yz(lcmesh%refelem%np)
216 real(rp) :: zx(lcmesh%refelem%np), zy(lcmesh%refelem%np), zz(lcmesh%refelem%np)
217 !-----------------------------------------------------------------------------
218
219 refelem => lcmesh%refElem3D
220
221 allocate( lcmesh%pos_en(refelem%Np,lcmesh%Ne,3) )
222 !allocate( mesh%fx(refElem%Nfaces*refElem%Nfp,mesh%Ne) )
223 !allocate( mesh%fy(refElem%Nfaces*refElem%Nfp,mesh%Ne) )
224 allocate( lcmesh%normal_fn(refelem%NfpTot,lcmesh%Ne,3) )
225 allocate( lcmesh%sJ(refelem%NfpTot,lcmesh%Ne) )
226 allocate( lcmesh%J(refelem%Np,lcmesh%Ne) )
227 allocate( lcmesh%Fscale(refelem%NfpTot,lcmesh%Ne) )
228 allocate( lcmesh%Escale(refelem%Np,lcmesh%Ne,3,3) )
229 allocate( lcmesh%zS(refelem%Np,lcmesh%Ne) )
230 allocate( lcmesh%Sz(refelem%Np,lcmesh%Ne) )
231 allocate( lcmesh%Gsqrt(refelem%Np,lcmesh%NeA) )
232 allocate( lcmesh%GsqrtH(refelem%Nfp_v,lcmesh%Ne2D) )
233 allocate( lcmesh%G_ij(refelem%Nfp_v,lcmesh%Ne2D,2,2) )
234 allocate( lcmesh%GIJ (refelem%Nfp_v,lcmesh%Ne2D,2,2) )
235 allocate( lcmesh%GI3 (refelem%Np,lcmesh%NeA,2) )
236 allocate( lcmesh%zlev(refelem%Np,lcmesh%Ne) )
237 allocate( lcmesh%gam(refelem%Np,lcmesh%NeA) )
238
239 allocate( lcmesh%lon2D(refelem%Nfp_v,lcmesh%Ne2D) )
240 allocate( lcmesh%lat2D(refelem%Nfp_v,lcmesh%Ne2D) )
241
242 do f=1, refelem%Nfaces_h
243 do i=1, refelem%Nfp_h
244 fid_h(i,f) = i + (f-1)*refelem%Nfp_h
245 fmask(fid_h(i,f)) = refelem%Fmask_h(i,f)
246 end do
247 end do
248 do f=1, refelem%Nfaces_v
249 do i=1, refelem%Nfp_v
250 fid_v(i,f) = i + refelem%Nfaces_h*refelem%Nfp_h + (f-1)*refelem%Nfp_v
251 fmask(fid_v(i,f)) = refelem%Fmask_v(i,f)
252 end do
253 end do
254
255 !$omp parallel private( ke, node_ids, vx, vy, vz, &
256 !$omp xX, xY, xZ, yX, yY, yZ, zX, zY, zZ, &
257 !$omp i, j, Escale_f, d )
258
259 !$omp do
260 do ke=1, lcmesh%Ne
261 node_ids(:) = lcmesh%EToV(ke,:)
262 vx(:) = lcmesh%pos_ev(node_ids(:),1)
263 vy(:) = lcmesh%pos_ev(node_ids(:),2)
264 vz(:) = lcmesh%pos_ev(node_ids(:),3)
265 call coord_conv( &
266 lcmesh%pos_en(:,ke,1), lcmesh%pos_en(:,ke,2), lcmesh%pos_en(:,ke,3), & ! (out)
267 xx, xy, xz, yx, yy, yz, zx, zy, zz, & ! (out)
268 vx, vy, vz, refelem ) ! (in)
269
270 lcmesh%J(:,ke) = xx(:)*(yy(:)*zz(:) - zy(:)*yz) &
271 - yx(:)*(xy(:)*zz(:) - zy(:)*xz) &
272 + zx(:)*(xy(:)*yz(:) - yy(:)*xz)
273
274 lcmesh%Escale(:,ke,1,1) = (yy(:)*zz(:) - zy(:)*yz(:))/lcmesh%J(:,ke)
275 lcmesh%Escale(:,ke,1,2) = - (xy(:)*zz(:) - zy(:)*xz(:))/lcmesh%J(:,ke)
276 lcmesh%Escale(:,ke,1,3) = (xy(:)*yz(:) - yy(:)*xz(:))/lcmesh%J(:,ke)
277
278 lcmesh%Escale(:,ke,2,1) = - (yx(:)*zz(:) - zx(:)*yz(:))/lcmesh%J(:,ke)
279 lcmesh%Escale(:,ke,2,2) = (xx(:)*zz(:) - zx(:)*xz(:))/lcmesh%J(:,ke)
280 lcmesh%Escale(:,ke,2,3) = - (xx(:)*yz(:) - yx(:)*xz(:))/lcmesh%J(:,ke)
281
282 lcmesh%Escale(:,ke,3,1) = (yx(:)*zy(:) - zx(:)*yy(:))/lcmesh%J(:,ke)
283 lcmesh%Escale(:,ke,3,2) = - (xx(:)*zy(:) - zx(:)*xy(:))/lcmesh%J(:,ke)
284 lcmesh%Escale(:,ke,3,3) = (xx(:)*yy(:) - yx(:)*xy(:))/lcmesh%J(:,ke)
285
286 !* Face
287
288 !
289 !mesh%fx(:,n) = mesh%x(fmask(:),n)
290 !mesh%fy(:,n) = mesh%y(fmask(:),n)
291
292 ! Calculate normal vectors
293 do j=1, 3
294 do i=1, 3
295 escale_f(:,i,j) = lcmesh%Escale(fmask(:),ke,i,j)
296 end do
297 end do
298 call calc_normal( lcmesh%normal_fn(:,ke,:), & ! (out)
299 escale_f, fid_h, fid_v, refelem ) ! (in)
300
301 lcmesh%sJ(:,ke) = sqrt( &
302 lcmesh%normal_fn(:,ke,1)**2 + lcmesh%normal_fn(:,ke,2)**2 + lcmesh%normal_fn(:,ke,3)**2 )
303 do d=1, 3
304 lcmesh%normal_fn(:,ke,d) = lcmesh%normal_fn(:,ke,d)/lcmesh%sJ(:,ke)
305 end do
306 lcmesh%sJ(:,ke) = lcmesh%sJ(:,ke)*lcmesh%J(fmask(:),ke)
307
308 lcmesh%Fscale(:,ke) = lcmesh%sJ(:,ke)/lcmesh%J(fmask(:),ke)
309 lcmesh%zlev(:,ke) = lcmesh%pos_en(:,ke,3)
310 end do
311 !$omp end do
312
313 !$omp workshare
314 lcmesh%Gsqrt (:,:) = 1.0_rp
315 lcmesh%GsqrtH(:,:) = 1.0_rp
316 lcmesh%GIJ (:,:,1,1) = 1.0_rp
317 lcmesh%GIJ (:,:,2,1) = 0.0_rp
318 lcmesh%GIJ (:,:,1,2) = 0.0_rp
319 lcmesh%GIJ (:,:,2,2) = 1.0_rp
320 lcmesh%G_ij (:,:,1,1) = 1.0_rp
321 lcmesh%G_ij (:,:,2,1) = 0.0_rp
322 lcmesh%G_ij (:,:,1,2) = 0.0_rp
323 lcmesh%G_ij (:,:,2,2) = 1.0_rp
324 lcmesh%GI3 (:,:,1) = 0.0_rp
325 lcmesh%GI3 (:,:,2) = 0.0_rp
326 lcmesh%gam (:,:) = 1.0_rp
327 !$omp end workshare
328
329 !$omp end parallel
330
331 return
332 end subroutine meshbase3d_setgeometricinfo
333
334end module scale_mesh_base3d
module FElib / Element / Base
module FElib / Mesh / Local 2D
module FElib / Mesh / Local 3D
subroutine, public localmesh3d_final(this, is_generated)
subroutine, public localmesh3d_init(this, lcdomid, refelem, myrank)
module FElib / Mesh / Local, Base
module FElib / Mesh / Base 2D
module FElib / Mesh / Base 3D
subroutine, public meshbase3d_init(this, refelem, nlocalmeshperprc, nsidetile, nproc, myrank)
integer, public meshbase3d_dimtypeid_y
integer, public meshbase3d_dimtypeid_zt
integer, public meshbase3d_dimtypeid_z
integer, public meshbase3d_dimtype_num
subroutine, public meshbase3d_final(this)
integer, public meshbase3d_dimtypeid_xy
subroutine, public meshbase3d_setgeometricinfo(lcmesh, coord_conv, calc_normal)
integer, public meshbase3d_dimtypeid_xyt
integer, public meshbase3d_dimtypeid_xyz
integer, public meshbase3d_dimtypeid_x
integer, public meshbase3d_dimtypeid_xyzt
module FElib / Mesh / Base
subroutine, public meshbase_final(this)
subroutine, public meshbase_init(this, ndimtype, refelem, nlocalmeshperprc, nsidetile, nprocs)