FE-Project
All Classes Namespaces Files Functions Variables Pages
scale_mesh_base2d.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_2d, only: &
20
21 use scale_mesh_base, only: &
23
25
26 !-----------------------------------------------------------------------------
27 implicit none
28 private
29
30 !-----------------------------------------------------------------------------
31 !
32 !++ Public type & procedure
33 !
34 type, abstract, public, extends(meshbase) :: meshbase2d
35 type(localmesh2d), allocatable :: lcmesh_list(:)
36 type(elementbase2d), pointer :: refelem2d
37 contains
38 procedure(meshbase2d_generate), deferred :: generate
39 procedure :: getlocalmesh => meshbase2d_get_localmesh
40 end type meshbase2d
41
42 interface
43 subroutine meshbase2d_generate(this)
44 import meshbase2d
45 class(meshbase2d), intent(inout), target :: this
46 end subroutine meshbase2d_generate
47 end interface
48
51
52 !-----------------------------------------------------------------------------
53 !
54 !++ Public parameters & variables
55 !
56 integer, public :: meshbase2d_dimtype_num = 4
57 integer, public :: meshbase2d_dimtypeid_x = 1
58 integer, public :: meshbase2d_dimtypeid_y = 2
59 integer, public :: meshbase2d_dimtypeid_xy = 3
60 integer, public :: meshbase2d_dimtypeid_xyt = 4
61
62 !-----------------------------------------------------------------------------
63 !
64 !++ Private procedure
65 !
66
67 !-----------------------------------------------------------------------------
68 !
69 !++ Private parameters & variables
70 !
71
72contains
73!OCL SERIAL
74 subroutine meshbase2d_init(this, &
75 refElem, NLocalMeshPerPrc, &
76 nprocs, myrank )
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
106 end subroutine meshbase2d_init
107
108!OCL SERIAL
109 subroutine meshbase2d_final( this )
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
127 end subroutine meshbase2d_final
128
129!OCL SERIAL
130 subroutine meshbase2d_get_localmesh( this, id, ptr_lcmesh )
132 implicit none
133
134 class(meshbase2d), target, intent(in) :: this
135 integer, intent(in) :: id
136 class(localmeshbase), pointer, intent(out) :: ptr_lcmesh
137 !-------------------------------------------------------------
138
139 ptr_lcmesh => this%lcmesh_list(id)
140 return
141 end subroutine meshbase2d_get_localmesh
142
143!OCL SERIAL
144 subroutine meshbase2d_setgeometricinfo( lcmesh, coord_conv, calc_normal )
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
263 end subroutine meshbase2d_setgeometricinfo
264
265end module scale_mesh_base2d
module FElib / Element / Base
module FElib / Mesh / Local 2D
subroutine, public localmesh2d_final(this, is_generated)
subroutine, public localmesh2d_init(this, lcdomid, refelem, myrank)
module FElib / Mesh / Local, Base
module FElib / Mesh / Base 2D
subroutine, public meshbase2d_final(this)
subroutine, public meshbase2d_init(this, refelem, nlocalmeshperprc, nprocs, myrank)
subroutine, public meshbase2d_setgeometricinfo(lcmesh, coord_conv, calc_normal)
integer, public meshbase2d_dimtypeid_xy
integer, public meshbase2d_dimtypeid_xyt
integer, public meshbase2d_dimtypeid_x
integer, public meshbase2d_dimtype_num
integer, public meshbase2d_dimtypeid_y
module FElib / Mesh / Base
subroutine, public meshbase_final(this)
subroutine, public meshbase_init(this, ndimtype, refelem, nlocalmeshperprc, nsidetile, nprocs)