v0.15.4
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
OpAdJointTestOp< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp > Struct Template Reference

Adjoint test operator for elastic problems. More...

Inheritance diagram for OpAdJointTestOp< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >:
[legend]
Collaboration diagram for OpAdJointTestOp< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >:
[legend]

Public Types

using OP = DomainBaseOp
 

Public Member Functions

 OpAdJointTestOp (const std::string field_name, boost::shared_ptr< HookeOps::CommonData > comm_ptr)
 

Protected Member Functions

MoFEMErrorCode iNtegrate (EntitiesFieldData::EntData &data)
 Integration implementation for adjoint test operator.
 

Protected Attributes

boost::shared_ptr< HookeOps::CommonDatacommPtr
 

Detailed Description

template<int SPACE_DIM>
struct OpAdJointTestOp< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >

Adjoint test operator for elastic problems.

This operator implements the adjoint (transpose) of the elastic stiffness operator. It's used in the adjoint equation K^T * λ = ∂f/∂u to compute sensitivities.

The operator computes: ∫_Ω ∇v : σ dΩ where v are test functions and σ is the Cauchy stress tensor.

This is the transpose of the standard elastic operator ∫_Ω ∇u : C : ∇v dΩ

Template Parameters
SPACE_DIMSpatial dimension (2 or 3)

Definition at line 1553 of file adjoint.cpp.

Member Typedef Documentation

◆ OP

template<int SPACE_DIM>
using OpAdJointTestOp< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >::OP = DomainBaseOp

Definition at line 1556 of file adjoint.cpp.

Constructor & Destructor Documentation

◆ OpAdJointTestOp()

template<int SPACE_DIM>
OpAdJointTestOp< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >::OpAdJointTestOp ( const std::string  field_name,
boost::shared_ptr< HookeOps::CommonData comm_ptr 
)
inline

Definition at line 1558 of file adjoint.cpp.

1560 : OP(field_name, field_name, DomainEleOp::OPROW), commPtr(comm_ptr) {}
constexpr auto field_name

Member Function Documentation

◆ iNtegrate()

template<int SPACE_DIM>
MoFEMErrorCode OpAdJointTestOp< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >::iNtegrate ( EntitiesFieldData::EntData row_data)
protected

Integration implementation for adjoint test operator.

Computes the integral: ∫_Ω ∇v : σ dΩ where v are test functions and σ is Cauchy stress tensor

This forms the right-hand side of the adjoint equation K^T * λ = ∂f/∂u

Parameters
row_dataTest function data for current element
Returns
MoFEMErrorCode Success or error code

Definition at line 1580 of file adjoint.cpp.

1581 {
1583 // Define tensor indices for Einstein notation
1590
1591 // Get element geometric information
1592 const double vol = OP::getMeasure();
1593 // Get integration weights for Gauss quadrature
1594 auto t_w = OP::getFTensor0IntegrationWeight();
1595 // Get test function gradients ∇v
1596 auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
1597 // Get Cauchy stress tensor σ from forward solution
1598 auto t_cauchy_stress =
1599 getFTensor2SymmetricFromMat<SPACE_DIM>(*(commPtr->getMatCauchyStress()));
1600
1601 // Loop over integration points for numerical integration
1602 for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
1603 // take into account Jacobian
1604 const double alpha = t_w * vol;
1605
1606 // get rhs vector
1607 auto t_nf = OP::template getNf<SPACE_DIM>();
1608 // loop over rows base functions
1609 int rr = 0;
1610 for (; rr != OP::nbRows / SPACE_DIM; rr++) {
1611 // Variation due to change in geometry (diff base)
1612 t_nf(j) += alpha * t_row_grad(i) * t_cauchy_stress(i, j);
1613
1614 ++t_row_grad;
1615 ++t_nf;
1616 }
1617 for (; rr < OP::nbRowBaseFunctions; ++rr) {
1618 ++t_row_grad;
1619 }
1620
1621 ++t_cauchy_stress;
1622 ++t_w; // move to another integration weight
1623 }
1625}
#define FTENSOR_INDEX(DIM, I)
constexpr int SPACE_DIM
[Define dimension]
Definition adjoint.cpp:28
constexpr IntegrationType I
Use Gauss quadrature for integration.
Definition adjoint.cpp:33
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.

Member Data Documentation

◆ commPtr

template<int SPACE_DIM>
boost::shared_ptr<HookeOps::CommonData> OpAdJointTestOp< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >::commPtr
protected

Definition at line 1563 of file adjoint.cpp.


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