v0.13.1
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 EntityHandle * getConn ()
 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 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)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is 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 call this function to loop over parent elements. This function calls finite element with is 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 call this function to loop over parent elements. This function calls finite element with is 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...
 
- 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

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 296 of file SmallStrainPlasticity.cpp.

299 :
301commonData(common_data) {
302 // sYmm = false;
303}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
constexpr auto field_name

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 
)
Examples
SmallStrainPlasticity.hpp.

Definition at line 305 of file SmallStrainPlasticity.cpp.

310 {
311 PetscFunctionBegin;
312 try {
313
314 int nb_dofs_row = row_data.getFieldData().size();
315 if(nb_dofs_row==0) PetscFunctionReturn(0);
316 int nb_dofs_col = col_data.getFieldData().size();
317 if(nb_dofs_col==0) PetscFunctionReturn(0);
318
319 K.resize(nb_dofs_row,nb_dofs_col,false);
320 K.clear();
321
322 int nb_gauss_pts = row_data.getN().size1();
323
324 if(commonData.bBar) {
325 double v = 0;
326 rowVolB.resize(6,nb_dofs_row,false);
327 rowVolB.clear();
328 colVolB.resize(6,nb_dofs_col,false);
329 colVolB.clear();
330 for(int gg = 0;gg<nb_gauss_pts;gg++) {
331 double val = getVolume()*getGaussPts()(3,gg);
332 v += val;
333 const MatrixAdaptor &diffN_row = row_data.getDiffN(gg,nb_dofs_row/3);
334 const MatrixAdaptor &diffN_col = col_data.getDiffN(gg,nb_dofs_col/3);
335 ierr = addVolumetricB(diffN_row,rowVolB,val); CHKERRQ(ierr);
336 ierr = addVolumetricB(diffN_col,colVolB,val); CHKERRQ(ierr);
337 }
338 rowVolB /= v;
339 colVolB /= v;
340 }
341
342 for(int gg = 0;gg<nb_gauss_pts;gg++) {
343
344 double val = getVolume()*getGaussPts()(3,gg);
345
346 const MatrixAdaptor &diffN_row = row_data.getDiffN(gg,nb_dofs_row/3);
347 const MatrixAdaptor &diffN_col = col_data.getDiffN(gg,nb_dofs_col/3);
348 ierr = makeB(diffN_row,rowB); CHKERRQ(ierr);
349 ierr = makeB(diffN_col,colB); CHKERRQ(ierr);
350 if(commonData.bBar) {
351 ierr = addVolumetricB(diffN_row,rowB,-1); CHKERRQ(ierr); // B-bar
352 ierr = addVolumetricB(diffN_col,colB,-1); CHKERRQ(ierr);
353 rowB += rowVolB;
354 colB += colVolB;
355 }
356
357 CB.resize(6,nb_dofs_col,false);
358 cblas_dgemm(
359 CblasRowMajor,CblasNoTrans,CblasNoTrans,
360 6,nb_dofs_col,6,
361 1.,
363 &colB(0,0),nb_dofs_col,
364 0,
365 &CB(0,0),nb_dofs_col
366 );
367 cblas_dgemm(
368 CblasRowMajor,CblasTrans,CblasNoTrans,
369 nb_dofs_row,nb_dofs_col,6,
370 val,
371 &rowB(0,0),nb_dofs_row,
372 &CB(0,0),nb_dofs_col,
373 1,
374 &K(0,0),nb_dofs_col
375 );
376 //noalias(CB) = prod(commonData.materialTangentOperator[gg],colB);
377 //noalias(K) += val*prod(trans(rowB),CB);
378
379 }
380
381 // cerr << K << endl;
382 // cerr << row_data.getIndices() << endl;
383 // cerr << col_data.getIndices() << endl;
384 // cerr << getFEMethod()->snes_A << endl;
385 // cerr << getFEMethod()->snes_B << endl;
386
387 int *row_indices_ptr,*col_indices_ptr;
388 row_indices_ptr = &row_data.getIndices()[0];
389 col_indices_ptr = &col_data.getIndices()[0];
390
392 getFEMethod()->snes_B,
393 nb_dofs_row,row_indices_ptr,
394 nb_dofs_col,col_indices_ptr,
395 &K(0,0),ADD_VALUES
396 ); CHKERRQ(ierr);
397
398 //is symmetric
399 if(row_side != col_side || row_type != col_type) {
400
401 transK.resize(nb_dofs_col,nb_dofs_row,false);
402 noalias(transK) = trans(K);
404 getFEMethod()->snes_B,
405 nb_dofs_col,col_indices_ptr,
406 nb_dofs_row,row_indices_ptr,
407 &transK(0,0),ADD_VALUES
408 ); CHKERRQ(ierr);
409
410 }
411
412 } catch (const std::exception& ex) {
413 ostringstream ss;
414 ss << "throw in method: " << ex.what() << endl;
415 SETERRQ(PETSC_COMM_SELF,1,ss.str().c_str());
416 }
417 PetscFunctionReturn(0);
418}
const double v
phase velocity of light in medium (cm/ns)
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
Definition: Types.hpp:132
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
vector< MatrixDouble > materialTangentOperator
PetscErrorCode makeB(const MatrixAdaptor &diffN, MatrixDouble &B)
PetscErrorCode addVolumetricB(const MatrixAdaptor &diffN, MatrixDouble &volB, double alpha)

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: