v0.15.0
Loading...
Searching...
No Matches
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.
 
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.
 
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.
 
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.
 
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.
 
template<class Tag >
MoFEMErrorCode createMPIAdjWithArrays (const std::string name, Mat *Adj, int verb=QUIET)
 Creates a sparse matrix representing an adjacency list.
 
template<class Tag >
MoFEMErrorCode createSeqAIJWithArrays (const std::string name, Mat *Aij, int verb=QUIET)
 Create sequencial matrix.
 
template<class Tag >
MoFEMErrorCode createSeqAIJWithArrays (const std::string name, SmartPetscObj< Mat > &aij_ptr, int verb=QUIET)
 Create sequencial matrix.
 
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
 
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
 
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
 
template<class Tag >
MoFEMErrorCode createHybridL2MPIAIJ (const std::string problem_name, SmartPetscObj< Mat > &aij_ptr, int verb=QUIET)
 Create a Hybrid MPIAij object.
 
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 createHybridL2MPIAIJ (const std::string problem_name, SmartPetscObj< Mat > &aij_ptr, 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)
 
template<>
MoFEMErrorCode createHybridL2MPIAIJ (const std::string problem_name, SmartPetscObj< Mat > &aij_ptr, int verb)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
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.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 

Detailed Description

Matrix manager is used to build and partition problems.

Examples
mofem/atom_tests/build_composite_problem.cpp, mofem/atom_tests/build_large_problem.cpp, mofem/atom_tests/build_problems.cpp, mofem/atom_tests/dm_build_partitioned_mesh.cpp, mofem/atom_tests/dm_create_subdm.cpp, mofem/atom_tests/dm_partitioned_no_field.cpp, mofem/atom_tests/prism_polynomial_approximation.cpp, mofem/atom_tests/quad_polynomial_approximation.cpp, mofem/atom_tests/remove_entities_from_problem.cpp, mofem/atom_tests/remove_entities_from_problem_not_partitioned.cpp, and mofem/users_modules/basic_finite_elements/nonlinear_elasticity/nonlinear_dynamics.cpp.

Definition at line 21 of file MatrixManager.hpp.

Constructor & Destructor Documentation

◆ MatrixManager()

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

Definition at line 606 of file MatrixManager.cpp.

