static char help[] =
"...\n\n";
constexpr
bool debug =
true;
private:
using MatSideArray = std::array<MatrixDouble, 2>;
struct SideData {
std::array<EntityHandle, 2> feSideHandle;
std::array<VectorInt, 2>
std::array<VectorInt, 2>
MatSideArray lVec;
MatSideArray velMat;
int currentFESide;
};
template <int FE_DIM>
static double get_velocity_potential(double x, double y, double z);
static double get_level_set(const double x, const double y, const double z);
std::tuple<double, Tag> evaluateError();
getZeroLevelVelOp(boost::shared_ptr<MatrixDouble> vel_ptr);
boost::shared_ptr<FaceSideEle>
getSideFE(boost::shared_ptr<SideData> side_data_ptr);
boost::function<double(double, double, double)> level_fun =
get_level_set);
boost::function<double(double, double, double)> vel_fun =
get_velocity_potential<FE_DIM>);
struct WrapperClass {
WrapperClass() = default;
virtual double getThreshold(const double max) = 0;
};
struct WrapperClassInitalSolution : public WrapperClass {
WrapperClassInitalSolution(boost::shared_ptr<double> max_ptr)
: WrapperClass(), maxPtr(max_ptr) {}
};
}
->synchroniseEntities(level);
}
double getThreshold(const double max) {
*maxPtr = std::max(*maxPtr, max);
return 0.05 * (*maxPtr);
}
private:
boost::shared_ptr<double> maxPtr;
};
struct WrapperClassErrorProjection : public WrapperClass {
WrapperClassErrorProjection(boost::shared_ptr<double> max_ptr)
: maxPtr(max_ptr) {}
->synchroniseEntities(level);
}
double getThreshold(const double max) { return 0.05 * (*maxPtr); }
private:
boost::shared_ptr<double> maxPtr;
};
struct OpRhsDomain;
struct OpLhsDomain;
struct OpRhsSkeleton;
struct OpLhsSkeleton;
G>::OpSource<1, DIM1 * DIM2>;
G>::OpSource<potential_velocity_field_dim, potential_velocity_field_dim>;
G>::OpBaseTimesVector<1, DIM1 * DIM2, 1>;
private:
boost::shared_ptr<double> maxPtr;
};
template <>
double LevelSet::get_velocity_potential<2>(double x, double y, double z) {
return (x * x - 0.25) * (y * y - 0.25);
}
constexpr
double xc = 0.1;
constexpr
double yc = 0.;
constexpr
double zc = 0.;
constexpr
double r = 0.2;
return std::sqrt(pow(x -
xc, 2) + pow(y -
yc, 2) + pow(z -
zc, 2)) -
r;
}
}
CHKERR initialiseFieldVelocity();
maxPtr = boost::make_shared<double>(0);
CHKERR refineMesh(WrapperClassInitalSolution(maxPtr));
}
boost::shared_ptr<MatrixDouble> l_ptr,
boost::shared_ptr<MatrixDouble> l_dot_ptr,
boost::shared_ptr<MatrixDouble> vel_ptr);
private:
boost::shared_ptr<MatrixDouble> lPtr;
boost::shared_ptr<MatrixDouble> lDotPtr;
boost::shared_ptr<MatrixDouble> velPtr;
};
boost::shared_ptr<MatrixDouble> vel_ptr);
private:
boost::shared_ptr<MatrixDouble> velPtr;
};
OpRhsSkeleton(boost::shared_ptr<SideData> side_data_ptr,
boost::shared_ptr<FaceSideEle> side_fe_ptr);
private:
boost::shared_ptr<SideData> sideDataPtr;
boost::shared_ptr<FaceSideEle>
sideFEPtr;
};
OpLhsSkeleton(boost::shared_ptr<SideData> side_data_ptr,
boost::shared_ptr<FaceSideEle> side_fe_ptr);
private:
boost::shared_ptr<SideData> sideDataPtr;
boost::shared_ptr<FaceSideEle>
sideFEPtr;
};
int main(
int argc,
char *argv[]) {
const char param_file[] = "param_file.petsc";
try {
DMType dm_name = "DMMOFEM";
auto core_log = logging::core::get();
core_log->add_sink(
}
return 0;
}
simple->getAddSkeletonFE() =
true;
simple->getAddBoundaryFE() =
true;
auto set_problem_bit = [&]() {
start_mask[s] = true;
#ifndef NDEBUG
auto proc_str = boost::lexical_cast<std::string>(mField.get_comm_rank());
CHKERR bit_mng->writeBitLevelByDim(
(proc_str + "level_base.vtk").c_str(), "VTK", "");
}
#endif
};
}
}
boost::shared_ptr<MatrixDouble> l_ptr,
boost::shared_ptr<MatrixDouble> l_dot_ptr,
boost::shared_ptr<MatrixDouble> vel_ptr)
lPtr(l_ptr), lDotPtr(l_dot_ptr), velPtr(vel_ptr) {}
boost::shared_ptr<MatrixDouble> vel_ptr)
velPtr(vel_ptr) {
this->sYmm = false;
}
boost::shared_ptr<SideData> side_data_ptr,
boost::shared_ptr<FaceSideEle> side_fe_ptr)
sideDataPtr(side_data_ptr), sideFEPtr(side_fe_ptr) {}
boost::shared_ptr<SideData> side_data_ptr,
boost::shared_ptr<FaceSideEle> side_fe_ptr)
sideDataPtr(side_data_ptr), sideFEPtr(side_fe_ptr) {}
const auto nb_int_points = getGaussPts().size2();
const auto nb_base_func = data.
getN().size2();
auto t_l = getFTensor2FromMat<DIM1, DIM2>(*
lPtr);
auto t_l_dot = getFTensor2FromMat<DIM1, DIM2>(*
lDotPtr);
auto t_vel = getFTensor1FromMat<SPACE_DIM>(*
velPtr);
auto t_w = getFTensor0IntegrationWeight();
for (auto gg = 0; gg != nb_int_points; ++gg) {
const auto alpha = t_w * getMeasure();
t_res0(
I,
J) = alpha * t_l_dot(
I,
J);
t_res1(
i,
I,
J) = (alpha * t_l(
I,
J)) * t_vel(
i);
++t_w;
++t_l;
++t_l_dot;
++t_vel;
auto t_nf = getFTensor2FromPtr<DIM1, DIM2>(&*nf.begin());
int rr = 0;
for (; rr != nb_dofs; ++rr) {
t_nf(
I,
J) += t_res0(
I,
J) * t_base - t_res1(
i,
I,
J) * t_diff_base(
i);
++t_base;
++t_diff_base;
++t_nf;
}
for (; rr < nb_base_func; ++rr) {
++t_base;
++t_diff_base;
}
}
}
const auto nb_int_points = getGaussPts().size2();
const auto nb_base_func = row_data.
getN().size2();
const auto nb_row_dofs = row_data.
getIndices().size();
const auto nb_col_dofs = col_data.
getIndices().size();
auto t_vel = getFTensor1FromMat<SPACE_DIM>(*velPtr);
auto t_w = getFTensor0IntegrationWeight();
for (auto gg = 0; gg != nb_int_points; ++gg) {
const auto alpha = t_w * getMeasure();
const auto beta = alpha * getTSa();
++t_w;
auto &mat = this->locMat;
int rr = 0;
for (; rr != nb_row_dofs; ++rr) {
auto t_mat = getFTensor2FromPtr<DIM1, DIM2>(&mat(rr *
DIM1, 0));
for (int cc = 0; cc != nb_col_dofs; ++cc) {
(beta * t_row_base - alpha * (t_row_diff_base(
i) * t_vel(
i))) *
t_col_base;
++t_col_base;
++t_mat;
}
++t_row_base;
++t_row_diff_base;
}
for (; rr < nb_base_func; ++rr) {
++t_row_base;
++t_row_diff_base;
}
++t_vel;
}
}
CHKERR loopSideFaces(
"dFE", *sideFEPtr);
const auto in_the_loop =
sideFEPtr->nInTheLoop;
auto not_side = [](auto s) {
};
&*base_mat.data().begin());
};
if (in_the_loop > 0) {
auto t_normal = getFTensor1Normal();
const auto nb_gauss_pts = getGaussPts().size2();
const auto nb_rows = sideDataPtr->indicesRowSideMap[s0].size();
if (nb_rows) {
resSkelton.resize(nb_rows, false);
resSkelton.clear();
const auto opposite_s0 = not_side(s0);
const auto sense_row = sideDataPtr->senseMap[s0];
#ifndef NDEBUG
const auto opposite_sense_row = sideDataPtr->senseMap[opposite_s0];
if (sense_row * opposite_sense_row > 0)
"Should be opposite sign");
#endif
const auto nb_row_base_functions =
sideDataPtr->rowBaseSideMap[s0].size2();
auto t_w = getFTensor0IntegrationWeight();
getFTensor2FromMat<DIM1, DIM2>(sideDataPtr->lVec[
LEFT_SIDE]),
getFTensor2FromMat<DIM1, DIM2>(sideDataPtr->lVec[
RIGHT_SIDE]));
getFTensor1FromMat<SPACE_DIM>(sideDataPtr->velMat[
LEFT_SIDE]),
getFTensor1FromMat<SPACE_DIM>(sideDataPtr->velMat[
RIGHT_SIDE]));
auto next = [&]() {
for (auto &t_l : arr_t_l)
++t_l;
for (auto &t_vel : arr_t_vel)
++t_vel;
};
#ifndef NDEBUG
if (nb_gauss_pts != sideDataPtr->rowBaseSideMap[s0].size1())
"Inconsistent number of DOFs");
#endif
auto t_row_base =
get_ntensor(sideDataPtr->rowBaseSideMap[s0]);
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
const auto dot = sense_row * (t_normal(
i) * t_vel(
i));
const auto l_upwind_side = (dot > 0) ? s0 : opposite_s0;
const auto l_upwind = arr_t_l[l_upwind_side];
t_res(
I,
J) = t_w * dot * l_upwind(
I,
J);
next();
++t_w;
auto t_res_skeleton =
getFTensor2FromPtr<DIM1, DIM2>(&*resSkelton.data().begin());
auto rr = 0;
for (; rr != nb_rows; ++rr) {
t_res_skeleton(
I,
J) += t_row_base * t_res(
I,
J);
++t_row_base;
++t_res_skeleton;
}
for (; rr < nb_row_base_functions; ++rr) {
++t_row_base;
}
}
sideDataPtr->indicesRowSideMap[s0].size(),
&*sideDataPtr->indicesRowSideMap[s0].begin(),
&*resSkelton.begin(), ADD_VALUES);
}
}
}
}
CHKERR loopSideFaces(
"dFE", *sideFEPtr);
const auto in_the_loop =
sideFEPtr->nInTheLoop;
auto not_side = [](auto s) {
};
&*base_mat.data().begin());
};
if (in_the_loop > 0) {
auto t_normal = getFTensor1Normal();
const auto nb_gauss_pts = getGaussPts().size2();
const auto nb_rows = sideDataPtr->indicesRowSideMap[s0].size();
if (nb_rows) {
const auto opposite_s0 = not_side(s0);
const auto sense_row = sideDataPtr->senseMap[s0];
const auto nb_row_base_functions =
sideDataPtr->rowBaseSideMap[s0].size2();
const auto nb_cols = sideDataPtr->indicesColSideMap[s1].size();
const auto sense_col = sideDataPtr->senseMap[s1];
matSkeleton.resize(nb_rows, nb_cols, false);
matSkeleton.clear();
auto t_w = getFTensor0IntegrationWeight();
getFTensor1FromMat<SPACE_DIM>(sideDataPtr->velMat[
LEFT_SIDE]),
getFTensor1FromMat<SPACE_DIM>(sideDataPtr->velMat[
RIGHT_SIDE]));
auto next = [&]() {
for (auto &t_vel : arr_t_vel)
++t_vel;
};
auto t_row_base =
get_ntensor(sideDataPtr->rowBaseSideMap[s0]);
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
const auto dot = sense_row * (t_normal(
i) * t_vel(
i));
const auto l_upwind_side = (dot > 0) ? s0 : opposite_s0;
const auto sense_upwind = sideDataPtr->senseMap[l_upwind_side];
next();
++t_w;
auto rr = 0;
if (s1 == l_upwind_side) {
for (; rr != nb_rows; ++rr) {
auto get_ntensor = [](
auto &base_mat,
auto gg,
auto bb) {
double *ptr = &base_mat(gg, bb);
};
auto t_col_base =
auto t_mat_skeleton =
getFTensor2FromPtr<DIM1, DIM2>(&matSkeleton(rr *
DIM1, 0));
t_res_row(
I,
J) = t_res(
I,
J) * t_row_base;
++t_row_base;
for (size_t cc = 0; cc != nb_cols; ++cc) {
t_mat_skeleton(
I,
J) += t_res_row(
I,
J) * t_col_base;
++t_col_base;
++t_mat_skeleton;
}
}
}
for (; rr < nb_row_base_functions; ++rr) {
++t_row_base;
}
}
sideDataPtr->indicesRowSideMap[s0].size(),
&*sideDataPtr->indicesRowSideMap[s0].begin(),
sideDataPtr->indicesColSideMap[s1].size(),
&*sideDataPtr->indicesColSideMap[s1].begin(),
&*matSkeleton.data().begin(), ADD_VALUES);
}
}
}
}
}
auto get_parent_vel_this = [&]() {
auto parent_fe_ptr = boost::make_shared<DomianParentEle>(
mField);
parent_fe_ptr->getOpPtrVector(), {potential_velocity_space});
parent_fe_ptr->getOpPtrVector().push_back(
"V", vel_ptr));
return parent_fe_ptr;
};
auto get_parents_vel_fe_ptr = [&](auto this_fe_ptr) {
std::vector<boost::shared_ptr<DomianParentEle>> parents_elems_ptr_vec;
parents_elems_ptr_vec.emplace_back(
boost::make_shared<DomianParentEle>(
mField));
parents_elems_ptr_vec[
l - 1]->getOpPtrVector().push_back(
}
return parents_elems_ptr_vec[0];
};
auto this_fe_ptr = get_parent_vel_this();
auto parent_fe_ptr = get_parents_vel_fe_ptr(this_fe_ptr);
}
pip->getOpDomainRhsPipeline().clear();
pip->setDomainLhsIntegrationRule([](int, int, int o) { return 3 * o; });
pip->setDomainRhsIntegrationRule([](int, int, int o) { return 3 * o; });
pip->getDomainLhsFE()->exeTestHook = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
pip->getDomainRhsFE()->exeTestHook = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
auto l_ptr = boost::make_shared<MatrixDouble>();
auto l_dot_ptr = boost::make_shared<MatrixDouble>();
auto vel_ptr = boost::make_shared<MatrixDouble>();
{L2});
pip->getOpDomainRhsPipeline().push_back(
pip->getOpDomainRhsPipeline().push_back(
pip->getOpDomainRhsPipeline().push_back(
new OpRhsDomain("L", l_ptr, l_dot_ptr, vel_ptr));
{L2});
pip->getOpDomainLhsPipeline().push_back(new OpLhsDomain("L", vel_ptr));
}
boost::shared_ptr<FaceSideEle>
auto l_ptr = boost::make_shared<MatrixDouble>();
auto vel_ptr = boost::make_shared<MatrixDouble>();
OpSideData(boost::shared_ptr<SideData> side_data_ptr)
sideDataPtr(side_data_ptr) {
std::fill(&doEntities[MBVERTEX], &doEntities[MBMAXTYPE], false);
for (
auto t = moab::CN::TypeDimensionMap[
FE_DIM].first;
t <= moab::CN::TypeDimensionMap[
FE_DIM].second; ++
t)
sYmm = false;
}
MoFEMErrorCode doWork(
int row_side,
int col_side, EntityType row_type,
EntityType col_type,
EntData &row_data,
if ((CN::Dimension(row_type) ==
FE_DIM) &&
(CN::Dimension(col_type) ==
FE_DIM)) {
auto reset = [&](auto nb_in_loop) {
sideDataPtr->feSideHandle[nb_in_loop] = 0;
sideDataPtr->indicesRowSideMap[nb_in_loop].clear();
sideDataPtr->indicesColSideMap[nb_in_loop].clear();
sideDataPtr->rowBaseSideMap[nb_in_loop].clear();
sideDataPtr->colBaseSideMap[nb_in_loop].clear();
sideDataPtr->senseMap[nb_in_loop] = 0;
};
const auto nb_in_loop = getFEMethod()->nInTheLoop;
if (nb_in_loop == 0)
for (auto s : {0, 1})
reset(s);
sideDataPtr->currentFESide = nb_in_loop;
sideDataPtr->senseMap[nb_in_loop] = getSkeletonSense();
} else {
}
};
private:
boost::shared_ptr<SideData> sideDataPtr;
};
OpSideDataOnParent(boost::shared_ptr<SideData> side_data_ptr,
boost::shared_ptr<MatrixDouble> l_ptr,
boost::shared_ptr<MatrixDouble> vel_ptr)
sideDataPtr(side_data_ptr), lPtr(l_ptr), velPtr(vel_ptr) {
std::fill(&doEntities[MBVERTEX], &doEntities[MBMAXTYPE], false);
for (
auto t = moab::CN::TypeDimensionMap[
FE_DIM].first;
t <= moab::CN::TypeDimensionMap[
FE_DIM].second; ++
t)
sYmm = false;
}
MoFEMErrorCode doWork(
int row_side,
int col_side, EntityType row_type,
EntityType col_type,
EntData &row_data,
if ((CN::Dimension(row_type) ==
FE_DIM) &&
(CN::Dimension(col_type) ==
FE_DIM)) {
const auto nb_in_loop = sideDataPtr->currentFESide;
sideDataPtr->feSideHandle[nb_in_loop] = getFEEntityHandle();
sideDataPtr->indicesRowSideMap[nb_in_loop] = row_data.
getIndices();
sideDataPtr->indicesColSideMap[nb_in_loop] = col_data.
getIndices();
sideDataPtr->rowBaseSideMap[nb_in_loop] = row_data.
getN();
sideDataPtr->colBaseSideMap[nb_in_loop] = col_data.
getN();
(sideDataPtr->lVec)[nb_in_loop] = *lPtr;
(sideDataPtr->velMat)[nb_in_loop] = *velPtr;
#ifndef NDEBUG
if ((sideDataPtr->lVec)[nb_in_loop].size2() !=
(sideDataPtr->velMat)[nb_in_loop].size2())
"Wrong number of integaration pts %d != %d",
(sideDataPtr->lVec)[nb_in_loop].size2(),
(sideDataPtr->velMat)[nb_in_loop].size2());
if ((sideDataPtr->velMat)[nb_in_loop].size1() !=
SPACE_DIM)
"Wrong size of velocity vector size = %d",
(sideDataPtr->velMat)[nb_in_loop].size1());
#endif
if (!nb_in_loop) {
(sideDataPtr->lVec)[1] = sideDataPtr->lVec[0];
(sideDataPtr->velMat)[1] = (sideDataPtr->velMat)[0];
} else {
#ifndef NDEBUG
if (sideDataPtr->rowBaseSideMap[0].size1() !=
sideDataPtr->rowBaseSideMap[1].size1()) {
"Wrong number of integration pt %d != %d",
sideDataPtr->rowBaseSideMap[0].size1(),
sideDataPtr->rowBaseSideMap[1].size1());
}
if (sideDataPtr->colBaseSideMap[0].size1() !=
sideDataPtr->colBaseSideMap[1].size1()) {
"Wrong number of integration pt");
}
#endif
}
} else {
}
};
private:
boost::shared_ptr<SideData> sideDataPtr;
boost::shared_ptr<MatrixDouble> lPtr;
boost::shared_ptr<MatrixDouble> velPtr;
};
auto get_parent_this = [&]() {
auto parent_fe_ptr = boost::make_shared<DomianParentEle>(
mField);
parent_fe_ptr->getOpPtrVector(), {L2});
parent_fe_ptr->getOpPtrVector().push_back(
parent_fe_ptr->getOpPtrVector().push_back(
new OpSideDataOnParent(side_data_ptr, l_ptr, vel_ptr));
return parent_fe_ptr;
};
auto get_parents_fe_ptr = [&](auto this_fe_ptr) {
std::vector<boost::shared_ptr<DomianParentEle>> parents_elems_ptr_vec;
parents_elems_ptr_vec.emplace_back(
boost::make_shared<DomianParentEle>(
mField));
parents_elems_ptr_vec[
l - 1]->getOpPtrVector().push_back(
}
return parents_elems_ptr_vec[0];
};
auto get_side_fe_ptr = [&]() {
auto side_fe_ptr = boost::make_shared<FaceSideEle>(
mField);
auto this_fe_ptr = get_parent_this();
auto parent_fe_ptr = get_parents_fe_ptr(this_fe_ptr);
side_fe_ptr->getOpPtrVector().push_back(new OpSideData(side_data_ptr));
side_fe_ptr->getOpPtrVector().push_back(
return side_fe_ptr;
};
return get_side_fe_ptr();
};
pip->getOpSkeletonRhsPipeline().clear();
pip->setSkeletonLhsIntegrationRule([](int, int, int o) { return 18; });
pip->setSkeletonRhsIntegrationRule([](int, int, int o) { return 18; });
pip->getSkeletonLhsFE()->exeTestHook = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
pip->getSkeletonRhsFE()->exeTestHook = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
auto side_data_ptr = boost::make_shared<SideData>();
pip->getOpSkeletonRhsPipeline().push_back(
new OpRhsSkeleton(side_data_ptr, side_fe_ptr));
pip->getOpSkeletonLhsPipeline().push_back(
new OpLhsSkeleton(side_data_ptr, side_fe_ptr));
}
OpErrorSkel(boost::shared_ptr<FaceSideEle> side_fe_ptr,
boost::shared_ptr<SideData> side_data_ptr,
sideFEPtr(side_fe_ptr), sideDataPtr(side_data_ptr),
errorSumPtr(error_sum_ptr), thError(th_error) {}
CHKERR loopSideFaces(
"dFE", *sideFEPtr);
const auto in_the_loop =
sideFEPtr->nInTheLoop;
auto not_side = [](auto s) {
};
auto nb_gauss_pts = getGaussPts().size2();
getFTensor2FromMat<DIM1, DIM2>(sideDataPtr->lVec[
LEFT_SIDE]),
getFTensor2FromMat<DIM1, DIM2>(sideDataPtr->lVec[
RIGHT_SIDE]));
getFTensor1FromMat<SPACE_DIM>(sideDataPtr->velMat[
LEFT_SIDE]),
getFTensor1FromMat<SPACE_DIM>(sideDataPtr->velMat[
RIGHT_SIDE]));
auto next = [&]() {
for (auto &t_l : arr_t_l)
++t_l;
for (auto &t_vel : arr_t_vel)
++t_vel;
};
double e = 0;
auto t_w = getFTensor0IntegrationWeight();
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
e += t_w * getMeasure() * t_diff(
I,
J) * t_diff(
I,
J);
next();
++t_w;
}
e = std::sqrt(e);
getNumeredEntFiniteElementPtr()->getBasicDataPtr()->moab;
const void *tags_ptr[2];
CHKERR moab.tag_get_by_ptr(thError, sideDataPtr->feSideHandle.data(), 2,
tags_ptr);
for (auto ff : {0, 1}) {
*((double *)tags_ptr[ff]) += e;
}
CHKERR VecSetValue(errorSumPtr, 0, e, ADD_VALUES);
};
}
private:
boost::shared_ptr<FaceSideEle> sideFEPtr;
boost::shared_ptr<SideData> sideDataPtr;
Tag thError;
};
Tag th_error;
double def_val = 0;
MB_TAG_CREAT | MB_TAG_SPARSE,
&def_val);
auto clear_tags = [&]() {
double zero;
};
auto evaluate_error = [&]() {
auto skel_fe = boost::make_shared<BoundaryEle>(
mField);
skel_fe->getRuleHook = [](
int,
int,
int o) {
return 3 * o; };
auto side_data_ptr = boost::make_shared<SideData>();
skel_fe->getOpPtrVector().push_back(
new OpErrorSkel(side_fe_ptr, side_data_ptr, error_sum_ptr, th_error));
skel_fe->exeTestHook = [&](
FEMethod *fe_ptr) {
};
simple->getSkeletonFEName(), skel_fe);
};
auto assemble_and_sum = [](auto vec) {
double sum;
return sum;
};
auto propagate_error_to_parents = [&]() {
auto fe_ptr = boost::make_shared<FEMethod>();
fe_ptr->exeTestHook = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
fe_ptr->preProcessHook = []() { return 0; };
fe_ptr->postProcessHook = []() { return 0; };
fe_ptr->operatorHook = [&]() {
auto fe_ent = fe_ptr->numeredEntFiniteElementPtr->getEnt();
auto parent = fe_ptr->numeredEntFiniteElementPtr->getParentEnt();
auto th_parent = fe_ptr->numeredEntFiniteElementPtr->getBasicDataPtr()
->th_RefParentHandle;
double error;
CHKERR moab.tag_get_data(th_error, &fe_ent, 1, &error);
[&](auto fe_ent, auto error) {
double *e_ptr;
CHKERR moab.tag_get_by_ptr(th_error, &fe_ent, 1,
(const void **)&e_ptr);
(*e_ptr) += error;
CHKERR moab.tag_get_data(th_parent, &fe_ent, 1, &parent);
if (parent != fe_ent && parent)
CHKERR add_error(parent, *e_ptr);
};
CHKERR add_error(parent, error);
};
fe_ptr);
};
return std::make_tuple(assemble_and_sum(error_sum_ptr), th_error);
}
DivergenceVol(boost::shared_ptr<MatrixDouble> l_ptr,
boost::shared_ptr<MatrixDouble> vel_ptr,
divVec(div_vec) {}
if (nb_dofs) {
const auto nb_gauss_pts = getGaussPts().size2();
const auto t_w = getFTensor0IntegrationWeight();
auto t_l = getFTensor2FromMat<DIM1, DIM2>(*lPtr);
auto t_vel = getFTensor1FromMat<SPACE_DIM>(*velPtr);
double div = 0;
for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
for (int rr = 0; rr != nb_dofs; ++rr) {
div += getMeasure() * t_w * t_l(
I,
I) * (t_diff(
i) * t_vel(
i));
++t_diff;
}
++t_w;
++t_l;
++t_vel;
}
CHKERR VecSetValue(divVec, 0, div, ADD_VALUES);
}
}
private:
boost::shared_ptr<MatrixDouble> lPtr;
boost::shared_ptr<MatrixDouble> velPtr;
};
DivergenceSkeleton(boost::shared_ptr<SideData> side_data_ptr,
boost::shared_ptr<FaceSideEle> side_fe_ptr,
sideDataPtr(side_data_ptr), sideFEPtr(side_fe_ptr), divVec(div_vec) {}
&*base_mat.data().begin());
};
auto not_side = [](auto s) {
};
CHKERR loopSideFaces(
"dFE", *sideFEPtr);
const auto in_the_loop =
sideFEPtr->nInTheLoop;
auto t_normal = getFTensor1Normal();
const auto nb_gauss_pts = getGaussPts().size2();
const auto nb_dofs = sideDataPtr->indicesRowSideMap[s0].size();
if (nb_dofs) {
auto t_base =
get_ntensor(sideDataPtr->rowBaseSideMap[s0]);
auto nb_row_base_functions = sideDataPtr->rowBaseSideMap[s0].size2();
auto side_sense = sideDataPtr->senseMap[s0];
auto opposite_s0 = not_side(s0);
getFTensor2FromMat<DIM1, DIM2>(sideDataPtr->lVec[
LEFT_SIDE]),
getFTensor2FromMat<DIM1, DIM2>(sideDataPtr->lVec[
RIGHT_SIDE]));
getFTensor1FromMat<SPACE_DIM>(sideDataPtr->velMat[
LEFT_SIDE]),
getFTensor1FromMat<SPACE_DIM>(sideDataPtr->velMat[
RIGHT_SIDE]));
auto next = [&]() {
for (auto &t_l : arr_t_l)
++t_l;
for (auto &t_vel : arr_t_vel)
++t_vel;
};
double div = 0;
auto t_w = getFTensor0IntegrationWeight();
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
const auto dot = (t_normal(
i) * t_vel(
i)) * side_sense;
const auto l_upwind_side = (dot > 0) ? s0 : opposite_s0;
const auto l_upwind =
arr_t_l[l_upwind_side];
auto res = t_w * l_upwind(
I,
I) * dot;
++t_w;
next();
int rr = 0;
for (; rr != nb_dofs; ++rr) {
div += t_base * res;
++t_base;
}
for (; rr < nb_row_base_functions; ++rr) {
++t_base;
}
}
CHKERR VecSetValue(divVec, 0, div, ADD_VALUES);
}
if (!in_the_loop)
break;
}
}
private:
boost::shared_ptr<SideData> sideDataPtr;
boost::shared_ptr<FaceSideEle> sideFEPtr;
boost::shared_ptr<MatrixDouble> velPtr;
};
auto vol_fe = boost::make_shared<DomainEle>(
mField);
auto skel_fe = boost::make_shared<BoundaryEle>(
mField);
vol_fe->getRuleHook = [](
int,
int,
int o) {
return 3 * o; };
skel_fe->getRuleHook = [](
int,
int,
int o) {
return 3 * o; };
auto l_ptr = boost::make_shared<MatrixDouble>();
auto vel_ptr = boost::make_shared<MatrixDouble>();
auto side_data_ptr = boost::make_shared<SideData>();
{L2});
vol_fe->getOpPtrVector().push_back(
vol_fe->getOpPtrVector().push_back(
new DivergenceVol(l_ptr, vel_ptr, div_vol_vec));
skel_fe->getOpPtrVector().push_back(
new DivergenceSkeleton(side_data_ptr, side_fe_ptr, div_skel_vec));
[](double x, double y, double) { return x - y; });
[](double x, double y, double) { return x - y; });
vol_fe->exeTestHook = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
skel_fe->exeTestHook = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
auto assemble_and_sum = [](auto vec) {
double sum;
return sum;
};
auto div_vol = assemble_and_sum(div_vol_vec);
auto div_skel = assemble_and_sum(div_skel_vec);
auto eps = std::abs((div_vol - div_skel) / (div_vol + div_skel));
MOFEM_LOG(
"WORLD", Sev::inform) <<
"Testing divergence volume: " << div_vol;
MOFEM_LOG(
"WORLD", Sev::inform) <<
"Testing divergence skeleton: " << div_skel
<<
" relative difference: " <<
eps;
constexpr double eps_err = 1e-6;
"No consistency between skeleton integral and volume integral");
};
auto post_proc = [&](auto dm, auto f_res, auto out_name) {
auto post_proc_fe =
boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(
mField);
auto l_vec = boost::make_shared<VectorDouble>();
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
{{"L", l_vec}},
{},
{}, {})
);
}
post_proc_fe);
post_proc_fe->writeFile(out_name);
};
constexpr
double eps = 1e-4;
auto x =
opt->setRandomFields(
simple->getDM(), {{
"L", {-1, 1}}, {
"V", {-1, 1}}});
auto dot_x = opt->setRandomFields(
simple->getDM(), {{
"L", {-1, 1}}});
auto diff_x = opt->setRandomFields(
simple->getDM(), {{
"L", {-1, 1}}});
auto test_domain_ops = [&](auto fe_name, auto lhs_pipeline,
auto rhs_pipeline) {
auto diff_res = opt->checkCentralFiniteDifference(
simple->getDM(), fe_name, rhs_pipeline, lhs_pipeline, x, dot_x,
CHKERR post_proc(
simple->getDM(), diff_res,
"tangent_op_error.h5m");
}
double fnorm;
CHKERR VecNorm(diff_res, NORM_2, &fnorm);
"Test consistency of tangent matrix %3.4e", fnorm);
constexpr double err = 1e-9;
if (fnorm > err)
"Norm of directional derivative too large err = %3.4e", fnorm);
};
CHKERR test_domain_ops(
simple->getDomainFEName(), pip->getDomainLhsFE(),
pip->getDomainRhsFE());
CHKERR test_domain_ops(
simple->getSkeletonFEName(), pip->getSkeletonLhsFE(),
pip->getSkeletonRhsFE());
};
boost::function<double(double, double, double)> level_fun) {
boost::shared_ptr<FEMethod> lhs_fe = boost::make_shared<DomainEle>(
mField);
boost::shared_ptr<FEMethod> rhs_fe = boost::make_shared<DomainEle>(
mField);
auto swap_fe = [&]() {
lhs_fe.swap(pip->getDomainLhsFE());
rhs_fe.swap(pip->getDomainRhsFE());
};
swap_fe();
pip->setDomainLhsIntegrationRule([](int, int, int o) { return 3 * o; });
pip->setDomainRhsIntegrationRule([](int, int, int o) { return 3 * o; });
remove_mask.flip();
CHKERR prb_mng->removeDofsOnEntities(
"LEVELSET_POJECTION",
"L",
auto test_bit_ele = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
pip->getDomainLhsFE()->exeTestHook = test_bit_ele;
pip->getDomainRhsFE()->exeTestHook = test_bit_ele;
{L2});
{L2});
pip->getOpDomainLhsPipeline().push_back(
new OpMassLL(
"L",
"L"));
pip->getOpDomainRhsPipeline().push_back(
new OpSourceL(
"L", level_fun));
auto ksp = pip->createKSP(sub_dm);
CHKERR KSPSetFromOptions(ksp);
CHKERR VecGhostUpdateBegin(
L, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(
L, INSERT_VALUES, SCATTER_FORWARD);
MOFEM_LOG(
"LevelSet", Sev::inform) <<
"Error indicator " << error;
#ifndef NDEBUG
auto fe_meshset =
std::vector<Tag> tags{th_error};
"PARALLEL=WRITE_PART", &fe_meshset, 1,
&*tags.begin(), tags.size());
#endif
auto post_proc = [&](auto dm, auto out_name, auto th_error) {
auto post_proc_fe =
boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(
mField);
post_proc_fe->setTagsToTransfer({th_error});
post_proc_fe->exeTestHook = test_bit_ele;
auto l_vec = boost::make_shared<VectorDouble>();
auto l_grad_mat = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector(), {L2});
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
{{"L", l_vec}},
{{"GradL", l_grad_mat}},
{}, {})
);
}
post_proc_fe);
post_proc_fe->writeFile(out_name);
};
CHKERR post_proc(sub_dm,
"initial_level_set.h5m", th_error);
swap_fe();
}
boost::function<double(double, double, double)> vel_fun) {
boost::shared_ptr<FEMethod> lhs_fe = boost::make_shared<DomainEle>(
mField);
boost::shared_ptr<FEMethod> rhs_fe = boost::make_shared<DomainEle>(
mField);
auto swap_fe = [&]() {
lhs_fe.swap(pip->getDomainLhsFE());
rhs_fe.swap(pip->getDomainRhsFE());
};
swap_fe();
pip->setDomainLhsIntegrationRule([](int, int, int o) { return 3 * o; });
pip->setDomainRhsIntegrationRule([](int, int, int o) { return 3 * o; });
remove_mask.flip();
CHKERR prb_mng->removeDofsOnEntities(
"VELOCITY_PROJECTION",
"V",
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(0);
};
pip->getDomainLhsFE()->exeTestHook = test_bit;
pip->getDomainRhsFE()->exeTestHook = test_bit;
{potential_velocity_space});
{potential_velocity_space});
pip->getOpDomainLhsPipeline().push_back(
new OpMassVV(
"V",
"V"));
pip->getOpDomainRhsPipeline().push_back(
new OpSourceV(
"V", vel_fun));
auto ksp = pip->createKSP(sub_dm);
CHKERR KSPSetFromOptions(ksp);
CHKERR VecGhostUpdateBegin(
L, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(
L, INSERT_VALUES, SCATTER_FORWARD);
auto post_proc = [&](auto dm, auto out_name) {
auto post_proc_fe =
boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(
mField);
post_proc_fe->exeTestHook = test_bit;
post_proc_fe->getOpPtrVector(), {potential_velocity_space});
auto potential_vec = boost::make_shared<VectorDouble>();
auto velocity_mat = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
{{"VelocityPotential", potential_vec}},
{{"Velocity", velocity_mat}},
{}, {})
);
} else {
"3d case not implemented");
}
post_proc_fe);
post_proc_fe->writeFile(out_name);
};
CHKERR post_proc(sub_dm,
"initial_velocity_potential.h5m");
swap_fe();
}
remove_mask.flip();
remove_mask);
auto add_post_proc_fe = [&]() {
auto post_proc_fe = boost::make_shared<PostProcEle>(
mField);
Tag th_error;
double def_val = 0;
"Error", 1, MB_TYPE_DOUBLE, th_error, MB_TAG_CREAT | MB_TAG_SPARSE,
&def_val);
post_proc_fe->setTagsToTransfer({th_error});
post_proc_fe->exeTestHook = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
auto vel_ptr = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector(), {L2});
auto l_vec = boost::make_shared<VectorDouble>();
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getPostProcMesh(),
post_proc_fe->getMapGaussPts(),
{{"L", l_vec}},
{{"V", vel_ptr}},
{}, {})
);
}
return post_proc_fe;
};
auto post_proc_fe = add_post_proc_fe();
auto set_time_monitor = [&](auto dm, auto ts) {
auto monitor_ptr = boost::make_shared<FEMethod>();
monitor_ptr->preProcessHook = []() { return 0; };
monitor_ptr->operatorHook = []() { return 0; };
monitor_ptr->postProcessHook = [&]() {
if (!post_proc_fe)
"Null pointer for post proc element");
post_proc_fe);
CHKERR post_proc_fe->writeFile(
"level_set_" +
boost::lexical_cast<std::string>(monitor_ptr->ts_step) + ".h5m");
};
boost::shared_ptr<FEMethod> null;
null, null);
return monitor_ptr;
};
auto ts = pip->createTSIM(sub_dm);
auto set_solution = [&](auto ts) {
};
auto monitor_pt = set_time_monitor(sub_dm, ts);
auto ts_pre_step = [](TS ts) {
MOFEM_LOG(
"LevelSet", Sev::inform) <<
"Error indicator " << error;
auto get_norm = [&](auto x) {
double nrm;
CHKERR VecNorm(x, NORM_2, &nrm);
return nrm;
};
auto set_solution = [&](auto ts) {
DM dm;
CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
<< "Problem " << prb_ptr->getName() << " solution vector norm "
<< get_norm(x);
};
auto refine_and_project = [&](auto ts) {
DM dm;
remove_mask
.flip();
->removeDofsOnEntities(
"ADVECTION",
"L",
BitRefLevel().set(),
remove_mask);
};
auto ts_reset_theta = [&](auto ts) {
DM dm;
};
CHKERR refine_and_project(ts);
};
auto ts_post_step = [](TS ts) { return 0; };
CHKERR TSSetPreStep(ts, ts_pre_step);
CHKERR TSSetPostStep(ts, ts_post_step);
}
ParallelComm *pcomm =
meshset_ptr->get_ptr(), 1);
};
auto get_refined_elements_meshset = [&](
auto bit,
auto mask) {
Tag th_error;
"get error handle");
std::vector<double> errors(fe_ents.size());
"get tag data");
auto it = std::max_element(errors.begin(), errors.end());
double max;
MOFEM_LOG(
"LevelSet", Sev::inform) <<
"Max error: " << max;
auto threshold = wp.getThreshold(max);
std::vector<EntityHandle> fe_to_refine;
fe_to_refine.reserve(fe_ents.size());
auto fe_it = fe_ents.begin();
auto error_it = errors.begin();
for (
auto i = 0;
i != fe_ents.size(); ++
i) {
if (*error_it > threshold) {
fe_to_refine.push_back(*fe_it);
}
++fe_it;
++error_it;
}
ents.insert_list(fe_to_refine.begin(), fe_to_refine.end());
auto get_neighbours_by_bridge_vertices = [&](auto &&ents) {
moab::Interface::UNION);
return ents;
};
ents = get_neighbours_by_bridge_vertices(ents);
#ifndef NDEBUG
"add entities to meshset");
(proc_str + "_fe_to_refine.vtk").c_str(), "VTK", "",
meshset_ptr->get_ptr(), 1);
}
#endif
return ents;
};
auto refine_mesh = [&](
auto l,
auto &&fe_to_refine) {
CHKERR bit_mng->getEntitiesByDimAndRefLevel(
fe_to_refine = intersect(level_ents, fe_to_refine);
level_ents = subtract(level_ents, fe_to_refine);
Range fe_to_refine_children;
bit_mng->updateRangeByChildren(fe_to_refine, fe_to_refine_children);
fe_to_refine_children = fe_to_refine_children.subset_by_dimension(
FE_DIM);
level_ents.merge(fe_to_refine_children);
auto fix_neighbour_level = [&](auto ll) {
auto level_ll = level_ents;
level_ll);
CHKERR skin.find_skin(0, level_ll,
false, skin_edges);
for (auto lll = 0; lll <= ll; ++lll) {
CHKERR bit_mng->updateRangeByParent(skin_edges, skin_parents);
}
for (auto lll = 0; lll <= ll - 2; ++lll) {
bad_bit[lll] = true;
}
skin_parents,
FE_DIM,
false, skin_adj_ents, moab::Interface::UNION);
skin_adj_ents);
skin_adj_ents = intersect(skin_adj_ents, level_ents);
if (!skin_adj_ents.empty()) {
level_ents = subtract(level_ents, skin_adj_ents);
Range skin_adj_ents_children;
bit_mng->updateRangeByChildren(skin_adj_ents, skin_adj_ents_children);
level_ents.merge(skin_adj_ents_children);
}
};
level_ents);
level_ents.subset_by_dimension(
FE_DIM), level_ents,
true);
} else {
level_ents.subset_by_dimension(
FE_DIM),
d,
false, level_ents,
moab::Interface::UNION);
}
}
level_ents);
#ifndef NDEBUG
CHKERR bit_mng->writeBitLevelByDim(
(boost::lexical_cast<std::string>(
l) +
"_" + proc_str +
"_ref_mesh.vtk")
.c_str(),
"VTK", "");
#endif
};
auto set_skelton_bit = [&](
auto l) {
CHKERR bit_mng->getEntitiesByDimAndRefLevel(
Range level_edges_parents;
CHKERR bit_mng->updateRangeByParent(level_edges, level_edges_parents);
level_edges_parents = level_edges_parents.subset_by_dimension(
FE_DIM - 1);
CHKERR bit_mng->filterEntitiesByRefLevel(
auto parent_skeleton = intersect(level_edges, level_edges_parents);
auto skeleton = subtract(level_edges, level_edges_parents);
moab::Interface::UNION);
#ifndef NDEBUG
CHKERR bit_mng->writeBitLevel(
(boost::lexical_cast<std::string>(
l) +
"_" + proc_str +
"_skeleton.vtk")
.c_str(),
"VTK", "");
#endif
};
start_mask[s] = true;
CHKERR bit_mng->lambdaBitRefLevel(
true);
level0);
MOFEM_LOG(
"WORLD", Sev::inform) <<
"Process level: " <<
l;
CHKERR refine_mesh(
l + 1, get_refined_elements_meshset(
CHKERR wp.setAggregateBit(*
this,
l + 1);
}
}
auto lhs_fe = boost::make_shared<DomainEle>(
mField);
auto rhs_fe_prj = boost::make_shared<DomainEle>(
mField);
auto rhs_fe_current = boost::make_shared<DomainEle>(
mField);
lhs_fe->getRuleHook = [](
int,
int,
int o) {
return 3 * o; };
rhs_fe_prj->getRuleHook = [](
int,
int,
int o) {
return 3 * o; };
rhs_fe_current->getRuleHook = [](
int,
int,
int o) {
return 3 * o; };
current_ents);
CHKERR bit_mng->getEntitiesByDimAndRefLevel(
CHKERR bit_mng->updateRangeByParent(prj_ents, prj_ents);
}
current_ents = subtract(
current_ents, prj_ents);
auto test_mesh_bit = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
auto test_prj_bit = [&](
FEMethod *fe_ptr) {
return fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel().test(
};
auto test_current_bit = [&](
FEMethod *fe_ptr) {
return current_ents.find(fe_ptr->getFEEntityHandle()) != current_ents.end();
};
lhs_fe->exeTestHook =
test_mesh_bit;
rhs_fe_prj->exeTestHook =
test_prj_bit;
rhs_fe_current->exeTestHook =
test_current_bit;
remove_mask.flip();
CHKERR prb_mng->removeDofsOnEntities(
"DG_PROJECTION",
"L",
BitRefLevel().set(), remove_mask,
nullptr, 0,
true);
{L2});
lhs_fe->getOpPtrVector().push_back(
auto set_prj_from_child = [&](auto rhs_fe_prj) {
rhs_fe_prj->getOpPtrVector(), {L2});
auto l_vec = boost::make_shared<MatrixDouble>();
rhs_fe_prj->getOpPtrVector().push_back(
auto get_parent_this = [&]() {
auto fe_parent_this = boost::make_shared<DomianParentEle>(
mField);
fe_parent_this->getOpPtrVector().push_back(
return fe_parent_this;
};
auto get_parents_fe_ptr = [&](auto this_fe_ptr) {
std::vector<boost::shared_ptr<DomianParentEle>> parents_elems_ptr_vec;
parents_elems_ptr_vec.emplace_back(
boost::make_shared<DomianParentEle>(
mField));
parents_elems_ptr_vec[
l - 1]->getOpPtrVector().push_back(
}
return parents_elems_ptr_vec[0];
};
auto this_fe_ptr = get_parent_this();
auto parent_fe_ptr = get_parents_fe_ptr(this_fe_ptr);
rhs_fe_prj->getOpPtrVector().push_back(
};
auto set_prj_from_parent = [&](auto rhs_fe_current) {
auto l_vec = boost::make_shared<MatrixDouble>();
auto get_parent_this = [&]() {
auto fe_parent_this = boost::make_shared<DomianParentEle>(
mField);
fe_parent_this->getOpPtrVector().push_back(
return fe_parent_this;
};
auto get_parents_fe_ptr = [&](auto this_fe_ptr) {
std::vector<boost::shared_ptr<DomianParentEle>> parents_elems_ptr_vec;
parents_elems_ptr_vec.emplace_back(
boost::make_shared<DomianParentEle>(
mField));
parents_elems_ptr_vec[
l - 1]->getOpPtrVector().push_back(
}
return parents_elems_ptr_vec[0];
};
auto this_fe_ptr = get_parent_this();
auto parent_fe_ptr = get_parents_fe_ptr(this_fe_ptr);
l_vec->resize(
static_cast<DomainEleOp *
>(op_ptr)->getGaussPts().size2(),
false);
l_vec->clear();
return 0;
};
rhs_fe_current->getOpPtrVector().push_back(reset_op_ptr);
rhs_fe_current->getOpPtrVector().push_back(
rhs_fe_current->getOpPtrVector().push_back(
new OpScalarFieldL(
"L", l_vec));
};
set_prj_from_child(rhs_fe_prj);
set_prj_from_parent(rhs_fe_current);
boost::shared_ptr<FEMethod> null_fe;
lhs_fe, null_fe, null_fe);
null_fe, null_fe);
rhs_fe_current, null_fe, null_fe);
CHKERR KSPSetFromOptions(ksp);
CHKERR VecGhostUpdateBegin(
L, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(
L, INSERT_VALUES, SCATTER_FORWARD);
MOFEM_LOG(
"LevelSet", Sev::inform) <<
"Error indicator " << error;
auto post_proc = [&](auto dm, auto out_name, auto th_error) {
auto post_proc_fe =
boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(
mField);
post_proc_fe->setTagsToTransfer({th_error});
post_proc_fe->exeTestHook = test_mesh_bit;
auto l_vec = boost::make_shared<VectorDouble>();
auto l_grad_mat = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector(), {L2});
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
{{"L", l_vec}},
{{"GradL", l_grad_mat}},
{}, {})
);
}
post_proc_fe);
post_proc_fe->writeFile(out_name);
};
CHKERR post_proc(sub_dm,
"dg_projection.h5m", th_error);
}