Skip to content
Snippets Groups Projects
Commit 391c1eaa authored by Kutalmış Berçin's avatar Kutalmış Berçin Committed by Andrew Heather
Browse files

ENH: snGradSchemes: new scheme relaxedSnGrad

Surface gradient scheme with under-/over-relaxed
full or limited explicit non-orthogonal correction.

A minimal example by using system/fvSchemes:

  snGradSchemes
  {
      snGrad(<term>)       relaxed;
  }

and by using system/fvSolution:

  relaxationFactors
  {
      fields
      {
          snGrad(<term>)   <relaxation factor>;
      }
  }
parent 37a6fc8d
1 merge request!502ENH: New schemes for non-orthogonal meshes
......@@ -456,6 +456,7 @@ $(snGradSchemes)/orthogonalSnGrad/orthogonalSnGrads.C
$(snGradSchemes)/quadraticFitSnGrad/quadraticFitSnGrads.C
$(snGradSchemes)/linearFitSnGrad/linearFitSnGrads.C
$(snGradSchemes)/skewCorrectedSnGrad/skewCorrectedSnGrads.C
$(snGradSchemes)/relaxedSnGrad/relaxedSnGrads.C
convectionSchemes = finiteVolume/convectionSchemes
$(convectionSchemes)/convectionScheme/convectionSchemes.C
......
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 "fv.H"
#include "relaxedSnGrad.H"
#include "volFields.H"
#include "surfaceFields.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh>>
Foam::fv::relaxedSnGrad<Type>::correction
(
const GeometricField<Type, fvPatchField, volMesh>& vf
) const
{
typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfFieldType;
// Calculate explicit correction field
tmp<SurfFieldType> tcorrection = correctedScheme_().correction(vf);
// Retrieve relaxation factor value
const word fieldName(vf.name());
const word oldFieldName(fieldName + "_0");
const scalar relax =
vf.mesh().fieldRelaxationFactor("snGrad("+fieldName+")");
// Return explicit correction field if
// previous-time step correction is unavailable
const objectRegistry& obr = vf.db();
if (!obr.foundObject<SurfFieldType>(oldFieldName))
{
SurfFieldType* oldCorrection =
new SurfFieldType(oldFieldName, tcorrection());
oldCorrection->store();
}
// Return under/over-relaxed explicit correction field
tmp<SurfFieldType> trelaxedCorrection(new SurfFieldType(tcorrection()));
SurfFieldType& oldCorrection =
obr.lookupObjectRef<SurfFieldType>(oldFieldName);
trelaxedCorrection.ref() *= relax;
trelaxedCorrection.ref() += (scalar(1) - relax)*oldCorrection;
oldCorrection = tcorrection;
return trelaxedCorrection;
}
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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::relaxedSnGrad
Group
grpFvSnGradSchemes
Description
Surface gradient scheme with under-/over-relaxed
full or limited explicit non-orthogonal correction.
Usage
Minimal example by using \c system/fvSchemes:
\verbatim
snGradSchemes
{
snGrad(<term>) relaxed;
}
\endverbatim
and by using \c system/fvSolution:
\verbatim
relaxationFactors
{
fields
{
snGrad(<term>) <relaxation factor>;
}
}
\endverbatim
SourceFiles
relaxedSnGrad.C
\*---------------------------------------------------------------------------*/
#ifndef relaxedSnGrad_H
#define relaxedSnGrad_H
#include "correctedSnGrad.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace fv
{
/*---------------------------------------------------------------------------*\
Class relaxedSnGrad Declaration
\*---------------------------------------------------------------------------*/
template<class Type>
class relaxedSnGrad
:
public snGradScheme<Type>
{
// Private Data
//- Type of correction scheme
tmp<snGradScheme<Type>> correctedScheme_;
// Private Member Functions
//- No copy assignment
void operator=(const relaxedSnGrad&) = delete;
public:
//- Runtime type information
TypeName("relaxed");
// Constructors
//- Construct from mesh
relaxedSnGrad(const fvMesh& mesh)
:
snGradScheme<Type>(mesh),
correctedScheme_(new correctedSnGrad<Type>(this->mesh()))
{}
//- Construct from mesh and data stream
relaxedSnGrad(const fvMesh& mesh, Istream& schemeData)
:
snGradScheme<Type>(mesh),
correctedScheme_(new correctedSnGrad<Type>(this->mesh()))
{}
//- Destructor
virtual ~relaxedSnGrad() = default;
// Member Functions
//- Return the interpolation weighting factors for the given field
virtual tmp<surfaceScalarField> deltaCoeffs
(
const GeometricField<Type, fvPatchField, volMesh>&
) const
{
return this->mesh().nonOrthDeltaCoeffs();
}
//- Return true if this scheme uses an explicit correction
virtual bool corrected() const noexcept
{
return true;
}
//- Return the explicit correction to the relaxedSnGrad
//- for the given field using the gradients of the field components
virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
correction(const GeometricField<Type, fvPatchField, volMesh>&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace fv
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "relaxedSnGrad.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 "relaxedSnGrad.H"
#include "fvMesh.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makeSnGradScheme(relaxedSnGrad)
// ************************************************************************* //
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