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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "mapPolyMesh.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(componentDisplacementMotionSolver, 0);
36 }
37 
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 Foam::direction Foam::componentDisplacementMotionSolver::cmpt
42 (
43  const word& cmptName
44 ) const
45 {
46  if (cmptName == "x")
47  {
48  return vector::X;
49  }
50  else if (cmptName == "y")
51  {
52  return vector::Y;
53  }
54  else if (cmptName == "z")
55  {
56  return vector::Z;
57  }
58  else
59  {
61  << "Given component name " << cmptName << " should be x, y or z"
62  << exit(FatalError);
63 
64  return 0;
65  }
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
71 Foam::componentDisplacementMotionSolver::componentDisplacementMotionSolver
72 (
73  const polyMesh& mesh,
74  const IOdictionary& dict,
75  const word& type
76 )
77 :
79  cmptName_(coeffDict().lookup("component")),
80  cmpt_(cmpt(cmptName_)),
81  points0_
82  (
84  (
85  IOobject
86  (
87  "points",
88  time().constant(),
90  mesh,
93  false
94  )
95  ).component(cmpt_)
96  ),
97  pointDisplacement_
98  (
99  IOobject
100  (
101  "pointDisplacement" + cmptName_,
102  mesh.time().timeName(),
103  mesh,
106  ),
108  )
109 {
110  if (points0_.size() != mesh.nPoints())
111  {
113  << "Number of points in mesh " << mesh.nPoints()
114  << " differs from number of points " << points0_.size()
115  << " read from file "
116  << typeFilePath<pointIOField>
117  (
118  IOobject
119  (
120  "points",
121  mesh.time().constant(),
123  mesh,
126  false
127  )
128  )
129  << exit(FatalError);
130  }
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
144  // No local data to update
145 }
146 
147 
149 {
150  // pointMesh already updates pointFields.
151 
153 
154  // Map points0_. Bit special since we somehow have to come up with
155  // a sensible points0 position for introduced points.
156  // Find out scaling between points0 and current points
157 
158  // Get the new points either from the map or the mesh
159  const scalarField points
160  (
161  mpm.hasMotionPoints()
162  ? mpm.preMotionPoints().component(cmpt_)
163  : mesh().points().component(cmpt_)
164  );
165 
166  // Get extents of points0 and points and determine scale
167  const scalar scale =
168  (gMax(points0_)-gMin(points0_))
169  /(gMax(points)-gMin(points));
170 
171  scalarField newPoints0(mpm.pointMap().size());
172 
173  forAll(newPoints0, pointi)
174  {
175  label oldPointi = mpm.pointMap()[pointi];
176 
177  if (oldPointi >= 0)
178  {
179  label masterPointi = mpm.reversePointMap()[oldPointi];
180 
181  if (masterPointi == pointi)
182  {
183  newPoints0[pointi] = points0_[oldPointi];
184  }
185  else
186  {
187  // New point. Assume motion is scaling.
188  newPoints0[pointi] =
189  points0_[oldPointi]
190  + scale*(points[pointi]-points[masterPointi]);
191  }
192  }
193  else
194  {
196  << "Cannot work out coordinates of introduced vertices."
197  << " New vertex " << pointi << " at coordinate "
198  << points[pointi] << exit(FatalError);
199  }
200  }
201  points0_.transfer(newPoints0);
202 }
203 
204 
205 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:44
Foam::Vector< scalar >::Z
Definition: Vector.H:80
Foam::Vector< scalar >::Y
Definition: Vector.H:80
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
mapPolyMesh.H
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:315
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::mapPolyMesh::preMotionPoints
const pointField & preMotionPoints() const
Pre-motion point positions.
Definition: mapPolyMesh.H:612
Foam::componentDisplacementMotionSolver::movePoints
virtual void movePoints(const pointField &)
Update local data for geometry changes.
Definition: componentDisplacementMotionSolver.C:142
componentDisplacementMotionSolver.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::componentDisplacementMotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
Definition: componentDisplacementMotionSolver.C:148
Foam::label
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:62
Foam::Field< vector >
Foam::cellModeller::lookup
const cellModel * lookup(const word &modelName)
Deprecated(2017-11) equivalent to cellModel::ptr static method.
Definition: cellModeller.H:46
Foam::Vector< scalar >::X
Definition: Vector.H:80
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::componentDisplacementMotionSolver::~componentDisplacementMotionSolver
virtual ~componentDisplacementMotionSolver()
Destructor.
Definition: componentDisplacementMotionSolver.C:136
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:562
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::motionSolver
Virtual base class for mesh motion solver.
Definition: motionSolver.H:58
Foam::mapPolyMesh::reversePointMap
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:468
Foam::mapPolyMesh::pointMap
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:395
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::motionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)=0
Update local data for topology changes.
Definition: motionSolver.C:219
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::primitiveMesh::nPoints
label nPoints() const
Number of mesh points.
Definition: primitiveMeshI.H:37
Foam::direction
uint8_t direction
Definition: direction.H:47
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::mapPolyMesh::hasMotionPoints
bool hasMotionPoints() const
Has valid preMotionPoints?
Definition: mapPolyMesh.H:618
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:88
constant
constant condensation/saturation model.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::IOobject::MUST_READ
Definition: IOobject.H:120