displacementLaplacianFvMotionSolver.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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "motionInterpolation.H"
28 #include "motionDiffusivity.H"
29 #include "fvmLaplacian.H"
31 #include "OFstream.H"
32 #include "meshTools.H"
33 #include "mapPolyMesh.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(displacementLaplacianFvMotionSolver, 0);
40 
42  (
43  motionSolver,
44  displacementLaplacianFvMotionSolver,
45  dictionary
46  );
47 
49  (
50  displacementMotionSolver,
51  displacementLaplacianFvMotionSolver,
52  displacement
53  );
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
60 (
61  const polyMesh& mesh,
62  const IOdictionary& dict
63 )
64 :
67  cellDisplacement_
68  (
69  IOobject
70  (
71  "cellDisplacement",
72  mesh.time().timeName(),
73  mesh,
74  IOobject::READ_IF_PRESENT,
75  IOobject::AUTO_WRITE
76  ),
77  fvMesh_,
79  (
80  "cellDisplacement",
81  pointDisplacement_.dimensions(),
82  vector::zero
83  ),
84  cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
85  ),
86  pointLocation_(NULL),
87  interpolationPtr_
88  (
89  coeffDict().found("interpolation")
90  ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
91  : motionInterpolation::New(fvMesh_)
92  ),
93  diffusivityPtr_
94  (
95  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
96  ),
97  frozenPointsZone_
98  (
99  coeffDict().found("frozenPointsZone")
100  ? fvMesh_.pointZones().findZoneID(coeffDict().lookup("frozenPointsZone"))
101  : -1
102  )
103 {
104  IOobject io
105  (
106  "pointLocation",
107  fvMesh_.time().timeName(),
108  fvMesh_,
109  IOobject::MUST_READ,
110  IOobject::AUTO_WRITE
111  );
112 
113  if (debug)
114  {
115  Info<< "displacementLaplacianFvMotionSolver:" << nl
116  << " diffusivity : " << diffusivityPtr_().type() << nl
117  << " frozenPoints zone : " << frozenPointsZone_ << endl;
118  }
119 
120 
121  if (io.headerOk())
122  {
123  pointLocation_.reset
124  (
125  new pointVectorField
126  (
127  io,
128  pointMesh::New(fvMesh_)
129  )
130  );
131 
132  if (debug)
133  {
134  Info<< "displacementLaplacianFvMotionSolver :"
135  << " Read pointVectorField "
136  << io.name()
137  << " to be used for boundary conditions on points."
138  << nl
139  << "Boundary conditions:"
140  << pointLocation_().boundaryField().types() << endl;
141  }
142  }
143 }
144 
145 
148 (
149  const polyMesh& mesh,
150  const IOdictionary& dict,
151  const pointVectorField& pointDisplacement,
152  const pointIOField& points0
153 )
154 :
155  displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName),
157  cellDisplacement_
158  (
159  IOobject
160  (
161  "cellDisplacement",
162  mesh.time().timeName(),
163  mesh,
164  IOobject::READ_IF_PRESENT,
165  IOobject::AUTO_WRITE
166  ),
167  fvMesh_,
169  (
170  "cellDisplacement",
171  pointDisplacement_.dimensions(),
172  vector::zero
173  ),
174  cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
175  ),
176  pointLocation_(NULL),
177  interpolationPtr_
178  (
179  coeffDict().found("interpolation")
180  ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
181  : motionInterpolation::New(fvMesh_)
182  ),
183  diffusivityPtr_
184  (
185  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
186  ),
187  frozenPointsZone_
188  (
189  coeffDict().found("frozenPointsZone")
190  ? fvMesh_.pointZones().findZoneID(coeffDict().lookup("frozenPointsZone"))
191  : -1
192  )
193 {
194  IOobject io
195  (
196  "pointLocation",
197  fvMesh_.time().timeName(),
198  fvMesh_,
199  IOobject::MUST_READ,
200  IOobject::AUTO_WRITE
201  );
202 
203  if (debug)
204  {
205  Info<< "displacementLaplacianFvMotionSolver:" << nl
206  << " diffusivity : " << diffusivityPtr_().type() << nl
207  << " frozenPoints zone : " << frozenPointsZone_ << endl;
208  }
209 
210 
211  if (io.headerOk())
212  {
213  pointLocation_.reset
214  (
215  new pointVectorField
216  (
217  io,
218  pointMesh::New(fvMesh_)
219  )
220  );
221 
222  if (debug)
223  {
224  Info<< "displacementLaplacianFvMotionSolver :"
225  << " Read pointVectorField "
226  << io.name()
227  << " to be used for boundary conditions on points."
228  << nl
229  << "Boundary conditions:"
230  << pointLocation_().boundaryField().types() << endl;
231  }
232  }
233 }
234 
235 
236 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
237 
240 {}
241 
242 
243 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
244 
247 {
248  if (!diffusivityPtr_.valid())
249  {
250  diffusivityPtr_ = motionDiffusivity::New
251  (
252  fvMesh_,
253  coeffDict().lookup("diffusivity")
254  );
255  }
256  return diffusivityPtr_();
257 }
258 
259 
262 {
263  interpolationPtr_->interpolate
264  (
265  cellDisplacement_,
266  pointDisplacement_
267  );
268 
269  if (pointLocation_.valid())
270  {
271  if (debug)
272  {
273  Info<< "displacementLaplacianFvMotionSolver : applying "
274  << " boundary conditions on " << pointLocation_().name()
275  << " to new point location."
276  << endl;
277  }
278 
279  pointLocation_().internalField() =
280  points0()
281  + pointDisplacement_.internalField();
282 
283  pointLocation_().correctBoundaryConditions();
284 
285  // Implement frozen points
286  if (frozenPointsZone_ != -1)
287  {
288  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
289 
290  forAll(pz, i)
291  {
292  pointLocation_()[pz[i]] = points0()[pz[i]];
293  }
294  }
295 
296  twoDCorrectPoints(pointLocation_().internalField());
297 
298  return tmp<pointField>(pointLocation_().internalField());
299  }
300  else
301  {
302  tmp<pointField> tcurPoints
303  (
304  points0() + pointDisplacement_.internalField()
305  );
306 
307  // Implement frozen points
308  if (frozenPointsZone_ != -1)
309  {
310  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
311 
312  forAll(pz, i)
313  {
314  tcurPoints()[pz[i]] = points0()[pz[i]];
315  }
316  }
317 
318  twoDCorrectPoints(tcurPoints());
319 
320  return tcurPoints;
321  }
322 }
323 
324 
326 {
327  // The points have moved so before interpolation update
328  // the motionSolver accordingly
329  movePoints(fvMesh_.points());
330 
331  diffusivity().correct();
332  pointDisplacement_.boundaryField().updateCoeffs();
333 
335  (
337  (
338  diffusivity().operator()(),
339  cellDisplacement_,
340  "laplacian(diffusivity,cellDisplacement)"
341  )
342  );
343 }
344 
345 
347 (
348  const mapPolyMesh& mpm
349 )
350 {
352 
353  // Update diffusivity. Note two stage to make sure old one is de-registered
354  // before creating/registering new one.
355  diffusivityPtr_.clear();
356 }
357 
358 
359 // ************************************************************************* //
Foam::displacementLaplacianFvMotionSolver::~displacementLaplacianFvMotionSolver
~displacementLaplacianFvMotionSolver()
Destructor.
Definition: displacementLaplacianFvMotionSolver.C:239
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
meshTools.H
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::IOField< vector >
Foam::displacementLaplacianFvMotionSolver::solve
virtual void solve()
Solve for motion.
Definition: displacementLaplacianFvMotionSolver.C:325
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::motionDiffusivity
Abstract base class for cell-centre mesh motion diffusivity.
Definition: motionDiffusivity.H:49
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
mapPolyMesh.H
Foam::pointZone
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list.
Definition: pointZone.H:62
Foam::displacementLaplacianFvMotionSolver::displacementLaplacianFvMotionSolver
displacementLaplacianFvMotionSolver(const displacementLaplacianFvMotionSolver &)
Disallow default bitwise copy construct.
displacementLaplacianFvMotionSolver.H
Foam::displacementMotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
Definition: displacementMotionSolver.C:248
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::displacementLaplacianFvMotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update topology.
Definition: displacementLaplacianFvMotionSolver.C:347
Foam::displacementLaplacianFvMotionSolver::diffusivity
motionDiffusivity & diffusivity()
Return reference to the diffusivity field.
Definition: displacementLaplacianFvMotionSolver.C:246
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
OFstream.H
Foam::IOobject::headerOk
bool headerOk()
Read and check header info.
Definition: IOobject.C:439
Foam::fvMotionSolverCore
Base class for fvMesh based motionSolvers.
Definition: fvMotionSolverCore.H:48
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::motionDiffusivity::New
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Definition: motionDiffusivity.C:48
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::displacementMotionSolver
Virtual base class for displacement motion solver.
Definition: displacementMotionSolver.H:55
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::displacementLaplacianFvMotionSolver::curPoints
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Definition: displacementLaplacianFvMotionSolver.C:261
motionDiffusivity.H
internalField
conserve internalField()+
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::solve
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
motionInterpolation.H
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress