layeredEngineMesh.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 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
28#include "layeredEngineMesh.H"
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{
49 engineDB_.engineDict().readIfPresent("pistonLayers", pistonLayers_);
50}
51
52
53// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
54
56{}
57
58
59// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60
62{
63 scalar deltaZ = engineDB_.pistonDisplacement().value();
64 Info<< "deltaZ = " << deltaZ << endl;
65
66 // Position of the top of the static mesh layers above the piston
67 scalar pistonPlusLayers = pistonPosition_.value() + pistonLayers_.value();
68
69 pointField newPoints(points());
70
71 forAll(newPoints, pointi)
72 {
73 point& p = newPoints[pointi];
74
75 if (p.z() < pistonPlusLayers) // In piston bowl
76 {
77 p.z() += deltaZ;
78 }
79 else if (p.z() < deckHeight_.value()) // In liner region
80 {
81 p.z() +=
82 deltaZ
83 *(deckHeight_.value() - p.z())
84 /(deckHeight_.value() - pistonPlusLayers);
85 }
86 }
87
88 if (engineDB_.foundObject<surfaceScalarField>("phi"))
89 {
91 engineDB_.lookupObjectRef<surfaceScalarField>("phi");
92
93 const volScalarField& rho =
94 engineDB_.lookupObject<volScalarField>("rho");
95
96 const volVectorField& U =
97 engineDB_.lookupObject<volVectorField>("U");
98
99 bool absolutePhi = false;
100 if (moving())
101 {
103 absolutePhi = true;
104 }
105
106 movePoints(newPoints);
107
108 if (absolutePhi)
109 {
111 }
112 }
113 else
114 {
115 movePoints(newPoints);
116 }
117
118 pistonPosition_.value() += deltaZ;
119 scalar pistonSpeed = deltaZ/engineDB_.deltaTValue();
120
121 Info<< "clearance: " << deckHeight_.value() - pistonPosition_.value() << nl
122 << "Piston speed = " << pistonSpeed << " m/s" << endl;
123}
124
125
126// ************************************************************************* //
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::layeredEngineMesh.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
volScalarField & p
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
const pointField & points
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.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
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
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333