ADflow  v1.0
ADflow is a finite volume RANS solver tailored for gradient-based aerodynamic design optimization.
oversetPackingRoutines.F90
Go to the documentation of this file.
2 contains
3 
4  subroutine getoblockbuffersizes(il, jl, kl, iSize, rSize)
5 
6  ! Subroutine to get the required buffer sizes. This only uses the
7  ! block dimensions and technically doesn't have anything to do with
8  ! oBlock. This allows us to figure out the sizes and perform a
9  ! global communication before actually building the ADTrees, which
10  ! may not be that well load-balanced.
11 
12  use constants
13  implicit none
14 
15  ! Input/OUtput
16  integer(kind=intType) :: il, jl, kl
17  integer(kind=intType) :: rSize, iSize
18 
19  ! Working paramters
20  integer(kind=intType) :: ie, je, ke, nBBox, nLeaves
21 
22  ! Initializeation
23  isize = 0
24  rsize = 0
25 
26  ! Create ie sizes as well
27  ie = il + 1; je = jl + 1; ke = kl + 1
28 
29  ! Count up the integers we want to send:
30 
31  isize = isize + 6 ! Blocks sizes + proc + nn + cluster
32 
33  isize = isize + 8 * il * jl * kl ! hexa conn
34 
35  isize = isize + (ie + 2) * (je + 2) * (ke + 2) ! global cell
36 
37  isize = isize + il * jl * kl ! nearWall
38 
39  isize = isize + ie * je * ke ! invalidDonor
40 
41  ! Number of boxes in the ADT is the same as the number of elements
42  ! (like hexa conn (without the *8 obviously)
43  nbbox = il * jl * kl
44  nleaves = nbbox - 1 ! See ADT/adtBuild.f90
45 
46  isize = isize + nleaves * 2 ! Two ints for the children in each leaf
47 
48  ! Count up the reals we ned to send:
49  rsize = rsize + ie * je * ke ! qualDonor
50 
51  rsize = rsize + 3 * ie * je * ke ! xADT
52 
53  rsize = rsize + nbbox * 6 ! Cell bounding boxes
54 
55  rsize = rsize + nleaves * 12 ! Bounding boxes for leaves
56 
57  rsize = rsize + 1 ! Min block volume
58  end subroutine getoblockbuffersizes
59 
60  subroutine packoblock(oBlock)
61 
62  use constants
63  use oversetdata, only: oversetblock
64  implicit none
65 
66  ! Pack up everything we need for this block into its own buffer
67  ! inlucding the data required for the ADTree
68 
69  ! Input/Output Parameters
70  type(oversetblock), intent(inout) :: oBlock
71 
72  ! Working paramters
73  integer(kind=intType) :: rSize, iSize, i, j, k, nHexa, nADT
74  integer(kind=intType) :: ie, je, ke, il, jl, kl
75 
76  ! If the buffer is already allocated, the block is packed and there
77  ! is nothing to do
78  if (allocated(oblock%rBuffer)) then
79  return
80  end if
81 
82  call getoblockbuffersizes(oblock%il, oblock%jl, oblock%kl, isize, rsize)
83 
84  ! Allocate the buffers
85  allocate (oblock%rBuffer(rsize), oblock%iBuffer(isize))
86 
87  ! Reset the integer counter and add all the integers on this pass
88  isize = 0
89 
90  oblock%iBuffer(1) = oblock%il
91  oblock%iBuffer(2) = oblock%jl
92  oblock%iBuffer(3) = oblock%kl
93  oblock%iBuffer(4) = oblock%proc
94  oblock%iBuffer(5) = oblock%block
95  oblock%iBuffer(6) = oblock%cluster
96 
97  isize = isize + 6
98 
99  il = oblock%il
100  jl = oblock%jl
101  kl = oblock%kl
102 
103  ie = il + 1
104  je = jl + 1
105  ke = kl + 1
106 
107  nhexa = oblock%il * oblock%jl * oblock%kl
108  nadt = ie * je * ke
109 
110  do j = 1, nhexa
111  do i = 1, 8
112  isize = isize + 1
113  oblock%iBuffer(isize) = oblock%hexaConn(i, j)
114  end do
115  end do
116 
117  do k = 0, ke + 1
118  do j = 0, je + 1
119  do i = 0, ie + 1
120  isize = isize + 1
121  oblock%iBuffer(isize) = oblock%globalCell(i, j, k)
122  end do
123  end do
124  end do
125 
126  do k = 1, kl
127  do j = 1, jl
128  do i = 1, il
129  isize = isize + 1
130  oblock%iBuffer(isize) = oblock%nearWall(i, j, k)
131  end do
132  end do
133  end do
134 
135  do k = 1, ke
136  do j = 1, je
137  do i = 1, ie
138  isize = isize + 1
139  oblock%iBuffer(isize) = oblock%invalidDonor(i, j, k)
140  end do
141  end do
142  end do
143 
144  do i = 1, oblock%ADT%nLeaves
145  isize = isize + 1
146  oblock%iBuffer(isize) = oblock%ADT%ADTree(i)%children(1)
147  isize = isize + 1
148  oblock%iBuffer(isize) = oblock%ADT%ADTree(i)%children(2)
149  end do
150 
151  ! Reset the real counter and add all the real values on this pass.
152  rsize = 0
153 
154  do i = 1, ie * je * ke
155  rsize = rsize + 1
156  oblock%rBuffer(rsize) = oblock%qualDonor(1, i)
157  end do
158 
159  do j = 1, ie * je * ke
160  do i = 1, 3
161  rsize = rsize + 1
162  oblock%rBuffer(rsize) = oblock%xADT(i, j)
163  end do
164  end do
165 
166  do i = 1, oblock%ADT%nBboxes
167  oblock%rBuffer(rsize + 1:rsize + 6) = oblock%ADT%xBBox(:, i)
168  rsize = rsize + 6
169  end do
170 
171  do i = 1, oblock%ADT%nLeaves
172  oblock%rBuffer(rsize + 1:rsize + 6) = oblock%ADT%ADTree(i)%xMin(:)
173  rsize = rsize + 6
174 
175  oblock%rBuffer(rsize + 1:rsize + 6) = oblock%ADT%ADTree(i)%xMax(:)
176  rsize = rsize + 6
177  end do
178 
179  rsize = rsize + 1
180  oblock%rBuffer(rsize) = oblock%minVol
181 
182  end subroutine packoblock
183 
184  subroutine unpackoblock(oBlock)
185 
186  use constants
187  use oversetdata, only: oversetblock
188  implicit none
189 
190  ! unPack everything we need for this block from its own buffer
191  ! and reconstitute the data required for the ADTree. It is assumed
192  ! the buffers are already allocated and the data is available. This
193  ! does the exact *OPPOSITE* operation as the packBlock() routine
194 
195  ! Input/Output Parameters
196  type(oversetblock), intent(inout) :: oBlock
197 
198  ! Working paramters
199  integer(kind=intType) :: rSize, iSize, i, j, k, nHexa, nADT
200  integer(kind=intType) :: ie, je, ke, il, jl, kl
201 
202  ! Reset the integer counter and add all the integers on this pass
203  isize = 0
204 
205  oblock%il = oblock%iBuffer(1)
206  oblock%jl = oblock%iBuffer(2)
207  oblock%kl = oblock%iBuffer(3)
208  oblock%proc = oblock%iBuffer(4)
209  oblock%block = oblock%iBuffer(5)
210  oblock%cluster = oblock%iBuffer(6)
211  isize = isize + 6
212 
213  il = oblock%il
214  jl = oblock%jl
215  kl = oblock%kl
216  ie = il + 1
217  je = jl + 1
218  ke = kl + 1
219 
220  nhexa = oblock%il * oblock%jl * oblock%kl
221  nadt = ie * je * ke
222 
223  ! Allocate the remainder of the arrays in oBlock.
224  allocate (oblock%hexaConn(8, nhexa))
225  allocate (oblock%globalCell(0:ie + 1, 0:je + 1, 0:ke + 1))
226  allocate (oblock%nearWall(1:il, 1:jl, 1:kl))
227  allocate (oblock%invalidDonor(1:ie, 1:je, 1:ke))
228  allocate (oblock%qualDonor(1, ie * je * ke))
229  allocate (oblock%xADT(3, nadt))
230 
231  ! -------------------------------------------------------------------
232  ! Once we know the sizes, allocate all the arrays in the
233  ! ADTree. Since we are not going to call the *actual* build routine
234  ! for the ADT, we need to set all the information ourselves. This
235  ! essentially does the same thing as buildSerialHex.
236  oblock%ADT%adtType = adtvolumeadt
237  oblock%ADT%nNodes = nadt
238  oblock%ADT%nTetra = 0
239  oblock%ADT%nPyra = 0
240  oblock%ADT%nPrisms = 0
241  oblock%ADT%nTria = 0
242  oblock%ADT%nQuads = 0
243  oblock%ADT%coor => oblock%xADT
244  oblock%ADT%hexaConn => oblock%hexaConn
245  nullify (oblock%ADT%tetraConn, oblock%ADT%pyraConn, oblock%ADT%prismsConn)
246  oblock%ADT%nBBoxes = nhexa
247  allocate (oblock%ADT%xBBOX(6, nhexa))
248  allocate (oblock%ADT%elementType(nhexa))
249  allocate (oblock%ADT%elementID(nhexa))
250  oblock%ADT%comm = mpi_comm_self
251  oblock%ADT%nProcs = 1
252  oblock%ADT%myID = 0
253 
254  ! All hexas
255  oblock%ADT%elementType = adthexahedron
256 
257  do i = 1, nhexa
258  oblock%ADT%elementID(i) = i
259  end do
260 
261  oblock%ADT%nLeaves = oblock%ADT%nBBoxes - 1
262  if (oblock%ADT%nBBoxes <= 1) oblock%ADT%nLeaves = oblock%ADT%nLeaves + 1
263  allocate (oblock%ADT%ADTree(oblock%ADT%nLeaves))
264 
265  ! -------------------------------------------------------------------
266 
267  ! Now continue copying out the integer values
268  do i = 1, nhexa
269  do j = 1, 8
270  isize = isize + 1
271  oblock%hexaConn(j, i) = oblock%iBuffer(isize)
272  end do
273  end do
274 
275  do k = 0, ke + 1
276  do j = 0, je + 1
277  do i = 0, ie + 1
278  isize = isize + 1
279  oblock%globalCell(i, j, k) = oblock%iBuffer(isize)
280  end do
281  end do
282  end do
283 
284  do k = 1, kl
285  do j = 1, jl
286  do i = 1, il
287  isize = isize + 1
288  oblock%nearWall(i, j, k) = oblock%iBuffer(isize)
289  end do
290  end do
291  end do
292 
293  do k = 1, ke
294  do j = 1, je
295  do i = 1, ie
296  isize = isize + 1
297  oblock%invalidDonor(i, j, k) = oblock%iBuffer(isize)
298  end do
299  end do
300  end do
301 
302  do i = 1, oblock%ADT%nLeaves
303  isize = isize + 1
304  oblock%ADT%ADTree(i)%children(1) = oblock%iBuffer(isize)
305  isize = isize + 1
306  oblock%ADT%ADTree(i)%children(2) = oblock%iBuffer(isize)
307  end do
308 
309  ! Now copy out the real values
310  rsize = 0
311 
312  do i = 1, ie * je * ke
313  rsize = rsize + 1
314  oblock%qualDonor(1, i) = oblock%rBuffer(rsize)
315  end do
316 
317  do j = 1, ie * je * ke
318  do i = 1, 3
319  rsize = rsize + 1
320  oblock%xADT(i, j) = oblock%rBuffer(rsize)
321  end do
322  end do
323 
324  do i = 1, oblock%ADT%nBboxes
325  oblock%ADT%xBBox(:, i) = oblock%rBuffer(rsize + 1:rsize + 6)
326  rsize = rsize + 6
327  end do
328 
329  do i = 1, oblock%ADT%nLeaves
330  oblock%ADT%ADTree(i)%xMin(:) = oblock%rBuffer(rsize + 1:rsize + 6)
331  rsize = rsize + 6
332 
333  oblock%ADT%ADTree(i)%xMax(:) = oblock%rBuffer(rsize + 1:rsize + 6)
334  rsize = rsize + 6
335  end do
336 
337  rsize = rsize + 1
338  oblock%minVol = oblock%rBuffer(rsize)
339 
340  ! Flag this oBlock as being allocated:
341  oblock%allocated = .true.
342  deallocate (oblock%iBuffer, oblock%rBuffer)
343 
344  end subroutine unpackoblock
345 
346  subroutine getofringebuffersizes(il, jl, kl, iSize, rSize)
347 
348  ! Subroutine to get the required buffer sizes. This one is pretty
349  ! easy, but we use a routine to make it look the same as for the
350  ! oBlock.
351 
352  use constants
353  implicit none
354 
355  ! Input/OUtput
356  integer(kind=intType), intent(in) :: il, jl, kl
357  integer(kind=intType), intent(out) :: rSize, iSize
358 
359  ! Working
360  integer(kind=intType) :: mm
361 
362  ! All arrays have the same size
363  mm = (il - 1) * (jl - 1) * (kl - 1) ! nx*ny*nz
364 
365  ! Initializeation
366  isize = mm * 3 + 5 ! We need wallInd, isWall, myIndex plus 5 for the sizes
367  rsize = mm * 6 ! Need to send x and xSeed (3 each)
368 
369  end subroutine getofringebuffersizes
370 
371  subroutine packofringe(oFringe)
372 
373  use constants
374  use oversetdata, only: oversetfringe
375 
376  implicit none
377 
378  ! Pack up the search coordines in this oFringe into its own buffer
379  ! so we are ready to send it.
380 
381  ! Input/Output Parameters
382  type(oversetfringe), intent(inout) :: oFringe
383 
384  ! Working paramters
385  integer(kind=intType) :: rSize, iSize, mm, i, ii
386 
387  ! If the buffer is already allocated, the block is packed and there
388  ! is nothing to do
389  if (allocated(ofringe%rBuffer)) then
390  return
391  end if
392 
393  call getofringebuffersizes(ofringe%il, ofringe%jl, ofringe%kl, &
394  isize, rsize)
395 
396  ! Allocate the buffers
397  allocate (ofringe%rBuffer(rsize), ofringe%iBuffer(isize))
398 
399  mm = (ofringe%nx) * (ofringe%ny) * (ofringe%nz)
400 
401  ofringe%iBuffer(1) = ofringe%il
402  ofringe%iBuffer(2) = ofringe%jl
403  ofringe%iBuffer(3) = ofringe%kl
404  ofringe%iBuffer(4) = ofringe%cluster
405  ofringe%iBuffer(5) = ofringe%block
406  ii = 5
407 
408  ! Copy the integers. Just wallInd and isWall
409  do i = 1, mm
410 
411  ii = ii + 1
412  ofringe%iBuffer(ii) = ofringe%wallInd(i)
413 
414  ii = ii + 1
415  ofringe%iBuffer(ii) = ofringe%isWall(i)
416 
417  ii = ii + 1
418  ofringe%iBuffer(ii) = ofringe%fringeIntBuffer(5, i) ! myIndex
419 
420  end do
421 
422  ! Copy the reals. Reset the buffer here.
423  ii = 0
424  do i = 1, mm
425  ofringe%rBuffer(ii + 1) = ofringe%x(1, i)
426  ofringe%rBuffer(ii + 2) = ofringe%x(2, i)
427  ofringe%rBuffer(ii + 3) = ofringe%x(3, i)
428  ofringe%rBuffer(ii + 4) = ofringe%xSeed(1, i)
429  ofringe%rBuffer(ii + 5) = ofringe%xSeed(2, i)
430  ofringe%rBuffer(ii + 6) = ofringe%xSeed(3, i)
431  ii = ii + 6
432  end do
433 
434  end subroutine packofringe
435 
436  subroutine unpackofringe(oFringe)
437 
438  use constants
439  use oversetdata, only: oversetfringe
440  implicit none
441 
442  ! Pack up the search coordines in this oFringe into its own buffer
443  ! so we are ready to send it.
444 
445  ! Input/Output Parameters
446  type(oversetfringe), intent(inout) :: oFringe
447 
448  ! Working paramters
449  integer(kind=intType) :: rSize, iSize, idom, i, ii, mm
450 
451  ! Set the sizes of this oFringe
452  ofringe%il = ofringe%iBuffer(1)
453  ofringe%jl = ofringe%iBuffer(2)
454  ofringe%kl = ofringe%iBuffer(3)
455  ofringe%nx = ofringe%il - 1
456  ofringe%ny = ofringe%jl - 1
457  ofringe%nz = ofringe%kl - 1
458  ofringe%cluster = ofringe%iBuffer(4)
459  ofringe%block = ofringe%ibuffer(5)
460 
461  mm = (ofringe%nx) * (ofringe%ny) * (ofringe%nz)
462 
463  allocate ( &
464  ofringe%x(3, mm), &
465  ofringe%xSeed(3, mm), &
466  ofringe%wallInd(mm), &
467  ofringe%isWall(mm))
468 
469  ! Assume each cell will get just one donor. It's just a guess, it
470  ! will be expanded if necessary so the exact value doesn't matter.
471  allocate (ofringe%fringeIntBuffer(5, mm), ofringe%fringeRealBuffer(4, mm))
472  ofringe%nDonor = 0
473 
474  ii = 5 ! Already copied out the sizes
475 
476  ! Copy out integers
477  do i = 1, mm
478  ii = ii + 1
479  ofringe%wallInd(i) = ofringe%iBuffer(ii)
480 
481  ii = ii + 1
482  ofringe%isWall(i) = ofringe%iBuffer(ii)
483 
484  ii = ii + 1
485  ofringe%fringeIntBuffer(5, i) = ofringe%iBuffer(ii)
486 
487  ofringe%fringeIntBuffer(4, i) = ofringe%block
488  end do
489 
490  ! Copy the reals. Reset the counter ii counter here.
491  ii = 0
492  do i = 1, mm
493  ofringe%x(1, i) = ofringe%rBuffer(ii + 1)
494  ofringe%x(2, i) = ofringe%rBuffer(ii + 2)
495  ofringe%x(3, i) = ofringe%rBuffer(ii + 3)
496 
497  ofringe%xSeed(1, i) = ofringe%rBuffer(ii + 4)
498  ofringe%xSeed(2, i) = ofringe%rBuffer(ii + 5)
499  ofringe%xSeed(3, i) = ofringe%rBuffer(ii + 6)
500 
501  ii = ii + 6
502  end do
503 
504  ! Flag this oFringe as being allocated:
505  ofringe%allocated = .true.
506  deallocate (ofringe%rBuffer, ofringe%iBuffer)
507 
508  end subroutine unpackofringe
509 
510  subroutine getwallsize(famList, nNodes, nCells, dualMesh)
511  ! Simple helper routine to return the number of wall nodes and cells
512  ! for the block pointed to by blockPointers.
513 
514  use constants
515  use blockpointers, only: bctype, nbocos, bcdata
516  use sorting, only: faminlist
517  implicit none
518 
519  ! Input
520  integer(kind=intType), intent(in), dimension(:) :: famList
521  logical :: dualMesh
522 
523  ! Output
524  integer(kind=intType), intent(out) :: nNodes, nCells
525 
526  ! Working:
527  integer(kind=intType) :: mm, iBeg, iEnd, jBeg, jEnd
528 
529  ! Figure out the size the wall is going to be.
530  nnodes = 0
531  ncells = 0
532  do mm = 1, nbocos
533  faminclude: if (faminlist(bcdata(mm)%famID, famlist)) then
534  if (dualmesh) then
535  jbeg = bcdata(mm)%jnBeg - 1; jend = bcdata(mm)%jnEnd
536  ibeg = bcdata(mm)%inBeg - 1; iend = bcdata(mm)%inEnd
537  else
538  jbeg = bcdata(mm)%jnBeg; jend = bcdata(mm)%jnEnd
539  ibeg = bcdata(mm)%inBeg; iend = bcdata(mm)%inEnd
540  end if
541 
542  nnodes = nnodes + (iend - ibeg + 1) * (jend - jbeg + 1)
543  ncells = ncells + (iend - ibeg) * (jend - jbeg)
544  end if faminclude
545  end do
546 
547  end subroutine getwallsize
548 
549  subroutine getosurfbuffersizes(famList, il, jl, kl, iSize, rSize, dualMesh)
550 
551  ! Subroutine to get the required buffer sizes. This one is pretty
552  ! easy, but we use a routine to make it look the same as for hte
553  ! oBlock. Note that these bufer sizes are over-estimates. They are
554  ! the maximum possible amount of data to send.
555 
556  use constants
557  implicit none
558 
559  ! Input/OUtput
560  integer(kind=intType), intent(in), dimension(:) :: famList
561  integer(kind=intType), intent(in) :: il, jl, kl
562  logical, intent(in) :: dualMesh
563  integer(kind=intType), intent(out) :: rSize, iSize
564 
565  ! Working
566  integer(kind=intType) :: mm, nNodes, nCells, nBBox, nLeaves
567 
568  ! Initalization
569  isize = 3 ! For the block sizes
570  isize = isize + 4 ! For the maxCells/nCells/nNodes variables
571  rsize = 0
572 
573  call getwallsize(famlist, nnodes, ncells, dualmesh)
574 
575  ! Note that nCells here is the maximum number size. This will result
576  ! in a slight overestimate of the buffer size. This is ok.
577 
578  if (nnodes > 0) then
579  ! Count up the integers we want to send:
580 
581  isize = isize + ncells * 4 ! This is for the connectivity
582 
583  isize = isize + ncells ! This is for the iblank array
584 
585  isize = isize + ncells ! This is for the cellPtr array
586 
587  ! Number of boxes in the ADT is the same as the number of elements
588  nbbox = ncells
589  nleaves = nbbox - 1 ! See ADT/adtBuild.f90
590 
591  isize = isize + nleaves * 2 ! Two ints for the children in each leaf
592 
593  ! Count up the reals we ned to send:
594  rsize = rsize + 3 * nnodes ! surface coordinates
595 
596  rsize = rsize + nnodes ! surface delta
597 
598  rsize = rsize + nbbox * 6 ! Cell bounding boxes
599 
600  rsize = rsize + nleaves * 12 ! Bounding boxes for leaves
601  end if
602  end subroutine getosurfbuffersizes
603 
604  subroutine packosurf(famList, oSurf, dualMesh)
605 
606  use constants
607  use oversetdata, only: oversetwall
608 
609  implicit none
610 
611  ! Pack up the search coordines in this oSurf into its own buffer
612  ! so we are ready to send it.
613 
614  ! Input/Output Parameters
615  integer(kind=intType), intent(in), dimension(:) :: famList
616  type(oversetwall), intent(inout) :: oSurf
617  logical, intent(in) :: dualMesh
618  ! Working paramters
619  integer(kind=intType) :: rSize, iSize, mm, i, j, nNodes, nCells
620 
621  call getosurfbuffersizes(famlist, osurf%il, osurf%kl, osurf%kl, isize, rsize, dualmesh)
622 
623  ! Allocate the buffers
624  allocate (osurf%rBuffer(rsize), osurf%iBuffer(isize))
625 
626  osurf%iBuffer(1) = osurf%il
627  osurf%iBuffer(2) = osurf%jl
628  osurf%iBuffer(3) = osurf%kl
629  osurf%iBuffer(4) = osurf%nNodes
630  osurf%iBuffer(5) = osurf%nCells
631  osurf%iBuffer(6) = osurf%maxCells
632  osurf%iBuffer(7) = osurf%cluster
633 
634  if (osurf%nNodes > 0) then
635  isize = 7
636  do j = 1, osurf%nCells
637  do i = 1, 4
638  isize = isize + 1
639  osurf%iBuffer(isize) = osurf%conn(i, j)
640  end do
641  end do
642 
643  do i = 1, osurf%maxCells
644  isize = isize + 1
645  osurf%iBuffer(isize) = osurf%iBlank(i)
646  end do
647 
648  do i = 1, osurf%nCells
649  isize = isize + 1
650  osurf%iBuffer(isize) = osurf%cellPtr(i)
651  end do
652 
653  do i = 1, osurf%ADT%nLeaves
654  isize = isize + 1
655  osurf%iBuffer(isize) = osurf%ADT%ADTree(i)%children(1)
656  isize = isize + 1
657  osurf%iBuffer(isize) = osurf%ADT%ADTree(i)%children(2)
658  end do
659 
660  ! Done with the integer values, do the real ones
661  rsize = 0
662 
663  do i = 1, osurf%nNodes
664  do j = 1, 3
665  rsize = rsize + 1
666  osurf%rBuffer(rsize) = osurf%x(j, i)
667  end do
668  end do
669 
670  do i = 1, osurf%nNodes
671  rsize = rsize + 1
672  osurf%rBuffer(rsize) = osurf%delta(i)
673  end do
674 
675  do i = 1, osurf%ADT%nBboxes
676  osurf%rBuffer(rsize + 1:rsize + 6) = osurf%ADT%xBBox(:, i)
677  rsize = rsize + 6
678  end do
679 
680  do i = 1, osurf%ADT%nLeaves
681  osurf%rBuffer(rsize + 1:rsize + 6) = osurf%ADT%ADTree(i)%xMin(:)
682  rsize = rsize + 6
683 
684  osurf%rBuffer(rsize + 1:rsize + 6) = osurf%ADT%ADTree(i)%xMax(:)
685  rsize = rsize + 6
686  end do
687  end if
688  end subroutine packosurf
689 
690  subroutine unpackosurf(oSurf)
691 
692  use constants
693  use oversetdata, only: oversetwall
694  use kdtree2_module, only: kdtree2_create
695  implicit none
696 
697  ! Input/Output Parameters
698  type(oversetwall), intent(inout) :: oSurf
699 
700  ! Working paramters
701  integer(kind=intType) :: rSize, iSize, idom, i, j, k, n, iNode
702 
703  ! Set the sizes of this oSurf
704  osurf%il = osurf%iBuffer(1)
705  osurf%jl = osurf%iBuffer(2)
706  osurf%kl = osurf%iBuffer(3)
707 
708  osurf%nNodes = osurf%iBuffer(4)
709  osurf%nCells = osurf%iBuffer(5)
710  osurf%maxCells = osurf%iBuffer(6)
711  osurf%cluster = osurf%iBuffer(7)
712 
713  isize = 7
714  rsize = 0
715 
716  ! Allocate the arrays now that we know the sizes
717  allocate (osurf%x(3, osurf%nNodes))
718  allocate (osurf%delta(osurf%nNodes))
719  allocate (osurf%conn(4, osurf%nCells))
720  allocate (osurf%iBlank(osurf%maxCells))
721  allocate (osurf%cellPtr(osurf%nCells))
722  allocate (osurf%nte(4, osurf%nNodes))
723  ! Once we know the sizes, allocate all the arrays in the
724  ! ADTree. Since we are not going to call the *actual* build routine
725  ! for the ADT, we need to set all the information ourselves. This
726  ! essentially does the same thing as buildSerialHex.
727  osurf%ADT%adtType = adtsurfaceadt
728  osurf%ADT%nNodes = osurf%nNodes
729  osurf%ADT%nTetra = 0
730  osurf%ADT%nPyra = 0
731  osurf%ADT%nPrisms = 0
732  osurf%ADT%nTria = 0
733  osurf%ADT%nQuads = osurf%nCells
734  osurf%ADT%coor => osurf%x
735  osurf%ADT%quadsConn => osurf%conn
736  nullify (osurf%ADT%triaConn)
737  osurf%ADT%nBBoxes = osurf%nCells
738  allocate (osurf%ADT%xBBOX(6, osurf%nCells))
739  allocate (osurf%ADT%elementType(osurf%nCells))
740  allocate (osurf%ADT%elementID(osurf%nCells))
741  osurf%ADT%comm = mpi_comm_self
742  osurf%ADT%nProcs = 1
743  osurf%ADT%myID = 0
744 
745  ! All hexas
746  osurf%ADT%elementType = adtquadrilateral
747 
748  do i = 1, osurf%nCells
749  osurf%ADT%elementID(i) = i
750  end do
751 
752  osurf%ADT%nLeaves = osurf%ADT%nBBoxes - 1
753  if (osurf%ADT%nBBoxes <= 1) osurf%ADT%nLeaves = osurf%ADT%nLeaves + 1
754 
755  allocate (osurf%ADT%ADTree(osurf%ADT%nLeaves))
756 
757  ! Now continue copying out the values if necessary:
758  if (osurf%nNodes > 0) then
759  do j = 1, osurf%nCells
760  do i = 1, 4
761  isize = isize + 1
762  osurf%conn(i, j) = osurf%iBuffer(isize)
763  end do
764  end do
765 
766  do i = 1, osurf%maxCells
767  isize = isize + 1
768  osurf%iBlank(i) = osurf%iBuffer(isize)
769  end do
770 
771  do i = 1, osurf%nCells
772  isize = isize + 1
773  osurf%cellPtr(i) = osurf%iBuffer(isize)
774  end do
775 
776  do i = 1, osurf%ADT%nLeaves
777  isize = isize + 1
778  osurf%ADT%ADTree(i)%children(1) = osurf%iBuffer(isize)
779  isize = isize + 1
780  osurf%ADT%ADTree(i)%children(2) = osurf%iBuffer(isize)
781  end do
782 
783  ! Done with the integer values, do the real ones
784  rsize = 0
785 
786  do i = 1, osurf%nNodes
787  do j = 1, 3
788  rsize = rsize + 1
789  osurf%x(j, i) = osurf%rBuffer(rsize)
790  end do
791  end do
792 
793  do i = 1, osurf%nNodes
794  rsize = rsize + 1
795  osurf%delta(i) = osurf%rBuffer(rsize)
796  end do
797 
798  do i = 1, osurf%ADT%nBboxes
799  osurf%ADT%xBBox(:, i) = osurf%rBuffer(rsize + 1:rsize + 6)
800  rsize = rsize + 6
801  end do
802 
803  do i = 1, osurf%ADT%nLeaves
804  osurf%ADT%ADTree(i)%xMin(:) = osurf%rBuffer(rsize + 1:rsize + 6)
805  rsize = rsize + 6
806 
807  osurf%ADT%ADTree(i)%xMax(:) = osurf%rBuffer(rsize + 1:rsize + 6)
808  rsize = rsize + 6
809  end do
810  end if
811 
812  ! Build the KDTree
813  if (osurf%nNodes > 0) then
814  osurf%tree => kdtree2_create(osurf%x)
815  end if
816 
817  ! Build the inverse of the connectivity, the nodeToElem array.
818  osurf%nte = 0
819  do i = 1, osurf%nCells
820  do j = 1, 4
821  n = osurf%conn(j, i)
822  inner: do k = 1, 4
823  if (osurf%nte(k, n) == 0) then
824  osurf%nte(k, n) = i
825  exit inner
826  end if
827  end do inner
828  end do
829  end do
830 
831  ! Flag this oSurf as being allocated:
832  osurf%allocated = .true.
833  deallocate (osurf%rBuffer, osurf%iBuffer)
834  end subroutine unpackosurf
835 end module oversetpackingroutines
Definition: BCData.F90:1
integer(kind=inttype) nbocos
integer(kind=inttype), dimension(:), pointer bctype
integer, parameter adtvolumeadt
Definition: constants.F90:246
integer(kind=adtelementtype), parameter adthexahedron
Definition: constants.F90:254
integer(kind=adtelementtype), parameter adtquadrilateral
Definition: constants.F90:250
integer, parameter adtsurfaceadt
Definition: constants.F90:245
type(kdtree2) function, pointer, public kdtree2_create(input_data, dim, sort, rearrange)
Definition: kd_tree.f90:588
subroutine getoblockbuffersizes(il, jl, kl, iSize, rSize)
subroutine getosurfbuffersizes(famList, il, jl, kl, iSize, rSize, dualMesh)
subroutine getwallsize(famList, nNodes, nCells, dualMesh)
subroutine getofringebuffersizes(il, jl, kl, iSize, rSize)
subroutine unpackofringe(oFringe)
subroutine packosurf(famList, oSurf, dualMesh)
logical function faminlist(famID, famList)
Definition: sorting.F90:7