v0.14.0
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
MoFEM::PostProcBrokenMeshInMoabBase< E > Struct Template Reference

#include <src/post_proc/PostProcBrokenMeshInMoabBase.hpp>

Inheritance diagram for MoFEM::PostProcBrokenMeshInMoabBase< E >:
[legend]
Collaboration diagram for MoFEM::PostProcBrokenMeshInMoabBase< E >:
[legend]

Public Member Functions

 PostProcBrokenMeshInMoabBase (MoFEM::Interface &m_field, std::string opts_prefix="")
 
 PostProcBrokenMeshInMoabBase (MoFEM::Interface &m_field, boost::shared_ptr< moab::Core > core_mesh_ptr, std::string opts_prefix="")
 
virtual ~PostProcBrokenMeshInMoabBase ()
 
auto & getMapGaussPts ()
 Get vector of vectors associated to integration points. More...
 
auto & getPostProcMesh ()
 Get postprocessing mesh. More...
 
auto & getPostProcElements ()
 Get postprocessing elements. More...
 
MoFEMErrorCode writeFile (const std::string file_name)
 wrote results in (MOAB) format, use "file_name.h5m" More...
 
MoFEMErrorCode setTagsToTransfer (std::vector< Tag > tags_to_transfer)
 Set tags to be transferred to post-processing mesh. More...
 

Protected Member Functions

MoFEMErrorCode setGaussPts (int order)
 
MoFEMErrorCode preProcess ()
 Generate vertices and elements. More...
 
MoFEMErrorCode postProcess ()
 
int getRule (int order)
 
virtual int getMaxLevel () const
 Determine refinement level based on fields approx ordre. More...
 
virtual MoFEMErrorCode preProcPostProc ()
 
virtual MoFEMErrorCode postProcPostProc ()
 
virtual MoFEMErrorCode transferTags ()
 

Protected Attributes

boost::shared_ptr< moab::Core > coreMeshPtr = boost::make_shared<moab::Core>()
 
std::vector< EntityHandlemapGaussPts
 
Range postProcElements
 
std::map< EntityType, PostProcGenerateRefMeshPtrrefElementsMap
 
std::vector< Tag > tagsToTransfer
 
std::string optionsPrefix = ""
 Prefix for options. More...
 

Detailed Description

template<typename E>
struct MoFEM::PostProcBrokenMeshInMoabBase< E >

Examples
adolc_plasticity.cpp, and nonlinear_elastic.cpp.

Definition at line 94 of file PostProcBrokenMeshInMoabBase.hpp.

Constructor & Destructor Documentation

◆ PostProcBrokenMeshInMoabBase() [1/2]

template<typename E >
MoFEM::PostProcBrokenMeshInMoabBase< E >::PostProcBrokenMeshInMoabBase ( MoFEM::Interface m_field,
std::string  opts_prefix = "" 
)

Definition at line 267 of file PostProcBrokenMeshInMoabBase.hpp.

269  : E(m_field), optionsPrefix(opts_prefix) {}

◆ PostProcBrokenMeshInMoabBase() [2/2]

template<typename E >
MoFEM::PostProcBrokenMeshInMoabBase< E >::PostProcBrokenMeshInMoabBase ( MoFEM::Interface m_field,
boost::shared_ptr< moab::Core >  core_mesh_ptr,
std::string  opts_prefix = "" 
)

Definition at line 272 of file PostProcBrokenMeshInMoabBase.hpp.

275  : PostProcBrokenMeshInMoabBase(m_field, opts_prefix) {
276  coreMeshPtr = core_mesh_ptr;
277 }

◆ ~PostProcBrokenMeshInMoabBase()

Definition at line 280 of file PostProcBrokenMeshInMoabBase.hpp.

280  {
281  ParallelComm *pcomm_post_proc_mesh =
282  ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
283  if (pcomm_post_proc_mesh != NULL)
284  delete pcomm_post_proc_mesh;
285 }

Member Function Documentation

◆ getMapGaussPts()

template<typename E >
auto & MoFEM::PostProcBrokenMeshInMoabBase< E >::getMapGaussPts
inline

Get vector of vectors associated to integration points.

Returns
std::vector<EntityHandle>&

Definition at line 637 of file PostProcBrokenMeshInMoabBase.hpp.

637  {
638  return mapGaussPts;
639 }

◆ getMaxLevel()

template<typename E >
int MoFEM::PostProcBrokenMeshInMoabBase< E >::getMaxLevel
protectedvirtual

Determine refinement level based on fields approx ordre.

level = (order - 1) / 2

Returns
int

Definition at line 193 of file PostProcBrokenMeshInMoabBase.hpp.

193  {
194  auto get_element_max_dofs_order = [&]() {
195  int max_order = 0;
196  auto dofs_vec = E::getDataVectorDofsPtr();
197  for (auto &dof : *dofs_vec) {
198  const int dof_order = dof->getDofOrder();
199  max_order = (max_order < dof_order) ? dof_order : max_order;
200  };
201  return max_order;
202  };
203  const auto dof_max_order = get_element_max_dofs_order();
204  return (dof_max_order > 0) ? (dof_max_order - 1) / 2 : 0;
205 };

◆ getPostProcElements()

template<typename E >
auto & MoFEM::PostProcBrokenMeshInMoabBase< E >::getPostProcElements
inline

Get postprocessing elements.

Returns
auto&

Definition at line 649 of file PostProcBrokenMeshInMoabBase.hpp.

649  {
650  return postProcElements;
651 }

◆ getPostProcMesh()

template<typename E >
auto & MoFEM::PostProcBrokenMeshInMoabBase< E >::getPostProcMesh
inline

Get postprocessing mesh.

Returns
moab::Interface&

Definition at line 641 of file PostProcBrokenMeshInMoabBase.hpp.

641  {
642  if (!coreMeshPtr)
643  CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Core mesh not set");
644  moab::Interface &post_proc_mesh_interface = *coreMeshPtr;
645  return post_proc_mesh_interface;
646 }

◆ getRule()

template<typename E >
int MoFEM::PostProcBrokenMeshInMoabBase< E >::getRule ( int  order)
protected

Definition at line 287 of file PostProcBrokenMeshInMoabBase.hpp.

287  {
288  return -1;
289 };

◆ postProcess()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::postProcess
protected

Definition at line 623 of file PostProcBrokenMeshInMoabBase.hpp.

623  {
625 
626  ParallelComm *pcomm_post_proc_mesh =
627  ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
628  if(pcomm_post_proc_mesh == nullptr)
629  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "PComm not allocated");
630 
632  CHKERR pcomm_post_proc_mesh->resolve_shared_ents(0);
633 
635 }

◆ postProcPostProc()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::postProcPostProc
protectedvirtual

Definition at line 534 of file PostProcBrokenMeshInMoabBase.hpp.

534  {
536 
537  auto update_elements = [&]() {
538  ReadUtilIface *iface;
539  CHKERR getPostProcMesh().query_interface(iface);
541 
542  Range ents;
543  for (auto &m : refElementsMap) {
544  if (m.second->nbEles) {
545  MOFEM_TAG_AND_LOG("SELF", Sev::noisy, "PostProc")
546  << "Update < " << moab::CN::EntityTypeName(m.first)
547  << " number of processed " << m.second->countEle;
548  CHKERR iface->update_adjacencies(
549  m.second->startingEleHandle, m.second->countEle,
550  m.second->levelRef[0].size2(), m.second->eleConn);
551  ents.merge(Range(m.second->startingEleHandle,
552  m.second->startingEleHandle + m.second->countEle - 1));
553  }
554  }
555 
557  };
558 
559  auto remove_obsolete_entities = [&]() {
561  Range ents, adj;
562  for (auto &m : refElementsMap) {
563  if (m.second->nbEles) {
564  ents.merge(Range(m.second->startingEleHandle,
565  m.second->startingEleHandle + m.second->countEle - 1));
566  const int dim = moab::CN::Dimension(m.first);
567  for (auto d = 1; d != dim; ++d) {
568  CHKERR getPostProcMesh().get_adjacencies(ents, d, false, adj,
569  moab::Interface::UNION);
570  }
571  }
572  }
573  CHKERR getPostProcMesh().delete_entities(adj);
575  };
576 
577  auto set_proc_tags = [&]() {
579  ParallelComm *pcomm_post_proc_mesh =
580  ParallelComm::get_pcomm(&(getPostProcMesh()), MYPCOMM_INDEX);
581  if (pcomm_post_proc_mesh) {
582  Range ents;
583  for (auto &m : refElementsMap) {
584  if (m.second->nbEles) {
585  ents.merge(
586  Range(m.second->startingEleHandle,
587  m.second->startingEleHandle + m.second->countEle - 1));
588  }
589  }
590 
591  int rank = E::mField.get_comm_rank();
592  CHKERR getPostProcMesh().tag_clear_data(pcomm_post_proc_mesh->part_tag(),
593  ents, &rank);
594  }
596  };
597 
598  CHKERR update_elements();
599  CHKERR remove_obsolete_entities();
600  CHKERR set_proc_tags();
601 
603 }

