From 49c79f073aa01dafa5f58acce6ab873f79491fc5 Mon Sep 17 00:00:00 2001 From: Henry Date: Wed, 2 Feb 2011 12:40:14 +0000 Subject: [PATCH] Higher-order upwinded: Reorganised linearUpwind and created an equivalent quadraticUpwind scheme from it --- src/finiteVolume/Make/files | 5 +- .../linearUpwind/linearUpwind.C | 15 +- .../linearUpwind/linearUpwind.H | 11 +- .../linearUpwind/linearUpwindV.C | 2 +- .../linearUpwind/linearUpwindV.H | 2 +- .../schemes/quadraticUpwind/quadraticUpwind.C | 38 +++++ .../schemes/quadraticUpwind/quadraticUpwind.H | 134 ++++++++++++++++++ 7 files changed, 189 insertions(+), 18 deletions(-) rename src/finiteVolume/interpolation/surfaceInterpolation/{limitedSchemes => schemes}/linearUpwind/linearUpwind.C (92%) rename src/finiteVolume/interpolation/surfaceInterpolation/{limitedSchemes => schemes}/linearUpwind/linearUpwind.H (92%) rename src/finiteVolume/interpolation/surfaceInterpolation/{limitedSchemes => schemes}/linearUpwind/linearUpwindV.C (98%) rename src/finiteVolume/interpolation/surfaceInterpolation/{limitedSchemes => schemes}/linearUpwind/linearUpwindV.H (98%) create mode 100644 src/finiteVolume/interpolation/surfaceInterpolation/schemes/quadraticUpwind/quadraticUpwind.C create mode 100644 src/finiteVolume/interpolation/surfaceInterpolation/schemes/quadraticUpwind/quadraticUpwind.H diff --git a/src/finiteVolume/Make/files b/src/finiteVolume/Make/files index 2e4dd38bb1..1e8508d860 100644 --- a/src/finiteVolume/Make/files +++ b/src/finiteVolume/Make/files @@ -243,13 +243,14 @@ $(schemes)/cubicUpwindFit/cubicUpwindFit.C $(schemes)/quadraticLinearPureUpwindFit/quadraticLinearPureUpwindFit.C */ $(schemes)/linearPureUpwindFit/linearPureUpwindFit.C +$(schemes)/linearUpwind/linearUpwind.C +$(schemes)/linearUpwind/linearUpwindV.C +$(schemes)/quadraticUpwind/quadraticUpwind.C limitedSchemes = $(surfaceInterpolation)/limitedSchemes $(limitedSchemes)/limitedSurfaceInterpolationScheme/limitedSurfaceInterpolationSchemes.C $(limitedSchemes)/upwind/upwind.C $(limitedSchemes)/blended/blended.C -$(limitedSchemes)/linearUpwind/linearUpwind.C -$(limitedSchemes)/linearUpwind/linearUpwindV.C $(limitedSchemes)/Gamma/Gamma.C $(limitedSchemes)/SFCD/SFCD.C $(limitedSchemes)/Minmod/Minmod.C diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwind.C b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwind.C similarity index 92% rename from src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwind.C rename to src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwind.C index 7f17534c94..77b6ae0d3b 100644 --- a/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwind.C +++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwind.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd. + \\ / A nd | Copyright (C) 2004-2011 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -25,7 +25,6 @@ License #include "linearUpwind.H" #include "fvMesh.H" -#include "zeroGradientFvPatchField.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // @@ -85,16 +84,8 @@ Foam::linearUpwind::correction forAll(faceFlux, facei) { - if (faceFlux[facei] > 0) - { - label own = owner[facei]; - sfCorr[facei] = (Cf[facei] - C[own]) & gradVf[own]; - } - else - { - label nei = neighbour[facei]; - sfCorr[facei] = (Cf[facei] - C[nei]) & gradVf[nei]; - } + label celli = (faceFlux[facei] > 0) ? owner[facei] : neighbour[facei]; + sfCorr[facei] = (Cf[facei] - C[celli]) & gradVf[celli]; } diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwind.H b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwind.H similarity index 92% rename from src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwind.H rename to src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwind.H index 603294c189..232172ba0d 100644 --- a/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwind.H +++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwind.H @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd. + \\ / A nd | Copyright (C) 2004-2011 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -26,7 +26,8 @@ Class Description linearUpwind interpolation scheme class derived from upwind and returns - upwind weighting factors but also applies an explicit correction. + upwind weighting factors and also applies a gradient-based explicit + correction. SourceFiles linearUpwind.C @@ -67,6 +68,12 @@ class linearUpwind //- Disallow default bitwise assignment void operator=(const linearUpwind&); + //- Hide the limiter because this is not formally a limited scheme + virtual tmp limiter + ( + const GeometricField& + ) const; + public: diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwindV.C b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwindV.C similarity index 98% rename from src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwindV.C rename to src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwindV.C index 695637e383..0747a3fb2b 100644 --- a/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwindV.C +++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwindV.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd. + \\ / A nd | Copyright (C) 2004-2011 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- License diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwindV.H b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwindV.H similarity index 98% rename from src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwindV.H rename to src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwindV.H index 6cd31b836a..262d54115d 100644 --- a/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/linearUpwind/linearUpwindV.H +++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwindV.H @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd. + \\ / A nd | Copyright (C) 2004-2011 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- License diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/quadraticUpwind/quadraticUpwind.C b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/quadraticUpwind/quadraticUpwind.C new file mode 100644 index 0000000000..d0779945d3 --- /dev/null +++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/quadraticUpwind/quadraticUpwind.C @@ -0,0 +1,38 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Copyright (C) 2004-2011 OpenCFD Ltd. + \\/ 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 2 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, write to the Free Software Foundation, + Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA + +\*---------------------------------------------------------------------------*/ + +#include "quadraticUpwind.H" + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +namespace Foam +{ + //makeSurfaceInterpolationScheme(quadraticUpwind); + makeSurfaceInterpolationTypeScheme(quadraticUpwind, scalar); + makeSurfaceInterpolationTypeScheme(quadraticUpwind, vector); +} + +// ************************************************************************* // diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/quadraticUpwind/quadraticUpwind.H b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/quadraticUpwind/quadraticUpwind.H new file mode 100644 index 0000000000..a2fb36f7fb --- /dev/null +++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/quadraticUpwind/quadraticUpwind.H @@ -0,0 +1,134 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Copyright (C) 2004-2011 OpenCFD Ltd. + \\/ 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 . + +Class + quadraticUpwind + +Description + quadraticUpwind interpolation scheme class derived from linearUpwind and + returns blended linear/upwind weighting factors and also applies a explicit + gradient-based correction obtained from the linearUpwind scheme. + +SourceFiles + quadraticUpwind.C + +\*---------------------------------------------------------------------------*/ + +#ifndef quadraticUpwind_H +#define quadraticUpwind_H + +#include "linearUpwind.H" + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +namespace Foam +{ + +/*---------------------------------------------------------------------------*\ + Class quadraticUpwind Declaration +\*---------------------------------------------------------------------------*/ + +template +class quadraticUpwind +: + public linearUpwind +{ + // Private Member Functions + + //- Disallow default bitwise copy construct + quadraticUpwind(const quadraticUpwind&); + + //- Disallow default bitwise assignment + void operator=(const quadraticUpwind&); + + +public: + + //- Runtime type information + TypeName("quadraticUpwind"); + + + // Constructors + + //- Construct from mesh and Istream + quadraticUpwind + ( + const fvMesh& mesh, + Istream& schemeData + ) + : + linearUpwind(mesh, schemeData) + {} + + //- Construct from mesh, faceFlux and Istream + quadraticUpwind + ( + const fvMesh& mesh, + const surfaceScalarField& faceFlux, + Istream& schemeData + ) + : + linearUpwind(mesh, faceFlux, schemeData) + {} + + + // Member Functions + + //- Return the interpolation weighting factors + virtual tmp weights + ( + const GeometricField& + ) const + { + return + 0.75*this->mesh().surfaceInterpolation::weights() + + 0.25*linearUpwind::weights(); + } + + //- Return true if this scheme uses an explicit correction + virtual bool corrected() const + { + return true; + } + + //- Return the explicit correction to the face-interpolate + virtual tmp > + correction + ( + const GeometricField& vf + ) const + { + return 0.25*linearUpwind::correction(vf); + } +}; + + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +} // End namespace Foam + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +#endif + +// ************************************************************************* //