From bb3c0c3593b21b62b85ad7c6c337f6c253993fb7 Mon Sep 17 00:00:00 2001 From: Chris Greenshields Date: Tue, 20 Sep 2016 18:38:15 +0100 Subject: [PATCH] Maxwell model for viscoelasticity using the upper-convected time derivative of the stress tensor. See http://en.wikipedia.org/wiki/Upper-convected_Maxwell_model The model includes an additional viscosity (nu) from the transport model from which it is instantiated, which makes it equivalent to the Oldroyd-B model for the case of an incompressible transport model. See https://en.wikipedia.org/wiki/Oldroyd-B_model --- .../turbulentFluidThermoModels.C | 3 + .../laminar/Maxwell/Maxwell.C | 240 ++++++++++++++++++ .../laminar/Maxwell/Maxwell.H | 174 +++++++++++++ .../pimpleFoam/laminar/planarContraction/0/U | 46 ++++ .../pimpleFoam/laminar/planarContraction/0/p | 42 +++ .../laminar/planarContraction/0/sigma | 38 +++ .../constant/transportProperties | 21 ++ .../constant/turbulenceProperties | 32 +++ .../planarContraction/system/blockMeshDict | 112 ++++++++ .../planarContraction/system/controlDict | 52 ++++ .../planarContraction/system/fvSchemes | 54 ++++ .../planarContraction/system/fvSolution | 67 +++++ .../laminar/planarContraction/system/graphs | 62 +++++ .../pimpleFoam/laminar/planarPoiseuille/0/U | 32 +++ .../pimpleFoam/laminar/planarPoiseuille/0/p | 31 +++ .../laminar/planarPoiseuille/0/sigma | 32 +++ .../laminar/planarPoiseuille/Allclean | 12 + .../laminar/planarPoiseuille/Allrun | 15 ++ .../planarPoiseuille/constant/fvOptions | 38 +++ .../constant/transportProperties | 21 ++ .../constant/turbulenceProperties | 32 +++ .../planarPoiseuille/system/blockMeshDict | 96 +++++++ .../planarPoiseuille/system/controlDict | 54 ++++ .../laminar/planarPoiseuille/system/fvSchemes | 54 ++++ .../planarPoiseuille/system/fvSolution | 67 +++++ .../laminar/planarPoiseuille/system/probes | 21 ++ .../laminar/planarPoiseuille/system/residuals | 19 ++ .../planarPoiseuille/system/singleGraph | 48 ++++ .../validation/WatersKing/Make/files | 3 + .../validation/WatersKing/Make/options | 14 + .../validation/WatersKing/WatersKing.C | 146 +++++++++++ .../validation/WatersKing/createFields.H | 64 +++++ .../planarPoiseuille/validation/createGraph | 27 ++ 33 files changed, 1769 insertions(+) create mode 100644 src/TurbulenceModels/turbulenceModels/laminar/Maxwell/Maxwell.C create mode 100644 src/TurbulenceModels/turbulenceModels/laminar/Maxwell/Maxwell.H create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/U create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/p create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/sigma create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarContraction/constant/transportProperties create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarContraction/constant/turbulenceProperties create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/blockMeshDict create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/controlDict create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/fvSchemes create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/fvSolution create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/graphs create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/U create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/p create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/sigma create mode 100755 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/Allclean create mode 100755 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/Allrun create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/fvOptions create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/transportProperties create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/turbulenceProperties create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/blockMeshDict create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/controlDict create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/fvSchemes create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/fvSolution create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/probes create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/residuals create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/singleGraph create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/Make/files create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/Make/options create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/WatersKing.C create mode 100644 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/createFields.H create mode 100755 tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/createGraph diff --git a/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C b/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C index 83c184f49a..6506c21387 100644 --- a/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C +++ b/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C @@ -45,6 +45,9 @@ makeBaseTurbulenceModel #include "Stokes.H" makeLaminarModel(Stokes); +#include "Maxwell.H" +makeLaminarModel(Maxwell); + // -------------------------------------------------------------------------- // // RAS models diff --git a/src/TurbulenceModels/turbulenceModels/laminar/Maxwell/Maxwell.C b/src/TurbulenceModels/turbulenceModels/laminar/Maxwell/Maxwell.C new file mode 100644 index 0000000000..1a7f409cd1 --- /dev/null +++ b/src/TurbulenceModels/turbulenceModels/laminar/Maxwell/Maxwell.C @@ -0,0 +1,240 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Copyright (C) 2016 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 . + +\*---------------------------------------------------------------------------*/ + +#include "Maxwell.H" +#include "fvOptions.H" + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +namespace Foam +{ +namespace laminarModels +{ + +// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // + +template +Maxwell::Maxwell +( + const alphaField& alpha, + const rhoField& rho, + const volVectorField& U, + const surfaceScalarField& alphaRhoPhi, + const surfaceScalarField& phi, + const transportModel& transport, + const word& propertiesName, + const word& type +) +: + laminarModel + ( + type, + alpha, + rho, + U, + alphaRhoPhi, + phi, + transport, + propertiesName + ), + + nuM_ + ( + dimensioned + ( + "nuM", + dimViscosity, + this->coeffDict_.lookup("nuM") + ) + ), + + lambda_ + ( + dimensioned + ( + "lambda", + dimTime, + this->coeffDict_.lookup("lambda") + ) + ), + + sigma_ + ( + IOobject + ( + IOobject::groupName("sigma", U.group()), + this->runTime_.timeName(), + this->mesh_, + IOobject::MUST_READ, + IOobject::AUTO_WRITE + ), + this->mesh_ + ) +{ + if (type == typeName) + { + this->printCoeffs(type); + } +} + + +// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // + +template +bool Maxwell::read() +{ + if (laminarModel::read()) + { + nuM_.readIfPresent(this->coeffDict()); + lambda_.readIfPresent(this->coeffDict()); + + return true; + } + else + { + return false; + } +} + +template +tmp +Maxwell::R() const +{ + return sigma_; +} + +template +tmp +Maxwell::devRhoReff() const +{ + return tmp + ( + new volSymmTensorField + ( + IOobject + ( + IOobject::groupName("devRhoReff", this->U_.group()), + this->runTime_.timeName(), + this->mesh_, + IOobject::NO_READ, + IOobject::NO_WRITE + ), + this->alpha_*this->rho_*sigma_ + - (this->alpha_*this->rho_*this->nu()) + *dev(twoSymm(fvc::grad(this->U_))) + ) + ); +} + + +template +tmp +Maxwell::divDevRhoReff +( + volVectorField& U +) const +{ + return + ( + fvc::div + ( + this->alpha_*this->rho_*this->nuM_*fvc::grad(U) + ) + + fvc::div(this->alpha_*this->rho_*sigma_) + - fvc::div(this->alpha_*this->rho_*this->nu()*dev2(T(fvc::grad(U)))) + - fvm::laplacian(this->alpha_*this->rho_*nu0(), U) + ); +} + + +template +tmp +Maxwell::divDevRhoReff +( + const volScalarField& rho, + volVectorField& U +) const +{ + return + ( + fvc::div + ( + this->alpha_*rho*this->nuM_*fvc::grad(U) + ) + + fvc::div(this->alpha_*rho*sigma_) + - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U)))) + - fvm::laplacian(this->alpha_*rho*nu0(), U) + ); +} + + +template +void Maxwell::correct() +{ + // Local references + const alphaField& alpha = this->alpha_; + const rhoField& rho = this->rho_; + const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_; + const volVectorField& U = this->U_; + volSymmTensorField& sigma = this->sigma_; + fv::options& fvOptions(fv::options::New(this->mesh_)); + + laminarModel::correct(); + + tmp tgradU(fvc::grad(U)); + const volTensorField& gradU = tgradU(); + dimensionedScalar rLambda = 1.0/(lambda_); + + // Note sigma is positive on lhs of momentum eqn + volSymmTensorField P + ( + twoSymm(sigma & gradU) + - nuM_*rLambda*twoSymm(gradU) + ); + + // Viscoelastic stress equation + tmp sigmaEqn + ( + fvm::ddt(alpha, rho, sigma) + + fvm::div(alphaRhoPhi, sigma) + + fvm::Sp(alpha*rho*rLambda, sigma) + == + alpha*rho*P + + fvOptions(alpha, rho, sigma) + ); + + sigmaEqn.ref().relax(); + fvOptions.constrain(sigmaEqn.ref()); + solve(sigmaEqn); + fvOptions.correct(sigma_); +} + + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +} // End namespace laminarModels +} // End namespace Foam + +// ************************************************************************* // diff --git a/src/TurbulenceModels/turbulenceModels/laminar/Maxwell/Maxwell.H b/src/TurbulenceModels/turbulenceModels/laminar/Maxwell/Maxwell.H new file mode 100644 index 0000000000..85c4d38d51 --- /dev/null +++ b/src/TurbulenceModels/turbulenceModels/laminar/Maxwell/Maxwell.H @@ -0,0 +1,174 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Copyright (C) 2016 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 . + +Class + Foam::laminarModels::Maxwell + +Group + grpLaminar + +Description + Maxwell model for viscoelasticity using the upper-convected time + derivative of the stress tensor. + See http://en.wikipedia.org/wiki/Upper-convected_Maxwell_model + + The model includes an additional viscosity (nu) from the transport + model from which it is instantiated, which makes it equivalent to + the Oldroyd-B model for the case of an incompressible transport + model (where nu is non-zero). + See https://en.wikipedia.org/wiki/Oldroyd-B_model + + Reference: + \verbatim + Amoreira, L. J., & Oliveira, P. J. (2010). + Comparison of different formulations for the numerical calculation + of unsteady incompressible viscoelastic fluid flow. + Adv. Appl. Math. Mech, 4, 483-502. + \endverbatim + +SourceFiles + Maxwell.C + +\*---------------------------------------------------------------------------*/ + +#ifndef Maxwell_H +#define Maxwell_H + +#include "laminarModel.H" + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +namespace Foam +{ +namespace laminarModels +{ + +/*---------------------------------------------------------------------------*\ + Class Maxwell Declaration +\*---------------------------------------------------------------------------*/ + +template +class Maxwell +: + public laminarModel +{ + +protected: + + // Protected data + + // Model coefficients + + dimensionedScalar nuM_; + dimensionedScalar lambda_; + + + // Fields + + volSymmTensorField sigma_; + + + // Protected Member Functions + + //- Return the turbulence viscosity + tmp nu0() const + { + return this->nu() + nuM_; + } + + +public: + + typedef typename BasicTurbulenceModel::alphaField alphaField; + typedef typename BasicTurbulenceModel::rhoField rhoField; + typedef typename BasicTurbulenceModel::transportModel transportModel; + + + //- Runtime type information + TypeName("Maxwell"); + + + // Constructors + + //- Construct from components + Maxwell + ( + const alphaField& alpha, + const rhoField& rho, + const volVectorField& U, + const surfaceScalarField& alphaRhoPhi, + const surfaceScalarField& phi, + const transportModel& transport, + const word& propertiesName = turbulenceModel::propertiesName, + const word& type = typeName + ); + + + //- Destructor + virtual ~Maxwell() + {} + + + // Member Functions + + //- Read model coefficients if they have changed + virtual bool read(); + + //- Return the Reynolds stress tensor + virtual tmp R() const; + + //- Return the effective stress tensor + virtual tmp devRhoReff() const; + + //- Return the source term for the momentum equation + virtual tmp divDevRhoReff(volVectorField& U) const; + + //- Return the source term for the momentum equation + virtual tmp divDevRhoReff + ( + const volScalarField& rho, + volVectorField& U + ) const; + + //- Solve the turbulence equations and correct eddy-Viscosity and + // related properties + virtual void correct(); +}; + + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +} // End namespace laminarModels +} // End namespace Foam + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +#ifdef NoRepository + #include "Maxwell.C" +#endif + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +#endif + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/U b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/U new file mode 100644 index 0000000000..ee40b539a6 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/U @@ -0,0 +1,46 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class volVectorField; + object U; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +Uinlet (0.03876 0 0); + +dimensions [0 1 -1 0 0 0 0]; + +internalField uniform (0 0 0); + +boundaryField +{ + inlet + { + type fixedValue; + value uniform $Uinlet; + } + + outlet + { + type zeroGradient; + value uniform (0 0 0); + } + + wall + { + type fixedValue; + value uniform (0 0 0); + } + + #includeEtc "caseDicts/setConstraintTypes" +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/p b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/p new file mode 100644 index 0000000000..f2d9ffb3b9 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/p @@ -0,0 +1,42 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class volScalarField; + object p; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +dimensions [0 2 -2 0 0 0 0]; + +internalField uniform 0; + +boundaryField +{ + inlet + { + type zeroGradient; + } + + outlet + { + type fixedValue; + value uniform 0; + } + + wall + { + type zeroGradient; + } + + #includeEtc "caseDicts/setConstraintTypes" +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/sigma b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/sigma new file mode 100644 index 0000000000..259ccfb7b7 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/0/sigma @@ -0,0 +1,38 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class volSymmTensorField; + object R; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +dimensions [0 2 -2 0 0 0 0]; + +internalField uniform (0 0 0 0 0 0); + +boundaryField +{ + inlet + { + type fixedValue; + value $internalField; + } + + ".*" + { + type zeroGradient; + } + + #includeEtc "caseDicts/setConstraintTypes" +} + + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarContraction/constant/transportProperties b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/constant/transportProperties new file mode 100644 index 0000000000..35bbd73e1e --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/constant/transportProperties @@ -0,0 +1,21 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: 3.0.1 | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object transportProperties; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +transportModel Newtonian; + +nu [0 2 -1 0 0 0 0] 1e-5; + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarContraction/constant/turbulenceProperties b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/constant/turbulenceProperties new file mode 100644 index 0000000000..6503cde5c3 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/constant/turbulenceProperties @@ -0,0 +1,32 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object turbulenceProperties; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +simulationType laminar; + +laminar +{ + laminarModel Maxwell; + + MaxwellCoeffs + { + nuM 0.002; + lambda 0.03; + } + + printCoeffs on; +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/blockMeshDict b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/blockMeshDict new file mode 100644 index 0000000000..b466803448 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/blockMeshDict @@ -0,0 +1,112 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object blockMeshDict; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +convertToMeters 0.0032; + +vertices +( + (-40 0 -1) + ( 0 0 -1) + ( 30 0 -1) + (-40 1 -1) + ( 0 1 -1) + ( 30 1 -1) + (-40 4 -1) + ( 0 4 -1) + + (-40 0 1) + ( 0 0 1) + ( 30 0 1) + (-40 1 1) + ( 0 1 1) + ( 30 1 1) + (-40 4 1) + ( 0 4 1) +); + +blocks +( + hex (0 1 4 3 8 9 12 11) (40 12 1) simpleGrading (0.02 0.4 1) + hex (1 2 5 4 9 10 13 12) (30 12 1) simpleGrading (50 0.4 1) + hex (3 4 7 6 11 12 15 14) (40 24 1) simpleGrading (0.02 ((0.5 0.5 4.0) (0.5 0.5 0.25)) 1) +); + +edges +( +); + +boundary +( + inlet + { + type patch; + faces + ( + (0 3 11 8) + (3 6 14 11) + ); + } + + walls + { + type wall; + faces + ( + (6 7 15 14) + (7 4 12 15) + (4 5 13 12) + ); + } + + outlet + { + type patch; + faces + ( + (2 5 13 10) + ); + } + + centreline + { + type symmetryPlane; + faces + ( + (0 1 9 8) + (1 2 10 9) + ); + } + + frontAndBack + { + type empty; + faces + ( + (0 1 4 3) + (3 4 7 6) + (1 2 5 4) + (8 9 12 11) + (11 12 15 14) + (9 10 13 12) + ); + } +); + +mergePatchPairs +( +); + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/controlDict b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/controlDict new file mode 100644 index 0000000000..66f6d6cfde --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/controlDict @@ -0,0 +1,52 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object controlDict; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +application pimpleFoam; + +startFrom latestTime; + +startTime 0; + +stopAt endTime; + +endTime 0.25; + +deltaT 2e-4; + +writeControl runTime; + +writeInterval 0.01; + +purgeWrite 0; + +writeFormat ascii; + +writePrecision 8; + +writeCompression off; + +timeFormat general; + +timePrecision 6; + +runTimeModifiable true; + +functions +{ + #includeFunc graphs +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/fvSchemes b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/fvSchemes new file mode 100644 index 0000000000..9f4397a084 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/fvSchemes @@ -0,0 +1,54 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object fvSchemes; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +ddtSchemes +{ + default backward; +} + +gradSchemes +{ + default Gauss linear; +} + +divSchemes +{ + default none; + + div(phi,U) Gauss linearUpwind grad(U); + div(phi,sigma) Gauss vanAlbada; + + div(sigma) Gauss linear; + div((nu*dev2(T(grad(U))))) Gauss linear; + div((nuM*grad(U))) Gauss linear; +} + +laplacianSchemes +{ + default Gauss linear corrected; +} + +interpolationSchemes +{ + default linear; +} + +snGradSchemes +{ + default corrected; +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/fvSolution b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/fvSolution new file mode 100644 index 0000000000..302add3af7 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/fvSolution @@ -0,0 +1,67 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object fvSolution; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +solvers +{ + p + { + solver GAMG; + smoother DIC; + tolerance 1e-6; + relTol 0.05; + } + + "(U|sigma)" + { + solver smoothSolver; + smoother symGaussSeidel; + tolerance 1e-6; + relTol 0.1; + } + + pFinal + { + $p; + relTol 0; + } + + "(U|sigma)Final" + { + $U; + relTol 0; + } +} + +PIMPLE +{ + momentumPredictor off; + nOuterCorrectors 15; + nCorrectors 1; + nNonOrthogonalCorrectors 0; + pRefCell 0; + pRefValue 0; + turbOnFinalIterOnly no; +} + +relaxationFactors +{ + equations + { + ".*" 1; + } +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/graphs b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/graphs new file mode 100644 index 0000000000..6c622052c3 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarContraction/system/graphs @@ -0,0 +1,62 @@ +/*--------------------------------*- C++ -*----------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Web: www.OpenFOAM.org + \\/ M anipulation | +------------------------------------------------------------------------------- +Description + Writes graph data for specified fields along a line, specified by start + and end points. + +\*---------------------------------------------------------------------------*/ + + +// Sampling and I/O settings +#includeEtc "caseDicts/postProcessing/graphs/sampleDict.cfg" + +// Override settings here, e.g. +// setConfig { type midPoint; } + +type sets; +libs ("libsampling.so"); + +writeControl writeTime; + +interpolationScheme cellPoint; + +setFormat raw; + +setConfig +{ + type midPoint; // midPoint + axis distance; // x, y, z, xyz +} + +sets +( + lineA + { + $setConfig; + start (-0.0016 0 0); + end (-0.0016 0.0128 0); + } + + lineB + { + $setConfig; + start (-0.0048 0 0); + end (-0.0048 0.0128 0); + } + + lineC + { + $setConfig; + start (-0.032 0 0); + end (-0.032 0.0128 0); + } +); + +fields (U); + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/U b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/U new file mode 100644 index 0000000000..6a258e941b --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/U @@ -0,0 +1,32 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class volVectorField; + object U; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +dimensions [0 1 -1 0 0 0 0]; + +internalField uniform (0 0 0); + +boundaryField +{ + wall + { + type fixedValue; + value uniform (0 0 0); + } + + #includeEtc "caseDicts/setConstraintTypes" +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/p b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/p new file mode 100644 index 0000000000..af35540f51 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/p @@ -0,0 +1,31 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class volScalarField; + object p; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +dimensions [0 2 -2 0 0 0 0]; + +internalField uniform 0; + +boundaryField +{ + walls + { + type zeroGradient; + } + + #includeEtc "caseDicts/setConstraintTypes" +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/sigma b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/sigma new file mode 100644 index 0000000000..491709d9ff --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/0/sigma @@ -0,0 +1,32 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class volSymmTensorField; + object R; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +dimensions [0 2 -2 0 0 0 0]; + +internalField uniform (0 0 0 0 0 0); + +boundaryField +{ + walls + { + type zeroGradient; + } + + #includeEtc "caseDicts/setConstraintTypes" +} + + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/Allclean b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/Allclean new file mode 100755 index 0000000000..0f2f7ba75f --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/Allclean @@ -0,0 +1,12 @@ +#!/bin/sh +cd ${0%/*} || exit 1 # Run from this directory + +# Source tutorial clean functions +. $WM_PROJECT_DIR/bin/tools/CleanFunctions + +cleanCase +rm -rf postProcessing *.dat validation/*.eps + +wclean validation/WatersKing + +#------------------------------------------------------------------------------ diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/Allrun b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/Allrun new file mode 100755 index 0000000000..1737ee272a --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/Allrun @@ -0,0 +1,15 @@ +#!/bin/sh +cd ${0%/*} || exit 1 # Run from this directory + +# Source tutorial run functions +. $WM_PROJECT_DIR/bin/tools/RunFunctions + +runApplication blockMesh +runApplication $(getApplication) + +wmake validation/WatersKing +runApplication WatersKing + +( cd validation && ./createGraph ) + +#------------------------------------------------------------------------------ diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/fvOptions b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/fvOptions new file mode 100644 index 0000000000..5a267f51d2 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/fvOptions @@ -0,0 +1,38 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + location "constant"; + object fvOptions; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +momentumSource +{ + type vectorSemiImplicitSource; + active yes; + + vectorSemiImplicitSourceCoeffs + { + timeStart 0.0; + duration 1000; + selectionMode all; + + volumeMode specific; + injectionRateSuSp + { + U ((5 0 0) 0); + } + } +} + + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/transportProperties b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/transportProperties new file mode 100644 index 0000000000..ca57ed7652 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/transportProperties @@ -0,0 +1,21 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object transportProperties; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +transportModel Newtonian; + +nu [0 2 -1 0 0 0 0] 0.1; // kinematic -> 0.002 dynamic + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/turbulenceProperties b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/turbulenceProperties new file mode 100644 index 0000000000..e7cf61543f --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/constant/turbulenceProperties @@ -0,0 +1,32 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object turbulenceProperties; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +simulationType laminar; + +laminar +{ + laminarModel Maxwell; + + MaxwellCoeffs + { + nuM 1; + lambda 5; + } + + printCoeffs on; +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/blockMeshDict b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/blockMeshDict new file mode 100644 index 0000000000..99174900c4 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/blockMeshDict @@ -0,0 +1,96 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object blockMeshDict; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +convertToMeters 0.1; + +vertices +( + (-1 0 -1) + ( 1 0 -1) + ( 1 10 -1) + (-1 10 -1) + + (-1 0 1) + ( 1 0 1) + ( 1 10 1) + (-1 10 1) +); + +blocks +( + hex (0 1 2 3 4 5 6 7) (1 40 1) simpleGrading (1 4 1) +); + +edges +( +); + +boundary +( + left + { + type cyclic; + neighbourPatch right; + faces + ( + (0 3 7 4) + ); + } + + right + { + type cyclic; + neighbourPatch left; + faces + ( + (1 2 6 5) + ); + } + + walls + { + type wall; + faces + ( + (0 1 5 4) + ); + } + + centreline + { + type symmetryPlane; + faces + ( + (2 3 7 6) + ); + } + + frontAndBack + { + type empty; + faces + ( + (0 1 2 3) + (4 5 6 7) + ); + } +); + +mergePatchPairs +( +); + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/controlDict b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/controlDict new file mode 100644 index 0000000000..a986b1ae75 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/controlDict @@ -0,0 +1,54 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object controlDict; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +application pimpleFoam; + +startFrom startTime; + +startTime 0; + +stopAt endTime; + +endTime 25; + +deltaT 5e-3; + +writeControl runTime; + +writeInterval 1; + +purgeWrite 0; + +writeFormat ascii; + +writePrecision 8; + +writeCompression off; + +timeFormat general; + +timePrecision 6; + +runTimeModifiable true; + +functions +{ + #includeFunc residuals + #includeFunc singleGraph + #includeFunc probes +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/fvSchemes b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/fvSchemes new file mode 100644 index 0000000000..c7802ef935 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/fvSchemes @@ -0,0 +1,54 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object fvSchemes; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +ddtSchemes +{ + default backward; +} + +gradSchemes +{ + default Gauss linear; +} + +divSchemes +{ + default none; + + div(phi,U) Gauss linearUpwind grad(U); + div(phi,sigma) Gauss vanAlbada; + + div(sigma) Gauss linear; + div((nu*dev2(T(grad(U))))) Gauss linear; + div((nuM*grad(U))) Gauss linear; +} + +laplacianSchemes +{ + default Gauss linear uncorrected; +} + +interpolationSchemes +{ + default linear; +} + +snGradSchemes +{ + default uncorrected; +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/fvSolution b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/fvSolution new file mode 100644 index 0000000000..c82c0fb010 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/fvSolution @@ -0,0 +1,67 @@ +/*--------------------------------*- C++ -*----------------------------------*\ +| ========= | | +| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox | +| \\ / O peration | Version: dev | +| \\ / A nd | Web: www.OpenFOAM.org | +| \\/ M anipulation | | +\*---------------------------------------------------------------------------*/ +FoamFile +{ + version 2.0; + format ascii; + class dictionary; + object fvSolution; +} +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +solvers +{ + p + { + solver GAMG; + smoother DIC; + tolerance 1e-6; + relTol 0.05; + } + + "(U|sigma)" + { + solver smoothSolver; + smoother symGaussSeidel; + tolerance 1e-6; + relTol 0.1; + } + + pFinal + { + $p; + relTol 0; + } + + "(U|sigma)Final" + { + $U; + relTol 0; + } +} + +PIMPLE +{ + momentumPredictor off; + nOuterCorrectors 15; + nCorrectors 3; + nNonOrthogonalCorrectors 0; + pRefCell 0; + pRefValue 0; + turbOnFinalIterOnly no; +} + +relaxationFactors +{ + equations + { + ".*" 1; + } +} + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/probes b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/probes new file mode 100644 index 0000000000..6bec7ba300 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/probes @@ -0,0 +1,21 @@ +/*--------------------------------*- C++ -*----------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Web: www.OpenFOAM.org + \\/ M anipulation | +------------------------------------------------------------------------------- +Description + Writes out values of fields from cells nearest to specified locations. + +\*---------------------------------------------------------------------------*/ + +#includeEtc "caseDicts/postProcessing/probes/probes.cfg" + +fields (U); +probeLocations +( + (0 1 0) +); + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/residuals b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/residuals new file mode 100644 index 0000000000..af061576bd --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/residuals @@ -0,0 +1,19 @@ +/*--------------------------------*- C++ -*----------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Web: www.OpenFOAM.org + \\/ M anipulation | +------------------------------------------------------------------------------- +Description + For specified fields, writes out the initial residuals for the first + solution of each time step; for non-scalar fields (e.g. vectors), writes + the largest of the residuals for each component (e.g. x, y, z). + +\*---------------------------------------------------------------------------*/ + +#includeEtc "caseDicts/postProcessing/numerical/residuals.cfg" + +fields (p sigma); + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/singleGraph b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/singleGraph new file mode 100644 index 0000000000..560dd1c51a --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/system/singleGraph @@ -0,0 +1,48 @@ +/*--------------------------------*- C++ -*----------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Web: www.OpenFOAM.org + \\/ M anipulation | +------------------------------------------------------------------------------- +Description + Writes graph data for specified fields along a line, specified by start + and end points. + +\*---------------------------------------------------------------------------*/ + + +// Sampling and I/O settings +#includeEtc "caseDicts/postProcessing/graphs/sampleDict.cfg" + +// Override settings here, e.g. +// setConfig { type midPoint; } + +type sets; +libs ("libsampling.so"); + +writeControl writeTime; + +interpolationScheme cellPoint; + +setFormat raw; + +setConfig +{ + type midPoint; // midPoint + axis distance; // x, y, z, xyz +} + +sets +( + line + { + $setConfig; + start (0 0 0); + end (0 1 0); + } +); + +fields (U); + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/Make/files b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/Make/files new file mode 100644 index 0000000000..af1834eae1 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/Make/files @@ -0,0 +1,3 @@ +WatersKing.C + +EXE = $(FOAM_USER_APPBIN)/WatersKing diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/Make/options b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/Make/options new file mode 100644 index 0000000000..780009c876 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/Make/options @@ -0,0 +1,14 @@ +EXE_INC = \ + -I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \ + -I$(LIB_SRC)/TurbulenceModels/incompressible/lnInclude \ + -I$(LIB_SRC)/transportModels \ + -I$(LIB_SRC)/transportModels/incompressible/singlePhaseTransportModel \ + -I$(LIB_SRC)/finiteVolume/lnInclude \ + -I$(LIB_SRC)/meshTools/lnInclude + +EXE_LIBS = \ + -lturbulenceModels \ + -lincompressibleTurbulenceModels \ + -lincompressibleTransportModels \ + -lfiniteVolume \ + -lmeshTools diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/WatersKing.C b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/WatersKing.C new file mode 100644 index 0000000000..4f9ebd896b --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/WatersKing.C @@ -0,0 +1,146 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Copyright (C) 2016 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 . + +Application + WatersKing + +Description + Analytical solution for the start-up planar Poiseuille flow of an + Oldroyd-B fluid. + + References: + \verbatim + Waters, N. D., & King, M. J. (1970). + Unsteady flow of an elasto-viscous liquid. + Rheologica Acta, 9, 345-355. + + Amoreira, L. J., & Oliveira, P. J. (2010). + Comparison of different formulations for the numerical + calculation of unsteady incompressible viscoelastic fluid + flow. Adv. Appl. Math. Mech, 4, 483-502. + \endverbatim + +\*---------------------------------------------------------------------------*/ + +#include "fvCFD.H" +#include "singlePhaseTransportModel.H" +#include "turbulentTransportModel.H" + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +int main(int argc, char *argv[]) +{ + #include "setRootCase.H" + #include "createTime.H" + + // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + + #include "createMesh.H" + #include "createFields.H" + + const scalar h = mesh.bounds().span().y(); + Info<< "Height from centreline to wall = " << h << endl; + + label centrelineID = mesh.boundary().findPatchID("centreline"); + const vector patchToCell = + mesh.boundary()[centrelineID].Cf()[0] + - mesh.C()[mesh.findNearestCell(location)]; + + const scalar y = patchToCell.y()/h; + Info<< "Normalised distance from centreline = " << y << nl << endl; + + const scalar nu0 = nu1 + nu2; + const scalar E = lambda*nu0/(rho*sqr(h)); + const scalar beta = nu2/nu0; + const scalar UInf = K*sqr(h)/3.0/nu0; + + Info<< "Waters and King parameters:" << nl + << "E = " << E << nl + << "beta = " << beta << nl + << "K = " << K << nl + << "UInf = " << UInf << nl << endl; + + label order = 8; + + scalarField ak(order, 0); + scalarField bk(order, 0); + scalarField ck(order, 0); + scalarField B(order, 0); + + forAll(ak, i) + { + scalar k = i + 1; + ak[i] = (2.0*k - 1)/2.0*constant::mathematical::pi*::sqrt(E); + bk[i] = (1.0 + beta*sqr(ak[i]))/2.0; + ck[i] = ::sqrt(mag(sqr(bk[i]) - sqr(ak[i]))); + B[i] = 48*::pow(-1, k)/::pow((2*k - 1)*constant::mathematical::pi, 3)* + ::cos((2*k - 1)*constant::mathematical::pi*y/2); + } + + scalarField A(order, 0); + OFstream file(runTime.path()/"WatersKing.dat"); + const scalar LOGVGREAT = ::log(VGREAT); + while (!runTime.end()) + { + scalar t = runTime.timeOutputValue()/lambda; + forAll(A, i) + { + if (bk[i]*t < LOGVGREAT) + { + if (bk[i] >= ak[i]) + { + A[i] = (bk[i] - sqr(ak[i]))/ck[i]*::sinh(ck[i]*t) + + ::cosh(ck[i]*t); + } + else + { + A[i] = (bk[i] - sqr(ak[i]))/ck[i]*::sin(ck[i]*t) + + ::cos(ck[i]*t); + } + A[i] *= ::exp(-bk[i]*t); + } + else + { + Info<< "Coefficient A[" << order << "] = 0" << endl; + order = i; + Info<< "Resizing A and B to " << order << endl; + A.resize(order); + B.resize(order); + } + } + scalar U = UInf*(1.5*(1 - sqr(y)) + sum(A*B)); + file<< runTime.timeName() << token::TAB << U << endl; + runTime++; + } + + Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s" + << " ClockTime = " << runTime.elapsedClockTime() << " s" + << nl << endl; + + Info<< "End\n" << endl; + + return 0; +} + + +// ************************************************************************* // diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/createFields.H b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/createFields.H new file mode 100644 index 0000000000..4ae80853b4 --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/WatersKing/createFields.H @@ -0,0 +1,64 @@ +Info<< "Reading transportProperties\n" << endl; +IOdictionary transportProperties +( + IOobject + ( + "transportProperties", + runTime.constant(), + mesh, + IOobject::MUST_READ, + IOobject::NO_WRITE + ) +); +const scalar nu2 = + dimensionedScalar + ( + "nu", + dimViscosity, + transportProperties.lookup("nu") + ).value(); + +Info<< "Reading viscoelastic properties\n" << endl; +IOdictionary turbulenceProperties +( + IOobject + ( + "turbulenceProperties", + runTime.constant(), + mesh, + IOobject::MUST_READ, + IOobject::NO_WRITE + ) +); +const dictionary& MaxwellCoeffs = + turbulenceProperties.subDict("laminar").subDict("MaxwellCoeffs"); +const scalar nu1 = readScalar(MaxwellCoeffs.lookup("nuM")); +const scalar lambda = readScalar(MaxwellCoeffs.lookup("lambda")); + +const scalar rho = 1; + +Info<< "Reading pressure gradient\n" << endl; +IOdictionary fvOptions +( + IOobject + ( + "fvOptions", + runTime.constant(), + mesh, + IOobject::MUST_READ, + IOobject::NO_WRITE + ) +); +const dictionary& gradPDict = + fvOptions.subDict("momentumSource").subDict + ( + "vectorSemiImplicitSourceCoeffs" + ).subDict + ( + "injectionRateSuSp" + ); +const scalar K = + Tuple2(gradPDict.lookup("U")).first().x(); + +dictionary probes(IFstream(runTime.system()/"probes")()); +const point location = pointField(probes.lookup("probeLocations"))[0]; diff --git a/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/createGraph b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/createGraph new file mode 100755 index 0000000000..a30402600e --- /dev/null +++ b/tutorials/incompressible/pimpleFoam/laminar/planarPoiseuille/validation/createGraph @@ -0,0 +1,27 @@ +#!/bin/sh + +tail -n +4 ../postProcessing/probes/0/U | \ + tr -s " " | tr -d '(' | cut -d " " -f2-3 > ../Numerical.dat + +if ! which gnuplot > /dev/null 2>&1 +then + echo "gnuplot not found - skipping graph creation" >&2 + exit 1 +fi + +gnuplot<