◆ preProcess()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::preProcess
protected

Generate vertices and elements.

Returns
MoFEMErrorCode

Definition at line 606 of file PostProcBrokenMeshInMoabBase.hpp.

606  {
608 
609  ParallelComm *pcomm_post_proc_mesh =
610  ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
611  if (pcomm_post_proc_mesh != NULL)
612  delete pcomm_post_proc_mesh;
613  CHKERR getPostProcMesh().delete_mesh();
614  pcomm_post_proc_mesh =
615  new ParallelComm(&(getPostProcMesh()), PETSC_COMM_WORLD);
616 
618 
620 }

◆ preProcPostProc()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::preProcPostProc
protectedvirtual

Definition at line 419 of file PostProcBrokenMeshInMoabBase.hpp.

419  {
421 
422 auto get_ref_ele = [&](const EntityType type) {
423  PostProcGenerateRefMeshPtr ref_ele_ptr;
424 
425  auto it = refElementsMap.find(type);
426  if (it != refElementsMap.end()) {
427  ref_ele_ptr = it->second;
428  } else {
429  switch (type) {
430  case MBTET:
431  ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBTET>>();
432  break;
433  case MBHEX:
434  ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBHEX>>();
435  break;
436  case MBTRI:
437  ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBTRI>>();
438  break;
439  case MBQUAD:
440  ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBQUAD>>();
441  break;
442  case MBEDGE:
443  ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBEDGE>>();
444  break;
445  default:
446  MOFEM_LOG("SELF", Sev::error)
447  << "Generation of reference elements for type < "
448  << moab::CN::EntityTypeName(type) << " > is not implemented";
449  CHK_THROW_MESSAGE(MOFEM_NOT_IMPLEMENTED, "Element not implemented");
450  }
451 
452  CHK_THROW_MESSAGE(ref_ele_ptr->getOptions(optionsPrefix), "getOptions");
453  CHK_THROW_MESSAGE(ref_ele_ptr->generateReferenceElementMesh(),
454  "Error when generating reference element");
455 
456  refElementsMap[type] = ref_ele_ptr;
457  }
458 
459  return ref_ele_ptr;
460 };
461 
462  auto fe_ptr = this->problemPtr->numeredFiniteElementsPtr;
463 
464  auto miit =
465  fe_ptr->template get<Composite_Name_And_Part_mi_tag>().lower_bound(
466  boost::make_tuple(this->getFEName(), this->getLoFERank()));
467  auto hi_miit =
468  fe_ptr->template get<Composite_Name_And_Part_mi_tag>().upper_bound(
469  boost::make_tuple(this->getFEName(), this->getHiFERank()));
470 
471  const int number_of_ents_in_the_loop = this->getLoopSize();
472  if (std::distance(miit, hi_miit) != number_of_ents_in_the_loop) {
473  SETERRQ(E::mField.get_comm(), MOFEM_DATA_INCONSISTENCY,
474  "Wrong size of indicies. Inconsistent size number of iterated "
475  "elements iterated by problem and from range.");
476  }
477 
478  for (auto &m : refElementsMap) {
479  m.second->nbVertices = 0;
480  m.second->nbEles = 0;
481  m.second->countEle = 0;
482  m.second->countVertEle = 0;
483  }
484 
485  for (; miit != hi_miit; ++miit) {
486  auto type = (*miit)->getEntType();
487  auto ref_ele = get_ref_ele(type);
488 
489  // Set pointer to element. So that getDataVectorDofsPtr in getMaxLevel
490  // can work
491  E::numeredEntFiniteElementPtr = *miit;
492  bool add = true;
493  if (E::exeTestHook) {
494  add = E::exeTestHook(this);
495  }
496 
497  if (add) {
498  size_t level = getMaxLevel();
499  level = std::min(level, ref_ele->levelGaussPtsOnRefMesh.size() - 1);
500  ref_ele->nbVertices += ref_ele->levelGaussPtsOnRefMesh[level].size2();
501  ref_ele->nbEles += ref_ele->levelRef[level].size1();
502  }
503  }
504 
505  auto alloc_vertices_and_elements_on_post_proc_mesh = [&]() {
507 
508  ReadUtilIface *iface;
509  CHKERR getPostProcMesh().query_interface(iface);
510 
511  for (auto &m : refElementsMap) {
512  if (m.second->nbEles) {
513  CHKERR iface->get_node_coords(3, m.second->nbVertices, 0,
514  m.second->startingVertEleHandle,
515  m.second->verticesOnEleArrays);
516  CHKERR iface->get_element_connect(
517  m.second->nbEles, m.second->levelRef[0].size2(), m.first, 0,
518  m.second->startingEleHandle, m.second->eleConn);
519 
520  m.second->countEle = 0;
521  m.second->countVertEle = 0;
522  }
523  }
524 
526  };
527 
528  CHKERR alloc_vertices_and_elements_on_post_proc_mesh();
529 
531 }

◆ setGaussPts()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::setGaussPts ( int  order)
protected

Definition at line 292 of file PostProcBrokenMeshInMoabBase.hpp.

292  {
294 
295  auto type = type_from_handle(this->getFEEntityHandle());
296 
298 
299  try {
300  ref_ele = refElementsMap.at(type);
301  } catch (const out_of_range &e) {
302  SETERRQ1(
303  PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
304  "Generation of reference elements for type <%s> is not implemented",
305  moab::CN::EntityTypeName(type));
306  }
307 
308  auto set_gauss_pts = [&](auto &level_gauss_pts_on_ref_mesh, auto &level_ref,
309  auto &level_shape_functions,
310 
311  auto start_vert_handle, auto start_ele_handle,
312  auto &verts_array, auto &conn, auto &ver_count,
313  auto &ele_count
314 
315  ) {
317 
318  size_t level = getMaxLevel();
319  level = std::min(level, level_gauss_pts_on_ref_mesh.size() - 1);
320 
321  auto &level_ref_gauss_pts = level_gauss_pts_on_ref_mesh[level];
322  auto &level_ref_ele = level_ref[level];
323  auto &shape_functions = level_shape_functions[level];
324  E::gaussPts.resize(level_ref_gauss_pts.size1(), level_ref_gauss_pts.size2(),
325  false);
326  noalias(E::gaussPts) = level_ref_gauss_pts;
327 
328  const auto fe_ent = E::numeredEntFiniteElementPtr->getEnt();
329  auto get_fe_coords = [&]() {
330  const EntityHandle *conn;
331  int num_nodes;
333  E::mField.get_moab().get_connectivity(fe_ent, conn, num_nodes, true),
334  "error get connectivity");
335  VectorDouble coords(num_nodes * 3);
337  E::mField.get_moab().get_coords(conn, num_nodes, &*coords.begin()),
338  "error get coordinates");
339  return coords;
340  };
341 
342  auto coords = get_fe_coords();
343 
344  const int num_nodes = level_ref_gauss_pts.size2();
345  mapGaussPts.resize(level_ref_gauss_pts.size2());
346 
349  &*shape_functions.data().begin());
351  &verts_array[0][ver_count], &verts_array[1][ver_count],
352  &verts_array[2][ver_count]);
353  for (int gg = 0; gg != num_nodes; ++gg, ++ver_count) {
354 
355  mapGaussPts[gg] = start_vert_handle + ver_count;
356 
357  auto set_float_precision = [](const double x) {
358  if (std::abs(x) < std::numeric_limits<float>::epsilon())
359  return 0.;
360  else
361  return x;
362  };
363 
364  t_coords(i) = 0;
365  auto t_ele_coords = getFTensor1FromArray<3, 3>(coords);
366  for (int nn = 0; nn != CN::VerticesPerEntity(type); ++nn) {
367  t_coords(i) += t_n * t_ele_coords(i);
368  ++t_ele_coords;
369  ++t_n;
370  }
371 
372  for (auto ii : {0, 1, 2})
373  t_coords(ii) = set_float_precision(t_coords(ii));
374 
375  ++t_coords;
376  }
377 
378  Tag th;
379  int def_in_the_loop = -1;
380  CHKERR getPostProcMesh().tag_get_handle(
381  "NB_IN_THE_LOOP", 1, MB_TYPE_INTEGER, th, MB_TAG_CREAT | MB_TAG_SPARSE,
382  &def_in_the_loop);
383 
384  postProcElements.clear();
385  const int num_el = level_ref_ele.size1();
386  const int num_nodes_on_ele = level_ref_ele.size2();
387  auto start_e = start_ele_handle + ele_count;
388  postProcElements = Range(start_e, start_e + num_el - 1);
389  for (auto tt = 0; tt != level_ref_ele.size1(); ++tt, ++ele_count) {
390  for (int nn = 0; nn != num_nodes_on_ele; ++nn) {
391  conn[num_nodes_on_ele * ele_count + nn] =
392  mapGaussPts[level_ref_ele(tt, nn)];
393  }
394  }
395 
396  const int n_in_the_loop = E::nInTheLoop;
397  CHKERR getPostProcMesh().tag_clear_data(th, postProcElements,
398  &n_in_the_loop);
400 
402  };
403 
404  CHKERR set_gauss_pts(
405 
406  ref_ele->levelGaussPtsOnRefMesh, ref_ele->levelRef,
407  ref_ele->levelShapeFunctions,
408 
409  ref_ele->startingVertEleHandle, ref_ele->startingEleHandle,
410  ref_ele->verticesOnEleArrays, ref_ele->eleConn, ref_ele->countVertEle,
411  ref_ele->countEle
412 
413  );
414 
416 };

