v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
OpElasticTools::OpSchurEnd Struct Reference

#include <users_modules/multifield_plasticity/src/ElasticOperators.hpp>

Inheritance diagram for OpElasticTools::OpSchurEnd:
[legend]
Collaboration diagram for OpElasticTools::OpSchurEnd:
[legend]

Public Member Functions

 OpSchurEnd (const std::string row_field, boost::shared_ptr< CommonData > common_data_ptr, const double eps=1e-12)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Private Attributes

const double eps
 
MatrixDouble invMat
 
MatrixDouble K
 
VectorInt iPIV
 
VectorDouble lapackWork
 
map< std::string, MatrixDouble > invBlockMat
 
map< std::string, BlockMatContainer > blockMat
 
boost::shared_ptr< CommonDatacommonDataPtr
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > >
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Definition at line 1113 of file ElasticOperators.hpp.

Constructor & Destructor Documentation

◆ OpSchurEnd()

OpElasticTools::OpSchurEnd::OpSchurEnd ( const std::string  row_field,
boost::shared_ptr< CommonData common_data_ptr,
const double  eps = 1e-12 
)
inline

Definition at line 1114 of file ElasticOperators.hpp.

1117 : DomainEleOp(row_field, row_field, OPROW),
1118 commonDataPtr(common_data_ptr), eps(eps) {
1119 sYmm = false;
1120 }
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
bool sYmm
If true assume that matrix is symmetric structure.
@ OPROW
operator doWork function is executed on FE rows
boost::shared_ptr< CommonData > commonDataPtr

Member Function Documentation

◆ doWork()

MoFEMErrorCode OpElasticTools::OpSchurEnd::doWork ( int  side,
EntityType  type,
EntData data 
)
virtual

Operator for linear form, usually to calculate values on right hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 964 of file ElasticOperators.cpp.

