v0.15.5
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Protected Attributes | List of all members
EshelbianPlasticity::OpTreeSearch Struct Reference
Inheritance diagram for EshelbianPlasticity::OpTreeSearch:
[legend]
Collaboration diagram for EshelbianPlasticity::OpTreeSearch:
[legend]

Public Types

using UOP = FaceElementForcesAndSourcesCore::UserDataOperator
 

Public Member Functions

 OpTreeSearch (boost::shared_ptr< ContactTree > contact_tree_ptr, boost::shared_ptr< MatrixDouble > u_h1_ptr, boost::shared_ptr< MatrixDouble > traction_ptr, Range r, moab::Interface *post_proc_mesh_ptr, std::vector< EntityHandle > *map_gauss_pts_ptr)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 

Protected Attributes

boost::shared_ptr< ContactTreecontactTreePtr
 
boost::shared_ptr< MatrixDouble > uH1Ptr
 
boost::shared_ptr< MatrixDouble > tractionPtr
 
moab::Interface * postProcMeshPtr = nullptr
 
std::vector< EntityHandle > * mapGaussPtsPtr = nullptr
 
Range contactRange
 

Detailed Description

Definition at line 1434 of file EshelbianContact.cpp.

Member Typedef Documentation

◆ UOP

using EshelbianPlasticity::OpTreeSearch::UOP = FaceElementForcesAndSourcesCore::UserDataOperator

Definition at line 1436 of file EshelbianContact.cpp.

Constructor & Destructor Documentation

◆ OpTreeSearch()

EshelbianPlasticity::OpTreeSearch::OpTreeSearch ( boost::shared_ptr< ContactTree contact_tree_ptr,
boost::shared_ptr< MatrixDouble >  u_h1_ptr,
boost::shared_ptr< MatrixDouble >  traction_ptr,
Range  r,
moab::Interface *  post_proc_mesh_ptr,
std::vector< EntityHandle > *  map_gauss_pts_ptr 
)

Definition at line 1459 of file EshelbianContact.cpp.

1468 : UOP(NOSPACE, UOP::OPSPACE), contactTreePtr(contact_tree_ptr),
1469 uH1Ptr(u_h1_ptr), tractionPtr(traction_ptr),
1470 postProcMeshPtr(post_proc_mesh_ptr), mapGaussPtsPtr(map_gauss_pts_ptr),
1471 contactRange(r) {}
@ NOSPACE
Definition definitions.h:83
boost::shared_ptr< ContactTree > contactTreePtr
std::vector< EntityHandle > * mapGaussPtsPtr
boost::shared_ptr< MatrixDouble > uH1Ptr
boost::shared_ptr< MatrixDouble > tractionPtr
FaceElementForcesAndSourcesCore::UserDataOperator UOP

Member Function Documentation

◆ doWork()

MoFEMErrorCode EshelbianPlasticity::OpTreeSearch::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)

Definition at line 1473 of file EshelbianContact.cpp.

