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 -------------------------------------------------------------------------------
11 License
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 
29 #include "pLaplacianMotionSolver.H"
30 #include "motionInterpolation.H"
32 #include "syncTools.H"
33 #include "fvmLaplacian.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(pLaplacianMotionSolver, 1);
40 
42  (
43  motionSolver,
44  pLaplacianMotionSolver,
45  dictionary
46  );
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 Foam::pLaplacianMotionSolver::pLaplacianMotionSolver
53 (
54  const polyMesh& mesh,
55  const IOdictionary& dict
56 )
57 :
58  motionSolver(mesh, dict, typeName),
60  useFixedValuePointMotionUBCs_
61  (coeffDict().getOrDefault<bool>("useFixedValuePointMotionUBCs", false)),
62  pointMotionU_
63  (
64  IOobject
65  (
66  "pointMotionU",
67  mesh.time().timeName(),
68  mesh,
69  IOobject::READ_IF_PRESENT,
70  IOobject::AUTO_WRITE
71  ),
72  pointMesh::New(mesh),
74  word
75  (
76  useFixedValuePointMotionUBCs_
77  ? fixedValuePointPatchVectorField::typeName
79  )
80  ),
81  cellMotionU_
82  (
83  IOobject
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];
122  if (isA<fixedValueFvPatchVectorField>(bField))
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_);
129  }
130  }
131  }
132 }
133 
134 
135 
137 {
138  interpolationPtr_->interpolate
139  (
140  cellMotionU_,
141  pointMotionU_
142  );
143 
145  (
146  fvMesh_,
147  pointMotionU_.primitiveFieldRef(),
148  maxEqOp<vector>(),
149  vector::zero
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 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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.
Calculate the matrix for the laplacian of the field.
virtual void updateMesh(const mapPolyMesh &)
Update the mesh corresponding to given map.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:54
virtual void solve()
Solve for motion.
const dimensionSet dimless
Dimensionless.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
A calculated boundary condition for pointField.
void setBoundaryConditions()
Set boundary conditions of cellMotionU based on pointMotionU.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
Macros for easy insertion into run-time selection tables.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Base class for fvMesh based motionSolvers.
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
dimensionedScalar exp(const dimensionedScalar &ds)
const pointField & points
A class for handling words, derived from Foam::string.
Definition: word.H:63
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
virtual void movePoints(const pointField &)
Update local data for geometry changes.
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers, to the points. This base class implements the default method which applies volPointInterpolation only.
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
const scalar gamma
Definition: EEqn.H:9
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: HashPtrTable.H:50
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
bool found
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127