607 : cOre(const_cast<MoFEM::Core &>(core)) {
608 PetscLogEventRegister("MatMngCrtMPIAIJ", 0, &MOFEM_EVENT_createMPIAIJ);
609 PetscLogEventRegister("MatMngCrtMPIAIJWthArr", 0,
611 PetscLogEventRegister("MatMngCrtMPIAdjWithArr", 0,
613 PetscLogEventRegister("MatMngCrtMPIAIJCUSPARSEWthArr", 0,
615 PetscLogEventRegister("MatMngCrtSeqAIJCUSPARSEWthArrs", 0,
617 PetscLogEventRegister("MatMngCrtSeqAIJWthArrs", 0,
619 PetscLogEventRegister("MatMngCrtCheckMPIAIJWthArrFillIn", 0,
621}
Core (interface) class.
Definition Core.hpp:82
PetscLogEvent MOFEM_EVENT_createSeqAIJCUSPARSEWithArrays
PetscLogEvent MOFEM_EVENT_createMPIAdjWithArrays
PetscLogEvent MOFEM_EVENT_checkMatrixFillIn
PetscLogEvent MOFEM_EVENT_createMPIAIJ
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.
PetscLogEvent MOFEM_EVENT_createMPIAIJWithArrays
PetscLogEvent MOFEM_EVENT_createSeqAIJWithArrays
PetscLogEvent MOFEM_EVENT_createMPIAIJCUSPARSEWithArrays

◆ ~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 899 of file MatrixManager.cpp.

901 {
902 MoFEM::CoreInterface &m_field = cOre;
904
906
907 struct TestMatrixFillIn : public FEMethod {
908 CoreInterface *mFieldPtr;
909
910 Mat A;
911
912 int rowPrint, colPrint;
913
914 TestMatrixFillIn(CoreInterface *m_field_ptr, Mat a, int row_print,
915 int col_print)
918
919 MoFEMErrorCode preProcess() {
922 }
923
924 MoFEMErrorCode operator()() {
926
927 if (refinedFiniteElementsPtr->find(
928 numeredEntFiniteElementPtr->getEnt()) ==
929 refinedFiniteElementsPtr->end()) {
931 "data inconsistency");
932 }
933
934 auto row_dofs = getRowDofsPtr();
935 auto col_dofs = getColDofsPtr();
936
937 for (auto cit = col_dofs->begin(); cit != col_dofs->end(); cit++) {
938
939 if (refinedEntitiesPtr->find((*cit)->getEnt()) ==
940 refinedEntitiesPtr->end()) {
942 "data inconsistency");
943 }
944 if (!(*cit)->getActive()) {
946 "data inconsistency");
947 }
948
949 FieldEntityEntFiniteElementAdjacencyMap_multiIndex::index<
950 Composite_Unique_mi_tag>::type::iterator ait;
951 ait = adjacenciesPtr->get<Composite_Unique_mi_tag>().find(
952 boost::make_tuple((*cit)->getFieldEntityPtr()->getLocalUniqueId(),
953 numeredEntFiniteElementPtr->getLocalUniqueId()));
954 if (ait == adjacenciesPtr->end()) {
956 "adjacencies data inconsistency");
957 } else {
958 UId uid = ait->getEntUniqueId();
959 if (entitiesPtr->find(uid) == entitiesPtr->end()) {
961 "data inconsistency");
962 }
963 if (dofsPtr->find((*cit)->getLocalUniqueId()) == dofsPtr->end()) {
965 "data inconsistency");
966 }
967 }
968
969 if ((*cit)->getEntType() != MBVERTEX) {
970
971 auto range =
972 col_dofs->get<Ent_mi_tag>().equal_range((*cit)->getEnt());
973 int nb_dofs_on_ent = std::distance(range.first, range.second);
974
975 int max_order = (*cit)->getMaxOrder();
976 if ((*cit)->getNbOfCoeffs() * (*cit)->getOrderNbDofs(max_order) !=
978
979 /* It could be that you have
980 removed DOFs from problem, and for example if this was vector field
981 with components {Ux,Uy,Uz}, you removed on Uz element.*/
982
983 MOFEM_LOG("SELF", Sev::warning)
984 << "Warning: Number of Dofs in Col different than number "
985 "of dofs for given entity order: "
986 << "Field name: " << (*cit)->getName() << " "
987 << (*cit)->getNbOfCoeffs() * (*cit)->getOrderNbDofs(max_order)
988 << " != " << nb_dofs_on_ent;
989 }
990 }
991 }
992
993 for (auto rit = row_dofs->begin(); rit != row_dofs->end(); rit++) {
994
995 if (refinedEntitiesPtr->find((*rit)->getEnt()) ==
996 refinedEntitiesPtr->end()) {
998 "data inconsistency");
999 }
1000 if (!(*rit)->getActive()) {
1002 "data inconsistency");
1003 }
1004
1005 FieldEntityEntFiniteElementAdjacencyMap_multiIndex::index<
1006 Composite_Unique_mi_tag>::type::iterator ait;
1007 ait = adjacenciesPtr->get<Composite_Unique_mi_tag>().find(
1008 boost::make_tuple((*rit)->getFieldEntityPtr()->getLocalUniqueId(),
1009 numeredEntFiniteElementPtr->getLocalUniqueId()));
1010 if (ait == adjacenciesPtr->end()) {
1012 MOFEM_LOG("SELF", Sev::error) << *(*rit);
1013 MOFEM_LOG("SELF", Sev::error) << *(*rit);
1014 MOFEM_LOG("SELF", Sev::error) << *numeredEntFiniteElementPtr;
1015 MOFEM_LOG("SELF", Sev::error) << "dof: " << (*rit)->getBitRefLevel();
1016 MOFEM_LOG("SELF", Sev::error)
1017 << "fe: " << numeredEntFiniteElementPtr->getBitRefLevel();
1018 MOFEM_LOG("SELF", Sev::error)
1019 << "problem: " << problemPtr->getBitRefLevel();
1020 MOFEM_LOG("SELF", Sev::error)
1021 << "problem mask: " << problemPtr->getBitRefLevelMask();
1023 "adjacencies data inconsistency");
1024 } else {
1025 UId uid = ait->getEntUniqueId();
1026 if (entitiesPtr->find(uid) == entitiesPtr->end()) {
1028 "data inconsistency");
1029 }
1030 if (dofsPtr->find((*rit)->getLocalUniqueId()) == dofsPtr->end()) {
1032 "data inconsistency");
1033 }
1034 }
1035 int row = (*rit)->getPetscGlobalDofIdx();
1036
1037 auto col_dofs = getColDofsPtr();
1038 for (auto cit = col_dofs->begin(); cit != col_dofs->end(); cit++) {
1039
1040 int col = (*cit)->getPetscGlobalDofIdx();
1041
1042 if (row == rowPrint && col == colPrint) {
1043 MOFEM_LOG("SELF", Sev::noisy) << "fe:\n"
1044 << *numeredEntFiniteElementPtr;
1045 MOFEM_LOG("SELF", Sev::noisy) << "row:\n" << *(*rit);
1046 MOFEM_LOG("SELF", Sev::noisy) << "col:\n" << *(*cit);
1047 MOFEM_LOG("SELF", Sev::noisy)
1048 << "fe:\n"
1049 << numeredEntFiniteElementPtr->getBitRefLevel();
1050 MOFEM_LOG("SELF", Sev::noisy) << "row:\n"
1051 << (*rit)->getBitRefLevel();
1052 MOFEM_LOG("SELF", Sev::noisy) << "col:\n"
1053 << (*cit)->getBitRefLevel();
1054 }
1055
1056 CHKERR MatSetValue(A, row, col, 1, INSERT_VALUES);
1057 }
1058
1059 if ((*rit)->getEntType() != MBVERTEX) {
1060
1061 auto range =
1062 row_dofs->get<Ent_mi_tag>().equal_range((*rit)->getEnt());
1063 int nb_dofs_on_ent = std::distance(range.first, range.second);
1064
1065 int max_order = (*rit)->getMaxOrder();
1066 if ((*rit)->getNbOfCoeffs() * (*rit)->getOrderNbDofs(max_order) !=
1068
1069 /* It could be that you have removed DOFs from problem, and for
1070 * example if this was vector field with components {Ux,Uy,Uz}, you
1071 * removed on Uz element. */
1072
1073 MOFEM_LOG("SELF", Sev::warning)
1074 << "Warning: Number of Dofs in Row different than number "
1075 "of dofs for given entity order "
1076 << "Field name: " << (*rit)->getName() << " "
1077 << (*rit)->getNbOfCoeffs() * (*rit)->getOrderNbDofs(max_order)
1078 << " != " << nb_dofs_on_ent;
1079 }
1080 }
1081 }
1082
1084 }
1085
1086 MoFEMErrorCode postProcess() {
1088
1089 // cerr << mFieldPtr->get_comm_rank() << endl;
1092
1094 }
1095 };
1096
1097 // create matrix
1099
1100 if (verb >= VERY_VERBOSE) {
1102 }
1103
1104 if (verb >= NOISY) {
1106 std::string wait;
1107 std::cin >> wait;
1108 }
1109
1111
1112 // get problem
1113 auto problems_ptr = m_field.get_problems();
1114 auto &prb_set = problems_ptr->get<Problem_mi_tag>();
1115 auto p_miit = prb_set.find(problem_name);
1116 if (p_miit == prb_set.end())
1118 "problem < %s > not found (top tip: check spelling)",
1119 problem_name.c_str());
1120 MOFEM_LOG_C("WORLD", Sev::inform, "check problem < %s >",
1121 problem_name.c_str());
1122
1123 // loop all elements in problem and check if assemble is without error
1124 auto fe_ptr = m_field.get_finite_elements();
1125 for (auto &fe : *fe_ptr) {
1126 MOFEM_LOG_C("WORLD", Sev::verbose, "\tcheck element %s",
1127 fe->getName().c_str());
1128 CHKERR m_field.loop_finite_elements(problem_name, fe->getName(), method,
1129 nullptr, MF_EXIST,
1130 CacheTupleSharedPtr(), verb);
1131 }
1132
1135
1137
1139}
#define MOFEM_LOG_C(channel, severity, format,...)
#define MatrixManagerFunctionBegin
Create adjacent matrices using different indices.
constexpr double a
@ VERY_VERBOSE
@ NOISY
@ MF_EXIST
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
virtual const Problem_multiIndex * get_problems() const =0
Get the problems object.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
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.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
uint128_t UId
Unique Id.
Definition Types.hpp:31
boost::shared_ptr< CacheTuple > CacheTupleSharedPtr
constexpr AssemblyType A
[Define dimension]
virtual MPI_Comm & get_comm() const =0

