ADflow  v1.0
ADflow is a finite volume RANS solver tailored for gradient-based aerodynamic design optimization.
adtData.F90
Go to the documentation of this file.
1 module adtdata
2 !
3 ! Module, which defines the derived data types and the arrays to
4 ! store multiple ADT's. An array is chosen to store multiple
5 ! ATD's rather than a linked list, because this is more
6 ! convenient during the search. When the ADT's are built there
7 ! is some additional work due to reallocation. However this is
8 ! negligible due to the usage of pointers.
9 !
10  use constants, only: inttype, realtype, adtelementtype, alwaysrealtype
11  implicit none
12  save
13 !
14 ! Define the functions needed for the sorting of the derived
15 ! data types to be private, i.e. they can only be accessed
16 ! within this module.
17 !
18  public
19  private :: adtbboxtargettypelessequal
20  private :: adtbboxtargettypeless
21  private :: adttypeassign
22 !
23 ! Definition of the derived data type store a leaf of an ADT.
24 ! The ADT itself is an array of these leaves.
25 !
27 
28  ! children(2): Children of the parent. If negative it means that
29  ! it is a terminal leaf and the absolute values
30  ! indicate the bounding box id's. Note that it is
31  ! allowed that 1 child is negative and the other
32  ! positive.
33  ! xMin(6): The minimum coordinates of the leaf.
34  ! xMax(6): The maximum coordinates of the leaf.
35 
36  integer(kind=intType), dimension(2) :: children
37  real(kind=realtype), dimension(6) :: xmin, xmax
38 
39  end type adtleaftype
40 !
41 ! The definition of adtBBoxTargetType, which stores the data of
42 ! a possible bounding box which minimizes the distances to the
43 ! given coordinate.
44 !
46 
47  ! ID: The id of the bounding box in the list.
48  ! posDist2: the possible minimum distance squared to the active
49  ! coordinate.
50 
51  integer(kind=intType) :: id
52  real(kind=realtype) :: posdist2
53 
54  end type adtbboxtargettype
55 
56  ! Interfaces for the extension of the operators <= and <.
57  ! These are needed for the sorting of BBoxTargetType. Note
58  ! that the = operator does not need to be defined, because
59  ! BBoxTargetType only contains primitive types.
60 
61  interface operator(<=)
62  module procedure adtbboxtargettypelessequal
63  end interface
64 
65  interface operator(<)
66  module procedure adtbboxtargettypeless
67  end interface
68 !
69 ! Definition of the derived data type to store an ADT.
70 !
71  type adttype
72 
73  ! comm: The communicator of this ADT.
74  ! nProcs: The number of processors which participate in this
75  ! ADT.
76  ! myID: My processor ID in the processor group of comm.
77 
78  integer :: comm, nprocs, myid
79 
80  ! adtType: Type of ADT. Possibilities are adtSurfaceADT and
81  ! adtVolumeADT.
82  ! adtID: The given ID of the ADT.
83  ! isActive: Whether or not the ADT is active. If not, this
84  ! entry could be used during a reallocation.
85 
86  integer :: adttype
87  character(len=64) :: adtid
88  logical :: isactive
89 
90  ! nNodes: Number of local nodes in the given grid.
91  ! nTria: Number of local triangles in the given grid.
92  ! nQuads: Idem for the quadrilaterals.
93  ! nTetra: Idem for the tetrahedra.
94  ! nPyra: Idem for the pyramids.
95  ! nPrisms: Idem for the prisms.
96  ! nHexa: Idem for the hexahedra.
97 
98  integer(kind=intType) :: nnodes, ntria, nquads
99  integer(kind=intType) :: ntetra, npyra, nprisms, nhexa
100 
101  ! coor(3,nNodes): Nodal coordinates of the local grid.
102  ! To save memory this pointer is not
103  ! allocated, but set to the data given.
104  real(kind=realtype), dimension(:, :), pointer :: coor
105 
106  ! triaConn(3,nTria): Local connectivity of the triangles.
107  ! To save memory this pointer is not
108  ! allocated, but set to the data given.
109  ! quadsConn(4,nQuads): Idem for the quadrilaterals.
110  ! tetraConn(4,nTetra): Idem for the tetrahedra.
111  ! pyraConn(5,nPyra): Idem for the pyramids.
112  ! prismsConn(6,nPrisms): Idem for the prisms.
113  ! hexaConn(8,nHexa): Idem for the hexahedra.
114 
115  integer(kind=intType), dimension(:, :), pointer :: triaconn
116  integer(kind=intType), dimension(:, :), pointer :: quadsconn
117  integer(kind=intType), dimension(:, :), pointer :: tetraconn
118  integer(kind=intType), dimension(:, :), pointer :: pyraconn
119  integer(kind=intType), dimension(:, :), pointer :: prismsconn
120  integer(kind=intType), dimension(:, :), pointer :: hexaconn
121 
122  ! nRootLeaves: Number of non-empty root leaves.
123  ! This number is of course less than or
124  ! equal to nProcs.
125  ! myEntryInRootProcs: If the local tree is not empty, this
126  ! contains the entry in rootLeavesProcs.
127  ! rootLeavesProcs(:): The corresponding processor ID's.
128  ! rootBBoxes(3,2,:): The 3D bounding boxes of the non-empty
129  ! root leaves.
130 
131  integer :: nrootleaves, myentryinrootprocs
132  integer, dimension(:), pointer :: rootleavesprocs
133  real(kind=realtype), dimension(:, :, :), pointer :: rootbboxes
134 
135  ! nBBoxes: Number of bounding boxes stored in
136  ! the ADT.
137  ! elementType(nBBoxes): The corresponding element type of the
138  ! bounding box.
139  ! elementID(nBBoxes): The corresponding entry in the element
140  ! connectivity of the bounding box.
141  ! xBBox(6,nBBoxes): The coordinates of the bounding boxes
142  ! of the elements stored in the ADT.
143 
144  integer(kind=intType) :: nbboxes
145 
146  integer(kind=adtElementType), dimension(:), pointer :: elementtype
147  integer(kind=intType), dimension(:), pointer :: elementid
148  real(kind=realtype), dimension(:, :), pointer :: xbbox
149 
150  ! nLeaves: Number of present in the ADT. Due to the
151  ! variable splitting the tree is optimally
152  ! balanced and therefore nLeaves = nBBoxes -1.
153  ! ADTree(nLeaves): The alternating digital tree.
154 
155  integer(kind=intType) :: nleaves
156  type(adtleaftype), dimension(:), pointer :: adtree
157 
158  end type adttype
159 
160  ! Interface for the extension of the operator =.
161 
162  interface assignment(=)
163  module procedure adttypeassign
164  end interface
165 !
166 ! Variables stored in this module.
167 !
168  ! ADTs(:): The array to store the different ADT's.
169 
170  type(adttype), dimension(:), allocatable, target :: adts
171 
172  ! nProcRecv: Number of processors from which I receive
173  ! coordinates that must be searched in my ADT.
174  ! nCoorMax: Maximum number of coordinates that can be
175  ! searched during an interpolation round.
176  ! nRounds: Number of rounds in the outer loop of the search
177  ! algorithm.
178  ! nLocalInterpol: Number of local coordinates that must be
179  ! searched in the locally stored tree.
180 
181  integer :: nprocrecv
182 
183  integer(kind=intType) :: ncoormax
184  integer(kind=intType) :: nrounds
185  integer(kind=intType) :: nlocalinterpol
186 
187  ! procRecv: Processor ID's from which I will receive
188  ! coordinates.
189  ! nCoorProcRecv: Number of coordinates I must receive from the
190  ! processors which send coordinates to me.
191  ! nCoorPerRootLeaf: Number of coordinates, which may be searched
192  ! in each of the local ADT's. The array is in
193  ! cumulative storage format.
194  ! mCoorPerRootLeaf: Idem, but its contents changes during the
195  ! iterative algorithm.
196  ! coorPerRootLeaf: The ID's of the corresponding coordinates.
197 
198  integer, dimension(:), allocatable :: procrecv
199 
200  integer(kind=intType), dimension(:), allocatable :: ncoorprocrecv
201  integer(kind=intType), dimension(:), allocatable :: ncoorperrootleaf
202  integer(kind=intType), dimension(:), allocatable :: mcoorperrootleaf
203  integer(kind=intType), dimension(:), allocatable :: coorperrootleaf
204 
205  !=================================================================
206 
207 contains
208 
209  !===============================================================
210 
211  logical function adtbboxtargettypelessequal(g1, g2)
212 !
213 ! This function returns .true. if g1 <= g2. The comparison is
214 ! firstly based on the possible minimum distance such that the
215 ! most likely candidates are treated first. In case of ties
216 ! the boundary box ID is considered.
217 ! Function intent(in) arguments.
218 ! ------------------------------
219 ! g1, g2: The two instances of the derived datatype that most
220 ! be compared.
221 !
222  implicit none
223 !
224 ! Function arguments.
225 !
226  type(adtbboxtargettype), intent(in) :: g1, g2
227 
228  ! Compare the possible minimum distances.
229 
230  if (g1%posDist2 < g2%posDist2) then
231  adtbboxtargettypelessequal = .true.
232  return
233  else if (g1%posDist2 > g2%posDist2) then
234  adtbboxtargettypelessequal = .false.
235  return
236  end if
237 
238  ! Compare the bounding box ID's.
239 
240  if (g1%ID < g2%ID) then
241  adtbboxtargettypelessequal = .true.
242  return
243  else if (g1%ID > g2%ID) then
244  adtbboxtargettypelessequal = .false.
245  return
246  end if
247 
248  ! g1 and g2 are identical. Return .true.
249 
250  adtbboxtargettypelessequal = .true.
251 
252  end function adtbboxtargettypelessequal
253 
254  !===============================================================
255 
256  logical function adtbboxtargettypeless(g1, g2)
257 !
258 ! This function returns .true. if g1 < g2. The comparison is
259 ! firstly based on the possible minimum distance such that the
260 ! most likely candidates are treated first. In case of ties
261 ! the boundary box ID is considered.
262 ! Function intent(in) arguments.
263 ! ------------------------------
264 ! g1, g2: The two instances of the derived datatype that most
265 ! be compared.
266 !
267  implicit none
268 !
269 ! Function arguments.
270 !
271  type(adtbboxtargettype), intent(in) :: g1, g2
272 
273  ! Compare the possible minimum distances.
274 
275  if (g1%posDist2 < g2%posDist2) then
276  adtbboxtargettypeless = .true.
277  return
278  else if (g1%posDist2 > g2%posDist2) then
279  adtbboxtargettypeless = .false.
280  return
281  end if
282 
283  ! Compare the bounding box ID's.
284 
285  if (g1%ID < g2%ID) then
286  adtbboxtargettypeless = .true.
287  return
288  else if (g1%ID > g2%ID) then
289  adtbboxtargettypeless = .false.
290  return
291  end if
292 
293  ! g1 and g2 are identical. Return .false.
294 
295  adtbboxtargettypeless = .false.
296 
297  end function adtbboxtargettypeless
298 
299  !===============================================================
300 
301  subroutine adttypeassign(g1, g2)
302 !
303 ! This subroutine defines the generic assignment operator for
304 ! the derived datatype adtType. The contents of g1 is copied
305 ! into g2, where pointers just point to the other pointers,
306 ! i.e. no additional allocation takes place.
307 ! Subroutine intent(in) arguments.
308 ! --------------------------------
309 ! g2: Entity whose data should be copied.
310 ! Subroutine intent(out) arguments.
311 ! ---------------------------------
312 ! g1: Entity whose data should be assigned.
313 !
314  implicit none
315 !
316 ! Subroutine arguments.
317 !
318  type(adttype), intent(in) :: g2
319  type(adttype), intent(out) :: g1
320 
321  g1%comm = g2%comm
322  g1%nProcs = g2%nProcs
323  g1%myID = g2%myID
324  g1%adtType = g2%adtType
325  g1%adtID = g2%adtID
326  g1%isActive = g2%isActive
327 
328  g1%nNodes = g2%nNodes
329  g1%nTria = g2%nTria
330  g1%nQuads = g2%nQuads
331  g1%nTetra = g2%nTetra
332  g1%nPyra = g2%nPyra
333  g1%nPrisms = g2%nPrisms
334  g1%nHexa = g2%nHexa
335 
336  g1%coor => g2%coor
337  g1%triaConn => g2%triaConn
338  g1%quadsConn => g2%quadsConn
339  g1%tetraConn => g2%tetraConn
340  g1%pyraConn => g2%pyraConn
341  g1%prismsConn => g2%prismsConn
342  g1%hexaConn => g2%hexaConn
343 
344  g1%nRootLeaves = g2%nRootLeaves
345  g1%myEntryInRootProcs = g2%myEntryInRootProcs
346  g1%rootLeavesProcs => g2%rootLeavesProcs
347  g1%rootBBoxes => g2%rootBBoxes
348 
349  g1%nBBoxes = g2%nBBoxes
350  g1%elementType => g2%elementType
351  g1%elementID => g2%elementID
352  g1%xBBox => g2%xBBox
353 
354  g1%nLeaves = g2%nLeaves
355  g1%ADTree => g2%ADTree
356 
357  end subroutine adttypeassign
358 
359 end module adtdata
integer(kind=inttype) nrounds
Definition: adtData.F90:184
integer nprocrecv
Definition: adtData.F90:181
integer(kind=inttype) ncoormax
Definition: adtData.F90:183
integer, dimension(:), allocatable procrecv
Definition: adtData.F90:198
integer(kind=inttype), dimension(:), allocatable ncoorprocrecv
Definition: adtData.F90:200
integer(kind=inttype), dimension(:), allocatable coorperrootleaf
Definition: adtData.F90:203
integer(kind=inttype) nlocalinterpol
Definition: adtData.F90:185
type(adttype), dimension(:), allocatable, target adts
Definition: adtData.F90:170
integer(kind=inttype), dimension(:), allocatable ncoorperrootleaf
Definition: adtData.F90:201
integer(kind=inttype), dimension(:), allocatable mcoorperrootleaf
Definition: adtData.F90:202