meshToMesh Class Reference

Class to calculate the cell-addressing between two overlapping meshes. More...

Collaboration diagram for meshToMesh:
[legend]

Public Types

enum class  interpolationMethod { imDirect , imMapNearest , imCellVolumeWeight , imCorrectedCellVolumeWeight }
 Enumeration specifying interpolation method. More...
 
enum class  procMapMethod { pmAABB , pmLOD }
 Enumeration specifying processor parallel map construction method. More...
 

Public Member Functions

 TypeName ("meshToMesh")
 Run-time type information. More...
 
 meshToMesh (const polyMesh &src, const polyMesh &tgt, const interpolationMethod &method, const procMapMethod &mapMethod=procMapMethod::pmAABB, const bool interpAllPatches=true)
 Construct from source and target meshes. More...
 
 meshToMesh (const polyMesh &src, const polyMesh &tgt, const word &methodName, const word &AMIMethodName, const procMapMethod &mapMethod=procMapMethod::pmAABB, const bool interpAllPatches=true)
 Construct from source and target meshes, generic mapping methods. More...
 
 meshToMesh (const polyMesh &src, const polyMesh &tgt, const interpolationMethod &method, const HashTable< word > &patchMap, const wordList &cuttingPatches, const procMapMethod &mapMethod=procMapMethod::pmAABB, const bool normalise=true)
 Construct from source and target meshes. More...
 
 meshToMesh (const polyMesh &src, const polyMesh &tgt, const word &methodName, const word &AMIMethodName, const HashTable< word > &patchMap, const wordList &cuttingPatches, const procMapMethod &mapMethod=procMapMethod::pmAABB, const bool normalise=true)
 Construct from source and target meshes, generic mapping methods. More...
 
virtual ~meshToMesh ()
 Destructor. More...
 
const polyMeshsrcRegion () const
 Return const access to the source mesh. More...
 
const polyMeshtgtRegion () const
 Return const access to the target mesh. More...
 
const labelListListsrcToTgtCellAddr () const
 Return const access to the source to target cell addressing. More...
 
const labelListListtgtToSrcCellAddr () const
 Return const access to the target to source cell addressing. More...
 
const scalarListListsrcToTgtCellWght () const
 Return const access to the source to target cell weights. More...
 
const scalarListListtgtToSrcCellWght () const
 Return const access to the target to source cell weights. More...
 
const pointListListsrcToTgtCellVec () const
 Return const access to the source to target offset vectors. More...
 
const pointListListtgtToSrcCellVec () const
 Return const access to the target to source offset vectors. More...
 
scalar V () const
 Return const access to the overlap volume. More...
 
const PtrList< AMIPatchToPatchInterpolation > & patchAMIs () const
 Return the list of AMIs between source and target patches. More...
 
const autoPtr< mapDistribute > & srcMap () const
 Source map pointer - valid if no singleMeshProc. More...
 
const autoPtr< mapDistribute > & tgtMap () const
 Target map pointer - valid if no singleMeshProc. More...
 
template<class Type , class CombineOp >
void mapSrcToTgt (const UList< Type > &srcFld, const CombineOp &cop, List< Type > &result) const
 Map field from src to tgt mesh with defined operation. More...
 
template<class Type , class CombineOp >
void mapSrcToTgt (const UList< Type > &srcField, const UList< typename outerProduct< vector, Type >::type > &, const CombineOp &cop, List< Type > &result) const
 Map extrapolated field (using gradient) from src to tgt. More...
 
template<class Type , class CombineOp >
tmp< Field< Type > > mapSrcToTgt (const Field< Type > &srcFld, const CombineOp &cop) const
 Return the src field mapped to the tgt mesh with a defined. More...
 
template<class Type , class CombineOp >
tmp< Field< Type > > mapSrcToTgt (const tmp< Field< Type > > &tsrcFld, const CombineOp &cop) const
 Convenience function to map a tmp field to the tgt mesh. More...
 
template<class Type >
tmp< Field< Type > > mapSrcToTgt (const Field< Type > &srcFld) const
 Convenience function to map a field to the tgt mesh with a. More...
 
template<class Type >
tmp< Field< Type > > mapSrcToTgt (const tmp< Field< Type > > &tsrcFld) const
 Convenience function to map a tmp field to the tgt mesh. More...
 
template<class Type , class CombineOp >
void mapTgtToSrc (const UList< Type > &tgtFld, const CombineOp &cop, List< Type > &result) const
 Map field from tgt to src mesh with defined operation. More...
 
template<class Type , class CombineOp >
void mapTgtToSrc (const UList< Type > &srcField, const UList< typename outerProduct< vector, Type >::type > &, const CombineOp &cop, List< Type > &result) const
 Map extrapolated field (using gradient) from tgt to src. More...
 
template<class Type , class CombineOp >
tmp< Field< Type > > mapTgtToSrc (const Field< Type > &tgtFld, const CombineOp &cop) const
 Return the tgt field mapped to the src mesh with a defined. More...
 
template<class Type , class CombineOp >
tmp< Field< Type > > mapTgtToSrc (const tmp< Field< Type > > &ttgtFld, const CombineOp &cop) const
 Convenience function to map a tmp field to the src mesh. More...
 
