v0.15.5
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
OpAdJointGradTimesSymTensor< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp > Struct Template Reference
Inheritance diagram for OpAdJointGradTimesSymTensor< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >:
[legend]
Collaboration diagram for OpAdJointGradTimesSymTensor< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >:
[legend]

Public Types

using OP = DomainBaseOp
 

Public Member Functions

 OpAdJointGradTimesSymTensor (const std::string field_name, boost::shared_ptr< HookeOps::CommonData > comm_ptr, boost::shared_ptr< MatrixDouble > jac, boost::shared_ptr< MatrixDouble > diff_jac, boost::shared_ptr< VectorDouble > cof_vals)
 

Protected Member Functions

MoFEMErrorCode iNtegrate (EntitiesFieldData::EntData &data)
 

Protected Attributes

boost::shared_ptr< HookeOps::CommonDatacommPtr
 
boost::shared_ptr< MatrixDoublejac
 
boost::shared_ptr< MatrixDoublediffJac
 
boost::shared_ptr< VectorDoublecofVals
 

Detailed Description

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

Definition at line 1323 of file adjoint.cpp.

Member Typedef Documentation

◆ OP

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

Definition at line 1326 of file adjoint.cpp.

Constructor & Destructor Documentation

◆ OpAdJointGradTimesSymTensor()

template<int SPACE_DIM>
OpAdJointGradTimesSymTensor< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >::OpAdJointGradTimesSymTensor ( const std::string  field_name,
boost::shared_ptr< HookeOps::CommonData comm_ptr,
boost::shared_ptr< MatrixDouble jac,
boost::shared_ptr< MatrixDouble diff_jac,
boost::shared_ptr< VectorDouble cof_vals 
)
inline

Member Function Documentation

◆ iNtegrate()

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

Definition at line 1382 of file adjoint.cpp.

1383 {
1391
1392 auto t_diff_symm = diff_symmetrize(FTensor::Number<SPACE_DIM>());
1393
1394 // get element volume
1395 const double vol = OP::getMeasure();
1396 // get integration weights
1397 auto t_w = OP::getFTensor0IntegrationWeight();
1398 // get Jacobian values
1399 auto t_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(jac));
1400 // get diff Jacobian values
1401 auto t_diff_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(diffJac));
1402 // get cofactor values
1403 auto t_cof = getFTensor0FromVec(*(cofVals));
1404 // get base function gradient on rows
1405 auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
1406 // get fradient of the field
1407 auto t_grad_u =
1408 getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(commPtr->matGradPtr));
1409 // get field gradient values
1410 auto t_cauchy_stress =
1411 getFTensor2SymmetricFromMat<SPACE_DIM>(*(commPtr->getMatCauchyStress()));
1412 // material stiffness tensor
1413 auto t_D =
1414 getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*(commPtr->matDPtr));
1415 // loop over integration points
1416 for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
1417 // take into account Jacobian
1418 const double alpha = t_w * vol;
1419
1420 auto t_det = determinantTensor(t_jac);
1422 CHKERR invertTensor(t_jac, t_det, t_inv_jac);
1423
1424 // Calculate the variation of the gradient due to geometry change
1426 t_diff_inv_jac(i, j) =
1427 -(t_inv_jac(i, I) * t_diff_jac(I, J)) * t_inv_jac(J, j);
1429 t_diff_grad(i, j) = t_grad_u(i, k) * t_diff_inv_jac(k, j);
1430
1431 // Calculate the variation of the strain tensor
1433 t_diff_strain(i, j) = t_diff_symm(i, j, k, l) * t_diff_grad(k, l);
1434
1435 // Calculate the variation of the stress tensor
1437 t_diff_stress(i, j) = t_D(i, j, k, l) * t_diff_strain(k, l);
1438
1439 // get rhs vector
1440 auto t_nf = OP::template getNf<SPACE_DIM>();
1441 // loop over rows base functions
1442 int rr = 0;
1443 for (; rr != OP::nbRows / SPACE_DIM; rr++) {
1444
1446 t_diff_row_grad(k) = t_row_grad(j) * t_diff_inv_jac(j, k);
1447
1448 // Variation due to change in geometry (diff base)
1449 t_nf(j) += alpha * t_diff_row_grad(i) * t_cauchy_stress(i, j);
1450
1451 // Variation due to change in domain (cofactor)
1452 t_nf(j) += (alpha * t_cof) * t_row_grad(i) * t_cauchy_stress(i, j);
1453
1454 // Variation due to change in stress (diff stress)
1455 t_nf(j) += alpha * t_row_grad(i) * t_diff_stress(i, j);
1456
1457 ++t_row_grad;
1458 ++t_nf;
1459 }
1460 for (; rr < OP::nbRowBaseFunctions; ++rr) {
1461 ++t_row_grad;
1462 }
1463
1464 ++t_grad_u;
1465 ++t_cauchy_stress;
1466 ++t_jac;
1467 ++t_diff_jac;
1468 ++t_cof;
1469 ++t_w; // move to another integration weight
1470 }
1472}
#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

template<int SPACE_DIM>
boost::shared_ptr<VectorDouble> OpAdJointGradTimesSymTensor< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >::cofVals
protected

Definition at line 1340 of file adjoint.cpp.

◆ commPtr

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

Definition at line 1337 of file adjoint.cpp.

◆ diffJac

template<int SPACE_DIM>
boost::shared_ptr<MatrixDouble> OpAdJointGradTimesSymTensor< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >::diffJac
protected

Definition at line 1339 of file adjoint.cpp.

◆ jac

template<int SPACE_DIM>
boost::shared_ptr<MatrixDouble> OpAdJointGradTimesSymTensor< SPACE_DIM, IntegrationType::GAUSS, DomainBaseOp >::jac
protected

Definition at line 1338 of file adjoint.cpp.


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