◆ checkMPIAIJMatrixFillIn() [1/3]

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

Definition at line 1155 of file MatrixManager.cpp.

1156 {
1158 // create matrix
1159 SmartPetscObj<Mat> A;
1164}
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

◆ checkMPIAIJMatrixFillIn() [2/3]

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

◆ 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 
)

Definition at line 1142 of file MatrixManager.cpp.

◆ checkMPIAIJWithArraysMatrixFillIn() [2/3]

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

◆ 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 }

◆ createHybridL2MPIAIJ() [1/3]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createHybridL2MPIAIJ ( const std::string  problem_name,
SmartPetscObj< Mat > &  aij_ptr,
int  verb 
)

Definition at line 1167 of file MatrixManager.cpp.

1168 {
1169 MoFEM::CoreInterface &m_field = cOre;
1171
1172 auto prb_ptr = m_field.get_problems();
1173 auto p_miit = prb_ptr->get<Problem_mi_tag>().find(problem_name);
1174 if (p_miit == prb_ptr->get<Problem_mi_tag>().end()) {
1175 SETERRQ(m_field.get_comm(), MOFEM_NOT_FOUND,
1176 "problem < %s > is not found (top tip: check spelling)",
1177 problem_name.c_str());
1178 }
1179
1180 auto nb_rows = p_miit->getNbDofsRow();
1181 auto nb_cols = p_miit->getNbDofsCol();
1182 auto nb_loc_rows = p_miit->getNbLocalDofsRow();
1183 auto nb_loc_cols = p_miit->getNbLocalDofsCol();
1184
1185 auto row_ptr = p_miit->getNumeredRowDofsPtr();
1186 auto col_ptr = p_miit->getNumeredColDofsPtr();
1187
1189 for (auto &c : *col_ptr) {
1190 fields_ids |= c->getId();
1191 }
1192 auto fields = m_field.get_fields();
1193 std::vector<int> fields_bit_numbers;
1194 for (auto &f : *fields) {
1195 if ((fields_ids & f->getId()).any()) {
1196 fields_bit_numbers.push_back(f->getBitNumber());
1197 }
1198 }
1199
1200 Range adj;
1202 int prev_dim = -1;
1203
1204 auto get_adj = [&](auto ent, auto dim) {
1205 if (prev_ent == ent && prev_dim == dim) {
1206 return adj;
1207 } else {
1208 adj.clear();
1209 Range bridge;
1210 CHKERR m_field.get_moab().get_adjacencies(&ent, 1, dim + 1, false,
1211 bridge);
1212 CHKERR m_field.get_moab().get_adjacencies(bridge, dim, false, adj,
1213 moab::Interface::UNION);
1214 prev_ent = ent;
1215 prev_dim = dim;
1216 return adj;
1217 }
1218 };
1219
1220 int prev_bit_number = -1;
1223 using IT = boost::multi_index::index<NumeredDofEntity_multiIndex,
1224 Unique_mi_tag>::type::iterator;
1225 std::pair<IT, IT> pair_lo_hi;
1226 auto get_col_it = [&](auto bit_number, auto first_ent, auto second_ent) {
1229 return pair_lo_hi;
1230 } else {
1231 auto lo_it = col_ptr->get<Unique_mi_tag>().lower_bound(
1233 auto hi_it = col_ptr->get<Unique_mi_tag>().upper_bound(
1235 pair_lo_hi = std::make_pair(lo_it, hi_it);
1239 return pair_lo_hi;
1240 }
1241 };
1242
1243 auto create_ghost_vec = [&]() {
1244 SmartPetscObj<Vec> v;
1245 CHKERR m_field.getInterface<VecManager>()->vecCreateGhost(problem_name, ROW,
1246 v);
1250 return v;
1251 };
1252
1253 auto v_o_nnz = create_ghost_vec();
1254 auto v_d_nnz = create_ghost_vec();
1255
1256
1257 double *o_nnz_real;
1259 double *d_nnz_real;
1261
1262 for (auto r = row_ptr->begin(); r != row_ptr->end(); ++r) {
1263
1264 auto row_loc_idx = (*r)->getPetscLocalDofIdx();
1265 if (row_loc_idx < 0)
1266 continue;
1267
1268 auto ent = (*r)->getEnt();
1269 auto dim = dimension_from_handle(ent);
1270 auto adj = get_adj(ent, dim);
1271
1272 for (auto p = adj.pair_begin(); p != adj.pair_end(); ++p) {
1273 auto first_ent = p->first;
1274 auto second_ent = p->second;
1275
1276 for (auto bit_number : fields_bit_numbers) {
1277
1279
1280 for (; lo_it != hi_it; ++lo_it) {
1281 auto col_loc_idx = (*lo_it)->getPetscLocalDofIdx();
1282 if (col_loc_idx < 0)
1283 continue;
1284
1285 if (
1286
1287 (*lo_it)->getOwnerProc() == (*r)->getOwnerProc()
1288
1289 ) {
1290 d_nnz_real[row_loc_idx] += 1;
1291 } else {
1292 o_nnz_real[row_loc_idx] += 1;
1293 }
1294 }
1295 }
1296 }
1297 }
1298
1301
1302 auto update_vec = [&](auto v) {
1309 };
1310
1313
1314 int o_nz = 0;
1315 int d_nz = 0;
1316
1319
1320 std::vector<int> d_nnz(nb_loc_rows, 0);
1321 std::vector<int> o_nnz(nb_loc_rows, 0);
1322 for (auto r = row_ptr->begin(); r != row_ptr->end(); ++r) {
1323
1324 auto row_loc_idx = (*r)->getPetscLocalDofIdx();
1325
1326 if (
1327
1329
1330 ) {
1331 auto row_loc_idx = (*r)->getPetscLocalDofIdx();
1336 }
1337
1338 }
1339
1342
1343 if (verb >= QUIET) {
1344 MOFEM_LOG("SYNC", Sev::verbose)
1345 << "Hybrid L2 matrix d_nz: " << d_nz << " o_nz: " << o_nz;
1346 MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::verbose);
1347 }
1348
1349 if (verb >= VERBOSE) {
1350 MOFEM_LOG("SYNC", Sev::noisy) << "Hybrid L2 matrix";
1351 int idx = 0;
1352 for (auto &d : d_nnz) {
1353 MOFEM_LOG("SYNC", Sev::noisy) << idx << ": " << d;
1354 ++idx;
1355 }
1356 MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::noisy);
1357 }
1358
1359 Mat a_raw;
1360 CHKERR MatCreateAIJ(m_field.get_comm(), nb_loc_rows, nb_loc_cols, nb_rows,
1361 nb_cols, 0, &*d_nnz.begin(), 0, &*o_nnz.begin(), &a_raw);
1363 aij_ptr = SmartPetscObj<Mat>(a_raw);
1364
1365 MOFEM_LOG_CHANNEL("WORLD");
1366 MOFEM_LOG_CHANNEL("SYNC");
1367
1369}
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
@ QUIET
@ VERBOSE
@ ROW
@ MOFEM_NOT_FOUND
Definition definitions.h:33
multi_index_container< boost::shared_ptr< NumeredDofEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< NumeredDofEntity::interface_type_DofEntity, UId, &NumeredDofEntity::getLocalUniqueId > >, ordered_non_unique< tag< Part_mi_tag >, member< NumeredDofEntity, unsigned int, &NumeredDofEntity::pArt > >, ordered_non_unique< tag< Idx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::dofIdx > >, ordered_non_unique< tag< PetscGlobalIdx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::petscGloablDofIdx > >, ordered_non_unique< tag< PetscLocalIdx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::petscLocalDofIdx > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< NumeredDofEntity::interface_type_DofEntity, EntityHandle, &NumeredDofEntity::getEnt > > > > NumeredDofEntity_multiIndex
MultiIndex container keeps NumeredDofEntity.
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition Types.hpp:42
auto dimension_from_handle(const EntityHandle h)
get entity dimension form handle
int r
Definition sdf.py:8
virtual moab::Interface & get_moab()=0
static UId getLoFieldEntityUId(const FieldBitNumber bit, const EntityHandle ent)
static UId getHiFieldEntityUId(const FieldBitNumber bit, const EntityHandle ent)
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
constexpr IntegrationType IT