template<class Type >
tmp< Field< Type > > mapTgtToSrc (const Field< Type > &tgtFld) const
 Convenience function to map a field to the src mesh with a. More...
 
template<class Type >
tmp< Field< Type > > mapTgtToSrc (const tmp< Field< Type > > &ttgtFld) const
 Convenience function to map a tmp field to the src mesh. More...
 
template<class Type , class CombineOp >
void mapSrcToTgt (const VolumeField< Type > &field, const CombineOp &cop, VolumeField< Type > &result, const bool secondOrder=true) const
 Interpolate a field with a defined operation. Values. More...
 
template<class Type , class CombineOp >
tmp< VolumeField< Type > > mapSrcToTgt (const VolumeField< Type > &field, const CombineOp &cop, const bool secondOrder=true) const
 Interpolate a field with a defined operation. The initial. More...
 
template<class Type , class CombineOp >
tmp< VolumeField< Type > > mapSrcToTgt (const tmp< VolumeField< Type > > &tfield, const CombineOp &cop, const bool secondOrder=true) const
 Interpolate a tmp field with a defined operation. The. More...
 
template<class Type >
tmp< VolumeField< Type > > mapSrcToTgt (const VolumeField< Type > &field, const bool secondOrder=true) const
 Convenience function to map a field with a default. More...
 
template<class Type >
tmp< VolumeField< Type > > mapSrcToTgt (const tmp< VolumeField< Type > > &tfield, const bool secondOrder=true) const
 Convenience function to map a tmp field with a default. More...
 
template<class Type , class CombineOp >
void mapTgtToSrc (const VolumeField< Type > &field, const CombineOp &cop, VolumeField< Type > &result, const bool secondOrder=true) const
 Interpolate a field with a defined operation. Values. More...
 
template<class Type , class CombineOp >
tmp< VolumeField< Type > > mapTgtToSrc (const VolumeField< Type > &field, const CombineOp &cop, const bool secondOrder=true) const
 Interpolate a field with a defined operation. The initial. More...
 
template<class Type , class CombineOp >
tmp< VolumeField< Type > > mapTgtToSrc (const tmp< VolumeField< Type > > &tfield, const CombineOp &cop, const bool secondOrder=true) const
 Interpolate a tmp field with a defined operation. The. More...
 
template<class Type >
tmp< VolumeField< Type > > mapTgtToSrc (const VolumeField< Type > &field, const bool secondOrder=true) const
 Convenience function to map a field with a default. More...
 
template<class Type >
tmp< VolumeField< Type > > mapTgtToSrc (const tmp< VolumeField< Type > > &tfield, const bool secondOrder=true) const
 Convenience function to map a tmp field with a default. More...
 
template<>
void mapInternalSrcToTgt (const VolumeField< sphericalTensor > &field, const plusEqOp< sphericalTensor > &cop, VolumeField< sphericalTensor > &result, const bool secondOrder) const
 
template<>
void mapInternalSrcToTgt (const VolumeField< sphericalTensor > &field, const minusEqOp< sphericalTensor > &cop, VolumeField< sphericalTensor > &result, const bool secondOrder) const
 
template<>
void mapInternalSrcToTgt (const VolumeField< symmTensor > &field, const plusEqOp< symmTensor > &cop, VolumeField< symmTensor > &result, const bool secondOrder) const
 
template<>
void mapInternalSrcToTgt (const VolumeField< symmTensor > &field, const minusEqOp< symmTensor > &cop, VolumeField< symmTensor > &result, const bool secondOrder) const
 
template<>
void mapInternalSrcToTgt (const VolumeField< tensor > &field, const plusEqOp< tensor > &cop, VolumeField< tensor > &result, const bool secondOrder) const
 
template<>
void mapInternalSrcToTgt (const VolumeField< tensor > &field, const minusEqOp< tensor > &cop, VolumeField< tensor > &result, const bool secondOrder) const
 
template<>
void mapInternalTgtToSrc (const VolumeField< sphericalTensor > &field, const plusEqOp< sphericalTensor > &cop, VolumeField< sphericalTensor > &result, const bool secondOrder) const
 
template<>
void mapInternalTgtToSrc (const VolumeField< sphericalTensor > &field, const minusEqOp< sphericalTensor > &cop, VolumeField< sphericalTensor > &result, const bool secondOrder) const
 
template<>
void mapInternalTgtToSrc (const VolumeField< symmTensor > &field, const plusEqOp< symmTensor > &cop, VolumeField< symmTensor > &result, const bool secondOrder) const
 
template<>
void mapInternalTgtToSrc (const VolumeField< symmTensor > &field, const minusEqOp< symmTensor > &cop, VolumeField< symmTensor > &result, const bool secondOrder) const
 
template<>
void mapInternalTgtToSrc (const VolumeField< tensor > &field, const plusEqOp< tensor > &cop, VolumeField< tensor > &result, const bool secondOrder) const
 
template<>
void mapInternalTgtToSrc (const VolumeField< tensor > &field, const minusEqOp< tensor > &cop, VolumeField< tensor > &result, const bool secondOrder) const
 
template<>
void mapAndOpSrcToTgt (const AMIPatchToPatchInterpolation &AMI, const Field< scalar > &srcField, Field< scalar > &tgtField, const plusEqOp< scalar > &cop) const
 
template<>
void mapAndOpSrcToTgt (const AMIPatchToPatchInterpolation &AMI, const Field< vector > &srcField, Field< vector > &tgtField, const plusEqOp< vector > &cop) const
 
template<>
void mapAndOpSrcToTgt (const AMIPatchToPatchInterpolation &AMI, const Field< sphericalTensor > &srcField, Field< sphericalTensor > &tgtField, const plusEqOp< sphericalTensor > &cop) const
 
template<>
void mapAndOpSrcToTgt (const AMIPatchToPatchInterpolation &AMI, const Field< symmTensor > &srcField, Field< symmTensor > &tgtField, const plusEqOp< symmTensor > &cop) const
 
template<>
void mapAndOpSrcToTgt (const AMIPatchToPatchInterpolation &AMI, const Field< tensor > &srcField, Field< tensor > &tgtField, const plusEqOp< tensor > &cop) const
 
template<>
void mapAndOpTgtToSrc (const AMIPatchToPatchInterpolation &AMI, Field< scalar > &srcField, const Field< scalar > &tgtField, const plusEqOp< scalar > &cop) const
 
template<>
void mapAndOpTgtToSrc (const AMIPatchToPatchInterpolation &AMI, Field< vector > &srcField, const Field< vector > &tgtField, const plusEqOp< vector > &cop) const
 
template<>
void mapAndOpTgtToSrc (const AMIPatchToPatchInterpolation &AMI, Field< sphericalTensor > &srcField, const Field< sphericalTensor > &tgtField, const plusEqOp< sphericalTensor > &cop) const
 
template<>
void mapAndOpTgtToSrc (const AMIPatchToPatchInterpolation &AMI, Field< symmTensor > &srcField, const Field< symmTensor > &tgtField, const plusEqOp< symmTensor > &cop) const
 
template<>
void mapAndOpTgtToSrc (const AMIPatchToPatchInterpolation &AMI, Field< tensor > &srcField, const Field< tensor > &tgtField, const plusEqOp< tensor > &cop) const
 
template<class Type , class CombineOp >
Foam::tmp< Foam::Field< Type > > mapSrcToTgt (const Field< Type > &srcField, const CombineOp &cop) const
 
template<class Type , class CombineOp >
Foam::tmp< Foam::Field< Type > > mapSrcToTgt (const tmp< Field< Type > > &tsrcField, const CombineOp &cop) const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > mapSrcToTgt (const Field< Type > &srcField) const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > mapSrcToTgt (const tmp< Field< Type > > &tsrcField) const
 
template<class Type , class CombineOp >
Foam::tmp< Foam::Field< Type > > mapTgtToSrc (const Field< Type > &tgtField, const CombineOp &cop) const
 
template<class Type , class CombineOp >
Foam::tmp< Foam::Field< Type > > mapTgtToSrc (const tmp< Field< Type > > &ttgtField, const CombineOp &cop) const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > mapTgtToSrc (const Field< Type > &tgtField) const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > mapTgtToSrc (const tmp< Field< Type > > &ttgtField) const
 
template<class Type , class CombineOp >
Foam::tmp< Foam::VolumeField< Type > > mapSrcToTgt (const VolumeField< Type > &field, const CombineOp &cop, const bool secondOrder) const
 
template<class Type , class CombineOp >
Foam::tmp< Foam::VolumeField< Type > > mapSrcToTgt (const tmp< VolumeField< Type > > &tfield, const CombineOp &cop, const bool secondOrder) const
 
template<class Type >
Foam::tmp< Foam::VolumeField< Type > > mapSrcToTgt (const VolumeField< Type > &field, const bool secondOrder) const
 
template<class Type >
Foam::tmp< Foam::VolumeField< Type > > mapSrcToTgt (const tmp< VolumeField< Type > > &tfield, const bool secondOrder) const
 
template<class Type , class CombineOp >
Foam::tmp< Foam::VolumeField< Type > > mapTgtToSrc (const VolumeField< Type > &field, const CombineOp &cop, const bool secondOrder) const
 
template<class Type , class CombineOp >
Foam::tmp< Foam::VolumeField< Type > > mapTgtToSrc (const tmp< VolumeField< Type > > &tfield, const CombineOp &cop, const bool secondOrder) const
 
template<class Type >
Foam::tmp< Foam::VolumeField< Type > > mapTgtToSrc (const VolumeField< Type > &field, const bool secondOrder) const
 
template<class Type >
Foam::tmp< Foam::VolumeField< Type > > mapTgtToSrc (const tmp< VolumeField< Type > > &tfield, const bool secondOrder) const
 

Static Public Member Functions

static word interpolationMethodAMI (const interpolationMethod method)
 Conversion between mesh and patch interpolation methods. More...
 

Static Public Attributes

static const Enum< interpolationMethodinterpolationMethodNames_
 
static const Enum< procMapMethodprocMapMethodNames_
 

Detailed Description

Class to calculate the cell-addressing between two overlapping meshes.

Mapping is performed using a run-time selectable interpolation mothod

See also
meshToMeshMethod
Source files

Definition at line 64 of file meshToMesh.H.

Member Enumeration Documentation

◆ interpolationMethod

enum class interpolationMethod
strong

Enumeration specifying interpolation method.

Enumerator
imDirect 
imMapNearest 
imCellVolumeWeight 
imCorrectedCellVolumeWeight 

