fvMotionSolverEngineMesh.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "fvcMeshPhi.H"
31#include "surfaceInterpolate.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39}
40
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
45:
47 pistonLayers_("pistonLayers", dimLength, Zero),
48 motionSolver_
49 (
50 *this,
51 engineDB_.engineDict()
52 )
53{
54 engineDB_.engineDict().readIfPresent("pistonLayers", pistonLayers_);
55}
56
57
58// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59
61{}
62
63
64// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65
67{
68 scalar deltaZ = engineDB_.pistonDisplacement().value();
69 Info<< "deltaZ = " << deltaZ << endl;
70
71 // Position of the top of the static mesh layers above the piston
72 scalar pistonPlusLayers = pistonPosition_.value() + pistonLayers_.value();
73
74 scalar pistonSpeed = deltaZ/engineDB_.deltaTValue();
75
76 motionSolver_.pointMotionU().boundaryFieldRef()[pistonIndex_] ==
77 pistonSpeed;
78
79 {
80 scalarField linerPoints
81 (
82 boundary()[linerIndex_].patch().localPoints().component(vector::Z)
83 );
84
85 motionSolver_.pointMotionU().boundaryFieldRef()[linerIndex_] ==
86 pistonSpeed*pos0(deckHeight_.value() - linerPoints)
87 *(deckHeight_.value() - linerPoints)
88 /(deckHeight_.value() - pistonPlusLayers);
89 }
90
91 motionSolver_.solve();
92
93 if (engineDB_.foundObject<surfaceScalarField>("phi"))
94 {
96 engineDB_.lookupObjectRef<surfaceScalarField>("phi");
97
98 const volScalarField& rho =
99 engineDB_.lookupObject<volScalarField>("rho");
100
101 const volVectorField& U =
102 engineDB_.lookupObject<volVectorField>("U");
103
104 bool absolutePhi = false;
105 if (moving())
106 {
108 absolutePhi = true;
109 }
110
111 movePoints(motionSolver_.curPoints());
112
113 if (absolutePhi)
114 {
116 }
117 }
118 else
119 {
120 movePoints(motionSolver_.curPoints());
121 }
122
123
124 pistonPosition_.value() += deltaZ;
125
126 Info<< "clearance: " << deckHeight_.value() - pistonPosition_.value() << nl
127 << "Piston speed = " << pistonSpeed << " m/s" << endl;
128}
129
130
131// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Foam::engineMesh.
Definition: engineMesh.H:56
const engineTime & engineDB_
Definition: engineMesh.H:70
const IOdictionary & engineDict() const
Return the engine geometry dictionary.
Definition: engineTime.H:130
Foam::fvMotionSolverEngineMesh.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
faceListList boundary
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:36
Namespace for OpenFOAM.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53