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

#include <users_modules/softmech/chemo_mech/src/ElecPhysOperators.hpp>

Inheritance diagram for ElectroPhysiology::OpAssembleSlowRhsV:
[legend]
Collaboration diagram for ElectroPhysiology::OpAssembleSlowRhsV:
[legend]

Public Types

typedef boost::function< double(const double, const double)> Feval_u
 
typedef boost::function< double(const double, const double)> FUVal
 
- 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 Member Functions

 OpAssembleSlowRhsV (std::string mass_field, boost::shared_ptr< PreviousData > &data_u, boost::shared_ptr< PreviousData > &data_v, Feval_u rhs_u)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 
 OpAssembleSlowRhsV (std::string mass_field, boost::shared_ptr< PreviousData > &common_datau, boost::shared_ptr< PreviousData > &common_datav, FUVal rhs_u)
 
MoFEMErrorCode doWork (int side, EntityType type, 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 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 Attributes

Feval_u rhsU
 
FTensor::Number< 0 > NX
 
FTensor::Number< 1 > NY
 
FTensor::Number< 2 > NZ
 
- 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...
 

Private Attributes

boost::shared_ptr< PreviousDatadataU
 
boost::shared_ptr< PreviousDatadataV
 
VectorDouble vecF
 
MatrixDouble mat
 
boost::shared_ptr< PreviousDatacommonDatau
 
boost::shared_ptr< PreviousDatacommonDatav
 
FUVal rhsU
 

Additional Inherited Members

- 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 285 of file ElecPhysOperators.hpp.

Member Typedef Documentation

◆ Feval_u

typedef boost::function<double(const double, const double)> ElectroPhysiology::OpAssembleSlowRhsV::Feval_u

Definition at line 288 of file ElecPhysOperators.hpp.

◆ FUVal

typedef boost::function<double(const double, const double)> ElectroPhysiology::OpAssembleSlowRhsV::FUVal

Definition at line 275 of file ElecPhysOperators2D.hpp.

Constructor & Destructor Documentation

◆ OpAssembleSlowRhsV() [1/2]

ElectroPhysiology::OpAssembleSlowRhsV::OpAssembleSlowRhsV ( std::string  mass_field,
boost::shared_ptr< PreviousData > &  data_u,
boost::shared_ptr< PreviousData > &  data_v,
Feval_u  rhs_u 
)
inline

Definition at line 289 of file ElecPhysOperators.hpp.

293 : OpVolEle(mass_field, OpVolEle::OPROW)
294 , dataU(data_u)
295 , dataV(data_v)
296 , rhsU(rhs_u) {}
VolEle::UserDataOperator OpVolEle
boost::shared_ptr< PreviousData > dataU
boost::shared_ptr< PreviousData > dataV
@ OPROW
operator doWork function is executed on FE rows

◆ OpAssembleSlowRhsV() [2/2]

ElectroPhysiology::OpAssembleSlowRhsV::OpAssembleSlowRhsV ( std::string  mass_field,
boost::shared_ptr< PreviousData > &  common_datau,
boost::shared_ptr< PreviousData > &  common_datav,
FUVal  rhs_u 
)
inline

Definition at line 276 of file ElecPhysOperators2D.hpp.

279 : OpFaceEle(mass_field, OpFaceEle::OPROW), commonDatau(common_datau),
280 commonDatav(common_datav), rhsU(rhs_u) {}
FaceEle::UserDataOperator OpFaceEle
boost::shared_ptr< PreviousData > commonDatav
boost::shared_ptr< PreviousData > commonDatau

Member Function Documentation

◆ doWork() [1/2]

MoFEMErrorCode ElectroPhysiology::OpAssembleSlowRhsV::doWork ( int  side,
EntityType  type,
EntData data 
)
inline

Definition at line 303 of file ElecPhysOperators.hpp.

303 {
305
306 const int nb_dofs = data.getIndices().size();
307 if (nb_dofs) {
308 if (nb_dofs != static_cast<int>(data.getN().size2()))
309 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
310 "wrong number of dofs");
311 vecF.resize(nb_dofs, false);
312 mat.resize(nb_dofs, nb_dofs, false);
313 vecF.clear();
314 mat.clear();
315 const int nb_integration_pts = getGaussPts().size2();
316 auto t_val_u = getFTensor0FromVec(dataU->mass_values);
317 auto t_val_v = getFTensor0FromVec(dataV->mass_values);
318
319
320 auto t_row_v_base = data.getFTensor0N();
321 auto t_w = getFTensor0IntegrationWeight();
322 const double vol = getMeasure();
323
324
325 for (int gg = 0; gg != nb_integration_pts; ++gg) {
326
327 double rhs = rhsU(t_val_u, t_val_v);
328 const double a = vol * t_w;
329
330 for (int rr = 0; rr != nb_dofs; ++rr) {
331 auto t_col_v_base = data.getFTensor0N(gg, 0);
332 vecF[rr] += a * rhs * t_row_v_base;
333 for (int cc = 0; cc != nb_dofs; ++cc) {
334 mat(rr, cc) += a * t_row_v_base * t_col_v_base;
335 ++t_col_v_base;
336 }
337 ++t_row_v_base;
338 }
339 ++t_val_u;
340 ++t_val_v;
341 ++t_w;
342
343 }
345 cholesky_solve(mat, vecF, ublas::lower());
346
347 CHKERR VecSetOption(getFEMethod()->ts_F, VEC_IGNORE_NEGATIVE_INDICES,
348 PETSC_TRUE);
349 CHKERR VecSetValues(getFEMethod()->ts_F, data, &*vecF.begin(),
350 ADD_VALUES);
351 }
353 }
constexpr double a
void cholesky_solve(const TRIA &L, VEC &x, ublas::lower)
solve system L L^T x = b inplace
Definition: cholesky.hpp:221
size_t cholesky_decompose(const MATRIX &A, TRIA &L)
decompose the symmetric positive definit matrix A into product L L^T.
Definition: cholesky.hpp:52
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
#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
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of dofs on entity.
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

