/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2016 OpenFOAM Foundation Copyright (C) 2023 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 "fixedGradientFvPatchField.H" #include "dictionary.H" // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * // template bool Foam::fixedGradientFvPatchField::readGradientEntry ( const dictionary& dict, IOobjectOption::readOption readOpt ) { if (!IOobjectOption::isAnyRead(readOpt)) return false; const auto& p = fvPatchFieldBase::patch(); const auto* eptr = dict.findEntry("gradient", keyType::LITERAL); if (eptr) { gradient_.assign(*eptr, p.size()); return true; } if (IOobjectOption::isReadRequired(readOpt)) { FatalIOErrorInFunction(dict) << "Required entry 'gradient' : missing for patch " << p.name() << " in dictionary " << dict.relativeName() << nl << exit(FatalIOError); } return false; } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // template Foam::fixedGradientFvPatchField::fixedGradientFvPatchField ( const fvPatch& p, const DimensionedField& iF ) : fvPatchField(p, iF), gradient_(p.size(), Zero) {} template Foam::fixedGradientFvPatchField::fixedGradientFvPatchField ( const fvPatch& p, const DimensionedField& iF, const dictionary& dict, IOobjectOption::readOption requireGrad ) : fvPatchField(p, iF, dict, IOobjectOption::NO_READ), gradient_(p.size()) { if (readGradientEntry(dict, requireGrad)) { evaluate(); } else { // Not read (eg, optional and missing): // - treat as zero-gradient, do not evaluate fvPatchField::extrapolateInternal(); gradient_ = Zero; } } template Foam::fixedGradientFvPatchField::fixedGradientFvPatchField ( const fixedGradientFvPatchField& ptf, const fvPatch& p, const DimensionedField& iF, const fvPatchFieldMapper& mapper ) : fvPatchField(ptf, p, iF, mapper), gradient_(ptf.gradient_, mapper) { if (notNull(iF) && mapper.hasUnmapped()) { WarningInFunction << "On field " << iF.name() << " patch " << p.name() << " patchField " << this->type() << " : mapper does not map all values." << nl << " To avoid this warning fully specify the mapping in derived" << " patch fields." << endl; } } template Foam::fixedGradientFvPatchField::fixedGradientFvPatchField ( const fixedGradientFvPatchField& ptf ) : fvPatchField(ptf), gradient_(ptf.gradient_) {} template Foam::fixedGradientFvPatchField::fixedGradientFvPatchField ( const fixedGradientFvPatchField& ptf, const DimensionedField& iF ) : fvPatchField(ptf, iF), gradient_(ptf.gradient_) {} // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // template void Foam::fixedGradientFvPatchField::autoMap ( const fvPatchFieldMapper& m ) { fvPatchField::autoMap(m); gradient_.autoMap(m); } template void Foam::fixedGradientFvPatchField::rmap ( const fvPatchField& ptf, const labelList& addr ) { fvPatchField::rmap(ptf, addr); const fixedGradientFvPatchField& fgptf = refCast>(ptf); gradient_.rmap(fgptf.gradient_, addr); } template void Foam::fixedGradientFvPatchField::evaluate(const Pstream::commsTypes) { if (!this->updated()) { this->updateCoeffs(); } Field::operator= ( this->patchInternalField() + gradient_/this->patch().deltaCoeffs() ); fvPatchField::evaluate(); } template Foam::tmp> Foam::fixedGradientFvPatchField::valueInternalCoeffs ( const tmp& ) const { return tmp>::New(this->size(), pTraits::one); } template Foam::tmp> Foam::fixedGradientFvPatchField::valueBoundaryCoeffs ( const tmp& ) const { return gradient()/this->patch().deltaCoeffs(); } template Foam::tmp> Foam::fixedGradientFvPatchField::gradientInternalCoeffs() const { return tmp>::New(this->size(), Zero); } template Foam::tmp> Foam::fixedGradientFvPatchField::gradientBoundaryCoeffs() const { return gradient(); } template void Foam::fixedGradientFvPatchField::write(Ostream& os) const { fvPatchField::write(os); gradient_.writeEntry("gradient", os); } // ************************************************************************* //