v0.14.0
Public Member Functions | Public Attributes | List of all members
SmallStrainPlasticity::OpAssembleLhs Struct Reference

#include <users_modules/small_strain_plasticity/src/SmallStrainPlasticity.hpp>

Inheritance diagram for SmallStrainPlasticity::OpAssembleLhs:
[legend]
Collaboration diagram for SmallStrainPlasticity::OpAssembleLhs:
[legend]

Public Member Functions

 OpAssembleLhs (string field_name, CommonData &common_data)
 
PetscErrorCode 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 Member Functions inherited from SmallStrainPlasticity::MakeB
PetscErrorCode makeB (const MatrixAdaptor &diffN, MatrixDouble &B)
 
PetscErrorCode addVolumetricB (const MatrixAdaptor &diffN, MatrixDouble &volB, double alpha)
 

Public Attributes

CommonDatacommonData
 
MatrixDouble K
 
MatrixDouble transK
 
MatrixDouble rowB
 
MatrixDouble colB
 
MatrixDouble CB
 
MatrixDouble rowVolB
 
MatrixDouble colVolB
 
- 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)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Definition at line 456 of file SmallStrainPlasticity.hpp.

Constructor & Destructor Documentation

◆ OpAssembleLhs()

SmallStrainPlasticity::OpAssembleLhs::OpAssembleLhs ( string  field_name,
CommonData common_data 
)
Examples
SmallStrainPlasticity.hpp.

Definition at line 302 of file SmallStrainPlasticity.cpp.

Member Function Documentation

◆ doWork()

PetscErrorCode SmallStrainPlasticity::OpAssembleLhs::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
DataForcesAndSourcesCore::EntData row_data,
DataForcesAndSourcesCore::EntData col_data 
)

< higher order geometry

< higher order geometry

Examples
SmallStrainPlasticity.hpp.

Definition at line 311 of file SmallStrainPlasticity.cpp.

316  {
317  PetscFunctionBegin;
318  try {
319 
320  int nb_dofs_row = row_data.getFieldData().size();
321  if(nb_dofs_row==0) PetscFunctionReturn(0);
322  int nb_dofs_col = col_data.getFieldData().size();
323  if(nb_dofs_col==0) PetscFunctionReturn(0);
324 
325  K.resize(nb_dofs_row,nb_dofs_col,false);
326  K.clear();
327 
328  int nb_gauss_pts = row_data.getN().size1();
329 
330  if(commonData.bBar) {
331  double v = 0;
332  rowVolB.resize(6,nb_dofs_row,false);
333  rowVolB.clear();
334  colVolB.resize(6,nb_dofs_col,false);
335  colVolB.clear();
336  for(int gg = 0;gg<nb_gauss_pts;gg++) {
337  double val = getVolume()*getGaussPts()(3,gg);
338  if(getHoGaussPtsDetJac().size()>0) {
339  val *= getHoGaussPtsDetJac()[gg]; ///< higher order geometry
340  }
341  v += val;
342  const MatrixAdaptor &diffN_row = row_data.getDiffN(gg,nb_dofs_row/3);
343  const MatrixAdaptor &diffN_col = col_data.getDiffN(gg,nb_dofs_col/3);
344  ierr = addVolumetricB(diffN_row,rowVolB,val); CHKERRQ(ierr);
345  ierr = addVolumetricB(diffN_col,colVolB,val); CHKERRQ(ierr);
346  }
347  rowVolB /= v;
348  colVolB /= v;
349  }
350 
351  for(int gg = 0;gg<nb_gauss_pts;gg++) {
352 
353  double val = getVolume()*getGaussPts()(3,gg);
354  if(getHoGaussPtsDetJac().size()>0) {
355  val *= getHoGaussPtsDetJac()[gg]; ///< higher order geometry
356  }
357 
358  const MatrixAdaptor &diffN_row = row_data.getDiffN(gg,nb_dofs_row/3);
359  const MatrixAdaptor &diffN_col = col_data.getDiffN(gg,nb_dofs_col/3);
360  ierr = makeB(diffN_row,rowB); CHKERRQ(ierr);
361  ierr = makeB(diffN_col,colB); CHKERRQ(ierr);
362  if(commonData.bBar) {
363  ierr = addVolumetricB(diffN_row,rowB,-1); CHKERRQ(ierr); // B-bar
364  ierr = addVolumetricB(diffN_col,colB,-1); CHKERRQ(ierr);
365  rowB += rowVolB;
366  colB += colVolB;
367  }
368 
369  CB.resize(6,nb_dofs_col,false);
370  cblas_dgemm(
371  CblasRowMajor,CblasNoTrans,CblasNoTrans,
372  6,nb_dofs_col,6,
373  1.,
375  &colB(0,0),nb_dofs_col,
376  0,
377  &CB(0,0),nb_dofs_col
378  );
379  cblas_dgemm(
380  CblasRowMajor,CblasTrans,CblasNoTrans,
381  nb_dofs_row,nb_dofs_col,6,
382  val,
383  &rowB(0,0),nb_dofs_row,
384  &CB(0,0),nb_dofs_col,
385  1,
386  &K(0,0),nb_dofs_col
387  );
388  //noalias(CB) = prod(commonData.materialTangentOperator[gg],colB);
389  //noalias(K) += val*prod(trans(rowB),CB);
390 
391  }
392 
393  // cerr << K << endl;
394  // cerr << row_data.getIndices() << endl;
395  // cerr << col_data.getIndices() << endl;
396  // cerr << getFEMethod()->snes_A << endl;
397  // cerr << getFEMethod()->snes_B << endl;
398 
399  int *row_indices_ptr,*col_indices_ptr;
400  row_indices_ptr = &row_data.getIndices()[0];
401  col_indices_ptr = &col_data.getIndices()[0];
402 
403  ierr = MatSetValues(
404  getFEMethod()->snes_B,
405  nb_dofs_row,row_indices_ptr,
406  nb_dofs_col,col_indices_ptr,
407  &K(0,0),ADD_VALUES
408  ); CHKERRQ(ierr);
409 
410  //is symmetric
411  if(row_side != col_side || row_type != col_type) {
412 
413  transK.resize(nb_dofs_col,nb_dofs_row,false);
414  noalias(transK) = trans(K);
415  ierr = MatSetValues(
416  getFEMethod()->snes_B,
417  nb_dofs_col,col_indices_ptr,
418  nb_dofs_row,row_indices_ptr,
419  &transK(0,0),ADD_VALUES
420  ); CHKERRQ(ierr);
421 
422  }
423 
424  } catch (const std::exception& ex) {
425  ostringstream ss;
426  ss << "throw in method: " << ex.what() << endl;
427  SETERRQ(PETSC_COMM_SELF,1,ss.str().c_str());
428  }
429  PetscFunctionReturn(0);
430 }

