Loading...
Searching...
No Matches
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"
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;
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 // Evaluate the bcs so they are consistent with the internal field
149 // Might fight the multi-patch behaviour inside volPointInterpolate
150 if
151 (
152 pointDisplacement_.boundaryField().size()
153 != cellDisplacement_.boundaryField().size()
154 )
155 {
156 pointDisplacement_.correctBoundaryConditions();
157 }
158
159
160 if (pointLocation_)
161 {
162 if (debug)
163 {
164 Info<< "displacementComponentLaplacianFvMotionSolver : applying "
165 << " boundary conditions on " << pointLocation_().name()
166 << " to new point location."
167 << endl;
168 }
169
170 // Apply pointLocation_ b.c. to mesh points.
171
172 pointLocation_().primitiveFieldRef() = fvMesh_.points();
173
174 pointLocation_().primitiveFieldRef().replace
175 (
176 cmpt_,
177 points0_ + pointDisplacement_.primitiveField()
178 );
179
180 pointLocation_().correctBoundaryConditions();
181
182 // Implement frozen points
183 if (frozenPointsZone_ != -1)
184 {
185 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
186
187 forAll(pz, i)
188 {
189 label pointi = pz[i];
190
191 pointLocation_()[pointi][cmpt_] = points0_[pointi];
192 }
193 }
194
195 twoDCorrectPoints(pointLocation_().primitiveFieldRef());
196
197 return tmp<pointField>(pointLocation_().primitiveField());
198 }
199 else
200 {
201 auto tcurPoints = tmp<pointField>::New(fvMesh_.points());
202 auto& curPoints = tcurPoints.ref();
203
204 curPoints.replace
205 (
206 cmpt_,
207 points0_ + pointDisplacement_.primitiveField()
208 );
209
210 // Implement frozen points
211 if (frozenPointsZone_ != -1)
212 {
213 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
214
215 forAll(pz, i)
216 {
217 label pointi = pz[i];
218
219 curPoints[pointi][cmpt_] = points0_[pointi];
220 }
221 }
222
224
225 return tcurPoints;
226 }
227}
228
229
231{
232 // The points have moved so before interpolation update
233 // the motionSolver accordingly
234 movePoints(fvMesh_.points());
235
236 diffusivityPtr_->correct();
237 pointDisplacement_.boundaryFieldRef().updateCoeffs();
238
239 fv::options& fvOptions(fv::options::New(fvMesh_));
240
241 // We explicitly do NOT want to interpolate the motion inbetween
242 // different regions so bypass all the matrix manipulation.
244 (
246 (
247 dimensionedScalar("viscosity", dimViscosity, 1.0)
248 *diffusivityPtr_->operator()(),
249 cellDisplacement_,
250 "laplacian(diffusivity,cellDisplacement)"
251 )
252 ==
253 fvOptions(cellDisplacement_)
254 );
256 fvOptions.constrain(TEqn);
257 TEqn.solveSegregatedOrCoupled();
258 fvOptions.correct(cellDisplacement_);
259}
260
261
263(
264 const mapPolyMesh& mpm
265)
266{
268
269 // Update diffusivity. Note two stage to make sure old one is de-registered
270 // before creating/registering new one.
271 diffusivityPtr_.reset(nullptr);
272 diffusivityPtr_ = motionDiffusivity::New
273 (
274 fvMesh_,
275 coeffDict().lookup("diffusivity")
276 );
277}
278
279
280// ************************************************************************* //
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
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ 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
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
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.
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.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
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.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Base class for fvMesh based motionSolvers.
const fvMesh & fvMesh_
The fvMesh to be moved.
wordList cellMotionBoundaryTypes(const typename GeometricField< Type, pointPatchField, pointMesh >::Boundary &pmUbf) const
Create the corresponding patch types for cellMotion from those.
fvMotionSolver(const polyMesh &)
Construct from polyMesh.
const fvMesh & mesh() const
Return reference to the fvMesh to be moved.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Abstract base class for cell-centre mesh motion diffusivity.
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers,...
Virtual base class for mesh motion solver.
static autoPtr< motionSolver > New(const polyMesh &)
Select constructed from polyMesh.
virtual void twoDCorrectPoints(pointField &) const
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
A subset of mesh points.
Definition pointZone.H:64
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Lookup type of boundary radiation properties.
Definition lookup.H:60
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
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
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
Namespace for handling debugging switches.
Definition debug.C:45
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Namespace for OpenFOAM.
const dimensionSet dimViscosity
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.
fvMatrix< scalar > fvScalarMatrix
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
conserve primitiveFieldRef()+