v0.15.5
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Attributes | List of all members
OpAdJointObjective Struct Reference
Inheritance diagram for OpAdJointObjective:
[legend]
Collaboration diagram for OpAdJointObjective:
[legend]

Public Types

using OP = ForcesAndSourcesCore::UserDataOperator
 

Public Member Functions

 OpAdJointObjective (boost::shared_ptr< ObjectiveFunctionData > python_ptr, boost::shared_ptr< HookeOps::CommonData > comm_ptr, boost::shared_ptr< MatrixDouble > jac_ptr, boost::shared_ptr< MatrixDouble > diff_jac, boost::shared_ptr< VectorDouble > cof_vals, boost::shared_ptr< MatrixDouble > d_grad_ptr, boost::shared_ptr< MatrixDouble > d_u_ptr, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< double > glob_objective_ptr, boost::shared_ptr< double > glob_objective_grad_ptr)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Compute objective function contributions at element level.
 

Private Attributes

boost::shared_ptr< ObjectiveFunctionDatapythonPtr
 
boost::shared_ptr< HookeOps::CommonDatacommPtr
 
boost::shared_ptr< MatrixDoublejacPtr
 
boost::shared_ptr< MatrixDoublediffJacPtr
 
boost::shared_ptr< VectorDoublecofVals
 
boost::shared_ptr< MatrixDoubledGradPtr
 
boost::shared_ptr< MatrixDoubledUPtr
 
boost::shared_ptr< MatrixDoubleuPtr
 
boost::shared_ptr< doubleglobObjectivePtr
 
boost::shared_ptr< doubleglobObjectiveGradPtr
 

Detailed Description

Examples
mofem/tutorials/vec-7_shape_optimisation/adjoint.cpp.

Definition at line 1573 of file adjoint.cpp.

Member Typedef Documentation

◆ OP

using OpAdJointObjective::OP = ForcesAndSourcesCore::UserDataOperator

Definition at line 1574 of file adjoint.cpp.

Constructor & Destructor Documentation

◆ OpAdJointObjective()

OpAdJointObjective::OpAdJointObjective ( boost::shared_ptr< ObjectiveFunctionData python_ptr,
boost::shared_ptr< HookeOps::CommonData comm_ptr,
boost::shared_ptr< MatrixDouble jac_ptr,
boost::shared_ptr< MatrixDouble diff_jac,
boost::shared_ptr< VectorDouble cof_vals,
boost::shared_ptr< MatrixDouble d_grad_ptr,
boost::shared_ptr< MatrixDouble d_u_ptr,
boost::shared_ptr< MatrixDouble u_ptr,
boost::shared_ptr< double glob_objective_ptr,
boost::shared_ptr< double glob_objective_grad_ptr 
)
inline

Definition at line 1576 of file adjoint.cpp.

1586 : OP(NOSPACE, OP::OPSPACE), pythonPtr(python_ptr), commPtr(comm_ptr),
1587 jacPtr(jac_ptr), diffJacPtr(diff_jac), cofVals(cof_vals),
1588 dGradPtr(d_grad_ptr), dUPtr(d_u_ptr), uPtr(u_ptr),
1589 globObjectivePtr(glob_objective_ptr),
1590 globObjectiveGradPtr(glob_objective_grad_ptr) {}
@ NOSPACE
Definition definitions.h:83
boost::shared_ptr< MatrixDouble > dGradPtr
Definition adjoint.cpp:1733
boost::shared_ptr< double > globObjectiveGradPtr
Definition adjoint.cpp:1738
boost::shared_ptr< ObjectiveFunctionData > pythonPtr
Definition adjoint.cpp:1728
boost::shared_ptr< double > globObjectivePtr
Definition adjoint.cpp:1737
ForcesAndSourcesCore::UserDataOperator OP
Definition adjoint.cpp:1574
boost::shared_ptr< HookeOps::CommonData > commPtr
Definition adjoint.cpp:1729
boost::shared_ptr< VectorDouble > cofVals
Definition adjoint.cpp:1732
boost::shared_ptr< MatrixDouble > uPtr
Definition adjoint.cpp:1735
boost::shared_ptr< MatrixDouble > diffJacPtr
Definition adjoint.cpp:1731
boost::shared_ptr< MatrixDouble > jacPtr
Definition adjoint.cpp:1730
boost::shared_ptr< MatrixDouble > dUPtr
Definition adjoint.cpp:1734

Member Function Documentation

◆ doWork()

MoFEMErrorCode OpAdJointObjective::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)
inline

Compute objective function contributions at element level.

Evaluates Python objective function with current displacement and stress state, and accumulates global objective value and gradients.

Examples
mofem/tutorials/vec-7_shape_optimisation/adjoint.cpp.

Definition at line 1598 of file adjoint.cpp.

