/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2015 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 "basicXiSubXiEq.H" #include "zeroGradientFvPatchFields.H" #include "addToRunTimeSelectionTable.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { namespace XiEqModels { defineTypeNameAndDebug(basicSubGrid, 0); addToRunTimeSelectionTable(XiEqModel, basicSubGrid, dictionary); } } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::XiEqModels::basicSubGrid::basicSubGrid ( const dictionary& XiEqProperties, const psiuReactionThermo& thermo, const compressible::RASModel& turbulence, const volScalarField& Su ) : XiEqModel(XiEqProperties, thermo, turbulence, Su), B_ ( IOobject ( "B", Su.mesh().facesInstance(), Su.mesh(), IOobject::MUST_READ, IOobject::NO_WRITE ), Su.mesh() ), XiEqModel_(XiEqModel::New(XiEqModelCoeffs_, thermo, turbulence, Su)) {} // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // Foam::XiEqModels::basicSubGrid::~basicSubGrid() {} // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // Foam::tmp Foam::XiEqModels::basicSubGrid::XiEq() const { const fvMesh& mesh = Su_.mesh(); const volVectorField& U = mesh.lookupObject("U"); const volScalarField& Nv = mesh.lookupObject("Nv"); const volSymmTensorField& nsv = mesh.lookupObject("nsv"); volScalarField magU(mag(U)); volVectorField Uhat ( U/(mag(U) + dimensionedScalar("Usmall", U.dimensions(), 1e-4)) ); const scalarField Cw = pow(mesh.V(), 2.0/3.0); tmp tN ( new volScalarField ( IOobject ( "tN", mesh.time().constant(), mesh, IOobject::NO_READ, IOobject::NO_WRITE ), mesh, dimensionedScalar("zero", Nv.dimensions(), 0.0), zeroGradientFvPatchVectorField::typeName ) ); volScalarField& N = tN(); N.internalField() = Nv.internalField()*Cw; tmp tns ( new volSymmTensorField ( IOobject ( "tns", U.mesh().time().timeName(), U.mesh(), IOobject::NO_READ, IOobject::NO_WRITE ), U.mesh(), dimensionedSymmTensor ( "zero", nsv.dimensions(), pTraits::zero ), zeroGradientFvPatchSymmTensorField::typeName ) ); volSymmTensorField& ns = tns(); ns.internalField() = nsv.internalField()*Cw; volScalarField n(max(N - (Uhat & ns & Uhat), scalar(1e-4))); volScalarField b((Uhat & B_ & Uhat)/sqrt(n)); volScalarField up(sqrt((2.0/3.0)*turbulence_.k())); volScalarField XiSubEq ( scalar(1) + max(2.2*sqrt(b), min(0.34*magU/up*sqrt(b), scalar(1.6))) * min(n, scalar(1)) ); return (XiSubEq*XiEqModel_->XiEq()); } bool Foam::XiEqModels::basicSubGrid::read(const dictionary& XiEqProperties) { XiEqModel::read(XiEqProperties); return XiEqModel_->read(XiEqModelCoeffs_); } // ************************************************************************* //