Turbomachinery addition. Oliver Borm and Fethi Tekin
This commit is contained in:
parent
db148bed71
commit
cd7a8ad234
9 changed files with 640 additions and 80 deletions
|
@ -27,9 +27,10 @@ Class
|
|||
|
||||
Author
|
||||
Hrvoje Jasak, Wikki Ltd. All rights reserved
|
||||
Fethi Tekin, All rights reserved,
|
||||
|
||||
Description
|
||||
Mass-conservative face interpolation: typedef for polyPatch to
|
||||
Mass-conservative face interpolation: typedef for stand-alone patch to
|
||||
stand-alone patch
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -49,7 +50,7 @@ namespace Foam
|
|||
{
|
||||
typedef GGIInterpolation
|
||||
<
|
||||
PrimitivePatch<face, SubList, const pointField&>,
|
||||
PrimitivePatch<face, List, pointField>,
|
||||
PrimitivePatch<face, List, pointField>
|
||||
> overlapGgiInterpolation;
|
||||
}
|
||||
|
|
|
@ -24,6 +24,7 @@ License
|
|||
|
||||
Author
|
||||
Hrvoje Jasak, Wikki Ltd. All rights reserved.
|
||||
Fethi Tekin, All rights reserved.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
|
@ -60,7 +61,8 @@ Foam::overlapGgiPolyPatch::overlapGgiPolyPatch
|
|||
shadowName_(word::null),
|
||||
shadowIndex_(-1),
|
||||
rotationAxis_(vector(0.0, 0.0, 1.0)),
|
||||
angle_(0),
|
||||
nCopies_(0),
|
||||
expandedMasterPtr_(NULL),
|
||||
expandedSlavePtr_(NULL),
|
||||
patchToPatchPtr_(NULL),
|
||||
reconFaceCellCentresPtr_(NULL)
|
||||
|
@ -77,14 +79,15 @@ Foam::overlapGgiPolyPatch::overlapGgiPolyPatch
|
|||
const polyBoundaryMesh& bm,
|
||||
const word& shadowName,
|
||||
const vector& axis,
|
||||
const scalar angle
|
||||
const scalar nCopies
|
||||
)
|
||||
:
|
||||
coupledPolyPatch(name, size, start, index, bm),
|
||||
shadowName_(shadowName),
|
||||
shadowIndex_(-1),
|
||||
rotationAxis_(axis),
|
||||
angle_(angle),
|
||||
nCopies_(nCopies),
|
||||
expandedMasterPtr_(NULL),
|
||||
expandedSlavePtr_(NULL),
|
||||
patchToPatchPtr_(NULL),
|
||||
reconFaceCellCentresPtr_(NULL)
|
||||
|
@ -104,7 +107,8 @@ Foam::overlapGgiPolyPatch::overlapGgiPolyPatch
|
|||
shadowName_(dict.lookup("shadowPatch")),
|
||||
shadowIndex_(-1),
|
||||
rotationAxis_(dict.lookup("rotationAxis")),
|
||||
angle_(readScalar(dict.lookup("angle"))),
|
||||
nCopies_(readScalar(dict.lookup("nCopies"))),
|
||||
expandedMasterPtr_(NULL),
|
||||
expandedSlavePtr_(NULL),
|
||||
patchToPatchPtr_(NULL),
|
||||
reconFaceCellCentresPtr_(NULL)
|
||||
|
@ -122,7 +126,8 @@ Foam::overlapGgiPolyPatch::overlapGgiPolyPatch
|
|||
shadowName_(pp.shadowName_),
|
||||
shadowIndex_(-1),
|
||||
rotationAxis_(pp.rotationAxis_),
|
||||
angle_(pp.angle_),
|
||||
nCopies_(pp.nCopies_),
|
||||
expandedMasterPtr_(NULL),
|
||||
expandedSlavePtr_(NULL),
|
||||
patchToPatchPtr_(NULL),
|
||||
reconFaceCellCentresPtr_(NULL)
|
||||
|
@ -143,7 +148,8 @@ Foam::overlapGgiPolyPatch::overlapGgiPolyPatch
|
|||
shadowName_(pp.shadowName_),
|
||||
shadowIndex_(-1),
|
||||
rotationAxis_(pp.rotationAxis_),
|
||||
angle_(pp.angle_),
|
||||
nCopies_(pp.nCopies_),
|
||||
expandedMasterPtr_(NULL),
|
||||
expandedSlavePtr_(NULL),
|
||||
patchToPatchPtr_(NULL),
|
||||
reconFaceCellCentresPtr_(NULL)
|
||||
|
@ -199,7 +205,6 @@ Foam::label Foam::overlapGgiPolyPatch::shadowIndex() const
|
|||
return shadowIndex_;
|
||||
}
|
||||
|
||||
|
||||
const Foam::overlapGgiPolyPatch&
|
||||
Foam::overlapGgiPolyPatch::shadow() const
|
||||
{
|
||||
|
@ -209,34 +214,21 @@ Foam::overlapGgiPolyPatch::shadow() const
|
|||
|
||||
Foam::label Foam::overlapGgiPolyPatch::nCopies() const
|
||||
{
|
||||
// Calculate number of copies to be made for the expanded slave
|
||||
// to completely cover the master
|
||||
if (!master())
|
||||
{
|
||||
FatalErrorIn("label overlapGgiPolyPatch::nCopies() const")
|
||||
<< "nCopies requested for a slave. Error in master-slave logic"
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
label ncp = 0;
|
||||
scalar remainder = angle_;
|
||||
|
||||
const scalar slaveAngle = shadow().angle();
|
||||
|
||||
while (remainder > SMALL)
|
||||
{
|
||||
remainder -= slaveAngle;
|
||||
ncp++;
|
||||
}
|
||||
|
||||
return ncp;
|
||||
// Read the number of copies to be made from the dictionary for the
|
||||
// expanded slave and expanded master to cover 360 degrees
|
||||
return nCopies_;
|
||||
}
|
||||
|
||||
|
||||
bool Foam::overlapGgiPolyPatch::master() const
|
||||
{
|
||||
// The first overlapggi interface is master,second one is slave
|
||||
if (angle() == shadow().angle())
|
||||
{
|
||||
return start() < shadow().start() ;
|
||||
}
|
||||
|
||||
// Master is the one with the larger angle
|
||||
return angle() >= shadow().angle();
|
||||
return angle() > shadow().angle();
|
||||
}
|
||||
|
||||
|
||||
|
@ -246,7 +238,9 @@ void Foam::overlapGgiPolyPatch::write(Ostream& os) const
|
|||
polyPatch::write(os);
|
||||
os.writeKeyword("rotationAxis") << rotationAxis_
|
||||
<< token::END_STATEMENT << nl;
|
||||
os.writeKeyword("angle") << angle_
|
||||
os.writeKeyword("nCopies") << nCopies_
|
||||
<< token::END_STATEMENT << nl;
|
||||
os.writeKeyword("shadowPatch") << shadowName_
|
||||
<< token::END_STATEMENT << nl;
|
||||
}
|
||||
|
||||
|
|
|
@ -26,17 +26,17 @@ Class
|
|||
overlapGgiPolyPatch
|
||||
|
||||
Description
|
||||
Partial overlap generalised grid interface (GGI) patch. Master side
|
||||
remains unchanged and the slave side will be copied as needed
|
||||
to pave the master patch surface.
|
||||
Partial overlap generalised grid interface (GGI) patch. Master and slave
|
||||
sides are copied as much as the given number to complete the 360 degree
|
||||
cicumferential surface.
|
||||
|
||||
This implies that the master patch has got a larger angular pitch than
|
||||
the slave and that master and slave are aligned at one edge.
|
||||
Master and slave will specify the pitch, based on which the expansion
|
||||
of the master side will be performed and used in interpolation.
|
||||
The data interpolation between master and slave patches do not depend on
|
||||
relative position of them, because of the full circumferential expansion
|
||||
for both sides.
|
||||
|
||||
Author
|
||||
Hrvoje Jasak, Wikki Ltd. All rights reserved
|
||||
Fethi Tekin, All rights reserved.
|
||||
|
||||
SourceFiles
|
||||
overlapGgiPolyPatch.C
|
||||
|
@ -78,12 +78,20 @@ class overlapGgiPolyPatch
|
|||
//- Rotation axis
|
||||
const vector rotationAxis_;
|
||||
|
||||
//- Wedge angle
|
||||
const scalar angle_;
|
||||
|
||||
// Number of copies in order to complete 360 degrees
|
||||
const scalar nCopies_;
|
||||
|
||||
// Interpolation
|
||||
|
||||
//- Pointer to expanded master patch
|
||||
mutable standAlonePatch* expandedMasterPtr_;
|
||||
|
||||
//- Pointer to expanded faceCentres of Master
|
||||
mutable vectorField* faceCentresPtr_;
|
||||
|
||||
//- Pointer to expanded faceareas of Master
|
||||
mutable vectorField* faceAreasPtr_;
|
||||
|
||||
//- Pointer to expanded slave patch
|
||||
mutable standAlonePatch* expandedSlavePtr_;
|
||||
|
||||
|
@ -97,17 +105,23 @@ class overlapGgiPolyPatch
|
|||
|
||||
// Private member functions
|
||||
|
||||
//- Return reference to patch-to-patch interpolation
|
||||
const overlapGgiInterpolation& patchToPatch() const;
|
||||
|
||||
//- Calculate expanded master patch
|
||||
void calcExpandedMaster() const;
|
||||
|
||||
//- Calculate expanded slave patch
|
||||
void calcExpandedSlave() const;
|
||||
void calcExpandedSlave() const;
|
||||
|
||||
//- Return reference to expanded master patch
|
||||
const standAlonePatch& expandedMaster() const;
|
||||
|
||||
//- Return reference to expanded slave patch
|
||||
const standAlonePatch& expandedSlave() const;
|
||||
|
||||
//- Calculate interpolation
|
||||
void calcPatchToPatch() const;
|
||||
|
||||
//- Return reference to patch-to-patch interpolation
|
||||
const overlapGgiInterpolation& patchToPatch() const;
|
||||
void calcPatchToPatch() const;
|
||||
|
||||
//- Calculate reconstructed cell centres
|
||||
void calcReconFaceCellCentres() const;
|
||||
|
@ -118,10 +132,6 @@ class overlapGgiPolyPatch
|
|||
//- Check definition: angles and offsets
|
||||
void checkDefinition() const;
|
||||
|
||||
//- Expand slave face field to full coverage
|
||||
template<class Type>
|
||||
tmp<Field<Type> > expandSlaveData(const Field<Type>& spf) const;
|
||||
|
||||
//- Clear out
|
||||
void clearOut();
|
||||
|
||||
|
@ -177,7 +187,7 @@ public:
|
|||
const polyBoundaryMesh& bm,
|
||||
const word& shadowName,
|
||||
const vector& axis,
|
||||
const scalar angle
|
||||
const scalar n
|
||||
);
|
||||
|
||||
//- Construct from dictionary
|
||||
|
@ -263,9 +273,9 @@ public:
|
|||
}
|
||||
|
||||
//- Return wedge angle
|
||||
const scalar& angle() const
|
||||
scalar angle() const
|
||||
{
|
||||
return angle_;
|
||||
return 360.0/nCopies();
|
||||
}
|
||||
|
||||
//- Return number of slave copies
|
||||
|
@ -280,6 +290,14 @@ public:
|
|||
return !master();
|
||||
}
|
||||
|
||||
//- Expand master face field to full for 360 degrees coverage
|
||||
template<class Type>
|
||||
tmp<Field<Type> > expandMasterData(const Field<Type>& spf) const;
|
||||
|
||||
//- Expand slave face field to full for 360 degrees coverage
|
||||
template<class Type>
|
||||
tmp<Field<Type> > expandSlaveData(const Field<Type>& spf) const;
|
||||
|
||||
//- Interpolate face field: given field on a the shadow side,
|
||||
// create an interpolated field on this side
|
||||
template<class Type>
|
||||
|
|
|
@ -24,11 +24,14 @@ License
|
|||
|
||||
Author
|
||||
Hrvoje Jasak, Wikki Ltd. All rights reserved.
|
||||
Fethi Tekin, All rights reserved.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "overlapGgiPolyPatch.H"
|
||||
#include "polyMesh.H"
|
||||
#include "polyPatch.H"
|
||||
#include "primitiveMesh.H"
|
||||
#include "demandDrivenData.H"
|
||||
#include "polyPatchID.H"
|
||||
#include "OFstream.H"
|
||||
|
@ -36,6 +39,108 @@ Author
|
|||
|
||||
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||
|
||||
void Foam::overlapGgiPolyPatch::calcExpandedMaster() const
|
||||
{
|
||||
// Create expanded master patch interpolation
|
||||
if (expandedMasterPtr_)
|
||||
{
|
||||
FatalErrorIn("void overlapGgiPolyPatch::calcExpandedMaster() const")
|
||||
<< "Expanded master already calculated"
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
if (master())
|
||||
{
|
||||
// Create expanded master patch
|
||||
const label ncpm = nCopies();
|
||||
|
||||
Info << "Number of master copies: " << ncpm << endl;
|
||||
|
||||
// Create expanded master points and faces
|
||||
const polyPatch& master = boundaryMesh()[index()];
|
||||
const pointField& masterLocalPoints = master.localPoints();
|
||||
|
||||
pointField MasterExpandedPoints(ncpm*masterLocalPoints.size());
|
||||
|
||||
const scalar masterAngle = angle();
|
||||
|
||||
Info << "Master Angle is: " << masterAngle << endl;
|
||||
|
||||
// Transform points
|
||||
label nPoints_master = 0;
|
||||
|
||||
for (label copyI = 0; copyI < ncpm; copyI++)
|
||||
{
|
||||
// Calculate transform
|
||||
const tensor curRotation =
|
||||
RodriguesRotation(rotationAxis_, copyI*(masterAngle));
|
||||
|
||||
forAll (masterLocalPoints, pointI)
|
||||
{
|
||||
MasterExpandedPoints[nPoints_master] =
|
||||
transform(curRotation, masterLocalPoints[pointI]);
|
||||
nPoints_master++;
|
||||
}
|
||||
}
|
||||
|
||||
// Transform faces
|
||||
const faceList& masterLocalFaces = master.localFaces();
|
||||
faceList MasterExpandedFaces(ncpm*masterLocalFaces.size());
|
||||
|
||||
label nFacesMaster = 0;
|
||||
|
||||
for (label copyI = 0; copyI < ncpm; copyI++)
|
||||
{
|
||||
const label copyOffsetMaster = copyI*masterLocalPoints.size();
|
||||
|
||||
forAll (masterLocalFaces, faceI)
|
||||
{
|
||||
const face& curMasterFace = masterLocalFaces[faceI];
|
||||
|
||||
face& MastercurExpandedFace =
|
||||
MasterExpandedFaces[nFacesMaster];
|
||||
|
||||
// Copy face with offsets
|
||||
MastercurExpandedFace.setSize(curMasterFace.size());
|
||||
|
||||
forAll (curMasterFace, fpI)
|
||||
{
|
||||
MastercurExpandedFace[fpI] =
|
||||
curMasterFace[fpI] + copyOffsetMaster;
|
||||
}
|
||||
|
||||
nFacesMaster++;
|
||||
}
|
||||
}
|
||||
|
||||
expandedMasterPtr_ =
|
||||
new standAlonePatch(MasterExpandedFaces, MasterExpandedPoints);
|
||||
|
||||
if (debug > 1)
|
||||
{
|
||||
Info << "Writing expanded master patch as VTK" << endl;
|
||||
|
||||
const polyMesh& mesh = boundaryMesh().mesh();
|
||||
|
||||
fileName fvPath(mesh.time().path()/"VTK");
|
||||
mkDir(fvPath);
|
||||
|
||||
standAlonePatch::writeVTK
|
||||
(
|
||||
fvPath/fileName("expandedMaster" + name() + shadow().name()),
|
||||
MasterExpandedFaces,
|
||||
MasterExpandedPoints
|
||||
);
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
FatalErrorIn("void overlapGgiPolyPatch::calcExpandedMaster() const")
|
||||
<< "Attempting to create expanded master on a shadow"
|
||||
<< abort(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
void Foam::overlapGgiPolyPatch::calcExpandedSlave() const
|
||||
{
|
||||
// Create expanded slave patch interpolation
|
||||
|
@ -49,7 +154,7 @@ void Foam::overlapGgiPolyPatch::calcExpandedSlave() const
|
|||
if (master())
|
||||
{
|
||||
// Create expanded patch
|
||||
const label ncp = nCopies();
|
||||
const label ncp = shadow().nCopies();
|
||||
|
||||
Info << "Number of slave copies: " << ncp << endl;
|
||||
|
||||
|
@ -135,11 +240,18 @@ void Foam::overlapGgiPolyPatch::calcExpandedSlave() const
|
|||
}
|
||||
|
||||
|
||||
const Foam::standAlonePatch& Foam::overlapGgiPolyPatch::expandedMaster() const
|
||||
{
|
||||
if (!expandedMasterPtr_)
|
||||
{
|
||||
calcExpandedMaster();
|
||||
}
|
||||
|
||||
return *expandedMasterPtr_;
|
||||
}
|
||||
|
||||
const Foam::standAlonePatch& Foam::overlapGgiPolyPatch::expandedSlave() const
|
||||
{
|
||||
// Note: expanded slave only exists for the master side.
|
||||
// Slave patch will use master for interpolation.
|
||||
// HJ, 14/Jan/2009
|
||||
if (!expandedSlavePtr_)
|
||||
{
|
||||
calcExpandedSlave();
|
||||
|
@ -151,7 +263,8 @@ const Foam::standAlonePatch& Foam::overlapGgiPolyPatch::expandedSlave() const
|
|||
|
||||
void Foam::overlapGgiPolyPatch::calcPatchToPatch() const
|
||||
{
|
||||
// Create patch-to-patch interpolation
|
||||
// Create patch-to-patch interpolation between the expanded master
|
||||
// and slave patches
|
||||
if (patchToPatchPtr_)
|
||||
{
|
||||
FatalErrorIn("void overlapGgiPolyPatch::calcPatchToPatch() const")
|
||||
|
@ -164,7 +277,7 @@ void Foam::overlapGgiPolyPatch::calcPatchToPatch() const
|
|||
patchToPatchPtr_ =
|
||||
new overlapGgiInterpolation
|
||||
(
|
||||
*this,
|
||||
expandedMaster(),
|
||||
expandedSlave(),
|
||||
forwardT(),
|
||||
reverseT(),
|
||||
|
@ -214,8 +327,6 @@ void Foam::overlapGgiPolyPatch::calcReconFaceCellCentres() const
|
|||
// Create neighbouring face centres using interpolation
|
||||
if (master())
|
||||
{
|
||||
vectorField alpha(shadow().size(), vector(0, 0, 1));
|
||||
|
||||
const label shadowID = shadowIndex();
|
||||
|
||||
// Get the transformed and interpolated shadow face cell centers
|
||||
|
@ -265,6 +376,7 @@ void Foam::overlapGgiPolyPatch::checkDefinition() const
|
|||
|
||||
void Foam::overlapGgiPolyPatch::clearOut()
|
||||
{
|
||||
deleteDemandDrivenData(expandedMasterPtr_);
|
||||
deleteDemandDrivenData(expandedSlavePtr_);
|
||||
deleteDemandDrivenData(patchToPatchPtr_);
|
||||
deleteDemandDrivenData(reconFaceCellCentresPtr_);
|
||||
|
|
|
@ -24,6 +24,7 @@ License
|
|||
|
||||
Author
|
||||
Hrvoje Jasak, Wikki Ltd. All rights reserved.
|
||||
Fethi Tekin, All rights reserved.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
|
@ -36,20 +37,10 @@ template<class Type>
|
|||
Foam::tmp<Foam::Field<Type> >
|
||||
Foam::overlapGgiPolyPatch::expandSlaveData(const Field<Type>& spf) const
|
||||
{
|
||||
if (spf.size() != shadow().size())
|
||||
{
|
||||
FatalErrorIn
|
||||
(
|
||||
"tmp<Field<Type> > overlapGgiPolyPatch::interpolate"
|
||||
"(const Field<Type>& spf) const"
|
||||
) << " Incorrect field size for expansion. Field size: "
|
||||
<< spf.size() << " patch size: " << shadow().size()
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
const scalar slaveAngle = shadow().angle();
|
||||
|
||||
const label ncp = nCopies();
|
||||
const label ncp = shadow().nCopies();
|
||||
|
||||
tmp<Field<Type> > tef(new Field<Type>(ncp*spf.size()));
|
||||
Field<Type>& ef = tef();
|
||||
|
@ -72,6 +63,35 @@ Foam::overlapGgiPolyPatch::expandSlaveData(const Field<Type>& spf) const
|
|||
return tef;
|
||||
}
|
||||
|
||||
template<class Type>
|
||||
Foam::tmp<Foam::Field<Type> >
|
||||
Foam::overlapGgiPolyPatch::expandMasterData(const Field<Type>& spf) const
|
||||
{
|
||||
const scalar masterAngle = angle();
|
||||
|
||||
const label ncpm = nCopies();
|
||||
|
||||
tmp<Field<Type> > tef(new Field<Type>(ncpm*spf.size()));
|
||||
Field<Type>& ef = tef();
|
||||
|
||||
label nFaces = 0;
|
||||
|
||||
for (label copyI = 0; copyI < ncpm; copyI++)
|
||||
{
|
||||
// Calculate transform
|
||||
const tensor curRotation =
|
||||
RodriguesRotation(rotationAxis_, copyI*(masterAngle));
|
||||
|
||||
forAll (spf, faceI)
|
||||
{
|
||||
ef[nFaces] = transform(curRotation, spf[faceI]);
|
||||
nFaces++;
|
||||
}
|
||||
}
|
||||
|
||||
return tef;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
|
@ -82,13 +102,20 @@ Foam::overlapGgiPolyPatch::interpolate(const Field<Type>& pf) const
|
|||
if (master())
|
||||
{
|
||||
// Expand slave data
|
||||
tmp<Field<Type> > expand = expandSlaveData(pf);
|
||||
tmp<Field<Type> > expandslave = expandSlaveData(pf);
|
||||
|
||||
return patchToPatch().slaveToMaster(expand);
|
||||
tmp<Field<Type> > tresult= patchToPatch().slaveToMaster(expandslave);
|
||||
// Truncate to size
|
||||
tresult().setSize(size());
|
||||
|
||||
return tresult;
|
||||
}
|
||||
else
|
||||
{
|
||||
tmp<Field<Type> > tresult = patchToPatch().masterToSlave(pf);
|
||||
// Expand master data
|
||||
tmp<Field<Type> > expandmaster = expandMasterData(pf);
|
||||
|
||||
tmp<Field<Type> > tresult = patchToPatch().masterToSlave(expandmaster);
|
||||
|
||||
// Truncate to size
|
||||
tresult().setSize(size());
|
||||
|
@ -107,11 +134,19 @@ Foam::overlapGgiPolyPatch::interpolate(const tmp<Field<Type> >& tpf) const
|
|||
// Expand slave data
|
||||
tmp<Field<Type> > expand = expandSlaveData(tpf());
|
||||
|
||||
return patchToPatch().slaveToMaster(expand);
|
||||
tmp<Field<Type> > tresult = patchToPatch().slaveToMaster(expand);
|
||||
|
||||
// Truncate to size
|
||||
tresult().setSize(size());
|
||||
|
||||
return tresult;
|
||||
}
|
||||
else
|
||||
{
|
||||
tmp<Field<Type> > tresult = patchToPatch().masterToSlave(tpf);
|
||||
// Expand master data
|
||||
tmp<Field<Type> > expandmaster = expandMasterData(tpf());
|
||||
|
||||
tmp<Field<Type> > tresult = patchToPatch().masterToSlave(expandmaster);
|
||||
|
||||
// Truncate to size
|
||||
tresult().setSize(size());
|
||||
|
|
|
@ -18,6 +18,7 @@ $(solidBodyMotionFunctions)/SKA/SKA.C
|
|||
$(solidBodyMotionFunctions)/translation/translation.C
|
||||
|
||||
mixerGgiFvMesh/mixerGgiFvMesh.C
|
||||
turboFvMesh/turboFvMesh.C
|
||||
|
||||
tetMetrics/tetMetric.C
|
||||
tetMetrics/tetMetrics.C
|
||||
|
|
43
src/dynamicFvMesh/turboFvMesh/dynamicMeshDict
Normal file
43
src/dynamicFvMesh/turboFvMesh/dynamicMeshDict
Normal file
|
@ -0,0 +1,43 @@
|
|||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: 1.5.dev |
|
||||
| \\ / A nd | Web: http://www.openfoam.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
FoamFile
|
||||
{
|
||||
version 2.0;
|
||||
format ascii;
|
||||
class dictionary;
|
||||
object dynamicMeshDict;
|
||||
}
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
dynamicFvMeshLib "libtopoChangerFvMesh.so";
|
||||
dynamicFvMesh turboFvMesh;
|
||||
|
||||
turboFvMeshCoeffs
|
||||
{
|
||||
coordinateSystem
|
||||
{
|
||||
type cylindrical;
|
||||
origin (0 0 0);
|
||||
axis (0 0 1);
|
||||
direction (1 0 0);
|
||||
}
|
||||
|
||||
rpm
|
||||
{
|
||||
Rotor1 60;
|
||||
Rotor2 -30;
|
||||
Stator 0;
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
|
||||
|
225
src/dynamicFvMesh/turboFvMesh/turboFvMesh.C
Normal file
225
src/dynamicFvMesh/turboFvMesh/turboFvMesh.C
Normal file
|
@ -0,0 +1,225 @@
|
|||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright held by original author
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of OpenFOAM
|
||||
|
||||
OpenFOAM is free software; you can redistribute it and/or modify it
|
||||
under the terms of the GNU General Public License as published by the
|
||||
Free Software Foundation; either version 2 of the License, or (at your
|
||||
option) any later version.
|
||||
|
||||
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
||||
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
||||
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
||||
for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OpenFOAM; if not, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
Instructions:
|
||||
This tool is used to have multiple rotating regions around the same origin
|
||||
with different rpms.
|
||||
Creating the cellZones is not implemented in this tool.
|
||||
The steps to obtain the cellZones are :
|
||||
|
||||
1) use regionCellSets utility. With this command you can have different
|
||||
cellSets for each region.
|
||||
|
||||
2) Change the name of the regions eg. CellRegion0 to Rotor1 ,CellRegion1
|
||||
to Stator and vice versa.
|
||||
|
||||
3) run command " setsToZones -noFlipMap ". After this command the
|
||||
cellSets are transformed to cellZones.
|
||||
|
||||
4) in dynamicMeshDict rpm section should be added.
|
||||
|
||||
5) The case is ready to be run.
|
||||
|
||||
Implemented by Fethi Tekin, 24.03.2010
|
||||
|
||||
Author
|
||||
Hrvoje Jasak, Wikki Ltd. All rights reserved.
|
||||
Fethi Tekin, All rights reserved.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
#include "turboFvMesh.H"
|
||||
#include "Time.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
defineTypeNameAndDebug(turboFvMesh, 0);
|
||||
addToRunTimeSelectionTable(dynamicFvMesh, turboFvMesh, IOobject);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||
|
||||
void Foam::turboFvMesh::addZonesAndModifiers()
|
||||
{
|
||||
// Add zones and modifiers for motion action
|
||||
// No functionality is implemented
|
||||
|
||||
if (cellZones().size() > 0)
|
||||
{
|
||||
Info<< "void turboFvMesh::addZonesAndModifiers() : "
|
||||
<< "Zones and modifiers already present. Skipping."
|
||||
<< endl;
|
||||
|
||||
return;
|
||||
}
|
||||
else
|
||||
{
|
||||
FatalErrorIn("turboFvMesh")
|
||||
<< "Cell Regions have to be created"
|
||||
<< abort(FatalError);
|
||||
}
|
||||
}
|
||||
void Foam::turboFvMesh::calcMovingPoints() const
|
||||
{
|
||||
if (debug)
|
||||
{
|
||||
Info<< "void turboFvMesh::calcMovingMasks() const : "
|
||||
<< "Calculating point and cell masks"
|
||||
<< endl;
|
||||
}
|
||||
|
||||
if (movingPointsPtr_)
|
||||
{
|
||||
FatalErrorIn("void turboFvMesh::calcMovingMasks() const")
|
||||
<< "point mask already calculated"
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
// Retrieve the zone Names
|
||||
const wordList zoneNames = cellZones().names();
|
||||
|
||||
// Set the points
|
||||
movingPointsPtr_ = new vectorField(allPoints().size(),vector::zero);
|
||||
|
||||
vectorField& movingPoints = *movingPointsPtr_;
|
||||
|
||||
const cellList& c = cells();
|
||||
const faceList& f = allFaces();
|
||||
|
||||
forAll(zoneNames,zoneI)
|
||||
{
|
||||
Info<< "Moving Region Zone Name:" << zoneNames[zoneI]<< nl<<endl;
|
||||
|
||||
const labelList& cellAddr =
|
||||
cellZones()[cellZones().findZoneID(zoneNames[zoneI])];
|
||||
|
||||
rpm_ = readScalar(dict_.subDict("rpm").lookup(zoneNames[zoneI]));
|
||||
|
||||
Info<< "rpm:" << rpm_<<nl<<endl;
|
||||
|
||||
forAll (cellAddr, cellI)
|
||||
{
|
||||
const cell& curCell = c[cellAddr[cellI]];
|
||||
|
||||
forAll (curCell, faceI)
|
||||
{
|
||||
const face& curFace = f[curCell[faceI]];
|
||||
|
||||
forAll (curFace, pointI)
|
||||
{
|
||||
//The rotation data is saved within the cell data. For
|
||||
//non-rotating regions rpm is zero,so mesh movement is
|
||||
//also zero. The conversion of rotational speed
|
||||
movingPoints[curFace[pointI]] =
|
||||
vector(0,rpm_*360.0*time().deltaT().value()/60.0, 0);
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
// Construct from components
|
||||
Foam::turboFvMesh::turboFvMesh
|
||||
(
|
||||
const IOobject& io
|
||||
)
|
||||
:
|
||||
dynamicFvMesh(io),
|
||||
dict_
|
||||
(
|
||||
IOdictionary
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"dynamicMeshDict",
|
||||
time().constant(),
|
||||
*this,
|
||||
IOobject::MUST_READ,
|
||||
IOobject::NO_WRITE
|
||||
)
|
||||
).subDict(typeName + "Coeffs")
|
||||
),
|
||||
csPtr_
|
||||
(
|
||||
coordinateSystem::New
|
||||
(
|
||||
"coordinateSystem",
|
||||
dict_.subDict("coordinateSystem")
|
||||
)
|
||||
),
|
||||
|
||||
movingPointsPtr_(NULL)
|
||||
{
|
||||
addZonesAndModifiers();
|
||||
|
||||
Info<< "Turbomachine Mixer mesh:" << nl
|
||||
<< " origin: " << cs().origin() << nl
|
||||
<< " axis : " << cs().axis() << endl;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::turboFvMesh::~turboFvMesh()
|
||||
{
|
||||
deleteDemandDrivenData(movingPointsPtr_);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
// Return moving points mask
|
||||
const Foam::vectorField& Foam::turboFvMesh::movingPoints() const
|
||||
{
|
||||
if (!movingPointsPtr_)
|
||||
{
|
||||
calcMovingPoints();
|
||||
}
|
||||
|
||||
return *movingPointsPtr_;
|
||||
}
|
||||
|
||||
|
||||
bool Foam::turboFvMesh::update()
|
||||
{
|
||||
movePoints
|
||||
(
|
||||
csPtr_->globalPosition
|
||||
(
|
||||
csPtr_->localPosition(allPoints()) + movingPoints()
|
||||
)
|
||||
);
|
||||
|
||||
// The mesh is not morphing
|
||||
return false;
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
131
src/dynamicFvMesh/turboFvMesh/turboFvMesh.H
Normal file
131
src/dynamicFvMesh/turboFvMesh/turboFvMesh.H
Normal file
|
@ -0,0 +1,131 @@
|
|||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2006-7 H. Jasak All rights reserved
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of OpenFOAM.
|
||||
|
||||
OpenFOAM is free software; you can redistribute it and/or modify it
|
||||
under the terms of the GNU General Public License as published by the
|
||||
Free Software Foundation; either version 2 of the License, or (at your
|
||||
option) any later version.
|
||||
|
||||
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
||||
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
||||
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
||||
for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OpenFOAM; if not, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
Class
|
||||
turboFvMesh
|
||||
|
||||
Description
|
||||
Simple mixer mesh using an overlapGGI interface
|
||||
|
||||
Author
|
||||
Hrvoje Jasak, Wikki Ltd. All rights reserved.
|
||||
Fethi Tekin, All rights reserved.
|
||||
|
||||
SourceFiles
|
||||
turboFvMesh.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef turboFvMesh_H
|
||||
#define turboFvMesh_H
|
||||
|
||||
#include "dynamicFvMesh.H"
|
||||
#include "cylindricalCS.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class turboFvMesh Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class turboFvMesh
|
||||
:
|
||||
public dynamicFvMesh
|
||||
{
|
||||
// Private data
|
||||
|
||||
//- Motion dictionary
|
||||
dictionary dict_;
|
||||
|
||||
//- Coordinate system
|
||||
autoPtr<coordinateSystem> csPtr_;
|
||||
|
||||
// - Rotational speed in rotations per minute (rpm)
|
||||
mutable scalar rpm_;
|
||||
|
||||
//- Markup field for points.
|
||||
mutable vectorField* movingPointsPtr_;
|
||||
|
||||
|
||||
// Private Member Functions
|
||||
|
||||
//- Disallow default bitwise copy construct
|
||||
turboFvMesh(const turboFvMesh&);
|
||||
|
||||
//- Disallow default bitwise assignment
|
||||
void operator=(const turboFvMesh&);
|
||||
|
||||
|
||||
//- Add mixer zones
|
||||
void addZonesAndModifiers();
|
||||
|
||||
//- Calculate moving Points
|
||||
void calcMovingPoints() const;
|
||||
|
||||
//- Return moving points
|
||||
const vectorField& movingPoints() const;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("turboFvMesh");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from IOobject
|
||||
explicit turboFvMesh(const IOobject& io);
|
||||
|
||||
|
||||
// Destructor
|
||||
|
||||
virtual ~turboFvMesh();
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Return coordinate system
|
||||
const coordinateSystem& cs() const
|
||||
{
|
||||
return csPtr_();
|
||||
}
|
||||
|
||||
//- Update the mesh for both mesh motion
|
||||
virtual bool update();
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
Reference in a new issue