◆ doWork() [2/2]

MoFEMErrorCode ElectroPhysiology::OpAssembleSlowRhsV::doWork ( int  side,
EntityType  type,
EntData data 
)
inline

Definition at line 282 of file ElecPhysOperators2D.hpp.

282 {
284 // cerr << "In OpAssembleSlowRhsV...." << endl;
285 const int nb_dofs = data.getIndices().size();
286 if (nb_dofs) {
287 // cerr << "In SlowRhsV..." << endl;
288 if (nb_dofs != static_cast<int>(data.getN().size2()))
289 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
290 "wrong number of dofs");
291 vecF.resize(nb_dofs, false);
292 mat.resize(nb_dofs, nb_dofs, false);
293 vecF.clear();
294 mat.clear();
295 const int nb_integration_pts = getGaussPts().size2();
296 auto t_u_value = getFTensor0FromVec(commonDatau->mass_values);
297 auto t_v_value = getFTensor0FromVec(commonDatav->mass_values);
298 auto t_row_v_base = data.getFTensor0N();
299 auto t_w = getFTensor0IntegrationWeight();
300 const double vol = getMeasure();
301
302 const double ct = getFEMethod()->ts_t - 0.01;
303 auto t_coords = getFTensor1CoordsAtGaussPts();
304 for (int gg = 0; gg != nb_integration_pts; ++gg) {
305 const double a = vol * t_w;
306
307 // double u_dot = exactDot(t_coords(NX), t_coords(NY), ct);
308 // double u_lap = exactLap(t_coords(NX), t_coords(NY), ct);
309
310 // double f = u_dot - u_lap;
311 double f = t_u_value * (1.0 - t_u_value);
312 for (int rr = 0; rr != nb_dofs; ++rr) {
313 double rhs = rhsU(t_u_value, t_v_value);
314 auto t_col_v_base = data.getFTensor0N(gg, 0);
315 vecF[rr] += a * rhs * t_row_v_base;
316 // vecF[rr] += a * f * t_row_v_base;
317 for (int cc = 0; cc != nb_dofs; ++cc) {
318 mat(rr, cc) += a * t_row_v_base * t_col_v_base;
319 ++t_col_v_base;
320 }
321 ++t_row_v_base;
322 }
323 ++t_u_value;
324 ++t_v_value;
325 ++t_w;
326 ++t_coords;
327 }
329 cholesky_solve(mat, vecF, ublas::lower());
330
331 CHKERR VecSetOption(getFEMethod()->ts_F, VEC_IGNORE_NEGATIVE_INDICES,
332 PETSC_TRUE);
333 CHKERR VecSetValues(getFEMethod()->ts_F, data, &*vecF.begin(),
334 ADD_VALUES);
335 }
337 }
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
PetscReal ts_t
time

Member Data Documentation

◆ commonDatau

boost::shared_ptr<PreviousData> ElectroPhysiology::OpAssembleSlowRhsV::commonDatau
private

Definition at line 340 of file ElecPhysOperators2D.hpp.

◆ commonDatav

boost::shared_ptr<PreviousData> ElectroPhysiology::OpAssembleSlowRhsV::commonDatav
private

Definition at line 341 of file ElecPhysOperators2D.hpp.

◆ dataU

boost::shared_ptr<PreviousData> ElectroPhysiology::OpAssembleSlowRhsV::dataU
private

Definition at line 356 of file ElecPhysOperators.hpp.

◆ dataV

boost::shared_ptr<PreviousData> ElectroPhysiology::OpAssembleSlowRhsV::dataV
private

Definition at line 357 of file ElecPhysOperators.hpp.

◆ mat

MatrixDouble ElectroPhysiology::OpAssembleSlowRhsV::mat
private

Definition at line 359 of file ElecPhysOperators.hpp.

◆ NX

FTensor::Number< 0 > ElectroPhysiology::OpAssembleSlowRhsV::NX

Definition at line 299 of file ElecPhysOperators.hpp.

◆ NY

FTensor::Number< 1 > ElectroPhysiology::OpAssembleSlowRhsV::NY

Definition at line 300 of file ElecPhysOperators.hpp.

◆ NZ

FTensor::Number< 2 > ElectroPhysiology::OpAssembleSlowRhsV::NZ

Definition at line 301 of file ElecPhysOperators.hpp.

◆ rhsU [1/2]

Feval_u ElectroPhysiology::OpAssembleSlowRhsV::rhsU

Definition at line 298 of file ElecPhysOperators.hpp.

◆ rhsU [2/2]

FUVal ElectroPhysiology::OpAssembleSlowRhsV::rhsU
private

Definition at line 344 of file ElecPhysOperators2D.hpp.

◆ vecF

VectorDouble ElectroPhysiology::OpAssembleSlowRhsV::vecF
private

Definition at line 358 of file ElecPhysOperators.hpp.


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