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.
ConservativeVolumeTransfer Class Reference

Detailed Description

Definition at line 44 of file ConservativeVolumeTransfer.H.

Public Member Functions

 ConservativeVolumeTransfer (meshBase *_source, meshBase *_target)
 
 ~ConservativeVolumeTransfer () override
 
int transferPointData (const std::vector< int > &arrayIDs, const std::vector< std::string > &newnames=std::vector< std::string >())
 Transfer point data with given ids from source to target. More...
 
int transferPointData (int i, vtkSmartPointer< vtkGenericCell > genCell, std::vector< vtkSmartPointer< vtkDoubleArray >> &dasSource, std::vector< vtkSmartPointer< vtkDoubleArray >> &dasTarget, bool flip)
 
int transferCellData (const std::vector< int > &arrayIDs, const std::vector< std::string > &newnames=std::vector< std::string >())
 Transfer cell data with given ids from source to target. More...
 
int transferCellData (int i, vtkSmartPointer< vtkGenericCell > genCell, std::vector< vtkSmartPointer< vtkDoubleArray >> &dasSourceToPoint, std::vector< vtkSmartPointer< vtkDoubleArray >> &dasTarget)
 
int run (const std::vector< std::string > &newnames=std::vector< std::string >()) override
 Transfer all fields. More...
 
int initialize ()
 
int transferSurfacePointData (const std::vector< int > &arrayIDs)
 
void enableSurfaceTransfer ()
 Enable surface transfer. More...
 
void disableSurfaceTransfer ()
 Disable surface transfer. More...
 
int constructSupermesh ()
 Construct supermesh. More...
 
int constructMassMatrix ()
 Construct mass matrix corresponding to target. More...
 
int constructMixedMassMatrix ()
 Construct mixed mass matrix corresponding to the source and target. More...
 
int constructLoadVector (int arrayId, int componentId)
 Construct load vector corresponding to a given array and component. More...
 
int transfer (int arrayId)
 Transfer array given its id. More...
 
int convertSupermeshToUnstructuredGrid ()
 Convertes supermesh to an unstructured grid. More...
 
vtkSmartPointer< vtkUnstructuredGrid > getSourceGrid ()
 Get the source grid in the transfer. More...
 
vtkSmartPointer< vtkUnstructuredGrid > getTargetGrid ()
 Get the target grid in the transfer. More...
 
vtkSmartPointer< vtkUnstructuredGrid > getSupermeshGrid ()
 Get the supermesh grid in the transfer. 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...
 
void setCheckQual (bool x)
 
void setContBool (bool x)
 

Static Public Member Functions

static ConservativeVolumeTransferCreate (meshBase *_source, meshBase *_target)
 
static std::shared_ptr< ConservativeVolumeTransferCreateShared (meshBase *_source, meshBase *_target)
 

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

void getLibSupermeshData (vtkDataSet *data, long &nnodes, int &dim, long &nelements, int &loc, double *&positions, long *&enlist, long &initTetId)
 
int interpolateCellDataToPoints (int cellArrayId)
 
vtkSmartPointer< vtkPolyData > extractSurface (vtkUnstructuredGrid *grid)
 
void volumeCheck ()
 

Private Attributes

double TET8 [32]
 
double TET8W [8]
 
double TET4 [16]
 
double TET4W [4]
 
double TET10 [40]
 
double TET10W [10]
 
long numSourceNodes
 
long numTargetNodes
 
long initSourceTetId
 
long initTargetTetId
 
vtkSmartPointer< vtkUnstructuredGrid > sourceGrid
 
vtkSmartPointer< vtkUnstructuredGrid > targetGrid
 
vtkSmartPointer< vtkUnstructuredGrid > supermeshGrid
 
vtkSmartPointer< vtkQuadratureSchemeDefinition > quadrature
 
Eigen::SparseMatrix< double, Eigen::ColMajor, long > massMatrix
 
Eigen::SparseMatrix< double, Eigen::ColMajor, long > mixedMassMatrix
 
Eigen::VectorXd loadVector
 
std::vector< std::vector< double > > children_a
 
std::vector< std::vector< double > > children_b
 
std::vector< std::vector< double > > tets_c
 
std::vector< long > parents_a
 
std::vector< long > parents_b
 
bool surfaceTransferEnabled = false
 

Inherits TransferBase.

Constructor & Destructor Documentation

◆ ConservativeVolumeTransfer()

ConservativeVolumeTransfer::ConservativeVolumeTransfer ( meshBase _source,
meshBase _target 
)

Definition at line 47 of file ConservativeVolumeTransfer.C.

References meshBase::getDataSet(), massMatrix, mixedMassMatrix, NEM::MSH::New(), numSourceNodes, numTargetNodes, quadrature, TransferBase::source, sourceGrid, TransferBase::target, targetGrid, TET4, and TET4W.

