adding pressure gradient volumetric source (kinematic only)

This commit is contained in:
andy
2008-09-10 10:31:01 +01:00
parent 44f647d5ad
commit e60638fd4f
3 changed files with 305 additions and 0 deletions

View File

@ -0,0 +1,170 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2008 OpenCFD Ltd.
\\/ 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
\*---------------------------------------------------------------------------*/
#include "pressureGradientExplicitSource.H"
#include "volFields.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::pressureGradientExplicitSource::pressureGradientExplicitSource
(
const word& sourceName,
const fvMesh& mesh,
volVectorField& U
)
:
dict_
(
IOobject
(
sourceName + "Properties",
mesh.time().constant(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
)
),
mesh_(mesh),
U_(U),
Ubar_(dict_.lookup("Ubar")),
gradPini_(readScalar(dict_.lookup("gradPini"))),
gradP_(gradPini_),
flowDir_(Ubar_/mag(Ubar_)),
cellSource_(dict_.lookup("cellSource")),
cellSelector_
(
topoSetSource::New
(
cellSource_,
mesh,
dict_.subDict(cellSource_ + "Coeffs")
)
),
selectedCellSet_
(
mesh,
"pressureGradientExplicitSourceCellSet",
mesh.nCells()/10 + 1 // Reasonable size estimate.
)
{
// Create the cell set
cellSelector_->applyToSet
(
topoSetSource::NEW,
selectedCellSet_
);
// Give some feedback
Info<< "pressureGradientExplicitSource(" << sourceName << ")" << nl
<< "Selected " << returnReduce(selectedCellSet_.size(), sumOp<label>())
<< " cells." << endl;
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp<Foam::DimensionedField<Foam::vector, Foam::volMesh> >
Foam::pressureGradientExplicitSource::Su() const
{
tmp<DimensionedField<vector, volMesh> > tSource
(
new DimensionedField<vector, volMesh>
(
IOobject
(
"pressureGradientExplicitSource",
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh_,
dimensionedVector("zero", dimVelocity/dimTime, vector::zero)
)
);
DimensionedField<vector, volMesh>& sourceField = tSource();
forAllConstIter(cellSet, selectedCellSet_, iter)
{
label cellI = iter.key();
sourceField[cellI] = flowDir_*gradP_;
}
return tSource;
}
void Foam::pressureGradientExplicitSource::update()
{
const volScalarField& rUA =
mesh_.lookupObject<volScalarField>("(1|A(" + U_.name() + "))");
// Integrate flow variables over cell set
scalar volTot = 0.0;
scalar magUbarAve = 0.0;
scalar rUAave = 0.0;
forAllConstIter(cellSet, selectedCellSet_, iter)
{
label cellI = iter.key();
scalar volCell = mesh_.V()[cellI];
volTot += volCell;
magUbarAve += (flowDir_ & U_[cellI])*volCell;
rUAave += rUA[cellI]*volCell;
}
// Collect across all processors
reduce(volTot, sumOp<scalar>());
reduce(magUbarAve, sumOp<scalar>());
reduce(rUAave, sumOp<scalar>());
// Volume averages
magUbarAve /= volTot;
rUAave /= volTot;
// Calculate the pressure gradient increment needed to adjust the average
// flow-rate to the desired value
scalar gradPplus = (mag(Ubar_) - magUbarAve)/rUAave;
// Apply correction to velocity field
forAllConstIter(cellSet, selectedCellSet_, iter)
{
label cellI = iter.key();
U_[cellI] += flowDir_*rUA[cellI]*gradPplus;
}
// Update pressure gradient
gradP_ += gradPplus;
Info<< "Uncorrected Ubar = " << magUbarAve << tab
<< "Pressure gradient = " << gradP_ << endl;
}
// ************************************************************************* //

View File

@ -0,0 +1,134 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2008 OpenCFD Ltd.
\\/ 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
Foam::pressureGradientExplicitSource
Description
Creates a cell set pressure gradient source
Note: Currently only handles kinematic pressure
SourceFiles
pressureGradientExplicitSource.C
\*---------------------------------------------------------------------------*/
#ifndef pressureGradientExplicitSource_H
#define pressureGradientExplicitSource_H
#include "autoPtr.H"
#include "topoSetSource.H"
#include "cellSet.H"
#include "fvMesh.H"
#include "volFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class pressureGradientExplicitSource Declaration
\*---------------------------------------------------------------------------*/
class pressureGradientExplicitSource
{
// Private data
//- Properties dictionary
IOdictionary dict_;
//- Reference to the mesh
const fvMesh& mesh_;
//- Reference to the velocity field
volVectorField& U_;
//- Average velocity
vector Ubar_;
//- Initial pressure gradient
scalar gradPini_;
//- Pressure gradient
scalar gradP_;
//- Flow direction
vector flowDir_;
//- Name of cell source
word cellSource_;
//- The method by which the cells will be selecetd
autoPtr<topoSetSource> cellSelector_;
//- The set of selected cells
cellSet selectedCellSet_;
// Private Member Functions
//- Disallow default bitwise copy construct
pressureGradientExplicitSource(const pressureGradientExplicitSource&);
//- Disallow default bitwise assignment
void operator=(const pressureGradientExplicitSource&);
public:
// Constructors
//- Construct from explicit source name and mesh
pressureGradientExplicitSource
(
const word& sourceName,
const fvMesh& mesh,
volVectorField& U
);
// Member Functions
// Access
//- Return a tmp field of the source
tmp<DimensionedField<vector, volMesh> > Su() const;
//- Correct driving force for a constant mass flow rate
void update();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //