v0.14.0
Public Member Functions | Public Attributes | Private Attributes | List of all members
MoFEM::MatrixManager Struct Reference

Matrix manager is used to build and partition problems. More...

#include <src/interfaces/MatrixManager.hpp>

Inheritance diagram for MoFEM::MatrixManager:
[legend]
Collaboration diagram for MoFEM::MatrixManager:
[legend]

Public Member Functions

MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 MatrixManager (const MoFEM::Core &core)
 
virtual ~MatrixManager ()=default
 Destructor. More...
 
template<class Tag >
MoFEMErrorCode createMPIAIJWithArrays (const std::string name, Mat *Aij, int verb=QUIET)
 Creates a MPI AIJ matrix using arrays that contain in standard CSR format the local rows. More...
 
template<class Tag >
MoFEMErrorCode createMPIAIJWithArrays (const std::string name, SmartPetscObj< Mat > &aij_ptr, int verb=QUIET)
 Creates a MPI AIJ matrix using arrays that contain in standard CSR format the local rows. More...
 
template<class Tag >
MoFEMErrorCode createMPIAIJCUSPARSEWithArrays (const std::string name, Mat *Aij, int verb=QUIET)
 
template<class Tag >
MoFEMErrorCode createMPIAIJCUSPARSEWithArrays (const std::string name, SmartPetscObj< Mat > &aij_ptr, int verb=QUIET)
 
template<class Tag >
MoFEMErrorCode createSeqAIJCUSPARSEWithArrays (const std::string name, Mat *Aij, int verb=QUIET)
 
template<class Tag >
MoFEMErrorCode createSeqAIJCUSPARSEWithArrays (const std::string name, SmartPetscObj< Mat > &aij_ptr, int verb=QUIET)
 
template<class Tag >
MoFEMErrorCode createMPIAIJ (const std::string name, Mat *Aij, int verb=QUIET)
 Creates a MPI AIJ matrix using arrays that contain in standard CSR format the local rows. More...
 
template<class Tag >
MoFEMErrorCode createMPIAIJ (const std::string name, SmartPetscObj< Mat > &aij_ptr, int verb=QUIET)
 Creates a MPI AIJ matrix using arrays that contain in standard CSR format the local rows. More...
 
template<class Tag >
MoFEMErrorCode createMPIAdjWithArrays (const std::string name, Mat *Adj, int verb=QUIET)
 Creates a sparse matrix representing an adjacency list. More...
 
template<class Tag >
MoFEMErrorCode createSeqAIJWithArrays (const std::string name, Mat *Aij, int verb=QUIET)
 Create sequencial matrix. More...
 
template<class Tag >
MoFEMErrorCode createSeqAIJWithArrays (const std::string name, SmartPetscObj< Mat > &aij_ptr, int verb=QUIET)
 Create sequencial matrix. More...
 
MoFEMErrorCode checkMatrixFillIn (const std::string problem_name, int row_print, int col_print, Mat A, int verb=QUIET)
 check if matrix fill in correspond to finite element indices More...
 
template<class Tag >
MoFEMErrorCode checkMPIAIJWithArraysMatrixFillIn (const std::string problem_name, int row_print, int col_print, int verb=QUIET)
 check if matrix fill in correspond to finite element indices More...
 
template<class Tag >
MoFEMErrorCode checkMPIAIJMatrixFillIn (const std::string problem_name, int row_print, int col_print, int verb=QUIET)
 check if matrix fill in correspond to finite element indices More...
 
template<>
MoFEMErrorCode createMPIAIJWithArrays (const std::string name, Mat *Aij, int verb)
 
template<>
MoFEMErrorCode createMPIAIJCUSPARSEWithArrays (const std::string name, Mat *Aij, int verb)
 
template<>
MoFEMErrorCode createSeqAIJCUSPARSEWithArrays (const std::string name, Mat *Aij, int verb)
 
template<>
MoFEMErrorCode createMPIAIJ (const std::string name, Mat *Aij, int verb)
 
template<>
MoFEMErrorCode createMPIAdjWithArrays (const std::string name, Mat *Adj, int verb)
 
template<>
MoFEMErrorCode createSeqAIJWithArrays (const std::string name, Mat *Aij, int verb)
 
template<>
MoFEMErrorCode checkMPIAIJWithArraysMatrixFillIn (const std::string problem_name, int row_print, int col_print, int verb)
 
template<>
MoFEMErrorCode checkMPIAIJMatrixFillIn (const std::string problem_name, int row_print, int col_print, int verb)
 
template<>
MoFEMErrorCode createMPIAIJWithArrays (const std::string name, Mat *Aij, int verb)
 
template<>
MoFEMErrorCode createMPIAIJCUSPARSEWithArrays (const std::string name, Mat *Aij, int verb)
 
template<>
MoFEMErrorCode createMPIAIJ (const std::string name, Mat *Aij, int verb)
 
template<>
MoFEMErrorCode createMPIAdjWithArrays (const std::string name, Mat *Adj, int verb)
 
template<>
MoFEMErrorCode createSeqAIJWithArrays (const std::string name, Mat *Aij, int verb)
 
template<>
MoFEMErrorCode createSeqAIJCUSPARSEWithArrays (const std::string name, Mat *Aij, int verb)
 
template<>
MoFEMErrorCode checkMPIAIJWithArraysMatrixFillIn (const std::string problem_name, int row_print, int col_print, int verb)
 
template<>
MoFEMErrorCode checkMPIAIJMatrixFillIn (const std::string problem_name, int row_print, int col_print, int verb)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 

Public Attributes

MoFEM::CorecOre
 

Private Attributes

PetscLogEvent MOFEM_EVENT_createMPIAIJ
 
PetscLogEvent MOFEM_EVENT_createMPIAIJWithArrays
 
PetscLogEvent MOFEM_EVENT_createMPIAdjWithArrays
 
PetscLogEvent MOFEM_EVENT_createSeqAIJWithArrays
 
