v0.14.0
Public Member Functions | Public Attributes | List of all members
MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs Struct Reference

Integrate vector on lhs,. More...

#include <users_modules/minimal_surface_equation/src/MinimalSurfaceElement.hpp>

Inheritance diagram for MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs:
[legend]
Collaboration diagram for MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs:
[legend]

Public Member Functions

 OpAssmebleBcLhs (const string field_name, Mat m_a)
 
PetscErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 

Public Attributes

Mat mA
 
MatrixDouble K
 
MatrixDouble transK
 

Detailed Description

Integrate vector on lhs,.

\[ \mathbf{A} = \int_L \mathbf{N}^\textrm{T}\mathbf{N} \textrm{d}L \]

Examples
minimal_surface_area.cpp.

Definition at line 139 of file MinimalSurfaceElement.hpp.

Constructor & Destructor Documentation

◆ OpAssmebleBcLhs()

MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs::OpAssmebleBcLhs ( const string  field_name,
Mat  m_a 
)
inline

Member Function Documentation

◆ doWork()

PetscErrorCode MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
DataForcesAndSourcesCore::EntData row_data,
DataForcesAndSourcesCore::EntData col_data 
)
inline

Definition at line 147 of file MinimalSurfaceElement.hpp.

150  {
152 
153  int row_nb_dofs = row_data.getIndices().size();
154  if (row_nb_dofs == 0) {
156  }
157  int col_nb_dofs = col_data.getIndices().size();
158  if (col_nb_dofs == 0) {
160  }
161  K.resize(row_nb_dofs, col_nb_dofs, false);
162  K.clear();
163  int nb_gauss_pts = row_data.getN().size1();
164  for (int gg = 0; gg != nb_gauss_pts; gg++) {
165  double val = getLength() * getGaussPts()(1, gg);
166  VectorAdaptor row_n = row_data.getN(gg);
167  VectorAdaptor col_n = col_data.getN(gg);
168  noalias(K) += val * outer_prod(row_n, col_n);
169  }
170  CHKERR MatSetValues(mA, row_nb_dofs,
171  &*row_data.getIndices().data().begin(), col_nb_dofs,
172  &*col_data.getIndices().data().begin(),
173  &*K.data().begin(), ADD_VALUES);
174  // Matrix is symmetric, assemble off diagonal part
175  if (row_type != col_type || row_side != col_side) {
176  transK.resize(col_nb_dofs, row_nb_dofs, false);
177  noalias(transK) = trans(K);
178  CHKERR MatSetValues(mA, col_nb_dofs,
179  &*col_data.getIndices().data().begin(), row_nb_dofs,
180  &*row_data.getIndices().data().begin(),
181  &*transK.data().begin(), ADD_VALUES);
182  }
183 
185  }

Member Data Documentation

◆ K

MatrixDouble MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs::K

Definition at line 146 of file MinimalSurfaceElement.hpp.

◆ mA

Mat MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs::mA

Definition at line 141 of file MinimalSurfaceElement.hpp.

◆ transK

MatrixDouble MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs::transK

Definition at line 146 of file MinimalSurfaceElement.hpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
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
MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs::transK
MatrixDouble transK
Definition: MinimalSurfaceElement.hpp:146
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs::K
MatrixDouble K
Definition: MinimalSurfaceElement.hpp:146
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MinimalSurfaceEquation::MinimalSurfaceElement::OpAssmebleBcLhs::mA
Mat mA
Definition: MinimalSurfaceElement.hpp:141
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346