◆ setTagsToTransfer()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::setTagsToTransfer ( std::vector< Tag >  tags_to_transfer)

Set tags to be transferred to post-processing mesh.

Parameters
tags_to_transfer
Returns
MoFEMErrorCode

Definition at line 186 of file PostProcBrokenMeshInMoabBase.hpp.

187  {
189  tagsToTransfer.swap(tags_to_transfer);
191 }

◆ transferTags()

template<typename E >
MoFEMErrorCode MoFEM::PostProcBrokenMeshInMoabBase< E >::transferTags
protectedvirtual

transfer tags from mesh to post-process mesh

Definition at line 208 of file PostProcBrokenMeshInMoabBase.hpp.

208  {
210 
211  auto &calc_mesh = this->mField.get_moab();
212 
213  auto name = [&](auto tag) {
214  std::string name;
215  CHK_MOAB_THROW(calc_mesh.tag_get_name(tag, name), "get name");
216  return name;
217  };
218 
219  auto data_type = [&](auto tag) {
220  moab::DataType data_type;
221  CHK_MOAB_THROW(calc_mesh.tag_get_data_type(tag, data_type),
222  "get data type");
223  return data_type;
224  };
225 
226  auto type = [&](auto tag) {
227  moab::TagType type;
228  CHK_MOAB_THROW(calc_mesh.tag_get_type(tag, type), "get tag type");
229  return type;
230  };
231 
232  auto length = [&](auto tag) {
233  int length;
234  CHK_MOAB_THROW(calc_mesh.tag_get_length(tag, length), "get length ");
235  return length;
236  };
237 
238  auto default_value = [&](auto tag) {
239  const void *def_val;
240  int size;
241  CHK_MOAB_THROW(calc_mesh.tag_get_default_value(tag, def_val, size),
242  "get default tag value");
243  return def_val;
244  };
245 
246  auto data_ptr = [&](auto tag) {
247  const void *tag_data;
248  auto core_ent = this->getFEEntityHandle();
249  CHK_MOAB_THROW(calc_mesh.tag_get_by_ptr(tag, &core_ent, 1, &tag_data),
250  "get tag data");
251  return tag_data;
252  };
253 
254  for (auto tag : tagsToTransfer) {
255  Tag tag_postproc;
256  CHKERR getPostProcMesh().tag_get_handle(
257  name(tag).c_str(), length(tag), data_type(tag), tag_postproc,
258  type(tag) | MB_TAG_CREAT, default_value(tag));
259  CHKERR getPostProcMesh().tag_clear_data(tag_postproc, postProcElements,
260  data_ptr(tag));
261  }
262 
264 };

Member Data Documentation

◆ coreMeshPtr

template<typename E >
boost::shared_ptr<moab::Core> MoFEM::PostProcBrokenMeshInMoabBase< E >::coreMeshPtr = boost::make_shared<moab::Core>()
protected

Definition at line 165 of file PostProcBrokenMeshInMoabBase.hpp.

◆ mapGaussPts

template<typename E >
std::vector<EntityHandle> MoFEM::PostProcBrokenMeshInMoabBase< E >::mapGaussPts
protected

Definition at line 167 of file PostProcBrokenMeshInMoabBase.hpp.

◆ optionsPrefix

