displacementComponentLaplacianFvMotionSolver.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 "polyTopoChangeMap.H"
31 #include "volPointInterpolation.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
38 
40  (
44  );
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
52 (
53  const word& name,
54  const polyMesh& mesh,
55  const dictionary& dict
56 )
57 :
60  cellDisplacement_
61  (
62  IOobject
63  (
64  "cellDisplacement" + cmptName_,
65  mesh.time().name(),
66  mesh,
67  IOobject::READ_IF_PRESENT,
68  IOobject::AUTO_WRITE
69  ),
70  fvMesh_,
71  dimensionedScalar(pointDisplacement_.dimensions(), 0),
72  cellMotionBoundaryTypes<scalar>(pointDisplacement_.boundaryField())
73  ),
74  pointLocation_(nullptr),
75  diffusivityType_(dict.lookup("diffusivity")),
76  diffusivityPtr_
77  (
78  motionDiffusivity::New(fvMesh_, diffusivityType_)
79  ),
80  frozenPointsZone_
81  (
82  dict.found("frozenPointsZone")
83  ? fvMesh_.pointZones().findIndex(dict.lookup("frozenPointsZone"))
84  : -1
85  )
86 {
88  (
89  "pointLocation",
90  fvMesh_.time().name(),
91  fvMesh_,
94  );
95 
96  if (debug)
97  {
98  Info<< "displacementComponentLaplacianFvMotionSolver:" << nl
99  << " diffusivity : " << diffusivityPtr_().type() << nl
100  << " frozenPoints zone : " << frozenPointsZone_ << endl;
101  }
102 
103  if (io.headerOk())
104  {
105  pointLocation_.reset
106  (
107  new pointVectorField
108  (
109  IOobject
110  (
111  "pointLocation",
112  fvMesh_.time().name(),
113  fvMesh_,
116  ),
118  )
119  );
120 
121  if (debug)
122  {
123  Info<< "displacementComponentLaplacianFvMotionSolver :"
124  << " Read pointVectorField "
125  << pointLocation_().name()
126  << " to be used for boundary conditions on points."
127  << nl
128  << "Boundary conditions:"
129  << pointLocation_().boundaryField().types() << endl;
130  }
131  }
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
136 
139 {}
140 
141 
142 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
143 
146 {
148  (
149  cellDisplacement_,
150  pointDisplacement_
151  );
152 
153  if (pointLocation_.valid())
154  {
155  if (debug)
156  {
157  Info<< "displacementComponentLaplacianFvMotionSolver : applying "
158  << " boundary conditions on " << pointLocation_().name()
159  << " to new point location."
160  << endl;
161  }
162 
163  // Apply pointLocation_ b.c. to mesh points.
164 
165  pointLocation_().primitiveFieldRef() = fvMesh_.points();
166 
167  pointLocation_().primitiveFieldRef().replace
168  (
169  cmpt_,
170  points0_ + pointDisplacement_.primitiveField()
171  );
172 
173  pointLocation_().correctBoundaryConditions();
174 
175  // Implement frozen points
176  if (frozenPointsZone_ != -1)
177  {
178  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
179 
180  forAll(pz, i)
181  {
182  label pointi = pz[i];
183 
184  pointLocation_()[pointi][cmpt_] = points0_[pointi];
185  }
186  }
187 
188  twoDCorrectPoints(pointLocation_().primitiveFieldRef());
189 
190  return tmp<pointField>(pointLocation_().primitiveField());
191  }
192  else
193  {
194  tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
195  pointField& curPoints = tcurPoints.ref();
196 
197  curPoints.replace
198  (
199  cmpt_,
200  points0_ + pointDisplacement_.primitiveField()
201  );
202 
203  // Implement frozen points
204  if (frozenPointsZone_ != -1)
205  {
206  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
207 
208  forAll(pz, i)
209  {
210  label pointi = pz[i];
211 
212  curPoints[pointi][cmpt_] = points0_[pointi];
213  }
214  }
215 
216  twoDCorrectPoints(curPoints);
217 
218  return tcurPoints;
219  }
220 }
221 
222 
224 {
225  // The points have moved so before interpolation update
226  // the motionSolver accordingly
227  movePoints(fvMesh_.points());
228 
229  diffusivityPtr_->correct();
230  pointDisplacement_.boundaryFieldRef().updateCoeffs();
231 
233  (
235  (
236  diffusivityPtr_->operator()(),
237  cellDisplacement_,
238  "laplacian(diffusivity,cellDisplacement)"
239  )
240  );
241 }
242 
243 
245 (
246  const polyTopoChangeMap& map
247 )
248 {
250 
251  // Update diffusivity. Note two stage to make sure old one is de-registered
252  // before creating/registering new one.
253  diffusivityPtr_.reset(nullptr);
254  diffusivityType_.rewind();
255  diffusivityPtr_ = motionDiffusivity::New
256  (
257  fvMesh_,
258  diffusivityType_
259  );
260 }
261 
262 
264 (
265  const polyMeshMap& map
266 )
267 {
269 
270  // Update diffusivity. Note two stage to make sure old one is de-registered
271  // before creating/registering new one.
272  diffusivityPtr_.reset(nullptr);
273  diffusivityType_.rewind();
274  diffusivityPtr_ = motionDiffusivity::New
275  (
276  fvMesh_,
277  diffusivityType_
278  );
279 }
280 
281 
282 // ************************************************************************* //
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.
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:498
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Virtual base class for displacement motion solver.
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
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 given component ...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void topoChange(const polyTopoChangeMap &)
Update corresponding to the given map.
displacementComponentLaplacianFvMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
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
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.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
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
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dictionary dict
conserve primitiveFieldRef()+