PetscLogEvent MOFEM_EVENT_createMPIAIJCUSPARSEWithArrays
 
PetscLogEvent MOFEM_EVENT_createSeqAIJCUSPARSEWithArrays
 
PetscLogEvent MOFEM_EVENT_checkMatrixFillIn
 

Additional Inherited Members

- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 

Detailed Description

Matrix manager is used to build and partition problems.

Examples
build_large_problem.cpp, build_problems.cpp, dm_build_partitioned_mesh.cpp, dm_create_subdm.cpp, dm_partitioned_no_field.cpp, nonlinear_dynamics.cpp, prism_polynomial_approximation.cpp, quad_polynomial_approximation.cpp, remove_entities_from_problem.cpp, and remove_entities_from_problem_not_partitioned.cpp.

Definition at line 21 of file MatrixManager.hpp.

Constructor & Destructor Documentation

◆ MatrixManager()

MoFEM::MatrixManager::MatrixManager ( const MoFEM::Core core)

Definition at line 597 of file MatrixManager.cpp.

598  : cOre(const_cast<MoFEM::Core &>(core)) {
599  PetscLogEventRegister("MatrixManagerCreateMPIAIJ", 0,
601  PetscLogEventRegister("MatrixManagerCreateMPIAIJWithArrays", 0,
603  PetscLogEventRegister("MatrixManagerCreateMPIAdjWithArrays", 0,
605  PetscLogEventRegister("MatrixManagerCreateMPIAIJCUSPARSEWithArrays", 0,
607  PetscLogEventRegister("MatrixManagerCreateSeqAIJCUSPARSEWithArrays", 0,
609  PetscLogEventRegister("MatrixManagerCreateSeqAIJWithArrays", 0,
611  PetscLogEventRegister("MatrixManagerCheckMPIAIJWithArraysMatrixFillIn", 0,
613 }

◆ ~MatrixManager()

virtual MoFEM::MatrixManager::~MatrixManager ( )
virtualdefault

Destructor.

Member Function Documentation

◆ checkMatrixFillIn()

MoFEMErrorCode MoFEM::MatrixManager::checkMatrixFillIn ( const std::string  problem_name,
int  row_print,
int  col_print,
Mat  A,
int  verb = QUIET 
)

check if matrix fill in correspond to finite element indices

This is used to check consistency of code. If problem is notices with additional non-zero elements in matrix, this function can help detect problem. Should be used as a part of atom tests

Parameters
problem_name
rowprint info at particular row
colprint info at particular col
Returns
MoFEMErrorCode

Definition at line 890 of file MatrixManager.cpp.

892  {
893  MoFEM::CoreInterface &m_field = cOre;
895 
896  PetscLogEventBegin(MOFEM_EVENT_checkMatrixFillIn, 0, 0, 0, 0);
897 
898  struct TestMatrixFillIn : public FEMethod {
899  CoreInterface *mFieldPtr;
900 
901  Mat A;
902 
903  int rowPrint, colPrint;
904 
905  TestMatrixFillIn(CoreInterface *m_field_ptr, Mat a, int row_print,
906  int col_print)
907  : mFieldPtr(m_field_ptr), A(a), rowPrint(row_print),
908  colPrint(col_print){};
909 
910  MoFEMErrorCode preProcess() {
913  }
914 
915  MoFEMErrorCode operator()() {
917 
918  if (refinedFiniteElementsPtr->find(
919  numeredEntFiniteElementPtr->getEnt()) ==
920  refinedFiniteElementsPtr->end()) {
921  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
922  "data inconsistency");
923  }
924 
925  auto row_dofs = getRowDofsPtr();
926  auto col_dofs = getColDofsPtr();
927 
928  for (auto cit = col_dofs->begin(); cit != col_dofs->end(); cit++) {
929 
930  if (refinedEntitiesPtr->find((*cit)->getEnt()) ==
931  refinedEntitiesPtr->end()) {
932  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
933  "data inconsistency");
934  }
935  if (!(*cit)->getActive()) {
936  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
937  "data inconsistency");
938  }
939 
940  FieldEntityEntFiniteElementAdjacencyMap_multiIndex::index<
941  Composite_Unique_mi_tag>::type::iterator ait;
942  ait = adjacenciesPtr->get<Composite_Unique_mi_tag>().find(
943  boost::make_tuple((*cit)->getFieldEntityPtr()->getLocalUniqueId(),
944  numeredEntFiniteElementPtr->getLocalUniqueId()));
945  if (ait == adjacenciesPtr->end()) {
946  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
947  "adjacencies data inconsistency");
948  } else {
949  UId uid = ait->getEntUniqueId();
950  if (entitiesPtr->find(uid) == entitiesPtr->end()) {
951  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
952  "data inconsistency");
953  }
954  if (dofsPtr->find((*cit)->getLocalUniqueId()) == dofsPtr->end()) {
955  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
956  "data inconsistency");
957  }
958  }
959 
960  if ((*cit)->getEntType() != MBVERTEX) {
961 
962  auto range =
963  col_dofs->get<Ent_mi_tag>().equal_range((*cit)->getEnt());
964  int nb_dofs_on_ent = std::distance(range.first, range.second);
965 
966  int max_order = (*cit)->getMaxOrder();
967  if ((*cit)->getNbOfCoeffs() * (*cit)->getOrderNbDofs(max_order) !=
968  nb_dofs_on_ent) {
969 
970  /* It could be that you have
971  removed DOFs from problem, and for example if this was vector filed
972  with components {Ux,Uy,Uz}, you removed on Uz element.*/
973 
974  MOFEM_LOG("SELF", Sev::warning)
975  << "Warning: Number of Dofs in Col different than number "
976  "of dofs for given entity order "
977  << (*cit)->getNbOfCoeffs() * (*cit)->getOrderNbDofs(max_order)
978  << " " << nb_dofs_on_ent;
979  }
980  }
981  }
982 
983  for (auto rit = row_dofs->begin(); rit != row_dofs->end(); rit++) {
984 
985  if (refinedEntitiesPtr->find((*rit)->getEnt()) ==
986  refinedEntitiesPtr->end()) {
987  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
988  "data inconsistency");
989  }
990  if (!(*rit)->getActive()) {
991  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
992  "data inconsistency");
993  }
994 
995  FieldEntityEntFiniteElementAdjacencyMap_multiIndex::index<
996  Composite_Unique_mi_tag>::type::iterator ait;
997  ait = adjacenciesPtr->get<Composite_Unique_mi_tag>().find(
998  boost::make_tuple((*rit)->getFieldEntityPtr()->getLocalUniqueId(),
999  numeredEntFiniteElementPtr->getLocalUniqueId()));
1000  if (ait == adjacenciesPtr->end()) {
1002  MOFEM_LOG("SELF", Sev::error) << *(*rit);
1003  MOFEM_LOG("SELF", Sev::error) << *(*rit);
1004  MOFEM_LOG("SELF", Sev::error) << *numeredEntFiniteElementPtr;
1005  MOFEM_LOG("SELF", Sev::error) << "dof: " << (*rit)->getBitRefLevel();
1006  MOFEM_LOG("SELF", Sev::error)
1007  << "fe: " << numeredEntFiniteElementPtr->getBitRefLevel();
1008  MOFEM_LOG("SELF", Sev::error)
1009  << "problem: " << problemPtr->getBitRefLevel();
1010  MOFEM_LOG("SELF", Sev::error)
1011  << "problem mask: " << problemPtr->getBitRefLevelMask();
1012  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
1013  "adjacencies data inconsistency");
1014  } else {
1015  UId uid = ait->getEntUniqueId();
1016  if (entitiesPtr->find(uid) == entitiesPtr->end()) {
1017  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
1018  "data inconsistency");
1019  }
1020  if (dofsPtr->find((*rit)->getLocalUniqueId()) == dofsPtr->end()) {
1021  SETERRQ(mFieldPtr->get_comm(), MOFEM_DATA_INCONSISTENCY,
1022  "data inconsistency");
1023  }
1024  }
1025  int row = (*rit)->getPetscGlobalDofIdx();
1026 
1027  auto col_dofs = getColDofsPtr();
1028  for (auto cit = col_dofs->begin(); cit != col_dofs->end(); cit++) {
1029 
1030  int col = (*cit)->getPetscGlobalDofIdx();
1031 
1032  if (row == rowPrint && col == colPrint) {
1033  MOFEM_LOG("SELF", Sev::noisy) << "fe:\n"
1034  << *numeredEntFiniteElementPtr;
1035  MOFEM_LOG("SELF", Sev::noisy) << "row:\n" << *(*rit);
1036  MOFEM_LOG("SELF", Sev::noisy) << "col:\n" << *(*cit);
1037  MOFEM_LOG("SELF", Sev::noisy)
1038  << "fe:\n"
1039  << numeredEntFiniteElementPtr->getBitRefLevel();
1040  MOFEM_LOG("SELF", Sev::noisy) << "row:\n"
1041  << (*rit)->getBitRefLevel();
1042  MOFEM_LOG("SELF", Sev::noisy) << "col:\n"
1043  << (*cit)->getBitRefLevel();
1044  }
1045 
1046  CHKERR MatSetValue(A, row, col, 1, INSERT_VALUES);
1047  }
1048 
1049  if ((*rit)->getEntType() != MBVERTEX) {
1050 
1051  auto range =
1052  row_dofs->get<Ent_mi_tag>().equal_range((*rit)->getEnt());
1053  int nb_dofs_on_ent = std::distance(range.first, range.second);
1054 
1055  int max_order = (*rit)->getMaxOrder();
1056  if ((*rit)->getNbOfCoeffs() * (*rit)->getOrderNbDofs(max_order) !=
1057  nb_dofs_on_ent) {
1058 
1059  /* It could be that you have removed DOFs from problem, and for
1060  * example if this was vector filed with components {Ux,Uy,Uz}, you
1061  * removed on Uz element. */
1062 
1063  MOFEM_LOG("SELF", Sev::warning)
1064  << "Warning: Number of Dofs in Row different than number "
1065  "of dofs for given entity order "
1066  << (*rit)->getNbOfCoeffs() * (*rit)->getOrderNbDofs(max_order)
1067  << " " << nb_dofs_on_ent;
1068  }
1069  }
1070  }
1071 
1073  }
1074 
1075  MoFEMErrorCode postProcess() {
1077 
1078  // cerr << mFieldPtr->get_comm_rank() << endl;
1079  CHKERR MatAssemblyBegin(A, MAT_FLUSH_ASSEMBLY);
1080  CHKERR MatAssemblyEnd(A, MAT_FLUSH_ASSEMBLY);
1081 
1083  }
1084  };
1085 
1086  // create matrix
1087  CHKERR MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
1088 
1089  if (verb >= VERY_VERBOSE) {
1090  MatView(A, PETSC_VIEWER_STDOUT_WORLD);
1091  }
1092 
1093  if (verb >= NOISY) {
1094  MatView(A, PETSC_VIEWER_DRAW_WORLD);
1095  std::string wait;
1096  std::cin >> wait;
1097  }
1098 
1099  TestMatrixFillIn method(&m_field, A, row_print, col_print);
1100 
1101  // get problem
1102  auto problems_ptr = m_field.get_problems();
1103  auto &prb_set = problems_ptr->get<Problem_mi_tag>();
1104  auto p_miit = prb_set.find(problem_name);
1105  if (p_miit == prb_set.end())
1106  SETERRQ1(m_field.get_comm(), MOFEM_DATA_INCONSISTENCY,
1107  "problem < %s > not found (top tip: check spelling)",
1108  problem_name.c_str());
1109  MOFEM_LOG_C("WORLD", Sev::inform, "check problem < %s >",
1110  problem_name.c_str());
1111 
1112  // loop all elements in problem and check if assemble is without error
1113  auto fe_ptr = m_field.get_finite_elements();
1114  for (auto &fe : *fe_ptr) {
1115  MOFEM_LOG_C("WORLD", Sev::verbose, "\tcheck element %s",
1116  fe->getName().c_str());
1117  CHKERR m_field.loop_finite_elements(problem_name, fe->getName(), method,
1118  nullptr, MF_EXIST,
1119  CacheTupleSharedPtr(), verb);
1120  }
1121 
1122  CHKERR MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
1123  CHKERR MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
1124 
1125  PetscLogEventEnd(MOFEM_EVENT_checkMatrixFillIn, 0, 0, 0, 0);
1126 
1128 }

◆ checkMPIAIJMatrixFillIn() [1/3]

template<>
MoFEMErrorCode MoFEM::MatrixManager::checkMPIAIJMatrixFillIn ( const std::string  problem_name,
int  row_print,
int  col_print,
int  verb 
)

◆ checkMPIAIJMatrixFillIn() [2/3]

template<>
MoFEMErrorCode MoFEM::MatrixManager::checkMPIAIJMatrixFillIn ( const std::string  problem_name,
int  row_print,
int  col_print,
int  verb 
)

Definition at line 1144 of file MatrixManager.cpp.

1145  {
1147  // create matrix
1148  SmartPetscObj<Mat> A;
1149  CHKERR createMPIAIJ<PetscGlobalIdx_mi_tag>(problem_name, A, verb);
1150  CHKERR MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
1151  CHKERR checkMatrixFillIn(problem_name, row_print, col_print, A, verb);
1153 }

◆ checkMPIAIJMatrixFillIn() [3/3]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::checkMPIAIJMatrixFillIn ( const std::string  problem_name,
int  row_print,
int  col_print,
int  verb = QUIET 
)
inline

check if matrix fill in correspond to finite element indices

This is used to check consistency of code. If problem is notices with additional non-zero elements in matrix, this function can help detect problem. Should be used as a part of atom tests