965 {
967 if (row_type != MBTET && row_type != MBHEX)
969 if (commonDataPtr->blockMatContainerPtr) {
970
971 auto SEpEp = commonDataPtr->SEpEp;
972 auto aoSEpEp = commonDataPtr->aoSEpEp;
973 auto STauTau = commonDataPtr->STauTau;
974 auto aoSTauTau = commonDataPtr->aoSTauTau;
975
976 if (SEpEp) {
977
978 auto find_block = [&](BlockMatContainer &add_bmc, auto &row_name,
979 auto &col_name, auto row_side, auto col_side,
980 auto row_type, auto col_type) {
981 return add_bmc.get<0>().find(boost::make_tuple(
982 row_name, col_name, row_type, col_type, row_side, col_side));
983 };
984
985 auto set_block = [&](BlockMatContainer &add_bmc, auto &row_name,
986 auto &col_name, auto row_side, auto col_side,
987 auto row_type, auto col_type, const auto &m,
988 const auto &row_ind, const auto &col_ind) {
989 auto it = find_block(add_bmc, row_name, col_name, row_side, col_side,
990 row_type, col_type);
991 if (it != add_bmc.get<0>().end()) {
992 it->setMat(m);
993 it->setInd(row_ind, col_ind);
994 it->setSetAtElement();
995 return it;
996 } else {
997 auto p_it = add_bmc.insert(BlockMatData(row_name, col_name, row_type,
998 col_type, row_side, col_side,
999 m, row_ind, col_ind));
1000 return p_it.first;
1001 }
1002 };
1003
1004 auto add_block = [&](BlockMatContainer &add_bmc, auto &row_name,
1005 auto &col_name, auto row_side, auto col_side,
1006 auto row_type, auto col_type, const auto &m,
1007 const auto &row_ind, const auto &col_ind) {
1008 auto it = find_block(add_bmc, row_name, col_name, row_side, col_side,
1009 row_type, col_type);
1010 if (it != add_bmc.get<0>().end()) {
1011 it->addMat(m);
1012 it->setInd(row_ind, col_ind);
1013 it->setSetAtElement();
1014 return it;
1015 } else {
1016 auto p_it = add_bmc.insert(BlockMatData(row_name, col_name, row_type,
1017 col_type, row_side, col_side,
1018 m, row_ind, col_ind));
1019 return p_it.first;
1020 }
1021 };
1022
1023 auto assemble_block = [&](auto &bit, Mat S) {
1025 const VectorInt &rind = bit.rowInd;
1026 const VectorInt &cind = bit.colInd;
1027 const MatrixDouble &m = bit.M;
1028
1029 CHKERR MatSetValues(S, rind.size(), &*rind.begin(), cind.size(),
1030 &*cind.begin(), &*m.data().begin(), ADD_VALUES);
1031
1033 };
1034
1035 auto invert_symm_mat = [&](MatrixDouble &m, auto &inv) {
1037 const int nb = m.size1();
1038
1039 inv.resize(nb, nb, false);
1040 inv.clear();
1041 for (int cc = 0; cc != nb; ++cc)
1042 inv(cc, cc) = -1;
1043
1044 iPIV.resize(nb, false);
1045 lapackWork.resize(nb * nb, false);
1046 const auto info = lapack_dsysv('L', nb, nb, &*m.data().begin(), nb,
1047 &*iPIV.begin(), &*inv.data().begin(), nb,
1048 &*lapackWork.begin(), nb * nb);
1049 // if (info != 0)
1050 // SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
1051 // "Can not invert matrix info = %d", info);
1052
1054 };
1055
1056 auto invert_symm_schur = [&](BlockMatContainer &bmc, std::string field,
1057 auto &inv) {
1059
1060 auto bit =
1061 bmc.get<1>().find(boost::make_tuple(field, field, row_type, row_type));
1062 if (bit != bmc.get<1>().end()) {
1063
1064 auto &m = *const_cast<MatrixDouble *>(&(bit->M));
1065 CHKERR invert_symm_mat(m, inv);
1066
1067 } else
1068 SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
1069 "%s matrix not found", field.c_str());
1070
1072 };
1073
1074 // auto invert_nonsymm_mat = [&](MatrixDouble &m, auto &inv) {
1075 // MoFEMFunctionBeginHot;
1076
1077 // const int nb = m.size1();
1078
1079 // MatrixDouble trans_m = trans(m);
1080 // MatrixDouble trans_inv;
1081 // trans_inv.resize(nb, nb, false);
1082 // trans_inv.clear();
1083 // for (int c = 0; c != nb; ++c)
1084 // trans_inv(c, c) = -1;
1085
1086 // iPIV.resize(nb, false);
1087 // const auto info =
1088 // lapack_dgesv(nb, nb, &*trans_m.data().begin(), nb,
1089 // &*iPIV.begin(),
1090 // &*trans_inv.data().begin(), nb);
1091 // if (info != 0)
1092 // SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
1093 // "Can not invert matrix info = %d", info);
1094
1095 // inv.resize(nb, nb, false);
1096 // noalias(inv) = trans(trans_inv);
1097
1098 // MoFEMFunctionReturnHot(0);
1099 // };
1100
1101 // auto invert_nonsymm_schur = [&](BlockMatContainer &bmc,
1102 // std::string field, auto &inv,
1103 // const bool debug = false) {
1104 // MoFEMFunctionBeginHot;
1105
1106 // auto bit =
1107 // bmc.get<1>().find(boost::make_tuple(field, field, row_type, row_type));
1108 // if (bit != bmc.get<1>().end()) {
1109
1110 // auto &m = *const_cast<MatrixDouble *>(&(bit->M));
1111 // CHKERR invert_nonsymm_mat(m, inv);
1112
1113 // if (debug) {
1114 // std::cerr << prod(m, inv) << endl;
1115 // std::cerr << endl;
1116 // }
1117
1118 // } else
1119 // SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
1120 // "%s matrix not found", field.c_str());
1121
1122 // MoFEMFunctionReturnHot(0);
1123 // };
1124
1125 auto create_block_schur = [&](BlockMatContainer &bmc,
1126 BlockMatContainer &add_bmc,
1127 std::string field, AO ao,
1128 MatrixDouble &inv) {
1130 // AOView(ao, PETSC_VIEWER_STDOUT_SELF);
1131 for (auto &bit : add_bmc) {
1132 bit.unSetAtElement();
1133 bit.clearMat();
1134 }
1135
1136 for (auto &bit : bmc) {
1137 if (bit.setAtElement && bit.rowField != field &&
1138 bit.colField != field) {
1139 VectorInt rind = bit.rowInd;
1140 VectorInt cind = bit.colInd;
1141 const MatrixDouble &m = bit.M;
1142 if (ao) {
1143 CHKERR AOApplicationToPetsc(ao, rind.size(), &*rind.begin());
1144 CHKERR AOApplicationToPetsc(ao, cind.size(), &*cind.begin());
1145 }
1146 auto it =
1147 set_block(add_bmc, bit.rowField, bit.colField, bit.rowSide,
1148 bit.colSide, bit.rowType, bit.colType, m, rind, cind);
1149 }
1150 }
1151
1152 for (auto &bit_col : bmc) {
1153 if (bit_col.setAtElement && bit_col.rowField == field &&
1154 bit_col.colField != field) {
1155 const MatrixDouble &cm = bit_col.M;
1156 VectorInt cind = bit_col.colInd;
1157 invMat.resize(inv.size1(), cm.size2(), false);
1158 noalias(invMat) = prod(inv, cm);
1159 if (ao)
1160 CHKERR AOApplicationToPetsc(ao, cind.size(), &*cind.begin());
1161 for (auto &bit_row : bmc) {
1162 if (bit_row.setAtElement && bit_row.rowField != field &&
1163 bit_row.colField == field) {
1164 const MatrixDouble &rm = bit_row.M;
1165 VectorInt rind = bit_row.rowInd;
1166 K.resize(rind.size(), cind.size(), false);
1167 noalias(K) = prod(rm, invMat);
1168 if (ao)
1169 CHKERR AOApplicationToPetsc(ao, rind.size(), &*rind.begin());
1170 auto it =
1171 add_block(add_bmc, bit_row.rowField, bit_col.colField,
1172 bit_row.rowSide, bit_col.colSide, bit_row.rowType,
1173 bit_col.colType, K, rind, cind);
1174 }
1175 }
1176 }
1177 }
1178
1180 };
1181
1182 auto assemble_schur = [&](BlockMatContainer &add_bmc, Mat S,
1183 bool debug = false) {
1185 if (debug) {
1186 for (auto &bit : add_bmc) {
1187 if (bit.setAtElement) {
1188 std::cerr << "assemble vol: " << bit.rowField << " "
1189 << bit.colField << endl;
1190 // std::cerr << bit.M << endl;
1191 }
1192 }
1193 std::cerr << std::endl;
1194 }
1195 for (auto &bit : add_bmc) {
1196 if (bit.setAtElement)
1197 CHKERR assemble_block(bit, S);
1198 }
1200 };
1201
1202 auto precondition_schur = [&](BlockMatContainer &bmc,
1203 BlockMatContainer &add_bmc,
1204 const std::string field,
1205 const MatrixDouble &diag_mat,
1206 const double eps) {
1208
1209 for (auto &bit : add_bmc) {
1210 bit.unSetAtElement();
1211 bit.clearMat();
1212 }
1213
1214 for (auto &bit : bmc) {
1215 if (bit.setAtElement) {
1216 if (bit.rowField != field || bit.colField != field)
1217 auto it = set_block(add_bmc, bit.rowField, bit.colField,
1218 bit.rowSide, bit.colSide, bit.rowType,
1219 bit.colType, bit.M, bit.rowInd, bit.colInd);
1220 }
1221 }
1222
1223 auto bit =
1224 bmc.get<1>().find(boost::make_tuple(field, field, row_type, row_type));
1225 if (bit->setAtElement && bit != bmc.get<1>().end()) {
1226 auto it = set_block(add_bmc, bit->rowField, bit->colField,
1227 bit->rowSide, bit->colSide, bit->rowType,
1228 bit->colType, bit->M, bit->rowInd, bit->colInd);
1229 MatrixDouble &m = const_cast<MatrixDouble &>(it->M);
1230 m += eps * diag_mat;
1231 } else {
1232 auto row_it = bmc.get<3>().lower_bound(field);
1233 for (; row_it != bmc.get<3>().upper_bound(field); ++row_it) {
1234 if (row_it->setAtElement) {
1235 auto it = set_block(add_bmc, field, field, 0, 0, row_type, row_type,
1236 diag_mat, row_it->rowInd, row_it->rowInd);
1237 MatrixDouble &m = const_cast<MatrixDouble &>(it->M);
1238 m *= eps;
1239 break;
1240 }
1241 }
1242 if (row_it == bmc.get<3>().end())
1243 SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
1244 "row field not found %s", field.c_str());
1245 }
1246
1248 };
1249
1250 const bool debug = false;
1251 auto &mat_cont = *commonDataPtr->blockMatContainerPtr;
1252
1253 CHKERR invert_symm_schur(mat_cont, "EP",
1254 invBlockMat["EPEP"]);
1255 CHKERR create_block_schur(mat_cont,
1256 blockMat["EPEP"], "EP", aoSEpEp,
1257 invBlockMat["EPEP"]);
1258
1259 auto &mass_mat_scalar = commonDataPtr->massMatScalar;
1260 constexpr double eps = 1e-8; //FIXME: TODO: figure out this later
1261 CHKERR precondition_schur(blockMat["EPEP"], blockMat["precEPEP"], "TAU",
1262 mass_mat_scalar, eps);
1263 // if (getPtrFE()->mField.check_field("SIGMA") && false) {
1264 // auto &mass_mat_tensor = commonDataPtr->massMatTensor;
1265 // CHKERR precondition_schur(blockMat["precEPEP"], blockMat["precSIGMA"],
1266 // "SIGMA", mass_mat_tensor, eps);
1267 // CHKERR assemble_schur(blockMat["precSIGMA"], SEpEp, debug);
1268 // blockMat["precEPEP"] = blockMat["precSIGMA"];
1269 // } else
1270 CHKERR assemble_schur(blockMat["precEPEP"], SEpEp, false);
1271
1272 if(STauTau){
1273 CHKERR invert_symm_schur(blockMat["precEPEP"], "TAU",
1274 invBlockMat["TAUTAU"]);
1275 CHKERR create_block_schur(blockMat["precEPEP"], blockMat["TAUTAU"], "TAU",
1276 aoSTauTau, invBlockMat["TAUTAU"]);
1277 CHKERR assemble_schur(blockMat["TAUTAU"], STauTau);
1278 }
1279
1280 blockMat.clear();
1281
1282 }
1283 }
1284
1286}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
static const bool debug
FTensor::Index< 'm', SPACE_DIM > m
auto bit
set bit
static __CLPK_integer lapack_dsysv(char uplo, __CLPK_integer n, __CLPK_integer nrhs, __CLPK_doublereal *a, __CLPK_integer lda, __CLPK_integer *ipiv, __CLPK_doublereal *b, __CLPK_integer ldb, __CLPK_doublereal *work, __CLPK_integer lwork)
Definition: lapack_wrap.h:220
UBlasVector< int > VectorInt
Definition: Types.hpp:67
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
multi_index_container< BlockMatData, indexed_by< ordered_unique< composite_key< BlockMatData, member< BlockMatData, std::string, &BlockMatData::rowField >, member< BlockMatData, std::string, &BlockMatData::colField >, member< BlockMatData, EntityType, &BlockMatData::rowType >, member< BlockMatData, EntityType, &BlockMatData::colType >, member< BlockMatData, int, &BlockMatData::rowSide >, member< BlockMatData, int, &BlockMatData::colSide > > >, ordered_non_unique< composite_key< BlockMatData, member< BlockMatData, std::string, &BlockMatData::rowField >, member< BlockMatData, std::string, &BlockMatData::colField >, member< BlockMatData, EntityType, &BlockMatData::rowType >, member< BlockMatData, EntityType, &BlockMatData::colType > > >, ordered_non_unique< composite_key< BlockMatData, member< BlockMatData, std::string, &BlockMatData::rowField >, member< BlockMatData, std::string, &BlockMatData::colField > > >, ordered_non_unique< member< BlockMatData, std::string, &BlockMatData::rowField > >, ordered_non_unique< member< BlockMatData, std::string, &BlockMatData::colField > > > > BlockMatContainer
map< std::string, MatrixDouble > invBlockMat
map< std::string, BlockMatContainer > blockMat

