FE-Project
Loading...
Searching...
No Matches
scale_meshutil_3d.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
9!-------------------------------------------------------------------------------
10#include "scaleFElib.h"
12 !-----------------------------------------------------------------------------
13 !
14 !++ used modules
15 !
16 use scale_precision
17
18 !-----------------------------------------------------------------------------
19 implicit none
20 private
21
22 !-----------------------------------------------------------------------------
23 !
24 !++ Public procedure
25 !
31
32contains
33
34!OCL SERIAL
35 subroutine meshutil3d_gencubedomain( pos_v, EToV, &
36 Ke_x, xmin, xmax, Ke_y, ymin, ymax, Ke_z, zmin, zmax, &
37 Fz )
38
39 implicit none
40
41 integer, intent(in) :: ke_x
42 real(rp), intent(in) :: xmin, xmax
43 integer, intent(in) :: ke_y
44 real(rp), intent(in) :: ymin, ymax
45 integer, intent(in) :: ke_z
46 real(rp), intent(in) :: zmin, zmax
47 real(rp), intent(in), optional :: fz(ke_z+1)
48
49 real(rp), intent(out) :: pos_v((ke_x+1)*(ke_y+1)*(ke_z+1),3)
50 integer, intent(out) :: etov(ke_x*ke_y*ke_z,8)
51
52 integer :: i, j, k
53 integer :: n
54 integer :: k_
55 integer :: nvx, nvy, nvz
56 !-----------------------------------------------------------------------------
57
58 nvx = ke_x + 1
59 nvy = ke_y + 1
60 nvz = ke_z + 1
61
62 !$omp parallel private(i,j,k,n)
63
64 !$omp do collapse(2)
65 do k=1, nvz
66 do j=1, nvy
67 do i=1, nvx
68 n = i + (j-1)*nvx + (k-1)*nvx*nvy
69 pos_v(n,1) = (xmax - xmin)*dble(i - 1)/dble(nvx - 1) + xmin
70 pos_v(n,2) = (ymax - ymin)*dble(j - 1)/dble(nvy - 1) + ymin
71 if ( present(fz) ) then
72 pos_v(n,3) = fz(k)
73 else
74 pos_v(n,3) = (zmax - zmin)*dble(k - 1)/dble(nvz - 1) + zmin
75 end if
76 end do
77 end do
78 end do
79
80 !$omp do collapse(2)
81 do k=1, ke_z
82 do j=1, ke_y
83 do i=1, ke_x
84 n = i + (j-1)*ke_x + (k-1)*ke_x*ke_y
85 etov(n,1:4) = (k-1)*nvx*nvy + (j-1)*nvx + &
86 & (/ i, i+1, &
87 & i + nvx, i+1 + nvx /)
88 etov(n,5:8) = k *nvx*nvy + (j-1)*nvx + &
89 & (/ i, i+1, &
90 & i + nvx, i+1 + nvx /)
91 end do
92 end do
93 end do
94
95 !$omp end parallel
96
97 !---
98 !!$
99 !!$ write(*,*) "-- vx, vy --"
100 !!$ do j=1, NvY
101 !!$ do i=1, NvX
102 !!$ n = i + (j-1)*NvY
103 !!$ write(*,*) i, j, " :", mesh%vx(n), mesh%vy(n)
104 !!$ end do
105 !!$ end do
106 !!$
107 !!$ write(*,*) "-- EToV --"
108 !!$ do j=1, Ke_y
109 !!$ do i=1, Ke_x
110 !!$ n = i + (j-1)*Ke_y
111 !!$ write(*,*) i, j, " :", mesh%EToV(n,:)
112 !!$ end do
113 !!$ end do
114
115 return
116 end subroutine meshutil3d_gencubedomain
117
118!OCL SERIAL
119 subroutine meshutil3d_genconnectivity( EToE, EToF, &
120 EToV, Ne, Nfaces )
121
123 implicit none
124
125 integer, intent(in) :: ne
126 integer, intent(in) :: nfaces
127 integer, intent(out) :: etoe(ne, nfaces)
128 integer, intent(out) :: etof(ne, nfaces)
129 integer, intent(in) :: etov(ne,8)
130
131 integer :: nodes(ne*nfaces,4)
132 integer(kind=8) :: face_ids(ne*nfaces)
133 integer :: ke
134 integer :: f
135 integer :: n
136 integer :: n1, n2
137 integer :: nnodes
138 integer :: tmp
139 integer :: nnodes_row
140 integer :: matchl(2,3), matchr(2,3)
141
142 real(rp) :: etoe_1d(ne*nfaces)
143 real(rp) :: etof_1d(ne*nfaces)
144
145 integer :: spnodetonode(3,ne*nfaces)
146 integer :: spnodetonoderowtmp(3,ne*nfaces)
147 integer :: sort_indx(ne*nfaces)
148 integer(kind=8) :: sorted_faceid(ne*nfaces)
149
150 real(rp) :: vtmp(4)
151
152 !-----------------------------------------------------------------------------
153
154 nnodes = maxval( etov )
155 nnodes_row = size(nodes,1)
156
157 !---------
158 do ke=1, ne
159 nodes(ke ,:) = etov(ke,(/ 1, 2, 6, 5 /))
160 nodes(ke+ne ,:) = etov(ke,(/ 2, 4, 8, 6 /))
161 nodes(ke+2*ne,:) = etov(ke,(/ 4, 3, 7, 8 /))
162 nodes(ke+3*ne,:) = etov(ke,(/ 3, 1, 5, 7 /))
163 nodes(ke+4*ne,:) = etov(ke,(/ 1, 3, 4, 2 /))
164 nodes(ke+5*ne,:) = etov(ke,(/ 5, 6, 8, 7 /))
165 end do
166
167 ! Sort
168 do n=1, nnodes_row
169 vtmp(:) = nodes(n,:)
170 call bubblesort( vtmp )
171 nodes(n,:) = vtmp(:)
172 ! write(*,*) n, ":", nodes(n,:)
173 end do
174 nodes = nodes - 1
175 !---------
176
177 do ke=1, ne
178 etoe(ke,:) = ke
179 etof(ke,:) = (/ 1, 2, 3, 4, 5, 6 /)
180 end do
181
182 face_ids(:) = nodes(:,1)*nnodes**3 + nodes(:,2)*nnodes**2 &
183 + nodes(:,3)*nnodes + nodes(:,4) + 1
184
185 do f=1, nfaces
186 do ke=1, ne
187 n = ke + (f-1)*ne
188 spnodetonode(:,n) = (/ n, etoe(ke,f), etof(ke,f) /)
189 sorted_faceid(n) = face_ids(n)
190 sort_indx(n) = n
191 ! write(*,*) "face_id, n, EToE, EToF:", spNodeToNode(:,n)
192 end do
193 end do
194
195 !- sort row
196 call quicksort_exec_with_idx( ne*nfaces, sorted_faceid, sort_indx )
197 spnodetonoderowtmp(:,:) = spnodetonode(:,:)
198
199 do n=1, nnodes_row
200 spnodetonode(:,n) = spnodetonoderowtmp(:,sort_indx(n))
201 ! write(*,'(a,4i10)') "(sorted) face_id, n, EToE, EToF:", spNodeToNode(:,n)
202 end do
203
204 etoe_1d(:) = -1
205 etof_1d(:) = -1
206 do n=1, nnodes_row-1
207 if ( sorted_faceid(n) - sorted_faceid(n+1) == 0 ) then
208 matchl(:,:) = transpose( spnodetonode(:,(/ n, n+1 /)) )
209 matchr(:,:) = transpose( spnodetonode(:,(/ n+1, n /)) )
210
211 etoe_1d(matchl(:,1)) = matchr(:,2)
212 etof_1d(matchl(:,1)) = matchr(:,3)
213 end if
214 end do
215
216 do f=1, nfaces
217 do ke=1, ne
218 n = ke + (f-1)*ne
219 if ( etoe_1d(n) /= -1 ) then
220 etoe(ke,f) = etoe_1d(n)
221 etof(ke,f) = etof_1d(n)
222 end if
223 end do
224 end do
225
226 !-------------------
227
228! write(*,*) "EToE------"
229! do k=1, mesh%Ne
230! write(*,*) "k=", k, ":", mesh%EToE(k,:)
231! end do
232! write(*,*) "EToF------"
233! do k=1, mesh%Ne
234! write(*,*) "k=", k, ":", mesh%EToF(k,:)
235! end do
236
237 return
238 end subroutine meshutil3d_genconnectivity
239
240!OCL SERIAL
241 subroutine bubblesort( array )
242 implicit none
243 real(rp), intent(inout) :: array(:)
244
245 integer :: i, j, n
246 real(rp) :: t
247 !--------------------
248
249 n = size(array)
250 do i=1, n-1
251 do j=i+1,n
252 t = array(i)
253 if(t > array(j)) then
254 array(i) = array(j)
255 array(j) = t
256 end if
257 end do
258 end do
259
260 return
261 end subroutine bubblesort
262
263!OCL SERIAL
264 subroutine meshutil3d_buildinteriormap( VMapM, VMapP, MapM, MapP, &
265 pos_en, pos_ev, EtoE, EtoF, EtoV, Fmask_h, Fmask_v, &
266 Ne, Nv, Np, Nfp_h, Nfp_v, NfpTot, Nfaces_h, Nfaces_v, Nfaces)
267
268 implicit none
269
270 integer, intent(in) :: ne
271 integer, intent(in) :: nv
272 integer, intent(in) :: np
273 integer, intent(in) :: nfp_h
274 integer, intent(in) :: nfp_v
275 integer, intent(in) :: nfptot
276 integer, intent(in) :: nfaces_h
277 integer, intent(in) :: nfaces_v
278 integer, intent(in) :: nfaces
279
280 integer, intent(out) :: vmapm(nfptot,ne)
281 integer, intent(out) :: vmapp(nfptot,ne)
282 integer, intent(out) :: mapm(nfptot,ne)
283 integer, intent(out) :: mapp(nfptot,ne)
284
285 real(rp), intent(in) :: pos_en(np,ne,3)
286 real(rp), intent(in) :: pos_ev(nv,3)
287 integer, intent(in) :: etoe(ne,nfaces)
288 integer, intent(in) :: etof(ne,nfaces)
289 integer, intent(in) :: etov(ne,nv)
290 integer, intent(in) :: fmask_h(nfp_h,nfaces_h)
291 integer, intent(in) :: fmask_v(nfp_v,nfaces_v)
292
293 integer :: ke, ke1, ke2
294 integer :: f, f1, f2
295 integer :: p
296 integer :: n
297 integer :: i
298 integer :: idp, idm
299 integer :: v1, v2
300 integer :: nodeids(np,ne)
301 real(rp) :: r_h(nfp_h,nfp_h,3,2)
302 real(rp) :: r_v(nfp_v,nfp_v,3,2)
303 real(rp) :: dist_h(nfp_h,nfp_h)
304 real(rp) :: dist_v(nfp_v,nfp_v)
305 real(rp) :: x(np*ne), y(np*ne), z(np*ne)
306
307 integer :: vmapm_h(nfp_h,nfaces_h,ne)
308 integer :: vmapp_h(nfp_h,nfaces_h,ne)
309 integer :: vmapm_v(nfp_v,nfaces_v,ne)
310 integer :: vmapp_v(nfp_v,nfaces_v,ne)
311 integer :: mapm_h(nfp_h,nfaces_h,ne)
312 integer :: mapp_h(nfp_h,nfaces_h,ne)
313 integer :: mapm_v(nfp_v,nfaces_v,ne)
314 integer :: mapp_v(nfp_v,nfaces_v,ne)
315
316 integer :: mindist_indx(1)
317 !-----------------------------------------------------------------------------
318
319 !$omp parallel private(ke,f,p,n)
320
321 !$omp do
322 do ke=1, ne
323 do p=1, np
324 n = p + (ke-1)*np
325 nodeids(p,ke) = n
326 x(n) = pos_en(p,ke,1)
327 y(n) = pos_en(p,ke,2)
328 z(n) = pos_en(p,ke,3)
329 end do
330 end do
331 !$omp end do
332
333 !$omp do
334 do ke=1, ne
335 do f=1, nfaces_h
336 do p=1, nfp_h
337 n = p + (f-1)*nfp_h + (ke-1)*nfptot
338 mapm_h(p,f,ke) = n
339 mapp_h(p,f,ke) = n
340 vmapm_h(p,f,ke) = nodeids(fmask_h(p,f),ke)
341 end do
342 end do
343 do f=1, nfaces_v
344 do p=1, nfp_v
345 n = p + nfaces_h*nfp_h + (f-1)*nfp_v + (ke-1)*nfptot
346 mapm_v(p,f,ke) = n
347 mapp_v(p,f,ke) = n
348 vmapm_v(p,f,ke) = nodeids(fmask_v(p,f),ke)
349 end do
350 end do
351 end do
352 !$omp end do
353
354 !$omp workshare
355 vmapp_h(:,:,:) = -1
356 vmapp_v(:,:,:) = -1
357 !$omp end workshare
358 !$omp end parallel
359
360 !$omp parallel private( &
361 !$omp ke1, f1, ke2, f2, v1, v2, &
362 !$omp r_h, r_v, dist_h, dist_v, mindist_indx, idP, idM )
363
364 !$omp do
365 do ke1=1, ne
366 do f1=1, nfaces_h
367 ke2 = etoe(ke1,f1); f2 = etof(ke1,f1)
368
369 v1 = etov(ke1,f1); v2 = etov(ke1,1+mod(f1,nfaces_h))
370
371 r_h(:,:,1,1) = spread( x(vmapm_h(:,f1,ke1)), 2, nfp_h )
372 r_h(:,:,1,2) = spread( x(vmapm_h(:,f2,ke2)), 1, nfp_h )
373 r_h(:,:,2,1) = spread( y(vmapm_h(:,f1,ke1)), 2, nfp_h )
374 r_h(:,:,2,2) = spread( y(vmapm_h(:,f2,ke2)), 1, nfp_h )
375 r_h(:,:,3,1) = spread( z(vmapm_h(:,f1,ke1)), 2, nfp_h )
376 r_h(:,:,3,2) = spread( z(vmapm_h(:,f2,ke2)), 1, nfp_h )
377
378 dist_h(:,:) = (r_h(:,:,1,1) - r_h(:,:,1,2))**2 &
379 + (r_h(:,:,2,1) - r_h(:,:,2,2))**2 &
380 + (r_h(:,:,3,1) - r_h(:,:,3,2))**2
381 do idm=1, nfp_h
382 mindist_indx(:) = minloc(dist_h(idm,:))
383 idp = mindist_indx(1)
384 vmapp_h(idm,f1,ke1) = vmapm_h(idp,f2,ke2)
385 mapp_h(idm,f1,ke1) = idp + (f2-1)*nfp_h + (ke2-1)*nfptot
386 end do
387 end do
388 end do
389 !omp end do
390
391 !$omp do
392 do ke1=1, ne
393 do f1=1, nfaces_v
394 ke2 = etoe(ke1,nfaces_h+f1); f2 = etof(ke1,nfaces_h+f1) - nfaces_h
395
396 v1 = etov(ke1,1); v2 = etov(ke1,nfaces_h+1)
397
398 r_v(:,:,1,1) = spread( x(vmapm_v(:,f1,ke1)), 2, nfp_v )
399 r_v(:,:,1,2) = spread( x(vmapm_v(:,f2,ke2)), 1, nfp_v )
400 r_v(:,:,2,1) = spread( y(vmapm_v(:,f1,ke1)), 2, nfp_v )
401 r_v(:,:,2,2) = spread( y(vmapm_v(:,f2,ke2)), 1, nfp_v )
402 r_v(:,:,3,1) = spread( z(vmapm_v(:,f1,ke1)), 2, nfp_v )
403 r_v(:,:,3,2) = spread( z(vmapm_v(:,f2,ke2)), 1, nfp_v )
404
405 dist_v(:,:) = (r_v(:,:,1,1) - r_v(:,:,1,2))**2 &
406 + (r_v(:,:,2,1) - r_v(:,:,2,2))**2 &
407 + (r_v(:,:,3,1) - r_v(:,:,3,2))**2
408 do idm=1, nfp_v
409 mindist_indx(:) = minloc(dist_v(idm,:))
410 idp = mindist_indx(1)
411 vmapp_v(idm,f1,ke1) = vmapm_v(idp,f2,ke2)
412 mapp_v(idm,f1,ke1) = idp + nfaces_h*nfp_h + (f2-1)*nfp_v + (ke2-1)*nfptot
413 end do
414 end do
415 end do
416 !omp end do
417 !$omp end parallel
418
419 !$omp parallel do private(ke,f,n,i)
420 do ke=1, ne
421 do f=1, nfaces_h
422 do n=1, nfp_h
423 i = n + (f-1)*nfp_h
424 vmapm(i,ke) = vmapm_h(n,f,ke)
425 mapm(i,ke) = mapm_h(n,f,ke)
426 vmapp(i,ke) = vmapp_h(n,f,ke)
427 mapp(i,ke) = mapp_h(n,f,ke)
428 end do
429 end do
430 do f=1, nfaces_v
431 do n=1, nfp_v
432 i = n + nfaces_h*nfp_h + (f-1)*nfp_v
433 vmapm(i,ke) = vmapm_v(n,f,ke)
434 mapm(i,ke) = mapm_v(n,f,ke)
435 vmapp(i,ke) = vmapp_v(n,f,ke)
436 mapp(i,ke) = mapp_v(n,f,ke)
437 end do
438 end do
439 end do
440
441 !-----
442 ! mapB_counter = 0
443 ! do k=1,mesh%Ne
444 ! do f=1,elem%Nfaces
445 ! do p=1,elem%Nfp
446 ! n = p + (f-1)*elem%Nfp + (k-1)*elem%Nfp*elem%Nfaces
447 ! if (mesh%VMapM(p,f,k) == mesh%VMapP(p,f,k)) then
448 ! mapB_counter = mapB_counter + 1
449 ! mapB_tmp(mapB_counter) = n
450 !
451 ! vmapB_tmp(mapB_counter) = mesh%VMapM(p,f,k)
452 ! mesh%VMapP(p,f,k) = elem%Np*mesh%NeE + mapB_counter
453 ! end if
454 ! end do
455 ! end do
456 ! end do
457 !
458 ! allocate( mesh%mapB(mapB_counter) )
459 ! allocate( mesh%vmapB(mapB_counter) )
460 ! mesh%mapB(:) = mapB_tmp(1:mapB_counter)
461 ! mesh%vmapB(:) = vmapB_tmp(1:mapB_counter)
462
463 !-------
464 ! write(*,*) "Build MapInfo: "
465 ! do k=mesh%NeS,mesh%NeE
466 ! write(*,*) "k=", k, "---"
467 ! write(*,*) " VMapM:", mesh%VMapM(:,:,k)
468 ! write(*,*) " VMapP:", mesh%VMapP(:,:,k)
469 ! end do
470
471 ! write(*,*) "mapB:", mesh%mapB(:)
472 ! write(*,*) "vmapB:", mesh%vmapB(:)
473
474 return
475 end subroutine meshutil3d_buildinteriormap
476
477!OCL SERIAL
478 subroutine meshutil3d_genpatchboundarymap( VMapB, MapB, VMapP, &
479 pos_en, xmin, xmax, ymin, ymax, zmin, zmax, &
480 Fmask_h, Fmask_v, Ne, Nv, Np, Nfp_h, Nfp_v, NfpTot, Nfaces_h, Nfaces_v, Nfaces )
481
482 implicit none
483
484 integer, intent(in) :: ne
485 integer, intent(in) :: nv
486 integer, intent(in) :: np
487 integer, intent(in) :: nfp_h
488 integer, intent(in) :: nfp_v
489 integer, intent(in) :: nfptot
490 integer, intent(in) :: nfaces_h
491 integer, intent(in) :: nfaces_v
492 integer, intent(in) :: nfaces
493
494 integer, intent(inout), allocatable :: vmapb(:)
495 integer, intent(inout), allocatable :: mapb(:)
496 integer, intent(inout) :: vmapp(nfptot,ne)
497
498 real(rp), intent(in) :: pos_en(np,ne,3)
499 real(rp), intent(in) :: xmin, xmax
500 real(rp), intent(in) :: ymin, ymax
501 real(rp), intent(in) :: zmin, zmax
502 integer, intent(in) :: fmask_h(nfp_h,nfaces_h)
503 integer, intent(in) :: fmask_v(nfp_v,nfaces_v)
504
505 integer :: ke, n
506 integer :: b
507 integer :: f
508 integer :: i, j
509 real(rp) :: x, y, z
510
511 real(rp), parameter :: node_tol = 1.0e-12_rp
512
513 integer :: elemids_h(nfp_h*ne,nfaces_h)
514 real(rp) :: ordinfo_h(nfp_h*ne,nfaces_h)
515 integer :: faceids_h(nfp_h*ne,nfaces_h)
516 integer :: counterb_h(nfaces_h)
517
518 integer :: elemids_v(nfp_v*ne,nfaces_v)
519 real(rp) :: ordinfo_v(nfp_v*ne,nfaces_v)
520 integer :: faceids_v(nfp_v*ne,nfaces_v)
521 integer :: counterb_v(nfaces_v)
522
523 integer :: mapb_counter
524 real(rp) :: rdomx, rdomy, rdomz
525 !-----------------------------------------------------------------------------
526
527 counterb_h(:) = 0
528 counterb_v(:) = 0
529
530 rdomx = 1.0_rp/(xmax - xmin)
531 rdomy = 1.0_rp/(ymax - ymin)
532 rdomz = 1.0_rp/abs(zmax - zmin)
533
534 do ke=1, ne
535 do f=1, nfaces_h
536 x = sum(pos_en(fmask_h(:,f),ke,1)) / dble(nfp_h)
537 y = sum(pos_en(fmask_h(:,f),ke,2)) / dble(nfp_h)
538
539 call eval_domain_boundary( &
540 elemids_h, ordinfo_h, faceids_h, counterb_h, & ! (inout)
541 1, y, ymin, x, ke, f, rdomy ) ! (in)
542 call eval_domain_boundary( &
543 elemids_h, ordinfo_h, faceids_h, counterb_h, & ! (inout)
544 2, x, xmax, y, ke, f, rdomx ) ! (in)
545 call eval_domain_boundary( &
546 elemids_h, ordinfo_h, faceids_h, counterb_h, & ! (inout)
547 3, y, ymax, x, ke, f, rdomy ) ! (in)
548 call eval_domain_boundary( &
549 elemids_h, ordinfo_h, faceids_h, counterb_h, & ! (inout)
550 4, x, xmin, y, ke, f, rdomx ) ! (in)
551 end do
552 do f=1, nfaces_v
553 x = sum(pos_en(fmask_v(:,f),ke,1)) / dble(nfp_v)
554 z = sum(pos_en(fmask_v(:,f),ke,3)) / dble(nfp_v)
555
556 call eval_domain_boundary( &
557 elemids_v, ordinfo_v, faceids_v, counterb_v, & ! (inout)
558 1, z, zmin, x, ke, f, rdomz ) ! (in)
559 call eval_domain_boundary( &
560 elemids_v, ordinfo_v, faceids_v, counterb_v, & ! (inout)
561 2, z, zmax, x, ke, f, rdomz ) ! (in)
562 end do
563 end do
564
565
566 allocate( mapb(sum(counterb_h(:)*nfp_h)+sum(counterb_v(:)*nfp_v)) )
567 allocate( vmapb(size(mapb)) )
568
569 mapb_counter = 1
570 do b = 1, nfaces_h
571 ! write(*,*) "LocalMesh boundary ID:", b
572 ! write(*,*) counterB_h(b)
573 ! write(*,*) ordInfo_h(1:counterB_h(1),b)
574 ! write(*,*) elemIds_h(1:counterB_h(1),b)
575 ! write(*,*) faceIds_h(1:counterB_h(1),b)
576
577 do i=1, counterb_h(b)
578 ke = elemids_h(i,b); f = faceids_h(i,b)
579 do j=1, nfp_h
580 n = j + (f-1)*nfp_h
581 vmapp(n,ke) = np*ne + mapb_counter
582 vmapb(mapb_counter) = fmask_h(j,f) + (ke-1)*np
583 mapb_counter = mapb_counter + 1
584 end do
585 end do
586 end do
587
588 do b = 1, nfaces_v
589 do i=1, counterb_v(b)
590 ke = elemids_v(i,b); f = faceids_v(i,b)
591 do j=1, nfp_v
592 n = j + nfp_h*nfaces_h + nfp_v*(f-1)
593 vmapp(n,ke) = np*ne + mapb_counter
594 vmapb(mapb_counter) = fmask_v(j,f) + (ke-1)*np
595 mapb_counter = mapb_counter + 1
596 end do
597 end do
598 end do
599
600 ! write(*,*) "VMapP:-----"
601 ! do b=1, 4
602 ! do i=1, counterB(b)
603 ! k = elemIds(i,b); f = faceIDs(i,b)
604 ! write(*,*) "bid=", b, ":", mesh%VmapP(:,f,k)
605 ! end do
606 ! end do
607 ! write(*,*) "-----"
608 ! write(*,*) "VMapB:", mesh%VmapB(:)
609 ! write(*,*) "NeA=", mesh%NeA
610
611 return
612
613 contains
614!OCL SERIAL
615 subroutine eval_domain_boundary( &
616 elemIDs, ordInfo, faceIDs, counterB, &
617 domb_id, r, rbc, ord_info, k_, f_, normalized_fac )
618 implicit none
619
620 integer, intent(inout) :: elemids(:,:)
621 real(rp), intent(inout) :: ordinfo(:,:)
622 integer, intent(inout) :: counterb(:)
623 integer, intent(inout) :: faceids(:,:)
624 integer, intent(in) :: domb_id
625 real(rp), intent(in) :: r
626 real(rp), intent(in) :: rbc
627 real(rp), intent(in) :: ord_info
628 integer, intent(in) :: k_, f_
629 real(rp), intent(in) :: normalized_fac
630 !-------------------------------------------------------------
631
632 if ( abs(r - rbc)*normalized_fac < node_tol ) then
633 counterb(domb_id) = counterb(domb_id) + 1
634 ordinfo(counterb(domb_id),domb_id) = ord_info
635 elemids(counterb(domb_id),domb_id) = k_
636 faceids(counterb(domb_id),domb_id) = f_
637 end if
638
639 return
640 end subroutine eval_domain_boundary
641 end subroutine meshutil3d_genpatchboundarymap
642
643!OCL SERIAL
645 panelID_table, pi_table, pj_table, pk_table, &
646 tileID_map, tileFaceID_map, tilePanelID_map, &
647 Ntile, NtileFace, NtileVertex, &
648 isPeriodicX, isPeriodicY, isPeriodicZ, &
649 Ne_x, Ne_y, Ne_z )
650
651 ! use scale_prc, only: PRC_isMaster
652 implicit none
653
654 integer, intent(in) :: ntile
655 integer, intent(in) :: ntileface
656 integer, intent(in) :: ntilevertex
657 integer, intent(out) :: panelid_table(ntile)
658 integer, intent(out) :: pi_table(ntile)
659 integer, intent(out) :: pj_table(ntile)
660 integer, intent(out) :: pk_table(ntile)
661 integer, intent(out) :: tileid_map(ntileface,ntile)
662 integer, intent(out) :: tilefaceid_map(ntileface,ntile)
663 integer, intent(out) :: tilepanelid_map(ntileface,ntile)
664 logical, intent(in) :: isperiodicx
665 logical, intent(in) :: isperiodicy
666 logical, intent(in) :: isperiodicz
667 integer, intent(in) :: ne_x
668 integer, intent(in) :: ne_y
669 integer, intent(in) :: ne_z
670
671 integer :: ntileperpanel
672 integer :: nv_x, nv_y, nv_z
673 integer, allocatable :: nodesid_3d(:,:,:)
674 integer, allocatable :: etov(:,:)
675 integer, allocatable :: etoe(:,:)
676 integer, allocatable :: etof(:,:)
677 integer :: i, j, k, f
678 integer :: tileid, tileid_r
679 integer :: counter
680
681 !-----------------------------------------------------------------------------
682
683 ntileperpanel = ntile/1
684
685 nv_x = ne_x + 1
686 nv_y = ne_y + 1
687 nv_z = ne_z + 1
688
689 allocate( nodesid_3d(nv_x,nv_y,nv_z) )
690 allocate( etov(ntile,ntilevertex), etoe(ntile,ntileface), etof(ntile,ntileface) )
691
692 counter = 0
693 do k = 1, nv_z
694 do j = 1, nv_y
695 do i = 1, nv_x
696 counter = counter + 1
697 nodesid_3d(i,j,k) = counter
698 end do
699 end do
700 end do
701
702
703 !----
704
705 tileid = 0
706 do k = 1, ne_z
707 do j = 1, ne_y
708 do i = 1, ne_x
709 tileid = tileid + 1
710 panelid_table(tileid) = 1
711 pi_table(tileid) = i; pj_table(tileid) = j; pk_table(tileid) = k
712 etov(tileid,:) = (/ nodesid_3d(i,j ,k ), nodesid_3d(i+1,j ,k ), &
713 nodesid_3d(i,j+1,k ), nodesid_3d(i+1,j+1,k ), &
714 nodesid_3d(i,j ,k+1), nodesid_3d(i+1,j ,k+1), &
715 nodesid_3d(i,j+1,k+1), nodesid_3d(i+1,j+1,k+1) /)
716 end do
717 end do
718 end do
719
720 call meshutil3d_genconnectivity( etoe, etof, &
721 etov, ntile, ntileface )
722 tileid_map(:,:) = transpose(etoe)
723 tilefaceid_map(:,:) = transpose(etof)
724
725 do tileid=1, ntile
726 do f=1, ntileface
727 tileid_r = tileid_map(f,tileid)
728 tilepanelid_map(f,tileid) = panelid_table(tileid_r)
729 end do
730 end do
731
732 if (isperiodicx) then
733 do tileid=1, ntile
734 if (pi_table(tileid) == 1 .and. tilefaceid_map(4,tileid) == 4) then
735 tileid_map(4,tileid) = ne_x + (pj_table(tileid) - 1)*ne_x + (pk_table(tileid) - 1)*ne_x*ne_y
736 tilefaceid_map(4,tileid) = 2
737 end if
738 if (pi_table(tileid) == ne_x .and. tilefaceid_map(2,tileid) == 2) then
739 tileid_map(2,tileid) = 1 + (pj_table(tileid) - 1)*ne_x + (pk_table(tileid) - 1)*ne_x*ne_y
740 tilefaceid_map(2,tileid) = 4
741 end if
742 end do
743 end if
744
745 if (isperiodicy) then
746 do tileid=1, ntile
747 if (pj_table(tileid) == 1 .and. tilefaceid_map(1,tileid) == 1) then
748 tileid_map(1,tileid) = pi_table(tileid) + (ne_y - 1)*ne_x + (pk_table(tileid) - 1)*ne_x*ne_y
749 tilefaceid_map(1,tileid) = 3
750 end if
751 if (pj_table(tileid) == ne_y .and. tilefaceid_map(3,tileid) == 3) then
752 tileid_map(3,tileid) = pi_table(tileid) + (pk_table(tileid) - 1)*ne_x*ne_y
753 tilefaceid_map(3,tileid) = 1
754 end if
755 end do
756 end if
757
758 if (isperiodicz) then
759 do tileid=1, ntile
760 if (pk_table(tileid) == 1 .and. tilefaceid_map(5,tileid) == 5) then
761 tileid_map(5,tileid) = pi_table(tileid) + (pj_table(tileid) - 1)*ne_x + (ne_z - 1)*ne_x*ne_y
762 tilefaceid_map(5,tileid) = 6
763 end if
764 if (pk_table(tileid) == ne_z .and. tilefaceid_map(6,tileid) == 6) then
765 tileid_map(6,tileid) = pi_table(tileid) + (pj_table(tileid) - 1)*ne_x
766 tilefaceid_map(6,tileid) = 5
767 end if
768 end do
769 end if
770
771 return
772
773 !--
774 ! if (PRC_isMaster) then
775 ! write(*,*) "TotTile", Ntile
776 ! do tileID = 1, Ntile
777 ! write(*,*) "tileID:", tileID, ", EtoV:", EtoV(tileID,:)
778 ! end do
779 ! write(*,*) "-----------"
780 ! do tileID = 1, Ntile
781 ! write(*,*) "tileID:", tileID, ", EtoE:", EtoE(tileID,:)
782 ! end do
783 ! write(*,*) "-----------"
784 ! do tileID = 1, Ntile
785 ! write(*,*) "tileID:", tileID, ", EtoF:", EtoF(tileID,:)
786 ! end do
787 ! end if
788
789 end subroutine meshutil3d_buildglobalmap
790
791end module scale_meshutil_3d
module FElib / Mesh / utility for 3D mesh
subroutine, public meshutil3d_buildglobalmap(panelid_table, pi_table, pj_table, pk_table, tileid_map, tilefaceid_map, tilepanelid_map, ntile, ntileface, ntilevertex, isperiodicx, isperiodicy, isperiodicz, ne_x, ne_y, ne_z)
subroutine, public meshutil3d_buildinteriormap(vmapm, vmapp, mapm, mapp, pos_en, pos_ev, etoe, etof, etov, fmask_h, fmask_v, ne, nv, np, nfp_h, nfp_v, nfptot, nfaces_h, nfaces_v, nfaces)
subroutine, public meshutil3d_genconnectivity(etoe, etof, etov, ne, nfaces)
subroutine, public meshutil3d_gencubedomain(pos_v, etov, ke_x, xmin, xmax, ke_y, ymin, ymax, ke_z, zmin, zmax, fz)
subroutine, public meshutil3d_genpatchboundarymap(vmapb, mapb, vmapp, pos_en, xmin, xmax, ymin, ymax, zmin, zmax, fmask_h, fmask_v, ne, nv, np, nfp_h, nfp_v, nfptot, nfaces_h, nfaces_v, nfaces)
module common / sort algorithm