Template Parameters
Tag
Parameters
problem_name
rowprint info at particular row
colprint info at particular col
Returns
MoFEMErrorCode

Definition at line 286 of file MatrixManager.hpp.

288  {
289  static_assert(!std::is_same<Tag, Tag>::value, "not implemented");
290  return 0;
291  }

◆ checkMPIAIJWithArraysMatrixFillIn() [1/3]

template<>
MoFEMErrorCode MoFEM::MatrixManager::checkMPIAIJWithArraysMatrixFillIn ( const std::string  problem_name,
int  row_print,
int  col_print,
int  verb 
)

◆ checkMPIAIJWithArraysMatrixFillIn() [2/3]

template<>
MoFEMErrorCode MoFEM::MatrixManager::checkMPIAIJWithArraysMatrixFillIn ( const std::string  problem_name,
int  row_print,
int  col_print,
int  verb 
)

Definition at line 1132 of file MatrixManager.cpp.

1133  {
1135  // create matrix
1136  SmartPetscObj<Mat> A;
1137  CHKERR createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(problem_name, A, verb);
1138  CHKERR MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);
1139  CHKERR checkMatrixFillIn(problem_name, row_print, col_print, A, verb);
1141 }

◆ checkMPIAIJWithArraysMatrixFillIn() [3/3]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::checkMPIAIJWithArraysMatrixFillIn ( const std::string  problem_name,
int  row_print,
int  col_print,
int  verb = QUIET 
)
inline

check if matrix fill in correspond to finite element indices

This is used to check consistency of code. If problem is notices with additional non-zero elements in matrix, this function can help detect problem. Should be used as a part of atom tests

Template Parameters
Tag
Parameters
problem_name
rowprint info at particular row
colprint info at particular col
Returns
MoFEMErrorCode

Definition at line 265 of file MatrixManager.hpp.

267  {
268  static_assert(!std::is_same<Tag, Tag>::value, "not implemented");
269  return 0;
270  }

◆ createMPIAdjWithArrays() [1/3]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createMPIAdjWithArrays ( const std::string  name,
Mat *  Adj,
int  verb 
)

◆ createMPIAdjWithArrays() [2/3]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createMPIAdjWithArrays ( const std::string  name,
Mat *  Adj,
int  verb 
)

Definition at line 812 of file MatrixManager.cpp.

813  {
814  MoFEM::CoreInterface &m_field = cOre;
815  CreateRowComressedADJMatrix *core_ptr =
816  static_cast<CreateRowComressedADJMatrix *>(&cOre);
818  PetscLogEventBegin(MOFEM_EVENT_createMPIAdjWithArrays, 0, 0, 0, 0);
819 
820  auto problems_ptr = m_field.get_problems();
821  auto &prb = problems_ptr->get<Problem_mi_tag>();
822  auto p_miit = prb.find(name);
823  if (p_miit == prb.end()) {
824  SETERRQ1(m_field.get_comm(), MOFEM_NOT_FOUND,
825  "problem < %s > is not found (top tip: check spelling)",
826  name.c_str());
827  }
828 
829  std::vector<int> i_vec, j_vec;
830  j_vec.reserve(10000);
831  CHKERR core_ptr->createMatArrays<Idx_mi_tag>(p_miit, MATMPIADJ, i_vec, j_vec,
832  true, verb);
833  int *_i, *_j;
834  CHKERR PetscMalloc(i_vec.size() * sizeof(int), &_i);
835  CHKERR PetscMalloc(j_vec.size() * sizeof(int), &_j);
836  copy(i_vec.begin(), i_vec.end(), _i);
837  copy(j_vec.begin(), j_vec.end(), _j);
838 
839  int nb_col_dofs = p_miit->getNbDofsCol();
840  CHKERR MatCreateMPIAdj(m_field.get_comm(), i_vec.size() - 1, nb_col_dofs, _i,
841  _j, PETSC_NULL, Adj);
842  CHKERR MatSetOption(*Adj, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE);
843 
844  PetscLogEventEnd(MOFEM_EVENT_createMPIAdjWithArrays, 0, 0, 0, 0);
846 }

