lduPrimitiveMeshAssemblyTemplates.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2019 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29
30#include "cyclicFvPatch.H"
31#include "cyclicAMIFvPatch.H"
32#include "cyclicACMIFvPatch.H"
33
35#include "AssemblyFvPatch.H"
36
37// * * * * * * * * * * * * * Public Member Functions * * * * * * * * * * * //
38
39template<class Type>
41(
43)
44{
45 label newFaces(0);
46 label newFacesProc(0);
47 label newPatches(0);
48
49 const label nMeshes(meshes_.size());
50
51 for (label i=0; i < nMeshes; ++i)
52 {
53 forAll(meshes_[i].interfaces(), patchI)
54 {
55 const polyPatch& pp = psis[i].mesh().boundaryMesh()[patchI];
56 const fvPatchField<Type>& fvp = psis[i].boundaryField()[patchI];
57
58 if (fvp.useImplicit())
59 {
60 if (pp.masterImplicit())
61 {
62 label newFacesPatch(0);
63 label newFacesProcPatch(0);
64
65 pp.newInternalProcFaces(newFacesPatch, newFacesProcPatch);
66
67 newFaces += newFacesPatch;
68 newFacesProc += newFacesProcPatch;
69
70 if (newFacesProc > 0)
71 {
73 << "The number of faces on either side of the "
74 << "coupled patch " << patchI << " are not "
75 << "the same. "
76 << "This might be due to the decomposition used. "
77 << "Please use decomposition preserving implicit "
78 << "patches on a single processor."
79 << exit(FatalError);
80 }
81
82 cellBoundMap_[i][patchI].setSize(newFacesPatch, -1);
83 facePatchFaceMap_[i][patchI].setSize(newFacesPatch, -1);
84 faceBoundMap_[i][patchI].setSize(newFacesPatch, -1);
85
86 const label nbrId = pp.neighbPolyPatchID();
87 const label meshNrbId = findNbrMeshId(pp, i);
88
89 cellBoundMap_[meshNrbId][nbrId].setSize
90 (
91 newFacesPatch,
92 -1
93 );
94 facePatchFaceMap_[meshNrbId][nbrId].setSize
95 (
96 newFacesPatch,
97 -1
98 );
99 faceBoundMap_[meshNrbId][nbrId].setSize
100 (
101 newFacesPatch,
102 -1
103 );
104 }
105 }
106 else
107 {
108 patchMap_[i][patchI] = newPatches;
109 patchLocalToGlobalMap_[i][patchI] = newPatches;
110 newPatches++;
111 }
112 }
113 }
114
115 label virtualPatches = newPatches;
116
117 // patchLocalToGlobalMap: map from original to asembled + extra Ids
118
119 for (label i=0; i < nMeshes; ++i)
120 {
121 forAll(meshes_[i].interfaces(), patchI)
122 {
123 if (patchLocalToGlobalMap_[i][patchI] == -1)
124 {
125 patchLocalToGlobalMap_[i][patchI] = virtualPatches++;
126 }
127 }
128 }
129
130 DebugInfo << patchMap_ << endl;
131 DebugInfo << patchLocalToGlobalMap_ << endl;
132
133 label oldFaces(0);
134 // Add the internal faces for each mesh
135 for (label i=0; i < nMeshes; ++i)
136 {
137 newFaces += meshes_[i].lduAddr().upperAddr().size();
138 oldFaces += meshes_[i].lduAddr().upperAddr().size();
139 }
140
141 if (debug)
142 {
143 Info<< " old total faces : " << oldFaces
144 << " new total faces (internal) : " << newFaces
145 << " new faces (remote) : " << newFacesProc
146 << " new Faces : " << newFaces - oldFaces << endl;
147
148 Info<< " new patches : " << newPatches << endl;
149
150 Info<< "Local to Global patch Map : "
151 << patchLocalToGlobalMap_ << endl;
152 }
153
154 // This gives the global cellId given the local patchId for interfaces
155 patchAddr_.setSize(newPatches);
156
157 for (label i=0; i < nMeshes; ++i)
158 {
159 const lduInterfacePtrsList interfacesLst = meshes_[i].interfaces();
160
161 forAll(interfacesLst, patchI)
162 {
163 label globalPatchId = patchMap_[i][patchI];
164
165 if (globalPatchId != -1)
166 {
167 const labelUList& faceCells =
168 meshes_[i].lduAddr().patchAddr(patchI);
169
170 // Fill local patchAddr for standard patches
171 if (!faceCells.empty())
172 {
173 patchAddr_[globalPatchId].setSize(faceCells.size(), -1);
174
175 for (label celli = 0; celli < faceCells.size(); ++celli)
176 {
177 patchAddr_[globalPatchId][celli] =
178 cellOffsets_[i] + faceCells[celli];
179 }
180 }
181 }
182 }
183 }
184
185 // Interfaces
186 interfaces().setSize(newPatches);
187 // Primitive interfaces
188 primitiveInterfaces().setSize(newPatches);
189
190 // The interfaces are conserved (cyclics, proc, etc)
191 for (label i=0; i < nMeshes; ++i)
192 {
193 const lduInterfacePtrsList interfacesLst = meshes_[i].interfaces();
194
195 forAll(interfacesLst, patchI)
196 {
197 label globalPatchId = patchMap_[i][patchI];
198 if (globalPatchId != -1)
199 {
200 if (interfacesLst.set(patchI))
201 {
202 const polyPatch& pp =
203 psis[i].mesh().boundaryMesh()[patchI];
204
205 const fvBoundaryMesh& bm = psis[i].mesh().boundary();
206
207 if (isA<cyclicLduInterface>(interfacesLst[patchI]))
208 {
209 label nbrId = refCast
210 <const cyclicLduInterface>
211 (
212 interfacesLst[patchI]
213 ).neighbPatchID();
214
215 label globalNbr = patchMap()[i][nbrId];
216
218 (
219 globalPatchId,
221 (
222 pp,
223 bm,
224 patchAddr_[globalNbr],
225 patchAddr_[globalPatchId],
226 globalNbr
227 )
228 );
229
231 (
232 globalPatchId,
233 &primitiveInterfaces()[globalPatchId]
234 );
235 }
236 else if
237 (
238 isA<cyclicAMILduInterface>(interfacesLst[patchI])
239 )
240 {
241 label nbrId =
242 refCast
243 <
245 >(pp).neighbPatchID();
246
247 label globalNbr = patchMap()[i][nbrId];
248
250 (
251 globalPatchId,
253 (
254 pp,
255 bm,
256 patchAddr_[globalNbr],
257 patchAddr_[globalPatchId],
258 globalNbr
259 )
260 );
262 (
263 globalPatchId,
264 &primitiveInterfaces()[globalPatchId]
265 );
266 }
267 else if
268 (
269 isA<cyclicACMILduInterface>(interfacesLst[patchI])
270 )
271 {
272 label nbrId =
273 refCast
274 <
276 >(pp).neighbPatchID();
277
278 label globalNbr = patchMap()[i][nbrId];
279
280 label nonOverlId =
281 refCast
282 <
284 >(pp).nonOverlapPatchID();
285
286 label globalnonOverlId = patchMap()[i][nonOverlId];
287
289 (
290 globalPatchId,
292 (
293 pp,
294 bm,
295 patchAddr_[globalNbr],
296 patchAddr_[globalPatchId],
297 globalNbr,
298 globalnonOverlId
299 )
300 );
302 (
303 globalPatchId,
304 &primitiveInterfaces()[globalPatchId]
305 );
306 }
307 else
308 {
310 (
311 globalPatchId,
312 nullptr
313 );
315 (
316 globalPatchId,
317 interfacesLst(patchI)
318 );
319 }
320 }
321 }
322 }
323 }
324
325 // Create new lower/upper addressing
326 lowerAddr().setSize(newFaces, -1);
327 upperAddr().setSize(newFaces, -1);
328
329 label startIndex = 0;
330
331 for (label i=0; i < nMeshes; ++i)
332 {
333 faceMap_[i].setSize(meshes_[i].lduAddr().lowerAddr().size(), -1);
334
335 const label nFaces = meshes_[i].lduAddr().upperAddr().size();
336
337 // Add individual addresses
338 SubList<label>(lowerAddr(), nFaces, startIndex) =
339 meshes_[i].lduAddr().lowerAddr();
340
341 SubList<label>(upperAddr(), nFaces, startIndex) =
342 meshes_[i].lduAddr().upperAddr();
343
344 // Offset cellsIDs to global cell addressing
345 label localFacei = 0;
346
347 for (label facei=startIndex; facei < startIndex + nFaces; ++facei)
348 {
349 lowerAddr()[facei] += cellOffsets_[i];
350 upperAddr()[facei] += cellOffsets_[i];
351
352 faceMap_[i][localFacei++] = facei;
353 }
354
355 startIndex += nFaces;
356 }
357 // Add new lower/upper adressing for new internal faces corresponding
358 // to patch faces that has a correspondent on the slave patch
359 // (i.e map, AMI,etc)
360 // Don't include faces that are in different proc
361
362 label nFaces = startIndex;
363
364 for (label i=0; i < nMeshes; ++i)
365 {
366 const lduInterfacePtrsList interfacesLst = meshes_[i].interfaces();
367
368 forAll(interfacesLst, patchI)
369 {
370 const polyPatch& pp = psis[i].mesh().boundaryMesh()[patchI];
371
372 const fvPatchField<Type>& fvp = psis[i].boundaryField()[patchI];
373
374 if (fvp.useImplicit())
375 {
376 label meshNrbId = this->findNbrMeshId(pp, i);
377
378 if (pp.masterImplicit())
379 {
380 const labelUList& nbrFaceCells = pp.nbrCells();
381 const label nbrPatchId = pp.neighbPolyPatchID();
382 refPtr<labelListList> tlocalFaceToFace =
384
385 const labelListList& localFaceToFace = tlocalFaceToFace();
386
387 // Compact target map
388 // key() = local face in proci
389 // *iter = compactId
390
391 label subFaceI(0);
392 forAll(pp.faceCells(), faceI)
393 {
394 const label cellI =
395 pp.faceCells()[faceI] + cellOffsets_[i];
396
397 const labelList& facesIds = localFaceToFace[faceI];
398
399 forAll(facesIds, j)
400 {
401 label nbrFaceId = facesIds[j];
402
403 // local faces
404 const label nbrCellI =
405 nbrFaceCells[nbrFaceId]
406 + cellOffsets_[meshNrbId];
407
408 lowerAddr()[nFaces] = min(cellI, nbrCellI);
409 upperAddr()[nFaces] = max(cellI, nbrCellI);
410
411 cellBoundMap_[i][patchI][subFaceI] = nbrCellI;
412 cellBoundMap_[meshNrbId][nbrPatchId][subFaceI] =
413 cellI;
414
415 facePatchFaceMap_[i][patchI][subFaceI] = faceI;
416 facePatchFaceMap_[meshNrbId][nbrPatchId][subFaceI]
417 = nbrFaceId;
418
419 faceBoundMap_[i][patchI][subFaceI] = nFaces;
420 faceBoundMap_[meshNrbId][nbrPatchId][subFaceI] =
421 nFaces;
422
423
424 ++subFaceI;
425 ++nFaces;
426 }
427 }
428 }
429 }
430 }
431 }
432
433 if (newFaces != nFaces)
434 {
436 << "Incorrrect total number of faces in the assembled lduMatrix: "
437 << newFaces << " != " << nFaces << nl
438 << exit(FatalError);
439 }
440
441 // Sort upper-tri order
442 {
443 labelList oldToNew
444 (
446 (
447 lduAddr().size(),
448 lowerAddr(),
449 upperAddr()
450 )
451 );
452 inplaceReorder(oldToNew, lowerAddr());
453 inplaceReorder(oldToNew, upperAddr());
454
455 for (labelList& faceMap : faceMap_)
456 {
457 for (label& facei : faceMap)
458 {
459 facei = oldToNew[facei];
460 }
461 }
462
463 for (labelListList& bMap : faceBoundMap_)
464 {
465 for (labelList& faceMap : bMap)
466 {
467 for (label& facei : faceMap)
468 {
469 if (facei != -1)
470 {
471 facei = oldToNew[facei];
472 }
473 }
474 }
475 }
476 }
477
478 if (debug & 2)
479 {
480 DebugVar(faceBoundMap_);
481 DebugVar(cellBoundMap_);
484 DebugVar(patchAddr_);
485 DebugVar(cellOffsets_);
486 DebugVar(faceMap_);
489 }
490}
491
492
493// ************************************************************************* //
An assembly of lduMatrix that is specific inter-region coupling through mapped patches.
Generic GeometricField class.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
A List obtained as a section of another List.
Definition: SubList.H:70
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
const T * set(const label i) const
Definition: UPtrList.H:248
void setSize(const label n)
Alias for resize()
Definition: UPtrList.H:261
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
Cyclic patch for Arbitrary Mesh Interface (AMI)
An abstract base class for cyclic coupled interfaces.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
Foam::fvBoundaryMesh.
const fvMesh & mesh() const noexcept
Return the mesh reference.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
bool useImplicit() const noexcept
Use implicit formulation for coupled patches only.
Definition: fvPatchField.H:314
virtual bool update()
Update the mesh for both mesh motion and topology change.
label size() const
Return number of equations.
label findNbrMeshId(const polyPatch &pp, const label iMesh) const
Find nrb mesh Id for mapped patches.
const labelListList & patchMap() const
Return patchMap.
const labelListList & faceMap() const
Return faceMap.
virtual const labelUList & upperAddr() const
Return Upper addressing.
static void checkUpperTriangular(const label size, const labelUList &l, const labelUList &u)
Check if in upper-triangular ordering.
PtrList< const lduInterface > & primitiveInterfaces()
Return a non-const list of primitive interfaces.
static labelList upperTriOrder(const label nCells, const labelUList &lower, const labelUList &upper)
Calculate upper-triangular order.
virtual lduInterfacePtrsList interfaces() const
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
virtual const labelUList & lowerAddr() const
Return Lower addressing.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
virtual void newInternalProcFaces(label &, label &) const
Return number of new internal of this polyPatch faces.
Definition: polyPatch.H:319
virtual label neighbPolyPatchID() const
Return nbr patchID.
Definition: polyPatch.H:332
virtual bool masterImplicit() const
Return implicit master.
Definition: polyPatch.H:346
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:315
virtual const labelUList & nbrCells() const
Return nbrCells.
Definition: polyPatch.H:325
virtual refPtr< labelListList > mapCollocatedFaces() const
Return mapped collocated faces.
Definition: polyPatch.H:339
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:371
A class for managing references or pointers (no reference counting)
Definition: refPtr.H:58
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
#define DebugInfo
Report an information message using Foam::Info.
#define DebugVar(var)
Report a variable name and value.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333