componentDisplacementMotionSolver.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) 2012-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 
27 #include "mapPolyMesh.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(componentDisplacementMotionSolver, 0);
34 }
35 
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 Foam::direction Foam::componentDisplacementMotionSolver::cmpt
40 (
41  const word& cmptName
42 ) const
43 {
44  if (cmptName == "x")
45  {
46  return vector::X;
47  }
48  else if (cmptName == "y")
49  {
50  return vector::Y;
51  }
52  else if (cmptName == "z")
53  {
54  return vector::Z;
55  }
56  else
57  {
59  << "Given component name " << cmptName << " should be x, y or z"
60  << exit(FatalError);
61 
62  return 0;
63  }
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
70 (
71  const polyMesh& mesh,
72  const dictionary& dict,
73  const word& type
74 )
75 :
76  motionSolver(mesh, dict, type),
77  cmptName_(coeffDict().lookup("component")),
78  cmpt_(cmpt(cmptName_)),
79  points0_
80  (
82  (
83  IOobject
84  (
85  "points",
86  mesh.time().constant(),
88  mesh,
91  false
92  )
93  ).component(cmpt_)
94  ),
95  pointDisplacement_
96  (
97  IOobject
98  (
99  "pointDisplacement" + cmptName_,
100  mesh.time().timeName(),
101  mesh,
104  ),
105  pointMesh::New(mesh)
106  )
107 {
108  if (points0_.size() != mesh.nPoints())
109  {
111  << "Number of points in mesh " << mesh.nPoints()
112  << " differs from number of points " << points0_.size()
113  << " read from file "
114  << typeFilePath<pointIOField>
115  (
116  IOobject
117  (
118  "points",
119  mesh.time().constant(),
121  mesh,
124  false
125  )
126  )
127  << exit(FatalError);
128  }
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
133 
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
141 {
142  // No local data to update
143 }
144 
145 
147 {
148  // pointMesh already updates pointFields.
149 
151 
152  // Map points0_. Bit special since we somehow have to come up with
153  // a sensible points0 position for introduced points.
154  // Find out scaling between points0 and current points
155 
156  // Get the new points either from the map or the mesh
157  const scalarField points
158  (
159  mpm.hasMotionPoints()
160  ? mpm.preMotionPoints().component(cmpt_)
161  : mesh().points().component(cmpt_)
162  );
163 
164  // Get extents of points0 and points and determine scale
165  const scalar scale =
166  (gMax(points0_)-gMin(points0_))
167  /(gMax(points)-gMin(points));
168 
169  scalarField newPoints0(mpm.pointMap().size());
170 
171  forAll(newPoints0, pointi)
172  {
173  label oldPointi = mpm.pointMap()[pointi];
174 
175  if (oldPointi >= 0)
176  {
177  label masterPointi = mpm.reversePointMap()[oldPointi];
178 
179  if (masterPointi == pointi)
180  {
181  newPoints0[pointi] = points0_[oldPointi];
182  }
183  else
184  {
185  // New point. Assume motion is scaling.
186  newPoints0[pointi] =
187  points0_[oldPointi]
188  + scale*(points[pointi]-points[masterPointi]);
189  }
190  }
191  else
192  {
194  << "Cannot work out coordinates of introduced vertices."
195  << " New vertex " << pointi << " at coordinate "
196  << points[pointi] << exit(FatalError);
197  }
198  }
199  points0_.transfer(newPoints0);
200 }
201 
202 
203 // ************************************************************************* //
#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:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Type gMin(const FieldField< Field, Type > &f)
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
uint8_t direction
Definition: direction.H:45
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:55
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
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
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual void updateMesh(const mapPolyMesh &)=0
Update local data for topology changes.
Definition: motionSolver.C:152
dynamicFvMesh & mesh
const pointField & points
componentDisplacementMotionSolver(const polyMesh &, const dictionary &, const word &type)
Construct from polyMesh and dictionary and type.
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:124
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:460
Type gMax(const FieldField< Field, Type > &f)
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:385
label nPoints() const
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
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:602
Namespace for OpenFOAM.