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

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

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

Public Member Functions

 OpSchurEndBoundary (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::FaceElementForcesAndSourcesCore::UserDataOperator
double getArea ()
 get area of face More...
 
VectorDoublegetNormal ()
 get triangle normal More...
 
VectorDoublegetTangent1 ()
 get triangle tangent 1 More...
 
VectorDoublegetTangent2 ()
 get triangle tangent 2 More...
 
auto getFTensor1Normal ()
 get normal as tensor More...
 
auto getFTensor1Tangent1 ()
 get tangentOne as tensor More...
 
auto getFTensor1Tangent2 ()
 get tangentTwo as tensor More...
 
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
VectorDoublegetCoords ()
 get triangle coordinates More...
 
auto getFTensor1Coords ()
 get get coords at gauss points More...
 
MatrixDoublegetNormalsAtGaussPts ()
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
MatrixDoublegetTangent1AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
auto getFTensor1NormalsAtGaussPts ()
 get normal at integration points More...
 
auto getFTensor1Tangent1AtGaussPts ()
 get tangent 1 at integration points More...
 
auto getFTensor1Tangent2AtGaussPts ()
 get tangent 2 at integration points More...
 
FaceElementForcesAndSourcesCoregetFaceFE ()
 return pointer to Generic Triangle Finite Element object More...
 
MoFEMErrorCode loopSideVolumes (const string fe_name, VolumeElementForcesAndSourcesCoreOnSide &fe_method)
 
- 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 []
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Definition at line 1148 of file ElasticOperators.hpp.

Constructor & Destructor Documentation

◆ OpSchurEndBoundary()

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

Definition at line 1149 of file ElasticOperators.hpp.

1152 : BoundaryEleOp(row_field, row_field, OPROW),
1153 commonDataPtr(common_data_ptr), eps(eps) {
1154 sYmm = false;
1155 }
BoundaryEle::UserDataOperator BoundaryEleOp
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::OpSchurEndBoundary::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 713 of file ElasticOperators.cpp.

714 {
716 if (row_type != MBTRI && row_type != MBQUAD)
718
719 // MoFEMFunctionReturnHot(0); //FIXME:
720
721 if (commonDataPtr->blockMatContainerPtr) {
722
723 auto SEpEp = commonDataPtr->SEpEp;
724 auto aoSEpEp = commonDataPtr->aoSEpEp;
725 auto STauTau = commonDataPtr->STauTau;
726 auto aoSTauTau = commonDataPtr->aoSTauTau;
727
728 if (SEpEp) {
729
730 auto find_block = [&](BlockMatContainer &add_bmc, auto &row_name,
731 auto &col_name, auto row_side, auto col_side,
732 auto row_type, auto col_type) {
733 return add_bmc.get<0>().find(boost::make_tuple(
734 row_name, col_name, row_type, col_type, row_side, col_side));
735 };
736
737 auto set_block = [&](BlockMatContainer &add_bmc, auto &row_name,
738 auto &col_name, auto row_side, auto col_side,
739 auto row_type, auto col_type, const auto &m,
740 const auto &row_ind, const auto &col_ind) {
741 auto it = find_block(add_bmc, row_name, col_name, row_side, col_side,
742 row_type, col_type);
743 if (it != add_bmc.get<0>().end()) {
744 it->setMat(m);
745 it->setInd(row_ind, col_ind);
746 it->setSetAtElement();
747 return it;
748 } else {
749 auto p_it = add_bmc.insert(BlockMatData(row_name, col_name, row_type,
750 col_type, row_side, col_side,
751 m, row_ind, col_ind));
752 return p_it.first;
753 }
754 };
755
756 auto add_block = [&](BlockMatContainer &add_bmc, auto &row_name,
757 auto &col_name, auto row_side, auto col_side,
758 auto row_type, auto col_type, const auto &m,
759 const auto &row_ind, const auto &col_ind) {
760 auto it = find_block(add_bmc, row_name, col_name, row_side, col_side,
761 row_type, col_type);
762 if (it != add_bmc.get<0>().end()) {
763 it->addMat(m);
764 it->setInd(row_ind, col_ind);
765 it->setSetAtElement();
766 return it;
767 } else {
768 auto p_it = add_bmc.insert(BlockMatData(row_name, col_name, row_type,
769 col_type, row_side, col_side,
770 m, row_ind, col_ind));
771 return p_it.first;
772 }
773 };
774
775 auto assemble_block = [&](auto &bit, Mat S) {
777 const VectorInt &rind = bit.rowInd;
778 const VectorInt &cind = bit.colInd;
779 const MatrixDouble &m = bit.M;
780
781 CHKERR MatSetValues(S, rind.size(), &*rind.begin(), cind.size(),
782 &*cind.begin(), &*m.data().begin(), ADD_VALUES);
783
785 };
786
787 auto create_block_schur_bc = [&](BlockMatContainer &bmc,
788 BlockMatContainer &add_bmc,
789 std::string field, AO ao) {
791 for (auto &bit : add_bmc) {
792 bit.unSetAtElement();
793 bit.clearMat();
794 }
795
796 for (auto &bit : bmc) {
797 if (bit.setAtElement && bit.rowField != field &&
798 bit.colField != field) {
799 VectorInt rind = bit.rowInd;
800 VectorInt cind = bit.colInd;
801 const MatrixDouble &m = bit.M;
802
803 // for (int i = 0; i != m.size2() * m.size1(); i++)
804 // if (abs(m.data()[i]) > 1e-8) {
805 // string wait;
806 // cout << "what" << endl;
807 // }
808
809 if (ao) {
810 CHKERR AOApplicationToPetsc(ao, rind.size(), &*rind.begin());
811 CHKERR AOApplicationToPetsc(ao, cind.size(), &*cind.begin());
812 }
813 auto it =
814 set_block(add_bmc, bit.rowField, bit.colField, bit.rowSide,
815 bit.colSide, bit.rowType, bit.colType, m, rind, cind);
816 }
817 }
818
819 // for (auto &bit_col : bmc) {
820 // if (bit_col.setAtElement && bit_col.rowField == field &&
821 // bit_col.colField != field) {
822 // const MatrixDouble &cm = bit_col.M;
823 // VectorInt cind = bit_col.colInd;
824 // // invMat.resize(inv.size1(), cm.size2(), false);
825 // // noalias(invMat) = prod(inv, cm);
826 // noalias(invMat) = cm;
827 // if (ao)
828 // CHKERR AOApplicationToPetsc(ao, cind.size(), &*cind.begin());
829 // for (auto &bit_row : bmc) {
830 // if (bit_row.setAtElement && bit_row.rowField != field &&
831 // bit_row.colField == field) {
832 // const MatrixDouble &rm = bit_row.M;
833 // VectorInt rind = bit_row.rowInd;
834 // K.resize(rind.size(), cind.size(), false);
835 // // noalias(K) = prod(rm, invMat);
836 // noalias(K) = rm;
837 // if (ao)
838 // CHKERR AOApplicationToPetsc(ao, rind.size(), &*rind.begin());
839 // auto it =
840 // add_block(add_bmc, bit_row.rowField, bit_col.colField,
841 // bit_row.rowSide, bit_col.colSide, bit_row.rowType,
842 // bit_col.colType, K, rind, cind);
843 // }
844 // }
845 // }
846 // }
847
849 };
850
851 auto assemble_schur_bc = [&](BlockMatContainer &add_bmc, Mat S,
852 bool debug = false) {
854 for (auto &bit : add_bmc) {
855 if (bit.setAtElement)
856 CHKERR assemble_block(bit, S);
857 }
858 if (debug) {
859 for (auto &bit : add_bmc) {
860 if (bit.setAtElement) {
861 std::cerr << "assemble: " << bit.rowField << " " << bit.colField
862 << endl;
863 // std::cerr << bit.M << endl;
864 }
865 }
866 std::cerr << std::endl;
867 }
869 };
870
871 auto precondition_schur = [&](BlockMatContainer &bmc,
872 BlockMatContainer &add_bmc,
873 const std::string field,
874 const MatrixDouble &diag_mat,
875 const double eps) {
877
878 for (auto &bit : add_bmc) {
879 bit.unSetAtElement();
880 bit.clearMat();
881 }
882
883 for (auto &bit : bmc) {
884 if (bit.setAtElement) {
885 if (bit.rowField != field || bit.colField != field)
886 auto it = set_block(add_bmc, bit.rowField, bit.colField,
887 bit.rowSide, bit.colSide, bit.rowType,
888 bit.colType, bit.M, bit.rowInd, bit.colInd);
889 }
890 }
891
892 auto bit =
893 bmc.get<1>().find(boost::make_tuple(field, field, row_type, row_type));
894 if (bit->setAtElement && bit != bmc.get<1>().end()) {
895 auto it = set_block(add_bmc, bit->rowField, bit->colField,
896 bit->rowSide, bit->colSide, bit->rowType,
897 bit->colType, bit->M, bit->rowInd, bit->colInd);
898 MatrixDouble &m = const_cast<MatrixDouble &>(it->M);
899 m += eps * diag_mat;
900 } else {
901 auto row_it = bmc.get<3>().lower_bound(field);
902 for (; row_it != bmc.get<3>().upper_bound(field); ++row_it) {
903 if (row_it->setAtElement) {
904 auto it = set_block(add_bmc, field, field, 0, 0, row_type, row_type,
905 diag_mat, row_it->rowInd, row_it->rowInd);
906 MatrixDouble &m = const_cast<MatrixDouble &>(it->M);
907 m *= eps;
908 break;
909 }
910 }
911 if (row_it == bmc.get<3>().end())
912 SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
913 "row field not found %s", field.c_str());
914 }
915
917 };
918
919 const bool debug = false;
920 auto &mat_cont = *commonDataPtr->blockMatContainerPtr;
921 auto &mass_mat_tensor = commonDataPtr->massMatTensor;
922
923 if (SEpEp) {
924 CHKERR create_block_schur_bc(mat_cont, blockMat["BC"], "EP", aoSEpEp);
925 // if (getPtrFE()->mField.check_field("SIGMA") && false) {
926 // constexpr double eps = 1e-8;
927 // CHKERR precondition_schur(blockMat["BC"], blockMat["precBC"],
928 // "SIGMA",
929 // mass_mat_tensor, eps);
930 // CHKERR assemble_schur_bc(blockMat["precBC"], SEpEp, debug);
931 // blockMat["BC"] = blockMat["precBC"];
932 // // CHKERR assemble_schur_bc(blockMat["BC"], SEpEp, debug);
933 // } else
934 CHKERR assemble_schur_bc(blockMat["BC"], SEpEp, debug);
935
936 if (STauTau) {
937 CHKERR create_block_schur_bc(blockMat["BC"], blockMat["TAUTAU"],
938 "TAU", aoSTauTau);
939 CHKERR assemble_schur_bc(blockMat["TAUTAU"], STauTau, debug);
940 }
941 }
942
943 blockMat.clear();
944 }
945 }
946
948}
#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
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, BlockMatContainer > blockMat

Member Data Documentation

◆ blockMat

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

Definition at line 1168 of file ElasticOperators.hpp.

◆ commonDataPtr

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

Definition at line 1169 of file ElasticOperators.hpp.

◆ eps

const double OpElasticTools::OpSchurEndBoundary::eps
private

Definition at line 1160 of file ElasticOperators.hpp.

◆ invBlockMat

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

Definition at line 1167 of file ElasticOperators.hpp.

◆ invMat

MatrixDouble OpElasticTools::OpSchurEndBoundary::invMat
private

Definition at line 1162 of file ElasticOperators.hpp.

◆ iPIV

VectorInt OpElasticTools::OpSchurEndBoundary::iPIV
private

Definition at line 1164 of file ElasticOperators.hpp.

◆ K

MatrixDouble OpElasticTools::OpSchurEndBoundary::K
private

Definition at line 1163 of file ElasticOperators.hpp.

◆ lapackWork

VectorDouble OpElasticTools::OpSchurEndBoundary::lapackWork
private

Definition at line 1165 of file ElasticOperators.hpp.


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