◆ createMPIAdjWithArrays() [3/3]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createMPIAdjWithArrays ( const std::string  name,
Mat *  Adj,
int  verb = QUIET 
)
inline

Creates a sparse matrix representing an adjacency list.

The matrix does not have numerical values associated with it, but is intended for ordering (to reduce bandwidth etc) and partitioning.

See PETSc for details

Note
This matrix object does not support most matrix operations, include MatSetValues(). You must NOT free the ii, values and jj arrays yourself. PETSc will free them when the matrix is destroyed; you must allocate them with PetscMalloc(). If you call from Fortran you need not create the arrays with PetscMalloc(). Should not include the matrix diagonals.
Template Parameters
Tag
Parameters
name
Adj
verb
Returns
MoFEMErrorCode

Definition at line 166 of file MatrixManager.hpp.

167  {
168  static_assert(!std::is_same<Tag, Tag>::value, "not implemented");
169  return MOFEM_NOT_IMPLEMENTED;
170  }

◆ createMPIAIJ() [1/4]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJ ( const std::string  name,
Mat *  Aij,
int  verb 
)

◆ createMPIAIJ() [2/4]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJ ( const std::string  name,
Mat *  Aij,
int  verb 
)

Definition at line 742 of file MatrixManager.cpp.

743  {
744  MoFEM::CoreInterface &m_field = cOre;
745  CreateRowComressedADJMatrix *core_ptr =
746  static_cast<CreateRowComressedADJMatrix *>(&cOre);
748  PetscLogEventBegin(MOFEM_EVENT_createMPIAIJ, 0, 0, 0, 0);
749 
750  auto problems_ptr = m_field.get_problems();
751  auto &prb = problems_ptr->get<Problem_mi_tag>();
752  auto p_miit = prb.find(name);
753  if (p_miit == prb.end()) {
754  SETERRQ1(m_field.get_comm(), MOFEM_NOT_FOUND,
755  "problem < %s > is not found (top tip: check spelling)",
756  name.c_str());
757  }
758 
759  std::vector<int> i_vec, j_vec;
760  j_vec.reserve(10000);
761  CHKERR core_ptr->createMatArrays<PetscGlobalIdx_mi_tag>(
762  p_miit, MATMPIAIJ, i_vec, j_vec, false, verb);
763 
764  int nb_row_dofs = p_miit->getNbDofsRow();
765  int nb_col_dofs = p_miit->getNbDofsCol();
766  int nb_local_dofs_row = p_miit->getNbLocalDofsRow();
767  int nb_local_dofs_col = p_miit->getNbLocalDofsCol();
768 
769  auto get_layout = [&]() {
770  int start_ranges, end_ranges;
771  PetscLayout layout;
772  CHKERR PetscLayoutCreate(m_field.get_comm(), &layout);
773  CHKERR PetscLayoutSetBlockSize(layout, 1);
774  CHKERR PetscLayoutSetLocalSize(layout, nb_local_dofs_col);
775  CHKERR PetscLayoutSetUp(layout);
776  CHKERR PetscLayoutGetRange(layout, &start_ranges, &end_ranges);
777  CHKERR PetscLayoutDestroy(&layout);
778  return std::make_pair(start_ranges, end_ranges);
779  };
780 
781  auto get_nnz = [&](auto &d_nnz, auto &o_nnz) {
783  auto layout = get_layout();
784  int j = 0;
785  for (int i = 0; i != nb_local_dofs_row; ++i) {
786  for (; j != i_vec[i + 1]; ++j) {
787  if (j_vec[j] < layout.second && j_vec[j] >= layout.first)
788  ++(d_nnz[i]);
789  else
790  ++(o_nnz[i]);
791  }
792  }
794  };
795 
796  std::vector<int> d_nnz(nb_local_dofs_row, 0), o_nnz(nb_local_dofs_row, 0);
797  CHKERR get_nnz(d_nnz, o_nnz);
798 
799  CHKERR MatCreate(m_field.get_comm(), Aij);
800  CHKERR MatSetSizes(*Aij, nb_local_dofs_row, nb_local_dofs_col, nb_row_dofs,
801  nb_col_dofs);
802  CHKERR MatSetType(*Aij, MATMPIAIJ);
803  CHKERR MatMPIAIJSetPreallocation(*Aij, 0, &*d_nnz.begin(), 0,
804  &*o_nnz.begin());
805 
806  PetscLogEventEnd(MOFEM_EVENT_createMPIAIJ, 0, 0, 0, 0);
808 }

◆ createMPIAIJ() [3/4]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJ ( const std::string  name,
Mat *  Aij,
int  verb = QUIET 
)
inline

Creates a MPI AIJ matrix using arrays that contain in standard CSR format the local rows.

See PETSc for details

Template Parameters
Tag
Parameters
name
Aij
verb
Returns
MoFEMErrorCode

Definition at line 123 of file MatrixManager.hpp.

124  {
125  static_assert(!std::is_same<Tag, Tag>::value, "not implemented");
126  return 0;
127  }

◆ createMPIAIJ() [4/4]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJ ( const std::string  name,
SmartPetscObj< Mat > &  aij_ptr,
int  verb = QUIET 
)
inline

Creates a MPI AIJ matrix using arrays that contain in standard CSR format the local rows.

See PETSc for details

Template Parameters
Tag
Parameters
name
Aij
verb
Returns
MoFEMErrorCode

Definition at line 132 of file MatrixManager.hpp.

134  {
136  Mat aij;
137  CHKERR createMPIAIJ<Tag>(name, &aij, verb);
138  aij_ptr.reset(aij, false);
140  }

◆ createMPIAIJCUSPARSEWithArrays() [1/4]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJCUSPARSEWithArrays ( const std::string  name,
Mat *  Aij,
int  verb 
)

