dynamicOversetFvMesh.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) 2014-2019 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify i
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 
28 #include "dynamicOversetFvMesh.H"
30 #include "cellCellStencilObject.H"
33 #include "globalIndex.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(dynamicOversetFvMesh, 0);
40  addToRunTimeSelectionTable(dynamicFvMesh, dynamicOversetFvMesh, IOobject);
41 }
42 
43 
44 // * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
45 
47 {
49 
50  // The (processor-local part of the) stencil determines the local
51  // faces to add to the matrix. tbd: parallel
52  const labelListList& stencil = overlap.cellStencil();
53 
54  // Get the base addressing
55  //const lduAddressing& baseAddr = dynamicMotionSolverFvMesh::lduAddr();
56  const lduAddressing& baseAddr = dynamicFvMesh::lduAddr();
57 
58  // Add to the base addressing
59  labelList lowerAddr;
60  labelList upperAddr;
61  label nExtraFaces;
62 
63  const globalIndex globalNumbering(baseAddr.size());
64  labelListList localFaceCells;
65  labelListList remoteFaceCells;
66 
67  labelList globalCellIDs(overlap.cellInterpolationMap().constructSize());
68  forAll(baseAddr, cellI)
69  {
70  globalCellIDs[cellI] = globalNumbering.toGlobal(cellI);
71  }
72  overlap.cellInterpolationMap().distribute(globalCellIDs);
73 
74 
76  (
77  baseAddr,
78  stencil,
79  nExtraFaces,
80  lowerAddr,
81  upperAddr,
83  globalNumbering,
84  globalCellIDs,
85  localFaceCells,
86  remoteFaceCells
87  );
88 
89  if (debug)
90  {
91  Pout<< "dynamicOversetFvMesh::update() : extended addressing from"
92  << " nFaces:" << baseAddr.lowerAddr().size()
93  << " to nFaces:" << lowerAddr.size()
94  << " nExtraFaces:" << nExtraFaces << endl;
95  }
96 
97 
98  // Now for the tricky bits. We want to hand out processor faces according
99  // to the localFaceCells/remoteFaceCells. Ultimately we need
100  // per entry in stencil:
101  // - the patch (or -1 for internal faces)
102  // - the face (is either an internal face index or a patch face index)
103 
105 
106  // Per processor to owner (local)/neighbour (remote)
108  List<DynamicList<label>> dynProcNeighbour(Pstream::nProcs());
109  forAll(stencil, celli)
110  {
111  const labelList& nbrs = stencil[celli];
112  stencilPatches_[celli].setSize(nbrs.size());
113  stencilPatches_[celli] = -1;
114 
115  forAll(nbrs, nbri)
116  {
117  if (stencilFaces_[celli][nbri] == -1)
118  {
119  const label nbrCelli = nbrs[nbri];
120  label globalNbr = globalCellIDs[nbrCelli];
121  label proci = globalNumbering.whichProcID(globalNbr);
122  label remoteCelli = globalNumbering.toLocal(proci, globalNbr);
123 
124  // Overwrite the face to be a patch face
125  stencilFaces_[celli][nbri] = procOwner[proci].size();
126  stencilPatches_[celli][nbri] = proci;
127  procOwner[proci].append(celli);
128  dynProcNeighbour[proci].append(remoteCelli);
129 
130  //Pout<< "From neighbour proc:" << proci
131  // << " allocating patchFace:" << stencilFaces_[celli][nbri]
132  // << " to get remote cell " << remoteCelli
133  // << " onto local cell " << celli << endl;
134  }
135  }
136  }
137  labelListList procNeighbour(dynProcNeighbour.size());
138  forAll(procNeighbour, i)
139  {
140  procNeighbour[i] = std::move(dynProcNeighbour[i]);
141  }
142  labelListList mySendCells;
143  Pstream::exchange<labelList, label>(procNeighbour, mySendCells);
144 
145  label nbri = 0;
146  forAll(procOwner, proci)
147  {
148  if (procOwner[proci].size())
149  {
150  nbri++;
151  }
152  if (mySendCells[proci].size())
153  {
154  nbri++;
155  }
156  }
157  remoteStencilInterfaces_.setSize(nbri);
158  nbri = 0;
159 
160  // E.g. if proc1 needs some data from proc2 and proc2 needs some data from
161  // proc1. We first want the pair : proc1 receive and proc2 send
162  // and then the pair : proc1 send, proc2 receive
163 
164 
165  labelList procToInterface(Pstream::nProcs(), -1);
166 
167  forAll(procOwner, proci)
168  {
169  if (proci < Pstream::myProcNo() && procOwner[proci].size())
170  {
171  if (debug)
172  {
173  Pout<< "Adding interface " << nbri
174  << " to receive my " << procOwner[proci].size()
175  << " from " << proci << endl;
176  }
177  procToInterface[proci] = nbri;
179  (
180  nbri++,
182  (
183  procOwner[proci],
185  proci,
186  tensorField(0),
187  Pstream::msgType()+2
188  )
189  );
190  }
191  else if (proci > Pstream::myProcNo() && mySendCells[proci].size())
192  {
193  if (debug)
194  {
195  Pout<< "Adding interface " << nbri
196  << " to send my " << mySendCells[proci].size()
197  << " to " << proci << endl;
198  }
200  (
201  nbri++,
203  (
204  mySendCells[proci],
206  proci,
207  tensorField(0),
208  Pstream::msgType()+2
209  )
210  );
211  }
212  }
213  forAll(procOwner, proci)
214  {
215  if (proci > Pstream::myProcNo() && procOwner[proci].size())
216  {
217  if (debug)
218  {
219  Pout<< "Adding interface " << nbri
220  << " to receive my " << procOwner[proci].size()
221  << " from " << proci << endl;
222  }
223  procToInterface[proci] = nbri;
225  (
226  nbri++,
228  (
229  procOwner[proci],
231  proci,
232  tensorField(0),
233  Pstream::msgType()+3
234  )
235  );
236  }
237  else if (proci < Pstream::myProcNo() && mySendCells[proci].size())
238  {
239  if (debug)
240  {
241  Pout<< "Adding interface " << nbri
242  << " to send my " << mySendCells[proci].size()
243  << " to " << proci << endl;
244  }
246  (
247  nbri++,
249  (
250  mySendCells[proci],
252  proci,
253  tensorField(0),
254  Pstream::msgType()+3
255  )
256  );
257  }
258  }
259 
260 
261  // Rewrite stencilPatches now we know the actual interface (procToInterface)
262  for (auto& patches : stencilPatches_)
263  {
264  for (auto& interface : patches)
265  {
266  if (interface != -1)
267  {
268  interface = procToInterface[interface]+boundary().size();
269  }
270  }
271  }
272 
273 
274  // Get addressing and interfaces of all interfaces
275 
276 
277  List<const labelUList*> patchAddr;
278  {
279  const fvBoundaryMesh& fvp = boundary();
280 
281  patchAddr.setSize(fvp.size() + remoteStencilInterfaces_.size());
282 
283  //allInterfaces_ = dynamicMotionSolverFvMesh::interfaces();
285  allInterfaces_.setSize(patchAddr.size());
286 
287  forAll(fvp, patchI)
288  {
289  patchAddr[patchI] = &fvp[patchI].faceCells();
290  }
292  {
293  label patchI = fvp.size()+i;
296 
297  //Pout<< "at patch:" << patchI
298  // << " have procPatch:" << pp.type()
299  // << " from:" << pp.myProcNo()
300  // << " to:" << pp.neighbProcNo()
301  // << " with fc:" << pp.faceCells().size() << endl;
302 
303  patchAddr[patchI] = &pp.faceCells();
304  allInterfaces_.set(patchI, &pp);
305  }
306  }
307  const lduSchedule ps
308  (
309  lduPrimitiveMesh::nonBlockingSchedule<processorLduInterface>
310  (
312  )
313  );
314 
315  lduPtr_.reset
316  (
318  (
319  nCells(),
320  std::move(lowerAddr),
321  std::move(upperAddr),
322  patchAddr,
323  ps
324  )
325  );
326 
327 
328  // Check
329  if (debug)
330  {
331  const lduAddressing& addr = lduPtr_(); //this->lduAddr();
332 
333  Pout<< "Adapted addressing:"
334  << " lower:" << addr.lowerAddr().size()
335  << " upper:" << addr.upperAddr().size() << endl;
336 
337  // Using lduAddressing::patch
338  forAll(patchAddr, patchI)
339  {
340  Pout<< " " << patchI << "\tpatchAddr:"
341  << addr.patchAddr(patchI).size()
342  << endl;
343  }
344 
345  // Using interfaces
346  const lduInterfacePtrsList& iFaces = allInterfaces_;
347  Pout<< "Adapted interFaces:" << iFaces.size() << endl;
348  forAll(iFaces, patchI)
349  {
350  if (iFaces.set(patchI))
351  {
352  Pout<< " " << patchI << "\tinterface:"
353  << iFaces[patchI].type() << endl;
354  }
355  }
356  }
357 
358  return true;
359 }
360 
361 
363 (
364  const labelList& types,
365  const labelList& nbrTypes,
366  const scalarField& norm,
367  const scalarField& nbrNorm,
368  const label celli,
369  bitSet& isFront
370 ) const
371 {
372  const labelList& own = faceOwner();
373  const labelList& nei = faceNeighbour();
374  const cell& cFaces = cells()[celli];
375 
376  scalar avg = 0.0;
377  label n = 0;
378  label nFront = 0;
379  for (const label facei : cFaces)
380  {
381  if (isInternalFace(facei))
382  {
383  label nbrCelli = (own[facei] == celli ? nei[facei] : own[facei]);
384  if (norm[nbrCelli] == -GREAT)
385  {
386  // Invalid neighbour. Add to front
387  if (isFront.set(facei))
388  {
389  nFront++;
390  }
391  }
392  else
393  {
394  // Valid neighbour. Add to average
395  avg += norm[nbrCelli];
396  n++;
397  }
398  }
399  else
400  {
401  if (nbrNorm[facei-nInternalFaces()] == -GREAT)
402  {
403  if (isFront.set(facei))
404  {
405  nFront++;
406  }
407  }
408  else
409  {
410  avg += nbrNorm[facei-nInternalFaces()];
411  n++;
412  }
413  }
414  }
415 
416  if (n > 0)
417  {
418  return avg/n;
419  }
420  else
421  {
422  return norm[celli];
423  }
424 }
425 
426 
428 (
429  const GAMGAgglomeration& agglom
430 ) const
431 {
432  labelList cellToCoarse(identity(nCells()));
433  labelListList coarseToCell(invertOneToMany(nCells(), cellToCoarse));
434 
435  // Write initial agglomeration
436  {
437  volScalarField scalarAgglomeration
438  (
439  IOobject
440  (
441  "agglomeration",
442  this->time().timeName(),
443  *this,
446  false
447  ),
448  *this,
450  );
451  scalarField& fld = scalarAgglomeration.primitiveFieldRef();
452  forAll(fld, celli)
453  {
454  fld[celli] = cellToCoarse[celli];
455  }
456  fld /= max(fld);
458  <
461  >(scalarAgglomeration.boundaryFieldRef(), false);
462  scalarAgglomeration.write();
463 
464  Info<< "Writing initial cell distribution to "
465  << this->time().timeName() << endl;
466  }
467 
468 
469  for (label level = 0; level < agglom.size(); level++)
470  {
471  const labelList& addr = agglom.restrictAddressing(level);
472  label coarseSize = max(addr)+1;
473 
474  Info<< "Level : " << level << endl
475  << returnReduce(addr.size(), sumOp<label>()) << endl
476  << " current size : "
477  << returnReduce(addr.size(), sumOp<label>()) << endl
478  << " agglomerated size : "
479  << returnReduce(coarseSize, sumOp<label>()) << endl;
480 
481  forAll(addr, fineI)
482  {
483  const labelList& cellLabels = coarseToCell[fineI];
484  forAll(cellLabels, i)
485  {
486  cellToCoarse[cellLabels[i]] = addr[fineI];
487  }
488  }
489  coarseToCell = invertOneToMany(coarseSize, cellToCoarse);
490 
491  // Write agglomeration
492  {
493  volScalarField scalarAgglomeration
494  (
495  IOobject
496  (
497  "agglomeration_" + Foam::name(level),
498  this->time().timeName(),
499  *this,
502  false
503  ),
504  *this,
506  );
507  scalarField& fld = scalarAgglomeration.primitiveFieldRef();
508  forAll(fld, celli)
509  {
510  fld[celli] = cellToCoarse[celli];
511  }
512  //if (normalise)
513  //{
514  // fld /= max(fld);
515  //}
517  <
520  >(scalarAgglomeration.boundaryFieldRef(), false);
521  scalarAgglomeration.write();
522  }
523  }
524 }
525 
526 
527 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
528 
529 Foam::dynamicOversetFvMesh::dynamicOversetFvMesh(const IOobject& io)
530 :
532  active_(false)
533 {
534  // Load stencil (but do not update)
535  (void)Stencil::New(*this, false);
536 }
537 
538 
539 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
540 
542 {}
543 
544 
545 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
546 
548 {
549  if (!active_)
550  {
551  //return dynamicMotionSolverFvMesh::lduAddr();
552  return dynamicFvMesh::lduAddr();
553  }
554  if (lduPtr_.empty())
555  {
556  // Build extended addressing
557  updateAddressing();
558  }
559  return *lduPtr_;
560 }
561 
562 
564 {
565  if (!active_)
566  {
567  //return dynamicMotionSolverFvMesh::interfaces();
568  return dynamicFvMesh::interfaces();
569  }
570  if (lduPtr_.empty())
571  {
572  // Build extended addressing
573  updateAddressing();
574  }
575  return allInterfaces_;
576 }
577 
578 
581 {
582  if (lduPtr_.empty())
583  {
585  << "Extended addressing not allocated" << abort(FatalError);
586  }
587 
588  return *lduPtr_;
589 }
590 
591 
593 {
594  //if (dynamicMotionSolverFvMesh::update())
596  {
597  // Calculate the local extra faces for the interpolation. Note: could
598  // let demand-driven lduAddr() trigger it but just to make sure.
599  updateAddressing();
600 
601  // Addressing and/or weights have changed. Make interpolated cells
602  // up to date with donors
603  interpolateFields();
604 
605  return true;
606  }
607 
608  return false;
609 }
610 
611 
613 {
614  if (name.ends_with("_0"))
615  {
616  return baseName(name.substr(0, name.size()-2));
617  }
618 
619  return name;
620 }
621 
622 
624 {
625  // Add the stencil suppression list
626  wordHashSet suppressed(Stencil::New(*this).nonInterpolatedFields());
627 
628  // Use whatever the solver has set up as suppression list
629  const dictionary* dictPtr
630  (
631  this->schemesDict().findDict("oversetInterpolationSuppressed")
632  );
633  if (dictPtr)
634  {
635  suppressed.insert(dictPtr->toc());
636  }
637 
638  interpolate<volScalarField>(suppressed);
639  interpolate<volVectorField>(suppressed);
640  interpolate<volSphericalTensorField>(suppressed);
641  interpolate<volSymmTensorField>(suppressed);
642  interpolate<volTensorField>(suppressed);
643 
644  return true;
645 }
646 
647 
648 
650 (
654  const bool valid
655 ) const
656 {
657  //bool ok = dynamicMotionSolverFvMesh::writeObject(fmt, ver, cmp, valid);
658  bool ok = dynamicMotionSolverListFvMesh::writeObject(fmt, ver, cmp, valid);
659 
660  // For postprocessing : write cellTypes and zoneID
661  {
663 
665 
666  volScalarField volTypes
667  (
668  IOobject
669  (
670  "cellTypes",
671  this->time().timeName(),
672  *this,
675  false
676  ),
677  *this,
679  zeroGradientFvPatchScalarField::typeName
680  );
681 
682  forAll(volTypes.internalField(), cellI)
683  {
684  volTypes[cellI] = cellTypes[cellI];
685  }
686  volTypes.correctBoundaryConditions();
687  volTypes.writeObject(fmt, ver, cmp, valid);
688  }
689  {
690  volScalarField volZoneID
691  (
692  IOobject
693  (
694  "zoneID",
695  this->time().timeName(),
696  *this,
699  false
700  ),
701  *this,
703  zeroGradientFvPatchScalarField::typeName
704  );
705 
707  const labelIOList& zoneID = overlap.zoneID();
708 
709  forAll(zoneID, cellI)
710  {
711  volZoneID[cellI] = zoneID[cellI];
712  }
713  volZoneID.correctBoundaryConditions();
714  volZoneID.writeObject(fmt, ver, cmp, valid);
715  }
716  if (debug)
717  {
719  const labelIOList& zoneID = overlap.zoneID();
720  const labelListList& cellStencil = overlap.cellStencil();
721 
722  // Get remote zones
723  labelList donorZoneID(zoneID);
724  overlap.cellInterpolationMap().distribute(donorZoneID);
725 
726  // Get remote cellCentres
727  pointField cc(C());
729 
730  volScalarField volDonorZoneID
731  (
732  IOobject
733  (
734  "donorZoneID",
735  this->time().timeName(),
736  *this,
739  false
740  ),
741  *this,
742  dimensionedScalar("minOne", dimless, scalar(-1)),
743  zeroGradientFvPatchScalarField::typeName
744  );
745 
746  forAll(cellStencil, cellI)
747  {
748  const labelList& stencil = cellStencil[cellI];
749  if (stencil.size())
750  {
751  volDonorZoneID[cellI] = donorZoneID[stencil[0]];
752  for (label i = 1; i < stencil.size(); i++)
753  {
754  if (donorZoneID[stencil[i]] != volDonorZoneID[cellI])
755  {
756  WarningInFunction << "Mixed donor meshes for cell "
757  << cellI << " at " << C()[cellI]
758  << " donors:" << UIndirectList<point>(cc, stencil)
759  << endl;
760  volDonorZoneID[cellI] = -2;
761  }
762  }
763  }
764  }
765  //- Do not correctBoundaryConditions since re-interpolates!
766  //volDonorZoneID.correctBoundaryConditions();
767  volDonorZoneID.writeObject(fmt, ver, cmp, valid);
768  }
769 
770  return ok;
771 }
772 
773 
774 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::lduAddressing
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Definition: lduAddressing.H:114
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::HashTable< regIOobject * >::size
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Foam::dynamicOversetFvMesh::stencilPatches_
labelListList stencilPatches_
Corresponding patches (in above lduPtr) to the stencil.
Definition: dynamicOversetFvMesh.H:86
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::UPtrList::size
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:90
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::lduPrimitiveProcessorInterface
Concrete implementation of processor interface. Used to temporarily store settings.
Definition: lduPrimitiveProcessorInterface.H:53
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:64
Foam::GAMGAgglomeration
Geometric agglomerated algebraic multigrid agglomeration class.
Definition: GAMGAgglomeration.H:64
Foam::dynamicOversetFvMesh::interfaces
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: dynamicOversetFvMesh.C:563
Foam::dynamicOversetFvMesh::update
virtual bool update()
Update the mesh for both mesh motion and topology change.
Definition: dynamicOversetFvMesh.C:592
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dynamicOversetFvMesh::reverseFaceMap_
labelList reverseFaceMap_
From old to new face labels.
Definition: dynamicOversetFvMesh.H:89
Foam::cellCellStencil::zoneID
static const labelIOList & zoneID(const fvMesh &)
Helper: get reference to registered zoneID. Loads volScalarField.
Definition: cellCellStencil.C:99
Foam::fvMesh::writeObject
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool valid) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:914
interface
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
globalIndex.H
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:426
Foam::cellCellStencilObject::cellTypes
virtual const labelUList & cellTypes() const
Return the cell type list.
Definition: cellCellStencilObject.H:121
Foam::MeshObject::New
static const Type & New(const Mesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::dynamicOversetFvMesh::writeAgglomeration
void writeAgglomeration(const GAMGAgglomeration &agglom) const
Debug: dump agglomeration.
Definition: dynamicOversetFvMesh.C:428
Foam::tensorField
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Definition: primitiveFieldsFwd.H:57
Foam::bitSet::set
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:563
dynamicOversetFvMesh.H
Foam::UPtrList::set
const T * set(const label i) const
Return const pointer to element (if set) or nullptr.
Definition: UPtrListI.H:176
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::GeometricField::internalField
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
Definition: GeometricFieldI.H:43
Foam::HashSet< word >
Foam::dynamicOversetFvMesh::lduPtr_
autoPtr< fvMeshPrimitiveLduAddressing > lduPtr_
Extended addressing (extended with local interpolation stencils)
Definition: dynamicOversetFvMesh.H:71
Foam::lduAddressing::size
label size() const
Return number of equations.
Definition: lduAddressing.H:171
Foam::lduAddressing::upperAddr
virtual const labelUList & upperAddr() const =0
Return upper addressing.
Foam::GAMGAgglomeration::size
label size() const
Definition: GAMGAgglomeration.H:325
C
volScalarField & C
Definition: readThermalProperties.H:102
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::lduPrimitiveProcessorInterface::faceCells
virtual const labelUList & faceCells() const
Return faceCell addressing.
Definition: lduPrimitiveProcessorInterface.H:119
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
correctBoundaryConditions
cellMask correctBoundaryConditions()
Foam::dynamicOversetFvMesh::remoteStencilInterfaces_
PtrList< const lduPrimitiveProcessorInterface > remoteStencilInterfaces_
Added (processor)lduInterfaces for remote bits of stencil.
Definition: dynamicOversetFvMesh.H:76
Foam::fvSchemes::debug
static int debug
Debug switch.
Definition: fvSchemes.H:105
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::GAMGAgglomeration::restrictAddressing
const labelField & restrictAddressing(const label leveli) const
Return cell restrict addressing of given level.
Definition: GAMGAgglomeration.H:343
Foam::fvMeshPrimitiveLduAddressing::addAddressing
static labelList addAddressing(const lduAddressing &addr, const labelListList &nbrCells, label &nExtraFaces, labelList &lower, labelList &upper, labelListList &nbrCellFaces, const globalIndex &, const labelList &globalCellIDs, labelListList &localFaceCells, labelListList &remoteFaceCells)
Given additional addressing (in the form of additional neighbour.
Definition: fvMeshPrimitiveLduAddressing.C:109
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::label
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:62
Foam::IOstreamOption::versionNumber
Representation of a major/minor version number.
Definition: IOstreamOption.H:79
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fvMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:595
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fvMeshPrimitiveLduAddressing
Variant of fvMeshLduAddressing that contains addressing instead of slices.
Definition: fvMeshPrimitiveLduAddressing.H:57
Foam::UPtrList< const lduInterface >
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::mapDistribute::distribute
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Definition: mapDistributeTemplates.C:152
Foam::dynamicOversetFvMesh::writeObject
virtual bool writeObject(IOstream::streamFormat, IOstream::versionNumber, IOstream::compressionType, const bool valid) const
Write using given format, version and compression.
Definition: dynamicOversetFvMesh.C:650
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::dynamicOversetFvMesh::interpolateFields
virtual bool interpolateFields()
Update fields when mesh is updated.
Definition: dynamicOversetFvMesh.C:623
Foam::IOstreamOption::streamFormat
streamFormat
Data format (ascii | binary)
Definition: IOstreamOption.H:64
Foam::dynamicOversetFvMesh::~dynamicOversetFvMesh
virtual ~dynamicOversetFvMesh()
Destructor.
Definition: dynamicOversetFvMesh.C:541
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dynamicOversetFvMesh::updateAddressing
virtual bool updateAddressing() const
Calculate the extended lduAddressing.
Definition: dynamicOversetFvMesh.C:46
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::cellCellStencilObject::cellInterpolationMap
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
Definition: cellCellStencilObject.H:133
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:735
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:444
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:909
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:491
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::dynamicOversetFvMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing. If active: is (extended)
Definition: dynamicOversetFvMesh.C:547
cellTypes
const labelList & cellTypes
Definition: setCellMask.H:33
Foam::dynamicOversetFvMesh::allInterfaces_
lduInterfacePtrsList allInterfaces_
Interfaces for above mesh. Contains both original and.
Definition: dynamicOversetFvMesh.H:80
Foam::cellCellStencilObject::cellStencil
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
Definition: cellCellStencilObject.H:140
Foam::dynamicOversetFvMesh::stencilFaces_
labelListList stencilFaces_
Corresponding faces (in above lduPtr) to the stencil.
Definition: dynamicOversetFvMesh.H:83
Foam::lduAddressing::lowerAddr
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
Foam::List< labelList >
Foam::string::ends_with
bool ends_with(const std::string &s) const
True if string ends with the given suffix (cf. C++20)
Definition: string.H:301
Foam::fvMesh::interfaces
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.H:279
Foam::dynamicOversetFvMesh::primitiveLduAddr
const fvMeshPrimitiveLduAddressing & primitiveLduAddr() const
Return extended ldu addressing.
Definition: dynamicOversetFvMesh.C:580
Foam::UList< label >
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::IOList< label >
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::lduAddressing::patchAddr
virtual const labelUList & patchAddr(const label patchNo) const =0
Return patch to internal addressing given patch number.
Foam::oversetFvPatchField
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
Definition: oversetFvPatchField.H:56
Foam::UList::size
void size(const label n) noexcept
Override size to be inconsistent with allocated storage.
Definition: UListI.H:360
Foam::dynamicMotionSolverListFvMesh
Dynamic mesh able to handle multiple motion solvers. NOTE: If the word entry "solvers" is not found i...
Definition: dynamicMotionSolverListFvMesh.H:54
Foam::IOstreamOption::compressionType
compressionType
Compression treatment (UNCOMPRESSED | COMPRESSED)
Definition: IOstreamOption.H:71
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:109
lduPrimitiveProcessorInterface.H
Foam::dynamicMotionSolverListFvMesh::update
virtual bool update()
Dummy update function which does not change the mesh.
Definition: dynamicMotionSolverListFvMesh.C:118
overlap
const cellCellStencilObject & overlap
Definition: correctPhi.H:57
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::dynamicOversetFvMesh::baseName
static word baseName(const word &name)
Helper: strip off trailing _0.
Definition: dynamicOversetFvMesh.C:612
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:113
cellCellStencilObject.H
Foam::dynamicOversetFvMesh::cellAverage
scalar cellAverage(const labelList &types, const labelList &nbrTypes, const scalarField &norm, const scalarField &nbrNorm, const label celli, bitSet &isFront) const
Average norm of valid neighbours.
Definition: dynamicOversetFvMesh.C:363
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
zeroGradientFvPatchFields.H
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::cellCellStencilObject
Definition: cellCellStencilObject.H:57
Foam::UPtrList::setSize
void setSize(const label newLen)
Same as resize()
Definition: UPtrListI.H:160