1474 {
1476
1477 auto &m_field = getPtrFE()->mField;
1478 auto fe_ent = getNumeredEntFiniteElementPtr()->getEnt();
1479 auto fe_id = id_from_handle(fe_ent);
1480
1481 if (contactRange.find(fe_ent) == contactRange.end())
1483
1484 FTensor::Index<'i', 3> i;
1485 FTensor::Index<'j', 3> j;
1486
1487 const auto nb_gauss_pts = getGaussPts().size2();
1488
1489 auto t_disp_h1 = getFTensor1FromMat<3>(*uH1Ptr);
1490 auto t_coords = getFTensor1CoordsAtGaussPts();
1491 auto t_traction = getFTensor1FromMat<3>(*tractionPtr);
1492
1493 auto next = [&]() {
1494 ++t_disp_h1;
1495 ++t_traction;
1496 ++t_coords;
1497 };
1498
1499 auto get_ele_centre = [i](auto t_ele_coords) {
1500 FTensor::Tensor1<double, 3> t_ele_center;
1501 t_ele_center(i) = 0;
1502 for (int nn = 0; nn != 3; nn++) {
1503 t_ele_center(i) += t_ele_coords(i);
1504 ++t_ele_coords;
1505 }
1506 t_ele_center(i) /= 3;
1507 return t_ele_center;
1508 };
1509
1510 auto get_ele_radius = [i](auto t_ele_center, auto t_ele_coords) {
1512 t_n0(i) = t_ele_center(i) - t_ele_coords(i);
1513 return t_n0.l2();
1514 };
1515
1516 auto get_face_conn = [this](auto face) {
1517 const EntityHandle *conn;
1518 int num_nodes;
1519 CHK_MOAB_THROW(contactTreePtr->getPostProcMesh().get_connectivity(
1520 face, conn, num_nodes, true),
1521 "get conn");
1522 if (num_nodes != 3) {
1523 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "face is not a triangle");
1524 }
1525 return conn;
1526 };
1527
1528 auto get_face_coords = [this](auto conn) {
1529 std::array<double, 9> coords;
1530 CHKERR contactTreePtr->getPostProcMesh().get_coords(conn, 3, coords.data());
1531 return coords;
1532 };
1533
1534 auto get_closet_face = [this](auto *point_ptr, auto r) {
1535 FTensor::Tensor1<double, 3> t_point_out;
1536 std::vector<EntityHandle> faces_out;
1538 contactTreePtr->getTreeSurfPtr()->sphere_intersect_triangles(
1539 point_ptr, r / 8, contactTreePtr->getRootSetSurf(), faces_out),
1540 "get closest faces");
1541 return faces_out;
1542 };
1543
1544 auto get_faces_out = [this](auto *point_ptr, auto *unit_ray_ptr, auto radius,
1545 auto eps) {
1546 std::vector<double> distances_out;
1547 std::vector<EntityHandle> faces_out;
1549
1550 contactTreePtr->getTreeSurfPtr()->ray_intersect_triangles(
1551 distances_out, faces_out, contactTreePtr->getRootSetSurf(), eps,
1552 point_ptr, unit_ray_ptr, &radius),
1553
1554 "get closest faces");
1555 return std::make_pair(faces_out, distances_out);
1556 };
1557
1558 auto get_normal = [](auto &ele_coords) {
1560 Tools::getTriNormal(ele_coords.data(), &t_normal(0));
1561 return t_normal;
1562 };
1563
1564 auto make_map = [&](auto &face_out, auto &face_dist, auto &t_ray_point,
1565 auto &t_unit_ray, auto &t_master_coord) {
1566 FTensor::Index<'i', 3> i;
1567 FTensor::Index<'j', 3> j;
1568 std::map<double, EntityHandle> m;
1569 for (auto ii = 0; ii != face_out.size(); ++ii) {
1570 auto face_conn = get_face_conn(face_out[ii]);
1571 auto face_coords = get_face_coords(face_conn);
1572 auto t_face_normal = get_normal(face_coords);
1573 t_face_normal.normalize();
1575 t_x(i) = t_ray_point(i) + t_unit_ray(i) * face_dist[ii];
1577 t_chi(i, j) =
1578 t_x(i) * t_face_normal(j) - t_master_coord(i) * t_unit_ray(j);
1579 if (t_unit_ray(i) * t_face_normal(i) > std::cos(M_PI / 3)) {
1580 auto dot = std::sqrt(t_chi(i, j) * t_chi(i, j));
1581 m[dot] = face_out[ii];
1582 }
1583 }
1584 return m;
1585 };
1586
1587 auto get_tag_data = [this](auto tag, auto face, auto &vec) {
1589 int tag_length;
1590 CHKERR contactTreePtr->getPostProcMesh().tag_get_length(tag, tag_length);
1591 vec.resize(tag_length);
1592 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(tag, &face, 1,
1593 &*vec.begin());
1595 };
1596
1597 auto create_tag = [this](const std::string tag_name, const int size) {
1598 double def_VAL[] = {0, 0, 0, 0, 0, 0, 0, 0, 0};
1599 Tag th;
1600 if (postProcMeshPtr) {
1601 CHKERR postProcMeshPtr->tag_get_handle(
1602 tag_name.c_str(), size, MB_TYPE_DOUBLE, th,
1603 MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
1604 CHKERR postProcMeshPtr->tag_clear_data(th, &*mapGaussPtsPtr->begin(),
1605 mapGaussPtsPtr->size(), def_VAL);
1606 }
1607 return th;
1608 };
1609
1610 auto set_float_precision = [](const double x) {
1611 if (std::abs(x) < std::numeric_limits<float>::epsilon())
1612 return 0.;
1613 else
1614 return x;
1615 };
1616
1617 // scalars
1618 auto save_scal_tag = [&](auto &th, auto v, const int gg) {
1620 if (postProcMeshPtr) {
1621 v = set_float_precision(v);
1622 CHKERR postProcMeshPtr->tag_set_data(th, &(*mapGaussPtsPtr)[gg], 1, &v);
1623 }
1625 };
1626
1627 // adjacencies
1628 auto get_fe_adjacencies = [this](auto fe_ent) {
1629 Range adj_faces;
1630 CHK_MOAB_THROW(getFaceFE()->mField.get_moab().get_adjacencies(
1631 &fe_ent, 1, 2, false, adj_faces, moab::Interface::UNION),
1632 "get adj");
1633 std::set<int> adj_ids;
1634 for (auto f : adj_faces) {
1635 adj_ids.insert(id_from_handle(f));
1636 }
1637 return adj_ids;
1638 };
1639
1640 auto get_face_id = [this](auto face) {
1641 int id;
1642 if (contactTreePtr->getPostProcMesh().tag_get_data(
1643 contactTreePtr->thEleId, &face, 1, &id) == MB_SUCCESS) {
1644 return id;
1645 }
1646 return -1;
1647 };
1648
1649 auto get_body_id = [this](auto face) {
1650 int id;
1651 if (contactTreePtr->getPostProcMesh().tag_get_data(
1652 contactTreePtr->thBodyId, &face, 1, &id) == MB_SUCCESS) {
1653 return id;
1654 }
1655 return -1;
1656 };
1657
1658 auto get_face_part = [this](auto face) {
1659 ParallelComm *pcomm_post_proc_mesh = ParallelComm::get_pcomm(
1660 &(contactTreePtr->getPostProcMesh()), MYPCOMM_INDEX);
1661 int part;
1662 if (contactTreePtr->getPostProcMesh().tag_get_data(
1663 pcomm_post_proc_mesh->part_tag(), &face, 1, &part) == MB_SUCCESS) {
1664 return part;
1665 }
1666 return -1;
1667 };
1668
1669 auto check_face = [&](auto face, auto fe_id, auto part) {
1670 auto face_id = get_face_id(face);
1671 auto face_part = get_face_part(face);
1672 if (face_id == fe_id && face_part == part)
1673 return true;
1674 return false;
1675 };
1676
1677 // vectors
1678 VectorDouble3 v(3);
1679 FTensor::Tensor1<FTensor::PackPtr<double *, 0>, 3> t_v(&v[0], &v[1], &v[2]);
1680 auto save_vec_tag = [&](auto &th, auto &t_d, const int gg) {
1682 if (postProcMeshPtr) {
1683 t_v(i) = t_d(i);
1684 for (auto &a : v.data())
1685 a = set_float_precision(a);
1686 CHKERR postProcMeshPtr->tag_set_data(th, &(*mapGaussPtsPtr)[gg], 1,
1687 &*v.data().begin());
1688 }
1690 };
1691
1692 Tag th_mark = create_tag("contact_mark", 1);
1693 Tag th_mark_slave = create_tag("contact_mark_slave", 1);
1694 Tag th_body_id = create_tag("contact_body_id", 1);
1695 Tag th_gap = create_tag("contact_gap", 1);
1696 Tag th_tn_master = create_tag("contact_tn_master", 1);
1697 Tag th_tn_slave = create_tag("contact_tn_slave", 1);
1698 Tag th_contact_traction = create_tag("contact_traction", 3);
1699 Tag th_contact_traction_master = create_tag("contact_traction_master", 3);
1700 Tag th_contact_traction_slave = create_tag("contact_traction_slave", 3);
1701 Tag th_c = create_tag("contact_c", 1);
1702 Tag th_normal = create_tag("contact_normal", 3);
1703 Tag th_dist = create_tag("contact_dip", 3);
1704
1705 auto t_ele_centre = get_ele_centre(getFTensor1Coords());
1706 auto ele_radius = get_ele_radius(t_ele_centre, getFTensor1Coords());
1707
1708 contactTreePtr->shadowDataMap.clear();
1709 auto &shadow_vec = contactTreePtr->shadowDataMap[fe_ent];
1710 shadow_vec.clear();
1711
1712 auto adj_fe_ids = get_fe_adjacencies(fe_ent);
1713
1714 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1715
1716 FTensor::Tensor1<double, 3> t_spatial_coords;
1717 t_spatial_coords(i) = t_coords(i) + t_disp_h1(i);
1718
1719 if (postProcMeshPtr) {
1720 CHKERR save_vec_tag(th_contact_traction, t_traction, gg);
1721 }
1722
1723 auto faces_close = get_closet_face(&t_spatial_coords(0), ele_radius);
1724 for (auto face_close : faces_close) {
1725 if (check_face(face_close, fe_id, m_field.get_comm_rank())) {
1726
1727 auto body_id = get_body_id(face_close);
1728
1729 auto master_face_conn = get_face_conn(face_close);
1730 std::array<double, 9> master_coords;
1731 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1732 contactTreePtr->thSmallX, master_face_conn, 3,
1733 master_coords.data());
1734 std::array<double, 9> master_traction;
1735 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1736 contactTreePtr->thTraction, master_face_conn, 3,
1737 master_traction.data());
1738 auto t_normal_face_close = get_normal(master_coords);
1739 t_normal_face_close.normalize();
1740
1741 if (postProcMeshPtr) {
1742 double m = 1;
1743 CHKERR save_scal_tag(th_mark, m, gg);
1744 CHKERR save_scal_tag(th_body_id, static_cast<double>(body_id), gg);
1745 CHKERR save_vec_tag(th_normal, t_normal_face_close, gg);
1746 CHKERR save_vec_tag(th_contact_traction, t_traction, gg);
1747 }
1748
1749 FTensor::Tensor1<double, 3> t_unit_ray;
1750 t_unit_ray(i) = -t_normal_face_close(i);
1751 FTensor::Tensor1<double, 3> t_ray_point;
1752 t_ray_point(i) =
1753 t_spatial_coords(i) -
1754 t_unit_ray(i) * ContactOps::airplane_ray_distance * ele_radius;
1755
1756 constexpr double eps = 1e-3;
1757 auto [faces_out, faces_dist] =
1758 get_faces_out(&t_ray_point(0), &t_unit_ray(0),
1759 2 * ContactOps::airplane_ray_distance * ele_radius,
1760 eps * ele_radius);
1761
1762 auto m = make_map(faces_out, faces_dist, t_ray_point, t_unit_ray,
1763 t_spatial_coords);
1764 for (auto m_it = m.begin(); m_it != m.end(); ++m_it) {
1765 auto face = m_it->second;
1766 if (face != face_close) {
1767
1768 if (
1769
1770 (adj_fe_ids.find(get_face_id(face)) == adj_fe_ids.end() ||
1771 get_face_part(face) != m_field.get_comm_rank())
1772
1773 ) {
1774
1775 shadow_vec.push_back(ContactTree::FaceData());
1776 shadow_vec.back().gaussPtNb = gg;
1777
1778 auto slave_face_conn = get_face_conn(face);
1779 std::array<double, 9> slave_coords;
1780 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1781 contactTreePtr->thSmallX, slave_face_conn, 3,
1782 slave_coords.data());
1783 auto t_normal_face = get_normal(slave_coords);
1784 std::array<double, 9> slave_tractions;
1785 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1786 contactTreePtr->thTraction, slave_face_conn, 3,
1787 slave_tractions.data());
1788
1789 auto t_master_point =
1790 getFTensor1FromPtr<3>(shadow_vec.back().masterPoint.data());
1791 auto t_slave_point =
1792 getFTensor1FromPtr<3>(shadow_vec.back().slavePoint.data());
1793 auto t_ray_point_data =
1794 getFTensor1FromPtr<3>(shadow_vec.back().rayPoint.data());
1795 auto t_unit_ray_data =
1796 getFTensor1FromPtr<3>(shadow_vec.back().unitRay.data());
1797
1798 t_slave_point(i) = t_ray_point(i) + m_it->first * t_unit_ray(i);
1799
1800 auto eval_position = [&](auto &&t_elem_coords, auto &&t_point) {
1801 std::array<double, 2> loc_coords;
1804 &t_elem_coords(0, 0), &t_point(0), 1,
1805 loc_coords.data()),
1806 "get local coords");
1807 FTensor::Tensor1<double, 3> t_shape_fun;
1808 CHK_THROW_MESSAGE(Tools::shapeFunMBTRI<0>(&t_shape_fun(0),
1809 &loc_coords[0],
1810 &loc_coords[1], 1),
1811 "calc shape fun");
1812 FTensor::Index<'i', 3> i;
1813 FTensor::Index<'j', 3> j;
1814 FTensor::Tensor1<double, 3> t_point_out;
1815 t_point_out(i) = t_shape_fun(j) * t_elem_coords(j, i);
1816 return t_point_out;
1817 };
1818
1819 auto t_master_point_updated = eval_position(
1820 getFTensor2FromPtr<3, 3>(master_coords.data()),
1821 getFTensor1FromPtr<3>(shadow_vec.back().masterPoint.data()));
1822 t_master_point(i) = t_master_point_updated(i);
1823
1824 auto t_slave_point_updated = eval_position(
1825 getFTensor2FromPtr<3, 3>(slave_coords.data()),
1826 getFTensor1FromPtr<3>(shadow_vec.back().slavePoint.data()));
1827 t_slave_point(i) = t_slave_point_updated(i);
1828
1829 t_ray_point_data(i) = t_ray_point(i);
1830 t_unit_ray_data(i) = t_unit_ray(i);
1831
1832 std::copy(master_coords.begin(), master_coords.end(),
1833 shadow_vec.back().masterPointNodes.begin());
1834 std::copy(master_traction.begin(), master_traction.end(),
1835 shadow_vec.back().masterTractionNodes.begin());
1836 std::copy(slave_coords.begin(), slave_coords.end(),
1837 shadow_vec.back().slavePointNodes.begin());
1838 std::copy(slave_tractions.begin(), slave_tractions.end(),
1839 shadow_vec.back().slaveTractionNodes.begin());
1840
1841 shadow_vec.back().eleRadius = ele_radius;
1842
1843 // CHKERR get_tag_data(contactTreePtr->thIds, face,
1844 // shadow_vec.back().dofsSlaveIds);
1845 // CHKERR get_tag_data(contactTreePtr->thCoeff, face,
1846 // shadow_vec.back().dofsSlaveCoeff);
1847 // CHKERR get_tag_data(contactTreePtr->thBases, face,
1848 // shadow_vec.back().baseSlaveFuncs);
1849
1850 if (postProcMeshPtr) {
1851 auto [gap, tn_master, tn_slave, c, t_master_traction,
1852 t_slave_traction] =
1853 multiGetGap(&(shadow_vec.back()), t_spatial_coords);
1855 t_gap_vec(i) = t_slave_point(i) - t_spatial_coords(i);
1856 CHKERR save_scal_tag(th_gap, gap, gg);
1857 CHKERR save_scal_tag(th_tn_master, tn_master, gg);
1858 CHKERR save_scal_tag(th_tn_slave, tn_slave, gg);
1859 CHKERR save_scal_tag(th_c, c, gg);
1860 double m = 1;
1861 CHKERR save_scal_tag(th_mark_slave, m, gg);
1862 CHKERR save_vec_tag(th_dist, t_gap_vec, gg);
1863 CHKERR save_vec_tag(th_contact_traction_master,
1864 t_master_traction, gg);
1865 CHKERR save_vec_tag(th_contact_traction_slave, t_slave_traction,
1866 gg);
1867 }
1868
1869 break;
1870 }
1871 }
1872 }
1873 break;
1874 }
1875 }
1876 next();
1877 }
1878
1880}
constexpr double a
static const double eps
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MYPCOMM_INDEX
default communicator number PCOMM
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
static const double face_coords[4][9]
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'j', 3 > j
double airplane_ray_distance
auto multiGetGap(ContactTree::FaceData *face_data_ptr, FTensor::Tensor1< T1, 3 > &t_spatial_coords)
VectorBoundedArray< double, 3 > VectorDouble3
Definition Types.hpp:92
auto id_from_handle(const EntityHandle h)
int r
Definition sdf.py:205
FTensor::Index< 'm', 3 > m
static MoFEMErrorCode getLocalCoordinatesOnReferenceThreeNodeTri(const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords)
Get the local coordinates on reference three node tri object.
Definition Tools.cpp:188
static MoFEMErrorCode getTriNormal(const double *coords, double *normal, double *d_normal=nullptr)
Get the Tri Normal objectGet triangle normal.
Definition Tools.cpp:353

Member Data Documentation

◆ contactRange

Range EshelbianPlasticity::OpTreeSearch::contactRange
protected

Definition at line 1456 of file EshelbianContact.cpp.

◆ contactTreePtr

boost::shared_ptr<ContactTree> EshelbianPlasticity::OpTreeSearch::contactTreePtr
protected

Definition at line 1450 of file EshelbianContact.cpp.

◆ mapGaussPtsPtr

std::vector<EntityHandle>* EshelbianPlasticity::OpTreeSearch::mapGaussPtsPtr = nullptr
protected

Definition at line 1454 of file EshelbianContact.cpp.

◆ postProcMeshPtr

moab::Interface* EshelbianPlasticity::OpTreeSearch::postProcMeshPtr = nullptr
protected

Definition at line 1453 of file EshelbianContact.cpp.

◆ tractionPtr

boost::shared_ptr<MatrixDouble> EshelbianPlasticity::OpTreeSearch::tractionPtr
protected

Definition at line 1452 of file EshelbianContact.cpp.

◆ uH1Ptr

boost::shared_ptr<MatrixDouble> EshelbianPlasticity::OpTreeSearch::uH1Ptr
protected

Definition at line 1451 of file EshelbianContact.cpp.


The documentation for this struct was generated from the following file: