pointFieldDecomposer.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-2023 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 "pointFieldDecomposer.H"
27 #include "fvMesh.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 Foam::labelList Foam::pointFieldDecomposer::patchFieldDecomposer::addressing
32 (
33  const pointPatch& completePatch,
34  const pointPatch& procPatch,
35  const labelList& pointProcAddressing
36 )
37 {
38  const labelList& completePatchPoints = completePatch.meshPoints();
39  const labelList& procPatchPoints = procPatch.meshPoints();
40 
41  // Create a map from complete mesh point index to complete patch point index
42  labelList map(completePatch.boundaryMesh().mesh().size(), -1);
43  forAll(completePatchPoints, pointi)
44  {
45  map[completePatchPoints[pointi]] = pointi;
46  }
47 
48  // Determine the complete patch point for every proc patch point, going via
49  // the complete mesh point index and using the above map
50  labelList result(procPatch.size(), -1);
51  forAll(procPatchPoints, pointi)
52  {
53  result[pointi] = map[pointProcAddressing[procPatchPoints[pointi]]];
54  }
55 
56  // Check that all the patch point addresses are set
57  if (result.size() && min(result) < 0)
58  {
60  << "Incomplete patch point addressing"
61  << abort(FatalError);
62  }
63 
64  return result;
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
71 (
72  const pointPatch& completePatch,
73  const pointPatch& procPatch,
74  const labelList& pointProcAddressing
75 )
76 :
77  labelList(addressing(completePatch, procPatch, pointProcAddressing)),
78  forwardFieldMapper(static_cast<const labelList&>(*this))
79 {}
80 
81 
83 (
84  const pointMesh& completeMesh,
85  const PtrList<fvMesh>& procMeshes,
86  const labelListList& pointProcAddressing
87 )
88 :
89  completeMesh_(completeMesh),
90  procMeshes_(procMeshes),
91  pointProcAddressing_(pointProcAddressing),
92  patchFieldDecomposers_(procMeshes_.size())
93 {
94  forAll(procMeshes_, proci)
95  {
96  const pointMesh& procMesh = pointMesh::New(procMeshes_[proci]);
97 
98  patchFieldDecomposers_.set
99  (
100  proci,
101  new PtrList<patchFieldDecomposer>(procMesh.boundary().size())
102  );
103 
104  forAll(procMesh.boundary(), procPatchi)
105  {
106  if (procPatchi < completeMesh_.boundary().size())
107  {
108  patchFieldDecomposers_[proci].set
109  (
110  procPatchi,
111  new patchFieldDecomposer
112  (
113  completeMesh_.boundary()[procPatchi],
114  procMesh.boundary()[procPatchi],
115  pointProcAddressing_[proci]
116  )
117  );
118  }
119  }
120  }
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
125 
127 {}
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
132 bool Foam::pointFieldDecomposer::decomposes(const IOobjectList& objects)
133 {
134  bool result = false;
135 
136  #define DO_POINT_FIELDS_TYPE(Type, nullArg) \
137  result = result \
138  || !objects.lookupClass(PointField<Type>::typeName).empty();
139  FOR_ALL_FIELD_TYPES(DO_POINT_FIELDS_TYPE)
140  #undef DO_POINT_FIELDS_TYPE
141 
142  return result;
143 }
144 
145 
146 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Foam::tmp< Foam::Field< Type > > map(const Field< Type > &mapF) const
patchFieldDecomposer(const pointPatch &completePatch, const pointPatch &procPatch, const labelList &pointProcAddressing)
Construct given patches and addressing.
pointFieldDecomposer(const pointMesh &completeMesh, const PtrList< fvMesh > &procMeshes, const labelListList &pointAddressing)
Construct from components.
static bool decomposes(const IOobjectList &objects)
Return whether anything in the object list gets decomposed.
~pointFieldDecomposer()
Destructor.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
objects