Loading...
Searching...
No Matches
velocityDisplacementMotionSolver.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) 2015 OpenFOAM Foundation
9 Copyright (C) 2016 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
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
39
41 (
45 );
46}
47
48
49// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50
52Foam::velocityDisplacementMotionSolver::pointDisplacementBoundaryTypes() const
53{
54 const pointVectorField::Boundary& pmUbf(pointMotionU().boundaryField());
55
56 wordList cmUbf = pmUbf.types();
57
58 forAll(pmUbf, patchI)
59 {
60 if (isA<fixedValuePointPatchField<vector>>(pmUbf[patchI]))
61 {
63 }
64 }
66 return cmUbf;
67}
68
69
70// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71
72Foam::velocityDisplacementMotionSolver::velocityDisplacementMotionSolver
73(
74 const polyMesh& mesh,
75 const IOdictionary& dict
76)
77:
79 displacementMotionSolverPtr_()
80{
82
83 pointVectorField pointDisplacement
84 (
86 (
87 "pointVelocityDisplacement",
88 mesh.time().timeName(),
89 mesh
90 ),
93 pointDisplacementBoundaryTypes()
94 );
95
96 pointDisplacement.primitiveFieldRef() = mesh.points() - points0;
97
98 displacementMotionSolverPtr_.reset
99 (
100 dynamic_cast<displacementMotionSolver*>
101 (
103 (
104 coeffDict().get<word>("solver"),
105 mesh,
107 (
109 (
110 dict.name() + "Coeffs",
111 mesh.time().constant(),
112 mesh
113 ),
114 coeffDict()
115 ),
116 pointDisplacement,
117 points0
118 ).ptr()
120 );
121}
122
123
124// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
125
130// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131
134{
135 return displacementMotionSolverPtr_->curPoints();
136}
137
138
140{
141 movePoints(mesh().points());
142
143 const scalar deltaT(mesh().time().deltaTValue());
144
145 // Current and old point displacements
146 pointVectorField& displacement
147 (
148 displacementMotionSolverPtr_->pointDisplacement()
149 );
150 const vectorField displacementOld
151 (
152 mesh().points() - displacementMotionSolverPtr_->points0()
153 );
154
155 // Update the velocity boundary conditions
156 pointMotionU().correctBoundaryConditions();
157
158 pointVectorField::Boundary& dispBf = displacement.boundaryFieldRef();
159
160 // Update the displacement boundary conditions
161 forAll(pointMotionU().boundaryField(), patchI)
162 {
163 const pointPatchVectorField& patchField
164 (
165 pointMotionU().boundaryField()[patchI]
166 );
167
168 dispBf[patchI] ==
169 patchField.patchInternalField()*deltaT
170 + patchField.patchInternalField(displacementOld);
171 }
172
173 // Run the sub-solver
174 displacementMotionSolverPtr_->solve();
176 // Update the velocity
177 pointMotionU().primitiveFieldRef() =
178 (displacement.primitiveField() - displacementOld)/deltaT;
179}
180
181
185
186 displacementMotionSolverPtr_->movePoints(p);
187}
188
189
191(
192 const mapPolyMesh& mpm
193)
194{
196
197 displacementMotionSolverPtr_->updateMesh(mpm);
198}
199
200
201// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static const char *const typeName
Typename for Field.
Definition Field.H:93
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< vector, pointPatchField, pointMesh > Boundary
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
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
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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.
Virtual base class for displacement motion solver.
static autoPtr< displacementMotionSolver > New(const word &solverTypeName, const polyMesh &, const IOdictionary &, const pointVectorField &pointDisplacement, const pointIOField &points0)
Select constructed from polyMesh, dictionary and components.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
void patchInternalField(const UList< Type1 > &internalData, const labelUList &addressing, UList< Type1 > &pfld) const
Extract field using specified addressing.
static IOobject points0IO(const polyMesh &mesh)
Return IO object for points0.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A class for managing temporary objects.
Definition tmp.H:75
Mesh motion solver for a polyMesh. Wraps a displacement motion solver in a velocity motion solver.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void movePoints(const pointField &)
Update geometry.
virtual void updateMesh(const mapPolyMesh &)
Update topology.
Virtual base class for velocity motion solver.
velocityMotionSolver(const velocityMotionSolver &)=delete
No copy construct.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
pointVectorField & pointMotionU()
Return reference to the point motion velocity field.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
const pointField & points
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vectorIOField pointIOField
pointIOField is a vectorIOField.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
Field< vector > vectorField
Specialisation of Field<T> for vector.
pointPatchField< vector > pointPatchVectorField
vectorField pointField
pointField is a vectorField.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))