◆ createMPIAIJCUSPARSEWithArrays() [2/4]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJCUSPARSEWithArrays ( const std::string  name,
Mat *  Aij,
int  verb 
)

Definition at line 654 of file MatrixManager.cpp.

655  {
657 
658  MoFEM::CoreInterface &m_field = cOre;
659  CreateRowComressedADJMatrix *core_ptr =
660  static_cast<CreateRowComressedADJMatrix *>(&cOre);
661  PetscLogEventBegin(MOFEM_EVENT_createMPIAIJCUSPARSEWithArrays, 0, 0, 0, 0);
662 
663  auto problems_ptr = m_field.get_problems();
664  auto &prb = problems_ptr->get<Problem_mi_tag>();
665  auto p_miit = prb.find(name);
666  if (p_miit == prb.end()) {
667  SETERRQ1(m_field.get_comm(), MOFEM_NOT_FOUND,
668  "problem < %s > is not found (top tip: check spelling)",
669  name.c_str());
670  }
671 
672  std::vector<int> i_vec, j_vec;
673  j_vec.reserve(10000);
674  CHKERR core_ptr->createMatArrays<PetscGlobalIdx_mi_tag>(
675  p_miit, MATAIJCUSPARSE, i_vec, j_vec, false, verb);
676 
677  int nb_local_dofs_row = p_miit->getNbLocalDofsRow();
678  int nb_local_dofs_col = p_miit->getNbLocalDofsCol();
679 
680  auto get_layout = [&]() {
681  int start_ranges, end_ranges;
682  PetscLayout layout;
683  CHKERR PetscLayoutCreate(m_field.get_comm(), &layout);
684  CHKERR PetscLayoutSetBlockSize(layout, 1);
685  CHKERR PetscLayoutSetLocalSize(layout, nb_local_dofs_col);
686  CHKERR PetscLayoutSetUp(layout);
687  CHKERR PetscLayoutGetRange(layout, &start_ranges, &end_ranges);
688  CHKERR PetscLayoutDestroy(&layout);
689  return std::make_pair(start_ranges, end_ranges);
690  };
691 
692  auto get_nnz = [&](auto &d_nnz, auto &o_nnz) {
694  auto layout = get_layout();
695  int j = 0;
696  for (int i = 0; i != nb_local_dofs_row; ++i) {
697  for (; j != i_vec[i + 1]; ++j) {
698  if (j_vec[j] < layout.second && j_vec[j] >= layout.first)
699  ++(d_nnz[i]);
700  else
701  ++(o_nnz[i]);
702  }
703  }
705  };
706 
707  std::vector<int> d_nnz(nb_local_dofs_row, 0), o_nnz(nb_local_dofs_row, 0);
708  CHKERR get_nnz(d_nnz, o_nnz);
709 
710 #ifdef PETSC_HAVE_CUDA
711  CHKERR ::MatCreateAIJCUSPARSE(m_field.get_comm(), nb_local_dofs_row,
712  nb_local_dofs_col, nb_row_dofs, nb_col_dofs, 0,
713  &*d_nnz.begin(), 0, &*o_nnz.begin(), Aij);
714 #else
715  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
716  "Error: To use this matrix type compile PETSc with CUDA.");
717 #endif
718 
719  PetscLogEventEnd(MOFEM_EVENT_createMPIAIJCUSPARSEWithArrays, 0, 0, 0, 0);
720 
722 }

◆ createMPIAIJCUSPARSEWithArrays() [3/4]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJCUSPARSEWithArrays ( const std::string  name,
Mat *  Aij,
int  verb = QUIET 
)
inline

Definition at line 69 of file MatrixManager.hpp.

70  {
71  static_assert(!std::is_same<Tag, Tag>::value, "not implemented");
72  return 0;
73  }

◆ createMPIAIJCUSPARSEWithArrays() [4/4]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJCUSPARSEWithArrays ( const std::string  name,
SmartPetscObj< Mat > &  aij_ptr,
int  verb = QUIET 
)
inline

Definition at line 78 of file MatrixManager.hpp.

80  {
82  Mat aij;
83  CHKERR createMPIAIJCUSPARSEWithArrays<Tag>(name, &aij, verb);
84  aij_ptr.reset(aij, false);
86  }

◆ createMPIAIJWithArrays() [1/4]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJWithArrays ( const std::string  name,
Mat *  Aij,
int  verb 
)

◆ createMPIAIJWithArrays() [2/4]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJWithArrays ( const std::string  name,
Mat *  Aij,
int  verb 
)

Definition at line 616 of file MatrixManager.cpp.

617  {
619 
620  MoFEM::CoreInterface &m_field = cOre;
621  CreateRowComressedADJMatrix *core_ptr =
622  static_cast<CreateRowComressedADJMatrix *>(&cOre);
623  PetscLogEventBegin(MOFEM_EVENT_createMPIAIJWithArrays, 0, 0, 0, 0);
624 
625  auto problems_ptr = m_field.get_problems();
626  auto &prb = problems_ptr->get<Problem_mi_tag>();
627  auto p_miit = prb.find(name);
628  if (p_miit == prb.end()) {
629  SETERRQ1(m_field.get_comm(), MOFEM_NOT_FOUND,
630  "problem < %s > is not found (top tip: check spelling)",
631  name.c_str());
632  }
633 
634  std::vector<int> i_vec, j_vec;
635  j_vec.reserve(10000);
636  CHKERR core_ptr->createMatArrays<PetscGlobalIdx_mi_tag>(
637  p_miit, MATMPIAIJ, i_vec, j_vec, false, verb);
638 
639  int nb_row_dofs = p_miit->getNbDofsRow();
640  int nb_col_dofs = p_miit->getNbDofsCol();
641  int nb_local_dofs_row = p_miit->getNbLocalDofsRow();
642  int nb_local_dofs_col = p_miit->getNbLocalDofsCol();
643 
644  CHKERR ::MatCreateMPIAIJWithArrays(
645  m_field.get_comm(), nb_local_dofs_row, nb_local_dofs_col, nb_row_dofs,
646  nb_col_dofs, &*i_vec.begin(), &*j_vec.begin(), PETSC_NULL, Aij);
647 
648  PetscLogEventEnd(MOFEM_EVENT_createMPIAIJWithArrays, 0, 0, 0, 0);
650 }

