9#ifndef __UNSATURATD_FLOW_HPP__
10#define __UNSATURATD_FLOW_HPP__
14using PostProcEle = PostProcBrokenMeshInMoab<VolumeElementForcesAndSourcesCore>;
52 const std::string &prefix)
const = 0;
54 virtual MoFEMErrorCode
calK() {
57 "Not implemented how to calculate hydraulic conductivity");
64 "Not implemented how to calculate derivative of hydraulic "
69 virtual MoFEMErrorCode
calC() {
72 "Not implemented how to calculate capacity");
79 "Not implemented how to calculate capacity");
86 "Not implemented how to calculate capacity");
90 virtual MoFEMErrorCode
calSe() {
93 "Not implemented how to calculate capacity");
112 if (
dM != PETSC_NULLPTR) {
114 CHKERRABORT(PETSC_COMM_WORLD,
ierr);
128 int &block_id)
const {
130 for (MaterialsDoubleMap::const_iterator mit =
dMatMap.begin();
132 if (mit->second->tEts.find(ent) != mit->second->tEts.end()) {
133 block_id = mit->first;
138 "Element not found, no material data");
148 boost::function<
double(
const double x,
const double y,
const double z)>
152 typedef map<int, boost::shared_ptr<BcData>>
BcMap;
169 const double x,
const double y,
const double z,
178 if (it->second->eNts.find(ent) != it->second->eNts.end()) {
179 block_id = it->first;
187 value =
bcValueMap.at(block_id)->hookFun(x, y, z);
211 const double y,
const double z,
double &flux) {
219 if (it->second->eNts.find(ent) != it->second->eNts.end()) {
220 block_id = it->first;
228 flux =
bcFluxMap.at(block_id)->hookFun(x, y, z);
251 boost::shared_ptr<MethodForForceScaling> &value_scale)
265 MoFEMErrorCode
doWork(
int side, EntityType type,
266 EntitiesFieldData::EntData &data) {
268 if (data.getFieldData().size() == 0)
273 nF.resize(data.getIndices().size());
276 int nb_gauss_pts = data.getN().size1();
277 for (
int gg = 0; gg < nb_gauss_pts; gg++) {
286 const double beta = w * (value - z);
287 noalias(
nF) += beta * prod(data.getVectorN<3>(gg),
getNormal());
295 CHKERR VecSetValues(f, data.getIndices().size(), &data.getIndices()[0],
317 MoFEMErrorCode
doWork(
int side, EntityType type,
318 EntitiesFieldData::EntData &data) {
320 const int nb_dofs = data.getIndices().size();
323 nF.resize(nb_dofs,
false);
331 boost::shared_ptr<GenericMaterial> &block_data =
cTx.
dMatMap.at(block_id);
334 auto t_n_hdiv = data.getFTensor1N<3>();
347 auto t_base_diff_hdiv = data.getFTensor2DiffN<3, 3>();
348 divVec.resize(nb_dofs,
false);
351 int nb_gauss_pts = data.getN().size1();
352 for (
int gg = 0; gg != nb_gauss_pts; gg++) {
355 v = t_base_diff_hdiv(
i,
i);
359 const double alpha = t_w * vol;
361 block_data->x = t_coords(0);
362 block_data->y = t_coords(1);
363 block_data->z = t_coords(2);
364 CHKERR block_data->calK();
365 const double K = block_data->K;
366 const double scaleZ = block_data->scaleZ;
367 const double z = t_coords(2);
369 noalias(
nF) -= alpha * (t_h - z * scaleZ) *
divVec;
372 for (
int rr = 0; rr != nb_dofs; rr++) {
373 t_nf += alpha * (1 / K) * (t_n_hdiv(
i) * t_flux(
i));
384 &*data.getIndices().begin(), &*
nF.begin(),
405 MoFEMErrorCode
doWork(
int side, EntityType type,
406 EntitiesFieldData::EntData &data) {
408 const int nb_dofs = data.getIndices().size();
412 nF.resize(nb_dofs,
false);
420 boost::shared_ptr<GenericMaterial> &block_data =
cTx.
dMatMap.at(block_id);
432 const double scale = block_data->sCale;
436 int nb_gauss_pts = data.getN().size1();
437 for (
int gg = 0; gg != nb_gauss_pts; gg++) {
438 const double alpha = t_w * vol *
scale;
440 block_data->x = t_coords(0);
441 block_data->y = t_coords(1);
442 block_data->z = t_coords(2);
443 CHKERR block_data->calC();
444 const double C = block_data->C;
446 noalias(
nF) += (alpha * (t_div_flux + C * t_h_t)) * data.getN(gg);
455 CHKERR VecSetValues(f, nb_dofs, &*data.getIndices().begin(), &*
nF.begin(),
467 const std::string flux_name)
488 MoFEMErrorCode
doWork(
int row_side,
int col_side, EntityType row_type,
490 EntitiesFieldData::EntData &row_data,
491 EntitiesFieldData::EntData &col_data) {
493 const int nb_row = row_data.getIndices().size();
494 const int nb_col = col_data.getIndices().size();
499 nN.resize(nb_row, nb_col,
false);
507 boost::shared_ptr<GenericMaterial> &block_data =
cTx.
dMatMap.at(block_id);
513 auto t_n_hdiv_row = row_data.getFTensor1N<3>();
518 int nb_gauss_pts = row_data.getN().size1();
519 for (
int gg = 0; gg != nb_gauss_pts; gg++) {
521 block_data->x = t_coords(0);
522 block_data->y = t_coords(1);
523 block_data->z = t_coords(2);
524 CHKERR block_data->calK();
525 const double K = block_data->K;
527 const double alpha = t_w * vol;
528 const double beta = alpha * (1 / K);
530 for (
int kk = 0; kk != nb_row; kk++) {
531 auto t_n_hdiv_col = col_data.getFTensor1N<3>(gg, 0);
532 for (
int ll = 0; ll != nb_col; ll++) {
533 t_a += beta * (t_n_hdiv_row(
j) * t_n_hdiv_col(
j));
544 CHKERR MatSetValues(
a, nb_row, &*row_data.getIndices().begin(), nb_col,
545 &*col_data.getIndices().begin(), &*
nN.data().begin(),
548 if (row_side != col_side || row_type != col_type) {
549 transNN.resize(nb_col, nb_row);
551 CHKERR MatSetValues(
a, nb_col, &*col_data.getIndices().begin(), nb_row,
552 &*row_data.getIndices().begin(),
553 &*
transNN.data().begin(), ADD_VALUES);
583 MoFEMErrorCode
doWork(
int row_side,
int col_side, EntityType row_type,
585 EntitiesFieldData::EntData &row_data,
586 EntitiesFieldData::EntData &col_data) {
588 int nb_row = row_data.getIndices().size();
589 int nb_col = col_data.getIndices().size();
594 nN.resize(nb_row, nb_col,
false);
602 boost::shared_ptr<GenericMaterial> &block_data =
cTx.
dMatMap.at(block_id);
614 const double scale = block_data->sCale;
619 int nb_gauss_pts = row_data.getN().size1();
621 auto t_n_row = row_data.getFTensor0N();
622 for (
int gg = 0; gg != nb_gauss_pts; gg++) {
624 double alpha = t_w * vol *
scale;
628 block_data->h_t = t_h_t;
629 block_data->x = t_coords(0);
630 block_data->y = t_coords(1);
631 block_data->z = t_coords(2);
632 CHKERR block_data->calC();
633 CHKERR block_data->calDiffC();
634 const double C = block_data->C;
635 const double diffC = block_data->diffC;
638 &*
nN.data().begin());
640 for (
int kk = 0; kk != nb_row; kk++) {
642 auto t_n_col = col_data.getFTensor0N(gg, 0);
644 for (
int ll = 0; ll != nb_col; ll++) {
646 t_a += (alpha * (C * ts_a + diffC * t_h_t)) * t_n_row * t_n_col;
659 CHKERR MatSetValues(
a, nb_row, &row_data.getIndices()[0], nb_col,
660 &col_data.getIndices()[0], &*
nN.data().begin(),
675 const std::string &val_name_row,
676 const std::string &flux_name_col)
696 MoFEMErrorCode
doWork(
int row_side,
int col_side, EntityType row_type,
698 EntitiesFieldData::EntData &row_data,
699 EntitiesFieldData::EntData &col_data) {
701 int nb_row = row_data.getFieldData().size();
702 int nb_col = col_data.getFieldData().size();
713 boost::shared_ptr<GenericMaterial> &block_data =
cTx.
dMatMap.at(block_id);
714 nN.resize(nb_row, nb_col,
false);
715 divVec.resize(nb_col,
false);
719 auto t_base_diff_hdiv = col_data.getFTensor2DiffN<3, 3>();
721 const double scale = block_data->sCale;
722 int nb_gauss_pts = row_data.getN().size1();
723 for (
int gg = 0; gg < nb_gauss_pts; gg++) {
726 v = t_base_diff_hdiv(
i,
i);
729 noalias(
nN) += alpha * outer_prod(row_data.getN(gg),
divVec);
732 &row_data.getIndices()[0], nb_col,
733 &col_data.getIndices()[0], &
nN(0, 0), ADD_VALUES);
747 const std::string &flux_name_col,
748 const std::string &val_name_row)
769 MoFEMErrorCode
doWork(
int row_side,
int col_side, EntityType row_type,
771 EntitiesFieldData::EntData &row_data,
772 EntitiesFieldData::EntData &col_data) {
774 int nb_row = row_data.getFieldData().size();
775 int nb_col = col_data.getFieldData().size();
780 nN.resize(nb_row, nb_col,
false);
781 divVec.resize(nb_row,
false);
789 boost::shared_ptr<GenericMaterial> &block_data =
cTx.
dMatMap.at(block_id);
799 auto t_n_hdiv_row = row_data.getFTensor1N<3>();
802 auto t_base_diff_hdiv = row_data.getFTensor2DiffN<3, 3>();
805 int nb_gauss_pts = row_data.getN().size1();
806 for (
int gg = 0; gg != nb_gauss_pts; gg++) {
808 block_data->x = t_coords(0);
809 block_data->y = t_coords(1);
810 block_data->z = t_coords(2);
811 CHKERR block_data->calK();
812 CHKERR block_data->calDiffK();
813 const double K = block_data->K;
815 const double KK = K * K;
816 const double diffK = block_data->diffK;
817 double alpha = t_w * vol;
819 v = t_base_diff_hdiv(
i,
i);
822 noalias(
nN) -= alpha * outer_prod(
divVec, col_data.getN(gg));
824 for (
int rr = 0; rr != nb_row; rr++) {
825 double beta = alpha * (-diffK / KK) * (t_n_hdiv_row(
i) * t_flux(
i));
826 auto t_n_col = col_data.getFTensor0N(gg, 0);
827 for (
int cc = 0; cc != nb_col; cc++) {
828 t_a += beta * t_n_col;
840 &row_data.getIndices()[0], nb_col,
841 &col_data.getIndices()[0], &
nN(0, 0), ADD_VALUES);
850 const std::string &val_name)
858 MoFEMErrorCode
doWork(
int side, EntityType type,
859 EntitiesFieldData::EntData &data) {
861 if (data.getFieldData().size() == 0)
863 auto nb_dofs = data.getFieldData().size();
864 if (nb_dofs != data.getN().size2()) {
866 "wrong number of dofs");
868 nN.resize(nb_dofs, nb_dofs,
false);
869 nF.resize(nb_dofs,
false);
881 boost::shared_ptr<GenericMaterial> &block_data =
cTx.
dMatMap.at(block_id);
885 for (
auto gg = 0; gg != nb_gauss_pts; gg++) {
892 nN += alpha * outer_prod(data.getN(gg), data.getN(gg));
893 nF += alpha * block_data->initialPcEval() * data.getN(gg);
902 CHKERR VecSetValues(
cTx.
D1, nb_dofs, &*data.getIndices().begin(),
903 &*
nF.begin(), INSERT_VALUES);
917 const std::string fluxes_name)
931 MoFEMErrorCode
doWork(
int side, EntityType type,
932 EntitiesFieldData::EntData &data) {
934 int nb_dofs = data.getFieldData().size();
938 auto t_n_hdiv = data.getFTensor1N<3>();
943 double flux_on_entity = 0;
944 int nb_gauss_pts = data.getN().size1();
945 for (
int gg = 0; gg < nb_gauss_pts; gg++) {
946 auto t_data = data.getFTensor0FieldData();
947 for (
int rr = 0; rr != nb_dofs; rr++) {
948 flux_on_entity -= (0.5 * t_data * t_w) * (t_n_hdiv(
i) * t_normal(
i));
972 moab::Interface &post_proc_mesh,
973 std::vector<EntityHandle> &map_gauss_pts,
979 MoFEMErrorCode
doWork(
int side, EntityType type,
980 EntitiesFieldData::EntData &data) {
982 int nb_dofs = data.getFieldData().size();
992 boost::shared_ptr<GenericMaterial> &block_data =
cTx.
dMatMap.at(block_id);
996 int def_block_id = -1;
998 MB_TAG_CREAT | MB_TAG_SPARSE,
1006 MB_TAG_CREAT | MB_TAG_SPARSE, &zero);
1009 MB_TAG_CREAT | MB_TAG_SPARSE, &zero);
1018 MB_TAG_CREAT | MB_TAG_SPARSE, &zero);
1021 MB_TAG_CREAT | MB_TAG_SPARSE, &zero);
1028 int nb_gauss_pts = data.getN().size1();
1029 for (
int gg = 0; gg < nb_gauss_pts; gg++) {
1030 block_data->h = t_h;
1031 block_data->x = t_coords(0);
1032 block_data->y = t_coords(1);
1033 block_data->z = t_coords(2);
1035 CHKERR block_data->calTheta();
1036 double theta = block_data->tHeta;
1038 CHKERR block_data->calSe();
1040 double Se = block_data->Se;
1043 CHKERR block_data->calK();
1044 double K = block_data->K;
1047 CHKERR block_data->calC();
1048 double C = block_data->C;
1079 boost::shared_ptr<PostProcEle> &post_proc,
1080 boost::shared_ptr<ForcesAndSourcesCore> flux_Integrate,
1081 const int frequency)
1104 PetscPrintf(PETSC_COMM_WORLD,
"Output results %d - %d\n", step,
1108 string(
"out_") + boost::lexical_cast<std::string>(step) +
".h5m");
1127 CHKERR PetscPrintf(PETSC_COMM_WORLD,
"Flux at time %6.4g %6.4g\n",
ts_t,
1136 MoFEMErrorCode
addFields(
const std::string &values,
const std::string &fluxes,
1150 if (it->getName().compare(0, 4,
"SOIL") != 0)
1157 MBTET, values +
"_t");
1173 const std::string &values_name) {
1186 values_name +
"_t");
1191 if (it->getName().compare(0, 4,
"SOIL") != 0)
1194 dMatMap[it->getMeshsetId()]->tEts, MBTET,
"MIX");
1209 if (it->getName().compare(0, 4,
"HEAD") != 0)
1212 bcValueMap[it->getMeshsetId()]->eNts, MBTRI,
"MIX_BCVALUE");
1227 if (it->getName().compare(0, 4,
"FLUX") != 0)
1230 bcFluxMap[it->getMeshsetId()]->eNts, MBTRI,
"MIX_BCFLUX");
1242 BitRefLevel ref_level = BitRefLevel().set(0)) {
1255 CHKERR DMCreate(PETSC_COMM_WORLD, &
dM);
1259 CHKERR DMMoFEMSetIsPartitioned(
dM, PETSC_TRUE);
1263 CHKERR DMMoFEMSetIsPartitioned(
dM, PETSC_TRUE);
1267 CHKERR DMMoFEMAddElement(
dM,
"MIX");
1268 CHKERR DMMoFEMAddElement(
dM,
"MIX_BCFLUX");
1269 CHKERR DMMoFEMAddElement(
dM,
"MIX_BCVALUE");
1275 "MIX",
"FLUXES", zero_flux_ents);
1277 PetscSection section;
1281 CHKERR PetscSectionDestroy(§ion);
1286 boost::shared_ptr<ForcesAndSourcesCore>
1288 boost::shared_ptr<ForcesAndSourcesCore>
1290 boost::shared_ptr<ForcesAndSourcesCore>
1292 boost::shared_ptr<ForcesAndSourcesCore>
1295 boost::shared_ptr<MethodForForceScaling>
1297 boost::shared_ptr<MethodForForceScaling>
1302 boost::shared_ptr<VectorDouble>
1311 int operator()(
int,
int,
int p_data)
const {
return 3 * p_data; }
1332 boost::shared_ptr<ForcesAndSourcesCore>
fePtr;
1337 boost::shared_ptr<ForcesAndSourcesCore> &fe_ptr
1343 CHKERR fePtr->mField.getInterface<VecManager>()->setOtherLocalGhostVector(
1344 fePtr->problemPtr,
"VALUES",
string(
"VALUES") +
"_t",
ROW,
1345 fePtr->ts_u_t, INSERT_VALUES, SCATTER_REVERSE);
1346 switch (
fePtr->ts_ctx) {
1347 case TSMethod::CTX_TSSETIFUNCTION:
1348 CHKERR VecSetOption(
fePtr->ts_F, VEC_IGNORE_NEGATIVE_INDICES,
1366 for (std::vector<int>::iterator it =
cTx.
bcVecIds.begin();
1369 fePtr->problemPtr->getColDofsByPetscGlobalDofIdx(*it)
1371 dof_ptr->getFieldData() = *vit;
1393 boost::shared_ptr<ForcesAndSourcesCore>
fePtr;
1397 boost::shared_ptr<ForcesAndSourcesCore> &fe_ptr
1402 switch (
fePtr->ts_ctx) {
1403 case TSMethod::CTX_TSSETIJACOBIAN: {
1404 CHKERR MatAssemblyBegin(
fePtr->ts_B, MAT_FINAL_ASSEMBLY);
1405 CHKERR MatAssemblyEnd(
fePtr->ts_B, MAT_FINAL_ASSEMBLY);
1412 CHKERR MatAssemblyBegin(
fePtr->ts_B, MAT_FINAL_ASSEMBLY);
1413 CHKERR MatAssemblyEnd(
fePtr->ts_B, MAT_FINAL_ASSEMBLY);
1418 case TSMethod::CTX_TSSETIFUNCTION: {
1446 ForcesAndSourcesCore::RuleHookFun face_rule =
FaceRule()) {
1450 feFaceBc = boost::shared_ptr<ForcesAndSourcesCore>(
1452 feFaceRhs = boost::shared_ptr<ForcesAndSourcesCore>(
1455 new VolumeElementForcesAndSourcesCore(
mField));
1456 feVolLhs = boost::shared_ptr<ForcesAndSourcesCore>(
1457 new VolumeElementForcesAndSourcesCore(
mField));
1458 feVolRhs = boost::shared_ptr<ForcesAndSourcesCore>(
1459 new VolumeElementForcesAndSourcesCore(
mField));
1485 feFaceBc->getOpPtrVector().push_back(
1502 feVolRhs->getOpPtrVector().push_back(
new OpCalculateScalarFieldValues(
1504 feVolRhs->getOpPtrVector().push_back(
1506 feVolRhs->getOpPtrVector().push_back(
1510 feVolRhs->getOpPtrVector().back().opType =
1511 ForcesAndSourcesCore::UserDataOperator::OPROW;
1514 feVolLhs->getOpPtrVector().push_back(
new OpCalculateScalarFieldValues(
1520 feVolLhs->getOpPtrVector().push_back(
1522 feVolLhs->getOpPtrVector().push_back(
1524 feVolLhs->getOpPtrVector().push_back(
1527 feVolLhs->getOpPtrVector().push_back(
1529 feVolLhs->getOpPtrVector().push_back(
1534 boost::shared_ptr<FEMethod> null;
1542 auto get_post_process_ele = [&]() {
1543 auto post_process = boost::make_shared<PostProcEle>(
mField);
1547 auto values_ptr = boost::make_shared<VectorDouble>();
1548 auto grad_ptr = boost::make_shared<MatrixDouble>();
1549 auto flux_ptr = boost::make_shared<MatrixDouble>();
1551 post_process->getOpPtrVector().push_back(
1552 new OpCalculateScalarFieldValues(
"VALUES", values_ptr));
1553 post_process->getOpPtrVector().push_back(
1554 new OpCalculateScalarFieldGradient<3>(
"VALUES", grad_ptr));
1555 post_process->getOpPtrVector().push_back(
1556 new OpCalculateHVecVectorField<3>(
"FLUXES", flux_ptr));
1558 using OpPPMap = OpPostProcMapInMoab<3, 3>;
1560 post_process->getOpPtrVector().push_back(
1562 new OpPPMap(post_process->getPostProcMesh(),
1563 post_process->getMapGaussPts(),
1565 {{
"VALUES", values_ptr}},
1567 {{
"GRAD", grad_ptr}, {
"FLUXES", flux_ptr}},
1573 return post_process;
1576 auto post_process = get_post_process_ele();
1578 post_process->getOpPtrVector().push_back(
1579 new OpValuesAtGaussPts(*
this,
"VALUES"));
1580 post_process->getOpPtrVector().push_back(
1581 new OpPostProcMaterial(*
this, post_process->getPostProcMesh(),
1582 post_process->getMapGaussPts(),
"VALUES"));
1585 boost::shared_ptr<ForcesAndSourcesCore> flux_integrate;
1586 flux_integrate = boost::shared_ptr<ForcesAndSourcesCore>(
1589 flux_integrate->getOpPtrVector().push_back(
1590 new OpIntegrateFluxes(*
this,
"FLUXES"));
1592 PetscOptionsBegin(PETSC_COMM_WORLD,
"",
"Monitor post proc",
"none");
1594 "-how_often_output",
1595 "frequency how often results are dumped on hard disk",
"", frequency,
1599 tsMonitor = boost::shared_ptr<FEMethod>(
1616 CHKERR DMCreateMatrix(dM, &Aij);
1617 CHKERR DMCreateGlobalVector(dM, &D0);
1618 CHKERR VecDuplicate(D0, &D1);
1622 int nb_locals = mField.get_comm_rank() == 0 ? 1 : 0;
1623 int nb_ghosts = mField.get_comm_rank() > 0 ? 1 : 0;
1624 CHKERR VecCreateGhost(PETSC_COMM_WORLD, nb_locals, 1, nb_ghosts, ghosts,
1640 CHKERR VecDestroy(&ghostFlux);
1651 CHKERR VecZeroEntries(D0);
1652 CHKERR VecGhostUpdateBegin(D0, INSERT_VALUES, SCATTER_FORWARD);
1653 CHKERR VecGhostUpdateEnd(D0, INSERT_VALUES, SCATTER_FORWARD);
1657 CHKERR DMoFEMLoopFiniteElements(dM,
"MIX_BCFLUX", feFaceBc);
1658 CHKERR VecAssemblyBegin(D0);
1659 CHKERR VecAssemblyEnd(D0);
1660 CHKERR VecGhostUpdateBegin(D0, INSERT_VALUES, SCATTER_FORWARD);
1661 CHKERR VecGhostUpdateEnd(D0, INSERT_VALUES, SCATTER_FORWARD);
1663 CHKERR VecNorm(D0, NORM_2, &norm2D0);
1665 PetscPrintf(PETSC_COMM_WORLD,
"norm2D0 = %6.4e\n", norm2D0);
1676 CHKERR VecZeroEntries(D1);
1677 CHKERR VecGhostUpdateBegin(D1, INSERT_VALUES, SCATTER_FORWARD);
1678 CHKERR VecGhostUpdateEnd(D1, INSERT_VALUES, SCATTER_FORWARD);
1680 CHKERR DMoFEMLoopFiniteElements(dM,
"MIX", feVolInitialPc);
1682 CHKERR VecAssemblyBegin(D1);
1683 CHKERR VecAssemblyEnd(D1);
1684 CHKERR VecGhostUpdateBegin(D1, INSERT_VALUES, SCATTER_FORWARD);
1685 CHKERR VecGhostUpdateEnd(D1, INSERT_VALUES, SCATTER_FORWARD);
1688 CHKERR VecNorm(D1, NORM_2, &norm2D1);
1690 PetscPrintf(PETSC_COMM_WORLD,
"norm2D1 = %6.4e\n", norm2D1);
1700 if (set_initial_pc) {
1702 CHKERR DMoFEMMeshToLocalVector(dM, D1, INSERT_VALUES, SCATTER_REVERSE);
1706 CHKERR DMoFEMMeshToLocalVector(dM,
D, INSERT_VALUES, SCATTER_FORWARD);
1710 CHKERR TSCreate(PETSC_COMM_WORLD, &ts);
1712 CHKERR TSSetType(ts, TSBEULER);
1715 CHKERR TSSetMaxTime(ts, ftime);
1717 CHKERR TSSetFromOptions(ts);
1720#if PETSC_VERSION_GE(3, 7, 0)
1721 CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);
1725 DMMoFEMGetTsCtx(dM, &
ts_ctx);
1726 CHKERR TSMonitorSet(ts, TsMonitorSet,
ts_ctx, PETSC_NULLPTR);
1733 CHKERR TSGetSNES(ts, &snes);
1735#if PETSC_VERSION_GE(3, 7, 0)
1737 PetscViewerAndFormat *vf;
1738 CHKERR PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,
1739 PETSC_VIEWER_DEFAULT, &vf);
1742 (MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
1743 void *))SNESMonitorFields,
1744 vf, (MoFEMErrorCode(*)(
void **))PetscViewerAndFormatDestroy);
1748 CHKERR SNESMonitorSet(snes,
1749 (MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
1750 void *))SNESMonitorFields,
1758 CHKERR TSGetTime(ts, &ftime);
1759 PetscInt steps, snesfails, rejects, nonlinits, linits;
1760 CHKERR TSGetStepNumber(ts, &steps);
1761 CHKERR TSGetSNESFailures(ts, &snesfails);
1762 CHKERR TSGetStepRejections(ts, &rejects);
1763 CHKERR TSGetSNESIterations(ts, &nonlinits);
1764 CHKERR TSGetKSPIterations(ts, &linits);
1765 PetscPrintf(PETSC_COMM_WORLD,
1766 "steps %d (%d rejected, %d SNES fails), ftime %g, nonlinits "
1768 steps, rejects, snesfails, ftime, nonlinits, linits);
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static PetscErrorCode ierr
void cholesky_solve(const TRIA &L, VEC &x, ublas::lower)
solve system L L^T x = b inplace
size_t cholesky_decompose(const MATRIX &A, TRIA &L)
decompose the symmetric positive definit matrix A into product L L^T.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ L2
field with C-1 continuity
@ HDIV
field with continuous normal traction
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
PostProcBrokenMeshInMoab< VolumeElementForcesAndSourcesCore > PostProcEle
implementation of Data Operators for Forces and Sources
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
constexpr auto field_name
Generic material model for unsaturated water transport.
double tHeta
Water content.
virtual MoFEMErrorCode calC()
virtual MoFEMErrorCode calSe()
virtual void printMatParameters(const int id, const std::string &prefix) const =0
virtual MoFEMErrorCode calTheta()
double diffC
Derivative of capacity [S^2/L^2 * L^2/F ].
static double ePsilon0
Regularization parameter.
virtual MoFEMErrorCode calK()
virtual MoFEMErrorCode calDiffC()
double sCale
Scale time dependent eq.
double K
Hydraulic conductivity [L/s].
Range tEts
Elements with this material.
double C
Capacity [S^2/L^2].
double Se
Effective saturation.
virtual ~GenericMaterial()
static double scaleZ
Scale z direction.
virtual double initialPcEval() const =0
Initialize head.
double h_t
rate of hydraulic head
static double ePsilon1
Regularization parameter.
virtual MoFEMErrorCode calDiffK()
double diffK
Derivative of hydraulic conductivity [L/s * L^2/F].
Evaluate boundary conditions on fluxes.
calculate flux at integration point
Calculate values at integration points.
MoFEM::Interface & mField
VectorDouble valuesAtGaussPts
values at integration points on element
VectorDouble divergenceAtGaussPts
divergence at integration points on element
MatrixDouble fluxesAtGaussPts
fluxes at integration points on element
Class storing information about boundary condition.
boost::function< double(const double x, const double y, const double z)> hookFun
Set integration rule to boundary elements.
int operator()(int p_row, int p_col, int p_data) const
MoFEMErrorCode operator()()
Main operator function executed for each loop iteration.
MonitorPostProc(UnsaturatedFlowElement &ctx, boost::shared_ptr< PostProcEle > &post_proc, boost::shared_ptr< ForcesAndSourcesCore > flux_Integrate, const int frequency)
boost::shared_ptr< PostProcEle > postProc
boost::shared_ptr< ForcesAndSourcesCore > fluxIntegrate
MoFEMErrorCode preProcess()
Pre-processing function executed at loop initialization.
UnsaturatedFlowElement & cTx
MoFEMErrorCode postProcess()
Post-processing function executed at loop completion.
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Do calculations.
FTensor::Index< 'i', 3 > i
OpDivTauU_HdivL2(UnsaturatedFlowElement &ctx, const std::string &flux_name_col, const std::string &val_name_row)
Constructor.
UnsaturatedFlowElement & cTx
UnsaturatedFlowElement & cTx
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
OpEvaluateInitiallHead(UnsaturatedFlowElement &ctx, const std::string &val_name)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Integrate boundary condition.
OpIntegrateFluxes(UnsaturatedFlowElement &ctx, const std::string fluxes_name)
Constructor.
FTensor::Index< 'i', 3 > i
UnsaturatedFlowElement & cTx
OpPostProcMaterial(UnsaturatedFlowElement &ctx, moab::Interface &post_proc_mesh, std::vector< EntityHandle > &map_gauss_pts, const std::string field_name)
moab::Interface & postProcMesh
std::vector< EntityHandle > & mapGaussPts
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
UnsaturatedFlowElement & cTx
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
UnsaturatedFlowElement & cTx
OpResidualFlux(UnsaturatedFlowElement &ctx, const std::string &flux_name)
FTensor::Index< 'i', 3 > i
OpResidualMass(UnsaturatedFlowElement &ctx, const std::string &val_name)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
UnsaturatedFlowElement & cTx
Evaluate boundary condition at the boundary.
OpRhsBcOnValues(UnsaturatedFlowElement &ctx, const std::string fluxes_name, boost::shared_ptr< MethodForForceScaling > &value_scale)
Constructor.
UnsaturatedFlowElement & cTx
boost::shared_ptr< MethodForForceScaling > valueScale
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Integrate boundary condition.
VectorDouble nF
Vector of residuals.
FTensor::Index< 'j', 3 > j
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Assemble matrix.
UnsaturatedFlowElement & cTx
OpTauDotSigma_HdivHdiv(UnsaturatedFlowElement &ctx, const std::string flux_name)
OpVDivSigma_L2Hdiv(UnsaturatedFlowElement &ctx, const std::string &val_name_row, const std::string &flux_name_col)
Constructor.
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Do calculations.
UnsaturatedFlowElement & cTx
OpVU_L2L2(UnsaturatedFlowElement &ctx, const std::string value_name)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Assemble matrix.
UnsaturatedFlowElement & cTx
Set integration rule to volume elements.
int operator()(int, int, int p_data) const
Post proces method for volume element Assemble vectors and matrices and apply essential boundary cond...
MoFEMErrorCode operator()()
postProcessVol(UnsaturatedFlowElement &ctx, boost::shared_ptr< ForcesAndSourcesCore > &fe_ptr)
UnsaturatedFlowElement & cTx
boost::shared_ptr< ForcesAndSourcesCore > fePtr
Pre-peprocessing Set head pressure rate and get inital essential boundary conditions.
boost::shared_ptr< ForcesAndSourcesCore > fePtr
MoFEMErrorCode operator()()
preProcessVol(UnsaturatedFlowElement &ctx, boost::shared_ptr< ForcesAndSourcesCore > &fe_ptr)
UnsaturatedFlowElement & cTx
Implementation of operators, problem and finite elements for unsaturated flow.
MoFEMErrorCode calculateEssentialBc()
Calculate boundary conditions for fluxes.
boost::shared_ptr< ForcesAndSourcesCore > feVolLhs
Assemble tangent matrix.
MoFEMErrorCode getBcOnValues(const EntityHandle ent, const int gg, const double x, const double y, const double z, double &value)
Get value on boundary.
std::map< int, boost::shared_ptr< GenericMaterial > > MaterialsDoubleMap
MoFEMErrorCode getBcOnFluxes(const EntityHandle ent, const double x, const double y, const double z, double &flux)
essential (Neumann) boundary condition (set fluxes)
MoFEMErrorCode setFiniteElements(ForcesAndSourcesCore::RuleHookFun vol_rule=VolRule(), ForcesAndSourcesCore::RuleHookFun face_rule=FaceRule())
Create finite element instances.
boost::shared_ptr< VectorDouble > headRateAtGaussPts
Vector keeps head rate.
boost::shared_ptr< MethodForForceScaling > scaleMethodFlux
Method scaling fluxes.
MoFEMErrorCode buildProblem(Range zero_flux_ents, BitRefLevel ref_level=BitRefLevel().set(0))
Build problem.
Vec D1
Vector with inital head capillary pressure.
map< int, boost::shared_ptr< BcData > > BcMap
MaterialsDoubleMap dMatMap
materials database
boost::shared_ptr< ForcesAndSourcesCore > feFaceRhs
Face element apply natural bc.
Vec ghostFlux
Ghost Vector of integrated fluxes.
MoFEMErrorCode addFields(const std::string &values, const std::string &fluxes, const int order)
add fields
std::vector< int > bcVecIds
EntityHandle lastEvalBcValEnt
boost::shared_ptr< ForcesAndSourcesCore > feVolRhs
Assemble residual vector.
EntityHandle lastEvalBcFluxEnt
UnsaturatedFlowElement(MoFEM::Interface &m_field)
boost::shared_ptr< FEMethod > tsMonitor
BcMap bcValueMap
Store boundary condition for head capillary pressure.
MoFEMErrorCode solveProblem(bool set_initial_pc=true)
solve problem
int lastEvalBcBlockFluxId
~UnsaturatedFlowElement()
MoFEMErrorCode destroyMatrices()
Delete matrices and vector when no longer needed.
boost::shared_ptr< ForcesAndSourcesCore > feVolInitialPc
Calculate inital boundary conditions.
MoFEMErrorCode calculateInitialPc()
Calculate inital pressure head distribution.
DM dM
Discrete manager for unsaturated flow problem.
MoFEMErrorCode createMatrices()
Create vectors and matrices.
boost::shared_ptr< ForcesAndSourcesCore > feFaceBc
Elemnet to calculate essential bc.
MoFEMErrorCode addFiniteElements(const std::string &fluxes_name, const std::string &values_name)
add finite elements
virtual MoFEMErrorCode getMaterial(const EntityHandle ent, int &block_id) const
For given element handle get material block Id.
boost::shared_ptr< MethodForForceScaling > scaleMethodValue
Method scaling values.
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MoFEMErrorCode add_field(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
virtual MPI_Comm & get_comm() const =0
bool sYmm
If true assume that matrix is symmetric structure.
Deprecated interface functions.
Structure for user loop methods on finite elements.
default operator for TRI element
VectorDouble & getNormal()
get triangle normal
auto getFTensor1Normal()
get normal as tensor
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)
auto getFTensor0IntegrationWeight()
Get integration weights.
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
@ OPROW
operator doWork function is executed on FE rows
@ OPROWCOL
operator doWork is executed on FE rows &columns
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Post post-proc data at points from hash maps.
TS ts
PETSc time stepping solver object.
PetscReal ts_t
Current time value.
Vec & ts_F
Reference to residual vector F(t,U,U_t,U_tt)
PetscReal ts_a
Shift parameter for U_t (see PETSc Time Solver documentation)
Mat & ts_B
Reference to preconditioner matrix for ts_A.
BasicMethodsSequence & getPostProcessMonitor()
Get the postProcess to do Monitor object.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
double getVolume() const
element volume (linear geometry)
Force scale operator for reading two columns.