fvcConsistentDdtPhiCorr for ddtSchemes in dbns library

This commit is contained in:
Vuko Vukcevic 2017-01-23 08:32:32 +01:00
parent 42306f97d6
commit b2d3172e52
4 changed files with 86 additions and 0 deletions

View file

@ -624,6 +624,30 @@ EulerLocalDdtScheme<Type>::fvcDdtPhiCorr
}
template<class Type>
tmp<typename EulerLocalDdtScheme<Type>::fluxFieldType>
EulerLocalDdtScheme<Type>::fvcDdtConsistentPhiCorr
(
const GeometricField<Type, fvsPatchField, surfaceMesh>& faceU,
const GeometricField<Type, fvPatchField, volMesh>& U,
const surfaceScalarField& rAUf
)
{
const objectRegistry& registry = this->mesh();
// Get access to the scalar beta[i]
const scalarField& beta =
registry.lookupObject<scalarField>(deltaTName_);
const surfaceScalarField rDeltaTf = fvc::interpolate
(
1.0/(beta[0]*registry.lookupObject<volScalarField>(deltaTauName_))
);
return (mesh().Sf() & faceU.oldTime())*rAUf*rDeltaTf;
}
template<class Type>
tmp<surfaceScalarField> EulerLocalDdtScheme<Type>::meshPhi
(

View file

@ -165,6 +165,16 @@ public:
const fluxFieldType& phi
);
// Member functions for the new time consistent formulation
tmp<fluxFieldType> fvcDdtConsistentPhiCorr
(
const GeometricField<Type, fvsPatchField, surfaceMesh>& faceU,
const GeometricField<Type, fvPatchField, volMesh>& U,
const surfaceScalarField& rAUf
);
tmp<surfaceScalarField> meshPhi
(
const GeometricField<Type, fvPatchField, volMesh>&
@ -191,6 +201,15 @@ tmp<surfaceScalarField> EulerLocalDdtScheme<scalar>::fvcDdtPhiCorr
);
template<>
tmp<surfaceScalarField> EulerLocalDdtScheme<scalar>::fvcDdtConsistentPhiCorr
(
const surfaceScalarField& faceU,
const volScalarField& U,
const surfaceScalarField& rAUf
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace fv

View file

@ -840,6 +840,31 @@ backwardDualDdtScheme<Type>::fvcDdtPhiCorr
}
template<class Type>
tmp<typename backwardDualDdtScheme<Type>::fluxFieldType>
backwardDualDdtScheme<Type>::fvcDdtConsistentPhiCorr
(
const GeometricField<Type, fvsPatchField, surfaceMesh>& faceU,
const GeometricField<Type, fvPatchField, volMesh>& U,
const surfaceScalarField& rAUf
)
{
notImplemented
(
type()
+ "::fvcDdtConsistentPhiCorr"
+ "\n("
+ "\n const GeometricField<Type, fvsPatchField, surfaceMesh>& faceU,"
+ "\n const GeometricField<Type, fvPatchField, volMesh>& U"
+ "\n const surfaceScalarField rAUf"
+ "\n)"
);
// Dummy return
return tmp<fluxFieldType>(NULL);
}
template<class Type>
tmp<surfaceScalarField> backwardDualDdtScheme<Type>::meshPhi
(

View file

@ -198,6 +198,15 @@ public:
);
// Member functions for the new time consistent formulation
tmp<fluxFieldType> fvcDdtConsistentPhiCorr
(
const GeometricField<Type, fvsPatchField, surfaceMesh>& faceU,
const GeometricField<Type, fvPatchField, volMesh>& U,
const surfaceScalarField& rAUf
);
tmp<surfaceScalarField> meshPhi
(
const GeometricField<Type, fvPatchField, volMesh>&
@ -224,6 +233,15 @@ tmp<surfaceScalarField> backwardDualDdtScheme<scalar>::fvcDdtPhiCorr
);
template<>
tmp<surfaceScalarField> backwardDualDdtScheme<scalar>::fvcDdtConsistentPhiCorr
(
const surfaceScalarField& faceU,
const volScalarField& U,
const surfaceScalarField& rAUf
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace fv