patchToPatch Class Referenceabstract

Class to generate coupling geometry between two primitive patches. More...

Inheritance diagram for patchToPatch:
Collaboration diagram for patchToPatch:

Public Member Functions

 TypeName ("patchToPatch")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, patchToPatch, bool,(const bool reverse),(reverse))
 Declare runtime constructor selection table. More...
 
 patchToPatch (const bool reverse)
 Construct from components. More...
 
 patchToPatch (const patchToPatch &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~patchToPatch ()
 Destructor. More...
 
bool reverse () const
 Flag to indicate that the two patches are co-directional and. More...
 
label singleProcess () const
 Index of the processor holding all faces of the patchToPatch,. More...
 
bool isSingleProcess () const
 Is this intersection on a single process? More...
 
PackedBoolList srcCoupled () const
 Return a list indicating which source faces are coupled. More...
 
PackedBoolList tgtCoupled () const
 Return a list indicating which target faces are coupled. More...
 
List< List< remote > > srcTgtProcFaces () const
 For each source face, the coupled target procs and faces. More...
 
List< List< remote > > tgtSrcProcFaces () const
 For each target face, the coupled source procs and faces. More...
 
template<class Type >
tmp< Field< Type > > srcToTgt (const Field< Type > &srcFld) const
 Interpolate a source patch field to the target with no left. More...
 
template<class Type >
tmp< Field< Type > > srcToTgt (const Field< Type > &srcFld, const Field< Type > &leftOverTgtFld) const
 Interpolate a source patch field to the target with left over. More...
 
template<class Type >
tmp< Field< Type > > tgtToSrc (const Field< Type > &tgtFld) const
 Interpolate a target patch field to the source with no left. More...
 
template<class Type >
tmp< Field< Type > > tgtToSrc (const Field< Type > &tgtFld, const Field< Type > &leftOverSrcFld) const
 Interpolate a target patch field to the source with left. More...
 
void update (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const transformer &tgtToSrc=NullObjectRef< transformer >())
 Update addressing and weights for the given patches. More...
 
void update (const primitivePatch &srcPatch, const vectorField &srcPointNormals, const primitivePatch &tgtPatch, const transformer &tgtToSrc=NullObjectRef< transformer >())
 Update addressing and weights for the given patches. More...
 
void operator= (const patchToPatch &)=delete
 Disallow default bitwise assignment. More...
 
template<class Type >
Foam::tmp< Foam::Field< Type > > srcToTgt (const Field< Type > &srcFld) const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > srcToTgt (const Field< Type > &srcFld, const Field< Type > &leftOverTgtFld) const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > tgtToSrc (const Field< Type > &tgtFld) const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > tgtToSrc (const Field< Type > &tgtFld, const Field< Type > &leftOverSrcFld) const
 

Static Public Member Functions

static autoPtr< patchToPatchNew (const word &patchToPatchType, const bool reverse)
 Select from name. More...
 

Protected Member Functions

virtual treeBoundBox srcBox (const face &srcFace, const pointField &srcPoints, const vectorField &srcPointNormals) const =0
 Get the bound box for a source face. More...
 
treeBoundBox srcBox (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const label srcFacei) const
 Get the bound box for a source face. More...
 
treeBoundBox srcBox (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0) const
 Get the bound box for the source patch. More...
 
treeBoundBox tgtBox (const primitiveOldTimePatch &tgtPatch, const label tgtFacei) const
 Get the bound box for a target face. More...
 
treeBoundBox tgtBox (const primitiveOldTimePatch &tgtPatch) const
 Get the bound box for the target patch. More...
 
virtual bool intersectFaces (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const label srcFacei, const label tgtFacei)=0
 Intersect two faces. More...
 
bool findOrIntersectFaces (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const label srcFacei, const label tgtFacei)
 Intersect two faces. More...
 
label intersectPatchQueue (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const bool isSrc, const DynamicList< labelPair > &queue, labelList &faceComplete, DynamicList< labelPair > &otherQueue, const labelList &otherFaceComplete, boolList &otherFaceQueued, boolList &otherFaceVisited)
 Intersect a queue of source-target face pairs. Update completion. More...
 
void intersectPatches (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch)
 Intersect the patches. More...
 
labelListList tgtPatchSendFaces (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch) const
 Determine which target faces need to be sent to the source. More...
 
virtual void initialise (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch)
 Initialisation. More...
 
virtual labelList finaliseLocal (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch)
 Finalise the intersection locally. Trims the local target patch. More...
 
virtual void distributeSrc (const primitiveOldTimePatch &srcPatch)
 Distribute the source patch so that any information needed by. More...
 
virtual void rDistributeTgt (const primitiveOldTimePatch &tgtPatch)
 Send data that resulted from an intersection between the source. More...
 
virtual label finalise (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const transformer &tgtToSrc)
 Finalising. More...
 
virtual tmpNrc< List< DynamicList< scalar > > > srcWeights () const =0
 For each source face, the coupled target weights. More...
 
virtual tmpNrc< List< DynamicList< scalar > > > tgtWeights () const =0
 For each target face, the coupled source weights. More...
 

Static Protected Member Functions

static List< remotedistributePatch (const distributionMap &map, const primitiveOldTimePatch &patch, autoPtr< PrimitiveOldTimePatch< faceList, pointField >> &localPatchPtr)
 Distribute a patch given its distribution map. More...
 

Protected Attributes

const bool reverse_
 Flag to indicate that the two patches are co-directional and. More...
 
label singleProcess_
 Index of the processor holding all faces of the patchToPatch, or -1. More...
 
List< DynamicList< label > > srcLocalTgtFaces_
 For each source face, the coupled local target faces. More...
 
List< DynamicList< label > > tgtLocalSrcFaces_
 For each target face, the coupled local source faces. More...
 
autoPtr< distributionMapsrcMapPtr_
 Map from source patch faces to target-local source patch faces. More...
 
autoPtr< distributionMaptgtMapPtr_
 Map from target patch faces to source-local target patch faces. More...
 
autoPtr< List< remote > > localSrcProcFacesPtr_
 When running in parallel, a map from local source face index to. More...
 
autoPtr< List< remote > > localTgtProcFacesPtr_
 When running in parallel, a map from local target face index to. More...
 

Detailed Description

Class to generate coupling geometry between two primitive patches.

Source files

Definition at line 55 of file patchToPatch.H.

Constructor & Destructor Documentation

◆ patchToPatch() [1/2]

patchToPatch ( const bool  reverse)

Construct from components.

Definition at line 764 of file patchToPatch.C.

◆ patchToPatch() [2/2]

patchToPatch ( const patchToPatch )
delete

Disallow default bitwise copy construction.

◆ ~patchToPatch()

~patchToPatch ( )
virtual

Destructor.

Definition at line 779 of file patchToPatch.C.

Member Function Documentation

◆ srcBox() [1/3]

virtual treeBoundBox srcBox ( const face srcFace,
const pointField srcPoints,
const vectorField srcPointNormals 
) const
protectedpure virtual

Get the bound box for a source face.

Implemented in nearby.

Referenced by patchToPatch::srcBox(), and patchToPatch::tgtPatchSendFaces().

Here is the caller graph for this function:

◆ srcBox() [2/3]

Foam::treeBoundBox srcBox ( const primitiveOldTimePatch srcPatch,
const vectorField srcPointNormals,
const vectorField srcPointNormals0,
const label  srcFacei 
) const
protected

◆ srcBox() [3/3]

Foam::treeBoundBox srcBox ( const primitiveOldTimePatch srcPatch,
const vectorField srcPointNormals,
const vectorField srcPointNormals0 
) const
protected

Get the bound box for the source patch.

Definition at line 118 of file patchToPatch.C.

References Foam::combine(), forAll, boundBox::inflate(), and treeBoundBox::invertedBox.

Here is the call graph for this function:

◆ tgtBox() [1/2]

Foam::treeBoundBox tgtBox ( const primitiveOldTimePatch tgtPatch,
const label  tgtFacei 
) const
protected

Get the bound box for a target face.

Definition at line 143 of file patchToPatch.C.

References Foam::combine(), PrimitiveOldTimePatch< FaceList, PointField >::has0(), PrimitivePatch< FaceList, PointField >::localFaces(), and PrimitivePatch< FaceList, PointField >::localPoints().

Referenced by patchToPatch::tgtPatchSendFaces().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ tgtBox() [2/2]

Foam::treeBoundBox tgtBox ( const primitiveOldTimePatch tgtPatch) const
protected

Get the bound box for the target patch.

Definition at line 163 of file patchToPatch.C.

References Foam::combine(), forAll, boundBox::inflate(), and treeBoundBox::invertedBox.

Here is the call graph for this function:

◆ intersectFaces()

virtual bool intersectFaces ( const primitiveOldTimePatch srcPatch,
const vectorField srcPointNormals,
const vectorField srcPointNormals0,
const primitiveOldTimePatch tgtPatch,
const label  srcFacei,
const label  tgtFacei 
)
protectedpure virtual

Intersect two faces.

Implemented in nearest, and nearby.

◆ findOrIntersectFaces()

bool findOrIntersectFaces ( const primitiveOldTimePatch srcPatch,
const vectorField srcPointNormals,
const vectorField srcPointNormals0,
const primitiveOldTimePatch tgtPatch,
const label  srcFacei,
const label  tgtFacei 
)
protected

Intersect two faces.

Definition at line 181 of file patchToPatch.C.

References forAll.

◆ intersectPatchQueue()

Foam::label intersectPatchQueue ( const primitiveOldTimePatch srcPatch,
const vectorField srcPointNormals,
const vectorField srcPointNormals0,
const primitiveOldTimePatch tgtPatch,
const bool  isSrc,
const DynamicList< labelPair > &  queue,
labelList faceComplete,
DynamicList< labelPair > &  otherQueue,
const labelList otherFaceComplete,
boolList otherFaceQueued,
boolList otherFaceVisited 
)
protected

◆ intersectPatches()

◆ tgtPatchSendFaces()

Foam::labelListList tgtPatchSendFaces ( const primitiveOldTimePatch srcPatch,
const vectorField srcPointNormals,
const vectorField srcPointNormals0,
const primitiveOldTimePatch tgtPatch 
) const
protected

Determine which target faces need to be sent to the source.

This is done before intersection. Bound boxes are used to estimate what faces will intersect.

Definition at line 32 of file patchToPatchParallelOps.C.

References List< T >::append(), forAll, Pstream::gatherList(), UPstream::myProcNo(), UPstream::nProcs(), Pstream::scatterList(), patchToPatch::srcBox(), patchToPatch::tgtBox(), and List< T >::transfer().

Here is the call graph for this function:

◆ distributePatch()

◆ initialise()

void initialise ( const primitiveOldTimePatch srcPatch,
const vectorField srcPointNormals,
const vectorField srcPointNormals0,
const primitiveOldTimePatch tgtPatch 
)
protectedvirtual

Initialisation.

Reimplemented in nearest, and nearby.

Definition at line 639 of file patchToPatch.C.

References forAll.

Referenced by nearby::initialise().

Here is the caller graph for this function:

◆ finaliseLocal()

Foam::labelList finaliseLocal ( const primitiveOldTimePatch srcPatch,
const vectorField srcPointNormals,
const vectorField srcPointNormals0,
const primitiveOldTimePatch tgtPatch 
)
protectedvirtual

Finalise the intersection locally. Trims the local target patch.

so that only parts that actually intersect the source remain. Returns new-to-old map for trimming target-associated data.

Reimplemented in nearest.

Definition at line 661 of file patchToPatch.C.

References forAll, and Foam::patchToPatchTools::trimDistributionMap().

Referenced by nearest::finaliseLocal().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ distributeSrc()

void distributeSrc ( const primitiveOldTimePatch srcPatch)
protectedvirtual

Distribute the source patch so that any information needed by.

the target is locally available

Definition at line 709 of file patchToPatch.C.

References Foam::patchToPatchTools::distributeAddressing().

Here is the call graph for this function:

◆ rDistributeTgt()

void rDistributeTgt ( const primitiveOldTimePatch tgtPatch)
protectedvirtual

Send data that resulted from an intersection between the source.

patch and a distributed source-local-target patch back to the original target processes

Reimplemented in nearest.

Definition at line 724 of file patchToPatch.C.

References Foam::patchToPatchTools::rDistributeTgtAddressing().

Referenced by nearest::rDistributeTgt().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalise()

Foam::label finalise ( const primitiveOldTimePatch srcPatch,
const vectorField srcPointNormals,
const vectorField srcPointNormals0,
const primitiveOldTimePatch tgtPatch,
const transformer tgtToSrc 
)
protectedvirtual

Finalising.

Reimplemented in nearest.

Definition at line 736 of file patchToPatch.C.

References forAll, and Foam::reduce().

Referenced by nearest::finalise().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ srcWeights()

virtual tmpNrc<List<DynamicList<scalar> > > srcWeights ( ) const
protectedpure virtual

For each source face, the coupled target weights.

Implemented in nearest.

◆ tgtWeights()

virtual tmpNrc<List<DynamicList<scalar> > > tgtWeights ( ) const
protectedpure virtual

For each target face, the coupled source weights.

Implemented in nearest.

Referenced by patchToPatch::srcToTgt().

Here is the caller graph for this function:

◆ TypeName()

TypeName ( "patchToPatch"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
patchToPatch  ,
bool  ,
(const bool reverse ,
(reverse  
)

Declare runtime constructor selection table.

◆ New()

Foam::autoPtr< Foam::patchToPatch > New ( const word patchToPatchType,
const bool  reverse 
)
static

Select from name.

Definition at line 785 of file patchToPatch.C.

References Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, and Foam::reverse().

Here is the call graph for this function:

◆ reverse()

bool reverse ( ) const
inline

Flag to indicate that the two patches are co-directional and.

that the orientation of one should therefore be reversed

Definition at line 31 of file patchToPatchI.H.

References patchToPatch::reverse_.

◆ singleProcess()

Foam::label singleProcess ( ) const
inline

Index of the processor holding all faces of the patchToPatch,.

or -1 if spread across multiple processors

Definition at line 37 of file patchToPatchI.H.

◆ isSingleProcess()

bool isSingleProcess ( ) const
inline

Is this intersection on a single process?

Definition at line 43 of file patchToPatchI.H.

◆ srcCoupled()

Foam::PackedBoolList srcCoupled ( ) const

Return a list indicating which source faces are coupled.

Definition at line 810 of file patchToPatch.C.

References PackedList< nBits >::empty(), and forAll.

Here is the call graph for this function:

◆ tgtCoupled()

Foam::PackedBoolList tgtCoupled ( ) const

Return a list indicating which target faces are coupled.

Definition at line 821 of file patchToPatch.C.

References PackedList< nBits >::empty(), and forAll.

Here is the call graph for this function:

◆ srcTgtProcFaces()

Foam::List< Foam::List< Foam::remote > > srcTgtProcFaces ( ) const

For each source face, the coupled target procs and faces.

Definition at line 833 of file patchToPatch.C.

References Foam::patchToPatchTools::localToRemote().

Here is the call graph for this function:

◆ tgtSrcProcFaces()

Foam::List< Foam::List< Foam::remote > > tgtSrcProcFaces ( ) const

For each target face, the coupled source procs and faces.

Definition at line 847 of file patchToPatch.C.

References Foam::patchToPatchTools::localToRemote().

Here is the call graph for this function:

◆ srcToTgt() [1/4]

tmp<Field<Type> > srcToTgt ( const Field< Type > &  srcFld) const

Interpolate a source patch field to the target with no left.

over values specified. If the interpolation weight sum is less than one for a face then they will be normalised. If the interpolation weight sum is zero for a face then that face's value will be NaN.

◆ srcToTgt() [2/4]

tmp<Field<Type> > srcToTgt ( const Field< Type > &  srcFld,
const Field< Type > &  leftOverTgtFld 
) const

Interpolate a source patch field to the target with left over.

values specified. If the interpolation weight sum is less than one for a face then the average will include the left over value multiplied by one minus the weight sum.

◆ tgtToSrc() [1/4]

tmp<Field<Type> > tgtToSrc ( const Field< Type > &  tgtFld) const

Interpolate a target patch field to the source with no left.

over values specified. As the corresponding srcToTgt.

◆ tgtToSrc() [2/4]

tmp<Field<Type> > tgtToSrc ( const Field< Type > &  tgtFld,
const Field< Type > &  leftOverSrcFld 
) const

Interpolate a target patch field to the source with left.

over values specified. As the corresponding srcToTgt.

◆ update() [1/2]

◆ update() [2/2]

void update ( const primitivePatch srcPatch,
const vectorField srcPointNormals,
const primitivePatch tgtPatch,
const transformer tgtToSrc = NullObjectRef<transformer>() 
)

Update addressing and weights for the given patches.

Definition at line 1067 of file patchToPatch.C.

◆ operator=()

void operator= ( const patchToPatch )
delete

Disallow default bitwise assignment.

◆ srcToTgt() [3/4]

Foam::tmp<Foam::Field<Type> > srcToTgt ( const Field< Type > &  srcFld) const

Definition at line 33 of file patchToPatchTemplates.C.

References Foam::patchToPatchTools::interpolate(), patchToPatch::srcMapPtr_, patchToPatch::tgtLocalSrcFaces_, and patchToPatch::tgtWeights().

Here is the call graph for this function:

◆ srcToTgt() [4/4]

Foam::tmp<Foam::Field<Type> > srcToTgt ( const Field< Type > &  srcFld,
const Field< Type > &  leftOverTgtFld 
) const

Definition at line 47 of file patchToPatchTemplates.C.

References Foam::patchToPatchTools::interpolate().

Here is the call graph for this function:

◆ tgtToSrc() [3/4]

Foam::tmp<Foam::Field<Type> > tgtToSrc ( const Field< Type > &  tgtFld) const

Definition at line 67 of file patchToPatchTemplates.C.

References Foam::patchToPatchTools::interpolate().

Here is the call graph for this function:

◆ tgtToSrc() [4/4]

Foam::tmp<Foam::Field<Type> > tgtToSrc ( const Field< Type > &  tgtFld,
const Field< Type > &  leftOverSrcFld 
) const

Definition at line 81 of file patchToPatchTemplates.C.

References Foam::patchToPatchTools::interpolate().

Here is the call graph for this function:

Member Data Documentation

◆ reverse_

const bool reverse_
protected

Flag to indicate that the two patches are co-directional and.

that the orientation of one should therefore be reversed

Definition at line 63 of file patchToPatch.H.

Referenced by patchToPatch::reverse().

◆ singleProcess_

label singleProcess_
protected

Index of the processor holding all faces of the patchToPatch, or -1.

if spread across multiple processors

Definition at line 67 of file patchToPatch.H.

◆ srcLocalTgtFaces_

List<DynamicList<label> > srcLocalTgtFaces_
protected

For each source face, the coupled local target faces.

Definition at line 70 of file patchToPatch.H.

Referenced by nearest::intersectFaces().

◆ tgtLocalSrcFaces_

List<DynamicList<label> > tgtLocalSrcFaces_
protected

For each target face, the coupled local source faces.

Definition at line 73 of file patchToPatch.H.

Referenced by nearest::intersectFaces(), and patchToPatch::srcToTgt().

◆ srcMapPtr_

autoPtr<distributionMap> srcMapPtr_
protected

Map from source patch faces to target-local source patch faces.

Definition at line 76 of file patchToPatch.H.

Referenced by patchToPatch::srcToTgt().

◆ tgtMapPtr_

autoPtr<distributionMap> tgtMapPtr_
protected

Map from target patch faces to source-local target patch faces.

Definition at line 79 of file patchToPatch.H.

◆ localSrcProcFacesPtr_

autoPtr<List<remote> > localSrcProcFacesPtr_
protected

When running in parallel, a map from local source face index to.

source processor and face index

Definition at line 83 of file patchToPatch.H.

◆ localTgtProcFacesPtr_

autoPtr<List<remote> > localTgtProcFacesPtr_
protected

When running in parallel, a map from local target face index to.

target processor and face index

Definition at line 87 of file patchToPatch.H.


The documentation for this class was generated from the following files: