v0.15.5
Loading...
Searching...
No Matches
EshelbianADOL-C.cpp
Go to the documentation of this file.
1/**
2 * \file EshelbianPlasticity.cpp
3 * \brief Implementation of automatic differentiation
4 */
5
6#include <MoFEM.hpp>
7using namespace MoFEM;
8
9#include <MatrixFunction.hpp>
10
12
13#include <EshelbianAux.hpp>
14
16
17#include <boost/math/constants/constants.hpp>
18
19constexpr double third = boost::math::constants::third<double>();
20
21namespace EshelbianPlasticity {
22
24
25 OpSpatialPhysical(const std::string &field_name,
26 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
27 const double alpha_u)
28 : OpAssembleVolume(field_name, data_ptr, OPROW), alphaU(alpha_u) {}
29
31
32private:
33 const double alphaU;
34};
35
38
40 auto t_L = symm_L_tensor();
41
42 int nb_dofs = data.getIndices().size();
43 int nb_integration_pts = data.getN().size1();
44 auto v = getVolume();
46 auto t_approx_P_adjoint_log_du =
47 getFTensor1FromMat<size_symm>(dataAtPts->adjointPdUAtPts);
48 auto t_P = getFTensor2FromMat<3, 3>(dataAtPts->PAtPts);
49 auto t_dot_log_u =
50 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchDotTensorAtPts);
51 auto t_diff_u =
52 getFTensor4DdgFromMat<3, 3, 1>(dataAtPts->diffStretchTensorAtPts);
53
54 FTensor::Index<'i', 3> i;
55 FTensor::Index<'j', 3> j;
56 FTensor::Index<'k', 3> k;
57 FTensor::Index<'l', 3> l;
58 auto get_ftensor2 = [](auto &v) {
60 &v[0], &v[1], &v[2], &v[3], &v[4], &v[5]);
61 };
62
63 int nb_base_functions = data.getN().size2();
64 auto t_row_base_fun = data.getFTensor0N();
65 for (int gg = 0; gg != nb_integration_pts; ++gg) {
66 double a = v * t_w;
67 auto t_nf = get_ftensor2(nF);
68
70 t_Ldiff_u(i, j, L) = t_diff_u(i, j, k, l) * t_L(k, l, L);
71
73 t_viscous_P(i, j) =
74 alphaU *
75 t_dot_log_u(i,
76 j); // That is cheat, should be split on axial and deviator
77
79 t_residual(L) =
80 a * (t_approx_P_adjoint_log_du(L) - t_Ldiff_u(i, j, L) * t_P(i, j) -
81 t_L(i, j, L) * t_viscous_P(i, j));
82
83 int bb = 0;
84 for (; bb != nb_dofs / 6; ++bb) {
85 t_nf(L) -= t_row_base_fun * t_residual(L);
86 ++t_nf;
87 ++t_row_base_fun;
88 }
89 for (; bb != nb_base_functions; ++bb)
90 ++t_row_base_fun;
91
92 ++t_w;
93 ++t_approx_P_adjoint_log_du;
94 ++t_P;
95 ++t_dot_log_u;
96 ++t_diff_u;
97 }
99}
100
102 const std::string &field_name,
103 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha_u) {
104 return new OpSpatialPhysical(field_name, data_ptr, alpha_u);
105}
106
108 const std::string &field_name,
109 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
110 boost::shared_ptr<ExternalStrainVec> external_strain_vec_ptr,
111 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv) {
112 return nullptr;
113}
114
116 OpSpatialPhysical_du_du(std::string row_field, std::string col_field,
117 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
118 const double alpha)
119 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false),
120 alphaU(alpha) {
121 sYmm = false;
122 }
123 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
124 MoFEMErrorCode integratePiola(EntData &row_data, EntData &col_data);
125
126private:
127 const double alphaU;
129
130};
131
133 EntData &col_data) {
134 return integratePiola(row_data, col_data);
135}
136
138 EntData &col_data) {
140
143 auto t_L = symm_L_tensor();
144 auto t_diff = diff_tensor();
145
146 int nb_integration_pts = row_data.getN().size1();
147 int row_nb_dofs = row_data.getIndices().size();
148 int col_nb_dofs = col_data.getIndices().size();
149
150 auto get_ftensor2 = [](MatrixDouble &m, const int r, const int c) {
152 size_symm>(
153
154 &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 0, c + 3),
155 &m(r + 0, c + 4), &m(r + 0, c + 5),
156
157 &m(r + 1, c + 0), &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 1, c + 3),
158 &m(r + 1, c + 4), &m(r + 1, c + 5),
159
160 &m(r + 2, c + 0), &m(r + 2, c + 1), &m(r + 2, c + 2), &m(r + 2, c + 3),
161 &m(r + 2, c + 4), &m(r + 2, c + 5),
162
163 &m(r + 3, c + 0), &m(r + 3, c + 1), &m(r + 3, c + 2), &m(r + 3, c + 3),
164 &m(r + 3, c + 4), &m(r + 3, c + 5),
165
166 &m(r + 4, c + 0), &m(r + 4, c + 1), &m(r + 4, c + 2), &m(r + 4, c + 3),
167 &m(r + 4, c + 4), &m(r + 4, c + 5),
168
169 &m(r + 5, c + 0), &m(r + 5, c + 1), &m(r + 5, c + 2), &m(r + 5, c + 3),
170 &m(r + 5, c + 4), &m(r + 5, c + 5)
171
172 );
173 };
174 FTensor::Index<'i', 3> i;
175 FTensor::Index<'j', 3> j;
176 FTensor::Index<'k', 3> k;
177 FTensor::Index<'l', 3> l;
178 FTensor::Index<'m', 3> m;
179 FTensor::Index<'n', 3> n;
180
181 auto v = getVolume();
182 auto t_w = getFTensor0IntegrationWeight();
183
184 auto get_dP = [&]() {
185 dP.resize(size_symm * size_symm, nb_integration_pts, false);
186 auto ts_a = getTSa();
187
188 auto t_P = getFTensor2FromMat<3, 3>(dataAtPts->PAtPts);
189 auto r_P_du = getFTensor4FromMat<3, 3, 3, 3>(dataAtPts->P_du);
190 auto t_approx_P_adjoint__dstretch =
191 getFTensor2FromMat<3, 3>(dataAtPts->adjointPdstretchAtPts);
192 auto t_dot_log_u =
193 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchDotTensorAtPts);
194 auto t_u = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
195 auto t_diff_u =
196 getFTensor4DdgFromMat<3, 3, 1>(dataAtPts->diffStretchTensorAtPts);
197 auto t_eigen_vals = getFTensor1FromMat<3>(dataAtPts->eigenVals);
198 auto t_eigen_vecs = getFTensor2FromMat<3, 3>(dataAtPts->eigenVecs);
199 auto &nbUniq = dataAtPts->nbUniq;
200
201 auto t_dP = getFTensor2FromMat<size_symm, size_symm>(dP);
202 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
203
205 t_deltaP(i, j) = t_approx_P_adjoint__dstretch(i, j) - t_P(i, j);
206
208 t_Ldiff_u(i, j, L) = t_diff_u(i, j, m, n) * t_L(m, n, L);
209 t_dP(L, J) =
210 -t_Ldiff_u(i, j, L) * (r_P_du(i, j, k, l) * t_Ldiff_u(k, l, J));
211 // viscous stress derivative
212 t_dP(L, J) -= (alphaU * ts_a) *
213 (t_L(i, j, L) * (t_diff(i, j, k, l) * t_L(k, l, J)));
214
218 t_deltaP_sym(i, j) = (t_deltaP(i, j) || t_deltaP(j, i));
219 t_deltaP_sym(i, j) /= 2.0;
220 auto t_diff2_uP2 = EigenMatrix::getDiffDiffMat(
221 t_eigen_vals, t_eigen_vecs, EshelbianCore::f, EshelbianCore::d_f,
222 EshelbianCore::dd_f, t_deltaP_sym, nbUniq[gg]);
223 t_dP(L, J) += t_L(i, j, L) * (t_diff2_uP2(i, j, k, l) * t_L(k, l, J));
224 }
225
226 ++t_P;
227 ++t_dP;
228 ++r_P_du;
229 ++t_approx_P_adjoint__dstretch;
230 ++t_dot_log_u;
231 ++t_u;
232 ++t_diff_u;
233 ++t_eigen_vals;
234 ++t_eigen_vecs;
235 }
236
237 return getFTensor2FromMat<size_symm, size_symm>(dP);
238 };
239
240 int row_nb_base_functions = row_data.getN().size2();
241 auto t_row_base_fun = row_data.getFTensor0N();
242
243 auto t_dP = get_dP();
244 for (int gg = 0; gg != nb_integration_pts; ++gg) {
245 double a = v * t_w;
246
247 int rr = 0;
248 for (; rr != row_nb_dofs / 6; ++rr) {
249 auto t_col_base_fun = col_data.getFTensor0N(gg, 0);
250 auto t_m = get_ftensor2(K, 6 * rr, 0);
251 for (int cc = 0; cc != col_nb_dofs / 6; ++cc) {
252 const double b = a * t_row_base_fun * t_col_base_fun;
253 t_m(L, J) -= b * t_dP(L, J);
254 ++t_m;
255 ++t_col_base_fun;
256 }
257 ++t_row_base_fun;
258 }
259
260 for (; rr != row_nb_base_functions; ++rr) {
261 ++t_row_base_fun;
262 }
263
264 ++t_w;
265 ++t_dP;
266 }
268}
269
271 std::string row_field, std::string col_field,
272 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha) {
273 return new OpSpatialPhysical_du_du(row_field, col_field, data_ptr, alpha);
274}
275
277 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
278 boost::shared_ptr<double> total_energy_ptr) {
279 return nullptr;
280}
281
283 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
284 boost::shared_ptr<PhysicalEquations> physics_ptr) {
285 return nullptr;
286}
287
289 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
290 boost::shared_ptr<PhysicalEquations> physics_ptr) {
291 return nullptr;
292};
293
295 OpGetScale(boost::shared_ptr<double> scale_ptr)
297 scalePtr(scale_ptr) {}
298
299 MoFEMErrorCode doWork(int side, EntityType type,
302 *(scalePtr) = 1.;
304 }
305
306private:
307 boost::shared_ptr<double> scalePtr;
308};
309
311 boost::shared_ptr<double> scale_ptr,
312 boost::shared_ptr<PhysicalEquations> physics_ptr) {
313 return new OpGetScale(scale_ptr);
314};
315
316struct OpHMHH : public OpJacobian {
317 OpHMHH(const int tag, const bool eval_rhs, const bool eval_lhs,
318 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
319 boost::shared_ptr<PhysicalEquations> physics_ptr)
320 : OpJacobian(tag, eval_rhs, eval_lhs, data_ptr, physics_ptr) {}
321
324};
325
328 FTensor::Index<'i', 3> i;
329 FTensor::Index<'j', 3> j;
330 FTensor::Index<'k', 3> k;
331 FTensor::Index<'l', 3> l;
332 FTensor::Index<'m', 3> m;
333 FTensor::Index<'n', 3> n;
334
335 const auto nb_integration_pts = getGaussPts().size2();
336 auto iu = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
337 auto t_grad_h1 = getFTensor2FromMat<3, 3>(dataAtPts->wGradH1AtPts);
338
339 auto create_data_vec = [nb_integration_pts](auto &v) {
340 v.resize(nb_integration_pts, false);
341 };
342 auto create_data_mat = [nb_integration_pts](auto &m) {
343 m.resize(9, nb_integration_pts, false);
344 };
345
346 create_data_mat(dataAtPts->PAtPts);
347
348 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
349 auto r_P = getFTensor2FromMat<3, 3>(dataAtPts->PAtPts);
350 for (unsigned int gg = 0; gg != nb_integration_pts; ++gg) {
351
353
356 t_h1(i, j) = t_kd(i, j);
357 physicsPtr->get_h()(i, j) = iu(i, j);
358 break;
359 case LARGE_ROT:
360 case MODERATE_ROT:
361 t_h1(i, j) = t_grad_h1(i, j) + t_kd(i, j);
362 physicsPtr->get_h()(i, j) = iu(i, k) * t_h1(k, j);
363 break;
364 case SMALL_ROT:
365 physicsPtr->get_h()(i, j) = iu(i, j);
366 break;
367 };
368
369 int r = ::function(tAg, physicsPtr->dependentVariablesPiola.size(),
370 physicsPtr->activeVariables.size(),
371 &physicsPtr->activeVariables[0],
372 &physicsPtr->dependentVariablesPiola[0]);
373 if (r < 0) { // function is locally analytic
374 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
375 "ADOL-C function evaluation with error r = %d", r);
376 }
377
380 r_P(i, j) = physicsPtr->get_P()(i, j);
381 break;
382 case LARGE_ROT:
383 case MODERATE_ROT: {
385 t_h1_du(i, j, m, n) = t_kd(i, m) * (t_kd(k, n) * t_h1(k, j));
386 r_P(k, l) = physicsPtr->get_P()(i, j) * t_h1_du(i, j, k, l);
387 }; break;
388 case SMALL_ROT:
389 r_P(i, j) = physicsPtr->get_P()(i, j);
390 break;
391 };
392
393 ++iu;
394 ++t_grad_h1;
395 ++r_P;
396 }
398}
399
402 FTensor::Index<'i', 3> i;
403 FTensor::Index<'j', 3> j;
404 FTensor::Index<'k', 3> k;
405 FTensor::Index<'l', 3> l;
406 FTensor::Index<'m', 3> m;
407 FTensor::Index<'n', 3> n;
408 FTensor::Index<'o', 3> o;
409 FTensor::Index<'p', 3> p;
410
411 const int number_of_active_variables = physicsPtr->activeVariables.size();
412 const int number_of_dependent_variables =
413 physicsPtr->dependentVariablesPiola.size();
414 std::vector<double *> jac_ptr(number_of_dependent_variables);
415 for (unsigned int n = 0; n != number_of_dependent_variables; ++n) {
416 jac_ptr[n] =
417 &(physicsPtr
418 ->dependentVariablesPiolaDirevatives[n *
419 number_of_active_variables]);
420 }
421
422 const auto nb_integration_pts = getGaussPts().size2();
423
424 auto create_data_mat = [nb_integration_pts](auto &m) {
425 m.resize(9, nb_integration_pts, false);
426 };
427
428 dataAtPts->P_du.resize(81, nb_integration_pts, false);
429
430 auto iu = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
431 auto t_grad_h1 = getFTensor2FromMat<3, 3>(dataAtPts->wGradH1AtPts);
432 auto r_P_du = getFTensor4FromMat<3, 3, 3, 3>(dataAtPts->P_du);
433
437
438 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
439
440 for (unsigned int gg = 0; gg != nb_integration_pts; ++gg) {
441
443
446 t_h1(i, j) = t_kd(i, j);
447 physicsPtr->get_h()(i, j) = iu(i, k) * t_h1(k, j);
448 break;
449 case LARGE_ROT:
450 case MODERATE_ROT:
451 t_h1(i, j) = t_grad_h1(i, j) + t_kd(i, j);
452 physicsPtr->get_h()(i, j) = iu(i, k) * t_h1(k, j);
453 break;
454 case SMALL_ROT:
455 physicsPtr->get_h()(i, j) = iu(i, j);
456 break;
457 };
458
459 // play recorder for jacobians
460 int r = ::jacobian(tAg, number_of_dependent_variables,
461 number_of_active_variables,
462 &physicsPtr->activeVariables[0], &jac_ptr[0]);
463 if (r < 0) {
464 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
465 "ADOL-C function evaluation with error");
466 }
467
469 t_P_dh_tmp(i, j, N0, k) = physicsPtr->get_P_dh0()(i, j, k);
470 t_P_dh_tmp(i, j, N1, k) = physicsPtr->get_P_dh1()(i, j, k);
471 t_P_dh_tmp(i, j, N2, k) = physicsPtr->get_P_dh2()(i, j, k);
472
474 case NO_H1_CONFIGURATION: {
476 t_h1_du(i, j, m, n) = t_kd(i, m) * t_kd(j, n);
477 r_P_du(k, l, m, n) =
478 (t_P_dh_tmp(i, j, o, p) * t_h1_du(o, p, m, n)) * t_h1_du(i, j, k, l);
479 } break;
480 case LARGE_ROT: {
482 t_h1_du(i, j, m, n) = t_kd(i, m) * (t_kd(k, n) * t_h1(k, j));
483 r_P_du(k, l, m, n) =
484 (t_P_dh_tmp(i, j, o, p) * t_h1_du(o, p, m, n)) * t_h1_du(i, j, k, l);
485 } break;
486 case MODERATE_ROT:
487 case SMALL_ROT:
488 r_P_du(i, j, m, n) = t_P_dh_tmp(i, j, m, n);
489 break;
490 };
491
492 ++iu;
493 ++t_grad_h1;
494 ++r_P_du;
495 }
497}
498
500 const int tag, const bool eval_rhs, const bool eval_lhs,
501 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
502 boost::shared_ptr<PhysicalEquations> physics_ptr) {
503 return (new OpHMHH(tag, eval_rhs, eval_lhs, data_ptr, physics_ptr));
504}
505
506} // namespace EshelbianPlasticity
507
509#include <impl/HMHNeohookean.cpp>
511#include <impl/HMHHencky.cpp>
512
513namespace EshelbianPlasticity {
514
516 const int tape, const double lambda, const double mu,
517 const double sigma_y) {
519 physicalEquations = boost::make_shared<HMHStVenantKirchhoff>(lambda, mu);
520 CHKERR physicalEquations->recordTape(tape, nullptr);
522}
523
525 const int tape, const double alpha, const double beta, const double lambda,
526 const double sigma_y) {
528 physicalEquations =
529 boost::make_shared<HMHPMooneyRivlinWriggersEq63>(alpha, beta, lambda);
530 CHKERR physicalEquations->recordTape(tape, nullptr);
532}
534 const double c10,
535 const double K) {
537 physicalEquations = boost::make_shared<HMHNeohookean>(mField, c10, K);
538 CHKERR physicalEquations->recordTape(tape, nullptr);
540}
541
544 physicalEquations = boost::make_shared<HMHHencky>(mField, E, nu);
546}
547
548}; // namespace EshelbianPlasticity
constexpr double third
Auxilary functions for Eshelbian plasticity.
Eshelbian plasticity interface.
Implementation StVenantKirchhoff.
Implementation of NeoHookean material.
Implement of MooneyRivlinWriggersEq63.
constexpr double a
Kronecker Delta class.
@ NOSPACE
Definition definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr auto t_kd
FTensor::Index< 'i', SPACE_DIM > i
static double lambda
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
const double n
refractive index of diffusive medium
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
auto getDiffDiffMat(A &&t_val, B &&t_vec, Fun< double > f, Fun< double > d_f, Fun< double > dd_f, C &&t_S, const int nb)
Get the Diff Diff Mat object.
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static constexpr auto size_symm
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
constexpr auto field_name
FTensor::Index< 'm', 3 > m
static enum StretchSelector stretchSelector
MoFEMErrorCode addMaterial_Hencky(double E, double nu)
MoFEMErrorCode addMaterial_HMHHStVenantKirchhoff(const int tape, const double lambda, const double mu, const double sigma_y)
MoFEMErrorCode addMaterial_HMHMooneyRivlin(const int tape, const double alpha, const double beta, const double lambda, const double sigma_y)
static enum RotSelector gradApproximator
static boost::function< double(const double)> f
static boost::function< double(const double)> dd_f
static boost::function< double(const double)> d_f
MoFEMErrorCode addMaterial_HMHNeohookean(const int tape, const double c10, const double K)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< double > scalePtr
OpGetScale(boost::shared_ptr< double > scale_ptr)
OpHMHH(const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
MoFEMErrorCode evaluateLhs(EntData &data)
MoFEMErrorCode evaluateRhs(EntData &data)
MoFEMErrorCode integratePiola(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
OpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
MoFEMErrorCode integrate(EntData &data)
virtual VolUserDataOperator * returnOpCalculateVarStretchFromStress(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
virtual VolUserDataOperator * returnOpCalculateEnergy(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< double > total_energy_ptr)
virtual VolUserDataOperator * returnOpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
virtual VolUserDataOperator * returnOpSpatialPhysicalExternalStrain(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< ExternalStrainVec > external_strain_vec_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv)
virtual VolUserDataOperator * returnOpCalculateStretchFromStress(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
virtual VolUserDataOperator * returnOpSetScale(boost::shared_ptr< double > scale_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
virtual VolUserDataOperator * returnOpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
virtual UserDataOperator * returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
bool sYmm
If true assume that matrix is symmetric structure.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
auto getFTensor0IntegrationWeight()
Get integration weights.
@ OPROW
operator doWork function is executed on FE rows
@ OPROWCOL
operator doWork is executed on FE rows &columns
@ OPSPACE
operator do Work is execute on space data
MatrixDouble K
local tangent matrix
VectorDouble nF
local right hand side vector
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts