displacementComponentLaplacianFvMotionSolver.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-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 "motionInterpolation.H"
31#include "motionDiffusivity.H"
32#include "fvmLaplacian.H"
34#include "mapPolyMesh.H"
35#include "fvOptions.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
42
44 (
48 );
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
54Foam::displacementComponentLaplacianFvMotionSolver::
55displacementComponentLaplacianFvMotionSolver
56(
57 const polyMesh& mesh,
58 const IOdictionary& dict
59)
60:
63 cellDisplacement_
64 (
66 (
67 "cellDisplacement" + cmptName_,
68 mesh.time().timeName(),
69 mesh,
70 IOobject::READ_IF_PRESENT,
71 IOobject::AUTO_WRITE
72 ),
73 fvMesh_,
74 dimensionedScalar(pointDisplacement_.dimensions(), Zero),
75 cellMotionBoundaryTypes<scalar>(pointDisplacement_.boundaryField())
76 ),
77 pointLocation_(nullptr),
78 interpolationPtr_
79 (
80 coeffDict().found("interpolation")
81 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
82 : motionInterpolation::New(fvMesh_)
83 ),
84 diffusivityPtr_
85 (
86 motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
87 ),
88 frozenPointsZone_
89 (
90 coeffDict().found("frozenPointsZone")
91 ? fvMesh_.pointZones().findZoneID
92 (
93 coeffDict().get<word>("frozenPointsZone")
94 )
95 : -1
96 )
97{
98 if (coeffDict().getOrDefault("applyPointLocation", true))
99 {
100 pointLocation_.reset
101 (
103 (
105 (
106 "pointLocation",
108 fvMesh_,
111 ),
113 )
114 );
115
116 //if (debug)
117 {
118 Info<< "displacementComponentLaplacianFvMotionSolver :"
119 << " Read pointVectorField "
120 << pointLocation_().name()
121 << " to be used for boundary conditions on points."
122 << nl
123 << "Boundary conditions:"
124 << pointLocation_().boundaryField().types() << endl;
125 }
126 }
127}
128
129
130// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
131
134{}
135
136
137// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138
141{
142 interpolationPtr_->interpolate
143 (
144 cellDisplacement_,
145 pointDisplacement_
146 );
147
148 if (pointLocation_)
149 {
150 if (debug)
151 {
152 Info<< "displacementComponentLaplacianFvMotionSolver : applying "
153 << " boundary conditions on " << pointLocation_().name()
154 << " to new point location."
155 << endl;
156 }
157
158 // Apply pointLocation_ b.c. to mesh points.
159
160 pointLocation_().primitiveFieldRef() = fvMesh_.points();
161
162 pointLocation_().primitiveFieldRef().replace
163 (
164 cmpt_,
165 points0_ + pointDisplacement_.primitiveField()
166 );
167
168 pointLocation_().correctBoundaryConditions();
169
170 // Implement frozen points
171 if (frozenPointsZone_ != -1)
172 {
173 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
174
175 forAll(pz, i)
176 {
177 label pointi = pz[i];
178
179 pointLocation_()[pointi][cmpt_] = points0_[pointi];
180 }
181 }
182
183 twoDCorrectPoints(pointLocation_().primitiveFieldRef());
184
185 return tmp<pointField>(pointLocation_().primitiveField());
186 }
187 else
188 {
189 tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
190 pointField& curPoints = tcurPoints.ref();
191
192 curPoints.replace
193 (
194 cmpt_,
195 points0_ + pointDisplacement_.primitiveField()
196 );
197
198 // Implement frozen points
199 if (frozenPointsZone_ != -1)
200 {
201 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
202
203 forAll(pz, i)
204 {
205 label pointi = pz[i];
206
207 curPoints[pointi][cmpt_] = points0_[pointi];
208 }
209 }
210
211 twoDCorrectPoints(curPoints);
212
213 return tcurPoints;
214 }
215}
216
217
219{
220 // The points have moved so before interpolation update
221 // the motionSolver accordingly
222 movePoints(fvMesh_.points());
223
224 diffusivityPtr_->correct();
225 pointDisplacement_.boundaryFieldRef().updateCoeffs();
226
228
229 // We explicitly do NOT want to interpolate the motion inbetween
230 // different regions so bypass all the matrix manipulation.
232 (
234 (
235 dimensionedScalar("viscosity", dimViscosity, 1.0)
236 *diffusivityPtr_->operator()(),
237 cellDisplacement_,
238 "laplacian(diffusivity,cellDisplacement)"
239 )
240 ==
241 fvOptions(cellDisplacement_)
242 );
243
246 fvOptions.correct(cellDisplacement_);
247}
248
249
251(
252 const mapPolyMesh& mpm
253)
254{
256
257 // Update diffusivity. Note two stage to make sure old one is de-registered
258 // before creating/registering new one.
259 diffusivityPtr_.reset(nullptr);
260 diffusivityPtr_ = motionDiffusivity::New
261 (
262 fvMesh_,
263 coeffDict().lookup("diffusivity")
264 );
265}
266
267
268// ************************************************************************* //
bool found
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
fv::options & fvOptions
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:557
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Virtual base class for displacement motion solver.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Mesh motion solver for an fvMesh. Based on solving the cell-centre Laplacian for the given component ...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:62
const dictionary & solverDict() const
Return the solver dictionary taking into account finalIteration.
Definition: fvMatrix.C:1558
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Base class for fvMesh based motionSolvers.
const fvMesh & fvMesh_
The fvMesh to be moved.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Finite-volume options.
Definition: fvOptions.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Abstract base class for cell-centre mesh motion diffusivity.
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers,...
void updateMesh()
Update for new mesh topology.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:61
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: motionSolver.H:150
A subset of mesh points.
Definition: pointZone.H:68
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
Lookup type of boundary radiation properties.
Definition: lookup.H:66
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
Calculate the matrix for the laplacian of the field.
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
word timeName
Definition: getTimeIndex.H:3
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Namespace for OpenFOAM.
const dimensionSet dimViscosity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
conserve primitiveFieldRef()+