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 | Copyright (C) 2012-2016 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 
69 Foam::componentDisplacementMotionSolver::componentDisplacementMotionSolver
70 (
71  const polyMesh& mesh,
72  const IOdictionary& 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  time().constant(),
88  mesh,
89  IOobject::MUST_READ,
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  <<
115  IOobject
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  // pointMesh already updates pointFields.
147 
149 
150  // Map points0_. Bit special since we somehow have to come up with
151  // a sensible points0 position for introduced points.
152  // Find out scaling between points0 and current points
153 
154  // Get the new points either from the map or the mesh
155  const scalarField points
156  (
157  mpm.hasMotionPoints()
158  ? mpm.preMotionPoints().component(cmpt_)
159  : mesh().points().component(cmpt_)
160  );
161 
162  // Get extents of points0 and points and determine scale
163  const scalar scale =
164  (gMax(points0_)-gMin(points0_))
165  /(gMax(points)-gMin(points));
166 
167  scalarField newPoints0(mpm.pointMap().size());
168 
169  forAll(newPoints0, pointi)
170  {
171  label oldPointi = mpm.pointMap()[pointi];
172 
173  if (oldPointi >= 0)
174  {
175  label masterPointi = mpm.reversePointMap()[oldPointi];
176 
177  if (masterPointi == pointi)
178  {
179  newPoints0[pointi] = points0_[oldPointi];
180  }
181  else
182  {
183  // New point. Assume motion is scaling.
184  newPoints0[pointi] =
185  points0_[oldPointi]
186  + scale*(points[pointi]-points[masterPointi]);
187  }
188  }
189  else
190  {
192  << "Cannot work out coordinates of introduced vertices."
193  << " New vertex " << pointi << " at coordinate "
194  << points[pointi] << exit(FatalError);
195  }
196  }
197  points0_.transfer(newPoints0);
198 }
199 
200 
201 // ************************************************************************* //
const Time & time() const
Return time.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#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.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
bool hasMotionPoints() const
Has valid preMotionPoints?
Definition: mapPolyMesh.H:613
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:390
Virtual base class for mesh motion solver.
Definition: motionSolver.H:53
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
static const pointMesh & New(const polyMesh &mesh)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual void updateMesh(const mapPolyMesh &)=0
Update local data for topology changes.
Definition: motionSolver.C:174
dynamicFvMesh & mesh
const pointField & points
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
Type gMax(const FieldField< Field, Type > &f)
defineTypeNameAndDebug(combustionModel, 0)
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:650
Constant dispersed-phase particle diameter model.
label nPoints() const
const pointField & preMotionPoints() const
Pre-motion point positions.
Definition: mapPolyMesh.H:607
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:91
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Namespace for OpenFOAM.
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:463