A single transformer object is now maintained within cyclic patches and returned from a single virtual functions massively simplifying the interface and allowing for further rationalisation of the calculation of the transformation.
194 lines
4.8 KiB
C
194 lines
4.8 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration | Website: https://openfoam.org
|
|
\\ / A nd | Copyright (C) 2011-2020 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 <http://www.gnu.org/licenses/>.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "jumpCyclicFvPatchField.H"
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
template<class Type>
|
|
Foam::jumpCyclicFvPatchField<Type>::jumpCyclicFvPatchField
|
|
(
|
|
const fvPatch& p,
|
|
const DimensionedField<Type, volMesh>& iF
|
|
)
|
|
:
|
|
cyclicFvPatchField<Type>(p, iF)
|
|
{}
|
|
|
|
|
|
template<class Type>
|
|
Foam::jumpCyclicFvPatchField<Type>::jumpCyclicFvPatchField
|
|
(
|
|
const jumpCyclicFvPatchField<Type>& ptf,
|
|
const fvPatch& p,
|
|
const DimensionedField<Type, volMesh>& iF,
|
|
const fvPatchFieldMapper& mapper
|
|
)
|
|
:
|
|
cyclicFvPatchField<Type>(ptf, p, iF, mapper)
|
|
{}
|
|
|
|
|
|
template<class Type>
|
|
Foam::jumpCyclicFvPatchField<Type>::jumpCyclicFvPatchField
|
|
(
|
|
const fvPatch& p,
|
|
const DimensionedField<Type, volMesh>& iF,
|
|
const dictionary& dict
|
|
)
|
|
:
|
|
cyclicFvPatchField<Type>(p, iF, dict)
|
|
{
|
|
// Call this evaluation in derived classes
|
|
// this->evaluate(Pstream::commsTypes::blocking);
|
|
}
|
|
|
|
|
|
template<class Type>
|
|
Foam::jumpCyclicFvPatchField<Type>::jumpCyclicFvPatchField
|
|
(
|
|
const jumpCyclicFvPatchField<Type>& ptf
|
|
)
|
|
:
|
|
cyclicFvPatchField<Type>(ptf)
|
|
{}
|
|
|
|
|
|
template<class Type>
|
|
Foam::jumpCyclicFvPatchField<Type>::jumpCyclicFvPatchField
|
|
(
|
|
const jumpCyclicFvPatchField<Type>& ptf,
|
|
const DimensionedField<Type, volMesh>& iF
|
|
)
|
|
:
|
|
cyclicFvPatchField<Type>(ptf, iF)
|
|
{}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
template<class Type>
|
|
Foam::tmp<Foam::Field<Type>>
|
|
Foam::jumpCyclicFvPatchField<Type>::patchNeighbourField() const
|
|
{
|
|
const Field<Type>& iField = this->primitiveField();
|
|
const labelUList& nbrFaceCells =
|
|
this->cyclicPatch().neighbFvPatch().faceCells();
|
|
|
|
tmp<Field<Type>> tpnf(new Field<Type>(this->size()));
|
|
Field<Type>& pnf = tpnf.ref();
|
|
|
|
Field<Type> jf(this->jump());
|
|
if (!this->cyclicPatch().owner())
|
|
{
|
|
jf *= -1.0;
|
|
}
|
|
|
|
if (this->doTransform())
|
|
{
|
|
forAll(*this, facei)
|
|
{
|
|
pnf[facei] =
|
|
this->transform().transform(iField[nbrFaceCells[facei]])
|
|
- jf[facei];
|
|
}
|
|
}
|
|
else
|
|
{
|
|
forAll(*this, facei)
|
|
{
|
|
pnf[facei] = iField[nbrFaceCells[facei]] - jf[facei];
|
|
}
|
|
}
|
|
|
|
return tpnf;
|
|
}
|
|
|
|
|
|
template<class Type>
|
|
void Foam::jumpCyclicFvPatchField<Type>::updateInterfaceMatrix
|
|
(
|
|
scalarField& result,
|
|
const scalarField& psiInternal,
|
|
const scalarField& coeffs,
|
|
const direction cmpt,
|
|
const Pstream::commsTypes
|
|
) const
|
|
{
|
|
NotImplemented;
|
|
}
|
|
|
|
|
|
template<class Type>
|
|
void Foam::jumpCyclicFvPatchField<Type>::updateInterfaceMatrix
|
|
(
|
|
Field<Type>& result,
|
|
const Field<Type>& psiInternal,
|
|
const scalarField& coeffs,
|
|
const Pstream::commsTypes
|
|
) const
|
|
{
|
|
Field<Type> pnf(this->size());
|
|
|
|
const labelUList& nbrFaceCells =
|
|
this->cyclicPatch().neighbFvPatch().faceCells();
|
|
|
|
// only apply jump to original field
|
|
if (&psiInternal == &this->primitiveField())
|
|
{
|
|
Field<Type> jf(this->jump());
|
|
|
|
if (!this->cyclicPatch().owner())
|
|
{
|
|
jf *= -1.0;
|
|
}
|
|
|
|
forAll(*this, facei)
|
|
{
|
|
pnf[facei] = psiInternal[nbrFaceCells[facei]] - jf[facei];
|
|
}
|
|
}
|
|
else
|
|
{
|
|
forAll(*this, facei)
|
|
{
|
|
pnf[facei] = psiInternal[nbrFaceCells[facei]];
|
|
}
|
|
}
|
|
|
|
// Transform according to the transformation tensors
|
|
this->transformCoupleField(pnf);
|
|
|
|
// Multiply the field by coefficients and add into the result
|
|
const labelUList& faceCells = this->cyclicPatch().faceCells();
|
|
forAll(faceCells, elemI)
|
|
{
|
|
result[faceCells[elemI]] -= coeffs[elemI]*pnf[elemI];
|
|
}
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|