12#ifndef __ADOLC_HUHU_HPP__
13#define __ADOLC_HUHU_HPP__
37 template <AssemblyType A, IntegrationType I,
typename DomainEleOp>
40 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
42 boost::shared_ptr<PhysicalEquations> physical_equations_ptr,
43 Sev sev = Sev::noisy) {
45 constexpr int DIM = (MODEL_TYPE ==
MODEL_3D) ? 3 : 2;
46 auto op_this = getPipThis<DomainEleOp>(m_field, pip, fe_name,
field_name,
47 physical_equations_ptr, sev);
48 auto m_grad_grad = boost::make_shared<MatrixDouble>();
49 op_this->getOpPtrVector().push_back(
50 new OpCalculateVectorFieldHessian<DIM, DIM>(
field_name, m_grad_grad));
51 op_this->getOpPtrVector().push_back(
52 physical_equations_ptr->createOp(physical_equations_ptr,
true,
false));
53 auto m_k = physical_equations_ptr->adolcDataPtr->getCommonDataPtr(
"k");
54 op_this->getOpPtrVector().push_back(
59 template <AssemblyType A, IntegrationType I,
typename DomainEleOp>
62 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
64 boost::shared_ptr<PhysicalEquations> physical_equations_ptr,
65 Sev sev = Sev::noisy) {
67 constexpr int DIM = (MODEL_TYPE ==
MODEL_3D) ? 3 : 2;
68 auto op_this = getPipThis<DomainEleOp>(m_field, pip, fe_name,
field_name,
69 physical_equations_ptr, sev);
70 auto m_grad_grad = boost::make_shared<MatrixDouble>();
71 op_this->getOpPtrVector().push_back(
72 new OpCalculateVectorFieldHessian<DIM, DIM>(
field_name, m_grad_grad));
73 op_this->getOpPtrVector().push_back(
74 physical_equations_ptr->createOp(physical_equations_ptr,
true,
true));
75 auto m_k = physical_equations_ptr->adolcDataPtr->getCommonDataPtr(
"k");
77 physical_equations_ptr->adolcDataPtr->getCommonDataPtr(
"k_dF");
78 op_this->getOpPtrVector().push_back(
80 op_this->getOpPtrVector().push_back(
87 template <
typename DomainEleOp>
90 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
92 boost::shared_ptr<PhysicalEquations> physical_equations_ptr,
95 auto ¶m_vec_by_range = physical_equations_ptr->paramVecByRange;
97 auto r = boost::make_shared<Range>();
98 for (
auto &p : param_vec_by_range) {
101 MOFEM_LOG(
"WORLD", Sev::inform) <<
"HuHu number of entities " << r->size();
103 constexpr int DIM = (MODEL_TYPE ==
MODEL_3D) ? 3 : 2;
107 auto space = field_structure->getSpace();
109 using DomainEle =
typename parent_of<DomainEleOp>::type;
110 auto op_this =
new OpLoopThis<DomainEle>(m_field, fe_name, r, sev);
111 pip.push_back(op_this);
112 auto this_fe_ptr = op_this->getThisFEPtr();
113 auto &this_pip = op_this->getOpPtrVector();
114 this_fe_ptr->getRuleHook = [](
int order_row,
int order_col,
116 return 2 * (order_data - 1);
119 auto base_mass = boost::make_shared<MatrixDouble>();
120 auto data_l2 = boost::make_shared<EntitiesFieldData>(MBENTITYSET);
121 auto jac_ptr = boost::make_shared<MatrixDouble>();
122 auto det_ptr = boost::make_shared<VectorDouble>();
123 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
126 this_pip.push_back(
new OpCalculateHOJac<DIM>(jac_ptr));
128 this_pip.push_back(
new OpInvertMatrix<DIM>(jac_ptr, det_ptr, inv_jac_ptr));
131 new OpBaseDerivativesMass<1>(base_mass, data_l2, base,
L2));
133 this_pip.push_back(
new OpBaseDerivativesNext<1>(
134 BaseDerivatives::SecondDerivative, base_mass, data_l2, base, space));
140 new OpSetHOInvJacToScalarBases<DIM, 1>(space, inv_jac_ptr));
143 new OpSetHOInvJacToScalarBases<DIM, 2>(space, inv_jac_ptr));
147 <<
"Unsupported space: " << space <<
" for field: " <<
field_name;
152 physical_equations_ptr->adolcDataPtr->getCommonDataPtr(
"grad");
154 new OpCalculateVectorFieldGradient<DIM, DIM>(
field_name, m_grad));
161boost::shared_ptr<PhysicalEquations>
163 boost::shared_ptr<ADolCData> adolc_data_ptr,
int tag);
166boost::shared_ptr<PhysicalEquations>
168 boost::shared_ptr<ADolCData> adolc_data_ptr,
int tag);
170template <
int FIELD_DIM,
int SPACE_DIM, AssemblyType A,
typename OpBase>
172 :
public FormsIntegrators<OpBase>::template Assembly<A>
::OpBase {
174 using OP =
typename FormsIntegrators<OpBase>::template Assembly<A>::OpBase;
177 boost::shared_ptr<MatrixDouble> mat_vals,
178 boost::shared_ptr<MatrixDouble> mat_K,
179 boost::shared_ptr<Range> ents_ptr =
nullptr)
184 boost::shared_ptr<MatrixDouble>
matK;
185 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &data);
188template <
int FIELD_DIM,
int SPACE_DIM, AssemblyType A,
typename OpBase>
190 :
public FormsIntegrators<OpBase>::template Assembly<A>
::OpBase {
192 using OP =
typename FormsIntegrators<OpBase>::template Assembly<A>::OpBase;
194 boost::shared_ptr<Range> ents_ptr =
nullptr)
198 boost::shared_ptr<MatrixDouble>
matK;
199 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
200 EntitiesFieldData::EntData &col_data);
203template <
int FIELD_DIM,
int SPACE_DIM, AssemblyType A,
typename OpBase>
205 :
public FormsIntegrators<OpBase>::template Assembly<A>
::OpBase {
206 using OP =
typename FormsIntegrators<OpBase>::template Assembly<A>::OpBase;
208 boost::shared_ptr<MatrixDouble> mat_vals,
209 boost::shared_ptr<MatrixDouble> mat_diff_K,
210 boost::shared_ptr<Range> ents_ptr =
nullptr)
212 matDiffK(mat_diff_K) {
219 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
220 EntitiesFieldData::EntData &col_data);
223template <
int FIELD_DIM,
int SPACE_DIM, AssemblyType A,
typename OpBase>
225 EntitiesFieldData::EntData &row_data) {
233 const double vol = OP::getMeasure();
235 auto t_w = OP::getFTensor0IntegrationWeight();
237 auto t_row_grad = row_data.getFTensor2DiffN2<
SPACE_DIM>();
239 auto t_val_grad = getFTensor3FromMat<FIELD_DIM, SPACE_DIM, SPACE_DIM>(
242 auto t_K = getFTensor1FromMat<1>(*(matK));
244 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
246 for (
int gg = 0; gg != OP::nbIntegrationPts; gg++) {
248 const double alpha = t_w * vol * t_K(0);
250 auto t_nf = getFTensor1FromPtr<FIELD_DIM>(OP::locF.data().data());
253 for (; rr != OP::nbRows /
FIELD_DIM; rr++) {
255 t_nf(
i) += alpha * (t_row_grad(
J, K) * t_val_grad(
i,
J, K));
259 for (; rr < OP::nbRowBaseFunctions; ++rr)
270template <
int FIELD_DIM,
int SPACE_DIM, AssemblyType A,
typename OpBase>
272 EntitiesFieldData::EntData &row_data,
273 EntitiesFieldData::EntData &col_data) {
284 const double vol = OP::getMeasure();
286 auto t_w = OP::getFTensor0IntegrationWeight();
288 auto t_row_grad = row_data.getFTensor2DiffN2<
SPACE_DIM>();
290 auto t_K = getFTensor1FromMat<1>(*(matK));
292 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
294 for (
int gg = 0; gg != OP::nbIntegrationPts; gg++) {
296 const double alpha = t_w * vol * t_K(0);
299 for (; rr != OP::nbRows /
FIELD_DIM; rr++) {
301 auto t_mat = getFTensor2FromPtr<FIELD_DIM, FIELD_DIM, FIELD_DIM>(
303 auto t_col_grad = col_data.getFTensor2DiffN2<
SPACE_DIM>(gg, 0);
306 for (
int cc = 0; cc != OP::nbCols /
FIELD_DIM; cc++) {
308 alpha * (t_row_grad(
J, K) * t_col_grad(
J, K)) *
t_kd(
i,
j);
316 for (; rr < OP::nbRowBaseFunctions; ++rr)
326template <
int FIELD_DIM,
int SPACE_DIM, AssemblyType A,
typename OpBase>
328 EntitiesFieldData::EntData &row_data,
329 EntitiesFieldData::EntData &col_data) {
339 const double vol = OP::getMeasure();
341 auto t_w = OP::getFTensor0IntegrationWeight();
343 auto t_row_grad = row_data.getFTensor2DiffN2<
SPACE_DIM>();
345 auto t_val_grad = getFTensor3FromMat<FIELD_DIM, SPACE_DIM, SPACE_DIM>(
348 auto t_diff_K = getFTensor2FromMat<FIELD_DIM, SPACE_DIM>(
351 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
353 for (
int gg = 0; gg != OP::nbIntegrationPts; gg++) {
355 const double alpha = t_w * vol;
358 for (; rr != OP::nbRows /
FIELD_DIM; rr++) {
360 auto t_mat = getFTensor2FromPtr<FIELD_DIM, FIELD_DIM, FIELD_DIM>(
362 auto t_col_grad = col_data.getFTensor1DiffN<
SPACE_DIM>(gg, 0);
363 for (
int bb = 0; bb != OP::nbCols /
FIELD_DIM; ++bb) {
364 t_mat(
i,
j) += alpha * (t_row_grad(
J, K) * t_val_grad(
i,
J, K)) *
365 (t_diff_K(
j, M) * t_col_grad(M));
375 for (; rr < OP::nbRowBaseFunctions; ++rr)
#define FTENSOR_INDEX(DIM, I)
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
@ L2
field with C-1 continuity
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_OPERATION_UNSUCCESSFUL
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
#define MOFEM_LOG(channel, severity)
Log.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'J', DIM1 > J
FTensor::Index< 'j', 3 > j
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< HUHU, MODEL_2D_PLANE_STRAIN >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< HUHU, MODEL_3D >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
constexpr IntegrationType I
constexpr auto field_name
boost::shared_ptr< MatrixDouble > matVals
OpLhsHuGrad(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, boost::shared_ptr< MatrixDouble > mat_diff_K, boost::shared_ptr< Range > ents_ptr=nullptr)
boost::shared_ptr< MatrixDouble > matDiffK
boost::shared_ptr< MatrixDouble > matK
OpLhsHuHu(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_K, boost::shared_ptr< Range > ents_ptr=nullptr)
OpMaterialFactory()=delete
static MoFEMErrorCode opRhsFactory(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string fe_name, std::string field_name, boost::shared_ptr< PhysicalEquations > physical_equations_ptr, Sev sev=Sev::noisy)
static auto getPipThis(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string fe_name, std::string field_name, boost::shared_ptr< PhysicalEquations > physical_equations_ptr, Sev sev)
static MoFEMErrorCode opLhsFactory(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string fe_name, std::string field_name, boost::shared_ptr< PhysicalEquations > physical_equations_ptr, Sev sev=Sev::noisy)
boost::shared_ptr< MatrixDouble > matVals
OpRhsHu(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, boost::shared_ptr< MatrixDouble > mat_K, boost::shared_ptr< Range > ents_ptr=nullptr)
boost::shared_ptr< MatrixDouble > matK
Deprecated interface functions.
FieldApproximationBase getApproxBase() const
Get approximation basis type.