ADflow  v1.0
ADflow is a finite volume RANS solver tailored for gradient-based aerodynamic design optimization.
adtUtils.F90
Go to the documentation of this file.
1 module adtutils
2  !
3  ! Module, which contains small subroutines which perform useful
4  ! tasks.
5  !
6  use constants
7  use adtdata
8  implicit none
9  !
10  ! Variables stored in this module.
11  !
12  ! nStack: Number of elements allocated in the stack array;
13  ! needed for a more efficient implementation of the
14  ! local qsort routines.
15  ! stack(:): The corresponding array to store the stack.
16  ! This is a pointer, such that the reallocation
17  ! is easier.
18 
19  integer(kind=intType) :: nstack
20  integer(kind=intType), dimension(:), pointer :: stack
21 
22  !=================================================================
23 
24 contains
25 
26  !===============================================================
27 
28  subroutine adtterminate(ADT, routineName, errorMessage)
29  !
30  ! This routine writes the given error message to standard
31  ! output and terminates the executation of the program.
32  ! Subroutine intent(in) arguments.
33  ! --------------------------------
34  ! routineName: Name of the routine where the error occured.
35  ! ADT: Currently active ADT.
36  ! Subroutine intent(inout) arguments.
37  ! -----------------------------------
38  ! errorMessage: On input it contains the error message to be
39  ! written. It is modified in this routine, such
40  ! that it fits on one line. On output its
41  ! contents is undefined, which does not matter
42  ! a whole lot.
43  !
44  implicit none
45  !
46  ! Subroutine arguments
47  !
48  type(adttype), intent(in) :: ADT
49 
50  character(len=*), intent(in) :: routineName
51  character(len=*), intent(in) :: errorMessage
52  !
53  ! Local parameter
54  !
55  integer, parameter :: maxCharLine = 55
56  !
57  ! Local variables
58  !
59  integer :: ierr, len, i2
60  logical :: firstTime
61 
62  character(len=len_trim(errorMessage)) :: message
63  character(len=8) :: integerString
64 
65  ! Copy the errorMessage into message. It is not possible to work
66  ! with errorMessage directly, because it is modified in this
67  ! routine. Sometimes a constant string is passed to this routine
68  ! and some compilers simply fail then.
69 
70  message = errormessage
71 
72  ! Print a nice error message. In case of a parallel executable
73  ! also the processor ID is printed.
74 
75  print "(a)", "#"
76  print "(a)", "#=========================== !!! Error !!! &
77  &============================"
78 
79 #ifndef SEQUENTIAL_MODE
80  write (integerstring, "(i8)") adt%myID
81  integerstring = adjustl(integerstring)
82 
83  print "(2a)", "#* adtTerminate called by processor ", &
84  trim(integerstring)
85 #endif
86 
87  print "(2a)", "#* Run-time error in procedure ", &
88  trim(routinename)
89 
90  ! Loop to write the error message. If the message is too long it
91  ! is split over several lines.
92 
93  firsttime = .true.
94  do
95  ! Determine the remaining error message to be written.
96  ! If longer than the maximum number of characters allowed
97  ! on a line, it is attempted to split the message.
98 
99  message = adjustl(message)
100  len = len_trim(message)
101  i2 = min(maxcharline, len)
102 
103  if (i2 < len) i2 = index(message(:i2), " ", .true.) - 1
104  if (i2 < 0) i2 = index(message, " ") - 1
105  if (i2 < 0) i2 = len
106 
107  ! Write this part of the error message. If it is the first
108  ! line of the message some additional stuff is printed.
109 
110  if (firsttime) then
111  print "(2a)", "#* Error message: ", trim(message(:i2))
112  firsttime = .false.
113  else
114  print "(2a)", "#* ", trim(message(:i2))
115  end if
116 
117  ! Exit the loop if the entire message has been written.
118 
119  if (i2 == len) exit
120 
121  ! Adapt the string for the next part to be written.
122 
123  message = message(i2 + 1:)
124 
125  end do
126 
127  ! Write the trailing message.
128 
129  print "(a)", "#*"
130  print "(a)", "#* Now exiting"
131  print "(a)", "#==========================================&
132  &============================"
133  print "(a)", "#"
134 
135  ! Call abort and stop the program. This stop should be done in
136  ! abort, but just to be sure.
137 
138  call mpi_abort(adt%comm, 1, ierr)
139  stop
140 
141  end subroutine adtterminate
142 
143  subroutine allocateadts
144  !
145  ! This routine allocates the memory for the first ADT and is
146  ! only called when no other ADT's are present.
147  !
148  implicit none
149  !
150  ! Local variables.
151  !
152  integer :: ierr
153 
154  ! Allocate the memory for 1 ADT. Note that adtTerminate is not
155  ! called when the memory allocation fails. The reason is that
156  ! the processor ID for the current tree is used in this routine
157  ! and that value has not been set yet.
158 
159  allocate (adts(1), stat=ierr)
160  if (ierr /= 0) stop "Allocation failure for ADTs"
161 
162  ! Nullify the pointers of ADTs(1).
163 
164  nullify (adts(1)%coor)
165 
166  nullify (adts(1)%triaConn)
167  nullify (adts(1)%quadsConn)
168  nullify (adts(1)%tetraConn)
169  nullify (adts(1)%pyraConn)
170  nullify (adts(1)%prismsConn)
171  nullify (adts(1)%hexaConn)
172 
173  nullify (adts(1)%rootLeavesProcs)
174  nullify (adts(1)%rootBBoxes)
175 
176  nullify (adts(1)%elementType)
177  nullify (adts(1)%elementID)
178  nullify (adts(1)%xBBox)
179 
180  nullify (adts(1)%ADTree)
181 
182  end subroutine allocateadts
183 
184  subroutine deallocateadts(adtID)
185  !
186  ! This routine deallocates the memory for the given entry in
187  ! the array ADTs and it tries to reallocate ADTs itself
188  ! accordingly.
189  ! Subroutine intent(in) arguments.
190  ! --------------------------------
191  ! adtID: The entry in ADTs to be deallocated.
192  !
193  implicit none
194  !
195  ! Subroutine arguments.
196  !
197  character(len=*), intent(in) :: adtID
198  !
199  ! Local variables.
200  !
201  integer :: ierr
202 
203  integer(kind=intType) :: jj, nn, nAlloc, nAllocNew
204 
205  type(adttype), dimension(:), allocatable :: tmpADTs
206 
207  ! Determine the index in the array ADTs, which stores the given
208  ! ID. As the number of trees stored is limited, a linear search
209  ! algorithm is okay.
210 
211  if (allocated(adts)) then
212  nalloc = ubound(adts, 1)
213  else
214  nalloc = 0
215  end if
216 
217  do jj = 1, nalloc
218  if (adtid == adts(jj)%adtID) exit
219  end do
220 
221  ! Return immediately if the ID is not present.
222 
223  if (jj > nalloc) return
224 
225  ! Deallocate the data for this ADT entry. Note that the memory
226  ! for the nodal coordinates and the connectivity should not be
227  ! deallocated, because these pointers are just set to the given
228  ! input. The deallocation only takes place if the tree is active.
229 
230  if (adts(jj)%isActive) then
231 
232  deallocate (adts(jj)%rootLeavesProcs, adts(jj)%rootBBoxes, &
233  adts(jj)%elementType, adts(jj)%elementID, &
234  adts(jj)%xBBox, adts(jj)%ADTree, &
235  stat=ierr)
236  if (ierr /= 0) &
237  call adtterminate(adts(jj), "deallocateADTs", &
238  "Deallocation failure for the ADT data")
239  end if
240 
241  ! Make sure the ADT is inactive and nullify the pointers.
242 
243  adts(jj)%isActive = .false.
244 
245  nullify (adts(jj)%coor)
246 
247  nullify (adts(jj)%triaConn)
248  nullify (adts(jj)%quadsConn)
249  nullify (adts(jj)%tetraConn)
250  nullify (adts(jj)%pyraConn)
251  nullify (adts(jj)%prismsConn)
252  nullify (adts(jj)%hexaConn)
253 
254  nullify (adts(jj)%rootLeavesProcs)
255  nullify (adts(jj)%rootBBoxes)
256 
257  nullify (adts(jj)%elementType)
258  nullify (adts(jj)%elementID)
259  nullify (adts(jj)%xBBox)
260 
261  nullify (adts(jj)%ADTree)
262 
263  ! Determine the highest entry in ADTs which is still valid.
264 
265  do nn = nalloc, 1, -1
266  if (adts(nn)%isActive) exit
267  end do
268 
269  ! Determine the situation we are having here.
270 
271  if (nn == 0) then
272 
273  ! No active ADT's anymore. Deallocte the entire array.
274  ! Note that adtTerminate cannot be called when something
275  ! goes wrong.
276 
277  deallocate (adts, stat=ierr)
278  if (ierr /= 0) stop "Deallocation failure for ADTs"
279 
280  else if (nn < nalloc) then
281 
282  ! There are still some active ADT's, but the highest ones
283  ! are inactive. Therefore ADTs is reallocated. First allocate
284  ! the memory for tmpADTs to be able to retrieve the currently
285  ! stored data later on.
286 
287  nallocnew = nn
288  allocate (tmpadts(nallocnew), stat=ierr)
289  if (ierr /= 0) &
290  call adtterminate(adts(jj), "adtDeallocateADTs", &
291  "Memory allocation failure for tmpADTs")
292 
293  ! Copy the data from ADTs to tmpADTs.
294 
295  do nn = 1, nallocnew
296  tmpadts(nn) = adts(nn)
297  end do
298 
299  ! Deallocate and allocate the memory for ADTs. Note that
300  ! adtTerminate is not called when the memory allocation fails.
301  ! The reason is that the processor ID for the current tree is
302  ! used in this routine and that value may not be available
303  ! anymore.
304 
305  deallocate (adts, stat=ierr)
306  if (ierr /= 0) stop "Deallocation failure for ADTs"
307 
308  allocate (adts(nallocnew), stat=ierr)
309  if (ierr /= 0) stop "Allocation failure for ADTs"
310 
311  ! Copy the data back into ADTs and release the memory of
312  ! tmpADTs afterwards.
313 
314  do nn = 1, nallocnew
315  adts(nn) = tmpadts(nn)
316  end do
317 
318  deallocate (tmpadts, stat=ierr)
319  if (ierr /= 0) stop "Deallocation failure for tmpADTs"
320 
321  end if
322 
323  end subroutine deallocateadts
324 
325  subroutine qsortbboxes(arr, nn, ADT, dir)
326  !
327  ! This routine sorts the integer array arr, such that the
328  ! coordinate of the corresponding bounding box in the
329  ! direction dir is in increasing order. Note that the array to
330  ! store the stack is stored in this module. The reason is that
331  ! this routine is called quite often and in this way constant
332  ! allocation, deallocation and reallocation of stack is
333  ! avoided.
334  ! Subroutine intent(in) arguments.
335  ! --------------------------------
336  ! nn: Size of the array to be sorted.
337  ! ADT: The ADTfrom which the coordinate of
338  ! the bounding box must be taken.
339  ! dir: Index of the coordinate, which must be sorted.
340  ! Subroutine intent(inout) arguments.
341  ! -----------------------------------
342  ! arr(:): On input it contains the bounding box ID's which
343  ! must be sorted. On output these ID's are sorted,
344  ! such that the given coordinate is in increasing
345  ! order.
346  !
347  implicit none
348  !
349  ! Subroutine arguments.
350  !
351  type(adttype), intent(in) :: ADT
352  integer(kind=intType), intent(in) :: nn, dir
353 
354  integer(kind=intType), dimension(:), intent(inout) :: arr
355  !
356  ! Local parameters.
357  !
358  integer(kind=intType), parameter :: m = 7
359  !
360  ! Local variables.
361  !
362  integer(kind=intType) :: i, j, k, r, l, jStack
363  integer(kind=intType) :: a, tmp
364 
365  real(kind=realtype) :: ra
366  real(kind=realtype), dimension(:, :), pointer :: xbbox
367 
368  ! Set the pointer for the coordinates of the bounding boxes.
369 
370  xbbox => adt%xBBox
371 
372  ! Initialize the variables that control the sorting.
373 
374  jstack = 0
375  l = 1
376  r = nn
377 
378  ! Start of the algorithm.
379 
380  sortloop: do
381 
382  ! Check for the size of the subarray.
383 
384  testinsertion: if ((r - l) < m) then
385 
386  ! Perform the insertion sort.
387 
388  do j = (l + 1), r
389  a = arr(j)
390  ra = xbbox(dir, a)
391  do i = (j - 1), l, -1
392  if (xbbox(dir, arr(i)) <= ra) exit
393  arr(i + 1) = arr(i)
394  end do
395  arr(i + 1) = a
396  end do
397 
398  ! In case there are no more elements on the stack, exit from
399  ! the outermost do-loop. Algorithm has finished.
400 
401  if (jstack == 0) exit sortloop
402 
403  ! Pop stack and begin a new round of partitioning.
404 
405  r = stack(jstack)
406  l = stack(jstack - 1)
407  jstack = jstack - 2
408 
409  else testinsertion
410 
411  ! Subarray is larger than the threshold for a linear sort.
412  ! Choose median of left, center and right elements as
413  ! partitioning element a. Also rearrange so that
414  ! (l) <= (l+1) <= (r).
415 
416  k = (l + r) / 2
417  tmp = arr(k) ! Swap the elements
418  arr(k) = arr(l + 1) ! k and l+1.
419  arr(l + 1) = tmp
420 
421  if (xbbox(dir, arr(r)) < xbbox(dir, arr(l))) then
422  tmp = arr(l) ! Swap the elements
423  arr(l) = arr(r) ! r and l.
424  arr(r) = tmp
425  end if
426 
427  if (xbbox(dir, arr(r)) < xbbox(dir, arr(l + 1))) then
428  tmp = arr(l + 1) ! Swap the elements
429  arr(l + 1) = arr(r) ! r and l+1.
430  arr(r) = tmp
431  end if
432 
433  if (xbbox(dir, arr(l + 1)) < xbbox(dir, arr(l))) then
434  tmp = arr(l + 1) ! Swap the elements
435  arr(l + 1) = arr(l) ! l and l+1.
436  arr(l) = tmp
437  end if
438 
439  ! Initialize the pointers for partitioning.
440 
441  i = l + 1
442  j = r
443  a = arr(l + 1)
444  ra = xbbox(dir, a)
445 
446  ! The innermost loop.
447 
448  innerloop: do
449 
450  ! Scan up to find element >= a.
451  do
452  i = i + 1
453  if (ra <= xbbox(dir, arr(i))) exit
454  end do
455 
456  ! Scan down to find element <= a.
457  do
458  j = j - 1
459  if (xbbox(dir, arr(j)) <= ra) exit
460  end do
461 
462  ! Exit the loop in case the pointers i and j crossed.
463 
464  if (j < i) exit innerloop
465 
466  ! Swap the element i and j.
467 
468  tmp = arr(i)
469  arr(i) = arr(j)
470  arr(j) = tmp
471 
472  end do innerloop
473 
474  ! Swap the entries j and l+1. Remember that a equals
475  ! arr(l+1).
476 
477  arr(l + 1) = arr(j)
478  arr(j) = a
479 
480  ! Push pointers to larger subarray on stack; process smaller
481  ! subarray immediately. Check if enough memory is available.
482  ! If not reallocate it.
483 
484  jstack = jstack + 2
485 
486  if (jstack > nstack) call reallocplus(stack, nstack, 100, adt)
487 
488  if ((r - i + 1) >= (j - l)) then
489  stack(jstack) = r
490  r = j - 1
491  stack(jstack - 1) = j
492  else
493  stack(jstack) = j - 1
494  stack(jstack - 1) = l
495  l = j
496  end if
497 
498  end if testinsertion
499  end do sortloop
500 
501  ! Check in debug mode if the sort has been done correctly.
502 
503  if (debug) then
504  do i = 1, (nn - 1)
505  if (xbbox(dir, arr(i + 1)) < xbbox(dir, arr(i))) then
506  call adtterminate(adt, "qsortBBoxes", &
507  "Array is not sorted correctly")
508  end if
509  end do
510  end if
511 
512  end subroutine qsortbboxes
513 
514  subroutine qsortbboxtargets(arr, nn, ADT)
515  !
516  ! This routine sorts the given number of bounding box targets
517  ! in increasing order, based on the generalized < operator.
518  !
519  implicit none
520  !
521  ! Subroutine arguments
522  !
523  type(adttype), intent(in) :: adt
524  integer(kind=intType), intent(in) :: nn
525 
526  type(adtbboxtargettype), dimension(:), pointer :: arr
527  !
528  ! Local variables
529  !
530  integer(kind=intType), parameter :: m = 7
531 
532  integer(kind=intType) :: i, j, k, r, l, jStack
533 
534  type(adtbboxtargettype) :: a, tmp
535 
536  ! Initialize the variables that control the sorting.
537 
538  jstack = 0
539  l = 1
540  r = nn
541 
542  ! Start of the algorithm
543 
544  sortloop: do
545 
546  ! Check for the size of the subarray.
547 
548  testinsertion: if ((r - l) < m) then
549 
550  ! Perform insertion sort
551 
552  do j = l + 1, r
553  a = arr(j)
554  do i = (j - 1), l, -1
555  if (arr(i) <= a) exit
556  arr(i + 1) = arr(i)
557  end do
558  arr(i + 1) = a
559  end do
560 
561  ! In case there are no more elements on the stack, exit from
562  ! the outermost do-loop. Algorithm has finished.
563 
564  if (jstack == 0) exit sortloop
565 
566  ! Pop stack and begin a new round of partitioning.
567 
568  r = stack(jstack)
569  l = stack(jstack - 1)
570  jstack = jstack - 2
571 
572  else testinsertion
573 
574  ! Subarray is larger than the threshold for a linear sort.
575  ! Choose median of left, center and right elements as
576  ! partitioning element a. Also rearrange so that
577  ! (l) <= (l+1) <= (r).
578 
579  k = (l + r) / 2
580  tmp = arr(k) ! Swap the elements
581  arr(k) = arr(l + 1) ! k and l+1.
582  arr(l + 1) = tmp
583 
584  if (arr(r) < arr(l)) then
585  tmp = arr(l) ! Swap the elements
586  arr(l) = arr(r) ! r and l.
587  arr(r) = tmp
588  end if
589 
590  if (arr(r) < arr(l + 1)) then
591  tmp = arr(l + 1) ! Swap the elements
592  arr(l + 1) = arr(r) ! r and l+1.
593  arr(r) = tmp
594  end if
595 
596  if (arr(l + 1) < arr(l)) then
597  tmp = arr(l + 1) ! Swap the elements
598  arr(l + 1) = arr(l) ! l and l+1.
599  arr(l) = tmp
600  end if
601 
602  ! Initialize the pointers for partitioning.
603 
604  i = l + 1
605  j = r
606  a = arr(l + 1)
607 
608  ! The innermost loop
609 
610  innerloop: do
611 
612  ! Scan up to find element >= a.
613  do
614  i = i + 1
615  if (a <= arr(i)) exit
616  end do
617 
618  ! Scan down to find element <= a.
619  do
620  j = j - 1
621  if (arr(j) <= a) exit
622  end do
623 
624  ! Exit the loop in case the pointers i and j crossed.
625 
626  if (j < i) exit innerloop
627 
628  ! Swap the element i and j.
629 
630  tmp = arr(i)
631  arr(i) = arr(j)
632  arr(j) = tmp
633 
634  end do innerloop
635 
636  ! Swap the entries j and l+1. Remember that a equals
637  ! arr(l+1).
638 
639  arr(l + 1) = arr(j)
640  arr(j) = a
641 
642  ! Push pointers to larger subarray on stack; process smaller
643  ! subarray immediately. Check if enough memory is available.
644  ! If not reallocate it.
645 
646  jstack = jstack + 2
647 
648  if (jstack > nstack) call reallocplus(stack, nstack, 100, adt)
649 
650  if ((r - i + 1) >= (j - l)) then
651  stack(jstack) = r
652  r = j - 1
653  stack(jstack - 1) = j
654  else
655  stack(jstack) = j - 1
656  stack(jstack - 1) = l
657  l = j
658  end if
659 
660  end if testinsertion
661  end do sortloop
662 
663  ! Check in debug mode whether the array is really sorted.
664 
665  if (debug) then
666  do i = 1, (nn - 1)
667  if (arr(i + 1) < arr(i)) &
668  call adtterminate(adt, "qsortBBoxTargets", &
669  "Array is not sorted correctly")
670  end do
671  end if
672 
673  end subroutine qsortbboxtargets
674 
675  subroutine reallocateadts(adtID, jj)
676  !
677  ! This routine reallocates the memory for the ADTs array, such
678  ! that it is possible to store a new ADT. First it is tried to
679  ! find an empty spot in the currently allocated array. If this
680  ! is not present a true reallocation takes place.
681  ! Subroutine intent(in) arguments.
682  ! --------------------------------
683  ! adtID: The ID of the ADT.
684  ! Subroutine intent(out) arguments.
685  ! ---------------------------------
686  ! jj: The index in the array ADTs, where the new entry will be
687  ! stored.
688  !
689  implicit none
690  !
691  ! Subroutine arguments.
692  !
693  character(len=*), intent(in) :: adtID
694  integer(kind=intType), intent(out) :: jj
695  !
696  ! Local variables.
697  !
698  integer :: ierr
699 
700  integer(kind=intType) :: nn, nAlloc
701 
702  type(adttype), dimension(:), allocatable :: tmpADTs
703 
704  ! Determine the current size of ADTs and look for an empty spot
705  ! in the currently allocated array. Also check if an ADT with
706  ! the given ID is not already active. A linear search algorithm
707  ! is used, because the number of ADT's stored is limited.
708 
709  nalloc = ubound(adts, 1)
710  jj = nalloc + 1
711  do nn = 1, nalloc
712  if (adts(nn)%isActive) then
713  if (adtid == adts(nn)%adtID) exit
714  else if (jj > nalloc) then
715  jj = nn
716  end if
717  end do
718 
719  ! If the given ID corresponds to an already active tree,
720  ! terminate. To avoid a messy output only processor 0 prints
721  ! an error message while the other ones wait to get killed.
722 
723  if (nn <= nalloc) then
724  if (adts(nn)%myID == 0) &
725  call adtterminate(adts(nn), "reallocateADTs", &
726  "Given ID corresponds to an already &
727  &active ADT")
728  call mpi_barrier(adts(nn)%comm, ierr)
729  end if
730 
731  ! Check if a reallocate must be done.
732 
733  checkreallocate: if (jj > nalloc) then
734 
735  ! No empty spot present in ADTs. A true reallocation must be
736  ! performed. First allocate the memory for tmpADTs to be able
737  ! to retrieve the currently stored data later on. Note that
738  ! adtTerminate is not called when the memory allocation fails.
739  ! The reason is that the processor ID for the current tree is
740  ! used in this routine and that value has not been set yet.
741 
742  allocate (tmpadts(nalloc), stat=ierr)
743  if (ierr /= 0) stop "Allocation failure for tmpADTs"
744 
745  ! Copy the data from ADTs to tmpADTs.
746 
747  do nn = 1, nalloc
748  tmpadts(nn) = adts(nn)
749  end do
750 
751  ! Release the memory of ADTs and allocate it again with
752  ! increased size.
753 
754  deallocate (adts, stat=ierr)
755  if (ierr /= 0) stop "Deallocation failure for ADTs"
756 
757  allocate (adts(jj), stat=ierr)
758  if (ierr /= 0) stop "Allocation failure for ADTs"
759 
760  ! Copy the data back from tmpADTs.
761 
762  do nn = 1, nalloc
763  adts(nn) = tmpadts(nn)
764  end do
765 
766  ! Release the memory of tmpADTs.
767 
768  deallocate (tmpadts, stat=ierr)
769  if (ierr /= 0) stop "Deallocation failure for tmpADTs"
770 
771  ! Nullify the pointers of the new entry.
772 
773  nullify (adts(jj)%coor)
774 
775  nullify (adts(jj)%triaConn)
776  nullify (adts(jj)%quadsConn)
777  nullify (adts(jj)%tetraConn)
778  nullify (adts(jj)%pyraConn)
779  nullify (adts(jj)%prismsConn)
780  nullify (adts(jj)%hexaConn)
781 
782  nullify (adts(jj)%rootLeavesProcs)
783  nullify (adts(jj)%rootBBoxes)
784 
785  nullify (adts(jj)%elementType)
786  nullify (adts(jj)%elementID)
787  nullify (adts(jj)%xBBox)
788 
789  nullify (adts(jj)%ADTree)
790 
791  end if checkreallocate
792 
793  end subroutine reallocateadts
794 
795  subroutine reallocbboxtargettypeplus(arr, nSize, nInc, ADT)
796  !
797  ! This routine reallocates the memory of the given
798  ! adtBBoxTargetType pointer array.
799  ! Subroutine intent(in) arguments.
800  ! --------------------------------
801  ! ADT: Currently active ADT.
802  ! nInc: Increment of the size of the array.
803  ! Subroutine intent(inout) arguments.
804  ! -----------------------------------
805  ! nSize: On input it contains the size of the given array.
806  ! On output this value is incremented by nInc.
807  ! Subroutine pointer arguments.
808  ! -----------------------------
809  ! arr: Array to be reallocated.
810  !
811  implicit none
812  !
813  ! Subroutine arguments.
814  !
815  type(adttype), intent(in) :: ADT
816  integer, intent(in) :: nInc
817  integer(kind=intType), intent(inout) :: nSize
818 
819  type(adtbboxtargettype), dimension(:), pointer :: arr
820  !
821  ! Local variables.
822  !
823  integer :: ierr
824  integer(kind=intType) :: i, nOld
825 
826  type(adtbboxtargettype), dimension(:), pointer :: tmp
827 
828  ! Store the input value of nSize and set the pointer tmp to the
829  ! original array.
830 
831  nold = nsize
832  tmp => arr
833 
834  ! Allocate the new memory for the array.
835 
836  nsize = nsize + ninc
837  allocate (arr(nsize), stat=ierr)
838  if (ierr /= 0) &
839  call adtterminate(adt, "reallocBBoxTargetTypePlus", &
840  "Memory allocation failure for arr.")
841 
842  ! Copy the data from the original array into arr.
843 
844  nold = min(nold, nsize)
845  do i = 1, nold
846  arr(i) = tmp(i)
847  end do
848 
849  ! Release the memory of tmp, which points to the original
850  ! memory of the given array.
851 
852  deallocate (tmp, stat=ierr)
853  if (ierr /= 0) &
854  call adtterminate(adt, "reallocBBoxTargetTypePlus", &
855  "Deallocation failure for tmp.")
856 
857  end subroutine reallocbboxtargettypeplus
858 
859  subroutine reallocplus(arr, nSize, nInc, ADT)
860  !
861  ! This internal routine reallocates the memory of the given
862  ! pointer array.
863  ! Subroutine intent(in) arguments.
864  ! --------------------------------
865  ! ADT: Currently active ADT.
866  ! nInc: Increment of the size of the array.
867  ! Subroutine intent(inout) arguments.
868  ! -----------------------------------
869  ! nSize: On input it contains the size of the given array.
870  ! On output this value is incremented by nInc.
871  ! Subroutine pointer arguments.
872  ! -----------------------------
873  ! arr: Array to be reallocated.
874  !
875  implicit none
876  !
877  ! Subroutine arguments.
878  !
879  type(adttype), intent(in) :: ADT
880  integer, intent(in) :: nInc
881  integer(kind=intType), intent(inout) :: nSize
882 
883  integer(kind=intType), dimension(:), pointer :: arr
884  !
885  ! Local variables.
886  !
887  integer :: ierr
888  integer(kind=intType) :: i, nOld
889 
890  integer(kind=intType), dimension(:), pointer :: tmp
891 
892  ! Store the input value of nSize and set the pointer tmp to the
893  ! original array.
894 
895  nold = nsize
896  tmp => arr
897 
898  ! Allocate the new memory for the array.
899 
900  nsize = nsize + ninc
901  allocate (arr(nsize), stat=ierr)
902  if (ierr /= 0) &
903  call adtterminate(adt, "reallocPlus", &
904  "Memory allocation failure for arr.")
905 
906  ! Copy the data from the original array into arr.
907 
908  nold = min(nold, nsize)
909  do i = 1, nold
910  arr(i) = tmp(i)
911  end do
912 
913  ! Release the memory of tmp, which points to the original
914  ! memory of the given array.
915 
916  deallocate (tmp, stat=ierr)
917  if (ierr /= 0) &
918  call adtterminate(adt, "reallocPlus", &
919  "Deallocation failure for tmp.")
920 
921  end subroutine reallocplus
922 
923 end module adtutils
type(adttype), dimension(:), allocatable, target adts
Definition: adtData.F90:170
subroutine deallocateadts(adtID)
Definition: adtUtils.F90:185
subroutine qsortbboxes(arr, nn, ADT, dir)
Definition: adtUtils.F90:326
subroutine qsortbboxtargets(arr, nn, ADT)
Definition: adtUtils.F90:515
subroutine reallocbboxtargettypeplus(arr, nSize, nInc, ADT)
Definition: adtUtils.F90:796
subroutine reallocplus(arr, nSize, nInc, ADT)
Definition: adtUtils.F90:860
integer(kind=inttype) nstack
Definition: adtUtils.F90:19
subroutine reallocateadts(adtID, jj)
Definition: adtUtils.F90:676
integer(kind=inttype), dimension(:), pointer stack
Definition: adtUtils.F90:20
subroutine adtterminate(ADT, routineName, errorMessage)
Definition: adtUtils.F90:29
subroutine allocateadts
Definition: adtUtils.F90:144