◆ createHybridL2MPIAIJ() [2/3]

template<>
MoFEMErrorCode MoFEM::MatrixManager::createHybridL2MPIAIJ ( const std::string  problem_name,
SmartPetscObj< Mat > &  aij_ptr,
int  verb 
)

◆ createHybridL2MPIAIJ() [3/3]

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

Create a Hybrid MPIAij object.

It assumes that L2 field is on skeleton, and dim+1 is bridge adjacency.

Template Parameters
Tag
Parameters
problem_name
aij_ptr
verb
Returns
MoFEMErrorCode

Definition at line 305 of file MatrixManager.hpp.

307 {
308 static_assert(!std::is_same<Tag, Tag>::value, "not implemented");
309 return 0;
310 }

◆ createMPIAdjWithArrays() [1/3]

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

Definition at line 820 of file MatrixManager.cpp.

822 {
823 MoFEM::CoreInterface &m_field = cOre;
824 CreateRowCompressedADJMatrix *core_ptr =
825 static_cast<CreateRowCompressedADJMatrix *>(&cOre);
828
829 auto problems_ptr = m_field.get_problems();
830 auto &prb = problems_ptr->get<Problem_mi_tag>();
831 auto p_miit = prb.find(name);
832 if (p_miit == prb.end()) {
834 "problem < %s > is not found (top tip: check spelling)",
835 name.c_str());
836 }
837
838 std::vector<int> i_vec, j_vec;
839 j_vec.reserve(10000);
840 CHKERR core_ptr->createMatArrays<Idx_mi_tag>(p_miit, MATMPIADJ, i_vec, j_vec,
841 true, verb);
842 int *_i, *_j;
843 CHKERR PetscMalloc(i_vec.size() * sizeof(int), &_i);
844 CHKERR PetscMalloc(j_vec.size() * sizeof(int), &_j);
845 copy(i_vec.begin(), i_vec.end(), _i);
846 copy(j_vec.begin(), j_vec.end(), _j);
847
848 int nb_col_dofs = p_miit->getNbDofsCol();
849 CHKERR MatCreateMPIAdj(m_field.get_comm(), i_vec.size() - 1, nb_col_dofs, _i,
852
855}

◆ createMPIAdjWithArrays() [2/3]

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

◆ 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");
170 }
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32

◆ createMPIAIJ() [1/4]

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

Definition at line 750 of file MatrixManager.cpp.

752 {
753 MoFEM::CoreInterface &m_field = cOre;
754 CreateRowCompressedADJMatrix *core_ptr =
755 static_cast<CreateRowCompressedADJMatrix *>(&cOre);
758
759 auto problems_ptr = m_field.get_problems();
760 auto &prb = problems_ptr->get<Problem_mi_tag>();
761 auto p_miit = prb.find(name);
762 if (p_miit == prb.end()) {
764 "problem < %s > is not found (top tip: check spelling)",
765 name.c_str());
766 }
767
768 std::vector<int> i_vec, j_vec;
769 j_vec.reserve(10000);
770 CHKERR core_ptr->createMatArrays<PetscGlobalIdx_mi_tag>(
771 p_miit, MATMPIAIJ, i_vec, j_vec, false, verb);
772
773 int nb_row_dofs = p_miit->getNbDofsRow();
774 int nb_col_dofs = p_miit->getNbDofsCol();
775 int nb_local_dofs_row = p_miit->getNbLocalDofsRow();
776 int nb_local_dofs_col = p_miit->getNbLocalDofsCol();
777
778 auto get_layout = [&]() {
787 return std::make_pair(start_ranges, end_ranges);
788 };
789
790 auto get_nnz = [&](auto &d_nnz, auto &o_nnz) {
792 auto layout = get_layout();
793 int j = 0;
794 for (int i = 0; i != nb_local_dofs_row; ++i) {
795 for (; j != i_vec[i + 1]; ++j) {
796 if (j_vec[j] < layout.second && j_vec[j] >= layout.first)
797 ++(d_nnz[i]);
798 else
799 ++(o_nnz[i]);
800 }
801 }
803 };
804
805 std::vector<int> d_nnz(nb_local_dofs_row, 0), o_nnz(nb_local_dofs_row, 0);
807
808 CHKERR MatCreate(m_field.get_comm(), Aij);
812 CHKERR MatMPIAIJSetPreallocation(*Aij, 0, &*d_nnz.begin(), 0,
813 &*o_nnz.begin());
814
817}
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j

