Loading...
Searching...
No Matches
pLaplacianMotionSolver.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) 2021-2023 PCOpt/NTUA
9 Copyright (C) 2021-2023 FOSS GP
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"
32#include "syncTools.H"
33#include "fvmLaplacian.H"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
40
42 (
46 );
47}
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
52Foam::pLaplacianMotionSolver::pLaplacianMotionSolver
53(
54 const polyMesh& mesh,
55 const IOdictionary& dict
56)
57:
60 useFixedValuePointMotionUBCs_
61 (coeffDict().getOrDefault<bool>("useFixedValuePointMotionUBCs", false)),
62 pointMotionU_
63 (
65 (
66 "pointMotionU",
67 mesh.time().timeName(),
68 mesh,
69 IOobject::READ_IF_PRESENT,
70 IOobject::AUTO_WRITE
71 ),
74 word
75 (
76 useFixedValuePointMotionUBCs_
77 ? fixedValuePointPatchVectorField::typeName
79 )
80 ),
81 cellMotionU_
82 (
84 (
85 "cellMotionU",
86 mesh.time().timeName(),
87 mesh,
88 IOobject::READ_IF_PRESENT,
89 IOobject::AUTO_WRITE
90 ),
91 fvMesh_,
92 dimensionedVector(pointMotionU_.dimensions(), Zero),
93 pointMotionU_.boundaryField().types()
94 ),
95 interpolationPtr_
96 (
97 coeffDict().found("interpolation")
98 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
99 : motionInterpolation::New(fvMesh_)
100 ),
101 nIters_(this->coeffDict().get<label>("iters")),
102 tolerance_(this->coeffDict().get<scalar>("tolerance")),
103 toleranceIntermediate_
104 (
105 this->coeffDict().
106 getOrDefault<scalar>("toleranceIntermediate", 100*tolerance_)
107 ),
108 exponent_(this->coeffDict().get<label>("exponent"))
109{}
110
111
112// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113
115{
116 pointMotionU_.boundaryFieldRef().updateCoeffs();
117 auto& cellMotionUbf = cellMotionU_.boundaryFieldRef();
118
119 forAll(cellMotionU_.boundaryField(), pI)
120 {
121 fvPatchVectorField& bField = cellMotionUbf[pI];
123 {
124 const pointField& points = fvMesh_.points();
125 const polyPatch& patch = fvMesh_.boundaryMesh()[pI];
126 forAll(bField, fI)
127 {
128 bField[fI] = patch[fI].average(points, pointMotionU_);
130 }
131 }
132}
133
134
135
136Foam::tmp<Foam::pointField> Foam::pLaplacianMotionSolver::curPoints() const
137{
138 interpolationPtr_->interpolate
139 (
140 cellMotionU_,
141 pointMotionU_
142 );
143
145 (
146 fvMesh_,
147 pointMotionU_.primitiveFieldRef(),
150 );
151
152 tmp<vectorField> tcurPoints
153 (
154 fvMesh_.points() + pointMotionU_.internalField()
155 );
157 twoDCorrectPoints(tcurPoints.ref());
158
159 return tcurPoints;
160}
161
162
164{
165// setBoundaryConditions();
166
167 for (label exp = 2; exp < exponent_ + 1; ++exp)
168 {
169 scalar tolerance
170 (exp == exponent_ ? tolerance_ : toleranceIntermediate_);
171 Info<< "Solving for exponent " << exp << endl;
172
173 for (label iter = 0; iter < nIters_; ++iter)
174 {
175 Info<< "Iteration " << iter << endl;
176 cellMotionU_.storePrevIter();
177 volScalarField gamma(pow(mag(fvc::grad(cellMotionU_)), exp - 2));
178 gamma.correctBoundaryConditions();
179 fvVectorMatrix dEqn
180 (
181 fvm::laplacian(gamma, cellMotionU_)
182 );
183
184 scalar residual = mag(dEqn.solve().initialResidual());
185
186 cellMotionU_.relax();
187
188 // Print execution time
189 fvMesh_.time().printExecutionTime(Info);
190
191 // Check convergence
192 if (residual < tolerance)
193 {
194 Info<< "\n***Reached mesh movement convergence limit at"
195 << " iteration " << iter << "***\n\n";
196 if (debug)
197 {
198 gamma.write();
199 }
200 break;
201 }
202 }
203 }
204}
205
208{
209 // Do nothing
210}
211
212
214{
215 // Do nothing
216}
217
218
219// ************************************************************************* //
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.
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 calculated boundary condition for pointField.
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...
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.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
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.
Similar to velocityLaplacian but with a variable diffusivity, based on the gradient of the displaceme...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
scalar tolerance_
Residual threshold.
label nIters_
Number of pLaplacian iterations per solution step.
void setBoundaryConditions()
Set boundary conditions of cellMotionU based on pointMotionU.
label exponent_
Exponent defining the order or the p-Laplacian.
autoPtr< motionInterpolation > interpolationPtr_
Interpolation used to transfer cell displacement to the points.
scalar toleranceIntermediate_
Residual threshold for intermediate exponents.
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.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
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
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
const scalar gamma
Definition EEqn.H:9
dynamicFvMesh & mesh
Calculate the matrix for the laplacian of the field.
const pointField & points
word timeName
Definition getTimeIndex.H:3
Namespace for handling debugging switches.
Definition debug.C:45
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Namespace for OpenFOAM.
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.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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.
Vector< scalar > vector
Definition vector.H:57
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