velocityLaplacianFvMotionSolver.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 "fvOptions.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
41
43 (
47 );
48}
49
50
51// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52
54(
55 const polyMesh& mesh,
56 const IOdictionary& dict
57)
58:
59 velocityMotionSolver(mesh, dict, typeName),
61 cellMotionU_
62 (
64 (
65 "cellMotionU",
66 mesh.time().timeName(),
67 mesh,
68 IOobject::READ_IF_PRESENT,
69 IOobject::AUTO_WRITE
70 ),
71 fvMesh_,
72 dimensionedVector(pointMotionU_.dimensions(), Zero),
73 cellMotionBoundaryTypes<vector>(pointMotionU_.boundaryField())
74 ),
75 interpolationPtr_
76 (
77 coeffDict().found("interpolation")
78 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
79 : motionInterpolation::New(fvMesh_)
80 ),
81 diffusivityPtr_
82 (
83 motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
84 )
85{}
86
87
88// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89
91{}
92
93
94// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95
98{
99 interpolationPtr_->interpolate
100 (
101 cellMotionU_,
102 pointMotionU_
103 );
104
105 tmp<pointField> tcurPoints
106 (
107 fvMesh_.points()
108 + fvMesh_.time().deltaTValue()*pointMotionU_.primitiveField()
109 );
110
111 twoDCorrectPoints(tcurPoints.ref());
112
113 return tcurPoints;
114}
115
116
118{
119 // The points have moved so before interpolation update
120 // the fvMotionSolver accordingly
121 movePoints(fvMesh_.points());
122
123 diffusivityPtr_->correct();
124 pointMotionU_.boundaryFieldRef().updateCoeffs();
125
127
128 const label nNonOrthCorr
129 (
130 getOrDefault<label>("nNonOrthogonalCorrectors", 1)
131 );
132
133 for (label i=0; i<nNonOrthCorr; ++i)
134 {
136 (
138 (
139 dimensionedScalar("viscosity", dimViscosity, 1.0)
140 * diffusivityPtr_->operator()(),
141 cellMotionU_,
142 "laplacian(diffusivity,cellMotionU)"
143 )
144 ==
145 fvOptions(cellMotionU_)
146 );
147
150 fvOptions.correct(cellMotionU_);
151 }
152}
153
154
156(
157 const mapPolyMesh& mpm
158)
159{
161
162 // Update diffusivity. Note two stage to make sure old one is de-registered
163 // before creating/registering new one.
164 diffusivityPtr_.reset(nullptr);
165 diffusivityPtr_ = motionDiffusivity::New
166 (
167 fvMesh_,
168 coeffDict().lookup("diffusivity")
169 );
170}
171
172
173// ************************************************************************* //
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...
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
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
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
Base class for fvMesh based motionSolvers.
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
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
Mesh motion solver for an fvMesh. Based on solving the cell-centre Laplacian for the motion velocity.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Virtual base class for velocity motion solver.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
fvVectorMatrix & UEqn
Definition: UEqn.H:13
dynamicFvMesh & mesh
Calculate the matrix for the laplacian of the field.
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.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const int nNonOrthCorr
dictionary dict