/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2017 OpenFOAM Foundation
Copyright (C) 2019-2025 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 "cyclicAMIPolyPatch.H"
#include "mapDistributeBase.H"
#include "AMIInterpolation.H"
#include "fvMatrix.H"
#include "volFields.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template
Foam::cyclicAMIFvPatchField::cyclicAMIFvPatchField
(
const fvPatch& p,
const DimensionedField& iF
)
:
cyclicAMILduInterfaceField(),
coupledFvPatchField(p, iF),
cyclicAMIPatch_(refCast(p)),
patchNeighbourFieldPtr_(nullptr)
{}
template
Foam::cyclicAMIFvPatchField::cyclicAMIFvPatchField
(
const fvPatch& p,
const DimensionedField& iF,
const dictionary& dict
)
:
cyclicAMILduInterfaceField(),
coupledFvPatchField(p, iF, dict, IOobjectOption::NO_READ),
cyclicAMIPatch_(refCast(p, dict)),
patchNeighbourFieldPtr_(nullptr)
{
if (!isA(p))
{
FatalIOErrorInFunction(dict)
<< "\n patch type '" << p.type()
<< "' not constraint type '" << typeName << "'"
<< "\n for patch " << p.name()
<< " of field " << this->internalField().name()
<< " in file " << this->internalField().objectPath()
<< exit(FatalIOError);
}
if (cacheNeighbourField())
{
// Handle neighbour value first, before any evaluate()
const auto* hasNeighbValue =
dict.findEntry("neighbourValue", keyType::LITERAL);
if (hasNeighbValue)
{
patchNeighbourFieldPtr_.reset
(
new Field(*hasNeighbValue, p.size())
);
}
}
// Use 'value' supplied, or evaluate (if coupled) or set to internal field
if (!this->readValueEntry(dict))
{
if (this->coupled())
{
// Tricky: avoid call to evaluate without call to initEvaluate.
// For now just disable the localConsistency to make it use the
// old logic (ultimately calls the fully self contained
// patchNeighbourField)
const auto oldConsistency = FieldBase::localBoundaryConsistency(0);
this->evaluate(UPstream::commsTypes::nonBlocking);
FieldBase::localBoundaryConsistency(oldConsistency);
}
else
{
this->extrapolateInternal(); // Zero-gradient patch values
}
}
}
template
Foam::cyclicAMIFvPatchField::cyclicAMIFvPatchField
(
const cyclicAMIFvPatchField& ptf,
const fvPatch& p,
const DimensionedField& iF,
const fvPatchFieldMapper& mapper
)
:
cyclicAMILduInterfaceField(),
coupledFvPatchField(ptf, p, iF, mapper),
cyclicAMIPatch_(refCast(p)),
patchNeighbourFieldPtr_(nullptr)
{
//if (ptf.patchNeighbourFieldPtr_ && cacheNeighbourField())
//{
// patchNeighbourFieldPtr_.reset
// (
// new Field(ptf.patchNeighbourFieldPtr_(), mapper)
// );
//}
if (!isA(this->patch()))
{
FatalErrorInFunction
<< "\n patch type '" << p.type()
<< "' not constraint type '" << typeName << "'"
<< "\n for patch " << p.name()
<< " of field " << this->internalField().name()
<< " in file " << this->internalField().objectPath()
<< exit(FatalError);
}
if (debug && !ptf.all_ready())
{
FatalErrorInFunction
<< "Outstanding request(s) on patch " << cyclicAMIPatch_.name()
<< abort(FatalError);
}
}
template
Foam::cyclicAMIFvPatchField::cyclicAMIFvPatchField
(
const cyclicAMIFvPatchField& ptf
)
:
cyclicAMILduInterfaceField(),
coupledFvPatchField(ptf),
cyclicAMIPatch_(ptf.cyclicAMIPatch_),
patchNeighbourFieldPtr_(nullptr)
{
if (debug && !ptf.all_ready())
{
FatalErrorInFunction
<< "Outstanding request(s) on patch " << cyclicAMIPatch_.name()
<< abort(FatalError);
}
}
template
Foam::cyclicAMIFvPatchField::cyclicAMIFvPatchField
(
const cyclicAMIFvPatchField& ptf,
const DimensionedField& iF
)
:
cyclicAMILduInterfaceField(),
coupledFvPatchField(ptf, iF),
cyclicAMIPatch_(ptf.cyclicAMIPatch_),
patchNeighbourFieldPtr_(nullptr)
{
if (debug && !ptf.all_ready())
{
FatalErrorInFunction
<< "Outstanding request(s) on patch " << cyclicAMIPatch_.name()
<< abort(FatalError);
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template
bool Foam::cyclicAMIFvPatchField::all_ready() const
{
int done = 0;
if
(
UPstream::finishedRequests
(
recvRequests_.start(),
recvRequests_.size()
)
)
{
recvRequests_.clear();
++done;
}
if
(
UPstream::finishedRequests
(
sendRequests_.start(),
sendRequests_.size()
)
)
{
sendRequests_.clear();
++done;
}
return (done == 2);
}
template
bool Foam::cyclicAMIFvPatchField::ready() const
{
if
(
UPstream::finishedRequests
(
recvRequests_.start(),
recvRequests_.size()
)
)
{
recvRequests_.clear();
if
(
UPstream::finishedRequests
(
sendRequests_.start(),
sendRequests_.size()
)
)
{
sendRequests_.clear();
}
return true;
}
return false;
}
template
void Foam::cyclicAMIFvPatchField::autoMap
(
const fvPatchFieldMapper& mapper
)
{
coupledFvPatchField::autoMap(mapper);
patchNeighbourFieldPtr_.reset(nullptr);
}
template
void Foam::cyclicAMIFvPatchField::rmap
(
const fvPatchField& ptf,
const labelList& addr
)
{
coupledFvPatchField::rmap(ptf, addr);
patchNeighbourFieldPtr_.reset(nullptr);
}
template
Foam::tmp>
Foam::cyclicAMIFvPatchField::getNeighbourField
(
const UList& internalData
) const
{
// By pass polyPatch to get nbrId. Instead use cyclicAMIFvPatch virtual
// neighbPatch()
const auto& neighbPatch = cyclicAMIPatch_.neighbPatch();
const labelUList& nbrFaceCells = neighbPatch.faceCells();
Field pnf(internalData, nbrFaceCells);
Field defaultValues;
if (cyclicAMIPatch_.applyLowWeightCorrection())
{
defaultValues = Field(internalData, cyclicAMIPatch_.faceCells());
}
tmp> tpnf = cyclicAMIPatch_.interpolate(pnf, defaultValues);
if (doTransform())
{
transform(tpnf.ref(), forwardT(), tpnf());
}
return tpnf;
}
template
bool Foam::cyclicAMIFvPatchField::cacheNeighbourField()
{
return (FieldBase::localBoundaryConsistency() != 0);
}
template
Foam::tmp>
Foam::cyclicAMIFvPatchField::getPatchNeighbourField
(
bool checkCommunicator
) const
{
const auto& AMI = this->ownerAMI();
if
(
AMI.distributed() && cacheNeighbourField()
&& (!checkCommunicator || AMI.comm() != -1)
)
{
if (!this->ready())
{
FatalErrorInFunction
<< "Outstanding recv request(s) on patch "
<< cyclicAMIPatch_.name()
<< " field " << this->internalField().name()
<< abort(FatalError);
}
const auto& fvp = this->patch();
if
(
patchNeighbourFieldPtr_
&& !fvp.boundaryMesh().mesh().upToDatePoints(this->internalField())
)
{
//DebugPout
// << "cyclicAMIFvPatchField::patchNeighbourField() :"
// << " field:" << this->internalField().name()
// << " patch:" << fvp.name()
// << " CLEARING patchNeighbourField"
// << endl;
patchNeighbourFieldPtr_.reset(nullptr);
}
// Initialise if not done in construct-from-dictionary
if (!patchNeighbourFieldPtr_)
{
//DebugPout
// << "cyclicAMIFvPatchField::patchNeighbourField() :"
// << " field:" << this->internalField().name()
// << " patch:" << fvp.name()
// << " caching patchNeighbourField"
// << endl;
// Do interpolation and store result
patchNeighbourFieldPtr_.reset
(
getNeighbourField(this->primitiveField()).ptr()
);
}
else
{
// Have cached value. Check
//if (debug)
//{
// tmp> tpnf
// (
// getNeighbourField(this->primitiveField())
// );
// if (tpnf() != patchNeighbourFieldPtr_())
// {
// FatalErrorInFunction
// << "On field " << this->internalField().name()
// << " patch " << fvp.name() << endl
// << "Cached patchNeighbourField :"
// << flatOutput(patchNeighbourFieldPtr_()) << endl
// << "Calculated patchNeighbourField:"
// << flatOutput(tpnf()) << exit(FatalError);
// }
//}
}
return patchNeighbourFieldPtr_();
}
else
{
// Do interpolation
return getNeighbourField(this->primitiveField());
}
}
template
Foam::tmp>
Foam::cyclicAMIFvPatchField::patchNeighbourField() const
{
return this->getPatchNeighbourField(true); // checkCommunicator = true
}
template
void Foam::cyclicAMIFvPatchField::patchNeighbourField
(
UList& pnf
) const
{
// checkCommunicator = false
auto tpnf = this->getPatchNeighbourField(false);
pnf.deepCopy(tpnf());
}
template
const Foam::cyclicAMIFvPatchField&
Foam::cyclicAMIFvPatchField::neighbourPatchField() const
{
const auto& fld =
static_cast&>
(
this->primitiveField()
);
return refCast>
(
fld.boundaryField()[cyclicAMIPatch_.neighbPatchID()]
);
}
template
void Foam::cyclicAMIFvPatchField::initEvaluate
(
const Pstream::commsTypes commsType
)
{
if (!this->updated())
{
this->updateCoeffs();
}
const auto& AMI = this->ownerAMI();
if (AMI.distributed() && cacheNeighbourField() && AMI.comm() != -1)
{
//DebugPout
// << "*** cyclicAMIFvPatchField::initEvaluate() :"
// << " field:" << this->internalField().name()
// << " patch:" << this->patch().name()
// << " sending patchNeighbourField"
// << endl;
if (commsType != UPstream::commsTypes::nonBlocking)
{
// Invalidate old field - or flag as fatal?
patchNeighbourFieldPtr_.reset(nullptr);
return;
}
// Start sending
// Bypass polyPatch to get nbrId.
// - use cyclicACMIFvPatch::neighbPatch() virtual instead
const cyclicAMIFvPatch& neighbPatch = cyclicAMIPatch_.neighbPatch();
const labelUList& nbrFaceCells = neighbPatch.faceCells();
const Field pnf(this->primitiveField(), nbrFaceCells);
const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
// Assert that all receives are known to have finished
if (!recvRequests_.empty())
{
FatalErrorInFunction
<< "Outstanding recv request(s) on patch "
<< cyclicAMIPatch_.name()
<< " field " << this->internalField().name()
<< abort(FatalError);
}
// Assume that sends are also OK
sendRequests_.clear();
cpp.initInterpolate
(
pnf,
sendRequests_,
sendBufs_,
recvRequests_,
recvBufs_
);
}
}
template
void Foam::cyclicAMIFvPatchField::evaluate
(
const Pstream::commsTypes commsType
)
{
if (!this->updated())
{
this->updateCoeffs();
}
const auto& AMI = this->ownerAMI();
if (AMI.distributed() && cacheNeighbourField() && AMI.comm() != -1)
{
// Calculate patchNeighbourField
if (commsType != UPstream::commsTypes::nonBlocking)
{
FatalErrorInFunction
<< "Can only evaluate distributed AMI with nonBlocking"
<< exit(FatalError);
}
patchNeighbourFieldPtr_.reset(nullptr);
const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
Field defaultValues;
if (AMI.applyLowWeightCorrection())
{
defaultValues = this->patchInternalField();
}
//DebugPout
// << "*** cyclicAMIFvPatchField::evaluate() :"
// << " field:" << this->internalField().name()
// << " patch:" << this->patch().name()
// << " receiving&caching patchNeighbourField"
// << endl;
patchNeighbourFieldPtr_.reset
(
cpp.interpolate
(
Field::null(), // Not used for distributed
recvRequests_,
recvBufs_,
defaultValues
).ptr()
);
// Receive requests all handled by last function call
recvRequests_.clear();
if (doTransform())
{
// In-place transform
auto& pnf = *patchNeighbourFieldPtr_;
transform(pnf, forwardT(), pnf);
}
}
// Use patchNeighbourField() and patchInternalField() to obtain face value
coupledFvPatchField::evaluate(commsType);
}
template
void Foam::cyclicAMIFvPatchField::initInterfaceMatrixUpdate
(
solveScalarField& result,
const bool add,
const lduAddressing& lduAddr,
const label patchId,
const solveScalarField& psiInternal,
const scalarField& coeffs,
const direction cmpt,
const Pstream::commsTypes commsType
) const
{
const auto& AMI = this->ownerAMI();
if (AMI.distributed() && AMI.comm() != -1)
{
// Start sending
if (commsType != UPstream::commsTypes::nonBlocking)
{
FatalErrorInFunction
<< "Can only evaluate distributed AMI with nonBlocking"
<< exit(FatalError);
}
const labelUList& nbrFaceCells =
lduAddr.patchAddr(cyclicAMIPatch_.neighbPatchID());
solveScalarField pnf(psiInternal, nbrFaceCells);
// Transform according to the transformation tensors
transformCoupleField(pnf, cmpt);
const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
// Assert that all receives are known to have finished
if (!recvRequests_.empty())
{
FatalErrorInFunction
<< "Outstanding recv request(s) on patch "
<< cyclicAMIPatch_.name()
<< " field " << this->internalField().name()
<< abort(FatalError);
}
// Assume that sends are also OK
sendRequests_.clear();
cpp.initInterpolate
(
pnf,
sendRequests_,
scalarSendBufs_,
recvRequests_,
scalarRecvBufs_
);
}
this->updatedMatrix(false);
}
template
void Foam::cyclicAMIFvPatchField::updateInterfaceMatrix
(
solveScalarField& result,
const bool add,
const lduAddressing& lduAddr,
const label patchId,
const solveScalarField& psiInternal,
const scalarField& coeffs,
const direction cmpt,
const Pstream::commsTypes commsType
) const
{
//DebugPout<< "cyclicAMIFvPatchField::updateInterfaceMatrix() :"
// << " field:" << this->internalField().name()
// << " patch:" << this->patch().name()
// << endl;
const labelUList& faceCells = lduAddr.patchAddr(patchId);
const auto& AMI = this->ownerAMI();
solveScalarField pnf;
if (AMI.distributed() && AMI.comm() != -1)
{
if (commsType != UPstream::commsTypes::nonBlocking)
{
FatalErrorInFunction
<< "Can only evaluate distributed AMI with nonBlocking"
<< exit(FatalError);
}
solveScalarField defaultValues;
if (AMI.applyLowWeightCorrection())
{
defaultValues = solveScalarField(psiInternal, faceCells);
}
const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
pnf =
cpp.interpolate
(
solveScalarField::null(), // Not used for distributed
recvRequests_,
scalarRecvBufs_,
defaultValues
);
// Receive requests all handled by last function call
recvRequests_.clear();
}
else
{
solveScalarField defaultValues;
if (cyclicAMIPatch_.applyLowWeightCorrection())
{
defaultValues = solveScalarField(psiInternal, faceCells);
}
const labelUList& nbrFaceCells =
lduAddr.patchAddr(cyclicAMIPatch_.neighbPatchID());
pnf = solveScalarField(psiInternal, nbrFaceCells);
// Transform according to the transformation tensors
transformCoupleField(pnf, cmpt);
pnf = cyclicAMIPatch_.interpolate(pnf, defaultValues);
}
// Multiply the field by coefficients and add into the result
this->addToInternalField(result, !add, faceCells, coeffs, pnf);
this->updatedMatrix(true);
}
template
void Foam::cyclicAMIFvPatchField::initInterfaceMatrixUpdate
(
Field& result,
const bool add,
const lduAddressing& lduAddr,
const label patchId,
const Field& psiInternal,
const scalarField& coeffs,
const Pstream::commsTypes commsType
) const
{
const auto& AMI = this->ownerAMI();
if (AMI.distributed() && AMI.comm() != -1)
{
if (commsType != UPstream::commsTypes::nonBlocking)
{
FatalErrorInFunction
<< "Can only evaluate distributed AMI with nonBlocking"
<< exit(FatalError);
}
const labelUList& nbrFaceCells =
lduAddr.patchAddr(cyclicAMIPatch_.neighbPatchID());
Field pnf(psiInternal, nbrFaceCells);
// Transform according to the transformation tensors
transformCoupleField(pnf);
const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
// Assert that all receives are known to have finished
if (!recvRequests_.empty())
{
FatalErrorInFunction
<< "Outstanding recv request(s) on patch "
<< cyclicAMIPatch_.name()
<< " field " << this->internalField().name()
<< abort(FatalError);
}
// Assume that sends are also OK
sendRequests_.clear();
cpp.initInterpolate
(
pnf,
sendRequests_,
sendBufs_,
recvRequests_,
recvBufs_
);
}
this->updatedMatrix(false);
}
template
void Foam::cyclicAMIFvPatchField::updateInterfaceMatrix
(
Field& result,
const bool add,
const lduAddressing& lduAddr,
const label patchId,
const Field& psiInternal,
const scalarField& coeffs,
const Pstream::commsTypes commsType
) const
{
//DebugPout<< "cyclicAMIFvPatchField::updateInterfaceMatrix() :"
// << " field:" << this->internalField().name()
// << " patch:" << this->patch().name()
// << endl;
const labelUList& faceCells = lduAddr.patchAddr(patchId);
const auto& AMI = this->ownerAMI();
Field pnf;
if (AMI.distributed() && AMI.comm() != -1)
{
if (commsType != UPstream::commsTypes::nonBlocking)
{
FatalErrorInFunction
<< "Can only evaluate distributed AMI with nonBlocking"
<< exit(FatalError);
}
const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
Field defaultValues;
if (AMI.applyLowWeightCorrection())
{
defaultValues = Field(psiInternal, faceCells);
}
pnf =
cpp.interpolate
(
Field::null(), // Not used for distributed
recvRequests_,
recvBufs_,
defaultValues
);
// Receive requests all handled by last function call
recvRequests_.clear();
}
else
{
const labelUList& nbrFaceCells =
lduAddr.patchAddr(cyclicAMIPatch_.neighbPatchID());
pnf = Field(psiInternal, nbrFaceCells);
// Transform according to the transformation tensors
transformCoupleField(pnf);
Field defaultValues;
if (cyclicAMIPatch_.applyLowWeightCorrection())
{
defaultValues = Field(psiInternal, faceCells);
}
pnf = cyclicAMIPatch_.interpolate(pnf, defaultValues);
}
// Multiply the field by coefficients and add into the result
this->addToInternalField(result, !add, faceCells, coeffs, pnf);
this->updatedMatrix(true);
}
template
void Foam::cyclicAMIFvPatchField::manipulateMatrix
(
fvMatrix& matrix,
const label mat,
const direction cmpt
)
{
if (this->cyclicAMIPatch().owner())
{
const label index = this->patch().index();
const label globalPatchID =
matrix.lduMeshAssembly().patchLocalToGlobalMap()[mat][index];
const Field intCoeffsCmpt
(
matrix.internalCoeffs()[globalPatchID].component(cmpt)
);
const Field boundCoeffsCmpt
(
matrix.boundaryCoeffs()[globalPatchID].component(cmpt)
);
tmp> tintCoeffs(coeffs(matrix, intCoeffsCmpt, mat));
tmp> tbndCoeffs(coeffs(matrix, boundCoeffsCmpt, mat));
const Field& intCoeffs = tintCoeffs.ref();
const Field& bndCoeffs = tbndCoeffs.ref();
const labelUList& u = matrix.lduAddr().upperAddr();
const labelUList& l = matrix.lduAddr().lowerAddr();
label subFaceI = 0;
const labelList& faceMap =
matrix.lduMeshAssembly().faceBoundMap()[mat][index];
forAll (faceMap, j)
{
label globalFaceI = faceMap[j];
const scalar boundCorr = -bndCoeffs[subFaceI];
const scalar intCorr = -intCoeffs[subFaceI];
matrix.upper()[globalFaceI] += boundCorr;
matrix.diag()[u[globalFaceI]] -= intCorr;
matrix.diag()[l[globalFaceI]] -= boundCorr;
if (matrix.asymmetric())
{
matrix.lower()[globalFaceI] += intCorr;
}
subFaceI++;
}
// Set internalCoeffs and boundaryCoeffs in the assembly matrix
// on clyclicAMI patches to be used in the individual matrix by
// matrix.flux()
if (matrix.psi(mat).mesh().fluxRequired(this->internalField().name()))
{
matrix.internalCoeffs().set
(
globalPatchID, intCoeffs*pTraits::one
);
matrix.boundaryCoeffs().set
(
globalPatchID, bndCoeffs*pTraits::one
);
const label nbrPathID =
cyclicAMIPatch_.cyclicAMIPatch().neighbPatchID();
const label nbrGlobalPatchID =
matrix.lduMeshAssembly().patchLocalToGlobalMap()
[mat][nbrPathID];
matrix.internalCoeffs().set
(
nbrGlobalPatchID, intCoeffs*pTraits::one
);
matrix.boundaryCoeffs().set
(
nbrGlobalPatchID, bndCoeffs*pTraits::one
);
}
}
}
template
Foam::tmp>
Foam::cyclicAMIFvPatchField::coeffs
(
fvMatrix& matrix,
const Field& coeffs,
const label mat
) const
{
const label index(this->patch().index());
const label nSubFaces
(
matrix.lduMeshAssembly().cellBoundMap()[mat][index].size()
);
auto tmapCoeffs = tmp>::New(nSubFaces, Zero);
auto& mapCoeffs = tmapCoeffs.ref();
const scalarListList& srcWeight =
cyclicAMIPatch_.cyclicAMIPatch().AMI().srcWeights();
label subFaceI = 0;
forAll(*this, faceI)
{
const scalarList& w = srcWeight[faceI];
for(label i=0; i
template
void Foam::cyclicAMIFvPatchField::collectStencilData
(
const refPtr& mapPtr,
const labelListList& stencil,
const Type2& data,
List& expandedData
)
{
expandedData.resize_nocopy(stencil.size());
if (mapPtr)
{
Type2 work(data);
mapPtr().distribute(work);
forAll(stencil, facei)
{
const auto& slots = stencil[facei];
expandedData[facei].push_back
(
UIndirectList(work, slots)
);
}
}
else
{
forAll(stencil, facei)
{
const auto& slots = stencil[facei];
expandedData[facei].push_back
(
UIndirectList(data, slots)
);
}
}
}
template
void Foam::cyclicAMIFvPatchField::write(Ostream& os) const
{
fvPatchField::write(os);
fvPatchField::writeValueEntry(os);
if (patchNeighbourFieldPtr_)
{
patchNeighbourFieldPtr_->writeEntry("neighbourValue", os);
}
}
// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
template
void Foam::cyclicAMIFvPatchField::operator=
(
const fvPatchField& ptf
)
{
fvPatchField::operator=(ptf);
//Pout<< "cyclicAMIFvPatchField::operator= :"
// << " field:" << this->internalField().name()
// << " patch:" << this->patch().name()
// << " copying from field:" << ptf.internalField().name()
// << endl;
const auto* cycPtr = isA>(ptf);
if
(
cycPtr
&& cycPtr->patchNeighbourFieldPtr_
&& cycPtr->patchNeighbourFieldPtr_->size() == this->size()
)
{
if (!patchNeighbourFieldPtr_)
{
patchNeighbourFieldPtr_ = autoPtr>::New();
}
// Copy values
*patchNeighbourFieldPtr_ = *(cycPtr->patchNeighbourFieldPtr_);
}
else
{
patchNeighbourFieldPtr_.reset(nullptr);
}
}
template
void Foam::cyclicAMIFvPatchField::operator==
(
const fvPatchField& ptf
)
{
fvPatchField::operator==(ptf);
//Pout<< "cyclicAMIFvPatchField::operator== :"
// << " field:" << this->internalField().name()
// << " patch:" << this->patch().name()
// << " copying from field:" << ptf.internalField().name()
// << endl;
const auto* cycPtr = isA>(ptf);
if
(
cycPtr
&& cycPtr->patchNeighbourFieldPtr_
&& cycPtr->patchNeighbourFieldPtr_->size() == this->size()
)
{
if (!patchNeighbourFieldPtr_)
{
patchNeighbourFieldPtr_ = autoPtr>::New();
}
// Copy values
*patchNeighbourFieldPtr_ = *(cycPtr->patchNeighbourFieldPtr_);
}
else
{
patchNeighbourFieldPtr_.reset(nullptr);
}
}
// ************************************************************************* //