Loading...
Searching...
No Matches
laplacianMotionSolver.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) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
31#include "motionInterpolation.H"
32#include "motionDiffusivity.H"
33#include "fvmLaplacian.H"
34#include "syncTools.H"
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
42
44 (
48 );
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
54Foam::laplacianMotionSolver::laplacianMotionSolver
55(
56 const polyMesh& mesh,
57 const IOdictionary& dict
58)
59:
62 pointMotionU_
63 (
65 (
66 "pointMotionU",
67 mesh.time().timeName(),
68 mesh,
69 IOobject::READ_IF_PRESENT,
70 IOobject::AUTO_WRITE
71 ),
74 fixedValuePointPatchVectorField::typeName
75 ),
76 cellMotionU_
77 (
79 (
80 "cellMotionU",
81 mesh.time().timeName(),
82 mesh,
83 IOobject::READ_IF_PRESENT,
84 IOobject::AUTO_WRITE
85 ),
86 fvMesh_,
87 dimensionedVector(pointMotionU_.dimensions(), Zero),
88 pointMotionU_.boundaryField().types()
89 ),
90 interpolationPtr_
91 (
92 coeffDict().found("interpolation")
93 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
94 : motionInterpolation::New(fvMesh_)
95 ),
96 diffusivityPtr_
97 (
98 motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
99 ),
100 nIters_(this->coeffDict().get<label>("iters")),
101 tolerance_(this->coeffDict().get<scalar>("tolerance"))
102{}
103
104
105// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106
108{
109 interpolationPtr_->interpolate
110 (
113 );
114
116 (
117 fvMesh_,
121 );
122
123 tmp<vectorField> tcurPoints
124 (
126 );
128 twoDCorrectPoints(tcurPoints.ref());
129
130 return tcurPoints;
131}
132
133
135{
136 diffusivityPtr_->correct();
137
138 // Iteratively solve the Laplace equation, to account for non-orthogonality
139 for (label iter = 0; iter < nIters_; ++iter)
140 {
141 Info<< "Iteration " << iter << endl;
142 fvVectorMatrix dEqn
143 (
145 (
146 dimensionedScalar("viscosity", dimViscosity, 1.0)
147 * diffusivityPtr_->operator()(),
148 cellMotionU_,
149 "laplacian(diffusivity,cellMotionU)"
150 )
151 );
152
153 scalar residual = mag(dEqn.solve().initialResidual());
154
155 // Print execution time
156 fvMesh_.time().printExecutionTime(Info);
157
158 // Check convergence
159 if (residual < tolerance_)
160 {
161 Info<< "\n***Reached mesh movement convergence limit at"
162 << " iteration " << iter << "***\n\n";
163 break;
164 }
165 }
166}
167
168
170{
171 pointMotionU_.boundaryFieldRef().updateCoeffs();
172 auto& cellMotionUbf = cellMotionU_.boundaryFieldRef();
173
174 forAll(cellMotionU_.boundaryField(), pI)
175 {
176 fvPatchVectorField& bField = cellMotionUbf[pI];
178 {
179 const pointField& points = fvMesh_.points();
180 const polyPatch& patch = fvMesh_.boundaryMesh()[pI];
181 forAll(bField, fI)
182 {
183 bField[fI] = patch[fI].average(points, pointMotionU_);
184 }
185 }
186 }
187}
188
191{
192 // Do nothing
193}
194
195
197{
198 // Update diffusivity. Note two stage to make sure old one is de-registered
199 // before creating/registering new one.
200 diffusivityPtr_.reset(nullptr);
201 diffusivityPtr_ = motionDiffusivity::New
202 (
203 fvMesh_,
204 coeffDict().lookup("diffusivity")
205 );
206}
207
208
209// ************************************************************************* //
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.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
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...
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ 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
const Type & initialResidual() const noexcept
Return initial residual.
static const Form zero
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.
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Base class for fvMesh based motionSolvers.
const fvMesh & fvMesh_
The fvMesh to be moved.
fvMotionSolver(const polyMesh &)
Construct from polyMesh.
Similar to velocityLaplacian but iteratively solves the mesh displacement PDEs to account for non-ort...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
autoPtr< motionDiffusivity > diffusivityPtr_
Diffusivity used to control the motion.
scalar tolerance_
Residual threshold.
label nIters_
Number of laplacian iterations per solution step.
void setBoundaryConditions()
Set boundary conditions of cellMotionU based on pointMotionU.
autoPtr< motionInterpolation > interpolationPtr_
Interpolation used to transfer cell displacement to the points.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual void updateMesh(const mapPolyMesh &)
Update the mesh corresponding to given map.
virtual void solve()
Solve for motion.
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.
const polyMesh & mesh() const
Return reference to mesh.
static autoPtr< motionSolver > New(const polyMesh &)
Select constructed from polyMesh.
motionSolver(const polyMesh &mesh)
Construct from polyMesh.
virtual void twoDCorrectPoints(pointField &) const
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Lookup type of boundary radiation properties.
Definition lookup.H:60
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
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
#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.
const pointField & points
word timeName
Definition getTimeIndex.H:3
const std::string patch
OpenFOAM patch number as a std::string.
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.
const dimensionSet dimless
Dimensionless.
messageStream Info
Information stream (stdout output on master, null elsewhere).
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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.
fvPatchField< vector > fvPatchVectorField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299