solidBodyMotionFvMesh.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 | Copyright (C) 2011-2015 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 "solidBodyMotionFvMesh.H"
28 #include "volFields.H"
29 #include "transformField.H"
30 #include "cellZoneMesh.H"
31 #include "cellSet.H"
32 #include "boolList.H"
33 #include "syncTools.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(solidBodyMotionFvMesh, 0);
40  addToRunTimeSelectionTable(dynamicFvMesh, solidBodyMotionFvMesh, IOobject);
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 Foam::solidBodyMotionFvMesh::solidBodyMotionFvMesh(const IOobject& io)
47 :
48  dynamicFvMesh(io),
49  dynamicMeshCoeffs_
50  (
52  (
53  IOobject
54  (
55  "dynamicMeshDict",
56  io.time().constant(),
57  *this,
58  IOobject::MUST_READ_IF_MODIFIED,
59  IOobject::NO_WRITE,
60  false
61  )
62  ).subDict(typeName + "Coeffs")
63  ),
64  SBMFPtr_(solidBodyMotionFunction::New(dynamicMeshCoeffs_, io.time())),
65  undisplacedPoints_
66  (
67  IOobject
68  (
69  "points",
70  io.time().constant(),
71  meshSubDir,
72  *this,
73  IOobject::MUST_READ,
74  IOobject::NO_WRITE,
75  false
76  )
77  ),
78  pointIDs_(),
79  moveAllCells_(false),
80  UName_(dynamicMeshCoeffs_.lookupOrDefault<word>("UName", "U"))
81 {
82  if (undisplacedPoints_.size() != nPoints())
83  {
85  (
86  "solidBodyMotionFvMesh::solidBodyMotionFvMesh(const IOobject&)",
87  dynamicMeshCoeffs_
88  ) << "Read " << undisplacedPoints_.size()
89  << " undisplaced points from " << undisplacedPoints_.objectPath()
90  << " but the current mesh has " << nPoints()
91  << exit(FatalIOError);
92  }
93 
94  word cellZoneName =
95  dynamicMeshCoeffs_.lookupOrDefault<word>("cellZone", "none");
96 
97  word cellSetName =
98  dynamicMeshCoeffs_.lookupOrDefault<word>("cellSet", "none");
99 
100  if ((cellZoneName != "none") && (cellSetName != "none"))
101  {
103  (
104  "solidBodyMotionFvMesh::solidBodyMotionFvMesh(const IOobject&)",
105  dynamicMeshCoeffs_
106  )
107  << "Either cellZone OR cellSet can be supplied, but not both. "
108  << "If neither is supplied, all cells will be included"
109  << exit(FatalIOError);
110  }
111 
112 
113  labelList cellIDs;
114  if (cellZoneName != "none")
115  {
116  Info<< "Applying solid body motion to cellZone " << cellZoneName
117  << endl;
118 
119  label zoneID = cellZones().findZoneID(cellZoneName);
120 
121  if (zoneID == -1)
122  {
124  (
125  "solidBodyMotionFvMesh::solidBodyMotionFvMesh(const IOobject&)"
126  ) << "Unable to find cellZone " << cellZoneName
127  << ". Valid cellZones are:"
128  << cellZones().names()
129  << exit(FatalError);
130  }
131 
132  cellIDs = cellZones()[zoneID];
133  }
134 
135  if (cellSetName != "none")
136  {
137  Info<< "Applying solid body motion to cellSet " << cellSetName
138  << endl;
139 
140  cellSet set(*this, cellSetName);
141 
142  cellIDs = set.toc();
143  }
144 
145  label nCells = returnReduce(cellIDs.size(), sumOp<label>());
146  moveAllCells_ = nCells == 0;
147 
148  if (moveAllCells_)
149  {
150  Info<< "Applying solid body motion to entire mesh" << endl;
151  }
152  else
153  {
154  // collect point IDs of points in cell zone
155 
156  boolList movePts(nPoints(), false);
157 
158  forAll(cellIDs, i)
159  {
160  label cellI = cellIDs[i];
161  const cell& c = cells()[cellI];
162  forAll(c, j)
163  {
164  const face& f = faces()[c[j]];
165  forAll(f, k)
166  {
167  label pointI = f[k];
168  movePts[pointI] = true;
169  }
170  }
171  }
172 
173  syncTools::syncPointList(*this, movePts, orEqOp<bool>(), false);
174 
175  DynamicList<label> ptIDs(nPoints());
176  forAll(movePts, i)
177  {
178  if (movePts[i])
179  {
180  ptIDs.append(i);
181  }
182  }
183 
184  pointIDs_.transfer(ptIDs);
185  }
186 }
187 
188 
189 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
190 
192 {}
193 
194 
195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 
198 {
199  static bool hasWarned = false;
200 
201  if (moveAllCells_)
202  {
204  (
205  transform
206  (
207  SBMFPtr_().transformation(),
208  undisplacedPoints_
209  )
210  );
211  }
212  else
213  {
214  pointField transformedPts(undisplacedPoints_);
215 
216  UIndirectList<point>(transformedPts, pointIDs_) =
217  transform
218  (
219  SBMFPtr_().transformation(),
220  pointField(transformedPts, pointIDs_)
221  );
222 
223  fvMesh::movePoints(transformedPts);
224  }
225 
226 
227  if (foundObject<volVectorField>(UName_))
228  {
229  const volVectorField& U = lookupObject<volVectorField>(UName_);
230 
232  }
233  else if (!hasWarned)
234  {
235  hasWarned = true;
236 
237  WarningIn("solidBodyMotionFvMesh::update()")
238  << "Did not find volVectorField " << UName_
239  << " Not updating " << UName_ << "boundary conditions."
240  << endl;
241  }
242 
243  return true;
244 }
245 
246 
247 // ************************************************************************* //
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
Base class for defining solid-body motions.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:390
Foam::cellZoneMesh.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
labelList f(nPoints)
A collection of cell labels.
Definition: cellSet.H:48
A class for handling words, derived from string.
Definition: word.H:59
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:310
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label nCells() const
const cellList & cells() const
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Namespace for OpenFOAM.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
virtual bool update()
Update the mesh for both mesh motion and topology change.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOerror FatalIOError
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
#define WarningIn(functionName)
Report a warning using Foam::Warning.
A List with indirect addressing.
Definition: fvMatrix.H:106
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
#define forAll(list, i)
Definition: UList.H:421
label k
Boltzmann constant.
Spatial transformation functions for primitive fields.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:363
Macros for easy insertion into run-time selection tables.
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
error FatalError
Constant dispersed-phase particle diameter model.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:354
const dimensionedScalar c
Speed of light in a vacuum.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1060
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:269
Abstract base class for geometry and/or topology changing fvMesh.
Definition: dynamicFvMesh.H:51
U correctBoundaryConditions()
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
label nPoints() const
bool movePoints()
Do what is neccessary if the mesh has moved.
U
Definition: pEqn.H:82
defineTypeNameAndDebug(combustionModel, 0)