48  {
49  source = _source;
50  target = _target;
51 
52  // Set grids, allocate vectors/matrices, quadrature scheme and construct
53  // transfer objects
54  sourceGrid = vtkUnstructuredGrid::SafeDownCast(source->getDataSet());
55  targetGrid = vtkUnstructuredGrid::SafeDownCast(target->getDataSet());
56 
57  numSourceNodes = sourceGrid->GetNumberOfPoints();
58  numTargetNodes = targetGrid->GetNumberOfPoints();
59 
60  massMatrix = Eigen::SparseMatrix<double, Eigen::ColMajor, long>(
62  mixedMassMatrix = Eigen::SparseMatrix<double, Eigen::ColMajor, long>(
64 
65  // set quadrature scheme
67  quadrature->Initialize(VTK_TETRA, 4, 4, TET4, TET4W);
68  // quadrature->Initialize(VTK_TETRA, 4, 8, TET8, TET8W);
69 }
vtkSmartPointer< vtkUnstructuredGrid > targetGrid
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
meshBase * target
Definition: TransferBase.H:106
Eigen::SparseMatrix< double, Eigen::ColMajor, long > massMatrix
vtkSmartPointer< vtkDataSet > getDataSet() const
get this meshes&#39; dataSet
Definition: meshBase.H:308
Eigen::SparseMatrix< double, Eigen::ColMajor, long > mixedMassMatrix
vtkSmartPointer< vtkUnstructuredGrid > sourceGrid
vtkSmartPointer< vtkQuadratureSchemeDefinition > quadrature
meshBase * source
Definition: TransferBase.H:105

◆ ~ConservativeVolumeTransfer()

ConservativeVolumeTransfer::~ConservativeVolumeTransfer ( )
inlineoverride

Definition at line 49 of file ConservativeVolumeTransfer.H.

49  {
50  std::cout << "Conservative Volume Transfer destroyed" << std::endl;
51  }

Member Function Documentation

◆ constructLoadVector()

int ConservativeVolumeTransfer::constructLoadVector ( int  arrayId,
int  componentId 
)

Definition at line 421 of file ConservativeVolumeTransfer.C.

References nemAux::Timer::elapsed(), loadVector, mixedMassMatrix, numSourceNodes, sourceGrid, nemAux::Timer::start(), and nemAux::Timer::stop().

Referenced by transfer().

422  {
423  nemAux::Timer timer;
424 
425  std::cout << "Constructing load vector ..." << std::endl;
426  timer.start();
427 
428  vtkSmartPointer<vtkDoubleArray> sourceArray = vtkDoubleArray::SafeDownCast(
429  sourceGrid->GetPointData()->GetArray(arrayId));
430 
431  // allocate and populate source values vector
432  Eigen::VectorXd sourceValuesVector = Eigen::VectorXd(numSourceNodes);
433  auto sourcePointValues = vtkDoubleArray::SafeDownCast(
434  sourceGrid->GetPointData()->GetArray(arrayId));
435  for (int i = 0; i < numSourceNodes; ++i) {
436  sourceValuesVector(i) = sourcePointValues->GetComponent(i, componentId);
437  }
438 
439  loadVector = mixedMassMatrix * sourceValuesVector;
440  timer.stop();
441 
442  std::cout << "Load vector constructed in " << timer.elapsed()
443  << " milliseconds." << std::endl;
444 
445  return 0;
446 }
Eigen::SparseMatrix< double, Eigen::ColMajor, long > mixedMassMatrix
vtkSmartPointer< vtkUnstructuredGrid > sourceGrid

◆ constructMassMatrix()

int ConservativeVolumeTransfer::constructMassMatrix ( )

Definition at line 280 of file ConservativeVolumeTransfer.C.

References nemAux::Timer::elapsed(), initTargetTetId, massMatrix, quadrature, nemAux::Timer::start(), nemAux::Timer::stop(), and targetGrid.

Referenced by initialize(), and run().

280  {
281  nemAux::Timer timer;
282  std::cout << "Constructing mass matrix..." << std::endl;
283  timer.start();
284 
285  const double *weights = quadrature->GetQuadratureWeights();
286 
287  using Tripletd = Eigen::Triplet<double>;
288  std::vector<Tripletd> tripletList;
289 
290  // TODO start from initTargetTetId
291  for (long cellId = initTargetTetId; cellId < targetGrid->GetNumberOfCells();
292  ++cellId) {
293  if (targetGrid->GetCellType(cellId) != VTK_TETRA)
294  continue;
295  vtkTetra *tet = vtkTetra::SafeDownCast(targetGrid->GetCell(cellId));
296 
297  double detJ = vtkMeshQuality::TetVolume(tet);
298 
299  for (int qi = 0; qi < quadrature->GetNumberOfQuadraturePoints(); ++qi) {
300  const double *psi = quadrature->GetShapeFunctionWeights(qi);
301  double wi = weights[qi];
302 
303  for (int r = 0; r < 4; ++r) {
304  for (int s = 0; s < 4; ++s) {
305  // TODO check that point ids are zero indexed
306  long rGlob = tet->GetPointId(r);
307  long sGlob = tet->GetPointId(s);
308  tripletList.push_back(
309  Tripletd(rGlob, sGlob, wi * psi[r] * psi[s] * detJ));
310  }
311  }
312  }
313  }
314 
315  massMatrix.setFromTriplets(tripletList.begin(), tripletList.end());
316 
317  timer.stop();
318  std::cout << "Mass matrix constructed in " << timer.elapsed()
319  << " milliseconds." << std::endl;
320 
321  return 0;
322 }
vtkSmartPointer< vtkUnstructuredGrid > targetGrid
Eigen::SparseMatrix< double, Eigen::ColMajor, long > massMatrix
vtkSmartPointer< vtkQuadratureSchemeDefinition > quadrature

◆ constructMixedMassMatrix()

int ConservativeVolumeTransfer::constructMixedMassMatrix ( )

Requires a constructed supermesh.

Definition at line 324 of file ConservativeVolumeTransfer.C.

References nemAux::Timer::elapsed(), initSourceTetId, initTargetTetId, mixedMassMatrix, parents_a, parents_b, quadrature, sourceGrid, nemAux::Timer::start(), nemAux::Timer::stop(), targetGrid, and tets_c.

Referenced by initialize(), and run().

324  {
325  nemAux::Timer timer;
326  timer.start();
327  std::cout << "Constructing mixed mass matrix ..." << std::endl;
328 
329  const double *weights = quadrature->GetQuadratureWeights();
330 
331  using Tripletd = Eigen::Triplet<double>;
332  std::vector<Tripletd> tripletList(tets_c.size());
333 
334  for (long subCellIdx = 0; subCellIdx < tets_c.size(); ++subCellIdx) {
335  std::vector<double> childTet = tets_c[subCellIdx];
336 
337  long sourceCellId = parents_a[subCellIdx];
338  long targetCellId = parents_b[subCellIdx];
339 
340  // TODO check indexing for initSource/initTargetTetIds
341  vtkTetra *sourceParentTet = vtkTetra::SafeDownCast(
342  sourceGrid->GetCell(initSourceTetId + sourceCellId));
343  vtkTetra *targetParentTet = vtkTetra::SafeDownCast(
344  targetGrid->GetCell(initTargetTetId + targetCellId));
345 
346  double childVertices[4][3];
347  for (int i = 0; i < 4; ++i) {
348  int offs = i * 3;
349  childVertices[i][0] = childTet[offs];
350  childVertices[i][1] = childTet[offs + 1];
351  childVertices[i][2] = childTet[offs + 2];
352  }
353 
354  double sourceParentVertices[4][3];
355  double targetParentVertices[4][3];
356  for (int i = 0; i < 4; ++i) {
357  double x[3];
358 
359  sourceParentTet->GetPoints()->GetPoint(i, x);
360  for (int j = 0; j < 3; ++j)
361  sourceParentVertices[i][j] = x[j];
362 
363  targetParentTet->GetPoints()->GetPoint(i, x);
364  for (int j = 0; j < 3; ++j)
365  targetParentVertices[i][j] = x[j];
366  }
367 
368  // TODO check vol = detJ
369  double detJ;
370  libsupermesh_tetrahedron_volume(&childTet[0], &detJ);
371  // if(detJ < 1e-20) continue;
372 
373  for (int qi = 0; qi < quadrature->GetNumberOfQuadraturePoints(); ++qi) {
374  const double *psi = quadrature->GetShapeFunctionWeights(qi);
375  double wi = weights[qi];
376 
377  double x[3] = {0., 0., 0.};
378 
379  // get geometric coordinate
380  // TODO make sure x is contained in subtet
381  for (int i = 0; i < 4; ++i) {
382  for (int j = 0; j < 3; ++j)
383  x[j] += psi[i] * childVertices[i][j];
384  }
385 
386  // TODO check source/target psis are positive and sum to 1
387  // source shape function val
388  double srcPsi[4] = {};
389  vtkTetra::BarycentricCoords(
390  x, sourceParentVertices[0], sourceParentVertices[1],
391  sourceParentVertices[2], sourceParentVertices[3], srcPsi);
392 
393  // target shape function val
394  double tgtPsi[4] = {};
395  vtkTetra::BarycentricCoords(
396  x, targetParentVertices[0], targetParentVertices[1],
397  targetParentVertices[2], targetParentVertices[3], tgtPsi);
398 
399  for (int r = 0; r < 4; ++r) {
400  for (int s = 0; s < 4; ++s) {
401  // TODO check that point ids are zero indexed
402  long rGlobal = targetParentTet->GetPointId(r);
403  long sGlobal = sourceParentTet->GetPointId(s);
404  tripletList.push_back(
405  Tripletd(rGlobal, sGlobal, wi * tgtPsi[r] * srcPsi[s] * detJ));
406  }
407  }
408  }
409  }
410 
411  mixedMassMatrix.setFromTriplets(tripletList.begin(), tripletList.end());
412 
413  timer.stop();
414 
415  std::cout << "Mixed mass matrix constructed in " << timer.elapsed()
416  << " milliseconds." << std::endl;
417 
418  return 0;
419 }
vtkSmartPointer< vtkUnstructuredGrid > targetGrid
std::vector< std::vector< double > > tets_c
Eigen::SparseMatrix< double, Eigen::ColMajor, long > mixedMassMatrix
vtkSmartPointer< vtkUnstructuredGrid > sourceGrid
vtkSmartPointer< vtkQuadratureSchemeDefinition > quadrature

◆ constructSupermesh()

int ConservativeVolumeTransfer::constructSupermesh ( )

Definition at line 131 of file ConservativeVolumeTransfer.C.

References nemAux::Timer::elapsed(), meshBase::getDataSet(), getLibSupermeshData(), initSourceTetId, initTargetTetId, parents_a, parents_b, TransferBase::source, nemAux::Timer::start(), nemAux::Timer::stop(), TransferBase::target, and tets_c.

Referenced by initialize(), and run().

131  {
132  nemAux::Timer supermeshTimer;
133  std::cout << "Constructing supermesh..." << std::endl;
134  supermeshTimer.start();
135 
136  /*
137  * A - source
138  * B - target
139  */
140 
141  long nnodes_a;
142  int dim_a;
143  long nelements_a;
144  int loc_a;
145  double *positions_a;
146  long *enlist_a;
147 
148  nemAux::Timer sourceConversionTimer;
149  sourceConversionTimer.start();
150 
151  getLibSupermeshData(source->getDataSet(), nnodes_a, dim_a, nelements_a, loc_a,
152  positions_a, enlist_a, initSourceTetId);
153 
154  sourceConversionTimer.stop();
155  std::cout << "source converted to libsupermesh format in "
156  << sourceConversionTimer.elapsed() << " milliseconds." << std::endl;
157 
158  long nnodes_b;
159  int dim_b;
160  long nelements_b;
161  int loc_b;
162  double *positions_b;
163  long *enlist_b;
164 
165  nemAux::Timer targetConversionTimer;
166  targetConversionTimer.start();
167 
168  getLibSupermeshData(target->getDataSet(), nnodes_b, dim_b, nelements_b, loc_b,
169  positions_b, enlist_b, initTargetTetId);
170 
171  targetConversionTimer.stop();
172  std::cout << "target converted to libsupermesh format in "
173  << targetConversionTimer.elapsed() << " milliseconds." << std::endl;
174 
175  // volumeCheck();
176 
177  nemAux::Timer intersectionSearchTimer;
178  intersectionSearchTimer.start();
179 
180  libsupermesh_tree_intersection_finder_set_input(
181  &nnodes_a, &dim_a, &nelements_a, &loc_a, &nnodes_b, &dim_b, &nelements_b,
182  &loc_b, positions_a, enlist_a, positions_b, enlist_b);
183 
184  intersectionSearchTimer.stop();
185  std::cout << "intersection search completed in "
186  << intersectionSearchTimer.elapsed() << std::endl;
187 
188  nemAux::Timer intersectionCalcTimer;
189  intersectionCalcTimer.start();
190 
191  long nindices;
192  libsupermesh_tree_intersection_finder_query_output(&nindices);
193 
194  long nelements = nelements_a;
195  long *indices = new long[nindices];
196  long *ind_ptr = new long[nelements_a + 1];
197  libsupermesh_tree_intersection_finder_get_output(&nelements, &nindices,
198  indices, ind_ptr);
199 
200  auto getTet = [](long i, long *enlist, double *positions) -> double * {
201  double *tet = new double[3 * 4];
202  for (int j = 0; j < 4; ++j) {
203  // i*4 + j : index of node j in tet i
204  // multiply this by 3 to account for offset in positions array
205  long ptIndex = enlist[i * 4 + j] * 3;
206  double pt[3] = {positions[ptIndex], positions[ptIndex + 1],
207  positions[ptIndex + 2]};
208  int offset = j * 3;
209  tet[offset] = pt[0];
210  tet[offset + 1] = pt[1];
211  tet[offset + 2] = pt[2];
212  }
213  return tet;
214  };
215 
216  double *tets_c_buf = new double[1000];
217 
218  auto pushData = [this](std::vector<double> tet_c, long parent_a,
219  long parent_b) {
220  this->tets_c.push_back(tet_c);
221  this->parents_a.push_back(parent_a);
222  this->parents_b.push_back(parent_b);
223  };
224 
225  for (long ai = 0; ai < nelements_a; ++ai) // cells in A
226  {
227  double *tet_a = getTet(ai, enlist_a, positions_a);
228 
229  for (long bptr = ind_ptr[ai]; bptr < ind_ptr[ai + 1];
230  ++bptr) // indices of intersection candidates in B
231  {
232  long bj = indices[bptr]; // get index in b
233  double *tet_b = getTet(bj, enlist_b, positions_b);
234 
235  int n_tets_c;
236  // TODO if n_tets_c > buffer size -> reallocate
237  libsupermesh_intersect_tets_real(tet_a, tet_b, tets_c_buf, &n_tets_c);
238  if (!n_tets_c)
239  continue;
240 
241  for (int ck = 0; ck < n_tets_c; ++ck) {
242  int offs = ck * 12;
243 
244  std::vector<double> tet_c;
245  tet_c.insert(tet_c.end(), &tets_c_buf[offs], &tets_c_buf[offs + 12]);
246 
247  // push tet and its parent indices in A and B (source and target)
248  pushData(tet_c, ai, bj);
249  }
250  }
251  }
252 
253  intersectionCalcTimer.stop();
254  std::cout << "intersection calculation completed in "
255  << intersectionCalcTimer.elapsed() << " milliseconds." << std::endl;
256 
257  double totalVol = 0.;
258  for (auto tet : tets_c) {
259  double vol;
260  libsupermesh_tetrahedron_volume(&tet[0], &vol);
261  if (vol < 0)
262  std::cout << "NEGATIVE VOL" << std::endl;
263  totalVol += vol;
264  }
265  std::cout << std::setprecision(15) << "supermesh volume : " << totalVol
266  << std::endl;
267  std::cout << "number of elems : " << tets_c.size() << std::endl;
268 
269  // free
270  delete[] indices;
271  delete[] ind_ptr;
272  delete[] tets_c_buf;
273 
274  supermeshTimer.stop();
275  std::cout << "Supermesh constructed in " << supermeshTimer.elapsed()
276  << " milliseconds." << std::endl;
277  return 0;
278 }
std::vector< std::vector< double > > tets_c
void getLibSupermeshData(vtkDataSet *data, long &nnodes, int &dim, long &nelements, int &loc, double *&positions, long *&enlist, long &initTetId)
meshBase * target
Definition: TransferBase.H:106
vtkSmartPointer< vtkDataSet > getDataSet() const
get this meshes&#39; dataSet
Definition: meshBase.H:308
meshBase * source
Definition: TransferBase.H:105

◆ convertSupermeshToUnstructuredGrid()

int ConservativeVolumeTransfer::convertSupermeshToUnstructuredGrid ( )
See also
getSupermeshGrid

Definition at line 578 of file ConservativeVolumeTransfer.C.

References NEM::MSH::New(), sourceGrid, supermeshGrid, and tets_c.

578  {
580 
581  auto mergePoints = vtkSmartPointer<vtkMergePoints>::New();
582  auto supermeshGridPoints = vtkSmartPointer<vtkPoints>::New();
583  double bounds[6];
584  sourceGrid->GetBounds(bounds);
585 
586  mergePoints->InitPointInsertion(supermeshGridPoints, bounds, tets_c.size());
587 
588  // insert points
589  for (auto tet : tets_c) {
590  for (int ptIdx = 0; ptIdx < 4; ++ptIdx) {
591  int offs = ptIdx * 3;
592 
593  double x[3] = {tet[offs], tet[offs + 1], tet[offs + 2]};
594  vtkIdType ptId;
595  int inserted = mergePoints->InsertUniquePoint(x, ptId);
596  }
597  }
598 
599  supermeshGrid->SetPoints(supermeshGridPoints);
600 
601  auto tetArray = vtkSmartPointer<vtkCellArray>::New();
602  // add tets
603  for (auto tet : tets_c) {
604  auto tetra = vtkSmartPointer<vtkTetra>::New();
605 
606  Eigen::Vector3d faceVertices[3];
607  Eigen::Vector3d topVertex;
608 
609  vtkIdType ptIds[4];
610 
611  // get points of face and "top" vertex
612  for (int ptIdx = 0; ptIdx < 4; ++ptIdx) {
613  int offs = ptIdx * 3;
614 
615  double x[3] = {tet[offs], tet[offs + 1], tet[offs + 2]};
616 
617  if (ptIdx < 3)
618  faceVertices[ptIdx] = Eigen::Vector3d(x[0], x[1], x[2]);
619  else
620  topVertex = Eigen::Vector3d(x[0], x[1], x[2]);
621 
622  ptIds[ptIdx] = mergePoints->IsInsertedPoint(x);
623 
624  if (ptIds[ptIdx] < 0)
625  std::cout << "INVALID POINT" << std::endl;
626  }
627 
628  Eigen::Vector3d cross = (faceVertices[2] - faceVertices[0])
629  .cross(faceVertices[1] - faceVertices[0]);
630  cross = cross / cross.norm();
631 
632  Eigen::Vector3d topDirection = topVertex - faceVertices[0];
633  topDirection = topDirection / topDirection.norm();
634 
635  // if negative, reverse first 3 point ids
636  if (cross.dot(topDirection) > 0)
637  std::swap(ptIds[0], ptIds[2]);
638 
639  tetra->GetPointIds()->SetId(0, ptIds[0]);
640  tetra->GetPointIds()->SetId(1, ptIds[1]);
641  tetra->GetPointIds()->SetId(2, ptIds[2]);
642  tetra->GetPointIds()->SetId(3, ptIds[3]);
643 
644  // take first 3 - orient so that corresponding normal points to 4th
645  tetArray->InsertNextCell(tetra);
646  }
647 
648  supermeshGrid->SetCells(VTK_TETRA, tetArray);
649 
650  std::cout << "Converted grid has : " << supermeshGrid->GetNumberOfCells()
651  << " cells." << std::endl;
652 
653  if (supermeshGrid->GetNumberOfCells() != tets_c.size()) {
654  std::cout << "Reference grid has " << tets_c.size() << " cells."
655  << std::endl;
656  }
657 
658  return 0;
659 }
std::vector< std::vector< double > > tets_c
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
vtkSmartPointer< vtkUnstructuredGrid > supermeshGrid
vtkSmartPointer< vtkUnstructuredGrid > sourceGrid

◆ Create()

static ConservativeVolumeTransfer* ConservativeVolumeTransfer::Create ( meshBase _source,
meshBase _target 
)
inlinestatic

Definition at line 54 of file ConservativeVolumeTransfer.H.

Referenced by CreateShared().

55  {
56  return new ConservativeVolumeTransfer(_source, _target);
57  }
ConservativeVolumeTransfer(meshBase *_source, meshBase *_target)

◆ CreateShared()

static std::shared_ptr<ConservativeVolumeTransfer> ConservativeVolumeTransfer::CreateShared ( meshBase _source,
meshBase _target 
)
inlinestatic

Definition at line 60 of file ConservativeVolumeTransfer.H.

References Create(), and TransferBase::transferPointData().

Referenced by NEM::DRV::TransferDriver::CreateTransferObject().

60  {
61  return std::shared_ptr<ConservativeVolumeTransfer>(
62  ConservativeVolumeTransfer::Create(_source, _target));
63  }
static ConservativeVolumeTransfer * Create(meshBase *_source, meshBase *_target)

◆ disableSurfaceTransfer()

void ConservativeVolumeTransfer::disableSurfaceTransfer ( )
inline

◆ enableSurfaceTransfer()

void ConservativeVolumeTransfer::enableSurfaceTransfer ( )
inline

This is necessary for meshes with a shared boundary.

Definition at line 103 of file ConservativeVolumeTransfer.H.

◆ extractSurface()

vtkSmartPointer< vtkPolyData > ConservativeVolumeTransfer::extractSurface ( vtkUnstructuredGrid *  grid)
private

Definition at line 677 of file ConservativeVolumeTransfer.C.

References NEM::MSH::New().

677  {
678  // extract surface
679  auto surfaceFilter = vtkSmartPointer<vtkDataSetSurfaceFilter>::New();
680  surfaceFilter->SetPassThroughPointIds(true);
681  surfaceFilter->SetOriginalPointIdsName("original_point_ids");
682 
683  surfaceFilter->SetInputData(grid);
684  surfaceFilter->Update();
685 
686  vtkPolyData *surface = surfaceFilter->GetOutput();
687 
688  // polydata appears to be freed once dataset filter is freed, make deep copy
689  auto surfaceCopy = vtkSmartPointer<vtkPolyData>::New();
690  surfaceCopy->DeepCopy(surface);
691 
692  return surfaceCopy;
693 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.

◆ getLibSupermeshData()

void ConservativeVolumeTransfer::getLibSupermeshData ( vtkDataSet *  data,
long &  nnodes,
int &  dim,
long &  nelements,
int &  loc,
double *&  positions,
long *&  enlist,
long &  initTetId 
)
private

Definition at line 530 of file ConservativeVolumeTransfer.C.

Referenced by constructSupermesh().

532  {
533  initTetId = 0;
534  for (int cellId = 0; cellId < data->GetNumberOfCells(); ++cellId) {
535  if (data->GetCellType(cellId) == VTK_TETRA) {
536  initTetId = cellId;
537  break;
538  }
539  }
540 
541  long numNodes = data->GetNumberOfPoints();
542  // assume remaining cells are tets
543  long numTets = data->GetNumberOfCells() - initTetId;
544 
545  nnodes = numNodes;
546  dim = 3;
547  nelements = numTets;
548  loc = 4;
549 
550  // coordinates of nodes
551  positions = new double[dim * nnodes];
552  // element-node graph
553  enlist = new long[loc * nelements];
554 
555  // get coordinates
556  for (int ptIdx = 0; ptIdx < numNodes; ++ptIdx) {
557  double pt[3];
558  data->GetPoint(ptIdx, pt);
559 
560  int offset = ptIdx * 3;
561  positions[offset] = pt[0];
562  positions[offset + 1] = pt[1];
563  positions[offset + 2] = pt[2];
564  }
565 
566  // get element-node graph
567  for (int cellId = initTetId; cellId < data->GetNumberOfCells(); ++cellId) {
568  vtkCell *tet = data->GetCell(cellId);
569 
570  int offset = (cellId - initTetId) * 4;
571  enlist[offset] = tet->GetPointId(0);
572  enlist[offset + 1] = tet->GetPointId(1);
573  enlist[offset + 2] = tet->GetPointId(2);
574  enlist[offset + 3] = tet->GetPointId(3);
575  }
576 }
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).

◆ getSourceGrid()

vtkSmartPointer<vtkUnstructuredGrid> ConservativeVolumeTransfer::getSourceGrid ( )
inline

Definition at line 141 of file ConservativeVolumeTransfer.H.

141 { return sourceGrid; }
vtkSmartPointer< vtkUnstructuredGrid > sourceGrid

◆ getSupermeshGrid()

vtkSmartPointer<vtkUnstructuredGrid> ConservativeVolumeTransfer::getSupermeshGrid ( )
inline

Call convertSupermeshToUnstructuredGrid before this query.

Definition at line 151 of file ConservativeVolumeTransfer.H.

151  {
152  return supermeshGrid;
153  }
vtkSmartPointer< vtkUnstructuredGrid > supermeshGrid

◆ getTargetGrid()

vtkSmartPointer<vtkUnstructuredGrid> ConservativeVolumeTransfer::getTargetGrid ( )
inline

Definition at line 145 of file ConservativeVolumeTransfer.H.

145 { return targetGrid; }
vtkSmartPointer< vtkUnstructuredGrid > targetGrid

◆ initialize()

int ConservativeVolumeTransfer::initialize ( )

Definition at line 94 of file ConservativeVolumeTransfer.C.

References constructMassMatrix(), constructMixedMassMatrix(), and constructSupermesh().

94  {
98 }
int constructMassMatrix()
Construct mass matrix corresponding to target.
int constructMixedMassMatrix()
Construct mixed mass matrix corresponding to the source and target.
int constructSupermesh()
Construct supermesh.

◆ interpolateCellDataToPoints()

int ConservativeVolumeTransfer::interpolateCellDataToPoints ( int  cellArrayId)
private

Definition at line 488 of file ConservativeVolumeTransfer.C.

References NEM::MSH::New(), and sourceGrid.

488  {
489  vtkSmartPointer<vtkDataArray> sourceCellValues =
490  sourceGrid->GetCellData()->GetArray(cellArrayId);
491 
492  auto sourcePointValues = vtkSmartPointer<vtkDoubleArray>::New();
493  sourcePointValues->SetName(sourceCellValues->GetName());
494  sourcePointValues->SetNumberOfComponents(
495  sourceCellValues->GetNumberOfComponents());
496 
497  for (int pointId = 0; pointId < sourceGrid->GetNumberOfPoints(); ++pointId) {
498  auto incidentCellIds = vtkSmartPointer<vtkIdList>::New();
499  sourceGrid->GetPointCells(pointId, incidentCellIds);
500 
501  double *tuple = new double[sourceCellValues->GetNumberOfComponents()]();
502 
503  double volSum = 0.;
504 
505  // get contribution from each incident cell
506  for (int i = 0; i < incidentCellIds->GetNumberOfIds(); ++i) {
507  int cellId = incidentCellIds->GetId(i);
508  vtkTetra *tet = vtkTetra::SafeDownCast(sourceGrid->GetCell(cellId));
509 
510  double vol = vtkMeshQuality::TetVolume(tet);
511  volSum += vol;
512  // ... at each component
513  for (int compIdx = 0; compIdx < sourceCellValues->GetNumberOfComponents();
514  ++compIdx) {
515  tuple[compIdx] += vol * sourceCellValues->GetComponent(cellId, compIdx);
516  }
517  }
518  double volAvg = volSum / incidentCellIds->GetNumberOfIds();
519 
520  // divide by volume average
521  for (int compIdx = 0; compIdx < sourceCellValues->GetNumberOfComponents();
522  ++compIdx) {
523  tuple[compIdx] = tuple[compIdx] / volAvg;
524  }
525  }
526 
527  return 0;
528 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
vtkSmartPointer< vtkUnstructuredGrid > sourceGrid

◆ run()

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

Implements TransferBase.

Definition at line 71 of file ConservativeVolumeTransfer.C.

References constructMassMatrix(), constructMixedMassMatrix(), constructSupermesh(), sourceGrid, and transfer().

71  {
75 
76  // transfer all point data
77  vtkSmartPointer<vtkPointData> sourcePointData = sourceGrid->GetPointData();
78  for (int arrayIdx = 0; arrayIdx < sourcePointData->GetNumberOfArrays();
79  ++arrayIdx) {
80  transfer(arrayIdx);
81  }
82 
83  /*
84  vtkSmartPointer<vtkCellData> sourceCellData = sourceGrid->GetCellData();
85  for(int arrayIdx = 0; arrayIdx < sourceCellData->GetNumberOfArrays();
86  ++arrayIdx)
87  {
88  // TODO implement transferCellData
89  }
90  */
91  return 0;
92 }
int constructMassMatrix()
Construct mass matrix corresponding to target.
vtkSmartPointer< vtkUnstructuredGrid > sourceGrid
int constructMixedMassMatrix()
Construct mixed mass matrix corresponding to the source and target.
int transfer(int arrayId)
Transfer array given its id.
int constructSupermesh()
Construct supermesh.

◆ setCheckQual()

void TransferBase::setCheckQual ( bool  x)
inlineinherited

Definition at line 100 of file TransferBase.H.

100 { checkQual = x; }

◆ setContBool()

void TransferBase::setContBool ( bool  x)
inlineinherited

Definition at line 102 of file TransferBase.H.

102 { continuous = x; }

◆ transfer()

int ConservativeVolumeTransfer::transfer ( int  arrayId)

Definition at line 448 of file ConservativeVolumeTransfer.C.

References constructLoadVector(), nemAux::Timer::elapsed(), loadVector, massMatrix, NEM::MSH::New(), numTargetNodes, sourceGrid, nemAux::Timer::start(), nemAux::Timer::stop(), and targetGrid.

Referenced by run(), and transferPointData().

448  {
449  nemAux::Timer timer;
450  timer.start();
451  std::cout << "Transferring ..." << std::endl;
452 
453  vtkSmartPointer<vtkDataArray> sourceArray =
454  sourceGrid->GetPointData()->GetArray(arrayId);
455 
456  auto targetPointValues = vtkSmartPointer<vtkDoubleArray>::New();
457  targetPointValues->SetName(sourceArray->GetName());
458  targetPointValues->SetNumberOfComponents(
459  sourceArray->GetNumberOfComponents());
460 
461  Eigen::SparseLU<Eigen::SparseMatrix<double, Eigen::ColMajor>,
462  Eigen::COLAMDOrdering<int>>
463  solver;
464  solver.analyzePattern(massMatrix);
465  solver.factorize(massMatrix);
466 
467  for (int componentId = 0; componentId < sourceArray->GetNumberOfComponents();
468  ++componentId) {
469  constructLoadVector(arrayId, componentId);
470 
471  Eigen::VectorXd targetValuesVector = solver.solve(loadVector);
472 
473  for (int targetNodeId = 0; targetNodeId < numTargetNodes; ++targetNodeId) {
474  targetPointValues->InsertComponent(targetNodeId, componentId,
475  targetValuesVector(targetNodeId));
476  }
477  }
478 
479  targetGrid->GetPointData()->AddArray(targetPointValues);
480 
481  timer.stop();
482  std::cout << "Transfer completed in " << timer.elapsed() << " milliseconds."
483  << std::endl;
484 
485  return 0;
486 }
vtkSmartPointer< vtkUnstructuredGrid > targetGrid
int constructLoadVector(int arrayId, int componentId)
Construct load vector corresponding to a given array and component.
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
Eigen::SparseMatrix< double, Eigen::ColMajor, long > massMatrix
vtkSmartPointer< vtkUnstructuredGrid > sourceGrid

◆ transferCellData() [1/3]

int ConservativeVolumeTransfer::transferCellData ( const std::vector< int > &  arrayIDs,
const std::vector< std::string > &  newnames = std::vector< std::string >() 
)
inlinevirtual
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

Implements TransferBase.

Definition at line 78 of file ConservativeVolumeTransfer.H.

80  {
81  return 0;
82  }

◆ transferCellData() [2/3]

int ConservativeVolumeTransfer::transferCellData ( int  i,
vtkSmartPointer< vtkGenericCell >  genCell,
std::vector< vtkSmartPointer< vtkDoubleArray >> &  dasSourceToPoint,
std::vector< vtkSmartPointer< vtkDoubleArray >> &  dasTarget 
)
inline

Definition at line 84 of file ConservativeVolumeTransfer.H.

References TransferBase::run().

87  {
88  return 0;
89  }

◆ transferCellData() [3/3]

int TransferBase::transferCellData ( const std::vector< std::string > &  arrayNames,
const std::vector< std::string > &  newnames = std::vector<std::string>() 
)
inherited
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 TransferBase::getArrayIDs(), meshBase::getDataSet(), TransferBase::source, and TransferBase::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/3]

int ConservativeVolumeTransfer::transferPointData ( const std::vector< int > &  arrayIDs,
const std::vector< std::string > &  newnames = std::vector< std::string >() 
)
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

Implements TransferBase.

Definition at line 100 of file ConservativeVolumeTransfer.C.

References transfer().

102  {
103  /*
104  if(surfaceTransferEnabled)
105  {
106  // extract source surface
107  vtkSmartPointer<vtkPolyData> sourceSurface = extractSurface(sourceGrid);
108  vtkSmartPointer<vtkPolyData> targetSurface = extractSurface(targetGrid);
109 
110  vtkMesh* sourceMesh = new vtkMesh(sourceSurface, "sourceSurface.vtp");
111  vtkMesh* targetMesh = new vtkMesh(targetSurface, "targetSurface.vtp");
112 
113  sourceMesh->write("source_surface_test.vtp");
114  targetMesh->write("target_surface_test.vtp");
115 
116  ConservativeSurfaceTransfer* surfaceTransfer = new
117  ConservativeSurfaceTransfer(sourceMesh, targetMesh);
118 
119  surfaceTransfer->transferPointData(arrayIDs);
120 
121  targetMesh->write("target_surface_test.vtp");
122  }
123  */
124 
125  for (const int &arrayId : arrayIDs) {
126  transfer(arrayId);
127  }
128  return 0;
129 }
int transfer(int arrayId)
Transfer array given its id.

◆ transferPointData() [2/3]

int ConservativeVolumeTransfer::transferPointData ( int  i,
vtkSmartPointer< vtkGenericCell >  genCell,
std::vector< vtkSmartPointer< vtkDoubleArray >> &  dasSource,
std::vector< vtkSmartPointer< vtkDoubleArray >> &  dasTarget,
bool  flip 
)
inline

Definition at line 70 of file ConservativeVolumeTransfer.H.

73  {
74  return 0;
75  }

◆ transferPointData() [3/3]

int TransferBase::transferPointData ( const std::vector< std::string > &  arrayNames,
const std::vector< std::string > &  newnames = std::vector<std::string>() 
)
inherited
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 TransferBase::getArrayIDs(), meshBase::getDataSet(), TransferBase::source, and TransferBase::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

◆ transferSurfacePointData()

int ConservativeVolumeTransfer::transferSurfacePointData ( const std::vector< int > &  arrayIDs)

◆ volumeCheck()

void ConservativeVolumeTransfer::volumeCheck ( )
private

Definition at line 661 of file ConservativeVolumeTransfer.C.

References initSourceTetId, initTargetTetId, sourceGrid, and targetGrid.

661  {
662  for (int i = initSourceTetId; i < sourceGrid->GetNumberOfCells(); ++i) {
663  vtkTetra *tet = vtkTetra::SafeDownCast(sourceGrid->GetCell(i));
664  if (vtkMeshQuality::TetVolume(tet) < 0.) {
665  std::cout << "source cell " << i << " has negative volume" << std::endl;
666  }
667  }
668  for (int i = initTargetTetId; i < targetGrid->GetNumberOfCells(); ++i) {
669  vtkTetra *tet = vtkTetra::SafeDownCast(targetGrid->GetCell(i));
670  if (vtkMeshQuality::TetVolume(tet) < 0.) {
671  std::cout << "target cell " << i << " has negative volume" << std::endl;
672  }
673  }
674 }
vtkSmartPointer< vtkUnstructuredGrid > targetGrid
vtkSmartPointer< vtkUnstructuredGrid > sourceGrid

Member Data Documentation

◆ c2cTrnsDistTol

double TransferBase::c2cTrnsDistTol
protectedinherited

Definition at line 116 of file TransferBase.H.

Referenced by FETransfer::transferCellData().

◆ checkQual

bool TransferBase::checkQual
protectedinherited

Definition at line 114 of file TransferBase.H.

Referenced by FETransfer::transferPointData().

◆ children_a

std::vector<std::vector<double> > ConservativeVolumeTransfer::children_a
private

Definition at line 227 of file ConservativeVolumeTransfer.H.

◆ children_b

std::vector<std::vector<double> > ConservativeVolumeTransfer::children_b
private

Definition at line 228 of file ConservativeVolumeTransfer.H.

◆ continuous

bool TransferBase::continuous
protectedinherited

Definition at line 115 of file TransferBase.H.

Referenced by FETransfer::transferCellData().

◆ initSourceTetId

long ConservativeVolumeTransfer::initSourceTetId
private

◆ initTargetTetId

long ConservativeVolumeTransfer::initTargetTetId
private

◆ loadVector

Eigen::VectorXd ConservativeVolumeTransfer::loadVector
private

Definition at line 225 of file ConservativeVolumeTransfer.H.

Referenced by constructLoadVector(), and transfer().

◆ massMatrix

Eigen::SparseMatrix<double, Eigen::ColMajor, long> ConservativeVolumeTransfer::massMatrix
private

◆ mixedMassMatrix

Eigen::SparseMatrix<double, Eigen::ColMajor, long> ConservativeVolumeTransfer::mixedMassMatrix
private

◆ numSourceNodes

long ConservativeVolumeTransfer::numSourceNodes
private

◆ numTargetNodes

long ConservativeVolumeTransfer::numTargetNodes
private

Definition at line 211 of file ConservativeVolumeTransfer.H.

Referenced by ConservativeVolumeTransfer(), and transfer().

◆ parents_a

std::vector<long> ConservativeVolumeTransfer::parents_a
private

Definition at line 231 of file ConservativeVolumeTransfer.H.

Referenced by constructMixedMassMatrix(), and constructSupermesh().

◆ parents_b

std::vector<long> ConservativeVolumeTransfer::parents_b
private

Definition at line 232 of file ConservativeVolumeTransfer.H.

Referenced by constructMixedMassMatrix(), and constructSupermesh().

◆ quadrature

vtkSmartPointer<vtkQuadratureSchemeDefinition> ConservativeVolumeTransfer::quadrature
private

◆ source

◆ sourceGrid

vtkSmartPointer<vtkUnstructuredGrid> ConservativeVolumeTransfer::sourceGrid
private

◆ srcCellLocator

vtkSmartPointer<vtkStaticCellLocator> TransferBase::srcCellLocator = nullptr
protectedinherited

Definition at line 108 of file TransferBase.H.

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

◆ srcPointLocator

vtkSmartPointer<vtkStaticPointLocator> TransferBase::srcPointLocator = nullptr
protectedinherited

Definition at line 111 of file TransferBase.H.

Referenced by FETransfer::FETransfer().

◆ supermeshGrid

vtkSmartPointer<vtkUnstructuredGrid> ConservativeVolumeTransfer::supermeshGrid
private

Definition at line 218 of file ConservativeVolumeTransfer.H.

Referenced by convertSupermeshToUnstructuredGrid().

◆ surfaceTransferEnabled

bool ConservativeVolumeTransfer::surfaceTransferEnabled = false
private

Definition at line 246 of file ConservativeVolumeTransfer.H.

◆ target

◆ targetGrid

vtkSmartPointer<vtkUnstructuredGrid> ConservativeVolumeTransfer::targetGrid
private

◆ TET10

double ConservativeVolumeTransfer::TET10[40]
private
Initial value:
= {
.7784952948213300, .0738349017262234, .0738349017262234,
.0738349017262234, .0738349017262234, .7784952948213300,
.0738349017262234, .0738349017262234, .0738349017262234,
.0738349017262234, .7784952948213300, .0738349017262234,
.0738349017262234, .0738349017262234, .0738349017262234,
.7784952948213300,
.4062443438840510, .4062443438840510, .0937556561159491,
.0937556561159491, .4062443438840510, .0937556561159491,
.4062443438840510, .0937556561159491, .4062443438840510,
.0937556561159491, .0937556561159491, .4062443438840510,
.0937556561159491, .4062443438840510, .4062443438840510,
.0937556561159491, .0937556561159491, .4062443438840510,
.0937556561159491, .4062443438840510, .0937556561159491,
.0937556561159491, .4062443438840510, .4062443438840510,
}

Definition at line 183 of file ConservativeVolumeTransfer.H.

◆ TET10W

double ConservativeVolumeTransfer::TET10W[10]
private
Initial value:
= {
.0476331348432089, .0476331348432089, .0476331348432089,
.0476331348432089,
.1349112434378610, .1349112434378610, .1349112434378610,
.1349112434378610, .1349112434378610, .1349112434378610,
}

Definition at line 202 of file ConservativeVolumeTransfer.H.

◆ TET4

double ConservativeVolumeTransfer::TET4[16]
private
Initial value:
= {.5854101966249684, .1381966011250105, .1381966011250105,
.1381966011250105, .1381966011250105, .5854101966249684,
.1381966011250105, .1381966011250105, .1381966011250105,
.1381966011250105, .5854101966249684, .1381966011250105,
.1381966011250105, .1381966011250105, .1381966011250105,
.5854101966249684}

Definition at line 173 of file ConservativeVolumeTransfer.H.

Referenced by ConservativeVolumeTransfer().

◆ TET4W

double ConservativeVolumeTransfer::TET4W[4]
private
Initial value:
= {0.250000000000000, 0.250000000000000, 0.250000000000000,
0.250000000000000}

Definition at line 180 of file ConservativeVolumeTransfer.H.

Referenced by ConservativeVolumeTransfer().

◆ TET8

double ConservativeVolumeTransfer::TET8[32]
private
Initial value:
= {
0.3281633025163817, 0.3281633025163817, 0.3281633025163817,
0.01551009245085488, 0.1080472498984286, 0.1080472498984286,
0.1080472498984286, 0.6758582503047142, 0.3281633025163817,
0.3281633025163817, 0.01551009245085488, 0.3281633025163817,
0.1080472498984286, 0.1080472498984286, 0.6758582503047142,
0.1080472498984286, 0.3281633025163817, 0.01551009245085488,
0.3281633025163817, 0.3281633025163817, 0.1080472498984286,
0.6758582503047142, 0.1080472498984286, 0.1080472498984286,
0.01551009245085488, 0.3281633025163817, 0.3281633025163817,
0.3281633025163817, 0.6758582503047142, 0.1080472498984286,
0.1080472498984286, 0.1080472498984286}

Definition at line 156 of file ConservativeVolumeTransfer.H.

◆ TET8W

double ConservativeVolumeTransfer::TET8W[8]
private
Initial value:
= {0.1362178425370874, 0.1137821574629126, 0.1362178425370874,
0.1137821574629126, 0.1362178425370874, 0.1137821574629126,
0.1362178425370874, 0.1137821574629126}

Definition at line 169 of file ConservativeVolumeTransfer.H.

◆ tets_c

std::vector<std::vector<double> > ConservativeVolumeTransfer::tets_c
private

◆ trgCellLocator

vtkSmartPointer<vtkStaticCellLocator> TransferBase::trgCellLocator = nullptr
protectedinherited

Definition at line 109 of file TransferBase.H.

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

◆ trgPointLocator

vtkSmartPointer<vtkStaticPointLocator> TransferBase::trgPointLocator = nullptr
protectedinherited

Definition at line 112 of file TransferBase.H.

Referenced by FETransfer::FETransfer().


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