processorTopology.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "processorTopology.H"
27 #include "polyBoundaryMesh.H"
28 #include "processorPolyPatch.H"
29 #include "commSchedule.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 Foam::labelList Foam::processorTopology::procNeighbours
34 (
35  const label nProcs,
36  const polyBoundaryMesh& patches
37 )
38 {
39  // Determine number of processor neighbours and max neighbour id.
40 
41  label nNeighbours = 0;
42 
43  label maxNb = 0;
44 
45  boolList isNeighbourProc(nProcs, false);
46 
48  {
49  const polyPatch& patch = patches[patchi];
50 
51  if (isA<processorPolyPatch>(patch))
52  {
53  const processorPolyPatch& procPatch =
54  refCast<const processorPolyPatch>(patch);
55 
56  label pNeighbProcNo = procPatch.neighbProcNo();
57 
58  if (!isNeighbourProc[pNeighbProcNo])
59  {
60  nNeighbours++;
61 
62  maxNb = max(maxNb, procPatch.neighbProcNo());
63 
64  isNeighbourProc[pNeighbProcNo] = true;
65  }
66  }
67  }
68 
69  labelList neighbours(nNeighbours, -1);
70 
71  nNeighbours = 0;
72 
73  forAll(isNeighbourProc, proci)
74  {
75  if (isNeighbourProc[proci])
76  {
77  neighbours[nNeighbours++] = proci;
78  }
79  }
80 
81  procPatchMap_.setSize(maxNb + 1);
82  procPatchMap_ = -1;
83 
85  {
86  const polyPatch& patch = patches[patchi];
87 
88  if (isA<processorPolyPatch>(patch))
89  {
90  const processorPolyPatch& procPatch =
91  refCast<const processorPolyPatch>(patch);
92 
93  // Construct reverse map
94  procPatchMap_[procPatch.neighbProcNo()] = patchi;
95  }
96  }
97 
98  return neighbours;
99 }
100 
101 
102 Foam::lduSchedule Foam::processorTopology::nonBlockingSchedule
103 (
104  const polyBoundaryMesh& patches
105 )
106 {
107  lduSchedule patchSchedule(2*patches.size());
108 
109  label patchEvali = 0;
110 
111  // 1. All non-processor patches
112  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
113 
114  // Have evaluate directly after initEvaluate. Could have them separated
115  // as long as they're not intermingled with processor patches since
116  // then e.g. any reduce parallel traffic would interfere with the
117  // processor swaps.
118 
120  {
121  if (!isA<processorPolyPatch>(patches[patchi]))
122  {
123  patchSchedule[patchEvali].patch = patchi;
124  patchSchedule[patchEvali++].init = true;
125  patchSchedule[patchEvali].patch = patchi;
126  patchSchedule[patchEvali++].init = false;
127  }
128  }
129 
130  // 2. All processor patches
131  // ~~~~~~~~~~~~~~~~~~~~~~~~
132 
133  // 2a. initEvaluate
135  {
136  if (isA<processorPolyPatch>(patches[patchi]))
137  {
138  patchSchedule[patchEvali].patch = patchi;
139  patchSchedule[patchEvali++].init = true;
140  }
141  }
142 
143  // 2b. evaluate
145  {
146  if (isA<processorPolyPatch>(patches[patchi]))
147  {
148  patchSchedule[patchEvali].patch = patchi;
149  patchSchedule[patchEvali++].init = false;
150  }
151  }
152 
153  return patchSchedule;
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
158 
160 (
161  const polyBoundaryMesh& patches,
162  const label comm
163 )
164 :
165  procNbrProcs_(Pstream::nProcs(comm)),
166  procPatchMap_(),
167  patchSchedule_(2*patches.size())
168 {
169  if (Pstream::parRun())
170  {
171  // Fill my 'slot' with my neighbours
172  procNbrProcs_[Pstream::myProcNo(comm)] =
173  procNeighbours(procNbrProcs_.size(), patches);
174 
175  // Distribute to all processors
176  Pstream::gatherList(procNbrProcs_, Pstream::msgType(), comm);
177  Pstream::scatterList(procNbrProcs_, Pstream::msgType(), comm);
178  }
179 
180  if
181  (
184  )
185  {
186  label patchEvali = 0;
187 
188  // 1. All non-processor patches
189  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
190 
192  {
193  if (!isA<processorPolyPatch>(patches[patchi]))
194  {
195  patchSchedule_[patchEvali].patch = patchi;
196  patchSchedule_[patchEvali++].init = true;
197  patchSchedule_[patchEvali].patch = patchi;
198  patchSchedule_[patchEvali++].init = false;
199  }
200  }
201 
202  // 2. All processor patches
203  // ~~~~~~~~~~~~~~~~~~~~~~~~
204 
205  // Determine the schedule for all. Insert processor pair once
206  // to determine the schedule. Each processor pair stands for both
207  // send and receive.
208  label nComms = 0;
209  forAll(procNbrProcs_, proci)
210  {
211  nComms += procNbrProcs_[proci].size();
212  }
213  DynamicList<labelPair> comms(nComms);
214 
215  forAll(procNbrProcs_, proci)
216  {
217  const labelList& nbrs = procNbrProcs_[proci];
218 
219  forAll(nbrs, i)
220  {
221  if (proci < nbrs[i])
222  {
223  comms.append(labelPair(proci, nbrs[i]));
224  }
225  }
226  }
227  comms.shrink();
228 
229  // Determine a schedule.
230  labelList mySchedule
231  (
233  (
234  Pstream::nProcs(comm),
235  comms
236  ).procSchedule()[Pstream::myProcNo(comm)]
237  );
238 
239  forAll(mySchedule, iter)
240  {
241  label commI = mySchedule[iter];
242 
243  // Get the other processor
244  label nb = comms[commI][0];
245  if (nb == Pstream::myProcNo(comm))
246  {
247  nb = comms[commI][1];
248  }
249  label patchi = procPatchMap_[nb];
250 
251  if (Pstream::myProcNo(comm) > nb)
252  {
253  patchSchedule_[patchEvali].patch = patchi;
254  patchSchedule_[patchEvali++].init = true;
255  patchSchedule_[patchEvali].patch = patchi;
256  patchSchedule_[patchEvali++].init = false;
257  }
258  else
259  {
260  patchSchedule_[patchEvali].patch = patchi;
261  patchSchedule_[patchEvali++].init = false;
262  patchSchedule_[patchEvali].patch = patchi;
263  patchSchedule_[patchEvali++].init = true;
264  }
265  }
266  }
267  else
268  {
269  patchSchedule_ = nonBlockingSchedule(patches);
270  }
271 }
272 
273 
274 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:252
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Inter-processor communications stream.
Definition: Pstream.H:56
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
Determines the order in which a set of processors should communicate with one another.
Definition: commSchedule.H:66
Foam::polyBoundaryMesh.
processorTopology(const polyBoundaryMesh &patches, const label comm)
Construct from boundaryMesh.
label patchi
const fvPatchList & patches
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
List< lduScheduleEntry > lduSchedule
Definition: lduSchedule.H:74
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)