◆ createMPIAIJWithArrays() [3/4]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJWithArrays ( const std::string  name,
Mat *  Aij,
int  verb = QUIET 
)
inline

Creates a MPI AIJ matrix using arrays that contain in standard CSR format the local rows.

See PETSc for details

Template Parameters
Tag
Parameters
name
Aij
verb
Returns
MoFEMErrorCode

Definition at line 49 of file MatrixManager.hpp.

50  {
51  static_assert(!std::is_same<Tag, Tag>::value, "not implemented");
52  return 0;
53  }

◆ createMPIAIJWithArrays() [4/4]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createMPIAIJWithArrays ( const std::string  name,
SmartPetscObj< Mat > &  aij_ptr,
int  verb = QUIET 
)
inline

Creates a MPI AIJ matrix using arrays that contain in standard CSR format the local rows.

See PETSc for details

Template Parameters
Tag
Parameters
name
Aij
verb
Returns
MoFEMErrorCode

Definition at line 58 of file MatrixManager.hpp.

60  {
62  Mat aij;
63  CHKERR createMPIAIJWithArrays<Tag>(name, &aij, verb);
64  aij_ptr.reset(aij, false);
66  }

◆ createSeqAIJCUSPARSEWithArrays() [1/4]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createSeqAIJCUSPARSEWithArrays ( const std::string  name,
Mat *  Aij,
int  verb 
)

◆ createSeqAIJCUSPARSEWithArrays() [2/4]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createSeqAIJCUSPARSEWithArrays ( const std::string  name,
Mat *  Aij,
int  verb 
)

Definition at line 726 of file MatrixManager.cpp.

727  {
729 
730 #ifdef PETSC_HAVE_CUDA
731  // CHKERR ::MatCreateSeqAIJCUSPARSE(MPI_Comm comm, PetscInt m, PetscInt n,
732  // PetscInt nz, const PetscInt nnz[], Mat *A);
733  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
734  "Not implemented type of matrix yet, try MPI version (aijcusparse)");
735 #endif
736 
738 }

◆ createSeqAIJCUSPARSEWithArrays() [3/4]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createSeqAIJCUSPARSEWithArrays ( const std::string  name,
Mat *  Aij,
int  verb = QUIET 
)
inline

Definition at line 89 of file MatrixManager.hpp.

90  {
91  static_assert(!std::is_same<Tag, Tag>::value, "not implemented");
92  return 0;
93  }

◆ createSeqAIJCUSPARSEWithArrays() [4/4]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createSeqAIJCUSPARSEWithArrays ( const std::string  name,
SmartPetscObj< Mat > &  aij_ptr,
int  verb = QUIET 
)
inline

Definition at line 98 of file MatrixManager.hpp.

100  {
102  Mat aij;
103  CHKERR createSeqAIJCUSPARSEWithArrays<Tag>(name, &aij, verb);
104  aij_ptr.reset(aij, false);
106  }

◆ createSeqAIJWithArrays() [1/4]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createSeqAIJWithArrays ( const std::string  name,
Mat *  Aij,
int  verb 
)

◆ createSeqAIJWithArrays() [2/4]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createSeqAIJWithArrays ( const std::string  name,
Mat *  Aij,
int  verb 
)

Definition at line 849 of file MatrixManager.cpp.

850  {
851  MoFEM::CoreInterface &m_field = cOre;
852  CreateRowComressedADJMatrix *core_ptr =
853  static_cast<CreateRowComressedADJMatrix *>(&cOre);
855  PetscLogEventBegin(MOFEM_EVENT_createMPIAIJWithArrays, 0, 0, 0, 0);
856 
857  auto problems_ptr = m_field.get_problems();
858  auto &prb = problems_ptr->get<Problem_mi_tag>();
859  auto p_miit = prb.find(name);
860  if (p_miit == prb.end()) {
861  SETERRQ1(m_field.get_comm(), MOFEM_NOT_FOUND,
862  "problem < %s > is not found (top tip: check spelling)",
863  name.c_str());
864  }
865 
866  std::vector<int> i_vec, j_vec;
867  j_vec.reserve(10000);
868  CHKERR core_ptr->createMatArrays<PetscGlobalIdx_mi_tag>(p_miit, MATAIJ, i_vec,
869  j_vec, false, verb);
870 
871  int nb_local_dofs_row = p_miit->getNbLocalDofsRow();
872  int nb_local_dofs_col = p_miit->getNbLocalDofsCol();
873 
874  double *_a;
875  CHKERR PetscMalloc(j_vec.size() * sizeof(double), &_a);
876 
877  Mat tmpMat;
878  CHKERR ::MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, nb_local_dofs_row,
879  nb_local_dofs_col, &*i_vec.begin(),
880  &*j_vec.begin(), _a, &tmpMat);
881  CHKERR MatDuplicate(tmpMat, MAT_SHARE_NONZERO_PATTERN, Aij);
882  CHKERR MatDestroy(&tmpMat);
883 
884  CHKERR PetscFree(_a);
885 
886  PetscLogEventEnd(MOFEM_EVENT_createMPIAIJWithArrays, 0, 0, 0, 0);
888 }

◆ createSeqAIJWithArrays() [3/4]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createSeqAIJWithArrays ( const std::string  name,
Mat *  Aij,
int  verb = QUIET 
)
inline

Create sequencial matrix.

Creates a sparse matrix in AIJ (compressed row) format (the default parallel PETSc format). For good matrix assembly performance the user should preallocate the matrix storage by setting the parameter nz (or the array nnz). By setting these parameters accurately, performance during matrix assembly can be increased by more than a factor of 50.

