displacementLinearMotionMotionSolver.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) 2018-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 
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
34 
36  (
40  );
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
48 (
49  const word& name,
50  const polyMesh& mesh,
51  const dictionary& dict
52 )
53 :
54  points0MotionSolver(name, mesh, dict, typeName),
55  axis_(normalised(coeffDict().lookup<vector>("axis", dimless))),
56  xFixed_(coeffDict().lookup<scalar>("xFixed", dimLength)),
57  xMoving_(coeffDict().lookup<scalar>("xMoving", dimLength)),
58  displacement_
59  (
60  Function1<scalar>::New
61  (
62  "displacement",
63  mesh.time().userUnits(),
64  dimLength,
65  coeffDict()
66  )
67  )
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
72 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
82 {
83  tmp<pointField> tcurPoints(new pointField(points0()));
84  pointField& curPoints = tcurPoints.ref();
85 
86  const scalar displacement = displacement_->value(mesh().time().value());
87 
88  forAll(curPoints, i)
89  {
90  const scalar lambda =
91  (xFixed_ - (axis_ & curPoints[i]))/(xFixed_ - xMoving_);
92 
93  if (lambda > 1)
94  {
95  curPoints[i] += axis_*displacement;
96  }
97  else if (lambda > 0)
98  {
99  curPoints[i] += axis_*lambda*displacement;
100  }
101  }
102 
103  return tcurPoints;
104 }
105 
106 
107 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Run-time selectable general function of one variable.
Definition: Function1.H:125
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh motion solver simple linear expansion and contraction of a mesh region defined by a motion axis ...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
displacementLinearMotionMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
Virtual base class for displacement motion solvers.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
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:181
A class for handling words, derived from string.
Definition: word.H:62
dimensionedScalar lambda(viscosity->lookup("lambda"))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
const dimensionSet dimLength
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:510
defineTypeNameAndDebug(combustionModel, 0)
dictionary dict