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-2022 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 {
33  defineTypeNameAndDebug(displacementLinearMotionMotionSolver, 0);
34 
36  (
37  motionSolver,
38  displacementLinearMotionMotionSolver,
39  dictionary
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(vector(coeffDict().lookup("axis")))),
56  xFixed_(coeffDict().lookup<scalar>("xFixed")),
57  xMoving_(coeffDict().lookup<scalar>("xMoving")),
58  displacement_(Function1<scalar>::New("displacement", coeffDict()))
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
73 {
74  tmp<pointField> tcurPoints(new pointField(points0()));
75  pointField& curPoints = tcurPoints.ref();
76 
77  const scalar t = mesh().time().userTimeValue();
78 
79  const scalar displacement = displacement_->value(t);
80 
81  forAll(curPoints, i)
82  {
83  const scalar lambda =
84  (xFixed_ - (axis_ & curPoints[i]))/(xFixed_ - xMoving_);
85 
86  if (lambda > 1)
87  {
88  curPoints[i] += axis_*displacement;
89  }
90  else if (lambda > 0)
91  {
92  curPoints[i] += axis_*lambda*displacement;
93  }
94  }
95 
96  return tcurPoints;
97 }
98 
99 
100 // ************************************************************************* //
Run-time selectable general function of one variable.
Definition: Function1.H:52
#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:156
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedScalar lambda(viscosity->lookup("lambda"))
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:845
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
stressControl lookup("compactNormalStress") >> compactNormalStress
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:413
A class for handling words, derived from string.
Definition: word.H:59
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
Virtual base class for displacement motion solvers.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
A class for managing temporary objects.
Definition: PtrList.H:53
displacementLinearMotionMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
Namespace for OpenFOAM.