solidBodyMotionSolver.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) 2016-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 "solidBodyMotionSolver.H"
28 #include "transformField.H"
29 #include "meshCellZones.H"
30 #include "cellSet.H"
31 #include "boolList.H"
32 #include "syncTools.H"
33 #include "polyTopoChangeMap.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(solidBodyMotionSolver, 0);
41  (
42  motionSolver,
43  solidBodyMotionSolver,
44  dictionary
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& name,
54  const polyMesh& mesh,
55  const dictionary& dict
56 )
57 :
58  points0MotionSolver(name, mesh, dict, typeName),
59  SBMFPtr_(solidBodyMotionFunction::New(coeffDict(), mesh.time())),
60  pointIDs_(),
61  moveAllCells_(false),
62  transform_(SBMFPtr_().transformation())
63 {
64  word cellZoneName =
65  coeffDict().lookupOrDefault<word>("cellZone", "none");
66 
67  word cellSetName =
68  coeffDict().lookupOrDefault<word>("cellSet", "none");
69 
70  if ((cellZoneName != "none") && (cellSetName != "none"))
71  {
72  FatalIOErrorInFunction(coeffDict())
73  << "Either cellZone OR cellSet can be supplied, but not both. "
74  << "If neither is supplied, all cells will be included"
75  << exit(FatalIOError);
76  }
77 
78  labelList cellIDs;
79  if (cellZoneName != "none")
80  {
81  Info<< "Applying solid body motion to cellZone " << cellZoneName
82  << endl;
83 
84  label zoneID = mesh.cellZones().findZoneID(cellZoneName);
85 
86  if (zoneID == -1)
87  {
89  << "Unable to find cellZone " << cellZoneName
90  << ". Valid cellZones are:"
91  << mesh.cellZones().names()
92  << exit(FatalError);
93  }
94 
95  cellIDs = mesh.cellZones()[zoneID];
96  }
97 
98  if (cellSetName != "none")
99  {
100  Info<< "Applying solid body motion to cellSet " << cellSetName
101  << endl;
102 
103  cellSet set(mesh, cellSetName);
104 
105  cellIDs = set.toc();
106  }
107 
108  label nCells = returnReduce(cellIDs.size(), sumOp<label>());
109  moveAllCells_ = nCells == 0;
110 
111  if (moveAllCells_)
112  {
113  Info<< "Applying solid body motion to entire mesh" << endl;
114  }
115  else
116  {
117  // collect point IDs of points in cell zone
118 
119  boolList movePts(mesh.nPoints(), false);
120 
121  forAll(cellIDs, i)
122  {
123  label celli = cellIDs[i];
124  const cell& c = mesh.cells()[celli];
125  forAll(c, j)
126  {
127  const face& f = mesh.faces()[c[j]];
128  forAll(f, k)
129  {
130  label pointi = f[k];
131  movePts[pointi] = true;
132  }
133  }
134  }
135 
136  syncTools::syncPointList(mesh, movePts, orEqOp<bool>(), false);
137 
138  DynamicList<label> ptIDs(mesh.nPoints());
139  forAll(movePts, i)
140  {
141  if (movePts[i])
142  {
143  ptIDs.append(i);
144  }
145  }
146 
147  pointIDs_.transfer(ptIDs);
148  }
149 }
150 
151 
152 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
153 
155 {}
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
161 {
162  transform_ = SBMFPtr_().transformation();
163 
164  if (moveAllCells_)
165  {
166  return transformPoints(transform_, points0_);
167  }
168  else
169  {
170  tmp<pointField> ttransformedPts(new pointField(mesh().points()));
171  pointField& transformedPts = ttransformedPts.ref();
172 
173  UIndirectList<point>(transformedPts, pointIDs_) = transformPoints
174  (
175  transform_,
176  pointField(points0_, pointIDs_)
177  );
178 
179  return ttransformedPts;
180  }
181 }
182 
183 
185 {
186  // pointMesh already updates pointFields
187 
188  // Get the new points either from the map or the mesh
189  const pointField& points =
190  (
191  map.hasMotionPoints()
192  ? map.preMotionPoints()
193  : mesh().points()
194  );
195 
196  pointField newPoints0(map.pointMap().size());
197 
198  forAll(newPoints0, pointi)
199  {
200  label oldPointi = map.pointMap()[pointi];
201 
202  if (oldPointi >= 0)
203  {
204  const label masterPointi = map.reversePointMap()[oldPointi];
205 
206  if (masterPointi == pointi)
207  {
208  newPoints0[pointi] = points0_[oldPointi];
209  }
210  else
211  {
212  newPoints0[pointi] =
213  transform_.invTransformPoint(points[pointi]);
214  }
215  }
216  else
217  {
219  << "Cannot determine co-ordinates of introduced vertices."
220  << " New vertex " << pointi << " at co-ordinate "
221  << points[pointi] << exit(FatalError);
222  }
223  }
224 
225  twoDCorrectPoints(newPoints0);
226 
227  points0_.transfer(newPoints0);
228 }
229 
230 
231 // ************************************************************************* //
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: MeshZones.C:341
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:501
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
bool hasMotionPoints() const
Has valid preMotionPoints?
const labelList & pointMap() const
Old point map.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const cellList & cells() const
label k
Boltzmann constant.
fvMesh & mesh
const dimensionedScalar c
Speed of light in a vacuum.
Macros for easy insertion into run-time selection tables.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
wordList names() const
Return a list of zone names.
Definition: MeshZones.C:256
Spatial transformation functions for primitive fields.
solidBodyMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from mesh and dictionary.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const pointField & points
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
static autoPtr< solidBodyMotionFunction > New(const dictionary &SBMFCoeffs, const Time &runTime)
Select constructed from the SBMFCoeffs dictionary and Time.
const labelList & reversePointMap() const
Reverse point map.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
const pointField & preMotionPoints() const
Pre-motion point positions.
labelList f(nPoints)
Foam::meshCellZones.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A collection of cell labels.
Definition: cellSet.H:48
Virtual base class for displacement motion solvers.
A List with indirect addressing.
Definition: fvMatrix.H:106
void transformPoints(vectorField &, const spatialTransform &, const vectorField &)
Transform given vectorField of coordinates with the given spatialTransform.
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
messageStream Info
label nPoints() const
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
IOerror FatalIOError
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.