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-2024 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 "polyTopoChangeMap.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
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 word& name,
72  const polyMesh& mesh,
73  const dictionary& dict,
74  const word& type
75 )
76 :
77  motionSolver(name, mesh, dict, type),
78  cmptName_(coeffDict().lookup("component")),
79  cmpt_(cmpt(cmptName_)),
80  points0_
81  (
83  (
84  IOobject
85  (
86  "points",
87  mesh.time().constant(),
88  polyMesh::meshSubDir,
89  mesh,
90  IOobject::MUST_READ,
91  IOobject::NO_WRITE,
92  false
93  )
94  ).component(cmpt_)
95  ),
96  pointDisplacement_
97  (
98  IOobject
99  (
100  "pointDisplacement" + cmptName_,
101  mesh.time().name(),
102  mesh,
103  IOobject::MUST_READ,
104  IOobject::AUTO_WRITE
105  ),
106  pointMesh::New(mesh)
107  )
108 {
109  if (points0_.size() != mesh.nPoints())
110  {
112  << "Number of points in mesh " << mesh.nPoints()
113  << " differs from number of points " << points0_.size()
114  << " read from file "
116  (
117  "points",
118  mesh.time().constant(),
120  mesh,
123  false
124  ).filePath()
125  << exit(FatalError);
126  }
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
131 
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
137 
139 {
140  // No local data to update
141 }
142 
143 
145 (
146  const polyTopoChangeMap& map
147 )
148 {
149  // pointMesh already updates pointFields.
150 
151  // Map points0_. Bit special since we somehow have to come up with
152  // a sensible points0 position for introduced points.
153  // Find out scaling between points0 and current points
154 
155  // Get the new points either from the map or the mesh
156  const scalarField points(mesh().points().component(cmpt_));
157 
158  // Get extents of points0 and points and determine scale
159  const scalar scale =
160  (gMax(points0_)-gMin(points0_))
161  /(gMax(points)-gMin(points));
162 
163  scalarField newPoints0(map.pointMap().size());
164 
165  forAll(newPoints0, pointi)
166  {
167  label oldPointi = map.pointMap()[pointi];
168 
169  if (oldPointi >= 0)
170  {
171  label masterPointi = map.reversePointMap()[oldPointi];
172 
173  if (masterPointi == pointi)
174  {
175  newPoints0[pointi] = points0_[oldPointi];
176  }
177  else
178  {
179  // New point. Assume motion is scaling.
180  newPoints0[pointi] =
181  points0_[oldPointi]
182  + scale*(points[pointi]-points[masterPointi]);
183  }
184  }
185  else
186  {
188  << "Cannot work out coordinates of introduced vertices."
189  << " New vertex " << pointi << " at coordinate "
190  << points[pointi] << exit(FatalError);
191  }
192  }
193  points0_.transfer(newPoints0);
194 }
195 
196 
198 {
199  points0_ = mesh().points().component(cmpt_);
200  pointDisplacement_ == Zero;
201 }
202 
203 
205 (
206  const polyDistributionMap& map
207 )
208 {}
209 
210 
211 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
Virtual base class for displacement motion solver.
scalarField points0_
Reference point field for this component.
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
virtual void distribute(const polyDistributionMap &)
Update corresponding to the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
componentDisplacementMotionSolver(const word &name, const polyMesh &, const dictionary &, const word &type)
Construct from polyMesh and dictionary and type.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:137
const Time & time() const
Return time.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:53
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:270
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reversePointMap() const
Reverse point map.
const labelList & pointMap() const
Old point map.
label nPoints() const
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
fileName filePath() const
Return the path for the file for this Type.
Definition: IOobject.H:563
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const pointField & points
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
defineTypeNameAndDebug(combustionModel, 0)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
error FatalError
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
uint8_t direction
Definition: direction.H:45
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
volScalarField & p