Skip to content
Snippets Groups Projects
Commit 2d1a22d6 authored by mattijs's avatar mattijs
Browse files

Merge branch 'master' of /home/dm4/OpenFOAM/OpenFOAM-dev

parents cbdba6b4 1ed0bb5e
Branches
Tags
No related merge requests found
......@@ -308,6 +308,7 @@ $(ddtSchemes)/backwardDdtScheme/backwardDdtSchemes.C
$(ddtSchemes)/boundedBackwardDdtScheme/boundedBackwardDdtScheme.C
$(ddtSchemes)/boundedBackwardDdtScheme/boundedBackwardDdtSchemes.C
$(ddtSchemes)/CrankNicholsonDdtScheme/CrankNicholsonDdtSchemes.C
$(ddtSchemes)/boundedDdtScheme/boundedDdtSchemes.C
d2dt2Schemes = finiteVolume/d2dt2Schemes
$(d2dt2Schemes)/d2dt2Scheme/d2dt2Schemes.C
......
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ 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 3 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, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "boundedDdtScheme.H"
#include "fvcDiv.H"
#include "fvcDdt.H"
#include "fvMatrices.H"
#include "fvmSup.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace fv
{
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class Type>
tmp<GeometricField<Type, fvPatchField, volMesh> >
boundedDdtScheme<Type>::fvcDdt
(
const dimensioned<Type>& dt
)
{
return scheme_().fvcDdt(dt);
}
template<class Type>
tmp<GeometricField<Type, fvPatchField, volMesh> >
boundedDdtScheme<Type>::fvcDdt
(
const GeometricField<Type, fvPatchField, volMesh>& vf
)
{
return scheme_().fvcDdt(vf);
}
template<class Type>
tmp<GeometricField<Type, fvPatchField, volMesh> >
boundedDdtScheme<Type>::fvcDdt
(
const dimensionedScalar& rho,
const GeometricField<Type, fvPatchField, volMesh>& vf
)
{
return scheme_().fvcDdt(rho, vf);
}
template<class Type>
tmp<GeometricField<Type, fvPatchField, volMesh> >
boundedDdtScheme<Type>::fvcDdt
(
const volScalarField& rho,
const GeometricField<Type, fvPatchField, volMesh>& vf
)
{
return scheme_().fvcDdt(rho, vf) - fvc::ddt(rho)*vf;
}
template<class Type>
tmp<fvMatrix<Type> >
boundedDdtScheme<Type>::fvmDdt
(
const GeometricField<Type, fvPatchField, volMesh>& vf
)
{
return scheme_().fvmDdt(vf);
}
template<class Type>
tmp<fvMatrix<Type> >
boundedDdtScheme<Type>::fvmDdt
(
const dimensionedScalar& rho,
const GeometricField<Type, fvPatchField, volMesh>& vf
)
{
return scheme_().fvmDdt(rho, vf);
}
template<class Type>
tmp<fvMatrix<Type> >
boundedDdtScheme<Type>::fvmDdt
(
const volScalarField& rho,
const GeometricField<Type, fvPatchField, volMesh>& vf
)
{
return scheme_().fvmDdt(rho, vf) - fvm::Sp(fvc::ddt(rho), vf);
}
template<class Type>
tmp<typename boundedDdtScheme<Type>::fluxFieldType>
boundedDdtScheme<Type>::fvcDdtPhiCorr
(
const volScalarField& rA,
const GeometricField<Type, fvPatchField, volMesh>& U,
const fluxFieldType& phi
)
{
return scheme_().fvcDdtPhiCorr(rA, U, phi);
}
template<class Type>
tmp<typename boundedDdtScheme<Type>::fluxFieldType>
boundedDdtScheme<Type>::fvcDdtPhiCorr
(
const volScalarField& rA,
const volScalarField& rho,
const GeometricField<Type, fvPatchField, volMesh>& U,
const fluxFieldType& phi
)
{
return scheme_().fvcDdtPhiCorr(rA, rho, U, phi);
}
template<class Type>
tmp<surfaceScalarField> boundedDdtScheme<Type>::meshPhi
(
const GeometricField<Type, fvPatchField, volMesh>& vf
)
{
return scheme_().meshPhi(vf);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace fv
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ 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 3 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, see <http://www.gnu.org/licenses/>.
Class
Foam::fv::boundedDdtScheme
Description
Bounded form of the selected ddt scheme.
Boundedness is achieved by subtracting ddt(phi)*vf or Sp(ddt(rho), vf)
which is non-conservative if ddt(rho) != 0 but conservative otherwise.
Can be used for the ddt of bounded scalar properties to improve stability
if insufficient convergence of the pressure equation causes temporary
divergence of the flux field.
SourceFiles
boundedDdtScheme.C
\*---------------------------------------------------------------------------*/
#ifndef boundedDdtScheme_H
#define boundedDdtScheme_H
#include "ddtScheme.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace fv
{
/*---------------------------------------------------------------------------*\
Class boundedDdtScheme Declaration
\*---------------------------------------------------------------------------*/
template<class Type>
class boundedDdtScheme
:
public fv::ddtScheme<Type>
{
// Private data
tmp<fv::ddtScheme<Type> > scheme_;
// Private Member Functions
//- Disallow default bitwise copy construct
boundedDdtScheme(const boundedDdtScheme&);
//- Disallow default bitwise assignment
void operator=(const boundedDdtScheme&);
public:
//- Runtime type information
TypeName("bounded");
// Constructors
//- Construct from mesh and Istream
boundedDdtScheme(const fvMesh& mesh, Istream& is)
:
ddtScheme<Type>(mesh, is),
scheme_
(
fv::ddtScheme<Type>::New(mesh, is)
)
{}
// Member Functions
//- Return mesh reference
const fvMesh& mesh() const
{
return fv::ddtScheme<Type>::mesh();
}
tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
(
const dimensioned<Type>&
);
tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
(
const GeometricField<Type, fvPatchField, volMesh>&
);
tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
(
const dimensionedScalar&,
const GeometricField<Type, fvPatchField, volMesh>&
);
tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
(
const volScalarField&,
const GeometricField<Type, fvPatchField, volMesh>&
);
tmp<fvMatrix<Type> > fvmDdt
(
const GeometricField<Type, fvPatchField, volMesh>&
);
tmp<fvMatrix<Type> > fvmDdt
(
const dimensionedScalar&,
const GeometricField<Type, fvPatchField, volMesh>&
);
tmp<fvMatrix<Type> > fvmDdt
(
const volScalarField&,
const GeometricField<Type, fvPatchField, volMesh>&
);
typedef typename ddtScheme<Type>::fluxFieldType fluxFieldType;
tmp<fluxFieldType> fvcDdtPhiCorr
(
const volScalarField& rA,
const GeometricField<Type, fvPatchField, volMesh>& U,
const fluxFieldType& phi
);
tmp<fluxFieldType> fvcDdtPhiCorr
(
const volScalarField& rA,
const volScalarField& rho,
const GeometricField<Type, fvPatchField, volMesh>& U,
const fluxFieldType& phi
);
tmp<surfaceScalarField> meshPhi
(
const GeometricField<Type, fvPatchField, volMesh>&
);
};
template<>
tmp<surfaceScalarField> boundedDdtScheme<scalar>::fvcDdtPhiCorr
(
const volScalarField& rA,
const volScalarField& U,
const surfaceScalarField& phi
);
template<>
tmp<surfaceScalarField> boundedDdtScheme<scalar>::fvcDdtPhiCorr
(
const volScalarField& rA,
const volScalarField& rho,
const volScalarField& U,
const surfaceScalarField& phi
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace fv
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
# include "boundedDdtScheme.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
\\/ 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 3 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, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "boundedDdtScheme.H"
#include "fvMesh.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace fv
{
makeFvDdtScheme(boundedDdtScheme)
}
}
// ************************************************************************* //
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment