NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
TransferBase Class Referenceabstract

Detailed Description

Definition at line 38 of file TransferBase.H.

Public Member Functions

 TransferBase ()
 
virtual ~TransferBase ()
 
virtual int transferPointData (const std::vector< int > &arrayIDs, const std::vector< std::string > &newnames=std::vector< std::string >())=0
 Transfer point data with given ids from source to target. More...
 
virtual int transferCellData (const std::vector< int > &arrayIDs, const std::vector< std::string > &newnames=std::vector< std::string >())=0
 Transfer cell data with given ids from source to target. More...
 
int transferPointData (const std::vector< std::string > &arrayNames, const std::vector< std::string > &newnames=std::vector< std::string >())
 Transfer point data with given field names from source to target. More...
 
int transferCellData (const std::vector< std::string > &arrayNames, const std::vector< std::string > &newnames=std::vector< std::string >())
 Transfer cell data with given field names from source to target. More...
 
virtual int run (const std::vector< std::string > &newnames=std::vector< std::string >())=0
 Transfer all fields. More...
 
void setCheckQual (bool x)
 
void setContBool (bool x)
 

Protected Attributes

meshBasesource
 
meshBasetarget
 
vtkSmartPointer< vtkStaticCellLocator > srcCellLocator = nullptr
 
vtkSmartPointer< vtkStaticCellLocator > trgCellLocator = nullptr
 
vtkSmartPointer< vtkStaticPointLocator > srcPointLocator = nullptr
 
vtkSmartPointer< vtkStaticPointLocator > trgPointLocator = nullptr
 
bool checkQual
 
bool continuous
 
double c2cTrnsDistTol
 

Private Member Functions

std::vector< int > getArrayIDs (const std::vector< std::string > &arrayNames, vtkFieldData *fieldData)
 given array names and field data, return vector of corresponding array ids in the field data More...
 
int getDataArrayIndex (const std::string &arrayName, vtkFieldData *data)
 given array name and field data, returns index of array with given name (-1 if not found) More...
 

Inherited by ConservativeSurfaceTransfer, ConservativeVolumeTransfer, and FETransfer.

Constructor & Destructor Documentation

◆ TransferBase()

TransferBase::TransferBase ( )
inline

Definition at line 40 of file TransferBase.H.

41  : source(nullptr),
42  target(nullptr),
43  checkQual(false),
44  continuous(false),
45  c2cTrnsDistTol(1.e-6) {
46  std::cout << "TransferBase constructed" << std::endl;
47  }
double c2cTrnsDistTol
Definition: TransferBase.H:116
meshBase * target
Definition: TransferBase.H:106
meshBase * source
Definition: TransferBase.H:105

◆ ~TransferBase()

virtual TransferBase::~TransferBase ( )
inlinevirtual

Definition at line 49 of file TransferBase.H.

49  {
50  std::cout << "TransferBase destroyed" << std::endl;
51  }

Member Function Documentation

◆ getArrayIDs()

std::vector< int > TransferBase::getArrayIDs ( const std::vector< std::string > &  arrayNames,
vtkFieldData *  fieldData 
)
private

Definition at line 56 of file TransferBase.C.

References getDataArrayIndex().

Referenced by transferCellData(), and transferPointData().

57 {
58  std::vector<int> arrayIDs;
59  for(auto arrayName : arrayNames)
60  {
61  int arrayIndex = getDataArrayIndex(arrayName, fieldData);
62  if(arrayIndex == -1)
63  {
64  std::cerr << "Array " << arrayName << " not found." << std::endl;
65  exit(1);
66  }
67  arrayIDs.push_back(arrayIndex);
68  }
69  return arrayIDs;
70 }
int getDataArrayIndex(const std::string &arrayName, vtkFieldData *data)
given array name and field data, returns index of array with given name (-1 if not found) ...
Definition: TransferBase.C:72

◆ getDataArrayIndex()

int TransferBase::getDataArrayIndex ( const std::string &  arrayName,
vtkFieldData *  data 
)
private
Parameters
arrayNamename of array
datamesh field data that includes the sought after arrays
Returns
array id

Definition at line 72 of file TransferBase.C.

Referenced by getArrayIDs().

73 {
74  for(int arrayIndex = 0; arrayIndex < data->GetNumberOfArrays(); ++arrayIndex)
75  {
76  if(arrayName == data->GetArrayName(arrayIndex)) return arrayIndex;
77  }
78  return -1;
79 }
data_type data
Edge/face with sorted point ids (a, b, c, ...) is located at some index i in data[b], with data[b][i].first == [a, c] (for edges, third point id treated as -1).

◆ run()

virtual int TransferBase::run ( const std::vector< std::string > &  newnames = std::vector< std::string >())
pure virtual
Parameters
newnamesoptional array of names to be applied to transferred fields
Returns
0 upon completion

Implemented in FETransfer, ConservativeVolumeTransfer, and ConservativeSurfaceTransfer.

Referenced by FETransfer::CreateShared(), and ConservativeVolumeTransfer::transferCellData().

◆ setCheckQual()

void TransferBase::setCheckQual ( bool  x)
inline

Definition at line 100 of file TransferBase.H.

100 { checkQual = x; }

◆ setContBool()

void TransferBase::setContBool ( bool  x)
inline

Definition at line 102 of file TransferBase.H.

