displacementLaplacianFvMotionSolver.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "motionDiffusivity.H"
28 #include "fvmLaplacian.H"
30 #include "OFstream.H"
31 #include "meshTools.H"
32 #include "polyTopoChangeMap.H"
33 #include "volPointInterpolation.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
40 
42  (
46  );
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const word& name,
55  const polyMesh& mesh,
56  const dictionary& dict
57 )
58 :
61  cellDisplacement_
62  (
63  IOobject
64  (
65  "cellDisplacement",
66  mesh.time().name(),
67  mesh,
68  IOobject::READ_IF_PRESENT,
69  IOobject::AUTO_WRITE
70  ),
71  fvMesh_,
73  (
74  "cellDisplacement",
75  pointDisplacement_.dimensions(),
76  Zero
77  ),
78  cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
79  ),
80  pointLocation_(nullptr),
81  diffusivityType_(dict.lookup("diffusivity")),
82  diffusivityPtr_
83  (
84  motionDiffusivity::New(fvMesh_, diffusivityType_)
85  ),
86  frozenPointsZone_
87  (
88  dict.found("frozenPointsZone")
89  ? fvMesh_.pointZones().findIndex(dict.lookup("frozenPointsZone"))
90  : -1
91  )
92 {
94  (
95  "pointLocation",
96  fvMesh_.time().name(),
97  fvMesh_,
100  );
101 
102  if (debug)
103  {
104  Info<< "displacementLaplacianFvMotionSolver:" << nl
105  << " diffusivity : " << diffusivityPtr_().type() << nl
106  << " frozenPoints zone : " << frozenPointsZone_ << endl;
107  }
108 
109  if (io.headerOk())
110  {
111  pointLocation_.reset
112  (
113  new pointVectorField
114  (
115  io,
117  )
118  );
119 
120  if (debug)
121  {
122  Info<< "displacementLaplacianFvMotionSolver :"
123  << " Read pointVectorField "
124  << io.name()
125  << " to be used for boundary conditions on points."
126  << nl
127  << "Boundary conditions:"
128  << pointLocation_().boundaryField().types() << endl;
129  }
130  }
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
145 {
146  if (!diffusivityPtr_.valid())
147  {
148  diffusivityType_.rewind();
149  diffusivityPtr_ = motionDiffusivity::New
150  (
151  fvMesh_,
152  diffusivityType_
153  );
154  }
155  return diffusivityPtr_();
156 }
157 
158 
161 {
163  (
164  cellDisplacement_,
165  pointDisplacement_
166  );
167 
168  if (pointLocation_.valid())
169  {
170  if (debug)
171  {
172  Info<< "displacementLaplacianFvMotionSolver : applying "
173  << " boundary conditions on " << pointLocation_().name()
174  << " to new point location."
175  << endl;
176  }
177 
178  pointLocation_().primitiveFieldRef() =
179  points0()
180  + pointDisplacement_.primitiveField();
181 
182  pointLocation_().correctBoundaryConditions();
183 
184  // Implement frozen points
185  if (frozenPointsZone_ != -1)
186  {
187  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
188 
189  forAll(pz, i)
190  {
191  pointLocation_()[pz[i]] = points0()[pz[i]];
192  }
193  }
194 
195  twoDCorrectPoints(pointLocation_().primitiveFieldRef());
196 
197  return tmp<pointField>(pointLocation_().primitiveField());
198  }
199  else
200  {
201  tmp<pointField> tcurPoints
202  (
203  points0() + pointDisplacement_.primitiveField()
204  );
205  pointField& curPoints = tcurPoints.ref();
206 
207  // Implement frozen points
208  if (frozenPointsZone_ != -1)
209  {
210  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
211 
212  forAll(pz, i)
213  {
214  curPoints[pz[i]] = points0()[pz[i]];
215  }
216  }
217 
218  twoDCorrectPoints(curPoints);
219 
220  return tcurPoints;
221  }
222 }
223 
224 
226 {
227  // The points have moved so before interpolation update
228  // the motionSolver accordingly
229  movePoints(fvMesh_.points());
230 
231  diffusivity().correct();
232  pointDisplacement_.boundaryFieldRef().updateCoeffs();
233 
235  (
237  (
238  diffusivity().operator()(),
239  cellDisplacement_,
240  "laplacian(diffusivity,cellDisplacement)"
241  )
242  );
243 }
244 
245 
247 (
248  const polyTopoChangeMap& map
249 )
250 {
252  diffusivityPtr_.clear();
253 }
254 
255 
257 (
258  const polyMeshMap& map
259 )
260 {
262  diffusivityPtr_.clear();
263 }
264 
265 
266 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:307
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const word & name() const
Return const reference to name.
Mesh motion solver for an fvMesh. Based on solving the cell-centre Laplacian for the motion displacem...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
displacementLaplacianFvMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
virtual void topoChange(const polyTopoChangeMap &)
Update corresponding to the given map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
motionDiffusivity & diffusivity()
Return reference to the diffusivity field.
Virtual base class for displacement motion solver.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
Base class for fvMesh based motionSolvers.
const fvMesh & fvMesh_
The fvMesh to be moved.
Abstract base class for cell-centre mesh motion diffusivity.
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
Named list of point indices representing a sub-set of the mesh faces.
Definition: pointZone.H:60
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:530
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
tmp< PointField< Type > > interpolate(const VolField< Type > &) const
Interpolate volField using inverse distance weighting.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the matrix for the laplacian of the field.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
defineTypeNameAndDebug(combustionModel, 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
static const char nl
Definition: Ostream.H:267
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dictionary dict
conserve primitiveFieldRef()+