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 -------------------------------------------------------------------------------
12 License
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 
30 #include "laplacianMotionSolver.H"
31 #include "motionInterpolation.H"
33 #include "fvmLaplacian.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(laplacianMotionSolver, 1);
40 
42  (
43  motionSolver,
44  laplacianMotionSolver,
45  dictionary
46  );
47 }
48 
49 
50 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51 
53 {
54  pointMotionU_.boundaryFieldRef().updateCoeffs();
55  auto& cellMotionUbf = cellMotionU_.boundaryFieldRef();
56 
58  {
59  fvPatchVectorField& bField = cellMotionUbf[pI];
60  if (isA<fixedValueFvPatchVectorField>(bField))
61  {
62  const pointField& points = fvMesh_.points();
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 
75 Foam::laplacianMotionSolver::laplacianMotionSolver
76 (
77  const polyMesh& mesh,
78  const IOdictionary& dict
79 )
80 :
81  motionSolver(mesh, dict, typeName),
83  pointMotionU_
84  (
85  IOobject
86  (
87  "pointMotionU",
88  mesh.time().timeName(),
89  mesh,
92  ),
95  fixedValuePointPatchVectorField::typeName
96  ),
97  cellMotionU_
98  (
99  IOobject
100  (
101  "cellMotionU",
102  mesh.time().timeName(),
103  mesh,
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 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::fvPatchField< vector >
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::laplacianMotionSolver::curPoints
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Definition: laplacianMotionSolver.C:124
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::laplacianMotionSolver::cellMotionU_
volVectorField cellMotionU_
Definition: laplacianMotionSolver.H:77
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
laplacianMotionSolver.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::fvMotionSolver
Base class for fvMesh based motionSolvers.
Definition: fvMotionSolver.H:50
Foam::laplacianMotionSolver::solve
virtual void solve()
Solve for motion.
Definition: laplacianMotionSolver.C:143
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::laplacianMotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update the mesh corresponding to given map.
Definition: laplacianMotionSolver.C:178
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::radiation::lookup
Lookup type of boundary radiation properties.
Definition: lookup.H:63
dict
dictionary dict
Definition: searchingEngine.H:14
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::laplacianMotionSolver::pointMotionU_
pointVectorField pointMotionU_
Definition: laplacianMotionSolver.H:76
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::motionInterpolation
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers,...
Definition: motionInterpolation.H:52
Foam::laplacianMotionSolver::setBoundaryConditions
void setBoundaryConditions()
Set boundary conditions of cellMotionU based on pointMotionU.
Definition: laplacianMotionSolver.C:52
Foam::motionSolver
Virtual base class for mesh motion solver.
Definition: motionSolver.H:58
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::motionInterpolation::New
static autoPtr< motionInterpolation > New(const fvMesh &mesh)
Select default.
Definition: motionInterpolation.C:66
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
motionInterpolation.H
Foam::laplacianMotionSolver::movePoints
virtual void movePoints(const pointField &)
Update local data for geometry changes.
Definition: laplacianMotionSolver.C:172
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::fvMotionSolver::fvMesh_
const fvMesh & fvMesh_
The fvMesh to be moved.
Definition: fvMotionSolver.H:58