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

#include <users_modules/basic_finite_elements/src/ConvectiveMassElement.hpp>

Inheritance diagram for ConvectiveMassElement::OpMassLhs_dM_dv:
[legend]
Collaboration diagram for ConvectiveMassElement::OpMassLhs_dM_dv:
[legend]

Public Member Functions

 OpMassLhs_dM_dv (const std::string vel_field, const std::string field_name, BlockData &data, CommonData &common_data, Range *forcesonlyonentities_ptr=NULL)
 
virtual MoFEMErrorCode getJac (EntitiesFieldData::EntData &col_data, int gg)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 

Public Attributes

BlockDatadAta
 
CommonDatacommonData
 
Range forcesOnlyOnEntities
 
MatrixDouble k
 
MatrixDouble jac
 

Detailed Description

Definition at line 222 of file ConvectiveMassElement.hpp.

Constructor & Destructor Documentation

◆ OpMassLhs_dM_dv()

ConvectiveMassElement::OpMassLhs_dM_dv::OpMassLhs_dM_dv ( const std::string  vel_field,
const std::string  field_name,
BlockData data,
CommonData common_data,
Range forcesonlyonentities_ptr = NULL 
)

Definition at line 439 of file ConvectiveMassElement.cpp.

442 : VolumeElementForcesAndSourcesCore::UserDataOperator(
443 vel_field, field_name,
444 ForcesAndSourcesCore::UserDataOperator::OPROWCOL),
445 dAta(data), commonData(common_data) {
446 sYmm = false;
447 if (forcesonlyonentities_ptr != NULL) {
448 forcesOnlyOnEntities = *forcesonlyonentities_ptr;
449 }
450}
constexpr auto field_name

Member Function Documentation

◆ doWork()

MoFEMErrorCode ConvectiveMassElement::OpMassLhs_dM_dv::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)

Definition at line 521 of file ConvectiveMassElement.cpp.

524 {
526
527 if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
528 dAta.tEts.end()) {
530 }
531
532 int nb_row = row_data.getIndices().size();
533 int nb_col = col_data.getIndices().size();
534 if (nb_row == 0)
536 if (nb_col == 0)
538
539 auto base = row_data.getFTensor0N();
540 int nb_base_functions = row_data.getN().size2();
541
542 {
543
544 k.resize(nb_row, nb_col);
545 k.clear();
546 jac.resize(3, nb_col);
547
548 for (unsigned int gg = 0; gg < row_data.getN().size1(); gg++) {
549
550 try {
551 CHKERR getJac(col_data, gg);
552 } catch (const std::exception &ex) {
553 std::ostringstream ss;
554 ss << "throw in method: " << ex.what() << std::endl;
555 SETERRQ(PETSC_COMM_SELF, 1, ss.str().c_str());
556 }
557
560
561 {
562 int dd1 = 0;
563 // integrate element stiffness matrix
564 for (; dd1 < nb_row / 3; dd1++) {
566 &jac(0, 0), &jac(0, 1), &jac(0, 2), &jac(1, 0), &jac(1, 1),
567 &jac(1, 2), &jac(2, 0), &jac(2, 1), &jac(2, 2), 3);
568 for (int dd2 = 0; dd2 < nb_col / 3; dd2++) {
570 &k(3 * dd1 + 0, 3 * dd2 + 0), &k(3 * dd1 + 0, 3 * dd2 + 1),
571 &k(3 * dd1 + 0, 3 * dd2 + 2), &k(3 * dd1 + 1, 3 * dd2 + 0),
572 &k(3 * dd1 + 1, 3 * dd2 + 1), &k(3 * dd1 + 1, 3 * dd2 + 2),
573 &k(3 * dd1 + 2, 3 * dd2 + 0), &k(3 * dd1 + 2, 3 * dd2 + 1),
574 &k(3 * dd1 + 2, 3 * dd2 + 2));
575 t_k(i, j) += base * t_jac(i, j);
576 ++t_jac;
577 }
578 ++base;
579 // for(int rr1 = 0;rr1<3;rr1++) {
580 // for(int dd2 = 0;dd2<nb_col;dd2++) {
581 // k(3*dd1+rr1,dd2) += row_data.getN()(gg,dd1)*jac(rr1,dd2);
582 // }
583 // }
584 }
585 for (; dd1 != nb_base_functions; dd1++) {
586 ++base;
587 }
588 }
589 }
590
591 if (!forcesOnlyOnEntities.empty()) {
592 VectorInt indices = row_data.getIndices();
593 VectorDofs &dofs = row_data.getFieldDofs();
594 VectorDofs::iterator dit = dofs.begin();
595 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
596 if (forcesOnlyOnEntities.find((*dit)->getEnt()) ==
597 forcesOnlyOnEntities.end()) {
598 indices[ii] = -1;
599 }
600 }
601 CHKERR MatSetValues(getFEMethod()->ts_B, nb_row, &indices[0], nb_col,
602 &col_data.getIndices()[0], &k(0, 0), ADD_VALUES);
603 } else {
604 CHKERR MatSetValues(getFEMethod()->ts_B, nb_row,
605 &row_data.getIndices()[0], nb_col,
606 &col_data.getIndices()[0], &k(0, 0), ADD_VALUES);
607 }
608 }
610}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
UBlasVector< int > VectorInt
Definition: Types.hpp:67
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
virtual MoFEMErrorCode getJac(EntitiesFieldData::EntData &col_data, int gg)
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 VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
const VectorInt & getIndices() const
Get global indices of dofs on entity.

◆ getJac()

MoFEMErrorCode ConvectiveMassElement::OpMassLhs_dM_dv::getJac ( EntitiesFieldData::EntData col_data,
int  gg 
)
virtual

