mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: Updated fvm field sources - added correct() function and updated grad(p) source
This commit is contained in:
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -69,53 +69,6 @@ void Foam::pressureGradientExplicitSource::writeGradP() const
|
||||
}
|
||||
|
||||
|
||||
void Foam::pressureGradientExplicitSource::update(fvMatrix<vector>& eqn)
|
||||
{
|
||||
volVectorField& U = const_cast<volVectorField&>(eqn.psi());
|
||||
|
||||
const volScalarField& rAU =
|
||||
mesh_.lookupObject<volScalarField>("(1|A(" + U.name() + "))");
|
||||
|
||||
// Integrate flow variables over cell set
|
||||
scalar magUbarAve = 0.0;
|
||||
scalar rAUave = 0.0;
|
||||
const scalarField& cv = mesh_.V();
|
||||
forAll(cells_, i)
|
||||
{
|
||||
label cellI = cells_[i];
|
||||
scalar volCell = cv[cellI];
|
||||
magUbarAve += (flowDir_ & U[cellI])*volCell;
|
||||
rAUave += rAU[cellI]*volCell;
|
||||
}
|
||||
|
||||
// Collect across all processors
|
||||
reduce(magUbarAve, sumOp<scalar>());
|
||||
|
||||
// Volume averages
|
||||
magUbarAve /= V_;
|
||||
rAUave /= V_;
|
||||
|
||||
// Calculate the pressure gradient increment needed to adjust the average
|
||||
// flow-rate to the desired value
|
||||
scalar gradPplus = (mag(Ubar_) - magUbarAve)/rAUave;
|
||||
|
||||
// Apply correction to velocity field
|
||||
forAll(cells_, i)
|
||||
{
|
||||
label cellI = cells_[i];
|
||||
U[cellI] += flowDir_*rAU[cellI]*gradPplus;
|
||||
}
|
||||
|
||||
// Update pressure gradient
|
||||
gradP_.value() += gradPplus;
|
||||
|
||||
Info<< "Uncorrected Ubar = " << magUbarAve << tab
|
||||
<< "Pressure gradient = " << gradP_.value() << endl;
|
||||
|
||||
writeGradP();
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::pressureGradientExplicitSource::pressureGradientExplicitSource
|
||||
@ -123,7 +76,7 @@ Foam::pressureGradientExplicitSource::pressureGradientExplicitSource
|
||||
const word& sourceName,
|
||||
const word& modelType,
|
||||
const dictionary& dict,
|
||||
const fvMesh& mesh
|
||||
const fvMesh& mesh
|
||||
)
|
||||
:
|
||||
basicSource(sourceName, modelType, dict, mesh),
|
||||
@ -133,6 +86,23 @@ Foam::pressureGradientExplicitSource::pressureGradientExplicitSource
|
||||
flowDir_(Ubar_/mag(Ubar_))
|
||||
{
|
||||
coeffs_.lookup("fieldNames") >> fieldNames_;
|
||||
|
||||
if (fieldNames_.size() != 1)
|
||||
{
|
||||
FatalErrorIn
|
||||
(
|
||||
"Foam::pressureGradientExplicitSource::"
|
||||
"pressureGradientExplicitSource"
|
||||
"("
|
||||
"onst word&, "
|
||||
"const word&, "
|
||||
"const dictionary&, "
|
||||
"const fvMesh&"
|
||||
")"
|
||||
) << "Source can only be applied to a single field. Current "
|
||||
<< "settings are:" << fieldNames_ << exit(FatalError);
|
||||
}
|
||||
|
||||
applied_.setSize(fieldNames_.size(), false);
|
||||
|
||||
// Read the initial pressure gradient from file if it exists
|
||||
@ -160,8 +130,6 @@ void Foam::pressureGradientExplicitSource::addSup
|
||||
const label fieldI
|
||||
)
|
||||
{
|
||||
update(eqn);
|
||||
|
||||
DimensionedField<vector, volMesh> Su
|
||||
(
|
||||
IOobject
|
||||
@ -178,7 +146,53 @@ void Foam::pressureGradientExplicitSource::addSup
|
||||
|
||||
UIndirectList<vector>(Su, cells_) = flowDir_*gradP_.value();
|
||||
|
||||
eqn -= Su;
|
||||
eqn += Su;
|
||||
}
|
||||
|
||||
|
||||
void Foam::pressureGradientExplicitSource::correct(volVectorField& U)
|
||||
{
|
||||
const volScalarField& rAU =
|
||||
mesh_.lookupObject<volScalarField>("(1|A(" + U.name() + "))");
|
||||
|
||||
// Integrate flow variables over cell set
|
||||
scalar magUbarAve = 0.0;
|
||||
scalar rAUave = 0.0;
|
||||
const scalarField& cv = mesh_.V();
|
||||
forAll(cells_, i)
|
||||
{
|
||||
label cellI = cells_[i];
|
||||
scalar volCell = cv[cellI];
|
||||
magUbarAve += (flowDir_ & U[cellI])*volCell;
|
||||
rAUave += rAU[cellI]*volCell;
|
||||
}
|
||||
|
||||
// Collect across all processors
|
||||
reduce(magUbarAve, sumOp<scalar>());
|
||||
reduce(rAUave, sumOp<scalar>());
|
||||
|
||||
// Volume averages
|
||||
magUbarAve /= V_;
|
||||
rAUave /= V_;
|
||||
|
||||
// Calculate the pressure gradient increment needed to adjust the average
|
||||
// flow-rate to the desired value
|
||||
scalar gradPplus = (mag(Ubar_) - magUbarAve)/rAUave;
|
||||
|
||||
// Apply correction to velocity field
|
||||
forAll(cells_, i)
|
||||
{
|
||||
label cellI = cells_[i];
|
||||
U[cellI] += flowDir_*rAU[cellI]*gradPplus;
|
||||
}
|
||||
|
||||
// Update pressure gradient
|
||||
gradP_.value() += gradPplus;
|
||||
|
||||
Info<< "Pressure gradient source: uncorrected Ubar = " << magUbarAve
|
||||
<< ", pressure gradient = " << gradP_.value() << endl;
|
||||
|
||||
writeGradP();
|
||||
}
|
||||
|
||||
|
||||
|
||||
Reference in New Issue
Block a user