◆ createMPIAIJ() [2/4]

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

◆ 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 
)

Definition at line 661 of file MatrixManager.cpp.

663 {
665
666 MoFEM::CoreInterface &m_field = cOre;
667 CreateRowCompressedADJMatrix *core_ptr =
668 static_cast<CreateRowCompressedADJMatrix *>(&cOre);
670
671 auto problems_ptr = m_field.get_problems();
672 auto &prb = problems_ptr->get<Problem_mi_tag>();
673 auto p_miit = prb.find(name);
674 if (p_miit == prb.end()) {
676 "problem < %s > is not found (top tip: check spelling)",
677 name.c_str());
678 }
679
680 std::vector<int> i_vec, j_vec;
681 j_vec.reserve(10000);
682 CHKERR core_ptr->createMatArrays<PetscGlobalIdx_mi_tag>(
683 p_miit, MATAIJCUSPARSE, i_vec, j_vec, false, verb);
684
685 int nb_local_dofs_row = p_miit->getNbLocalDofsRow();
686 int nb_local_dofs_col = p_miit->getNbLocalDofsCol();
687
688 auto get_layout = [&]() {
697 return std::make_pair(start_ranges, end_ranges);
698 };
699
700 auto get_nnz = [&](auto &d_nnz, auto &o_nnz) {
702 auto layout = get_layout();
703 int j = 0;
704 for (int i = 0; i != nb_local_dofs_row; ++i) {
705 for (; j != i_vec[i + 1]; ++j) {
706 if (j_vec[j] < layout.second && j_vec[j] >= layout.first)
707 ++(d_nnz[i]);
708 else
709 ++(o_nnz[i]);
710 }
711 }
713 };
714
715 std::vector<int> d_nnz(nb_local_dofs_row, 0), o_nnz(nb_local_dofs_row, 0);
717
718#ifdef PETSC_HAVE_CUDA
719 CHKERR ::MatCreateAIJCUSPARSE(m_field.get_comm(), nb_local_dofs_row,
720 nb_local_dofs_col, nb_row_dofs, nb_col_dofs, 0,
721 &*d_nnz.begin(), 0, &*o_nnz.begin(), Aij);
722#else
724 "Error: To use this matrix type compile PETSc with CUDA.");
725#endif
726
728
730}