1599 {
1601
1602 // Define tensor indices for calculations
1607
1610
1611 constexpr auto symm_size = (SPACE_DIM * (SPACE_DIM + 1)) / 2; // size of symmetric tensor in Voigt notation
1612
1613 auto t_diff_symm = diff_symmetrize(FTensor::Number<SPACE_DIM>()); // fourth order tensor for symmetrization to convert from gradient to strain
1614
1615 auto nb_gauss_pts = getGaussPts().size2(); // number of gauss points in the element
1616 auto objective_ptr = boost::make_shared<VectorDouble>(nb_gauss_pts); // objective function values at gauss points
1617 auto objective_dstress =
1618 boost::make_shared<MatrixDouble>(symm_size, nb_gauss_pts); // objective function gradient w.r.t. stress at gauss points
1619 auto objective_dstrain =
1620 boost::make_shared<MatrixDouble>(symm_size, nb_gauss_pts); // objective function gradient w.r.t. strain at gauss points
1621 auto objective_du = boost::make_shared<MatrixDouble>(SPACE_DIM, nb_gauss_pts);
1622// We have dJ/dsigma and dJ/depsilon at gauss points, we need to convert it to dJ/du
1623
1624 auto evaluate_python = [&]() {
1626 auto &coords = OP::getCoordsAtGaussPts();
1627 CHKERR pythonPtr->evalInteriorObjectiveFunction(
1628 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1629 objective_ptr);
1630 CHKERR pythonPtr->evalInteriorObjectiveGradientStress(
1631 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1632 objective_dstress);
1633 CHKERR pythonPtr->evalInteriorObjectiveGradientStrain(
1634 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1635 objective_dstrain);
1636 CHKERR pythonPtr->evalInteriorObjectiveGradientU(
1637 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1638 objective_du);
1639
1640 auto t_grad_u =
1641 getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(commPtr->matGradPtr));
1642 auto t_D =
1643 getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*(commPtr->matDPtr));
1644 auto t_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(jacPtr));
1645 auto t_diff_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(diffJacPtr));
1646 auto t_cof = getFTensor0FromVec(*(cofVals));
1647 auto t_d_grad = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(dGradPtr));
1648
1649 auto t_obj = getFTensor0FromVec(*objective_ptr);
1650 auto t_obj_dstress =
1651 getFTensor2SymmetricFromMat<SPACE_DIM>(*objective_dstress);
1652 auto t_obj_dstrain =
1653 getFTensor2SymmetricFromMat<SPACE_DIM>(*objective_dstrain);
1654 auto t_obj_du = getFTensor1FromMat<SPACE_DIM>(*objective_du);
1655 auto t_d_u = getFTensor1FromMat<SPACE_DIM>(*dUPtr);
1656
1657 auto vol = OP::getMeasure();
1658 auto t_w = getFTensor0IntegrationWeight();
1659 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1660
1661 auto t_det = determinantTensor(t_jac);
1663 CHKERR invertTensor(t_jac, t_det, t_inv_jac);
1664
1666 t_diff_inv_jac(i, j) =
1667 -(t_inv_jac(i, I) * t_diff_jac(I, J)) * t_inv_jac(J, j);
1669 t_diff_grad(i, j) = t_grad_u(i, k) * t_diff_inv_jac(k, j);
1670
1672 t_d_strain(i, j) = t_diff_symm(i, j, k, l) * (
1673
1674 t_d_grad(k, l)
1675
1676 +
1677
1678 t_diff_grad(k, l)
1679
1680 );
1681
1682 auto alpha = t_w * vol;
1683
1684 (*globObjectivePtr) += alpha * t_obj;
1685 (*globObjectiveGradPtr) +=
1686 alpha *
1687 (
1688
1689 t_obj_dstress(i, j) * (t_D(i, j, k, l) * t_d_strain(k, l))
1690
1691 +
1692
1693 t_obj_dstrain(i, j) * t_d_strain(i, j)
1694
1695 +
1696
1697 t_obj_du(i) * t_d_u(i)
1698
1699 +
1700
1701 t_obj * t_cof
1702
1703 );
1704
1705 ++t_w;
1706 ++t_jac;
1707 ++t_diff_jac;
1708 ++t_cof;
1709
1710 ++t_obj;
1711 ++t_obj_dstress;
1712 ++t_obj_dstrain;
1713 ++t_obj_du;
1714
1715 ++t_grad_u;
1716 ++t_d_grad;
1717 ++t_d_u;
1718 }
1720 };
1721
1722 CHKERR evaluate_python();
1723
1725 }
#define FTENSOR_INDEX(DIM, I)
auto diff_symmetrize(FTensor::Number< DIM >)
Definition adjoint.cpp:1344
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()
#define CHKERR
Inline error check.
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
static MoFEMErrorCode invertTensor(FTensor::Tensor2< T1, DIM, DIM > &t, T2 &det, FTensor::Tensor2< T3, DIM, DIM > &inv_t)
static auto determinantTensor(FTensor::Tensor2< T, DIM, DIM > &t)
Calculate the determinant of a tensor of rank DIM.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.

Member Data Documentation

◆ cofVals

boost::shared_ptr<VectorDouble> OpAdJointObjective::cofVals
private

◆ commPtr

boost::shared_ptr<HookeOps::CommonData> OpAdJointObjective::commPtr
private

◆ dGradPtr

boost::shared_ptr<MatrixDouble> OpAdJointObjective::dGradPtr
private

◆ diffJacPtr

boost::shared_ptr<MatrixDouble> OpAdJointObjective::diffJacPtr
private

◆ dUPtr

boost::shared_ptr<MatrixDouble> OpAdJointObjective::dUPtr
private

◆ globObjectiveGradPtr

boost::shared_ptr<double> OpAdJointObjective::globObjectiveGradPtr
private

◆ globObjectivePtr

boost::shared_ptr<double> OpAdJointObjective::globObjectivePtr
private

◆ jacPtr

boost::shared_ptr<MatrixDouble> OpAdJointObjective::jacPtr
private

◆ pythonPtr

boost::shared_ptr<ObjectiveFunctionData> OpAdJointObjective::pythonPtr
private

◆ uPtr

boost::shared_ptr<MatrixDouble> OpAdJointObjective::uPtr
private

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