102 { continuous = x; }

◆ transferCellData() [1/2]

virtual int TransferBase::transferCellData ( const std::vector< int > &  arrayIDs,
const std::vector< std::string > &  newnames = std::vector< std::string >() 
)
pure virtual
Parameters
arrayIDsarray of array ids to specify which fields to transfer
newnamesoptional array of names to be applied to transferred fields
Returns
0 upon completion

Implemented in FETransfer, ConservativeVolumeTransfer, and ConservativeSurfaceTransfer.

Referenced by FETransfer::CreateShared(), and transferCellData().

◆ transferCellData() [2/2]

int TransferBase::transferCellData ( const std::vector< std::string > &  arrayNames,
const std::vector< std::string > &  newnames = std::vector<std::string>() 
)
Parameters
arrayIDsarray of array names to specify which fields to transfer
newnamesoptional array of names to be applied to transferred fields
Returns
0 upon completion

Definition at line 47 of file TransferBase.C.

References getArrayIDs(), meshBase::getDataSet(), source, and transferCellData().

49 {
50  vtkCellData* cellData = source->getDataSet()->GetCellData();
51  std::vector<int> arrayIDs = getArrayIDs(arrayNames, cellData);
52 
53  return transferCellData(arrayIDs, newnames);
54 }
virtual int transferCellData(const std::vector< int > &arrayIDs, const std::vector< std::string > &newnames=std::vector< std::string >())=0
Transfer cell data with given ids from source to target.
vtkSmartPointer< vtkDataSet > getDataSet() const
get this meshes&#39; dataSet
Definition: meshBase.H:308
std::vector< int > getArrayIDs(const std::vector< std::string > &arrayNames, vtkFieldData *fieldData)
given array names and field data, return vector of corresponding array ids in the field data ...
Definition: TransferBase.C:56
meshBase * source
Definition: TransferBase.H:105

◆ transferPointData() [1/2]

virtual int TransferBase::transferPointData ( const std::vector< int > &  arrayIDs,
const std::vector< std::string > &  newnames = std::vector< std::string >() 
)
pure virtual
Parameters
arrayIDsarray of array ids to specify which fields to transfer
newnamesoptional array of names to be applied to transferred fields
Returns
0 upon completion

Implemented in FETransfer, ConservativeVolumeTransfer, and ConservativeSurfaceTransfer.

Referenced by FETransfer::CreateShared(), ConservativeSurfaceTransfer::CreateShared(), ConservativeVolumeTransfer::CreateShared(), and transferPointData().

◆ transferPointData() [2/2]

int TransferBase::transferPointData ( const std::vector< std::string > &  arrayNames,
const std::vector< std::string > &  newnames = std::vector<std::string>() 
)
Parameters
arrayIDsarray of array names to specify which fields to transfer
newnamesoptional array of names to be applied to transferred fields
Returns
0 upon completion

Definition at line 36 of file TransferBase.C.

References getArrayIDs(), meshBase::getDataSet(), source, and transferPointData().

38 {
39  vtkPointData* pointData = source->getDataSet()->GetPointData();
40  std::vector<int> arrayIDs = getArrayIDs(arrayNames, pointData);
41 
42  return transferPointData(arrayIDs, newnames);
43 }
vtkSmartPointer< vtkDataSet > getDataSet() const
get this meshes&#39; dataSet
Definition: meshBase.H:308
virtual int transferPointData(const std::vector< int > &arrayIDs, const std::vector< std::string > &newnames=std::vector< std::string >())=0
Transfer point data with given ids from source to target.
std::vector< int > getArrayIDs(const std::vector< std::string > &arrayNames, vtkFieldData *fieldData)
given array names and field data, return vector of corresponding array ids in the field data ...
Definition: TransferBase.C:56
meshBase * source
Definition: TransferBase.H:105

Member Data Documentation

◆ c2cTrnsDistTol

double TransferBase::c2cTrnsDistTol
protected

Definition at line 116 of file TransferBase.H.

Referenced by FETransfer::transferCellData().

◆ checkQual

bool TransferBase::checkQual
protected

Definition at line 114 of file TransferBase.H.

Referenced by FETransfer::transferPointData().

◆ continuous

bool TransferBase::continuous
protected

Definition at line 115 of file TransferBase.H.

Referenced by FETransfer::transferCellData().

◆ source

◆ srcCellLocator

vtkSmartPointer<vtkStaticCellLocator> TransferBase::srcCellLocator = nullptr
protected

Definition at line 108 of file TransferBase.H.

Referenced by FETransfer::FETransfer(), and FETransfer::getClosestSourceCell().

◆ srcPointLocator

vtkSmartPointer<vtkStaticPointLocator> TransferBase::srcPointLocator = nullptr
protected

Definition at line 111 of file TransferBase.H.

Referenced by FETransfer::FETransfer().

◆ target

◆ trgCellLocator

vtkSmartPointer<vtkStaticCellLocator> TransferBase::trgCellLocator = nullptr
protected

Definition at line 109 of file TransferBase.H.

Referenced by FETransfer::FETransfer(), and FETransfer::getClosestSourceCell().

◆ trgPointLocator

vtkSmartPointer<vtkStaticPointLocator> TransferBase::trgPointLocator = nullptr
protected

Definition at line 112 of file TransferBase.H.

Referenced by FETransfer::FETransfer().


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