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-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 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"
33#include "fvmLaplacian.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
40
42 (
46 );
47}
48
49
50// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51
53{
55 auto& cellMotionUbf = cellMotionU_.boundaryFieldRef();
56
58 {
59 fvPatchVectorField& bField = cellMotionUbf[pI];
60 if (isA<fixedValueFvPatchVectorField>(bField))
61 {
63 const polyPatch& patch = fvMesh_.boundaryMesh()[pI];
64 forAll(bField, fI)
65 {
66 bField[fI] = patch[fI].average(points, pointMotionU_);
67 }
68 }
69 }
70}
71
72
73// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
74
76(
77 const polyMesh& mesh,
78 const IOdictionary& dict
79)
80:
81 motionSolver(mesh, dict, typeName),
83 pointMotionU_
84 (
86 (
87 "pointMotionU",
88 mesh.time().timeName(),
89 mesh,
90 IOobject::READ_IF_PRESENT,
91 IOobject::AUTO_WRITE
92 ),
95 fixedValuePointPatchVectorField::typeName
96 ),
97 cellMotionU_
98 (
100 (
101 "cellMotionU",
102 mesh.time().timeName(),
103 mesh,
104 IOobject::READ_IF_PRESENT,
105 IOobject::AUTO_WRITE
106 ),
107 fvMesh_,
108 dimensionedVector(pointMotionU_.dimensions(), Zero),
109 pointMotionU_.boundaryField().types()
110 ),
111 interpolationPtr_
112 (
113 coeffDict().found("interpolation")
114 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
115 : motionInterpolation::New(fvMesh_)
116 ),
117 nIters_(this->coeffDict().get<label>("iters")),
118 tolerance_(this->coeffDict().get<scalar>("tolerance"))
119{}
120
121
122// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123
125{
126 interpolationPtr_->interpolate
127 (
128 cellMotionU_,
129 pointMotionU_
130 );
131
132 tmp<vectorField> tcurPoints
133 (
134 fvMesh_.points() + pointMotionU_.internalField()
135 );
136
137 twoDCorrectPoints(tcurPoints.ref());
138
139 return tcurPoints;
140}
141
142
144{
145 setBoundaryConditions();
146
147 // Iteratively solve the Laplace equation, to account for non-orthogonality
148 for (label iter = 0; iter < nIters_; ++iter)
149 {
150 Info<< "Iteration " << iter << endl;
151 fvVectorMatrix dEqn
152 (
153 fvm::laplacian(cellMotionU_)
154 );
155
156 scalar residual = mag(dEqn.solve().initialResidual());
157
158 // Print execution time
159 fvMesh_.time().printExecutionTime(Info);
160
161 // Check convergence
162 if (residual < tolerance_)
163 {
164 Info<< "\n***Reached mesh movement convergence limit at"
165 << " iteration " << iter << "***\n\n";
166 break;
167 }
168 }
169}
170
171
173{
174 // Do nothing
175}
176
177
179{
180 // Do nothing
181}
182
183
184// ************************************************************************* //
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.
void updateCoeffs()
Update the boundary condition coefficients.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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
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 > solve(const dictionary &)
Solve returning the solution statistics.
Base class for fvMesh based motionSolvers.
const fvMesh & fvMesh_
The fvMesh to be moved.
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.
void setBoundaryConditions()
Set boundary conditions of cellMotionU based on pointMotionU.
virtual void solve()
Solve for motion.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers,...
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:61
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:55
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
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
#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.
const pointField & points
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 dimless
Dimensionless.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333