patchToPatch Class Referenceabstract

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

Inheritance diagram for patchToPatch:
Collaboration diagram for patchToPatch:

Classes

struct  procFace
 Structure to conveniently store processor and face indices. More...
 

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...
 
List< List< procFace > > srcTgtProcFaces () const
 For each source face, the coupled target procs and faces. More...
 
List< List< procFace > > tgtSrcProcFaces () const
 For each target face, the coupled source procs and faces. More...
 
virtual tmpNrc< List< DynamicList< scalar > > > srcWeights (const primitivePatch &srcPatch) const =0
 For each source face, the coupled target weights. More...
 
virtual tmpNrc< List< DynamicList< scalar > > > tgtWeights (const primitivePatch &tgtPatch) const =0
 For each target face, the coupled source weights. 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...
 

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...
 
void calcSingleProcess (const primitiveOldTimePatch &srcPatch, const primitiveOldTimePatch &tgtPatch)
 Determine whether or not the intersection of the given 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...
 
labelListList srcPatchSendFaces () const
 Determine which source faces need to be sent to the target. More...
 
distributionMap patchDistributionMap (labelListList &&sendFaces) const
 Create a distribution map from a list-list of faces to be sent. More...
 
void distributePatch (const distributionMap &map, List< procFace > &localProcFaces) const
 Distribute a patch given its distribution map. Just generate. More...
 
PrimitiveOldTimePatch< faceList, pointFielddistributePatch (const distributionMap &map, const primitiveOldTimePatch &patch, List< procFace > &localProcFaces) const
 Distribute a patch given its distribution map. More...
 
virtual void initialise (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch)
 Initialisation. More...
 
virtual tmpNrc< PrimitiveOldTimePatch< faceList, pointField > > distributeTgt (const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, distributionMap &tgtMap)
 Distribute the target patch so that enough is locally available. More...
 
virtual tmpNrc< PrimitiveOldTimePatch< faceList, pointField > > distributeSrc (const primitiveOldTimePatch &srcPatch, distributionMap &srcMap)
 Distribute the source patch so that everything the target. More...
 
virtual void rDistributeTgt (const primitiveOldTimePatch &tgtPatch, const distributionMap &tgtMap)
 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...
 

Static Protected Member Functions

template<class SubListA , class SubListB >
static void transferListList (List< SubListA > &a, List< SubListB > &b)
 Transfer list-list b into list-list a. More...
 
template<class Type >
static void rDistributeListList (const label size, const distributionMap &map, List< List< Type >> &data)
 Reverse distribute a list-list given the map. More...
 
template<class Type >
static void rDistributeListList (const label size, const distributionMap &map, List< DynamicList< Type >> &data)
 Reverse distribute a dynamic list-list given the map. More...
 
static List< List< procFace > > localFacesToProcFaces (const List< DynamicList< label >> &localFaces, const List< procFace > &map=NullObjectRef< List< procFace >>())
 Map local faces to proc faces. More...
 
static List< DynamicList< label > > procFacesToLocalFaces (const List< List< procFace >> &procFaces, const HashTable< label, procFace, Hash< procFace >> &map)
 Map proc faces to local faces. 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...
 
autoPtr< List< procFace > > localSrcProcFacesPtr_
 When running in parallel, a map from local source face index to. More...
 
autoPtr< List< procFace > > localTgtProcFacesPtr_
 When running in parallel, a map from local target face index to. 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...
 

Detailed Description

Class to generate coupling geometry between two primitive patches.

Source files

Definition at line 53 of file patchToPatch.H.

Constructor & Destructor Documentation

◆ patchToPatch() [1/2]

patchToPatch ( const bool  reverse)

Construct from components.

Definition at line 831 of file patchToPatch.C.

◆ patchToPatch() [2/2]

patchToPatch ( const patchToPatch )
delete

Disallow default bitwise copy construction.

◆ ~patchToPatch()

~patchToPatch ( )
virtual

Destructor.

Definition at line 844 of file patchToPatch.C.

References patchToPatch::New().

Here is the call graph for this function:

Member Function Documentation

◆ transferListList()

void transferListList ( List< SubListA > &  a,
List< SubListB > &  b 
)
inlinestaticprotected

Transfer list-list b into list-list a.

Definition at line 32 of file patchToPatchI.H.

References forAll, patchToPatch::rDistributeListList(), List< T >::resize(), List< T >::size(), and List< T >::transfer().

Here is the call graph for this function:

◆ rDistributeListList() [1/2]

void rDistributeListList ( const label  size,
const distributionMap map,
List< List< Type >> &  data 
)
inlinestaticprotected

Reverse distribute a list-list given the map.

Definition at line 47 of file patchToPatchI.H.

References distributionMapBase::constructMap(), and distributionMapBase::subMap().

Referenced by patchToPatch::transferListList().

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

◆ rDistributeListList() [2/2]

void rDistributeListList ( const label  size,
const distributionMap map,
List< DynamicList< Type >> &  data 
)
inlinestaticprotected

Reverse distribute a dynamic list-list given the map.

Definition at line 72 of file patchToPatchI.H.

◆ localFacesToProcFaces()

Foam::List< Foam::List< Foam::patchToPatch::procFace > > localFacesToProcFaces ( const List< DynamicList< label >> &  localFaces,
const List< procFace > &  map = NullObjectRef<List<procFace>>() 
)
staticprotected

Map local faces to proc faces.

Definition at line 97 of file patchToPatch.C.

References forAll, Foam::isNull(), UPstream::myProcNo(), and patchToPatch::procFacesToLocalFaces().

Referenced by Foam::combine(), patchToPatch::srcTgtProcFaces(), and patchToPatch::tgtSrcProcFaces().

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

◆ procFacesToLocalFaces()

Foam::List< Foam::DynamicList< Foam::label > > procFacesToLocalFaces ( const List< List< procFace >> &  procFaces,
const HashTable< label, procFace, Hash< procFace >> &  map 
)
staticprotected

Map proc faces to local faces.

Definition at line 123 of file patchToPatch.C.

References forAll, Foam::isNull(), and patchToPatch::srcBox().

Referenced by patchToPatch::localFacesToProcFaces().

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

◆ 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.

Referenced by patchToPatch::procFacesToLocalFaces(), patchToPatch::srcBox(), and patchToPatch::update().

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 170 of file patchToPatch.C.

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

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 195 of file patchToPatch.C.

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

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

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 215 of file patchToPatch.C.

References Foam::combine(), patchToPatch::findOrIntersectFaces(), 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.

◆ 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 233 of file patchToPatch.C.

References forAll, and patchToPatch::intersectPatchQueue().

Referenced by patchToPatch::tgtBox().

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

◆ 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()

◆ calcSingleProcess()

void calcSingleProcess ( const primitiveOldTimePatch srcPatch,
const primitiveOldTimePatch tgtPatch 
)
protected

Determine whether or not the intersection of the given patches.

is on a single process. Set singleProcess_ to that process if so, and to -1 if not.

Definition at line 34 of file patchToPatchParallelOps.C.

References Foam::count(), Foam::findIndex(), and patchToPatch::tgtPatchSendFaces().

Referenced by patchToPatch::update().

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

◆ 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 72 of file patchToPatchParallelOps.C.

References List< T >::append(), forAll, and List< T >::transfer().

Referenced by patchToPatch::calcSingleProcess().

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

◆ srcPatchSendFaces()

Foam::labelListList srcPatchSendFaces ( ) const
protected

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

This is done after intersection, using the addressing generated by the intersection.

Definition at line 129 of file patchToPatchParallelOps.C.

References forAll, patchToPatch::localTgtProcFacesPtr_, UPstream::nProcs(), patchToPatch::patchDistributionMap(), and patchToPatch::tgtLocalSrcFaces_.

Here is the call graph for this function:

◆ patchDistributionMap()

Foam::distributionMap patchDistributionMap ( labelListList &&  sendFaces) const
protected

Create a distribution map from a list-list of faces to be sent.

(i.e., the result of calling one of the above methods).

Definition at line 157 of file patchToPatchParallelOps.C.

References patchToPatch::distributePatch(), forAll, Pstream::gatherList(), Foam::identity(), UPstream::myProcNo(), n, UPstream::nProcs(), Pstream::scatterList(), List< T >::setSize(), and List< T >::size().

Referenced by patchToPatch::srcPatchSendFaces().

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

◆ distributePatch() [1/2]

void distributePatch ( const distributionMap map,
List< procFace > &  localProcFaces 
) const
protected

Distribute a patch given its distribution map. Just generate.

local-proc-face addressing; not the distributed patch itself.

Definition at line 217 of file patchToPatchParallelOps.C.

References distributionMapBase::constructMap(), PstreamBuffers::finishedSends(), forAll, UPstream::myProcNo(), UPstream::nonBlocking, UPstream::nProcs(), List< T >::setSize(), List< T >::size(), and distributionMapBase::subMap().

Referenced by patchToPatch::patchDistributionMap().

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

◆ distributePatch() [2/2]

◆ initialise()

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

Initialisation.

Definition at line 690 of file patchToPatch.C.

References patchToPatch::distributeTgt(), and forAll.

Referenced by patchToPatch::intersectPatches(), intersection::srcBoxStatic(), and patchToPatch::update().

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

◆ distributeTgt()

Foam::tmpNrc< Foam::PrimitiveOldTimePatch< Foam::faceList, Foam::pointField > > distributeTgt ( const primitiveOldTimePatch srcPatch,
const vectorField srcPointNormals,
const vectorField srcPointNormals0,
const primitiveOldTimePatch tgtPatch,
distributionMap tgtMap 
)
protectedvirtual

Distribute the target patch so that enough is locally available.

for its intersection with the source patch can be computed. Happens before intersection. Bound boxes are used to determine what is needed where. Return the resulting local patch and the map from which it was calculated.

Definition at line 713 of file patchToPatch.C.

References patchToPatch::distributeSrc().

Referenced by patchToPatch::initialise(), and patchToPatch::update().

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

