ADflow  v1.0
ADflow is a finite volume RANS solver tailored for gradient-based aerodynamic design optimization.
preprocessingModules.F90
Go to the documentation of this file.
2  !
3  ! This local module contains the derived data type used to
4  ! determine the indirect halo's as well as an array of this type.
5  !
6  use precision
7  implicit none
8  save
9 
10  public
11  private :: lessequalindirecthalotype
12  private :: lessindirecthalotype
13  !
14  ! The definition of the derived data type indirectHaloType.
15  !
17 
18  ! myBlock : Local block ID of the halo.
19  ! myI, myJ, myK: i,j,k indices of the halo.
20  ! myDirectHalo : Index in the haloListType where the
21  ! corresponding direct halo is stored.
22  ! levOfInd : Level of indirectness.
23  ! donorProc : Processor where donor of the direct halo is
24  ! stored. In case this halo is a boundary
25  ! halo, donorProc is set to -1.
26 
27  integer(kind=intType) :: myblock
28  integer(kind=intType) :: myi, myj, myk
29  integer(kind=intType) :: mydirecthalo
30  integer(kind=intType) :: levofind
31  integer(kind=intType) :: donorproc
32 
33  end type indirecthalotype
34 
35  ! Interface for the extension of the operators <= and <.
36  ! These are needed for the sorting of indirectHaloType.
37 
38  interface operator(<=)
39  module procedure lessequalindirecthalotype
40  end interface operator(<=)
41 
42  interface operator(<)
43  module procedure lessindirecthalotype
44  end interface operator(<)
45 
46  ! nIndHalo : Number of indirect halo's to be treated.
47  ! indHalo(nIndHalo): The indirect halo's.
48 
49  integer(kind=intType) :: nindhalo
50  type(indirecthalotype), dimension(:), allocatable :: indhalo
51 
52  ! nLevOfInd : Number of levels of indirectness
53  ! nHaloPerLev(0:nLevOfInd): Number of indirect halo's per level
54  ! of indirectness; stored in
55  ! cumulative storage format.
56  ! nHaloPerProc(0:nProc) : Number of indirect halo's per
57  ! processor for a given level of
58  ! indirectness; cumulative storage.
59  ! nHaloPerProc(0) is not 0,
60  ! because of the presence of boundary
61  ! halo's, which get proc ID -1.
62 
63  integer(kind=intType) :: nlevofind
64  integer(kind=intType), dimension(:), allocatable :: nhaloperlev
65  integer(kind=intType), dimension(:), allocatable :: nhaloperproc
66 
67 contains
68  !
69  ! Functions to simulate the operators <= and <.
70  !
71  logical function lessequalindirecthalotype(g1, g2)
72  !
73  ! This function returns .true. if g1 <= g2 and .false.
74  ! otherwise. The comparison is firstly based on the level of
75  ! indirectness followed by the donor processor, the
76  ! corresponding direct halo, my block ID and finally the i, j
77  ! and k indices.
78  !
79  implicit none
80  !
81  ! Function arguments.
82  !
83  type(indirecthalotype), intent(in) :: g1, g2
84 
85  ! Compare the level of indirectness. If not equal, set
86  ! lessEqual appropriately and return.
87 
88  if (g1%levOfInd < g2%levOfInd) then
89  lessequalindirecthalotype = .true.
90  return
91  else if (g1%levOfInd > g2%levOfInd) then
92  lessequalindirecthalotype = .false.
93  return
94  end if
95 
96  ! Compare the donor processors.
97 
98  if (g1%donorProc < g2%donorProc) then
99  lessequalindirecthalotype = .true.
100  return
101  else if (g1%donorProc > g2%donorProc) then
102  lessequalindirecthalotype = .false.
103  return
104  end if
105 
106  ! Compare the direct halo.
107 
108  if (g1%myDirectHalo < g2%myDirectHalo) then
109  lessequalindirecthalotype = .true.
110  return
111  else if (g1%myDirectHalo > g2%myDirectHalo) then
112  lessequalindirecthalotype = .false.
113  return
114  end if
115 
116  ! Compare my block ID.
117 
118  if (g1%myBlock < g2%myBlock) then
119  lessequalindirecthalotype = .true.
120  return
121  else if (g1%myBlock > g2%myBlock) then
122  lessequalindirecthalotype = .false.
123  return
124  end if
125 
126  ! Finally compare the halo indices. Start with k.
127 
128  if (g1%myK < g2%myK) then
129  lessequalindirecthalotype = .true.
130  return
131  else if (g1%myK > g2%myK) then
132  lessequalindirecthalotype = .false.
133  return
134  end if
135 
136  ! The j index.
137 
138  if (g1%myJ < g2%myJ) then
139  lessequalindirecthalotype = .true.
140  return
141  else if (g1%myJ > g2%myJ) then
142  lessequalindirecthalotype = .false.
143  return
144  end if
145 
146  ! The i index.
147 
148  if (g1%myI < g2%myI) then
149  lessequalindirecthalotype = .true.
150  return
151  else if (g1%myI > g2%myI) then
152  lessequalindirecthalotype = .false.
153  return
154  end if
155 
156  ! Both entities are identical. So set lessEqual to .true.
157 
158  lessequalindirecthalotype = .true.
159 
160  end function lessequalindirecthalotype
161 
162  ! ================================================================
163 
164  logical function lessindirecthalotype(g1, g2)
165  !
166  ! This function returns .true. If g1 < g2 and .false.
167  ! otherwise. It is basically the same as the lessEqual
168  ! function, except that the equality is now considered as
169  ! .false.
170  !
171  implicit none
172  !
173  ! Function arguments.
174  !
175  type(indirecthalotype), intent(in) :: g1, g2
176 
177  ! Compare the level of indirectness. If not equal, set
178  ! lessIndirectHaloType appropriately and return.
179 
180  if (g1%levOfInd < g2%levOfInd) then
181  lessindirecthalotype = .true.
182  return
183  else if (g1%levOfInd > g2%levOfInd) then
184  lessindirecthalotype = .false.
185  return
186  end if
187 
188  ! Compare the donor processors.
189 
190  if (g1%donorProc < g2%donorProc) then
191  lessindirecthalotype = .true.
192  return
193  else if (g1%donorProc > g2%donorProc) then
194  lessindirecthalotype = .false.
195  return
196  end if
197 
198  ! Compare the direct halo.
199 
200  if (g1%myDirectHalo < g2%myDirectHalo) then
201  lessindirecthalotype = .true.
202  return
203  else if (g1%myDirectHalo > g2%myDirectHalo) then
204  lessindirecthalotype = .false.
205  return
206  end if
207 
208  ! Compare my block id.
209 
210  if (g1%myBlock < g2%myBlock) then
211  lessindirecthalotype = .true.
212  return
213  else if (g1%myBlock > g2%myBlock) then
214  lessindirecthalotype = .false.
215  return
216  end if
217 
218  ! Finally compare the halo indices. Start with k.
219 
220  if (g1%myK < g2%myK) then
221  lessindirecthalotype = .true.
222  return
223  else if (g1%myK > g2%myK) then
224  lessindirecthalotype = .false.
225  return
226  end if
227 
228  ! The j index.
229 
230  if (g1%myJ < g2%myJ) then
231  lessindirecthalotype = .true.
232  return
233  else if (g1%myJ > g2%myJ) then
234  lessindirecthalotype = .false.
235  return
236  end if
237 
238  ! The i index.
239 
240  if (g1%myI < g2%myI) then
241  lessindirecthalotype = .true.
242  return
243  else if (g1%myI > g2%myI) then
244  lessindirecthalotype = .false.
245  return
246  end if
247 
248  ! Both entities are identical.
249  ! So set lessIndirectHaloType to .false.
250 
251  lessindirecthalotype = .false.
252 
253  end function lessindirecthalotype
254 
255 end module indirecthalo
256 module halolist
257  !
258  ! This local module contains temporary variables to create the
259  ! list of halo cells and nodes.
260  !
261  use precision
262  implicit none
263  save
264 
265  public
266  private :: lessequalhalolisttype
267  private :: lesshalolisttype
268  !
269  ! The definition of the variables for the 3 lists.
270  !
272 
273  ! myBlock: local block ID of the halo.
274  ! myI, myJ, myK: i,j,k indices of the halo.
275  ! donorProc : processor where donor is stored. In case
276  ! the halo is a boundary halo, donorProc
277  ! is set to -1.
278  ! donorBlock: block ID of the donor. In case the halo
279  ! is a boundary halo donorBlock is set
280  ! to the corresponding boundary condition.
281  ! dI, dJ, dK: i,j,k indices of the donor.
282  ! levOfInd: level of indirectness.
283  ! interp(..): interpolants for the donor stencil; only
284  ! allocated for lists requiring this info.
285  ! nPeriodicSubfaces: Number of periodic subfaces that are
286  ! crossed when going from the halo to the
287  ! donor. This is at most the level of
288  ! indirectness of the halo.
289  ! periodicSubfaces: The corresponding subfaces ID's according
290  ! to the sequence defined in periodicGlobal.
291 
292  integer(kind=intType) :: myblock
293  integer(kind=intType) :: myi, myj, myk
294  integer(kind=intType) :: donorproc, donorblock
295  integer(kind=intType) :: di, dj, dk
296  integer(kind=intType) :: levofind
297  real(kind=realtype), dimension(:), pointer :: interp
298 
299  integer(kind=intType) :: nperiodicsubfaces
300  integer(kind=intType), dimension(:), pointer :: periodicsubfaces
301 
302  end type halolisttype
303 
304  ! Interface for the extension of the operators <= and <.
305  ! These are needed for the sorting of haloListType.
306 
307  interface operator(<=)
308  module procedure lessequalhalolisttype
309  end interface operator(<=)
310 
311  interface operator(<)
312  module procedure lesshalolisttype
313  end interface operator(<)
314 
315  ! nCellHalo1st: # of 1st level cell halo's
316  ! nCellHalo2nd: # of 2nd level cell halo's
317  ! nNodeHalo1st: # of 1st level node halo's
318 
319  integer(kind=intType) :: ncellhalo1st, ncellhalo2nd
320  integer(kind=intType) :: nnodehalo1st
321 
322  ! iiCell1st: Counter variable for the 1st level cell halo's
323  ! iiCell2nd: Counter variable for the 2nd level cell halo's
324  ! iiNode1st: Counter variable for the 1st level node halo's
325 
326  integer(kind=intType) :: iicell1st, iicell2nd, iinode1st
327 
328  ! cellHalo1st(nCellHalo1st) :: List of halo info for 1st
329  ! level cell halo's.
330  ! cellHalo2nd(nCellHalo2nd) :: Idem for 2nd level cell halo's.
331  ! nodeHalo1st(nNodeHalo1st) :: Idem for 1st level node halo's.
332 
333  type(halolisttype), dimension(:), allocatable :: cellhalo1st
334  type(halolisttype), dimension(:), allocatable :: cellhalo2nd
335  type(halolisttype), dimension(:), allocatable :: nodehalo1st
336 
337  ! transformCell(nCellHalo1st,3) :: Short hand for the transformation
338  ! matrix between the halo and
339  ! the donor for cell based halo's.
340  ! In principle the size equals the
341  ! number of faces (i.e. direct)
342  ! halo's, but the difference is
343  ! not so large.
344  ! transformNode(nNodeHalo1st,3) :: Idem for the nodes.
345 
346  integer(kind=intType), dimension(:, :), allocatable :: transformcell
347  integer(kind=intType), dimension(:, :), allocatable :: transformnode
348  !
349  ! The definition of the index variables, which store for each
350  ! i,j,k in the block the index in the corresponding list.
351  ! I know I'm wasting memory here (because only the halo's are
352  ! relevant), but that's not too much of a problem. The reason is
353  ! that neither the metrics nor the variables have been allocated
354  ! yet. So later on, much more memory is needed than the single
355  ! integer for each cell/node used here.
356  !
358 
359  ! entryList(:,:,:) :: Corresponding entry in the list.
360  ! Dimensions are either 0:ie,0:je,0:ke for
361  ! the node or 0:ib,0:jb,0:kb for the cell
362  ! based halo's. The latter is then suited
363  ! for the 2nd level halo's.
364 
365  integer(kind=intType), dimension(:, :, :), pointer :: entrylist
366 
367  end type indexlisttype
368 
369  ! nodeIndex(nDom) :: The node indices for every block.
370  ! cellIndex(nDom) :: Idem for the cells.
371 
372  type(indexlisttype), allocatable, dimension(:) :: nodeindex
373  type(indexlisttype), allocatable, dimension(:) :: cellindex
374 
375 contains
376  !
377  ! Functions to simulate the operators <= and <.
378  !
379  logical function lessequalhalolisttype(g1, g2)
380  !
381  ! lessEqual returns .true. if g1 <= g2 and .false. otherwise.
382  ! The comparison is firstly based on the processor ID of the
383  ! donor. After that it depends whether the halo is a boundary
384  ! halo or not. Note that boundary halo's have a donor processor
385  ! if of -1, such that they are always first in the list.
386  !
387  implicit none
388  !
389  ! Function arguments.
390  !
391  type(halolisttype), intent(in) :: g1, g2
392 
393  ! Compare the donor processors first. If not equal,
394  ! set lessEqual appropriately and return.
395 
396  if (g1%donorProc < g2%donorProc) then
397  lessequalhalolisttype = .true.
398  return
399  else if (g1%donorProc > g2%donorProc) then
400  lessequalhalolisttype = .false.
401  return
402  end if
403 
404  ! Donor processors are identical. Now it depends whether we are
405  ! dealing with boundary halo's or not.
406 
407  boundary: if (g1%donorProc == -1) then ! And thus
408  ! g2%donorProc == -1
409 
410  ! Both halo's are boundary halo's. Compare the block ID of
411  ! the halo's.
412 
413  if (g1%myBlock < g2%myBlock) then
414  lessequalhalolisttype = .true.
415  return
416  else if (g1%myBlock > g2%myBlock) then
417  lessequalhalolisttype = .false.
418  return
419  end if
420 
421  ! Compare the boundary conditions, which are stored in
422  ! donorBlock. Note that the sequence in BCTypes is such that
423  ! the most important BC has the highest number.
424 
425  if (g1%donorBlock < g2%donorBlock) then
426  lessequalhalolisttype = .true.
427  return
428  else if (g1%donorBlock > g2%donorBlock) then
429  lessequalhalolisttype = .false.
430  return
431  end if
432 
433  ! As it is possible that indirect halo's need donor info from
434  ! direct halo's or even indirect halo's with a smaller level
435  ! of indirectness, compare the level of indirectness.
436 
437  if (g1%levOfInd < g2%levOfInd) then
438  lessequalhalolisttype = .true.
439  return
440  else if (g1%levOfInd > g2%levOfInd) then
441  lessequalhalolisttype = .false.
442  return
443  end if
444 
445  ! Compare the indices of the halo. First k, then j and
446  ! finally i.
447 
448  if (g1%myK < g2%myK) then
449  lessequalhalolisttype = .true.
450  return
451  else if (g1%myK > g2%myK) then
452  lessequalhalolisttype = .false.
453  return
454  end if
455 
456  if (g1%myJ < g2%myJ) then
457  lessequalhalolisttype = .true.
458  return
459  else if (g1%myJ > g2%myJ) then
460  lessequalhalolisttype = .false.
461  return
462  end if
463 
464  if (g1%myI < g2%myI) then
465  lessequalhalolisttype = .true.
466  return
467  else if (g1%myI > g2%myI) then
468  lessequalhalolisttype = .false.
469  return
470  end if
471 
472  ! No need to compare anything else; g1 == g2.
473 
474  else boundary
475 
476  ! Both halo's are internal halo's, whose donor is stored on
477  ! the same processor. Compare the donor blocks.
478 
479  if (g1%donorBlock < g2%donorBlock) then
480  lessequalhalolisttype = .true.
481  return
482  else if (g1%donorBlock > g2%donorBlock) then
483  lessequalhalolisttype = .false.
484  return
485  end if
486 
487  ! Also the blocks are identical. Compare the donor indices.
488  ! First the k index.
489 
490  if (g1%dK < g2%dK) then
491  lessequalhalolisttype = .true.
492  return
493  else if (g1%dK > g2%dK) then
494  lessequalhalolisttype = .false.
495  return
496  end if
497 
498  ! The j index.
499 
500  if (g1%dJ < g2%dJ) then
501  lessequalhalolisttype = .true.
502  return
503  else if (g1%dJ > g2%dJ) then
504  lessequalhalolisttype = .false.
505  return
506  end if
507 
508  ! And the i index.
509 
510  if (g1%dI < g2%dI) then
511  lessequalhalolisttype = .true.
512  return
513  else if (g1%dI > g2%dI) then
514  lessequalhalolisttype = .false.
515  return
516  end if
517 
518  ! The donors are identical. Compare the halo's.
519  ! First the block id.
520 
521  if (g1%myBlock < g2%myBlock) then
522  lessequalhalolisttype = .true.
523  return
524  else if (g1%myBlock > g2%myBlock) then
525  lessequalhalolisttype = .false.
526  return
527  end if
528 
529  ! Halo blocks are also identical. Finally compare the
530  ! halo indices. Start with k.
531 
532  if (g1%myK < g2%myK) then
533  lessequalhalolisttype = .true.
534  return
535  else if (g1%myK > g2%myK) then
536  lessequalhalolisttype = .false.
537  return
538  end if
539 
540  ! The j index.
541 
542  if (g1%myJ < g2%myJ) then
543  lessequalhalolisttype = .true.
544  return
545  else if (g1%myJ > g2%myJ) then
546  lessequalhalolisttype = .false.
547  return
548  end if
549 
550  ! The i index.
551 
552  if (g1%myI < g2%myI) then
553  lessequalhalolisttype = .true.
554  return
555  else if (g1%myI > g2%myI) then
556  lessequalhalolisttype = .false.
557  return
558  end if
559 
560  end if boundary
561 
562  ! Both entities are identical. So set lessEqual to .true.
563 
564  lessequalhalolisttype = .true.
565 
566  end function lessequalhalolisttype
567 
568  ! ================================================================
569 
570  logical function lesshalolisttype(g1, g2)
571  !
572  ! This function returns .true. if g1 < g2 and .false.
573  ! otherwise. It is basically the same as the lessEqual
574  ! function, except that the equality is now considered as
575  ! .false.
576  !
577  implicit none
578  !
579  ! Function arguments.
580  !
581  type(halolisttype), intent(in) :: g1, g2
582 
583  ! Compare the donor processors first. If not equal,
584  ! set the function appropriately and return.
585 
586  if (g1%donorProc < g2%donorProc) then
587  lesshalolisttype = .true.
588  return
589  else if (g1%donorProc > g2%donorProc) then
590  lesshalolisttype = .false.
591  return
592  end if
593 
594  ! Donor processors are identical. Now it depends whether we are
595  ! dealing with boundary halo's or not.
596 
597  boundary: if (g1%donorProc == -1) then ! And thus
598  ! g2%donorProc == -1
599 
600  ! Both halo's are boundary halo's. Compare the block ID of
601  ! the halo's.
602 
603  if (g1%myBlock < g2%myBlock) then
604  lesshalolisttype = .true.
605  return
606  else if (g1%myBlock > g2%myBlock) then
607  lesshalolisttype = .false.
608  return
609  end if
610 
611  ! Compare the boundary conditions, which are stored in
612  ! donorBlock. Note that the sequence in BCTypes is such that
613  ! the most important bc has the highest number.
614 
615  if (g1%donorBlock < g2%donorBlock) then
616  lesshalolisttype = .true.
617  return
618  else if (g1%donorBlock > g2%donorBlock) then
619  lesshalolisttype = .false.
620  return
621  end if
622 
623  ! As it is possible that indirect halo's need donor info from
624  ! direct halo's or even indirect halo's with a smaller level
625  ! of indirectness, compare the level of indirectness.
626 
627  if (g1%levOfInd < g2%levOfInd) then
628  lesshalolisttype = .true.
629  return
630  else if (g1%levOfInd > g2%levOfInd) then
631  lesshalolisttype = .false.
632  return
633  end if
634 
635  ! Compare the indices of the halo. First k, then j and
636  ! finally i.
637 
638  if (g1%myK < g2%myK) then
639  lesshalolisttype = .true.
640  return
641  else if (g1%myK > g2%myK) then
642  lesshalolisttype = .false.
643  return
644  end if
645 
646  if (g1%myJ < g2%myJ) then
647  lesshalolisttype = .true.
648  return
649  else if (g1%myJ > g2%myJ) then
650  lesshalolisttype = .false.
651  return
652  end if
653 
654  if (g1%myI < g2%myI) then
655  lesshalolisttype = .true.
656  return
657  else if (g1%myI > g2%myI) then
658  lesshalolisttype = .false.
659  return
660  end if
661 
662  ! No need to compare anything else. G1 == g2.
663 
664  else boundary
665 
666  ! Both halo's are internal halo's, whose donor is stored on
667  ! the same processor. Compare the donor blocks.
668 
669  if (g1%donorBlock < g2%donorBlock) then
670  lesshalolisttype = .true.
671  return
672  else if (g1%donorBlock > g2%donorBlock) then
673  lesshalolisttype = .false.
674  return
675  end if
676 
677  ! Also the blocks are identical. Compare the donor indices.
678  ! First the k index.
679 
680  if (g1%dK < g2%dK) then
681  lesshalolisttype = .true.
682  return
683  else if (g1%dK > g2%dK) then
684  lesshalolisttype = .false.
685  return
686  end if
687 
688  ! The j index.
689 
690  if (g1%dJ < g2%dJ) then
691  lesshalolisttype = .true.
692  return
693  else if (g1%dJ > g2%dJ) then
694  lesshalolisttype = .false.
695  return
696  end if
697 
698  ! And the i index.
699 
700  if (g1%dI < g2%dI) then
701  lesshalolisttype = .true.
702  return
703  else if (g1%dI > g2%dI) then
704  lesshalolisttype = .false.
705  return
706  end if
707 
708  ! The donors are identical. Compare the halo's.
709  ! First the block id.
710 
711  if (g1%myBlock < g2%myBlock) then
712  lesshalolisttype = .true.
713  return
714  else if (g1%myBlock > g2%myBlock) then
715  lesshalolisttype = .false.
716  return
717  end if
718 
719  ! Halo blocks are also identical. Finally compare the
720  ! halo indices. Start with k.
721 
722  if (g1%myK < g2%myK) then
723  lesshalolisttype = .true.
724  return
725  else if (g1%myK > g2%myK) then
726  lesshalolisttype = .false.
727  return
728  end if
729 
730  ! The j index.
731 
732  if (g1%myJ < g2%myJ) then
733  lesshalolisttype = .true.
734  return
735  else if (g1%myJ > g2%myJ) then
736  lesshalolisttype = .false.
737  return
738  end if
739 
740  ! The i index.
741 
742  if (g1%myI < g2%myI) then
743  lesshalolisttype = .true.
744  return
745  else if (g1%myI > g2%myI) then
746  lesshalolisttype = .false.
747  return
748  end if
749 
750  end if boundary
751 
752  ! Both entities are identical.
753  ! So set lessHaloListType to .false.
754 
755  lesshalolisttype = .false.
756 
757  end function lesshalolisttype
758 
759 end module halolist
760 
762  !
763  ! Local module, which contains the definition of the derived
764  ! datatype used to test for negative volumes in the grid.
765  !
766  implicit none
767  save
768 
770 
771  ! blockHasNegVol: Whether or not the block
772  ! contains negative volumes.
773  ! volumeIsNeg(2:il,2:jl,2:kl): Whether or not the owned volumes
774  ! are negative.
775 
776  logical :: blockhasnegvol
777  logical, dimension(:, :, :), pointer :: volumeisneg
778  logical :: blockhasskewedvol
779  logical, dimension(:, :, :), pointer :: volumeisskewed
780 
781  end type checkvolblocktype
782 
783 end module checkvolblock
785  !
786  ! Local module that contains derived datatypes as well as arrays
787  ! of these derived datatypes to store information related to
788  ! periodicity.
789  !
790  use precision
791  implicit none
792  save
793 
794  public
795  private :: lesscgnsperiodictype
796  private :: equalcgnsperiodictype
797  private :: lessperiodicsubfaceshalot
798  private :: lessequalperiodicsubfaceshalot
799  private :: equalperiodicsubfaceshalot
800 
801  ! Definition of the derived data type for storing the periodic
802  ! faces of the cgns grid a bit easier.
803 
805 
806  ! cgnsBlock :: the block ID in the cgns grid.
807  ! cgnsSubface :: the suface ID in this block.
808 
809  integer(kind=intType) :: cgnsblock, cgnssubface
810 
811  end type cgnsperiodictype
812 
813  ! Interface for the extension of the operators < and ==.
814  ! These are needed for the sorting and searching of
815  ! cgnsPeriodicType.
816 
817  interface operator(<)
818  module procedure lesscgnsperiodictype
819  end interface operator(<)
820 
821  interface operator(==)
822  module procedure equalcgnsperiodictype
823  end interface operator(==)
824 
825  ! nPeriodicGlobal :: Total number of periodic faces in cgns grid.
826  ! periodicGlobal :: The corresponding faces.
827 
828  integer(kind=intType) :: nperiodicglobal
829  type(cgnsperiodictype), dimension(:), allocatable :: periodicglobal
830 
831  ! Definition of the derived data type to store the periodic
832  ! subfaces that are crossed when the halo and the donor are
833  ! connected. The direction is from the halo to the donor.
834 
836 
837  ! internalHalo: Whether or not the halo is an internal
838  ! halo, i.e. it is stored on the
839  ! same processor as the donor.
840  ! indexInHaloList: The corresponding index in either the
841  ! node or cell halo list.
842  ! nPeriodicSubfaces: Number of periodic subfaces that are
843  ! crossed. This is at most the level of
844  ! indirectness of the halo.
845  ! periodicSubfaces: The corresponding subfaces ID's according
846  ! to the sequence defined in periodicGlobal.
847 
848  logical :: internalhalo
849  integer(kind=intType) :: indexinhalolist
850  integer(kind=intType) :: nperiodicsubfaces
851  integer(kind=intType), dimension(:), pointer :: periodicsubfaces
852 
853  end type periodicsubfaceshalotype
854 
855  ! Interface for the extension of the operators <, <= and ==.
856  ! This is needed for the sorting and comparing of variables of
857  ! the type periodicSubfacesHaloType
858 
859  interface operator(<)
860  module procedure lessperiodicsubfaceshalot
861  end interface operator(<)
862 
863  interface operator(<=)
864  module procedure lessequalperiodicsubfaceshalot
865  end interface operator(<=)
866 
867  interface operator(==)
868  module procedure equalperiodicsubfaceshalot
869  end interface operator(==)
870 
871  !=================================================================
872 
873 contains
874 
875  !=================================================================
876  !
877  ! Functions to simulate the operators < and ==.
878  !
879  logical function lesscgnsperiodictype(g1, g2)
880  !
881  ! lessCGNSPeriodicType returns .true. if g1 is considered
882  ! smaller than g2. This comparison is first based on the block
883  ! ID followed by the subface id.
884  !
885  implicit none
886  !
887  ! Function arguments.
888  !
889  type(cgnsperiodictype), intent(in) :: g1, g2
890 
891  ! Compare the block ID. If not equal set lessCGNSPeriodicType
892  ! accordingly.
893 
894  if (g1%cgnsBlock < g2%cgnsBlock) then
895  lesscgnsperiodictype = .true.
896  return
897  else if (g1%cgnsBlock > g2%cgnsBlock) then
898  lesscgnsperiodictype = .false.
899  return
900  end if
901 
902  ! Block ID's are identical. Compare the subfaces.
903 
904  if (g1%cgnsSubface < g2%cgnsSubface) then
905  lesscgnsperiodictype = .true.
906  return
907  else if (g1%cgnsSubface > g2%cgnsSubface) then
908  lesscgnsperiodictype = .false.
909  return
910  end if
911 
912  ! Both objects are identical.
913  ! Set lessCGNSPeriodicType to .false.
914 
915  lesscgnsperiodictype = .false.
916 
917  end function lesscgnsperiodictype
918 
919  ! ================================================================
920 
921  logical function equalcgnsperiodictype(g1, g2)
922  !
923  ! equalCGNSPeriodicType returns .true. if g1 is considered
924  ! equal to g2, i.e. both the block and subface ID must match,
925  !
926  implicit none
927  !
928  ! Function arguments.
929  !
930  type(cgnsperiodictype), intent(in) :: g1, g2
931 
932  equalcgnsperiodictype = .false.
933  if (g1%cgnsBlock == g2%cgnsBlock .and. &
934  g1%cgnsSubface == g2%cgnsSubface) &
935  equalcgnsperiodictype = .true.
936 
937  end function equalcgnsperiodictype
938 
939  ! ================================================================
940 
941  logical function lessperiodicsubfaceshalot(g1, g2)
942  !
943  ! lessPeriodicSubfacesHaloT returns .true. if g1 is
944  ! considered smaller than g2.
945  !
946  implicit none
947  !
948  ! Function arguments.
949  !
950  type(periodicsubfaceshalotype), intent(in) :: g1, g2
951  !
952  ! Local variables.
953  !
954  integer(kind=intType) :: nn, i1, i2
955 
956  ! First compare whether or not both g1 and g2 are internal
957  ! halo's. Fortran does not allow a direct comparison of
958  ! logicals and therefore the integers i1 and i2 are used.
959 
960  i1 = 1; if (g1%internalHalo) i1 = 0
961  i2 = 1; if (g2%internalHalo) i2 = 0
962 
963  if (i1 < i2) then
964  lessperiodicsubfaceshalot = .true.
965  return
966  else if (i1 > i2) then
967  lessperiodicsubfaceshalot = .false.
968  return
969  end if
970 
971  ! Compare the number of periodic subfaces.
972 
973  if (g1%nPeriodicSubfaces < g2%nPeriodicSubfaces) then
974  lessperiodicsubfaceshalot = .true.
975  return
976  else if (g1%nPeriodicSubfaces > g2%nPeriodicSubfaces) then
977  lessperiodicsubfaceshalot = .false.
978  return
979  end if
980 
981  ! The number of periodic subfaces is the same. Compare the
982  ! subfaces themselves. It is assumed that the subfaces are
983  ! sorted in increading order. This can be done, because the
984  ! periodic transformations are commuting matrices.
985 
986  do nn = 1, g1%nPeriodicSubfaces
987  if (g1%periodicSubfaces(nn) < g2%periodicSubfaces(nn)) then
988  lessperiodicsubfaceshalot = .true.
989  return
990  else if (g1%periodicSubfaces(nn) > g2%periodicSubfaces(nn)) then
991  lessperiodicsubfaceshalot = .false.
992  return
993  end if
994  end do
995 
996  ! The periodic subfaces are identical as well. Compare the
997  ! indices in the list.
998 
999  if (g1%indexInHaloList < g2%indexInHaloList) then
1000  lessperiodicsubfaceshalot = .true.
1001  return
1002  else if (g1%indexInHaloList > g2%indexInHaloList) then
1003  lessperiodicsubfaceshalot = .false.
1004  return
1005  end if
1006 
1007  ! Both objects are the same. Return .false.
1008 
1009  lessperiodicsubfaceshalot = .false.
1010 
1011  end function lessperiodicsubfaceshalot
1012 
1013  ! ================================================================
1014 
1015  logical function lessequalperiodicsubfaceshalot(g1, g2)
1016  !
1017  ! lessEqualPeriodicSubfacesHaloT returns .true. if g1 is
1018  ! considered smaller than or equal to g2.
1019  !
1020  implicit none
1021  !
1022  ! Function arguments.
1023  !
1024  type(periodicsubfaceshalotype), intent(in) :: g1, g2
1025  !
1026  ! Local variables.
1027  !
1028  integer(kind=intType) :: nn, i1, i2
1029 
1030  ! First compare whether or not both g1 and g2 are internal
1031  ! halo's. Fortran does not allow a direct comparison of
1032  ! logicals and therefore the integers i1 and i2 are used.
1033 
1034  i1 = 1; if (g1%internalHalo) i1 = 0
1035  i2 = 1; if (g2%internalHalo) i2 = 0
1036 
1037  if (i1 < i2) then
1038  lessequalperiodicsubfaceshalot = .true.
1039  return
1040  else if (i1 > i2) then
1041  lessequalperiodicsubfaceshalot = .false.
1042  return
1043  end if
1044 
1045  ! Compare the number of periodic subfaces.
1046 
1047  if (g1%nPeriodicSubfaces < g2%nPeriodicSubfaces) then
1048  lessequalperiodicsubfaceshalot = .true.
1049  return
1050  else if (g1%nPeriodicSubfaces > g2%nPeriodicSubfaces) then
1051  lessequalperiodicsubfaceshalot = .false.
1052  return
1053  end if
1054 
1055  ! The number of periodic subfaces is the same. Compare the
1056  ! subfaces themselves. It is assumed that the subfaces are
1057  ! sorted in increading order. This can be done, because the
1058  ! periodic transformations are commuting matrices.
1059 
1060  do nn = 1, g1%nPeriodicSubfaces
1061  if (g1%periodicSubfaces(nn) < g2%periodicSubfaces(nn)) then
1062  lessequalperiodicsubfaceshalot = .true.
1063  return
1064  else if (g1%periodicSubfaces(nn) > g2%periodicSubfaces(nn)) then
1065  lessequalperiodicsubfaceshalot = .false.
1066  return
1067  end if
1068  end do
1069 
1070  ! The periodic subfaces are identical as well. Compare the
1071  ! indices in the list.
1072 
1073  if (g1%indexInHaloList < g2%indexInHaloList) then
1074  lessequalperiodicsubfaceshalot = .true.
1075  return
1076  else if (g1%indexInHaloList > g2%indexInHaloList) then
1077  lessequalperiodicsubfaceshalot = .false.
1078  return
1079  end if
1080 
1081  ! Both objects are the same. Return .true.
1082 
1083  lessequalperiodicsubfaceshalot = .true.
1084 
1085  end function lessequalperiodicsubfaceshalot
1086 
1087  ! ================================================================
1088 
1089  logical function equalperiodicsubfaceshalot(g1, g2)
1090  !
1091  ! equalPeriodicSubfacesHaloT returns .true. if g1 is
1092  ! considered equal to g2. The equal operator is only used to
1093  ! find the different number of periodic transformations in
1094  ! determinePeriodicData. Hence only the periodic subfaces of
1095  ! the halo's are compared and g1 and g2 are considered equal
1096  ! if the subfaces are equal, even if other member variables
1097  ! differ.
1098  !
1099  implicit none
1100  !
1101  ! Function arguments.
1102  !
1103  type(periodicsubfaceshalotype), intent(in) :: g1, g2
1104  !
1105  ! Local variables.
1106  !
1107  integer(kind=intType) :: nn
1108 
1109  if (g1%nPeriodicSubfaces /= g2%nPeriodicSubfaces) then
1110  equalperiodicsubfaceshalot = .false.
1111  return
1112  end if
1113 
1114  do nn = 1, g1%nPeriodicSubfaces
1115  if (g1%periodicSubfaces(nn) /= g2%periodicSubfaces(nn)) then
1116  equalperiodicsubfaceshalot = .false.
1117  return
1118  end if
1119  end do
1120 
1121  equalperiodicsubfaceshalot = .true.
1122 
1123  end function equalperiodicsubfaceshalot
1124 
1125 end module periodicinfo
1126 module bchalo
1127  !
1128  ! This local module contains the derived datatype bcHaloType,
1129  ! which is used to determine the boundary condition for an
1130  ! indirect halo when the nearest direct halo's are all boundary
1131  ! halo's.
1132  !
1133  use precision
1134  implicit none
1135  save
1136 
1137  public
1138  private :: lessequalbchalotype
1139  !
1140  ! The definition of the derived datatype.
1141  !
1143 
1144  ! directHalo: Index in the haloListType where the
1145  ! corresponding direct halo is stored.
1146  ! BC: Corresponding boundary condition.
1147 
1148  integer(kind=intType) :: directhalo, bc
1149 
1150  end type bchalotype
1151 
1152  ! Interface for the extension of the operator <= needed for the
1153  ! sorting of bcHaloType. Note that the = operator does not
1154  ! need to be defined, because bcHaloType only contains
1155  ! primitive types.
1156 
1157  interface operator(<=)
1158  module procedure lessequalbchalotype
1159  end interface operator(<=)
1160 
1161 contains
1162  !
1163  logical function lessequalbchalotype(g1, g2)
1164  !
1165  ! Function to simulate the operator <= for bcHaloType.
1166  ! It first compares the boundary condition. If equal the index
1167  ! of the direct halo is compared, although this is not really
1168  ! important.
1169  ! LessEqual returns .true. if g1 <= g2 and .false. otherwise.
1170  !
1171  implicit none
1172  !
1173  ! Function arguments.
1174  !
1175  type(bchalotype), intent(in) :: g1, g2
1176 
1177  ! First compare the boundary conditions. Note that the sequence
1178  ! in BCTypes is such that the most important BC has the
1179  ! highest number.
1180 
1181  if (g1%BC < g2%BC) then
1182  lessequalbchalotype = .true.
1183  return
1184  else if (g1%BC > g2%BC) then
1185  lessequalbchalotype = .false.
1186  return
1187  end if
1188 
1189  ! Boundary conditions are equal. Just compare the index.
1190 
1191  if (g1%directHalo < g2%directHalo) then
1192  lessequalbchalotype = .true.
1193  return
1194  else if (g1%directHalo > g2%directHalo) then
1195  lessequalbchalotype = .false.
1196  return
1197  end if
1198 
1199  ! g1 and g2 are equal. Return .true.
1200 
1201  lessequalbchalotype = .true.
1202 
1203  end function lessequalbchalotype
1204 
1205  ! ================================================================
1206 
1207  subroutine sortbchalotype(bcHaloArray, nn)
1208  !
1209  ! SortBCHaloType sorts the given number of BCHalo's in
1210  ! increasing order. Note that this routine is called sort and
1211  ! not qsort, because only an insertion sort is done here. The
1212  ! reason is that nn <= 3 and thus an insertion sort is okay.
1213  !
1214  implicit none
1215  !
1216  ! Subroutine arguments
1217  !
1218  integer(kind=intType), intent(in) :: nn
1219  type(bchalotype), dimension(*), intent(inout) :: bcHaloArray
1220  !
1221  ! Local variables.
1222  !
1223  integer(kind=intType) :: i, j
1224 
1225  type(bchalotype) :: a
1226 
1227  do j = 1, nn
1228  a = bchaloarray(j)
1229  do i = (j - 1), 1, -1
1230  if (bchaloarray(i) <= a) exit
1231  bchaloarray(i + 1) = bchaloarray(i)
1232  end do
1233  bchaloarray(i + 1) = a
1234  end do
1235 
1236  end subroutine sortbchalotype
1237 
1238 end module bchalo
1239 
1241  !
1242  ! This local module contains the derived datatype
1243  ! coarse1to1SubfaceType, which is used to determine the 1 to 1
1244  ! block boundaries for the coarser grids.
1245  !
1246  use precision
1247  implicit none
1248  save
1249  !
1250  ! The definition of the derived datatype.
1251  !
1253 
1254  ! Nodal range in the three coordinates directions for the
1255  ! coarse grid subface.
1256 
1257  integer(kind=intType) :: ibeg, jbeg, kbeg, iend, jend, kend
1258 
1259  ! Processor and block id of the neighboring block.
1260 
1261  integer(kind=intType) :: neighproc, neighblock
1262 
1263  ! Number of points in the three coordinate directions for the
1264  ! coarse grid donor subface.
1265 
1266  integer(kind=intType) :: ndi, ndj, ndk
1267 
1268  ! Corresponding i, j and k indices of the fine grid donor block
1269  ! for each of the coarse grid subface lines.
1270 
1271  integer(kind=intType), dimension(:), pointer :: idfine
1272  integer(kind=intType), dimension(:), pointer :: jdfine
1273  integer(kind=intType), dimension(:), pointer :: kdfine
1274 
1275  end type coarse1to1subfacetype
1276 
1277  ! Number of 1 to 1 fine grid subfaces on this processor.
1278 
1279  integer(kind=intType) :: nsubface1to1
1280 
1281  ! Array of 1 to 1 subfaces.
1282 
1283  type(coarse1to1subfacetype), dimension(:), allocatable :: subface1to1
1284 
1285 end module coarse1to1subface
1286 
1287 ! ==================================================================
1288 
1290  !
1291  ! This local module contains the derived datatype
1292  ! coarseningInfoType, which stores for a given block the grid
1293  ! lines to keep for the coarse grid.
1294  !
1296 
1297  ! Logical, which indicate whether or not a fine grid 1 to 1
1298  ! block boundary is still a 1 to 1 block boundary on the
1299  ! coarse grid.
1300 
1301  logical, dimension(:), pointer :: coarseis1to1
1302 
1303  end type coarseninginfotype
1304 
1305  ! Array to store the info for all the blocks.
1306 
1307  type(coarseninginfotype), dimension(:), allocatable :: coarseinfo
1308 
1309 end module coarseninginfo
subroutine sortbchalotype(bcHaloArray, nn)
type(coarse1to1subfacetype), dimension(:), allocatable subface1to1
integer(kind=inttype) nsubface1to1
type(coarseninginfotype), dimension(:), allocatable coarseinfo
integer(kind=inttype) nnodehalo1st
type(halolisttype), dimension(:), allocatable cellhalo1st
type(halolisttype), dimension(:), allocatable cellhalo2nd
integer(kind=inttype), dimension(:, :), allocatable transformcell
type(halolisttype), dimension(:), allocatable nodehalo1st
integer(kind=inttype) ncellhalo2nd
integer(kind=inttype) iinode1st
integer(kind=inttype) iicell1st
integer(kind=inttype), dimension(:, :), allocatable transformnode
type(indexlisttype), dimension(:), allocatable nodeindex
integer(kind=inttype) ncellhalo1st
type(indexlisttype), dimension(:), allocatable cellindex
integer(kind=inttype) iicell2nd
integer(kind=inttype), dimension(:), allocatable nhaloperproc
integer(kind=inttype) nindhalo
integer(kind=inttype) nlevofind
integer(kind=inttype), dimension(:), allocatable nhaloperlev
type(indirecthalotype), dimension(:), allocatable indhalo
type(cgnsperiodictype), dimension(:), allocatable periodicglobal
integer(kind=inttype) nperiodicglobal
integer, parameter realtype
Definition: precision.F90:112