multiSolidBodyMotionFvMesh.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-2016 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 
28 #include "volFields.H"
29 #include "transformField.H"
30 #include "cellZoneMesh.H"
31 #include "boolList.H"
32 #include "syncTools.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(multiSolidBodyMotionFvMesh, 0);
40  (
41  dynamicFvMesh,
42  multiSolidBodyMotionFvMesh,
43  IOobject
44  );
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 Foam::multiSolidBodyMotionFvMesh::multiSolidBodyMotionFvMesh(const IOobject& io)
51 :
52  dynamicFvMesh(io),
53  dynamicMeshCoeffs_
54  (
56  (
57  IOobject
58  (
59  "dynamicMeshDict",
60  io.time().constant(),
61  *this,
62  IOobject::MUST_READ_IF_MODIFIED,
63  IOobject::NO_WRITE,
64  false
65  )
66  ).subDict(typeName + "Coeffs")
67  ),
68  undisplacedPoints_
69  (
70  IOobject
71  (
72  "points",
73  io.time().constant(),
74  meshSubDir,
75  *this,
76  IOobject::MUST_READ,
77  IOobject::NO_WRITE,
78  false
79  )
80  )
81 {
82  if (undisplacedPoints_.size() != nPoints())
83  {
85  (
86  dynamicMeshCoeffs_
87  ) << "Read " << undisplacedPoints_.size()
88  << " undisplaced points from " << undisplacedPoints_.objectPath()
89  << " but the current mesh has " << nPoints()
90  << exit(FatalIOError);
91  }
92 
93 
94  zoneIDs_.setSize(dynamicMeshCoeffs_.size());
95  SBMFs_.setSize(dynamicMeshCoeffs_.size());
96  pointIDs_.setSize(dynamicMeshCoeffs_.size());
97  label zoneI = 0;
98 
99  forAllConstIter(dictionary, dynamicMeshCoeffs_, iter)
100  {
101  if (iter().isDict())
102  {
103  zoneIDs_[zoneI] = cellZones().findZoneID(iter().keyword());
104 
105  if (zoneIDs_[zoneI] == -1)
106  {
108  (
109  dynamicMeshCoeffs_
110  ) << "Cannot find cellZone named " << iter().keyword()
111  << ". Valid zones are " << cellZones().names()
112  << exit(FatalIOError);
113  }
114 
115  const dictionary& subDict = iter().dict();
116 
117  SBMFs_.set
118  (
119  zoneI,
120  solidBodyMotionFunction::New(subDict, io.time())
121  );
122 
123  // Collect points of cell zone.
124  const cellZone& cz = cellZones()[zoneIDs_[zoneI]];
125 
126  boolList movePts(nPoints(), false);
127 
128  forAll(cz, i)
129  {
130  label celli = cz[i];
131  const cell& c = cells()[celli];
132  forAll(c, j)
133  {
134  const face& f = faces()[c[j]];
135  forAll(f, k)
136  {
137  label pointi = f[k];
138  movePts[pointi] = true;
139  }
140  }
141  }
142 
143  syncTools::syncPointList(*this, movePts, orEqOp<bool>(), false);
144 
145  DynamicList<label> ptIDs(nPoints());
146  forAll(movePts, i)
147  {
148  if (movePts[i])
149  {
150  ptIDs.append(i);
151  }
152  }
153 
154  pointIDs_[zoneI].transfer(ptIDs);
155 
156  Info<< "Applying solid body motion " << SBMFs_[zoneI].type()
157  << " to " << pointIDs_[zoneI].size() << " points of cellZone "
158  << iter().keyword() << endl;
159 
160  zoneI++;
161  }
162  }
163  zoneIDs_.setSize(zoneI);
164  SBMFs_.setSize(zoneI);
165  pointIDs_.setSize(zoneI);
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
170 
172 {}
173 
174 
175 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176 
178 {
179  static bool hasWarned = false;
180 
181  pointField transformedPts(undisplacedPoints_);
182 
183  forAll(zoneIDs_, i)
184  {
185  const labelList& zonePoints = pointIDs_[i];
186 
187  UIndirectList<point>(transformedPts, zonePoints) =
189  (
190  SBMFs_[i].transformation(),
191  pointField(transformedPts, zonePoints)
192  );
193  }
194 
195  fvMesh::movePoints(transformedPts);
196 
197  if (foundObject<volVectorField>("U"))
198  {
199  const_cast<volVectorField&>(lookupObject<volVectorField>("U"))
201  }
202  else if (!hasWarned)
203  {
204  hasWarned = true;
205 
207  << "Did not find volVectorField U."
208  << " Not updating U boundary conditions." << endl;
209  }
210 
211  return true;
212 }
213 
214 
215 // ************************************************************************* //
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:363
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void transformPoints(vectorField &, const septernion &, const vectorField &)
Transform given vectorField of coordinates with the given septernion.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:602
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
bool movePoints()
Do what is neccessary if the mesh has moved.
label k
Boltzmann constant.
label size() const
Return number of elements in list.
Definition: DLListBaseI.H:77
const cellList & cells() const
Macros for easy insertion into run-time selection tables.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Spatial transformation functions for primitive fields.
Foam::cellZoneMesh.
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:341
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const Time & time() const
Return time.
Definition: IOobject.C:227
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:256
static autoPtr< solidBodyMotionFunction > New(const dictionary &SBMFCoeffs, const Time &runTime)
Select constructed from the SBMFCoeffs dictionary and Time.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
A subset of mesh cells.
Definition: cellZone.H:61
void setSize(const label)
Reset size of List.
Definition: List.C:295
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
#define WarningInFunction
Report a warning using Foam::Warning.
Constant dispersed-phase particle diameter model.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Abstract base class for geometry and/or topology changing fvMesh.
Definition: dynamicFvMesh.H:51
A List with indirect addressing.
Definition: fvMatrix.H:106
const dimensionedScalar c
Speed of light in a vacuum.
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
label nPoints() const
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:365
virtual bool update()
Update the mesh for both mesh motion and topology change.
Namespace for OpenFOAM.
IOerror FatalIOError