FE-Project
Loading...
Searching...
No Matches
scale_mesh_rectdom2d.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
9!-------------------------------------------------------------------------------
10#include "scaleFElib.h"
12
13 !-----------------------------------------------------------------------------
14 !
15 !++ used modules
16 !
17 use scale_precision
18
19 use scale_mesh_base2d, only: &
22
23 use scale_localmesh_2d, only: &
27
28 !-----------------------------------------------------------------------------
29 implicit none
30 private
31
32 !-----------------------------------------------------------------------------
33 !
34 !++ Public type & procedure
35 !
36 type, extends(meshbase2d), public :: meshrectdom2d
37 integer :: negx
38 integer :: negy
39
40 integer :: nprcx
41 integer :: nprcy
42
43 real(rp), public :: xmin_gl, xmax_gl
44 real(rp), public :: ymin_gl, ymax_gl
45 integer, allocatable :: rcdomij2lcmeshid(:,:)
46
47 logical :: isperiodicx
48 logical :: isperiodicy
49 contains
50 procedure :: init => meshrectdom2d_init
51 procedure :: final => meshrectdom2d_final
52 procedure :: generate => meshrectdom2d_generate
53 procedure :: assigndomid => meshrectdom2d_assigndomid
54 end type meshrectdom2d
55
58
59 !-----------------------------------------------------------------------------
60 !
61 !++ Public parameters & variables
62 !
63
64 !-----------------------------------------------------------------------------
65 !
66 !++ Private procedure
67 !
68
69 !-----------------------------------------------------------------------------
70 !
71 !++ Private parameters & variables
72 !
73
74contains
75 subroutine meshrectdom2d_init(this, &
76 NeGX, NeGY, &
77 dom_xmin, dom_xmax, dom_ymin, dom_ymax, &
78 isPeriodicX, isPeriodicY, &
79 refElem, NLocalMeshPerPrc, &
80 NprcX, NprcY, &
81 nproc, myrank )
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
122 end subroutine meshrectdom2d_init
123
124 subroutine meshrectdom2d_final( this )
125
126 class(meshrectdom2d), intent(inout) :: this
127
128 integer :: n
129
130 !-----------------------------------------------------------------------------
131
132 if (this%isGenerated) then
133 if ( allocated(this%rcdomIJ2LCMeshID) ) then
134 deallocate( this%rcdomIJ2LCMeshID )
135 end if
136 end if
137
138 call meshbase2d_final( this )
139
140 return
141 end subroutine meshrectdom2d_final
142
143 subroutine meshrectdom2d_generate( this )
144
145 implicit none
146
147 class(meshrectdom2d), intent(inout), target :: this
148
149 integer :: n
150 integer :: p
151 type(localmesh2d), pointer :: mesh
152
153 integer :: tileID_table(this%LOCAL_MESH_NUM, this%PRC_NUM)
154 integer :: panelID_table(this%LOCAL_MESH_NUM*this%PRC_NUM)
155 integer :: pi_table(this%LOCAL_MESH_NUM*this%PRC_NUM)
156 integer :: pj_table(this%LOCAL_MESH_NUM*this%PRC_NUM)
157
158 integer :: TILE_NUM_PER_PANEL
159 real(RP) :: delx, dely
160 integer :: tileID
161
162 !-----------------------------------------------------------------------------
163
164 tile_num_per_panel = this%LOCAL_MESH_NUM_global / 1
165
166
167 !--- Construct the connectivity of patches (only master node)
168
169 call this%AssignDomID( &
170 tileid_table, panelid_table, & ! (out)
171 pi_table, pj_table ) ! (out)
172
173 !--- Setup local meshes managed by my process
174
175 do n=1, this%LOCAL_MESH_NUM
176 mesh => this%lcmesh_list(n)
177 tileid = tileid_table(n, mesh%PRC_myrank+1)
178
179 call meshrectdom2d_setuplocaldom( mesh, &
180 tileid, panelid_table(tileid), &
181 pi_table(tileid), pj_table(tileid), this%NprcX, this%NprcY, &
182 this%xmin_gl, this%xmax_gl, this%ymin_gl, this%ymax_gl, &
183 this%NeGX/this%NprcX, this%NeGY/this%NprcY )
184
185 !---
186 ! write(*,*) "** my_rank=", mesh%PRC_myrank
187 ! write(*,*) " tileID:", mesh%tileID
188 ! write(*,*) " pnlID:", mesh%panelID, "-- i,j (within a panel)=", pi_table(tileID), pj_table(tileID)
189 ! write(*,*) " local mesh:", n, "( total", this%LOCAL_MESH_NUM, ")"
190 ! write(*,*) " panel_connect:", this%tilePanelID_globalMap(:,mesh%tileID)
191 ! write(*,*) " tile_connect:", this%tileID_globalMap(:,mesh%tileID)
192 ! write(*,*) " face_connect:", this%tileFaceID_globalMap(:,mesh%tileID)
193 ! write(*,*) " domain size"
194 ! write(*,*) " NeX, NeY:", mesh%NeX, mesh%NeY
195 ! write(*,*) " [X], [Y]:", mesh%xmin, mesh%xmax, ":", mesh%ymin, mesh%ymax
196 end do
197
198 this%isGenerated = .true.
199
200 return
201 end subroutine meshrectdom2d_generate
202
203 !- private ------------------------------------------------------
204
205 subroutine meshrectdom2d_setuplocaldom( lcmesh, &
206 tileID, panelID, &
207 i, j, NprcX, NprcY, &
208 dom_xmin, dom_xmax, dom_ymin, dom_ymax, &
209 NeX, NeY )
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
293 end subroutine meshrectdom2d_setuplocaldom
294
295 subroutine meshrectdom2d_assigndomid( this, &
296 tileID_table, panelID_table, &
297 pi_table, pj_table )
298
299 use scale_meshutil_2d, only: &
301 use scale_io
302 implicit none
303
304 class(meshrectdom2d), target, intent(inout) :: this
305 integer, intent(out) :: tileid_table(this%local_mesh_num, this%prc_num)
306 integer, intent(out) :: panelid_table(this%local_mesh_num*this%prc_num)
307 integer, intent(out) :: pi_table(this%local_mesh_num*this%prc_num)
308 integer, intent(out) :: pj_table(this%local_mesh_num*this%prc_num)
309
310 integer :: n
311 integer :: prc
312 integer :: tileid
313 integer :: is_lc, js_lc
314 integer :: ilc_count, jlc_count
315 integer :: ilc, jlc
316
317 type(localmesh2d), pointer :: lcmesh
318 !-----------------------------------------------------------------------------
319
321 panelid_table, pi_table, pj_table, & ! (out)
322 this%tileID_globalMap, this%tileFaceID_globalMap, this%tilePanelID_globalMap, & ! (out)
323 this%LOCAL_MESH_NUM_global, this%isPeriodicX, this%isPeriodicY, & ! (in)
324 this%NprcX, this%NprcY ) ! (in)
325
326 !----
327
328 do prc=1, this%PRC_NUM
329 do n=1, this%LOCAL_MESH_NUM
330 tileid = n + (prc-1)*this%LOCAL_MESH_NUM
331 lcmesh => this%lcmesh_list(n)
332 !-
333 tileid_table(n,prc) = tileid
334 this%tileID_global2localMap(tileid) = n
335 this%PRCRank_globalMap(tileid) = prc - 1
336
337 !-
338 if ( this%PRCRank_globalMap(tileid) == lcmesh%PRC_myrank ) then
339 if (n==1) then
340 is_lc = pi_table(tileid); ilc_count = 1
341 js_lc = pj_table(tileid); jlc_count = 1
342 end if
343 if(is_lc < pi_table(tileid)) ilc_count = ilc_count + 1
344 if(js_lc < pj_table(tileid)) jlc_count = jlc_count + 1
345 end if
346 end do
347 end do
348
349 allocate( this%rcdomIJ2LCMeshID(ilc_count,jlc_count) )
350 do jlc=1, jlc_count
351 do ilc=1, ilc_count
352 this%rcdomIJ2LCMeshID(ilc,jlc) = ilc + (jlc - 1)*ilc_count
353 end do
354 end do
355
356 return
357 end subroutine meshrectdom2d_assigndomid
358
359 subroutine meshrectdom2d_coord_conv( x, y, xr, xs, yr, ys, &
360 vx, vy, elem )
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
380 end subroutine meshrectdom2d_coord_conv
381
382 subroutine meshrectdom2d_calc_normal( normal_fn, &
383 Escale_f, fid, elem )
384
385 implicit none
386
387 type(elementbase2d), intent(in) :: elem
388 real(rp), intent(out) :: normal_fn(elem%nfptot,2)
389 integer, intent(in) :: fid(elem%nfp,elem%nfaces)
390 real(rp), intent(in) :: escale_f(elem%nfptot,2,2)
391
392 integer :: d
393 !-------------------------------------------------
394
395 do d=1, 2
396 normal_fn(fid(:,1),d) = - escale_f(fid(:,1),2,d)
397 normal_fn(fid(:,2),d) = + escale_f(fid(:,2),1,d)
398 normal_fn(fid(:,3),d) = + escale_f(fid(:,3),2,d)
399 normal_fn(fid(:,4),d) = - escale_f(fid(:,4),1,d)
400 end do
401
402 return
403 end subroutine meshrectdom2d_calc_normal
404
405end module scale_mesh_rectdom2d
module FElib / Element / Base
module FElib / Element / Quadrilateral
module FElib / Mesh / Local 2D
module FElib / Mesh / Local, Base
integer, parameter, public bctype_interior
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)
module FElib / Mesh / Rectangle 2D domain
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)
subroutine meshrectdom2d_init(this, negx, negy, dom_xmin, dom_xmax, dom_ymin, dom_ymax, isperiodicx, isperiodicy, refelem, nlocalmeshperprc, nprcx, nprcy, nproc, myrank)
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_buildglobalmap(panelid_table, pi_table, pj_table, tileid_map, tilefaceid_map, tilepanelid_map, ntile, isperiodicx, isperiodicy, ne_x, ne_y)
subroutine, public meshutil2d_genpatchboundarymap(vmapb, mapb, vmapp, pos_en, xmin, xmax, ymin, ymax, fmask, ne, np, nfp, nfaces, nv)