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-2019 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 "mapPolyMesh.H"
33 #include "volPointInterpolation.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 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const polyMesh& mesh,
55  const dictionary& dict
56 )
57 :
58  displacementMotionSolver(mesh, dict, typeName),
59  fvMotionSolver(mesh),
60  cellDisplacement_
61  (
62  IOobject
63  (
64  "cellDisplacement",
65  mesh.time().timeName(),
66  mesh,
69  ),
70  fvMesh_,
72  (
73  "cellDisplacement",
74  pointDisplacement_.dimensions(),
75  Zero
76  ),
77  cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
78  ),
79  pointLocation_(nullptr),
80  diffusivityPtr_
81  (
82  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
83  ),
84  frozenPointsZone_
85  (
86  coeffDict().found("frozenPointsZone")
87  ? fvMesh_.pointZones().findZoneID(coeffDict().lookup("frozenPointsZone"))
88  : -1
89  )
90 {
91  IOobject io
92  (
93  "pointLocation",
94  fvMesh_.time().timeName(),
95  fvMesh_,
98  );
99 
100  if (debug)
101  {
102  Info<< "displacementLaplacianFvMotionSolver:" << nl
103  << " diffusivity : " << diffusivityPtr_().type() << nl
104  << " frozenPoints zone : " << frozenPointsZone_ << endl;
105  }
106 
107 
108  if (io.typeHeaderOk<pointVectorField>(true))
109  {
110  pointLocation_.reset
111  (
112  new pointVectorField
113  (
114  io,
115  pointMesh::New(fvMesh_)
116  )
117  );
118 
119  if (debug)
120  {
121  Info<< "displacementLaplacianFvMotionSolver :"
122  << " Read pointVectorField "
123  << io.name()
124  << " to be used for boundary conditions on points."
125  << nl
126  << "Boundary conditions:"
127  << pointLocation_().boundaryField().types() << endl;
128  }
129  }
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
134 
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
144 {
145  if (!diffusivityPtr_.valid())
146  {
147  diffusivityPtr_ = motionDiffusivity::New
148  (
149  fvMesh_,
150  coeffDict().lookup("diffusivity")
151  );
152  }
153  return diffusivityPtr_();
154 }
155 
156 
159 {
161  (
162  cellDisplacement_,
163  pointDisplacement_
164  );
165 
166  if (pointLocation_.valid())
167  {
168  if (debug)
169  {
170  Info<< "displacementLaplacianFvMotionSolver : applying "
171  << " boundary conditions on " << pointLocation_().name()
172  << " to new point location."
173  << endl;
174  }
175 
176  pointLocation_().primitiveFieldRef() =
177  points0()
178  + pointDisplacement_.primitiveField();
179 
180  pointLocation_().correctBoundaryConditions();
181 
182  // Implement frozen points
183  if (frozenPointsZone_ != -1)
184  {
185  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
186 
187  forAll(pz, i)
188  {
189  pointLocation_()[pz[i]] = points0()[pz[i]];
190  }
191  }
192 
193  twoDCorrectPoints(pointLocation_().primitiveFieldRef());
194 
195  return tmp<pointField>(pointLocation_().primitiveField());
196  }
197  else
198  {
199  tmp<pointField> tcurPoints
200  (
201  points0() + pointDisplacement_.primitiveField()
202  );
203  pointField& curPoints = tcurPoints.ref();
204 
205  // Implement frozen points
206  if (frozenPointsZone_ != -1)
207  {
208  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
209 
210  forAll(pz, i)
211  {
212  curPoints[pz[i]] = points0()[pz[i]];
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  diffusivity().correct();
230  pointDisplacement_.boundaryFieldRef().updateCoeffs();
231 
233  (
235  (
236  diffusivity().operator()(),
237  cellDisplacement_,
238  "laplacian(diffusivity,cellDisplacement)"
239  )
240  );
241 }
242 
243 
245 (
246  const mapPolyMesh& mpm
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_.clear();
254 }
255 
256 
257 // ************************************************************************* //
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Virtual base class for displacement motion solver.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Calculate the matrix for the laplacian of the field.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
conserve primitiveFieldRef()+
Macros for easy insertion into run-time selection tables.
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
Base class for fvMesh based motionSolvers.
virtual void updateMesh(const mapPolyMesh &)
Update topology.
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Abstract base class for cell-centre mesh motion diffusivity.
static const char nl
Definition: Ostream.H:265
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
displacementLaplacianFvMotionSolver(const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list...
Definition: pointZone.H:62
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
messageStream Info
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:53
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
bool found
Namespace for OpenFOAM.
motionDiffusivity & diffusivity()
Return reference to the diffusivity field.