Files
openfoam/src/OpenFOAM/fields/pointPatchFields/basic/value/valuePointPatchField.C
Mark Olesen 3a6e427409 ENH: simplify dictionary search for value/refValue in BCs
- in expressions BCs in particular, there is various logic handling
  for if value/refValue/refGradient etc are found or not.

  Handle the lookups as findEntry and branch to use Field assign
  or other handling, depending on its existence.

STYLE: use wordList instead of wordRes for copy/filter dictionary
2022-10-04 15:51:26 +02:00

254 lines
5.5 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / 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) 2019-2022 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "valuePointPatchField.H"
#include "pointPatchFieldMapper.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
template<class Type>
Foam::valuePointPatchField<Type>::valuePointPatchField
(
const pointPatch& p,
const DimensionedField<Type, pointMesh>& iF
)
:
pointPatchField<Type>(p, iF),
Field<Type>(p.size())
{}
template<class Type>
Foam::valuePointPatchField<Type>::valuePointPatchField
(
const pointPatch& p,
const DimensionedField<Type, pointMesh>& iF,
const dictionary& dict,
const bool valueRequired
)
:
pointPatchField<Type>(p, iF, dict),
Field<Type>(p.size())
{
const auto* hasValue = dict.findEntry("value", keyType::LITERAL);
if (hasValue)
{
Field<Type>::assign(*hasValue, p.size());
}
else if (!valueRequired)
{
Field<Type>::operator=(Zero);
}
else
{
FatalIOErrorInFunction(dict)
<< "Essential entry 'value' missing on patch "
<< p.name() << endl
<< exit(FatalIOError);
}
}
template<class Type>
Foam::valuePointPatchField<Type>::valuePointPatchField
(
const valuePointPatchField<Type>& ptf,
const pointPatch& p,
const DimensionedField<Type, pointMesh>& iF,
const pointPatchFieldMapper& mapper
)
:
pointPatchField<Type>(ptf, p, iF, mapper),
Field<Type>(ptf, mapper)
{}
template<class Type>
Foam::valuePointPatchField<Type>::valuePointPatchField
(
const valuePointPatchField<Type>& ptf,
const DimensionedField<Type, pointMesh>& iF
)
:
pointPatchField<Type>(ptf, iF),
Field<Type>(ptf)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
void Foam::valuePointPatchField<Type>::autoMap
(
const pointPatchFieldMapper& m
)
{
Field<Type>::autoMap(m);
}
template<class Type>
void Foam::valuePointPatchField<Type>::rmap
(
const pointPatchField<Type>& ptf,
const labelList& addr
)
{
Field<Type>::rmap
(
refCast<const valuePointPatchField<Type>>
(
ptf
),
addr
);
}
template<class Type>
void Foam::valuePointPatchField<Type>::updateCoeffs()
{
if (this->updated())
{
return;
}
// Get internal field to insert values into
Field<Type>& iF = const_cast<Field<Type>&>(this->primitiveField());
this->setInInternalField(iF, *this);
pointPatchField<Type>::updateCoeffs();
}
template<class Type>
void Foam::valuePointPatchField<Type>::evaluate(const Pstream::commsTypes)
{
// Get internal field to insert values into
Field<Type>& iF = const_cast<Field<Type>&>(this->primitiveField());
this->setInInternalField(iF, *this);
pointPatchField<Type>::evaluate();
}
template<class Type>
void Foam::valuePointPatchField<Type>::write(Ostream& os) const
{
pointPatchField<Type>::write(os);
this->writeEntry("value", os);
}
// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
template<class Type>
void Foam::valuePointPatchField<Type>::operator=
(
const valuePointPatchField<Type>& ptf
)
{
Field<Type>::operator=(ptf);
}
template<class Type>
void Foam::valuePointPatchField<Type>::operator=
(
const pointPatchField<Type>& ptf
)
{
Field<Type>::operator=(this->patchInternalField());
}
template<class Type>
void Foam::valuePointPatchField<Type>::operator=
(
const Field<Type>& tf
)
{
Field<Type>::operator=(tf);
}
template<class Type>
void Foam::valuePointPatchField<Type>::operator=
(
const Type& t
)
{
Field<Type>::operator=(t);
}
template<class Type>
void Foam::valuePointPatchField<Type>::operator==
(
const valuePointPatchField<Type>& ptf
)
{
Field<Type>::operator=(ptf);
}
template<class Type>
void Foam::valuePointPatchField<Type>::operator==
(
const pointPatchField<Type>& ptf
)
{
Field<Type>::operator=(this->patchInternalField());
}
template<class Type>
void Foam::valuePointPatchField<Type>::operator==
(
const Field<Type>& tf
)
{
Field<Type>::operator=(tf);
}
template<class Type>
void Foam::valuePointPatchField<Type>::operator==
(
const Type& t
)
{
Field<Type>::operator=(t);
}
// ************************************************************************* //