ENH: add pointPatchField::patchInternalField with three-parameters
- allows extraction without constructing a tmp. The definition is similar to {faPatch,fvPatch} handling.
This commit is contained in:
parent
1a1322a2b0
commit
dde4b12687
@ -237,7 +237,8 @@ void Foam::valuePointPatchField<Type>::operator=
|
||||
const pointPatchField<Type>& ptf
|
||||
)
|
||||
{
|
||||
Field<Type>::operator=(this->patchInternalField());
|
||||
// pointPatchField has no values to copy, assign internal values
|
||||
this->extrapolateInternal();
|
||||
}
|
||||
|
||||
|
||||
@ -277,7 +278,8 @@ void Foam::valuePointPatchField<Type>::operator==
|
||||
const pointPatchField<Type>& ptf
|
||||
)
|
||||
{
|
||||
Field<Type>::operator=(this->patchInternalField());
|
||||
// pointPatchField has no values to copy, assign internal values
|
||||
this->extrapolateInternal();
|
||||
}
|
||||
|
||||
|
||||
|
@ -93,15 +93,16 @@ void Foam::processorCyclicPointPatchField<Type>::initSwapAddSeparated
|
||||
Field<Type>& pField
|
||||
) const
|
||||
{
|
||||
if (Pstream::parRun())
|
||||
if (UPstream::parRun())
|
||||
{
|
||||
// Get internal field into correct order for opposite side. Note use
|
||||
// of member data buffer since using non-blocking. Could be optimised
|
||||
// out if not using non-blocking...
|
||||
sendBuf_ = this->patchInternalField
|
||||
this->patchInternalField
|
||||
(
|
||||
pField,
|
||||
procPatch_.reverseMeshPoints()
|
||||
procPatch_.reverseMeshPoints(),
|
||||
sendBuf_
|
||||
);
|
||||
|
||||
if (commsType == Pstream::commsTypes::nonBlocking)
|
||||
@ -138,7 +139,7 @@ void Foam::processorCyclicPointPatchField<Type>::swapAddSeparated
|
||||
Field<Type>& pField
|
||||
) const
|
||||
{
|
||||
if (Pstream::parRun())
|
||||
if (UPstream::parRun())
|
||||
{
|
||||
// If nonblocking, data is already in receive buffer
|
||||
|
||||
|
@ -6,7 +6,7 @@
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
Copyright (C) 2011-2016 OpenFOAM Foundation
|
||||
Copyright (C) 2020-2022 OpenCFD Ltd.
|
||||
Copyright (C) 2020-2023 OpenCFD Ltd.
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of OpenFOAM.
|
||||
@ -108,6 +108,64 @@ void Foam::pointPatchField<Type>::write(Ostream& os) const
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
template<class Type1>
|
||||
void Foam::pointPatchField<Type>::patchInternalField
|
||||
(
|
||||
const UList<Type1>& internalData,
|
||||
const labelUList& addressing,
|
||||
Field<Type1>& pfld
|
||||
) const
|
||||
{
|
||||
// Check size
|
||||
if (internalData.size() != primitiveField().size())
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Internal field size: " << internalData.size()
|
||||
<< " != mesh size: " << primitiveField().size() << nl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
const label len = this->size();
|
||||
|
||||
pfld.resize_nocopy(len);
|
||||
|
||||
for (label i = 0; i < len; ++i)
|
||||
{
|
||||
pfld[i] = internalData[addressing[i]];
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
template<class Type1>
|
||||
Foam::tmp<Foam::Field<Type1>>
|
||||
Foam::pointPatchField<Type>::patchInternalField
|
||||
(
|
||||
const UList<Type1>& internalData,
|
||||
const labelUList& addressing
|
||||
) const
|
||||
{
|
||||
auto tpfld = tmp<Field<Type1>>::New();
|
||||
this->patchInternalField(internalData, addressing, tpfld.ref());
|
||||
return tpfld;
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
template<class Type1>
|
||||
Foam::tmp<Foam::Field<Type1>>
|
||||
Foam::pointPatchField<Type>::patchInternalField
|
||||
(
|
||||
const UList<Type1>& internalData
|
||||
) const
|
||||
{
|
||||
auto tpfld = tmp<Field<Type1>>::New();
|
||||
this->patchInternalField(internalData, patch().meshPoints(), tpfld.ref());
|
||||
return tpfld;
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
Foam::tmp<Foam::Field<Type>>
|
||||
Foam::pointPatchField<Type>::patchInternalField() const
|
||||
@ -116,41 +174,6 @@ Foam::pointPatchField<Type>::patchInternalField() const
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
template<class Type1>
|
||||
Foam::tmp<Foam::Field<Type1>>
|
||||
Foam::pointPatchField<Type>::patchInternalField
|
||||
(
|
||||
const Field<Type1>& iF,
|
||||
const labelUList& meshPoints
|
||||
) const
|
||||
{
|
||||
// Check size
|
||||
if (iF.size() != primitiveField().size())
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "given internal field does not correspond to the mesh. "
|
||||
<< "Field size: " << iF.size()
|
||||
<< " mesh size: " << primitiveField().size()
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
return tmp<Field<Type1>>::New(iF, meshPoints);
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
template<class Type1>
|
||||
Foam::tmp<Foam::Field<Type1>>
|
||||
Foam::pointPatchField<Type>::patchInternalField
|
||||
(
|
||||
const Field<Type1>& iF
|
||||
) const
|
||||
{
|
||||
return patchInternalField(iF, patch().meshPoints());
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
template<class Type1>
|
||||
void Foam::pointPatchField<Type>::addToInternalField
|
||||
@ -163,18 +186,16 @@ void Foam::pointPatchField<Type>::addToInternalField
|
||||
if (iF.size() != primitiveField().size())
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "given internal field does not correspond to the mesh. "
|
||||
<< "Field size: " << iF.size()
|
||||
<< " mesh size: " << primitiveField().size()
|
||||
<< "Internal field size: " << iF.size()
|
||||
<< " != mesh size: " << primitiveField().size() << nl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
if (pF.size() != size())
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "given patch field does not correspond to the mesh. "
|
||||
<< "Field size: " << pF.size()
|
||||
<< " mesh size: " << size()
|
||||
<< "Patch field size: " << pF.size()
|
||||
<< " != patch size: " << size() << nl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
@ -201,18 +222,16 @@ void Foam::pointPatchField<Type>::addToInternalField
|
||||
if (iF.size() != primitiveField().size())
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "given internal field does not correspond to the mesh. "
|
||||
<< "Field size: " << iF.size()
|
||||
<< " mesh size: " << primitiveField().size()
|
||||
<< "Internal field size: " << iF.size()
|
||||
<< " != mesh size: " << primitiveField().size() << nl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
if (pF.size() != size())
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "given patch field does not correspond to the mesh. "
|
||||
<< "Field size: " << pF.size()
|
||||
<< " mesh size: " << size()
|
||||
<< "Patch field size: " << pF.size()
|
||||
<< " != patch size: " << size() << nl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
@ -240,18 +259,16 @@ void Foam::pointPatchField<Type>::setInInternalField
|
||||
if (iF.size() != primitiveField().size())
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "given internal field does not correspond to the mesh. "
|
||||
<< "Field size: " << iF.size()
|
||||
<< " mesh size: " << primitiveField().size()
|
||||
<< "Internal field size: " << iF.size()
|
||||
<< " != mesh size: " << primitiveField().size() << nl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
if (pF.size() != meshPoints.size())
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "given patch field does not correspond to the meshPoints. "
|
||||
<< "Field size: " << pF.size()
|
||||
<< " meshPoints size: " << size()
|
||||
<< "Patch field size: " << pF.size()
|
||||
<< " != meshPoints size: " << meshPoints.size() << nl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
|
@ -431,26 +431,43 @@ public:
|
||||
|
||||
// Evaluation Functions
|
||||
|
||||
//- Return field created from appropriate internal field values
|
||||
tmp<Field<Type>> patchInternalField() const;
|
||||
|
||||
//- Return field created from appropriate internal field values
|
||||
//- given internal field reference
|
||||
//- Extract field using specified addressing
|
||||
// \param internalData The internal field to extract from
|
||||
// \param addressing Addressing (mesh-points) into internal field
|
||||
// \param [out] pfld The extracted patch field.
|
||||
// It is always resized according to the patch size(),
|
||||
// which can be smaller than the addressing size
|
||||
template<class Type1>
|
||||
tmp<Field<Type1>> patchInternalField
|
||||
void patchInternalField
|
||||
(
|
||||
const Field<Type1>& iF
|
||||
const UList<Type1>& internalData,
|
||||
const labelUList& addressing,
|
||||
Field<Type1>& pfld
|
||||
) const;
|
||||
|
||||
//- Return field created from selected internal field values
|
||||
//- given internal field reference
|
||||
// \param internalData The internal field to extract from
|
||||
// \param addressing Addressing (mesh-points) into internal field
|
||||
template<class Type1>
|
||||
tmp<Field<Type1>> patchInternalField
|
||||
(
|
||||
const Field<Type1>& iF,
|
||||
const labelUList& meshPoints
|
||||
const UList<Type1>& internalData,
|
||||
const labelUList& addressing
|
||||
) const;
|
||||
|
||||
//- Return field created from appropriate internal field values
|
||||
//- given internal field reference
|
||||
template<class Type1>
|
||||
tmp<Field<Type1>> patchInternalField
|
||||
(
|
||||
const UList<Type1>& internalData
|
||||
) const;
|
||||
|
||||
//- Return field created from appropriate internal field values
|
||||
tmp<Field<Type>> patchInternalField() const;
|
||||
|
||||
|
||||
//- Given the internal field and a patch field,
|
||||
//- add the patch field to the internal field
|
||||
template<class Type1>
|
||||
|
Loading…
Reference in New Issue
Block a user