points0MotionSolver.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-2019 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const polyMesh& mesh,
41  const dictionary& dict,
42  const word& type
43 )
44 :
45  motionSolver(mesh, dict, type),
46  points0_(pointIOField(polyMesh::points0IO(mesh)))
47 {
48  if (points0_.size() != mesh.nPoints())
49  {
51  << "Number of points in mesh " << mesh.nPoints()
52  << " differs from number of points " << points0_.size()
53  << " read from file "
54  << typeFilePath<pointIOField>
55  (
56  IOobject
57  (
58  "points",
59  mesh.time().constant(),
61  mesh,
64  false
65  )
66  )
67  << exit(FatalError);
68  }
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {}
82 
83 
85 {
86  // pointMesh already updates pointFields
87 
89 
90  // Map points0_. Bit special since we somehow have to come up with
91  // a sensible points0 position for introduced points.
92  // Find out scaling between points0 and current points
93 
94  // Get the new points either from the map or the mesh
95  const pointField& points =
96  (
97  mpm.hasMotionPoints()
98  ? mpm.preMotionPoints()
99  : mesh().points()
100  );
101 
102  // Note: boundBox does reduce
103  const vector span0 = boundBox(points0_).span();
104  const vector span = boundBox(points).span();
105 
106  vector scaleFactors(cmptDivide(span0, span));
107 
108  pointField newPoints0(mpm.pointMap().size());
109 
110  forAll(newPoints0, pointi)
111  {
112  label oldPointi = mpm.pointMap()[pointi];
113 
114  if (oldPointi >= 0)
115  {
116  label masterPointi = mpm.reversePointMap()[oldPointi];
117 
118  if (masterPointi == pointi)
119  {
120  newPoints0[pointi] = points0_[oldPointi];
121  }
122  else
123  {
124  // New point - assume motion is scaling
125  newPoints0[pointi] = points0_[oldPointi] + cmptMultiply
126  (
127  scaleFactors,
128  points[pointi] - points[masterPointi]
129  );
130  }
131  }
132  else
133  {
135  << "Cannot determine co-ordinates of introduced vertices."
136  << " New vertex " << pointi << " at co-ordinate "
137  << points[pointi] << exit(FatalError);
138  }
139  }
140 
141  twoDCorrectPoints(newPoints0);
142 
143  points0_.transfer(newPoints0);
144 
145  // points0 changed - set to write and check-in to database
146  points0_.rename("points0");
147  points0_.writeOpt() = IOobject::AUTO_WRITE;
148  points0_.instance() = mesh().time().timeName();
149  points0_.checkIn();
150 }
151 
152 
153 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:458
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
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:323
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:164
Virtual base class for mesh motion solver.
Definition: motionSolver.H:55
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
bool hasMotionPoints() const
Has valid preMotionPoints?
Definition: mapPolyMesh.H:608
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 > &)
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
virtual void updateMesh(const mapPolyMesh &)=0
Update local data for topology changes.
Definition: motionSolver.C:147
dynamicFvMesh & mesh
const pointField & points
points0MotionSolver(const polyMesh &, const dictionary &, const word &type)
Construct from mesh and dictionary.
A class for handling words, derived from string.
Definition: word.H:59
virtual void movePoints(const pointField &)
Update local data for geometry changes.
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:385
static IOobject points0IO(const polyMesh &mesh)
Return IO object for points0.
Definition: polyMesh.C:1220
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
const pointField & preMotionPoints() const
Pre-motion point positions.
Definition: mapPolyMesh.H:602
Namespace for OpenFOAM.