template<typename E >
std::string MoFEM::PostProcBrokenMeshInMoabBase< E >::optionsPrefix = ""
protected

Prefix for options.

Definition at line 176 of file PostProcBrokenMeshInMoabBase.hpp.

◆ postProcElements

template<typename E >
Range MoFEM::PostProcBrokenMeshInMoabBase< E >::postProcElements
protected

Definition at line 168 of file PostProcBrokenMeshInMoabBase.hpp.

◆ refElementsMap

template<typename E >
std::map<EntityType, PostProcGenerateRefMeshPtr> MoFEM::PostProcBrokenMeshInMoabBase< E >::refElementsMap
protected

Storing data about element types, and data for ReadUtilIface to create element entities on post-process mesh.

Definition at line 171 of file PostProcBrokenMeshInMoabBase.hpp.

◆ tagsToTransfer

template<typename E >
std::vector<Tag> MoFEM::PostProcBrokenMeshInMoabBase< E >::tagsToTransfer
protected

Set of tags on mesh to transfer to postprocessing mesh

Definition at line 174 of file PostProcBrokenMeshInMoabBase.hpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
CHK_MOAB_THROW
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:576
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
sdf_hertz.d
float d
Definition: sdf_hertz.py:5
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EntityHandle
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:596
E
MoFEM::th
Tag th
Definition: Projection10NodeCoordsOnField.cpp:122
MoFEM::PostProcBrokenMeshInMoabBase::postProcPostProc
virtual MoFEMErrorCode postProcPostProc()
Definition: PostProcBrokenMeshInMoabBase.hpp:534
MoFEM::PostProcBrokenMeshInMoabBase::mapGaussPts
std::vector< EntityHandle > mapGaussPts
Definition: PostProcBrokenMeshInMoabBase.hpp:167
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
MoFEM::PostProcBrokenMeshInMoabBase::refElementsMap
std::map< EntityType, PostProcGenerateRefMeshPtr > refElementsMap
Definition: PostProcBrokenMeshInMoabBase.hpp:171
MoFEM::PostProcBrokenMeshInMoabBase::transferTags
virtual MoFEMErrorCode transferTags()
Definition: PostProcBrokenMeshInMoabBase.hpp:208
MoFEM::PostProcBrokenMeshInMoabBase::coreMeshPtr
boost::shared_ptr< moab::Core > coreMeshPtr
Definition: PostProcBrokenMeshInMoabBase.hpp:165
MoFEM::PostProcBrokenMeshInMoabBase::tagsToTransfer
std::vector< Tag > tagsToTransfer
Definition: PostProcBrokenMeshInMoabBase.hpp:174
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::PostProcBrokenMeshInMoabBase::getMaxLevel
virtual int getMaxLevel() const
Determine refinement level based on fields approx ordre.
Definition: PostProcBrokenMeshInMoabBase.hpp:193
convert.type
type
Definition: convert.py:64
MoFEM::PostProcBrokenMeshInMoabBase::PostProcBrokenMeshInMoabBase
PostProcBrokenMeshInMoabBase(MoFEM::Interface &m_field, std::string opts_prefix="")
Definition: PostProcBrokenMeshInMoabBase.hpp:267
MoFEM::type_from_handle
auto type_from_handle(const EntityHandle h)
get type from entity handle
Definition: Templates.hpp:1869
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
FTensor::Index< 'i', 3 >
Range
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
FTensor::Tensor0
Definition: Tensor0.hpp:16
MoFEM::PostProcBrokenMeshInMoabBase::preProcPostProc
virtual MoFEMErrorCode preProcPostProc()
Definition: PostProcBrokenMeshInMoabBase.hpp:419
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
MoFEM::PostProcBrokenMeshInMoabBase::getPostProcMesh
auto & getPostProcMesh()
Get postprocessing mesh.
Definition: PostProcBrokenMeshInMoabBase.hpp:641
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MoFEM::PostProcBrokenMeshInMoabBase::postProcElements
Range postProcElements
Definition: PostProcBrokenMeshInMoabBase.hpp:168
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::PostProcBrokenMeshInMoabBase::optionsPrefix
std::string optionsPrefix
Prefix for options.
Definition: PostProcBrokenMeshInMoabBase.hpp:176
MoFEM::PostProcGenerateRefMeshPtr
boost::shared_ptr< PostProcGenerateRefMeshBase > PostProcGenerateRefMeshPtr
Definition: PostProcBrokenMeshInMoabBase.hpp:54