v0.15.0
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Attributes | List of all members
OpAdJointObjective Struct Reference

Operator for computing objective function and gradients in topology optimization. More...

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 > 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< MatrixDoubleuPtr
 
boost::shared_ptr< doubleglobObjectivePtr
 
boost::shared_ptr< doubleglobObjectiveGradPtr
 

Detailed Description

Operator for computing objective function and gradients in topology optimization.

This operator interfaces with Python-defined objective functions to:

  1. Evaluate objective function f(u, x) at current design point
  2. Compute gradients ∂f/∂u and ∂f/∂x for adjoint sensitivity analysis

The objective function typically depends on:

This enables flexible objective function definition (compliance, stress constraints, volume constraints, etc.) through Python scripting.

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

Definition at line 1794 of file adjoint.cpp.

Member Typedef Documentation

◆ OP

using OpAdJointObjective::OP = ForcesAndSourcesCore::UserDataOperator

Definition at line 1795 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 u_ptr,
boost::shared_ptr< double glob_objective_ptr,
boost::shared_ptr< double glob_objective_grad_ptr 
)
inline

Definition at line 1797 of file adjoint.cpp.

1806 : OP(NOSPACE, OP::OPSPACE), pythonPtr(python_ptr), commPtr(comm_ptr),
1807 jacPtr(jac_ptr), diffJacPtr(diff_jac), cofVals(cof_vals),
1808 dGradPtr(d_grad_ptr), uPtr(u_ptr), globObjectivePtr(glob_objective_ptr),
1809 globObjectiveGradPtr(glob_objective_grad_ptr) {}
@ NOSPACE
Definition definitions.h:83
boost::shared_ptr< MatrixDouble > dGradPtr
Definition adjoint.cpp:1940
boost::shared_ptr< double > globObjectiveGradPtr
Definition adjoint.cpp:1944
boost::shared_ptr< ObjectiveFunctionData > pythonPtr
Definition adjoint.cpp:1935
boost::shared_ptr< double > globObjectivePtr
Definition adjoint.cpp:1943
ForcesAndSourcesCore::UserDataOperator OP
Definition adjoint.cpp:1795
boost::shared_ptr< HookeOps::CommonData > commPtr
Definition adjoint.cpp:1936
boost::shared_ptr< VectorDouble > cofVals
Definition adjoint.cpp:1939
boost::shared_ptr< MatrixDouble > uPtr
Definition adjoint.cpp:1941
boost::shared_ptr< MatrixDouble > diffJacPtr
Definition adjoint.cpp:1938
boost::shared_ptr< MatrixDouble > jacPtr
Definition adjoint.cpp:1937

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/adjoint.cpp.

Definition at line 1817 of file adjoint.cpp.

