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

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

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

Public Member Functions

 OpAssembleRhs (string field_name, CommonData &common_data)
 
PetscErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &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 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
 
VectorDouble nF
 
MatrixDouble B
 
MatrixDouble volB
 
- 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

Definition at line 438 of file SmallStrainPlasticity.hpp.

Constructor & Destructor Documentation

◆ OpAssembleRhs()

SmallStrainPlasticity::OpAssembleRhs::OpAssembleRhs ( string  field_name,
CommonData common_data 
)

Definition at line 219 of file SmallStrainPlasticity.cpp.

222 :
223VolumeElementForcesAndSourcesCore::UserDataOperator(field_name,UserDataOperator::OPROW),
224commonData(common_data) {
225}
constexpr auto field_name

Member Function Documentation

◆ doWork()

PetscErrorCode SmallStrainPlasticity::OpAssembleRhs::doWork ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData &  data 
)
Examples
SmallStrainPlasticity.hpp.

Definition at line 227 of file SmallStrainPlasticity.cpp.

229 {
230 PetscFunctionBegin;
231
232 try {
233
234 int nb_dofs = data.getFieldData().size();
235 if(nb_dofs==0) PetscFunctionReturn(0);
236 int nb_gauss_pts = data.getN().size1();
237 nF.resize(nb_dofs,false);
238 nF.clear();
239
240 if(commonData.bBar) {
241 double v = 0;
242 volB.resize(6,nb_dofs,false);
243 volB.clear();
244 for(int gg = 0;gg<nb_gauss_pts;gg++) {
245 double val = getVolume()*getGaussPts()(3,gg);
246 v += val;
247 const MatrixAdaptor &diffN = data.getDiffN(gg,nb_dofs/3);
248 ierr = addVolumetricB(diffN,volB,val); CHKERRQ(ierr);
249 }
250 volB /= v;
251 }
252
253 for(int gg = 0;gg<nb_gauss_pts;gg++) {
254
255 double val = getVolume()*getGaussPts()(3,gg);
256 const MatrixAdaptor &diffN = data.getDiffN(gg,nb_dofs/3);
257 ierr = makeB(diffN,B); CHKERRQ(ierr);
258 if(commonData.bBar) {
259 ierr = addVolumetricB(diffN,B,-1); CHKERRQ(ierr);
260 B += volB;
261 }
262 // cerr << gg << endl;
263 // cerr << diffN << endl;
264 // cerr << B << endl;
265 // cerr << commonData.sTress[gg] << endl;
266 // cerr << prod(trans(B),commonData.sTress[gg]) << endl;
267
268 noalias(nF) += val*prod(trans(B),commonData.sTress[gg]);
269 // cerr << "nF " << nF << endl;
270
271 }
272
273 // cerr << side << " " << type << " " << data.getN().size1() << " "<< data.getN().size2() << endl;
274
275 int *indices_ptr;
276 indices_ptr = &data.getIndices()[0];
277
279 getFEMethod()->snes_f,
280 nb_dofs,
281 indices_ptr,
282 &nF[0],
283 ADD_VALUES
284 ); CHKERRQ(ierr);
285
286
287 } catch (const std::exception& ex) {
288 ostringstream ss;
289 ss << "throw in method: " << ex.what() << endl;
290 SETERRQ(PETSC_COMM_SELF,1,ss.str().c_str());
291 }
292
293 PetscFunctionReturn(0);
294}
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
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
PetscErrorCode makeB(const MatrixAdaptor &diffN, MatrixDouble &B)
PetscErrorCode addVolumetricB(const MatrixAdaptor &diffN, MatrixDouble &volB, double alpha)

Member Data Documentation

◆ B

MatrixDouble SmallStrainPlasticity::OpAssembleRhs::B
Examples
SmallStrainPlasticity.hpp.

Definition at line 442 of file SmallStrainPlasticity.hpp.

◆ commonData

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

Definition at line 440 of file SmallStrainPlasticity.hpp.

◆ nF

VectorDouble SmallStrainPlasticity::OpAssembleRhs::nF
Examples
SmallStrainPlasticity.hpp.

Definition at line 441 of file SmallStrainPlasticity.hpp.

◆ volB

MatrixDouble SmallStrainPlasticity::OpAssembleRhs::volB
Examples
SmallStrainPlasticity.hpp.

Definition at line 442 of file SmallStrainPlasticity.hpp.


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