◆ distributeSrc()

Foam::tmpNrc< Foam::PrimitiveOldTimePatch< Foam::faceList, Foam::pointField > > distributeSrc ( const primitiveOldTimePatch srcPatch,
distributionMap srcMap 
)
protectedvirtual

Distribute the source patch so that everything the target.

intersects is locally available. Happens after intersection. The intersection addressing is used to determine what is needed where. Return the resulting local patch and the map from which it was calculated.

Definition at line 751 of file patchToPatch.C.

References patchToPatch::rDistributeTgt().

Referenced by patchToPatch::distributeTgt(), and patchToPatch::update().

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

◆ rDistributeTgt()

void rDistributeTgt ( const primitiveOldTimePatch tgtPatch,
const distributionMap tgtMap 
)
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.

Definition at line 775 of file patchToPatch.C.

References patchToPatch::finalise(), forAll, and HashTable< T, Key, Hash >::insert().

Referenced by patchToPatch::distributeSrc(), intersection::srcBoxStatic(), and patchToPatch::update().

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.

Definition at line 804 of file patchToPatch.C.

References forAll, and Foam::reduce().

Referenced by patchToPatch::rDistributeTgt(), intersection::srcBoxStatic(), and patchToPatch::update().

Here is the call graph for this function:
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 851 of file patchToPatch.C.

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

Referenced by patchToPatch::~patchToPatch().

Here is the call graph for this function:
Here is the caller 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 87 of file patchToPatchI.H.

References patchToPatch::reverse_.

Referenced by patchToPatch::New().

Here is the caller graph for this function:

◆ 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 93 of file patchToPatchI.H.

References patchToPatch::singleProcess_.

◆ isSingleProcess()

bool isSingleProcess ( ) const
inline

Is this intersection on a single process?

Definition at line 99 of file patchToPatchI.H.

References patchToPatch::singleProcess_.

Referenced by patchToPatch::srcTgtProcFaces(), patchToPatch::tgtSrcProcFaces(), and patchToPatch::update().

Here is the caller graph for this function:

◆ srcTgtProcFaces()

Foam::List< Foam::List< Foam::patchToPatch::procFace > > srcTgtProcFaces ( ) const
inline

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

Definition at line 106 of file patchToPatchI.H.

References patchToPatch::isSingleProcess(), patchToPatch::localFacesToProcFaces(), patchToPatch::localTgtProcFacesPtr_, and patchToPatch::srcLocalTgtFaces_.

Referenced by Foam::meshEdge().

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

◆ tgtSrcProcFaces()

Foam::List< Foam::List< Foam::patchToPatch::procFace > > tgtSrcProcFaces ( ) const
inline

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

Definition at line 116 of file patchToPatchI.H.

References patchToPatch::isSingleProcess(), patchToPatch::localFacesToProcFaces(), patchToPatch::localSrcProcFacesPtr_, and patchToPatch::tgtLocalSrcFaces_.

Referenced by Foam::meshEdge().

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

◆ srcWeights()

virtual tmpNrc<List<DynamicList<scalar> > > srcWeights ( const primitivePatch srcPatch) const
pure virtual

For each source face, the coupled target weights.

Implemented in intersection, rays, inverseDistance, and nearest.

◆ tgtWeights()

virtual tmpNrc<List<DynamicList<scalar> > > tgtWeights ( const primitivePatch tgtPatch) const
pure virtual

For each target face, the coupled source weights.

Implemented in intersection, rays, inverseDistance, and nearest.

◆ 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 1045 of file patchToPatch.C.

References patchToPatch::update().

Here is the call graph for this function:

◆ operator=()

void operator= ( const patchToPatch )
delete

Disallow default bitwise assignment.

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 100 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 104 of file patchToPatch.H.

Referenced by patchToPatch::isSingleProcess(), and patchToPatch::singleProcess().

◆ localSrcProcFacesPtr_

autoPtr<List<procFace> > localSrcProcFacesPtr_
protected

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

source processor and face index

Definition at line 108 of file patchToPatch.H.

Referenced by patchToPatch::tgtSrcProcFaces().

◆ localTgtProcFacesPtr_

autoPtr<List<procFace> > localTgtProcFacesPtr_
protected

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

target processor and face index

Definition at line 112 of file patchToPatch.H.

Referenced by patchToPatch::srcPatchSendFaces(), and patchToPatch::srcTgtProcFaces().

◆ srcLocalTgtFaces_

List<DynamicList<label> > srcLocalTgtFaces_
protected

For each source face, the coupled local target faces.

Definition at line 115 of file patchToPatch.H.

Referenced by patchToPatch::srcTgtProcFaces(), and nearest::srcWeights().

◆ tgtLocalSrcFaces_

List<DynamicList<label> > tgtLocalSrcFaces_
protected

For each target face, the coupled local source faces.

Definition at line 118 of file patchToPatch.H.

Referenced by patchToPatch::srcPatchSendFaces(), patchToPatch::tgtSrcProcFaces(), and nearest::tgtWeights().


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