Member Data Documentation

◆ CB

MatrixDouble SmallStrainPlasticity::OpAssembleLhs::CB
Examples
SmallStrainPlasticity.hpp.

Definition at line 460 of file SmallStrainPlasticity.hpp.

◆ colB

MatrixDouble SmallStrainPlasticity::OpAssembleLhs::colB
Examples
SmallStrainPlasticity.hpp.

Definition at line 460 of file SmallStrainPlasticity.hpp.

◆ colVolB

MatrixDouble SmallStrainPlasticity::OpAssembleLhs::colVolB
Examples
SmallStrainPlasticity.hpp.

Definition at line 461 of file SmallStrainPlasticity.hpp.

◆ commonData

CommonData& SmallStrainPlasticity::OpAssembleLhs::commonData
Examples
SmallStrainPlasticity.hpp.

Definition at line 458 of file SmallStrainPlasticity.hpp.

◆ K

MatrixDouble SmallStrainPlasticity::OpAssembleLhs::K
Examples
SmallStrainPlasticity.hpp.

Definition at line 459 of file SmallStrainPlasticity.hpp.

◆ rowB

MatrixDouble SmallStrainPlasticity::OpAssembleLhs::rowB
Examples
SmallStrainPlasticity.hpp.

Definition at line 460 of file SmallStrainPlasticity.hpp.

◆ rowVolB

MatrixDouble SmallStrainPlasticity::OpAssembleLhs::rowVolB
Examples
SmallStrainPlasticity.hpp.

Definition at line 461 of file SmallStrainPlasticity.hpp.

◆ transK

MatrixDouble SmallStrainPlasticity::OpAssembleLhs::transK
Examples
SmallStrainPlasticity.hpp.

Definition at line 459 of file SmallStrainPlasticity.hpp.


The documentation for this struct was generated from the following files:
SmallStrainPlasticity::OpAssembleLhs::colVolB
MatrixDouble colVolB
Definition: SmallStrainPlasticity.hpp:461
MoFEM::MatSetValues
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
Definition: EntitiesFieldData.hpp:1631
SmallStrainPlasticity::OpAssembleLhs::transK
MatrixDouble transK
Definition: SmallStrainPlasticity.hpp:459
SmallStrainPlasticity::OpAssembleLhs::colB
MatrixDouble colB
Definition: SmallStrainPlasticity.hpp:460
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getVolume
double getVolume() const
element volume (linear geometry)
Definition: VolumeElementForcesAndSourcesCore.hpp:161
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
SmallStrainPlasticity::OpAssembleLhs::K
MatrixDouble K
Definition: SmallStrainPlasticity.hpp:459
MoFEM::ForcesAndSourcesCore::UserDataOperator::getGaussPts
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Definition: ForcesAndSourcesCore.hpp:1235
SmallStrainPlasticity::OpAssembleLhs::rowVolB
MatrixDouble rowVolB
Definition: SmallStrainPlasticity.hpp:461
SmallStrainPlasticity::OpAssembleLhs::rowB
MatrixDouble rowB
Definition: SmallStrainPlasticity.hpp:460
SmallStrainPlasticity::CommonData::materialTangentOperator
vector< MatrixDouble > materialTangentOperator
Definition: SmallStrainPlasticity.hpp:380
SmallStrainPlasticity::MakeB::addVolumetricB
PetscErrorCode addVolumetricB(const MatrixAdaptor &diffN, MatrixDouble &volB, double alpha)
Definition: SmallStrainPlasticity.cpp:200
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEMethod
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Definition: ForcesAndSourcesCore.hpp:1041
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEM::Types::MatrixAdaptor
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
Definition: Types.hpp:132
SmallStrainPlasticity::CommonData::bBar
bool bBar
Definition: SmallStrainPlasticity.hpp:389
SmallStrainPlasticity::MakeB::makeB
PetscErrorCode makeB(const MatrixAdaptor &diffN, MatrixDouble &B)
Definition: SmallStrainPlasticity.cpp:172
SmallStrainPlasticity::OpAssembleLhs::CB
MatrixDouble CB
Definition: SmallStrainPlasticity.hpp:460
SmallStrainPlasticity::OpAssembleLhs::commonData
CommonData & commonData
Definition: SmallStrainPlasticity.hpp:458