points0MotionSolver.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) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "points0MotionSolver.H"
30 #include "mapPolyMesh.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 }
38 
39 
40 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
41 
43 {
44  const word instance =
46  (
47  mesh.meshDir(),
48  "points0",
50  );
51  IOobject io
52  (
53  "points0",
54  instance,
56  mesh,
60  );
61 
62  // If points0 are located in constant directory, verify their existence
63  // or fallback to a copy of the original mesh points
64  if
65  (
66  instance == mesh.time().constant()
67  && !io.typeHeaderOk<pointIOField>()
68  )
69  {
70  io.rename("points");
71  }
72 
73  return io;
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
80 (
81  const polyMesh& mesh,
82  const IOdictionary& dict,
83  const word& type
84 )
85 :
87  zoneMotion(coeffDict(), mesh),
88  points0_(points0IO(mesh))
89 {
90  if
91  (
93  && (points0_.size() > mesh.nPoints())
94  )
95  {
96  // Allowed
97  }
98  else if (points0_.size() != mesh.nPoints())
99  {
100  const fileName fName
101  (
102  IOobject
103  (
104  "points",
105  time().constant(),
107  mesh,
111  ).typeFilePath<pointIOField>()
112  );
113 
115  << "Number of points in mesh " << mesh.nPoints()
116  << " differs from number of points " << points0_.size()
117  << " read from file " << fName << nl
118  << exit(FatalError);
119  }
120 }
121 
122 
124 (
125  const polyMesh& mesh,
126  const IOdictionary& dict,
127  const pointIOField& points0,
128  const word& type
129 )
130 :
131  motionSolver(mesh, dict, type),
132  zoneMotion(coeffDict(), mesh),
133  points0_(points0)
134 {
135  if (points0_.size() != mesh.nPoints())
136  {
138  << "Number of points in mesh " << mesh.nPoints()
139  << " differs from number of points " << points0_.size()
140  << " read from file " << points0.filePath()
142  }
143 }
144 
146 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147 
149 {}
150 
151 
153 {
154  // pointMesh already updates pointFields
155 
157 
158  // Map points0_. Bit special since we somehow have to come up with
159  // a sensible points0 position for introduced points.
160  // Find out scaling between points0 and current points
161 
162  // Get the new points either from the map or the mesh
163  const pointField& points =
164  (
165  mpm.hasMotionPoints()
166  ? mpm.preMotionPoints()
167  : mesh().points()
168  );
169 
170  // Note: boundBox does reduce
171  const vector span0 = boundBox(points0_).span();
172  const vector span = boundBox(points).span();
173 
174  vector scaleFactors(cmptDivide(span0, span));
175 
176  pointField newPoints0(mpm.pointMap().size());
177 
178  forAll(newPoints0, pointi)
179  {
180  label oldPointi = mpm.pointMap()[pointi];
181 
182  if (oldPointi >= 0)
183  {
184  label masterPointi = mpm.reversePointMap()[oldPointi];
185 
186  if (masterPointi == pointi)
187  {
188  newPoints0[pointi] = points0_[oldPointi];
189  }
190  else
191  {
192  // New point - assume motion is scaling
193  newPoints0[pointi] = points0_[oldPointi] + cmptMultiply
194  (
195  scaleFactors,
196  points[pointi] - points[masterPointi]
197  );
198  }
199  }
200  else
201  {
203  << "Cannot determine coordinates of introduced vertices."
204  << " New vertex " << pointi << " at coordinate "
205  << points[pointi] << exit(FatalError);
206  }
207  }
208 
209  twoDCorrectPoints(newPoints0);
210 
211  points0_.transfer(newPoints0);
212 
213  // points0 changed - set to write and check-in to database
214  points0_.rename("points0");
215  points0_.writeOpt(IOobject::AUTO_WRITE);
216  points0_.instance() = time().timeName();
217  points0_.checkIn();
218 }
219 
220 
221 // ************************************************************************* //
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:579
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
word findInstance(const fileName &dir, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null) const
Return time instance (location) of dir that contains the file name (eg, used in reading mesh data)...
Definition: Time.C:725
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
vectorIOField pointIOField
pointIOField is a vectorIOField.
Definition: pointIOField.H:38
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:410
Virtual base class for mesh motion solver.
Definition: motionSolver.H:54
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Ignore writing from objectRegistry::writeObject()
bool hasMotionPoints() const
Has valid preMotionPoints?
Definition: mapPolyMesh.H:764
Intermediate class for handling "zonified" motion. Can select motion for the entire mesh...
Definition: zoneMotion.H:66
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static IOobject points0IO(const polyMesh &mesh)
Return IO object for points0.
virtual void updateMesh(const mapPolyMesh &)=0
Update local data for topology changes.
Definition: motionSolver.C:206
dynamicFvMesh & mesh
const pointField & points
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
Reading is optional [identical to LAZY_READ].
points0MotionSolver(const points0MotionSolver &)=delete
No copy construct.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
defineTypeNameAndDebug(combustionModel, 0)
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:481
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:266
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))
Virtual base class for displacement motion solvers, where the point motion is relative to a set of fi...
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:841
pointIOField points0_
Starting points.
Automatically write from objectRegistry::writeObject()
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:192
pointField & points0() noexcept
Return reference to the reference (&#39;0&#39;) pointField.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:165
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
static bool allowConstructFromLargerSize
Permit read construct from a larger size.
Definition: Field.H:93
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
A primitive field of type <T> with automated input and output.
Do not request registration (bool: false)
const pointField & preMotionPoints() const
Pre-motion point positions.
Definition: mapPolyMesh.H:756
Namespace for OpenFOAM.