Loading...
Searching...
No Matches
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{
37 defineTypeNameAndDebug(layeredEngineMesh, 0);
38 addToRunTimeSelectionTable(engineMesh, layeredEngineMesh, IOobject);
39}
40
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
44Foam::layeredEngineMesh::layeredEngineMesh(const IOobject& io)
45:
47 pistonLayers_("pistonLayers", dimLength, Zero)
49 engineDB_.engineDict().readIfPresent("pistonLayers", pistonLayers_);
50}
51
52
53// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
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
89 auto* phiPtr = engineDB_.getObjectPtr<surfaceScalarField>("phi");
90
91 if (phiPtr)
92 {
93 auto& phi = *phiPtr;
94
95 const auto& rho = engineDB_.lookupObject<volScalarField>("rho");
96 const auto& U = engineDB_.lookupObject<volVectorField>("U");
97
98 const bool absolutePhi = moving();
99 if (absolutePhi)
100 {
101 // cf. fvc::makeAbsolute
103 }
104
105 movePoints(newPoints);
106
107 if (absolutePhi)
108 {
109 // cf. fvc::makeRelative
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.
propsDict readIfPresent("fields", acceptFields)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Foam::engineMesh.
Definition engineMesh.H:52
dimensionedScalar deckHeight_
Definition engineMesh.H:76
const engineTime & engineDB_
Definition engineMesh.H:70
dimensionedScalar pistonPosition_
Definition engineMesh.H:77
const surfaceScalarField & phi() const
Return cell face motion fluxes.
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition fvMesh.C:884
Foam::layeredEngineMesh.
bool moving() const noexcept
Is mesh moving.
Definition polyMesh.H:732
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
volScalarField & p
const auto & io
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
const pointField & points
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:29
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
vectorField pointField
pointField is a vectorField.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299