1818 {
1820
1821 // Define tensor indices for calculations
1826
1829
1830 constexpr auto symm_size = (SPACE_DIM * (SPACE_DIM + 1)) / 2;
1831
1832 auto t_diff_symm = diff_symmetrize(FTensor::Number<SPACE_DIM>());
1833
1834 auto nb_gauss_pts = getGaussPts().size2();
1835 auto objective_ptr = boost::make_shared<VectorDouble>(nb_gauss_pts);
1836 auto objective_dstress =
1837 boost::make_shared<MatrixDouble>(symm_size, nb_gauss_pts);
1838 auto objective_dstrain =
1839 boost::make_shared<MatrixDouble>(symm_size, nb_gauss_pts);
1840 auto obj_grad = boost::make_shared<MatrixDouble>(SPACE_DIM, nb_gauss_pts);
1841
1842 auto evaluate_python = [&]() {
1844 auto &coords = OP::getCoordsAtGaussPts();
1845 CHKERR pythonPtr->evalObjectiveFunction(
1846 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1847 objective_ptr);
1848 CHKERR pythonPtr->evalObjectiveGradientStress(
1849 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1850 objective_dstress);
1851 CHKERR pythonPtr->evalObjectiveGradientStrain(
1852 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1853 objective_dstrain);
1854
1855 auto t_grad_u =
1856 getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(commPtr->matGradPtr));
1857 auto t_D =
1858 getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*(commPtr->matDPtr));
1859 auto t_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(jacPtr));
1860 auto t_diff_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(diffJacPtr));
1861 auto t_cof = getFTensor0FromVec(*(cofVals));
1862 auto t_d_grad = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(dGradPtr));
1863
1864 auto t_obj = getFTensor0FromVec(*objective_ptr);
1865 auto t_obj_dstress =
1866 getFTensor2SymmetricFromMat<SPACE_DIM>(*objective_dstress);
1867 auto t_obj_dstrain =
1868 getFTensor2SymmetricFromMat<SPACE_DIM>(*objective_dstrain);
1869
1870 auto vol = OP::getMeasure();
1871 auto t_w = getFTensor0IntegrationWeight();
1872 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1873
1874 auto t_det = determinantTensor(t_jac);
1876 CHKERR invertTensor(t_jac, t_det, t_inv_jac);
1877
1879 t_diff_inv_jac(i, j) =
1880 -(t_inv_jac(i, I) * t_diff_jac(I, J)) * t_inv_jac(J, j);
1882 t_diff_grad(i, j) = t_grad_u(i, k) * t_diff_inv_jac(k, j);
1883
1885 t_d_strain(i, j) = t_diff_symm(i, j, k, l) * (
1886
1887 t_d_grad(k, l)
1888
1889 +
1890
1891 t_diff_grad(k, l)
1892
1893 );
1894
1895 auto alpha = t_w * vol;
1896
1897 (*globObjectivePtr) += alpha * t_obj;
1898 (*globObjectiveGradPtr) +=
1899 alpha *
1900 (
1901
1902 t_obj_dstress(i, j) * (t_D(i, j, k, l) * t_d_strain(k, l))
1903
1904 +
1905
1906 t_obj_dstrain(i, j) * t_d_strain(i, j)
1907
1908 +
1909
1910 t_obj * t_cof
1911
1912 );
1913
1914 ++t_w;
1915 ++t_jac;
1916 ++t_diff_jac;
1917 ++t_cof;
1918
1919 ++t_obj;
1920 ++t_obj_dstress;
1921 ++t_obj_dstrain;
1922
1923 ++t_grad_u;
1924 ++t_d_grad;
1925 }
1927 };
1928
1929 CHKERR evaluate_python();
1930
1932 }
#define FTENSOR_INDEX(DIM, I)
auto diff_symmetrize(FTensor::Number< DIM >)
Definition adjoint.cpp:1649
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
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 1939 of file adjoint.cpp.

◆ commPtr

boost::shared_ptr<HookeOps::CommonData> OpAdJointObjective::commPtr
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 1936 of file adjoint.cpp.

◆ dGradPtr

boost::shared_ptr<MatrixDouble> OpAdJointObjective::dGradPtr
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 1940 of file adjoint.cpp.

◆ diffJacPtr

boost::shared_ptr<MatrixDouble> OpAdJointObjective::diffJacPtr
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 1938 of file adjoint.cpp.

◆ globObjectiveGradPtr

boost::shared_ptr<double> OpAdJointObjective::globObjectiveGradPtr
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 1944 of file adjoint.cpp.

◆ globObjectivePtr

boost::shared_ptr<double> OpAdJointObjective::globObjectivePtr
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 1943 of file adjoint.cpp.

◆ jacPtr

boost::shared_ptr<MatrixDouble> OpAdJointObjective::jacPtr
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 1937 of file adjoint.cpp.

◆ pythonPtr

boost::shared_ptr<ObjectiveFunctionData> OpAdJointObjective::pythonPtr
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 1935 of file adjoint.cpp.

◆ uPtr

boost::shared_ptr<MatrixDouble> OpAdJointObjective::uPtr
private
Examples
mofem/tutorials/vec-7/adjoint.cpp.

Definition at line 1941 of file adjoint.cpp.


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