See PETSc for details

Template Parameters
Tag
Parameters
name
Aij
i
j
v
verb
Returns
MoFEMErrorCode

Definition at line 195 of file MatrixManager.hpp.

196  {
197  static_assert(!std::is_same<Tag, Tag>::value, "not implemented");
198  return 0;
199  }

◆ createSeqAIJWithArrays() [4/4]

template<class Tag >
MoFEMErrorCode MoFEM::MatrixManager::createSeqAIJWithArrays ( const std::string  name,
SmartPetscObj< Mat > &  aij_ptr,
int  verb = QUIET 
)
inline

Create sequencial matrix.

Creates a sparse matrix in AIJ (compressed row) format (the default parallel PETSc format). For good matrix assembly performance the user should preallocate the matrix storage by setting the parameter nz (or the array nnz). By setting these parameters accurately, performance during matrix assembly can be increased by more than a factor of 50.

See PETSc for details

Template Parameters
Tag
Parameters
name
Aij(SmartPetscObj)
i
j
v
verb
Returns
MoFEMErrorCode

Definition at line 224 of file MatrixManager.hpp.

226  {
228  Mat aij;
229  CHKERR createSeqAIJWithArrays<Tag>(name, &aij, verb);
230  aij_ptr.reset(aij, false);
232  }

◆ query_interface()

MoFEMErrorCode MoFEM::MatrixManager::query_interface ( boost::typeindex::type_index  type_index,
UnknownInterface **  iface 
) const
virtual

Implements MoFEM::UnknownInterface.

Definition at line 591 of file MatrixManager.cpp.

592  {
593  *iface = const_cast<MatrixManager *>(this);
594  return 0;
595 }

Member Data Documentation

◆ cOre

MoFEM::Core& MoFEM::MatrixManager::cOre

Definition at line 26 of file MatrixManager.hpp.

◆ MOFEM_EVENT_checkMatrixFillIn

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_checkMatrixFillIn
private

Definition at line 300 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createMPIAdjWithArrays

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createMPIAdjWithArrays
private

Definition at line 296 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createMPIAIJ

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createMPIAIJ
private

Definition at line 294 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createMPIAIJCUSPARSEWithArrays

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createMPIAIJCUSPARSEWithArrays
private

Definition at line 298 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createMPIAIJWithArrays

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createMPIAIJWithArrays
private

Definition at line 295 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createSeqAIJCUSPARSEWithArrays

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createSeqAIJCUSPARSEWithArrays
private

Definition at line 299 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createSeqAIJWithArrays

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createSeqAIJWithArrays
private

Definition at line 297 of file MatrixManager.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::CoreInterface::loop_finite_elements
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEM::MatrixManager::MOFEM_EVENT_createMPIAIJCUSPARSEWithArrays
PetscLogEvent MOFEM_EVENT_createMPIAIJCUSPARSEWithArrays
Definition: MatrixManager.hpp:298
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
NOISY
@ NOISY
Definition: definitions.h:211
MatrixManagerFunctionBegin
#define MatrixManagerFunctionBegin
Create adjacent matrices using different indices.
Definition: MatrixManager.cpp:6
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MOFEM_LOG_ATTRIBUTES
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
Definition: LogManager.hpp:296
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
MoFEM::MatrixManager::MOFEM_EVENT_createSeqAIJCUSPARSEWithArrays
PetscLogEvent MOFEM_EVENT_createSeqAIJCUSPARSEWithArrays
Definition: MatrixManager.hpp:299
MoFEM::CoreInterface::get_finite_elements
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
MoFEM::MatrixManager::MOFEM_EVENT_checkMatrixFillIn
PetscLogEvent MOFEM_EVENT_checkMatrixFillIn
Definition: MatrixManager.hpp:300
MoFEM::LogManager::BitScope
@ BitScope
Definition: LogManager.hpp:49
MoFEM::MatrixManager::MOFEM_EVENT_createMPIAIJWithArrays
PetscLogEvent MOFEM_EVENT_createMPIAIJWithArrays
Definition: MatrixManager.hpp:295
FEMethod
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
a
constexpr double a
Definition: approx_sphere.cpp:30
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
MoFEM::Types::UId
uint128_t UId
Unique Id.
Definition: Types.hpp:31
MoFEM::MatrixManager::MOFEM_EVENT_createMPIAdjWithArrays
PetscLogEvent MOFEM_EVENT_createMPIAdjWithArrays
Definition: MatrixManager.hpp:296
MoFEM::MatrixManager::MatrixManager
MatrixManager(const MoFEM::Core &core)
Definition: MatrixManager.cpp:597
MoFEM::CoreInterface::get_problems
virtual const Problem_multiIndex * get_problems() const =0
Get the problems object.
MoFEM::CacheTupleSharedPtr
boost::shared_ptr< CacheTuple > CacheTupleSharedPtr
Definition: FEMultiIndices.hpp:495
MoFEM::MatrixManager::MOFEM_EVENT_createMPIAIJ
PetscLogEvent MOFEM_EVENT_createMPIAIJ
Definition: MatrixManager.hpp:294
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MOFEM_NOT_FOUND
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::CoreInterface
Interface.
Definition: Interface.hpp:27
MF_EXIST
@ MF_EXIST
Definition: definitions.h:100
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEM::MatrixManager::cOre
MoFEM::Core & cOre
Definition: MatrixManager.hpp:26
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
VERY_VERBOSE
@ VERY_VERBOSE
Definition: definitions.h:210
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::MatrixManager::checkMatrixFillIn
MoFEMErrorCode checkMatrixFillIn(const std::string problem_name, int row_print, int col_print, Mat A, int verb=QUIET)
check if matrix fill in correspond to finite element indices
Definition: MatrixManager.cpp:890
MoFEM::MatrixManager::MOFEM_EVENT_createSeqAIJWithArrays
PetscLogEvent MOFEM_EVENT_createSeqAIJWithArrays
Definition: MatrixManager.hpp:297