7 #ifndef __POISSONOPERATORS_HPP__
8 #define __POISSONOPERATORS_HPP__
11 using namespace MoFEM;
60 nbIntegrationPts = getGaussPts().size2();
62 if (row_side == col_side && row_type == col_type) {
68 CHKERR iNtegrate(row_data, col_data);
70 CHKERR aSsemble(row_data, col_data);
95 locMat.resize(nbRows, nbCols,
false);
99 double vol = getVolume();
101 auto t_w = getFTensor0IntegrationWeight();
105 for (
int gg = 0; gg != nbIntegrationPts; gg++) {
107 const double alpha = t_w * vol;
113 for (
int rr = 0; rr != nbRows; rr++) {
117 for (
int cc = 0; cc != nbCols; cc++) {
119 a += alpha * (t_row_grad(
i) * t_col_grad(
i));
141 const int *row_indices = &*row_data.
getIndices().data().begin();
143 const int *col_indices = &*col_data.
getIndices().data().begin();
144 Mat B = getFEMethod()->ksp_B != PETSC_NULL ? getFEMethod()->ksp_B
145 : getFEMethod()->snes_B;
148 &*locMat.data().begin(), ADD_VALUES);
150 if (!isDiag && sYmm) {
153 locMat = trans(locMat);
155 &*locMat.data().begin(), ADD_VALUES);
183 nbIntegrationPts = OPBASE::getGaussPts().size2();
185 CHKERR iNtegrate(row_data);
187 CHKERR aSsemble(row_data);
220 :
public OpBaseRhs<VolumeElementForcesAndSourcesCore::UserDataOperator> {
222 typedef boost::function<
double(
const double,
const double,
const double)>
245 locVec.resize(nbRows,
false);
249 double vol = getVolume();
251 auto t_w = getFTensor0IntegrationWeight();
255 auto t_coords = getFTensor1CoordsAtGaussPts();
257 for (
int gg = 0; gg != nbIntegrationPts; gg++) {
260 vol * t_w * fSource(t_coords(NX), t_coords(NY), t_coords(NZ));
263 &*locVec.data().begin());
265 for (
int rr = 0; rr != nbRows; rr++) {
285 const int *indices = &*data.
getIndices().data().begin();
287 const double *vals = &*locVec.data().begin();
288 Vec f = getFEMethod()->ksp_f != PETSC_NULL ? getFEMethod()->ksp_f
289 : getFEMethod()->snes_f;
307 OpC(
const bool assemble_transpose)
310 assembleTranspose(assemble_transpose) {}
328 nbIntegrationPts = getGaussPts().size2();
330 CHKERR iNtegrate(row_data, col_data);
332 CHKERR aSsemble(row_data, col_data);
352 locMat.resize(nbRows, nbCols,
false);
356 const double area = getArea();
358 auto t_w = getFTensor0IntegrationWeight();
362 for (
int gg = 0; gg != nbIntegrationPts; gg++) {
363 const double alpha = area * t_w;
366 &*locMat.data().begin());
368 for (
int rr = 0; rr != nbRows; rr++) {
372 for (
int cc = 0; cc != nbCols; cc++) {
374 c += alpha * t_row * t_col;
392 const int *row_indices = &*row_data.
getIndices().data().begin();
394 const int *col_indices = &*col_data.
getIndices().data().begin();
395 Mat B = getFEMethod()->ksp_B != PETSC_NULL ? getFEMethod()->ksp_B
396 : getFEMethod()->snes_B;
399 &*locMat.data().begin(), ADD_VALUES);
401 if (assembleTranspose) {
403 locMat = trans(locMat);
405 &*locMat.data().begin(), ADD_VALUES);
420 :
public OpBaseRhs<FaceElementForcesAndSourcesCore::UserDataOperator> {
422 typedef boost::function<
double(
const double,
const double,
const double)>
428 fValue(f_value), bEta(beta) {}
445 locVec.resize(nbRows,
false);
449 const double area = getArea() * bEta;
451 auto t_w = getFTensor0IntegrationWeight();
455 auto t_coords = getFTensor1CoordsAtGaussPts();
457 for (
int gg = 0; gg != nbIntegrationPts; gg++) {
461 area * t_w * fValue(t_coords(NX), t_coords(NY), t_coords(NZ));
464 &*locVec.data().begin());
466 for (
int rr = 0; rr != nbRows; rr++) {
482 const int *indices = &*data.
getIndices().data().begin();
483 const double *vals = &*locVec.data().begin();
484 Vec f = getFEMethod()->ksp_f != PETSC_NULL ? getFEMethod()->ksp_f
485 : getFEMethod()->snes_f;
495 :
public OpBaseRhs<VolumeElementForcesAndSourcesCore::UserDataOperator> {
497 typedef boost::function<
double(
const double,
const double,
const double)>
499 typedef boost::function<FTensor::Tensor1<double, 3>(
500 const double,
const double,
const double)>
504 boost::shared_ptr<VectorDouble> &field_vals,
505 boost::shared_ptr<MatrixDouble> &grad_vals,
Vec global_error)
507 globalError(global_error), uValue(u_value), gValue(g_value),
508 fieldVals(field_vals), gradVals(grad_vals) {}
516 nbIntegrationPts = getGaussPts().size2();
517 CHKERR iNtegrate(row_data);
518 CHKERR aSsemble(row_data);
543 const double vol = getVolume();
545 auto t_w = getFTensor0IntegrationWeight();
549 auto t_grad = getFTensor1FromMat<3>(*gradVals);
551 auto t_coords = getFTensor1CoordsAtGaussPts();
555 for (
int gg = 0; gg != nbIntegrationPts; gg++) {
556 double alpha = vol * t_w;
558 double exact_u = uValue(t_coords(NX), t_coords(NY), t_coords(NZ));
560 t_exact_grad = gValue(t_coords(NX), t_coords(NY), t_coords(NZ));
562 t_error_grad(
i) = t_grad(
i) - t_exact_grad(
i);
564 double error = pow(t_u - exact_u, 2) + t_error_grad(
i) * t_error_grad(
i);
590 OpKt(boost::function<
double(
const double)>
a,
591 boost::function<
double(
const double)> diff_a,
592 boost::shared_ptr<VectorDouble> &field_vals,
593 boost::shared_ptr<MatrixDouble> &grad_vals)
594 :
OpK(false),
A(
a), diffA(diff_a), fieldVals(field_vals),
595 gradVals(grad_vals) {}
608 locMat.resize(nbRows, nbCols,
false);
612 double vol = getVolume();
614 auto t_w = getFTensor0IntegrationWeight();
618 auto t_grad = getFTensor1FromMat<3>(*gradVals);
622 for (
int gg = 0; gg != nbIntegrationPts; gg++) {
624 const double alpha = t_w * vol;
625 const double beta = alpha *
A(t_u);
627 t_gamma(
i) = (alpha * diffA(t_u)) * t_grad(
i);
630 &*locMat.data().begin());
632 for (
int rr = 0; rr != nbRows; rr++) {
638 for (
int cc = 0; cc != nbCols; cc++) {
640 a += (t_row_grad(
i) * beta) * t_col_grad(
i) +
641 t_row_grad(
i) * (t_gamma(
i) * t_col);
665 boost::shared_ptr<VectorDouble> &field_vals,
666 boost::shared_ptr<MatrixDouble> &grad_vals)
667 :
OpF(f_source),
A(
a), fieldVals(field_vals), gradVals(grad_vals) {}
680 locVec.resize(nbRows,
false);
684 double vol = getVolume();
686 auto t_w = getFTensor0IntegrationWeight();
690 auto t_grad = getFTensor1FromMat<3>(*gradVals);
696 auto t_coords = getFTensor1CoordsAtGaussPts();
698 for (
int gg = 0; gg != nbIntegrationPts; gg++) {
700 const double alpha = vol * t_w;
701 const double source_term =
702 alpha * fSource(t_coords(NX), t_coords(NY), t_coords(NZ));
704 grad_term(
i) = (alpha *
A(t_u)) * t_grad(
i);
707 &*locVec.data().begin());
709 for (
int rr = 0; rr != nbRows; rr++) {
711 t_a += t_v_grad(
i) * grad_term(
i) + t_v * source_term;
731 OpRes_g(
FVal f_value, boost::shared_ptr<VectorDouble> &field_vals)
732 :
Op_g(f_value,
"L", 1), fieldVals(field_vals) {}
743 locVec.resize(nbRows,
false);
747 const double area = getArea() * bEta;
749 auto t_w = getFTensor0IntegrationWeight();
755 auto t_coords = getFTensor1CoordsAtGaussPts();
757 for (
int gg = 0; gg != nbIntegrationPts; gg++) {
760 double alpha = area * t_w;
763 &*locVec.data().begin());
764 for (
int rr = 0; rr != nbRows; rr++) {
766 (t_u - fValue(t_coords(NX), t_coords(NY), t_coords(NZ)));
781 :
Op_g(
FVal(),
"U", 1), lambdaVals(lambda_vals) {}
792 locVec.resize(nbRows,
false);
796 const double area = getArea() * bEta;
798 auto t_w = getFTensor0IntegrationWeight();
804 for (
int gg = 0; gg != nbIntegrationPts; gg++) {
807 double alpha = area * t_w;
810 &*locVec.data().begin());
811 for (
int rr = 0; rr != nbRows; rr++) {
812 t_a += alpha * t_u * t_lambda;
835 int operator()(
int,
int,
int p)
const {
return 2 * (p - 1); }
851 return 2 * p_data + 1;
869 boost::function<
double(
const double,
const double,
const double)> f_u,
870 boost::function<
double(
const double,
const double,
const double)>
872 boost::shared_ptr<ForcesAndSourcesCore> &domain_lhs_fe,
873 boost::shared_ptr<ForcesAndSourcesCore> &boundary_lhs_fe,
874 boost::shared_ptr<ForcesAndSourcesCore> &domain_rhs_fe,
875 boost::shared_ptr<ForcesAndSourcesCore> &boundary_rhs_fe,
876 bool trans =
true)
const {
880 domain_lhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
882 boundary_lhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
884 domain_rhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
886 boundary_rhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
890 domain_lhs_fe->getRuleHook =
VolRule();
891 domain_rhs_fe->getRuleHook =
VolRule();
892 boundary_lhs_fe->getRuleHook =
FaceRule();
893 boundary_rhs_fe->getRuleHook =
FaceRule();
897 domain_lhs_fe->getOpPtrVector().push_back(
new OpK());
899 domain_rhs_fe->getOpPtrVector().push_back(
new OpF(f_source));
901 boundary_lhs_fe->getOpPtrVector().push_back(
new OpC(trans));
903 boundary_rhs_fe->getOpPtrVector().push_back(
new Op_g(f_u));
912 boost::function<
double(
const double,
const double,
const double)> f_u,
917 boost::shared_ptr<ForcesAndSourcesCore> &domain_error)
const {
920 domain_error = boost::shared_ptr<ForcesAndSourcesCore>(
922 domain_error->getRuleHook =
VolRule();
926 boost::shared_ptr<VectorDouble> values_at_integration_ptr =
927 boost::make_shared<VectorDouble>();
929 boost::shared_ptr<MatrixDouble> grad_at_integration_ptr =
930 boost::make_shared<MatrixDouble>();
932 domain_error->getOpPtrVector().push_back(
935 domain_error->getOpPtrVector().push_back(
938 domain_error->getOpPtrVector().push_back(
939 new OpError(f_u, g_u, values_at_integration_ptr,
940 grad_at_integration_ptr, global_error));
948 boost::shared_ptr<PostProcFE> &post_proc_volume)
const {
962 boost::shared_ptr<PostProcFE>(
new PostProcFE(mField));
965 auto det_ptr = boost::make_shared<VectorDouble>();
966 auto jac_ptr = boost::make_shared<MatrixDouble>();
967 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
968 post_proc_volume->getOpPtrVector().push_back(
970 post_proc_volume->getOpPtrVector().push_back(
972 post_proc_volume->getOpPtrVector().push_back(
975 auto u_ptr = boost::make_shared<VectorDouble>();
976 auto grad_ptr = boost::make_shared<MatrixDouble>();
977 auto e_ptr = boost::make_shared<VectorDouble>();
979 post_proc_volume->getOpPtrVector().push_back(
981 post_proc_volume->getOpPtrVector().push_back(
983 post_proc_volume->getOpPtrVector().push_back(
988 post_proc_volume->getOpPtrVector().push_back(
992 post_proc_volume->getPostProcMesh(),
993 post_proc_volume->getMapGaussPts(),
995 {{
"U", u_ptr}, {
"ERROR", e_ptr}},
997 {{
"GRAD", grad_ptr}},
1014 boost::function<
double(
const double,
const double,
const double)> f_u,
1015 boost::function<
double(
const double,
const double,
const double)>
1017 boost::function<
double(
const double)>
a,
1018 boost::function<
double(
const double)> diff_a,
1019 boost::shared_ptr<ForcesAndSourcesCore> &domain_lhs_fe,
1020 boost::shared_ptr<ForcesAndSourcesCore> &boundary_lhs_fe,
1021 boost::shared_ptr<ForcesAndSourcesCore> &domain_rhs_fe,
1022 boost::shared_ptr<ForcesAndSourcesCore> &boundary_rhs_fe,
1025 bool trans =
true)
const {
1029 domain_lhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
1031 boundary_lhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
1033 domain_rhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
1035 boundary_rhs_fe = boost::shared_ptr<ForcesAndSourcesCore>(
1039 domain_lhs_fe->getRuleHook = vol_rule;
1040 domain_rhs_fe->getRuleHook = vol_rule;
1041 boundary_lhs_fe->getRuleHook = face_rule;
1042 boundary_rhs_fe->getRuleHook = face_rule;
1047 boost::shared_ptr<VectorDouble> values_at_integration_ptr =
1048 boost::make_shared<VectorDouble>();
1050 boost::shared_ptr<MatrixDouble> grad_at_integration_ptr =
1051 boost::make_shared<MatrixDouble>();
1053 boost::shared_ptr<VectorDouble> multiplier_at_integration_ptr =
1054 boost::make_shared<VectorDouble>();
1058 domain_lhs_fe->getOpPtrVector().push_back(
1061 domain_lhs_fe->getOpPtrVector().push_back(
1064 domain_lhs_fe->getOpPtrVector().push_back(
new OpKt(
1065 a, diff_a, values_at_integration_ptr, grad_at_integration_ptr));
1068 domain_rhs_fe->getOpPtrVector().push_back(
1071 domain_rhs_fe->getOpPtrVector().push_back(
1074 domain_rhs_fe->getOpPtrVector().push_back(
new OpResF_Domain(
1075 f_source,
a, values_at_integration_ptr, grad_at_integration_ptr));
1078 boundary_lhs_fe->getOpPtrVector().push_back(
new OpC(trans));
1081 boundary_rhs_fe->getOpPtrVector().push_back(
1084 boundary_rhs_fe->getOpPtrVector().push_back(
1087 boundary_rhs_fe->getOpPtrVector().push_back(
1088 new OpRes_g(f_u, values_at_integration_ptr));
1089 boundary_rhs_fe->getOpPtrVector().push_back(
1101 #endif //__POISSONOPERATORS_HPP__