/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2013-2016 OpenFOAM Foundation \\/ M anipulation | Copyright (C) 2016 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 . \*---------------------------------------------------------------------------*/ #include "gaussConvectionScheme.H" #include "boundedConvectionScheme.H" #include "blendedSchemeBase.H" #include "fvcCellReduce.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // template void Foam::functionObjects::blendingFactor::calcBlendingFactor ( const GeometricField& field, const typename fv::convectionScheme& cs ) { if (!isA>(cs)) { WarningInFunction << "Scheme for field " << field.name() << " is not a " << fv::gaussConvectionScheme::typeName << " scheme. Not calculating " << resultName_ << endl; return; } const fv::gaussConvectionScheme& gcs = refCast>(cs); const surfaceInterpolationScheme& interpScheme = gcs.interpScheme(); if (!isA>(interpScheme)) { WarningInFunction << interpScheme.type() << " is not a blended scheme" << ". Not calculating " << resultName_ << endl; return; } // Retrieve the face-based blending factor const blendedSchemeBase& blendedScheme = refCast>(interpScheme); const surfaceScalarField factorf(blendedScheme.blendingFactor(field)); // Convert into vol field whose values represent the local face minima // Note: // - factor applied to 1st scheme, and (1-factor) to 2nd scheme // - not using the store(...) mechanism due to need to correct BCs volScalarField& indicator = const_cast ( lookupObject(resultName_) ); indicator = 1 - fvc::cellReduce(factorf, minEqOp(), GREAT); indicator.correctBoundaryConditions(); } template bool Foam::functionObjects::blendingFactor::calcScheme() { typedef GeometricField FieldType; if (!foundObject(fieldName_)) { return false; } const FieldType& field = lookupObject(fieldName_); const word divScheme("div(" + phiName_ + ',' + fieldName_ + ')'); ITstream& its = mesh_.divScheme(divScheme); const surfaceScalarField& phi = lookupObject(phiName_); tmp> tcs = fv::convectionScheme::New(mesh_, phi, its); if (isA>(tcs())) { const fv::boundedConvectionScheme& bcs = refCast>(tcs()); calcBlendingFactor(field, bcs.scheme()); } else { const fv::gaussConvectionScheme& gcs = refCast>(tcs()); calcBlendingFactor(field, gcs); } return true; } // ************************************************************************* //