◆ createMPIAIJCUSPARSEWithArrays() [2/4]

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

◆ 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;
84 aij_ptr.reset(aij, false);
86 }

◆ createMPIAIJWithArrays() [1/4]

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

Definition at line 624 of file MatrixManager.cpp.

625 {
627
628 MoFEM::CoreInterface &m_field = cOre;
629 CreateRowCompressedADJMatrix *core_ptr =
630 static_cast<CreateRowCompressedADJMatrix *>(&cOre);
632
633 auto problems_ptr = m_field.get_problems();
634 auto &prb = problems_ptr->get<Problem_mi_tag>();
635 auto p_miit = prb.find(name);
636 if (p_miit == prb.end()) {
638 "problem < %s > is not found (top tip: check spelling)",
639 name.c_str());
640 }
641
642 std::vector<int> i_vec, j_vec;
643 j_vec.reserve(10000);
644 CHKERR core_ptr->createMatArrays<PetscGlobalIdx_mi_tag>(
645 p_miit, MATMPIAIJ, i_vec, j_vec, false, verb);
646
647 int nb_row_dofs = p_miit->getNbDofsRow();
648 int nb_col_dofs = p_miit->getNbDofsCol();
649 int nb_local_dofs_row = p_miit->getNbLocalDofsRow();
650 int nb_local_dofs_col = p_miit->getNbLocalDofsCol();
651
652 CHKERR ::MatCreateMPIAIJWithArrays(
653 m_field.get_comm(), nb_local_dofs_row, nb_local_dofs_col, nb_row_dofs,
654 nb_col_dofs, &*i_vec.begin(), &*j_vec.begin(), PETSC_NULLPTR, Aij);
655
658}