Member Data Documentation

◆ blockMat

map<std::string, BlockMatContainer> OpElasticTools::OpSchurEnd::blockMat
private

Definition at line 1133 of file ElasticOperators.hpp.

◆ commonDataPtr

boost::shared_ptr<CommonData> OpElasticTools::OpSchurEnd::commonDataPtr
private

Definition at line 1134 of file ElasticOperators.hpp.

◆ eps

const double OpElasticTools::OpSchurEnd::eps
private

Definition at line 1125 of file ElasticOperators.hpp.

◆ invBlockMat

map<std::string, MatrixDouble> OpElasticTools::OpSchurEnd::invBlockMat
private

Definition at line 1132 of file ElasticOperators.hpp.

◆ invMat

MatrixDouble OpElasticTools::OpSchurEnd::invMat
private

Definition at line 1127 of file ElasticOperators.hpp.

◆ iPIV

VectorInt OpElasticTools::OpSchurEnd::iPIV
private

Definition at line 1129 of file ElasticOperators.hpp.

◆ K

MatrixDouble OpElasticTools::OpSchurEnd::K
private

Definition at line 1128 of file ElasticOperators.hpp.

◆ lapackWork

VectorDouble OpElasticTools::OpSchurEnd::lapackWork
private

Definition at line 1130 of file ElasticOperators.hpp.


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