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 | Copyright (C) 2016-2017 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 "points0MotionSolver.H"
27 #include "mapPolyMesh.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(points0MotionSolver, 0);
34 }
35 
36 
37 // * * * * * * * * * * * * * Protected Data Members * * * * * * * * * * * * * //
38 
40 (
41  const polyMesh& mesh
42 ) const
43 {
44  const word instance =
45  time().findInstance
46  (
47  mesh.meshDir(),
48  "points0",
50  );
51 
52  if (instance != time().constant())
53  {
54  // points0 written to a time folder
55 
56  return
57  IOobject
58  (
59  "points0",
60  instance,
62  mesh,
65  false
66  );
67  }
68  else
69  {
70  // check that points0 are actually in constant directory
71 
72  IOobject io
73  (
74  "points0",
75  instance,
77  mesh,
80  false
81  );
82 
83  if (io.typeHeaderOk<pointIOField>())
84  {
85  return io;
86  }
87  else
88  {
89  // copy of original mesh points
90 
91  return
92  IOobject
93  (
94  "points",
95  instance,
97  mesh,
100  false
101  );
102  }
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108 
109 Foam::points0MotionSolver::points0MotionSolver
110 (
111  const polyMesh& mesh,
112  const IOdictionary& dict,
113  const word& type
114 )
115 :
116  motionSolver(mesh, dict, type),
117  points0_(pointIOField(points0IO(mesh)))
118 {
119  if (points0_.size() != mesh.nPoints())
120  {
122  << "Number of points in mesh " << mesh.nPoints()
123  << " differs from number of points " << points0_.size()
124  << " read from file "
125  << typeFilePath<pointIOField>
126  (
127  IOobject
128  (
129  "points",
130  time().constant(),
132  mesh,
135  false
136  )
137  )
138  << exit(FatalError);
139  }
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
144 
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
152 {}
153 
154 
156 {
157  // pointMesh already updates pointFields
158 
160 
161  // Map points0_. Bit special since we somehow have to come up with
162  // a sensible points0 position for introduced points.
163  // Find out scaling between points0 and current points
164 
165  // Get the new points either from the map or the mesh
166  const pointField& points =
167  (
168  mpm.hasMotionPoints()
169  ? mpm.preMotionPoints()
170  : mesh().points()
171  );
172 
173  // Note: boundBox does reduce
174  const vector span0 = boundBox(points0_).span();
175  const vector span = boundBox(points).span();
176 
177  vector scaleFactors(cmptDivide(span0, span));
178 
179  pointField newPoints0(mpm.pointMap().size());
180 
181  forAll(newPoints0, pointi)
182  {
183  label oldPointi = mpm.pointMap()[pointi];
184 
185  if (oldPointi >= 0)
186  {
187  label masterPointi = mpm.reversePointMap()[oldPointi];
188 
189  if (masterPointi == pointi)
190  {
191  newPoints0[pointi] = points0_[oldPointi];
192  }
193  else
194  {
195  // New point - assume motion is scaling
196  newPoints0[pointi] = points0_[oldPointi] + cmptMultiply
197  (
198  scaleFactors,
199  points[pointi] - points[masterPointi]
200  );
201  }
202  }
203  else
204  {
206  << "Cannot determine co-ordinates of introduced vertices."
207  << " New vertex " << pointi << " at co-ordinate "
208  << points[pointi] << exit(FatalError);
209  }
210  }
211 
212  twoDCorrectPoints(newPoints0);
213 
214  points0_.transfer(newPoints0);
215 
216  // points0 changed - set to write and check-in to database
217  points0_.rename("points0");
218  points0_.writeOpt() = IOobject::AUTO_WRITE;
219  points0_.instance() = time().timeName();
220  points0_.checkIn();
221 }
222 
223 
224 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:463
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
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
vectorIOField pointIOField
pointIOField is a vectorIOField.
Definition: pointIOField.H:42
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
bool hasMotionPoints() const
Has valid preMotionPoints?
Definition: mapPolyMesh.H:613
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
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:52
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1011
virtual void updateMesh(const mapPolyMesh &)=0
Update local data for topology changes.
Definition: motionSolver.C:216
dynamicFvMesh & mesh
const pointField & points
IOobject points0IO(const polyMesh &mesh) const
Return IO object for points0.
A class for handling words, derived from string.
Definition: word.H:59
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:390
Constant dispersed-phase particle diameter model.
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:789
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
virtual ~points0MotionSolver()
Destructor.
label nPoints() const
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
const pointField & preMotionPoints() const
Pre-motion point positions.
Definition: mapPolyMesh.H:607
Namespace for OpenFOAM.