Loading...
Searching...
No Matches
displacementLaplacianFvMotionSolver.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-2017 OpenFOAM Foundation
9 Copyright (C) 2015-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 "OFstream.H"
35#include "meshTools.H"
36#include "mapPolyMesh.H"
37#include "fvOptions.H"
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
44
46 (
50 );
51
53 (
56 displacement
57 );
58}
59
60
61// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62
63Foam::displacementLaplacianFvMotionSolver::displacementLaplacianFvMotionSolver
64(
65 const polyMesh& mesh,
66 const IOdictionary& dict
67)
68:
71 cellDisplacement_
72 (
74 (
75 "cellDisplacement",
76 mesh.time().timeName(),
77 mesh,
78 IOobject::READ_IF_PRESENT,
79 IOobject::AUTO_WRITE
80 ),
81 fvMesh_,
82 dimensionedVector(pointDisplacement_.dimensions(), Zero),
83 cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
84 ),
85 pointLocation_(nullptr),
86 interpolationPtr_
87 (
88 coeffDict().found("interpolation")
89 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
90 : motionInterpolation::New(fvMesh_)
91 ),
92 diffusivityPtr_
93 (
94 motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
95 ),
96 frozenPointsZone_
97 (
98 coeffDict().found("frozenPointsZone")
99 ? fvMesh_.pointZones().findZoneID
100 (
101 coeffDict().get<word>("frozenPointsZone")
102 )
103 : -1
104 )
105{
107 (
108 "pointLocation",
110 fvMesh_,
113 );
114
115 if (debug)
116 {
117 Info<< "displacementLaplacianFvMotionSolver:" << nl
118 << " diffusivity : " << diffusivityPtr_().type() << nl
119 << " frozenPoints zone : " << frozenPointsZone_ << endl;
120 }
121
122
123 if (io.typeHeaderOk<pointVectorField>(true))
124 {
125 pointLocation_.reset
126 (
128 (
129 io,
131 )
132 );
133
134 if (debug)
135 {
136 Info<< "displacementLaplacianFvMotionSolver :"
137 << " Read pointVectorField "
138 << io.name()
139 << " to be used for boundary conditions on points."
140 << nl
141 << "Boundary conditions:"
142 << pointLocation_().boundaryField().types() << endl;
143 }
144 }
145}
146
147
148Foam::displacementLaplacianFvMotionSolver::
149displacementLaplacianFvMotionSolver
150(
151 const polyMesh& mesh,
152 const IOdictionary& dict,
153 const pointVectorField& pointDisplacement,
154 const pointIOField& points0
155)
156:
157 displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName),
159 cellDisplacement_
160 (
162 (
163 "cellDisplacement",
164 mesh.time().timeName(),
165 mesh,
166 IOobject::READ_IF_PRESENT,
167 IOobject::AUTO_WRITE
168 ),
169 fvMesh_,
170 dimensionedVector(pointDisplacement_.dimensions(), Zero),
171 cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
172 ),
173 pointLocation_(nullptr),
174 interpolationPtr_
175 (
176 coeffDict().found("interpolation")
177 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
178 : motionInterpolation::New(fvMesh_)
179 ),
180 diffusivityPtr_
181 (
182 motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
183 ),
184 frozenPointsZone_
185 (
186 coeffDict().found("frozenPointsZone")
187 ? fvMesh_.pointZones().findZoneID
188 (
189 coeffDict().get<word>("frozenPointsZone")
190 )
191 : -1
192 )
193{
194 IOobject io
195 (
196 "pointLocation",
197 fvMesh_.time().timeName(),
198 fvMesh_,
201 );
202
203 if (debug)
204 {
205 Info<< "displacementLaplacianFvMotionSolver:" << nl
206 << " diffusivity : " << diffusivityPtr_().type() << nl
207 << " frozenPoints zone : " << frozenPointsZone_ << endl;
208 }
209
210
211 if (io.typeHeaderOk<pointVectorField>(true))
212 {
213 pointLocation_.reset
214 (
216 (
217 io,
219 )
220 );
221
222 if (debug)
223 {
224 Info<< "displacementLaplacianFvMotionSolver :"
225 << " Read pointVectorField "
226 << io.name()
227 << " to be used for boundary conditions on points."
228 << nl
229 << "Boundary conditions:"
230 << pointLocation_().boundaryField().types() << endl;
232 }
233}
234
235
236// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
237
240{}
241
242
243// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
244
247{
248 if (!diffusivityPtr_)
249 {
250 diffusivityPtr_ = motionDiffusivity::New
251 (
252 fvMesh_,
253 coeffDict().lookup("diffusivity")
254 );
256
257 return *diffusivityPtr_;
258}
259
260
263{
264 interpolationPtr_->interpolate
265 (
266 cellDisplacement_,
267 pointDisplacement_
268 );
269
270 // Evaluate the bcs so they are consistent with the internal field
271 // Might fight the multi-patch behaviour inside volPointInterpolate
272 if
273 (
274 pointDisplacement_.boundaryField().size()
275 != cellDisplacement_.boundaryField().size()
276 )
277 {
278 pointDisplacement_.correctBoundaryConditions();
279 }
280
281 if (pointLocation_)
282 {
283 if (debug)
284 {
285 Info<< "displacementLaplacianFvMotionSolver : applying "
286 << " boundary conditions on " << pointLocation_().name()
287 << " to new point location."
288 << endl;
289 }
290
291 pointLocation_().primitiveFieldRef() =
292 points0()
293 + pointDisplacement_.primitiveField();
294
295 pointLocation_().correctBoundaryConditions();
296
297 // Implement frozen points
298 if (frozenPointsZone_ != -1)
299 {
300 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
301
302 forAll(pz, i)
303 {
304 pointLocation_()[pz[i]] = points0()[pz[i]];
305 }
306 }
307
308 twoDCorrectPoints(pointLocation_().primitiveFieldRef());
309
310 return tmp<pointField>(pointLocation_().primitiveField());
311 }
312 else
313 {
314 tmp<pointField> tcurPoints
315 (
316 points0() + pointDisplacement_.primitiveField()
317 );
318 pointField& curPoints = tcurPoints.ref();
319
320 // Implement frozen points
321 if (frozenPointsZone_ != -1)
322 {
323 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
324
325 forAll(pz, i)
326 {
327 curPoints[pz[i]] = points0()[pz[i]];
328 }
329 }
330
332
333 return tcurPoints;
334 }
335}
336
337
339{
340 // The points have moved so before interpolation update
341 // the motionSolver accordingly
342 movePoints(fvMesh_.points());
343
344 diffusivity().correct();
345 pointDisplacement_.boundaryFieldRef().updateCoeffs();
346
347 // Make sure the cellMotion bcs are consistent with the pointDisplacement.
348 // This makes sure the grad below uses more up-to-date values.
349 cellDisplacement_.boundaryFieldRef().updateCoeffs();
350
351 fv::options& fvOptions(fv::options::New(fvMesh_));
352
353 // We explicitly do NOT want to interpolate the motion inbetween
354 // different regions so bypass all the matrix manipulation.
356 (
358 (
359 dimensionedScalar("viscosity", dimViscosity, 1.0)
360 *diffusivity().operator()(),
361 cellDisplacement_,
362 "laplacian(diffusivity,cellDisplacement)"
363 )
364 ==
365 fvOptions(cellDisplacement_)
366 );
368 fvOptions.constrain(TEqn);
369 TEqn.solveSegregatedOrCoupled();
370 fvOptions.correct(cellDisplacement_);
371}
372
373
375(
376 const mapPolyMesh& mpm
377)
378{
380
381 // Update diffusivity. Note two stage to make sure old one is de-registered
382 // before creating/registering new one.
383 diffusivityPtr_.clear();
384}
385
386
387// ************************************************************************* //
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
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.
Mesh motion solver for an fvMesh. Based on solving the cell-centre Laplacian for the motion displacem...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
motionDiffusivity & diffusivity()
Return reference to the diffusivity field.
virtual void updateMesh(const mapPolyMesh &)
Update topology.
Virtual base class for displacement motion solver.
pointVectorField & pointDisplacement() noexcept
Return reference to the point motion displacement field.
displacementMotionSolver(const displacementMotionSolver &)=delete
No copy construct.
pointVectorField pointDisplacement_
Point motion field.
static autoPtr< displacementMotionSolver > New(const word &solverTypeName, const polyMesh &, const IOdictionary &, const pointVectorField &pointDisplacement, const pointIOField &points0)
Select constructed from polyMesh, dictionary and components.
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.
virtual void twoDCorrectPoints(pointField &) const
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
A subset of mesh points.
Definition pointZone.H:64
pointField & points0() noexcept
Return reference to the reference ('0') pointField.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
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
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
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
const auto & io
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
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Namespace for handling debugging switches.
Definition debug.C:45
int debug
Static debugging option.
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.
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vectorIOField pointIOField
pointIOField is a vectorIOField.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
fvMatrix< vector > fvVectorMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
dimensioned< vector > dimensionedVector
Dimensioned vector 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
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))
conserve primitiveFieldRef()+