rhoCentralFoam: Ensure fixed value boundary conditions are preserved

Resolves bug-report http://www.openfoam.org/mantisbt/view.php?id=1748
This commit is contained in:
Henry Weller
2015-06-17 08:46:46 +01:00
parent 280f0fbf31
commit bc493180d9
5 changed files with 49 additions and 26 deletions

View File

@ -32,7 +32,6 @@ volVectorField U
mesh
);
#include "rhoBoundaryTypes.H"
volScalarField rho
(
IOobject
@ -44,7 +43,7 @@ volScalarField rho
IOobject::AUTO_WRITE
),
thermo.rho(),
rhoBoundaryTypes
derivedPatchFieldTypes(p)
);
volVectorField rhoU
@ -57,7 +56,8 @@ volVectorField rhoU
IOobject::NO_READ,
IOobject::NO_WRITE
),
rho*U
rho*U,
derivedPatchFieldTypes(U)
);
volScalarField rhoE
@ -70,7 +70,8 @@ volScalarField rhoE
IOobject::NO_READ,
IOobject::NO_WRITE
),
rho*(e + 0.5*magSqr(U))
rho*(e + 0.5*magSqr(U)),
derivedPatchFieldTypes(T)
);
surfaceScalarField pos
@ -98,7 +99,18 @@ surfaceScalarField neg
);
surfaceScalarField phi("phi", mesh.Sf() & fvc::interpolate(rhoU));
surfaceScalarField phi
(
IOobject
(
"phi",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh.Sf() & fvc::interpolate(rhoU)
);
Info<< "Creating turbulence model\n" << endl;
autoPtr<compressible::turbulenceModel> turbulence

View File

@ -31,8 +31,9 @@ tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > interpolate
{
if
(
!sf.boundaryField()[patchi].coupled()
!sf.boundaryField()[patchi].coupled()
&& sf.boundaryField()[patchi].size()
&& !vf.boundaryField()[patchi].fixesValue()
&& dir.boundaryField()[patchi][0] > 0
)
{
@ -44,4 +45,28 @@ tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > interpolate
return tsf;
}
template<class Type>
wordList derivedPatchFieldTypes
(
const GeometricField<Type, fvPatchField, volMesh>& vf
)
{
wordList phiTypes
(
vf.boundaryField().size(),
calculatedFvPatchField<Type>::typeName
);
forAll(vf.boundaryField(), patchi)
{
if (vf.boundaryField()[patchi].fixesValue())
{
phiTypes[patchi] = fixedValueFvPatchField<Type>::typeName;
}
}
return phiTypes;
}
}

View File

@ -1,14 +0,0 @@
const volScalarField::GeometricBoundaryField& pbf = p.boundaryField();
wordList rhoBoundaryTypes = pbf.types();
forAll(rhoBoundaryTypes, patchi)
{
if (rhoBoundaryTypes[patchi] == "waveTransmissive")
{
rhoBoundaryTypes[patchi] = zeroGradientFvPatchScalarField::typeName;
}
else if (pbf[patchi].fixesValue())
{
rhoBoundaryTypes[patchi] = fixedRhoFvPatchScalarField::typeName;
}
}

View File

@ -187,7 +187,7 @@ int main(int argc, char *argv[])
rhoU.dimensionedInternalField()
/rho.dimensionedInternalField();
U.correctBoundaryConditions();
rhoU.boundaryField() = rho.boundaryField()*U.boundaryField();
rhoU.boundaryField() == rho.boundaryField()*U.boundaryField();
if (!inviscid)
{
@ -221,7 +221,7 @@ int main(int argc, char *argv[])
e = rhoE/rho - 0.5*magSqr(U);
e.correctBoundaryConditions();
thermo.correct();
rhoE.boundaryField() =
rhoE.boundaryField() ==
rho.boundaryField()*
(
e.boundaryField() + 0.5*magSqr(U.boundaryField())
@ -242,7 +242,7 @@ int main(int argc, char *argv[])
rho.dimensionedInternalField()
/psi.dimensionedInternalField();
p.correctBoundaryConditions();
rho.boundaryField() = psi.boundaryField()*p.boundaryField();
rho.boundaryField() == psi.boundaryField()*p.boundaryField();
turbulence->correct();

View File

@ -169,7 +169,7 @@ int main(int argc, char *argv[])
rhoU.dimensionedInternalField()
/rho.dimensionedInternalField();
U.correctBoundaryConditions();
rhoU.boundaryField() = rho.boundaryField()*U.boundaryField();
rhoU.boundaryField() == rho.boundaryField()*U.boundaryField();
if (!inviscid)
{
@ -203,7 +203,7 @@ int main(int argc, char *argv[])
e = rhoE/rho - 0.5*magSqr(U);
e.correctBoundaryConditions();
thermo.correct();
rhoE.boundaryField() =
rhoE.boundaryField() ==
rho.boundaryField()*
(
e.boundaryField() + 0.5*magSqr(U.boundaryField())
@ -224,7 +224,7 @@ int main(int argc, char *argv[])
rho.dimensionedInternalField()
/psi.dimensionedInternalField();
p.correctBoundaryConditions();
rho.boundaryField() = psi.boundaryField()*p.boundaryField();
rho.boundaryField() == psi.boundaryField()*p.boundaryField();
turbulence->correct();