componentDisplacementMotionSolver.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) 2020 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 
30 #include "mapPolyMesh.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(componentDisplacementMotionSolver, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 Foam::direction Foam::componentDisplacementMotionSolver::cmpt
43 (
44  const word& cmptName
45 ) const
46 {
47  if (cmptName == "x")
48  {
49  return vector::X;
50  }
51  else if (cmptName == "y")
52  {
53  return vector::Y;
54  }
55  else if (cmptName == "z")
56  {
57  return vector::Z;
58  }
59  else
60  {
62  << "Given component name " << cmptName << " should be x, y or z"
63  << exit(FatalError);
64 
65  return 0;
66  }
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71 
72 Foam::componentDisplacementMotionSolver::componentDisplacementMotionSolver
73 (
74  const polyMesh& mesh,
75  const IOdictionary& dict,
76  const word& type
77 )
78 :
79  motionSolver(mesh, dict, type),
80  cmptName_(coeffDict().get<word>("component")),
81  cmpt_(cmpt(cmptName_)),
82  points0_
83  (
85  (
86  IOobject
87  (
88  "points",
89  time().constant(),
90  polyMesh::meshSubDir,
91  mesh,
92  IOobject::MUST_READ,
93  IOobject::NO_WRITE,
94  IOobject::NO_REGISTER
95  )
96  ).component(cmpt_)
97  ),
98  pointDisplacement_
99  (
100  IOobject
101  (
102  "pointDisplacement" + cmptName_,
103  mesh.time().timeName(),
104  mesh,
105  IOobject::MUST_READ,
106  IOobject::AUTO_WRITE
107  ),
108  pointMesh::New(mesh)
109  )
110 {
111  if (points0_.size() != mesh.nPoints())
112  {
114  << "Number of points in mesh " << mesh.nPoints()
115  << " differs from number of points " << points0_.size()
116  << " read from file "
117  << typeFilePath<pointIOField>
118  (
119  IOobject
120  (
121  "points",
122  mesh.time().constant(),
124  mesh,
128  )
129  )
131  }
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144 {
145  // No local data to update
146 }
147 
148 
150 {
151  // pointMesh already updates pointFields.
152 
154 
155  // Map points0_. Bit special since we somehow have to come up with
156  // a sensible points0 position for introduced points.
157  // Find out scaling between points0 and current points
158 
159  // Get the new points either from the map or the mesh
160  const scalarField points
161  (
162  mpm.hasMotionPoints()
163  ? mpm.preMotionPoints().component(cmpt_)
164  : mesh().points().component(cmpt_)
165  );
166 
167  // Get extents of points0 and points and determine scale
168  const scalar scale =
169  (gMax(points0_)-gMin(points0_))
170  /(gMax(points)-gMin(points));
171 
172  scalarField newPoints0(mpm.pointMap().size());
173 
174  forAll(newPoints0, pointi)
175  {
176  label oldPointi = mpm.pointMap()[pointi];
177 
178  if (oldPointi >= 0)
179  {
180  label masterPointi = mpm.reversePointMap()[oldPointi];
181 
182  if (masterPointi == pointi)
183  {
184  newPoints0[pointi] = points0_[oldPointi];
185  }
186  else
187  {
188  // New point. Assume motion is scaling.
189  newPoints0[pointi] =
190  points0_[oldPointi]
191  + scale*(points[pointi]-points[masterPointi]);
192  }
193  }
194  else
195  {
197  << "Cannot work out coordinates of introduced vertices."
198  << " New vertex " << pointi << " at coordinate "
199  << points[pointi] << exit(FatalError);
200  }
201  }
202  points0_.transfer(newPoints0);
203 }
204 
205 
206 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
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
uint8_t direction
Definition: direction.H:48
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
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:578
Type gMin(const FieldField< Field, Type > &f)
vectorIOField pointIOField
pointIOField is a vectorIOField.
Definition: pointIOField.H:38
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:402
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Ignore writing from objectRegistry::writeObject()
bool hasMotionPoints() const
Has valid preMotionPoints?
Definition: mapPolyMesh.H:764
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1073
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
word timeName
Definition: getTimeIndex.H:3
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
virtual void updateMesh(const mapPolyMesh &)=0
Update local data for topology changes.
Definition: motionSolver.C:206
dynamicFvMesh & mesh
const pointField & points
const Time & time() const noexcept
Return time registry.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:608
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:89
Type gMax(const FieldField< Field, Type > &f)
defineTypeNameAndDebug(combustionModel, 0)
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:481
PtrList< volScalarField > & Y
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:165
volScalarField & p
scalarField points0_
Reference point field for this component.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Do not request registration (bool: false)
const pointField & preMotionPoints() const
Pre-motion point positions.
Definition: mapPolyMesh.H:756
Namespace for OpenFOAM.