Definition at line 71 of file meshToMesh.H.

◆ procMapMethod

enum class procMapMethod
strong

Enumeration specifying processor parallel map construction method.

Enumerator
pmAABB 
pmLOD 

Definition at line 82 of file meshToMesh.H.

Constructor & Destructor Documentation

◆ meshToMesh() [1/4]

meshToMesh ( const polyMesh src,
const polyMesh tgt,
const interpolationMethod method,
const procMapMethod mapMethod = procMapMethod::pmAABB,
const bool  interpAllPatches = true 
)

Construct from source and target meshes.

Definition at line 843 of file meshToMesh.C.

References meshToMesh::interpolationMethodAMI(), and meshToMesh::interpolationMethodNames_.

Here is the call graph for this function:

◆ meshToMesh() [2/4]

meshToMesh ( const polyMesh src,
const polyMesh tgt,
const word methodName,
const word AMIMethodName,
const procMapMethod mapMethod = procMapMethod::pmAABB,
const bool  interpAllPatches = true 
)

Construct from source and target meshes, generic mapping methods.

Definition at line 879 of file meshToMesh.C.

◆ meshToMesh() [3/4]

meshToMesh ( const polyMesh src,
const polyMesh tgt,
const interpolationMethod method,
const HashTable< word > &  patchMap,
const wordList cuttingPatches,
const procMapMethod mapMethod = procMapMethod::pmAABB,
const bool  normalise = true 
)

Construct from source and target meshes.

Definition at line 911 of file meshToMesh.C.

References meshToMesh::interpolationMethodAMI(), and meshToMesh::interpolationMethodNames_.

Here is the call graph for this function:

◆ meshToMesh() [4/4]

meshToMesh ( const polyMesh src,
const polyMesh tgt,
const word methodName,
const word AMIMethodName,
const HashTable< word > &  patchMap,
const wordList cuttingPatches,
const procMapMethod mapMethod = procMapMethod::pmAABB,
const bool  normalise = true 
)

Construct from source and target meshes, generic mapping methods.

Definition at line 949 of file meshToMesh.C.

◆ ~meshToMesh()

~meshToMesh ( )
virtual

Destructor.

Definition at line 990 of file meshToMesh.C.

Member Function Documentation

◆ TypeName()

TypeName ( "meshToMesh"  )

Run-time type information.

◆ srcRegion()

const Foam::polyMesh & srcRegion ( ) const
inline

Return const access to the source mesh.

Definition at line 33 of file meshToMeshI.H.

Referenced by Foam::MapMesh(), and Foam::MapVolFields().

Here is the caller graph for this function:

◆ tgtRegion()

const Foam::polyMesh & tgtRegion ( ) const
inline

Return const access to the target mesh.

Definition at line 39 of file meshToMeshI.H.

Referenced by Foam::MapMesh(), and Foam::MapVolFields().

Here is the caller graph for this function:

◆ srcToTgtCellAddr()

const Foam::labelListList & srcToTgtCellAddr ( ) const
inline

Return const access to the source to target cell addressing.

Definition at line 45 of file meshToMeshI.H.

Referenced by cellVolumeWeight::update().

Here is the caller graph for this function:

◆ tgtToSrcCellAddr()

const Foam::labelListList & tgtToSrcCellAddr ( ) const
inline

Return const access to the target to source cell addressing.

Definition at line 51 of file meshToMeshI.H.

Referenced by cellVolumeWeight::update().

Here is the caller graph for this function:

◆ srcToTgtCellWght()

const Foam::scalarListList & srcToTgtCellWght ( ) const
inline

Return const access to the source to target cell weights.

Definition at line 57 of file meshToMeshI.H.

Referenced by cellVolumeWeight::update().

Here is the caller graph for this function:

◆ tgtToSrcCellWght()

const Foam::scalarListList & tgtToSrcCellWght ( ) const
inline

Return const access to the target to source cell weights.

Definition at line 63 of file meshToMeshI.H.

Referenced by cellVolumeWeight::update().

Here is the caller graph for this function:

◆ srcToTgtCellVec()

const Foam::pointListList & srcToTgtCellVec ( ) const
inline

Return const access to the source to target offset vectors.

Definition at line 69 of file meshToMeshI.H.

◆ tgtToSrcCellVec()

const Foam::pointListList & tgtToSrcCellVec ( ) const
inline

Return const access to the target to source offset vectors.

Definition at line 75 of file meshToMeshI.H.

◆ V()

Foam::scalar V ( ) const
inline

Return const access to the overlap volume.

Definition at line 81 of file meshToMeshI.H.

◆ interpolationMethodAMI()

Foam::word interpolationMethodAMI ( const interpolationMethod  method)
static

Conversion between mesh and patch interpolation methods.

Definition at line 643 of file meshToMesh.C.

References Foam::abort(), Foam::FatalError, FatalErrorInFunction, and ensightPTraits< Type >::typeName.

Referenced by meshToMesh::meshToMesh().

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

◆ patchAMIs()

const Foam::PtrList< Foam::AMIPatchToPatchInterpolation > & patchAMIs ( ) const
inline

Return the list of AMIs between source and target patches.

Definition at line 102 of file meshToMeshI.H.

◆ srcMap()