◆ createMPIAIJWithArrays() [2/4]

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

◆ 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;
64 aij_ptr.reset(aij, false);
66 }

◆ createSeqAIJCUSPARSEWithArrays() [1/4]

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

Definition at line 733 of file MatrixManager.cpp.

735 {
737
738#ifdef PETSC_HAVE_CUDA
739 // CHKERR ::MatCreateSeqAIJCUSPARSE(MPI_Comm comm, PetscInt m, PetscInt n,
740 // PetscInt nz, const PetscInt nnz[], Mat
741 // *A);
743 "Not implemented type of matrix yet, try MPI version (aijcusparse)");
744#endif
745
747}

◆ createSeqAIJCUSPARSEWithArrays() [2/4]

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

◆ 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;
104 aij_ptr.reset(aij, false);
106 }

◆ createSeqAIJWithArrays() [1/4]

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

Definition at line 858 of file MatrixManager.cpp.

859 {
860 MoFEM::CoreInterface &m_field = cOre;
861 CreateRowCompressedADJMatrix *core_ptr =
862 static_cast<CreateRowCompressedADJMatrix *>(&cOre);
865
866 auto problems_ptr = m_field.get_problems();
867 auto &prb = problems_ptr->get<Problem_mi_tag>();
868 auto p_miit = prb.find(name);
869 if (p_miit == prb.end()) {
871 "problem < %s > is not found (top tip: check spelling)",
872 name.c_str());
873 }
874
875 std::vector<int> i_vec, j_vec;
876 j_vec.reserve(10000);
877 CHKERR core_ptr->createMatArrays<PetscGlobalIdx_mi_tag>(p_miit, MATAIJ, i_vec,
878 j_vec, false, verb);
879
880 int nb_local_dofs_row = p_miit->getNbLocalDofsRow();
881 int nb_local_dofs_col = p_miit->getNbLocalDofsCol();
882
883 double *_a;
884 CHKERR PetscMalloc(j_vec.size() * sizeof(double), &_a);
885
886 Mat tmpMat;
887 CHKERR ::MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, nb_local_dofs_row,
888 nb_local_dofs_col, &*i_vec.begin(),
889 &*j_vec.begin(), _a, &tmpMat);
892
894
897}

◆ createSeqAIJWithArrays() [2/4]

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

◆ 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;
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 600 of file MatrixManager.cpp.

601 {
602 *iface = const_cast<MatrixManager *>(this);
603 return 0;
604}
MatrixManager(const MoFEM::Core &core)

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 319 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createMPIAdjWithArrays

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createMPIAdjWithArrays
private

Definition at line 315 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createMPIAIJ

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createMPIAIJ
private

Definition at line 313 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createMPIAIJCUSPARSEWithArrays

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createMPIAIJCUSPARSEWithArrays
private

Definition at line 317 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createMPIAIJWithArrays

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createMPIAIJWithArrays
private

Definition at line 314 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createSeqAIJCUSPARSEWithArrays

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createSeqAIJCUSPARSEWithArrays
private

Definition at line 318 of file MatrixManager.hpp.

◆ MOFEM_EVENT_createSeqAIJWithArrays

PetscLogEvent MoFEM::MatrixManager::MOFEM_EVENT_createSeqAIJWithArrays
private

Definition at line 316 of file MatrixManager.hpp.


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