Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
34namespace Foam
35{
37}
38
39
40// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41
42Foam::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
72Foam::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 {
113 const fileName fName
114 (
115 IOobject
116 (
117 "points",
118 mesh.time().constant(),
119 polyMesh::meshSubDir,
120 mesh,
121 IOobject::MUST_READ,
122 IOobject::NO_WRITE,
123 IOobject::NO_REGISTER
124 ).typeFilePath<pointIOField>()
125 );
126
127 FatalErrorInFunction
128 << "Number of points in mesh " << mesh.nPoints()
129 << " differs from number of points " << points0_.size()
130 << " read from file " << fName << nl
131 << exit(FatalError);
132 }
133}
134
135
136// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
141
142// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145{
146 // No local data to update
147}
148
149
151{
152 // pointMesh already updates pointFields.
153
155
156 // Map points0_. Bit special since we somehow have to come up with
157 // a sensible points0 position for introduced points.
158 // Find out scaling between points0 and current points
159
160 // Get the new points either from the map or the mesh
161 const scalarField points
162 (
163 mpm.hasMotionPoints()
164 ? mpm.preMotionPoints().component(cmpt_)
165 : mesh().points().component(cmpt_)
166 );
167
168 // Get extents of points0 and points and determine scale
169 const scalar scale = gMinMax(points0_).span() / gMinMax(points).span();
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// ************************************************************************* //
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition Field.C:607
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Virtual base class for displacement motion solver.
scalarField points0_
Reference point field for this component.
pointScalarField pointDisplacement_
Point motion field.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const pointField & preMotionPoints() const noexcept
Pre-motion point positions.
bool hasMotionPoints() const noexcept
Has valid preMotionPoints?
const labelList & reversePointMap() const noexcept
Reverse point map.
const labelList & pointMap() const noexcept
Old point map.
virtual void updateMesh(const mapPolyMesh &)=0
Update local data for topology changes.
const polyMesh & mesh() const
Return reference to mesh.
static autoPtr< motionSolver > New(const polyMesh &)
Select constructed from polyMesh.
motionSolver(const polyMesh &mesh)
Construct from polyMesh.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
constant condensation/saturation model.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
word timeName
Definition getTimeIndex.H:3
Different types of constants.
Namespace for OpenFOAM.
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.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
vectorIOField pointIOField
pointIOField is a vectorIOField.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
uint8_t direction
Definition direction.H:49
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299