const Foam::autoPtr< Foam::mapDistribute > & srcMap ( ) const
inline

Source map pointer - valid if no singleMeshProc.

Definition at line 88 of file meshToMeshI.H.

Referenced by cellVolumeWeight::update().

Here is the caller graph for this function:

◆ tgtMap()

const Foam::autoPtr< Foam::mapDistribute > & tgtMap ( ) const
inline

Target map pointer - valid if no singleMeshProc.

Definition at line 95 of file meshToMeshI.H.

Referenced by cellVolumeWeight::update().

Here is the caller graph for this function:

◆ mapSrcToTgt() [1/19]

void mapSrcToTgt ( const UList< Type > &  srcFld,
const CombineOp &  cop,
List< Type > &  result 
) const

Map field from src to tgt mesh with defined operation.

Values passed in via 'result' are used to initialise the return value

Definition at line 53 of file meshToMeshTemplates.C.

References Foam::abort(), mapDistribute::distribute(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::nl, UList< T >::size(), and Foam::sum().

Referenced by interRegionHeatTransferModel::interpolate(), meshToMesh::mapInternalSrcToTgt(), and Foam::MapVolFields().

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

◆ mapSrcToTgt() [2/19]

void mapSrcToTgt ( const UList< Type > &  srcField,
const UList< typename outerProduct< vector, Type >::type > &  srcGradField,
const CombineOp &  cop,
List< Type > &  result 
) const

Map extrapolated field (using gradient) from src to tgt.

mesh with defined operation. Falls back to non-extrapolated mapping (above) if not constructed with method that supports getting offset vectors. Extrapolation only for internal values. Values passed in via 'result' are used to initialise the return value.

Definition at line 121 of file meshToMeshTemplates.C.

References Foam::abort(), mapDistribute::distribute(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::nl, Foam::returnReduce(), UList< T >::size(), and Foam::sum().

Here is the call graph for this function:

◆ mapSrcToTgt() [3/19]

tmp< Field< Type > > mapSrcToTgt ( const Field< Type > &  srcFld,
const CombineOp &  cop 
) const

Return the src field mapped to the tgt mesh with a defined.

operation. Initial values of the result are set to zero

◆ mapSrcToTgt() [4/19]

tmp< Field< Type > > mapSrcToTgt ( const tmp< Field< Type > > &  tsrcFld,
const CombineOp &  cop 
) const

Convenience function to map a tmp field to the tgt mesh.

with a defined operation

◆ mapSrcToTgt() [5/19]

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

Convenience function to map a field to the tgt mesh with a.

default operation (plusEqOp)

◆ mapSrcToTgt() [6/19]

tmp< Field< Type > > mapSrcToTgt ( const tmp< Field< Type > > &  tsrcFld) const

Convenience function to map a tmp field to the tgt mesh.

with a default operation (plusEqOp)

◆ mapTgtToSrc() [1/19]

void mapTgtToSrc ( const UList< Type > &  tgtFld,
const CombineOp &  cop,
List< Type > &  result 
) const

Map field from tgt to src mesh with defined operation.

Values passed in via 'result' are used to initialise the return value

Definition at line 268 of file meshToMeshTemplates.C.

References Foam::abort(), mapDistribute::distribute(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::nl, UList< T >::size(), and Foam::sum().

Referenced by interRegionHeatTransferModel::interpolate().

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

◆ mapTgtToSrc() [2/19]

void mapTgtToSrc ( const UList< Type > &  srcField,
const UList< typename outerProduct< vector, Type >::type > &  tgtGradField,
const CombineOp &  cop,
List< Type > &  result 
) const

Map extrapolated field (using gradient) from tgt to src.

mesh with defined operation. Falls back to non-extrapolated mapping (above) if not constructed with method that supports getting offset vectors. Extrapolation only for internal values. Values passed in via 'result' are used to initialise the return value

Definition at line 334 of file meshToMeshTemplates.C.

References Foam::abort(), mapDistribute::distribute(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::nl, Foam::returnReduce(), UList< T >::size(), and Foam::sum().

Here is the call graph for this function:

◆ mapTgtToSrc() [3/19]

tmp< Field< Type > > mapTgtToSrc ( const Field< Type > &  tgtFld,
const CombineOp &  cop 
) const

Return the tgt field mapped to the src mesh with a defined.

operation. Initial values of the result are set to zero

◆ mapTgtToSrc() [4/19]

tmp< Field< Type > > mapTgtToSrc ( const tmp< Field< Type > > &  ttgtFld,
const CombineOp &  cop 
) const

Convenience function to map a tmp field to the src mesh.

with a defined operation

◆ mapTgtToSrc() [5/19]

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

Convenience function to map a field to the src mesh with a.

default operation (plusEqOp)

◆ mapTgtToSrc() [6/19]

tmp< Field< Type > > mapTgtToSrc ( const tmp< Field< Type > > &  ttgtFld) const

Convenience function to map a tmp field to the src mesh.

with a default operation (plusEqOp)

◆ mapSrcToTgt() [7/19]

void mapSrcToTgt ( const VolumeField< Type > &  field,
const CombineOp &  cop,
VolumeField< Type > &  result,
const bool  secondOrder = true 
) const

Interpolate a field with a defined operation. Values.

passed in via 'result' are used to initialise the return value. Optionally uses gradient correction (internal field only) if interpolationMethod supports it

Definition at line 520 of file meshToMeshTemplates.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), field(), forAll, Foam::identity(), fvPatchField< Type >::patch(), fvPatchField< Type >::patchInternalField(), fvPatchField< Type >::rmap(), and UList< T >::size().

Here is the call graph for this function:

◆ mapSrcToTgt() [8/19]

tmp< VolumeField< Type > > mapSrcToTgt ( const VolumeField< Type > &  field,
const CombineOp &  cop,
const bool  secondOrder = true 
) const

Interpolate a field with a defined operation. The initial.

values of the result are set to zero

◆ mapSrcToTgt() [9/19]

tmp< VolumeField< Type > > mapSrcToTgt ( const tmp< VolumeField< Type > > &  tfield,
const CombineOp &  cop,
const bool  secondOrder = true 
) const

Interpolate a tmp field with a defined operation. The.

initial values of the result are set to zero

◆ mapSrcToTgt() [10/19]

tmp< VolumeField< Type > > mapSrcToTgt ( const VolumeField< Type > &  field,
const bool  secondOrder = true 
) const

Convenience function to map a field with a default.

operation (plusEqOp)

◆ mapSrcToTgt() [11/19]

tmp< VolumeField< Type > > mapSrcToTgt ( const tmp< VolumeField< Type > > &  tfield,
const bool  secondOrder = true 
) const

Convenience function to map a tmp field with a default.

operation (plusEqOp)

◆ mapTgtToSrc() [7/19]

void mapTgtToSrc ( const VolumeField< Type > &  field,
const CombineOp &  cop,
VolumeField< Type > &  result,
const bool  secondOrder = true 
) const

Interpolate a field with a defined operation. Values.

passed in via 'result' are used to initialise the return value. Optionally uses gradient correction (internal field only) if interpolationMethod supports it

Definition at line 753 of file meshToMeshTemplates.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), field(), forAll, Foam::identity(), fvPatchField< Type >::patch(), fvPatchField< Type >::patchInternalField(), fvPatchField< Type >::rmap(), and UList< T >::size().

Here is the call graph for this function:

◆ mapTgtToSrc() [8/19]

tmp< VolumeField< Type > > mapTgtToSrc ( const VolumeField< Type > &  field,
const CombineOp &  cop,
const bool  secondOrder = true 
) const

Interpolate a field with a defined operation. The initial.

values of the result are set to zero

◆ mapTgtToSrc() [9/19]

tmp< VolumeField< Type > > mapTgtToSrc ( const tmp< VolumeField< Type > > &  tfield,
const CombineOp &  cop,
const bool  secondOrder = true 
) const

Interpolate a tmp field with a defined operation. The.

initial values of the result are set to zero

◆ mapTgtToSrc() [10/19]

tmp< VolumeField< Type > > mapTgtToSrc ( const VolumeField< Type > &  field,
const bool  secondOrder = true 
) const

Convenience function to map a field with a default.

operation (plusEqOp)

◆ mapTgtToSrc() [11/19]

tmp< VolumeField< Type > > mapTgtToSrc ( const tmp< VolumeField< Type > > &  tfield,
const bool  secondOrder = true 
) const

Convenience function to map a tmp field with a default.

operation (plusEqOp)

◆ mapInternalSrcToTgt() [1/6]

void mapInternalSrcToTgt ( const VolumeField< sphericalTensor > &  field,
const plusEqOp< sphericalTensor > &  cop,
VolumeField< sphericalTensor > &  result,
const bool  secondOrder 
) const

Definition at line 75 of file meshToMesh.C.

References field(), meshToMesh::mapSrcToTgt(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalSrcToTgt() [2/6]

void mapInternalSrcToTgt ( const VolumeField< sphericalTensor > &  field,
const minusEqOp< sphericalTensor > &  cop,
VolumeField< sphericalTensor > &  result,
const bool  secondOrder 
) const

Definition at line 88 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalSrcToTgt() [3/6]

void mapInternalSrcToTgt ( const VolumeField< symmTensor > &  field,
const plusEqOp< symmTensor > &  cop,
VolumeField< symmTensor > &  result,
const bool  secondOrder 
) const

Definition at line 101 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalSrcToTgt() [4/6]

void mapInternalSrcToTgt ( const VolumeField< symmTensor > &  field,
const minusEqOp< symmTensor > &  cop,
VolumeField< symmTensor > &  result,
const bool  secondOrder 
) const

Definition at line 114 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalSrcToTgt() [5/6]

void mapInternalSrcToTgt ( const VolumeField< tensor > &  field,
const plusEqOp< tensor > &  cop,
VolumeField< tensor > &  result,
const bool  secondOrder 
) const

Definition at line 127 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalSrcToTgt() [6/6]

void mapInternalSrcToTgt ( const VolumeField< tensor > &  field,
const minusEqOp< tensor > &  cop,
VolumeField< tensor > &  result,
const bool  secondOrder 
) const

Definition at line 140 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalTgtToSrc() [1/6]

void mapInternalTgtToSrc ( const VolumeField< sphericalTensor > &  field,
const plusEqOp< sphericalTensor > &  cop,
VolumeField< sphericalTensor > &  result,
const bool  secondOrder 
) const

Definition at line 153 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalTgtToSrc() [2/6]

void mapInternalTgtToSrc ( const VolumeField< sphericalTensor > &  field,
const minusEqOp< sphericalTensor > &  cop,
VolumeField< sphericalTensor > &  result,
const bool  secondOrder 
) const

Definition at line 166 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalTgtToSrc() [3/6]

void mapInternalTgtToSrc ( const VolumeField< symmTensor > &  field,
const plusEqOp< symmTensor > &  cop,
VolumeField< symmTensor > &  result,
const bool  secondOrder 
) const

Definition at line 179 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalTgtToSrc() [4/6]

void mapInternalTgtToSrc ( const VolumeField< symmTensor > &  field,
const minusEqOp< symmTensor > &  cop,
VolumeField< symmTensor > &  result,
const bool  secondOrder 
) const

Definition at line 192 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalTgtToSrc() [5/6]

void mapInternalTgtToSrc ( const VolumeField< tensor > &  field,
const plusEqOp< tensor > &  cop,
VolumeField< tensor > &  result,
const bool  secondOrder 
) const

Definition at line 205 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapInternalTgtToSrc() [6/6]

void mapInternalTgtToSrc ( const VolumeField< tensor > &  field,
const minusEqOp< tensor > &  cop,
VolumeField< tensor > &  result,
const bool  secondOrder 
) const

Definition at line 218 of file meshToMesh.C.

References field(), and GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef().

Here is the call graph for this function:

◆ mapAndOpSrcToTgt() [1/5]

void mapAndOpSrcToTgt ( const AMIPatchToPatchInterpolation AMI,
const Field< scalar > &  srcField,
Field< scalar > &  tgtField,
const plusEqOp< scalar > &  cop 
) const

Definition at line 231 of file meshToMesh.C.

◆ mapAndOpSrcToTgt() [2/5]

void mapAndOpSrcToTgt ( const AMIPatchToPatchInterpolation AMI,
const Field< vector > &  srcField,
Field< vector > &  tgtField,
const plusEqOp< vector > &  cop 
) const

Definition at line 242 of file meshToMesh.C.

◆ mapAndOpSrcToTgt() [3/5]

void mapAndOpSrcToTgt ( const AMIPatchToPatchInterpolation AMI,
const Field< sphericalTensor > &  srcField,
Field< sphericalTensor > &  tgtField,
const plusEqOp< sphericalTensor > &  cop 
) const

Definition at line 253 of file meshToMesh.C.

◆ mapAndOpSrcToTgt() [4/5]

void mapAndOpSrcToTgt ( const AMIPatchToPatchInterpolation AMI,
const Field< symmTensor > &  srcField,
Field< symmTensor > &  tgtField,
const plusEqOp< symmTensor > &  cop 
) const

Definition at line 264 of file meshToMesh.C.

◆ mapAndOpSrcToTgt() [5/5]

void mapAndOpSrcToTgt ( const AMIPatchToPatchInterpolation AMI,
const Field< tensor > &  srcField,
Field< tensor > &  tgtField,
const plusEqOp< tensor > &  cop 
) const

Definition at line 275 of file meshToMesh.C.

◆ mapAndOpTgtToSrc() [1/5]

void mapAndOpTgtToSrc ( const AMIPatchToPatchInterpolation AMI,
Field< scalar > &  srcField,
const Field< scalar > &  tgtField,
const plusEqOp< scalar > &  cop 
) const

Definition at line 286 of file meshToMesh.C.

◆ mapAndOpTgtToSrc() [2/5]

void mapAndOpTgtToSrc ( const AMIPatchToPatchInterpolation AMI,
Field< vector > &  srcField,
const Field< vector > &  tgtField,
const plusEqOp< vector > &  cop 
) const

Definition at line 297 of file meshToMesh.C.

◆ mapAndOpTgtToSrc() [3/5]

void mapAndOpTgtToSrc ( const AMIPatchToPatchInterpolation AMI,
Field< sphericalTensor > &  srcField,
const Field< sphericalTensor > &  tgtField,
const plusEqOp< sphericalTensor > &  cop 
) const

Definition at line 308 of file meshToMesh.C.

◆ mapAndOpTgtToSrc() [4/5]

void mapAndOpTgtToSrc ( const AMIPatchToPatchInterpolation AMI,
Field< symmTensor > &  srcField,
const Field< symmTensor > &  tgtField,
const plusEqOp< symmTensor > &  cop 
) const

Definition at line 319 of file meshToMesh.C.

◆ mapAndOpTgtToSrc() [5/5]

void mapAndOpTgtToSrc ( const AMIPatchToPatchInterpolation AMI,
Field< tensor > &  srcField,
const Field< tensor > &  tgtField,
const plusEqOp< tensor > &  cop 
) const

Definition at line 330 of file meshToMesh.C.

◆ mapSrcToTgt() [12/19]

Foam::tmp< Foam::Field< Type > > mapSrcToTgt ( const Field< Type > &  srcField,
const CombineOp &  cop 
) const

Definition at line 215 of file meshToMeshTemplates.C.

References tmp< T >::ref(), and Foam::Zero.

Here is the call graph for this function:

◆ mapSrcToTgt() [13/19]

Foam::tmp< Foam::Field< Type > > mapSrcToTgt ( const tmp< Field< Type > > &  tsrcField,
const CombineOp &  cop 
) const

Definition at line 237 of file meshToMeshTemplates.C.

◆ mapSrcToTgt() [14/19]

Foam::tmp< Foam::Field< Type > > mapSrcToTgt ( const Field< Type > &  srcField) const

Definition at line 248 of file meshToMeshTemplates.C.

◆ mapSrcToTgt() [15/19]

Foam::tmp< Foam::Field< Type > > mapSrcToTgt ( const tmp< Field< Type > > &  tsrcField) const

Definition at line 258 of file meshToMeshTemplates.C.

◆ mapTgtToSrc() [12/19]

Foam::tmp< Foam::Field< Type > > mapTgtToSrc ( const Field< Type > &  tgtField,
const CombineOp &  cop 
) const

Definition at line 420 of file meshToMeshTemplates.C.

References tmp< T >::ref(), and Foam::Zero.

Here is the call graph for this function:

◆ mapTgtToSrc() [13/19]

Foam::tmp< Foam::Field< Type > > mapTgtToSrc ( const tmp< Field< Type > > &  ttgtField,
const CombineOp &  cop 
) const

Definition at line 442 of file meshToMeshTemplates.C.

◆ mapTgtToSrc() [14/19]

Foam::tmp< Foam::Field< Type > > mapTgtToSrc ( const Field< Type > &  tgtField) const

Definition at line 453 of file meshToMeshTemplates.C.

◆ mapTgtToSrc() [15/19]

Foam::tmp< Foam::Field< Type > > mapTgtToSrc ( const tmp< Field< Type > > &  ttgtField) const

Definition at line 463 of file meshToMeshTemplates.C.

◆ mapSrcToTgt() [16/19]

Foam::tmp< Foam::VolumeField< Type > > mapSrcToTgt ( const VolumeField< Type > &  field,
const CombineOp &  cop,
const bool  secondOrder 
) const

Definition at line 584 of file meshToMeshTemplates.C.

References fvMesh::boundary(), field(), forAll, primitiveMesh::nCells(), Foam::New(), IOobject::NO_READ, IOobject::NO_WRITE, PtrList< T >::set(), UPtrList< T >::size(), fvMesh::time(), Time::timeName(), Foam::type(), and Foam::Zero.

Here is the call graph for this function:

◆ mapSrcToTgt() [17/19]

Foam::tmp< Foam::VolumeField< Type > > mapSrcToTgt ( const tmp< VolumeField< Type > > &  tfield,
const CombineOp &  cop,
const bool  secondOrder 
) const

Definition at line 670 of file meshToMeshTemplates.C.

◆ mapSrcToTgt() [18/19]

Foam::tmp< Foam::VolumeField< Type > > mapSrcToTgt ( const VolumeField< Type > &  field,
const bool  secondOrder 
) const

Definition at line 683 of file meshToMeshTemplates.C.

References field().

Here is the call graph for this function:

◆ mapSrcToTgt() [19/19]

Foam::tmp< Foam::VolumeField< Type > > mapSrcToTgt ( const tmp< VolumeField< Type > > &  tfield,
const bool  secondOrder 
) const

Definition at line 695 of file meshToMeshTemplates.C.

◆ mapTgtToSrc() [16/19]

Foam::tmp< Foam::VolumeField< Type > > mapTgtToSrc ( const VolumeField< Type > &  field,
const CombineOp &  cop,
const bool  secondOrder 
) const

Definition at line 815 of file meshToMeshTemplates.C.

References fvMesh::boundary(), field(), forAll, primitiveMesh::nCells(), Foam::New(), IOobject::NO_READ, IOobject::NO_WRITE, PtrList< T >::set(), UPtrList< T >::size(), fvMesh::time(), Time::timeName(), Foam::type(), and Foam::Zero.

Here is the call graph for this function:

◆ mapTgtToSrc() [17/19]

Foam::tmp< Foam::VolumeField< Type > > mapTgtToSrc ( const tmp< VolumeField< Type > > &  tfield,
const CombineOp &  cop,
const bool  secondOrder 
) const

Definition at line 901 of file meshToMeshTemplates.C.

◆ mapTgtToSrc() [18/19]

Foam::tmp< Foam::VolumeField< Type > > mapTgtToSrc ( const VolumeField< Type > &  field,
const bool  secondOrder 
) const

Definition at line 914 of file meshToMeshTemplates.C.

References field().

Here is the call graph for this function:

◆ mapTgtToSrc() [19/19]

Foam::tmp< Foam::VolumeField< Type > > mapTgtToSrc ( const tmp< VolumeField< Type > > &  tfield,
const bool  secondOrder 
) const

Definition at line 926 of file meshToMeshTemplates.C.

Member Data Documentation

◆ interpolationMethodNames_

const Foam::Enum< Foam::meshToMesh::interpolationMethod > interpolationMethodNames_
static

Definition at line 79 of file meshToMesh.H.

Referenced by meshToMesh::meshToMesh(), and interRegionOption::setMapper().

◆ procMapMethodNames_

const Foam::Enum< Foam::meshToMesh::procMapMethod > procMapMethodNames_
static
Initial value:

Definition at line 88 of file meshToMesh.H.

Referenced by interRegionOption::setMapper().


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