v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
BoneRemodeling::OpAssmbleRhoLhs_dF Struct Reference

Off diagonal block of tangent matrix \(K_{\rho u}\). More...

#include <users_modules/bone_remodelling/src/Remodeling.hpp>

Inheritance diagram for BoneRemodeling::OpAssmbleRhoLhs_dF:
[legend]
Collaboration diagram for BoneRemodeling::OpAssmbleRhoLhs_dF:
[legend]

Public Member Functions

 OpAssmbleRhoLhs_dF (Remodeling::CommonData &common_data)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 
- 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...
 

Public Attributes

Remodeling::CommonDatacommonData
 
MatrixDouble locK_rho_F
 
- 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...
 

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

Off diagonal block of tangent matrix \(K_{\rho u}\).

\[ K_{\rho u}=-\intop_{V} c \left[\frac{\rho_{0}}{\rho_{0}^{\ast}}\right]^{n-m} \nabla N_j P_{ij} N_i \,dV \]

Definition at line 493 of file Remodeling.hpp.

Constructor & Destructor Documentation

◆ OpAssmbleRhoLhs_dF()

BoneRemodeling::OpAssmbleRhoLhs_dF::OpAssmbleRhoLhs_dF ( Remodeling::CommonData common_data)
Examples
Remodeling.cpp.

Definition at line 1066 of file Remodeling.cpp.

1067 : VolumeElementForcesAndSourcesCore::UserDataOperator(
1068 "RHO", "DISPLACEMENTS", UserDataOperator::OPROWCOL),
1069 commonData(common_data) {
1070 sYmm = false;
1071}
Remodeling::CommonData & commonData
Definition: Remodeling.hpp:496
bool sYmm
If true assume that matrix is symmetric structure.

Member Function Documentation

◆ doWork()

MoFEMErrorCode BoneRemodeling::OpAssmbleRhoLhs_dF::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
DataForcesAndSourcesCore::EntData &  row_data,
DataForcesAndSourcesCore::EntData &  col_data 
)
Examples
Remodeling.cpp, and Remodeling.hpp.

Definition at line 1074 of file Remodeling.cpp.

1077 {
1078
1080
1081 const int row_nb_dofs = row_data.getIndices().size();
1082 if (!row_nb_dofs)
1084 const int col_nb_dofs = col_data.getIndices().size();
1085 if (!col_nb_dofs)
1087
1088 // Set size can clear local tangent matrix
1089 locK_rho_F.resize(row_nb_dofs, col_nb_dofs, false);
1090 locK_rho_F.clear();
1091
1092 const int row_nb_gauss_pts = row_data.getN().size1();
1093 const int col_nb_base_functions = col_data.getN().size2();
1094
1095 const double c = commonData.c;
1096 const double n = commonData.n;
1097 const double m = commonData.m;
1098 const double rho_ref = commonData.rHo_ref;
1099
1102
1104 auto P = getFTensor2FromMat<3, 3>(commonData.data.matP);
1105
1107 auto col_diff_base_functions = col_data.getFTensor1DiffN<3>();
1108
1109 for (int gg = 0; gg != row_nb_gauss_pts; gg++) {
1110
1111 // Get volume and integration weight
1112 double w = getVolume() * getGaussPts()(3, gg);
1113 const double kp = commonData.getCFromDensity(rho);
1114 const double a = w * c * kp * pow(rho / rho_ref, n - m);
1115
1116 int col_bb = 0;
1117 for (; col_bb != col_nb_dofs / 3; col_bb++) {
1118 t1(i) = a * P(i, I) * col_diff_base_functions(I);
1119 auto row_base_functions = row_data.getFTensor0N(gg, 0);
1120 for (int row_bb = 0; row_bb != row_nb_dofs; row_bb++) {
1121 FTensor::Tensor1<double *, 3> k(&locK_rho_F(row_bb, 3 * col_bb + 0),
1122 &locK_rho_F(row_bb, 3 * col_bb + 1),
1123 &locK_rho_F(row_bb, 3 * col_bb + 2));
1124 k(i) -= row_base_functions * t1(i);
1125 ++row_base_functions;
1126 }
1127 ++col_diff_base_functions;
1128 }
1129 for (; col_bb != col_nb_base_functions; col_bb++) {
1130 ++col_diff_base_functions;
1131 }
1132
1133 ++P;
1134 ++rho;
1135 }
1136
1137 CHKERR MatSetValues(getFEMethod()->ts_B, row_nb_dofs,
1138 &*row_data.getIndices().begin(), col_nb_dofs,
1139 &*col_data.getIndices().begin(), &locK_rho_F(0, 0),
1140 ADD_VALUES);
1141
1143}
constexpr double a
#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
#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
FTensor::Index< 'n', SPACE_DIM > n
FTensor::Index< 'm', SPACE_DIM > m
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
FTensor::Index< 'k', 3 > k
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
double w(const double g, const double t)
constexpr IntegrationType I
double rho
Definition: plastic.cpp:191
boost::shared_ptr< VectorDouble > rhoPtr
Ptr to density matrix container.
Definition: Remodeling.hpp:187
double getCFromDensity(const double &rho)
Definition: Remodeling.cpp:218
double n
porosity exponent [-]
Definition: Remodeling.hpp:148
double m
algorithmic exponent [-]
Definition: Remodeling.hpp:147
double c
density evolution (growth) velocity [d/m^2]
Definition: Remodeling.hpp:146
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

Member Data Documentation

◆ commonData

Remodeling::CommonData& BoneRemodeling::OpAssmbleRhoLhs_dF::commonData
Examples
Remodeling.hpp.

Definition at line 496 of file Remodeling.hpp.

◆ locK_rho_F

MatrixDouble BoneRemodeling::OpAssmbleRhoLhs_dF::locK_rho_F
Examples
Remodeling.hpp.

Definition at line 497 of file Remodeling.hpp.


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