Reimplemented in ConvectiveMassElement::OpMassLhs_dM_dx, ConvectiveMassElement::OpMassLhs_dM_dX, ConvectiveMassElement::OpVelocityLhs_dV_dv, ConvectiveMassElement::OpVelocityLhs_dV_dx, ConvectiveMassElement::OpVelocityLhs_dV_dX, ConvectiveMassElement::OpEshelbyDynamicMaterialMomentumLhs_dv, ConvectiveMassElement::OpEshelbyDynamicMaterialMomentumLhs_dx, and ConvectiveMassElement::OpEshelbyDynamicMaterialMomentumLhs_dX.

Definition at line 452 of file ConvectiveMassElement.cpp.

453 {
455 int nb_col = col_data.getIndices().size();
456 jac.clear();
457 if (!nb_col)
462 FTensor::Tensor2<double *, 3, 3> t_jac(&jac(0, 0), &jac(0, 1), &jac(0, 2),
463 &jac(1, 0), &jac(1, 1), &jac(1, 2),
464 &jac(2, 0), &jac(2, 1), &jac(2, 2), 3);
466 &commonData.jacMass[gg](0, 0), &commonData.jacMass[gg](0, 1),
467 &commonData.jacMass[gg](0, 2), &commonData.jacMass[gg](1, 0),
468 &commonData.jacMass[gg](1, 1), &commonData.jacMass[gg](1, 2),
469 &commonData.jacMass[gg](2, 0), &commonData.jacMass[gg](2, 1),
470 &commonData.jacMass[gg](2, 2));
471 double *base_ptr = const_cast<double *>(&col_data.getN(gg)[0]);
472 FTensor::Tensor0<double *> base(base_ptr, 1);
473 if (commonData.dataAtGaussPts["DOT_" + commonData.meshPositions].size() ==
474 0) {
475 for (int dd = 0; dd < nb_col / 3; dd++) {
476 t_jac(i, j) += t_mass1(i, j) * base * getFEMethod()->ts_a;
477 ++base;
478 ++t_jac;
479 }
480 } else {
481 const int s = 3 + 9;
483 // T* d000, T* d001, T* d002,
484 // T* d010, T* d011, T* d012,
485 // T* d020, T* d021, T* d022,
486 // T* d100, T* d101, T* d102,
487 // T* d110, T* d111, T* d112,
488 // T* d120, T* d121, T* d122,
489 // T* d200, T* d201, T* d202,
490 // T* d210, T* d211, T* d212,
491 // T* d220, T* d221, T* d222,
492 &commonData.jacMass[gg](0, s + 0), &commonData.jacMass[gg](0, s + 1),
493 &commonData.jacMass[gg](0, s + 2), &commonData.jacMass[gg](0, s + 3),
494 &commonData.jacMass[gg](0, s + 4), &commonData.jacMass[gg](0, s + 5),
495 &commonData.jacMass[gg](0, s + 6), &commonData.jacMass[gg](0, s + 7),
496 &commonData.jacMass[gg](0, s + 8), &commonData.jacMass[gg](1, s + 0),
497 &commonData.jacMass[gg](1, s + 1), &commonData.jacMass[gg](1, s + 2),
498 &commonData.jacMass[gg](1, s + 3), &commonData.jacMass[gg](1, s + 4),
499 &commonData.jacMass[gg](1, s + 5), &commonData.jacMass[gg](1, s + 6),
500 &commonData.jacMass[gg](1, s + 7), &commonData.jacMass[gg](1, s + 8),
501 &commonData.jacMass[gg](2, s + 0), &commonData.jacMass[gg](2, s + 1),
502 &commonData.jacMass[gg](2, s + 2), &commonData.jacMass[gg](2, s + 3),
503 &commonData.jacMass[gg](2, s + 4), &commonData.jacMass[gg](2, s + 5),
504 &commonData.jacMass[gg](2, s + 6), &commonData.jacMass[gg](2, s + 7),
505 &commonData.jacMass[gg](2, s + 8));
506
507 double *diff_ptr =
508 const_cast<double *>(&(col_data.getDiffN(gg, nb_col / 3)(0, 0)));
509 FTensor::Tensor1<double *, 3> diff(diff_ptr, &diff_ptr[1], &diff_ptr[2], 3);
510 for (int dd = 0; dd < nb_col / 3; dd++) {
511 t_jac(i, j) += t_mass1(i, j) * base * getFEMethod()->ts_a;
512 t_jac(i, j) += t_mass3(i, j, k) * diff(k);
513 ++base;
514 ++diff;
515 ++t_jac;
516 }
517 }
519}
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
std::map< std::string, std::vector< VectorDouble > > dataAtGaussPts
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions

Member Data Documentation

◆ commonData

CommonData& ConvectiveMassElement::OpMassLhs_dM_dv::commonData

Definition at line 227 of file ConvectiveMassElement.hpp.

◆ dAta

BlockData& ConvectiveMassElement::OpMassLhs_dM_dv::dAta

Definition at line 226 of file ConvectiveMassElement.hpp.

◆ forcesOnlyOnEntities

Range ConvectiveMassElement::OpMassLhs_dM_dv::forcesOnlyOnEntities

Definition at line 228 of file ConvectiveMassElement.hpp.

◆ jac

MatrixDouble ConvectiveMassElement::OpMassLhs_dM_dv::jac

Definition at line 234 of file ConvectiveMassElement.hpp.

◆ k

MatrixDouble ConvectiveMassElement::OpMassLhs_dM_dv::k

Definition at line 234 of file ConvectiveMassElement.hpp.


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