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

Assemble tangent matrix. More...

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

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

Public Member Functions

 OpAssembleTangent (const 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 Attributes

CommonDatacommonData
 
MatrixDouble kTangent
 
VectorDouble auxVec
 
MatrixDouble auxMat
 

Detailed Description

Assemble tangent matrix.

Examples
minimal_surface_area.cpp.

Definition at line 328 of file MinimalSurfaceElement.hpp.

Constructor & Destructor Documentation

◆ OpAssembleTangent()

MinimalSurfaceEquation::MinimalSurfaceElement::OpAssembleTangent::OpAssembleTangent ( const string  field_name,
CommonData common_data 
)
inline

Definition at line 331 of file MinimalSurfaceElement.hpp.

334  commonData(common_data) {
335  sYmm = false; // matrix is not symmetric
336  }

Member Function Documentation

◆ doWork()

PetscErrorCode MinimalSurfaceEquation::MinimalSurfaceElement::OpAssembleTangent::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 340 of file MinimalSurfaceElement.hpp.

343  {
345 
346  int row_nb_dofs = row_data.getIndices().size();
347  if (row_nb_dofs == 0) {
349  }
350  int col_nb_dofs = col_data.getIndices().size();
351  if (col_nb_dofs == 0) {
353  }
354  kTangent.resize(row_nb_dofs, col_nb_dofs, false);
355  kTangent.clear();
356  auxVec.resize(col_nb_dofs, false);
357  auxMat.resize(col_nb_dofs, 2, false);
358  int nb_gauss_pts = row_data.getN().size1();
359  for (int gg = 0; gg != nb_gauss_pts; gg++) {
360  // cerr << "gg " << gg << endl;
361  double val = getGaussPts()(2, gg) * getArea();
362  MatrixAdaptor row_diff_n = row_data.getDiffN(gg);
363  MatrixAdaptor col_diff_n = col_data.getDiffN(gg);
364  double an = commonData.aN[gg];
365  ublas::matrix_row<MatrixDouble> grad_u(commonData.gradU, gg);
366  ublas::matrix_row<MatrixDouble> an_pow3_by_grad_u(
368  // cerr << "grad_u " << grad_u << endl;
369  // cerr << "an_pow3_by_grad_u " << an_pow3_by_grad_u << endl;
370  noalias(auxVec) = prod(col_diff_n, an_pow3_by_grad_u);
371  // cerr << "auxVec " << auxVec << endl;
372  noalias(auxMat) = outer_prod(auxVec, grad_u);
373  // cerr << "auxMat " << auxMat << endl;
374  // cerr << row_diff_n << endl;
375  // cerr << col_diff_n << endl;
376  // cerr << prod(row_diff_n,trans(val*(an*col_diff_n+auxMat))) << endl;
377  noalias(kTangent) +=
378  val * prod(row_diff_n, trans(an * col_diff_n - auxMat));
379  // cerr << "kTangent " << kTangent << endl;
380  }
381  CHKERR MatSetValues(getFEMethod()->snes_B, row_nb_dofs,
382  &*row_data.getIndices().data().begin(), col_nb_dofs,
383  &*col_data.getIndices().data().begin(),
384  &*kTangent.data().begin(), ADD_VALUES);
386  }

Member Data Documentation

◆ auxMat

MatrixDouble MinimalSurfaceEquation::MinimalSurfaceElement::OpAssembleTangent::auxMat

Definition at line 339 of file MinimalSurfaceElement.hpp.

◆ auxVec

VectorDouble MinimalSurfaceEquation::MinimalSurfaceElement::OpAssembleTangent::auxVec

Definition at line 338 of file MinimalSurfaceElement.hpp.

◆ commonData

CommonData& MinimalSurfaceEquation::MinimalSurfaceElement::OpAssembleTangent::commonData

Definition at line 330 of file MinimalSurfaceElement.hpp.

◆ kTangent

MatrixDouble MinimalSurfaceEquation::MinimalSurfaceElement::OpAssembleTangent::kTangent

Definition at line 337 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:460
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:1644
MinimalSurfaceEquation::MinimalSurfaceElement::CommonData::gradU
MatrixDouble gradU
nb_gauss_pts x 2 gradients at integration pts
Definition: MinimalSurfaceElement.hpp:76
MinimalSurfaceEquation::MinimalSurfaceElement::CommonData::aNpow3byGradU
MatrixDouble aNpow3byGradU
nb_gauss_pts x 2,
Definition: MinimalSurfaceElement.hpp:81
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:548
MinimalSurfaceEquation::MinimalSurfaceElement::OpAssembleTangent::kTangent
MatrixDouble kTangent
Definition: MinimalSurfaceElement.hpp:337
MinimalSurfaceEquation::MinimalSurfaceElement::OpAssembleTangent::commonData
CommonData & commonData
Definition: MinimalSurfaceElement.hpp:330
MinimalSurfaceEquation::MinimalSurfaceElement::OpAssembleTangent::auxVec
VectorDouble auxVec
Definition: MinimalSurfaceElement.hpp:338
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MinimalSurfaceEquation::MinimalSurfaceElement::OpAssembleTangent::auxMat
MatrixDouble auxMat
Definition: MinimalSurfaceElement.hpp:339
MinimalSurfaceEquation::MinimalSurfaceElement::CommonData::aN
VectorDouble aN
size of nb_gauss_pts,
Definition: MinimalSurfaceElement.hpp:79
MoFEM::Types::MatrixAdaptor
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
Definition: Types.hpp:132
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359