v0.15.5
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Public Attributes | Private Attributes | Static Private Attributes | List of all members
EshelbianPlasticity::SetIntegrationAtFrontFace Struct Reference
Collaboration diagram for EshelbianPlasticity::SetIntegrationAtFrontFace:
[legend]

Classes

struct  Fe
 

Public Types

using FunRule = boost::function< int(int)>
 

Public Member Functions

 SetIntegrationAtFrontFace (boost::shared_ptr< Range > front_nodes, boost::shared_ptr< Range > front_edges)
 
 SetIntegrationAtFrontFace (boost::shared_ptr< Range > front_nodes, boost::shared_ptr< Range > front_edges, FunRule fun_rule)
 
MoFEMErrorCode operator() (ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)
 

Public Attributes

FunRule funRule = [](int p) { return 2 * (p + 1) + 1; }
 

Private Attributes

boost::shared_ptr< RangefrontNodes
 
boost::shared_ptr< RangefrontEdges
 

Static Private Attributes

static std::map< long int, MatrixDouble > mapRefCoords
 

Detailed Description

Definition at line 661 of file EshelbianPlasticity.cpp.

Member Typedef Documentation

◆ FunRule

Constructor & Destructor Documentation

◆ SetIntegrationAtFrontFace() [1/2]

EshelbianPlasticity::SetIntegrationAtFrontFace::SetIntegrationAtFrontFace ( boost::shared_ptr< Range front_nodes,
boost::shared_ptr< Range front_edges 
)
inline

Definition at line 666 of file EshelbianPlasticity.cpp.

668 : frontNodes(front_nodes), frontEdges(front_edges){};

◆ SetIntegrationAtFrontFace() [2/2]

EshelbianPlasticity::SetIntegrationAtFrontFace::SetIntegrationAtFrontFace ( boost::shared_ptr< Range front_nodes,
boost::shared_ptr< Range front_edges,
FunRule  fun_rule 
)
inline

Definition at line 670 of file EshelbianPlasticity.cpp.

673 : frontNodes(front_nodes), frontEdges(front_edges), funRule(fun_rule){};

Member Function Documentation

◆ operator()()

MoFEMErrorCode EshelbianPlasticity::SetIntegrationAtFrontFace::operator() ( ForcesAndSourcesCore fe_raw_ptr,
int  order_row,
int  order_col,
int  order_data 
)
inline
Examples
/home/lk58p/mofem_install/vanilla_dev_release/mofem-cephas/mofem/users_modules/eshelbian_plasticity/src/impl/EshelbianPlasticity.cpp.

Definition at line 675 of file EshelbianPlasticity.cpp.

676 {
678
679 constexpr bool debug = false;
680
681 constexpr int numNodes = 3;
682 constexpr int numEdges = 3;
683 constexpr int refinementLevels = 6;
684
685 auto &m_field = fe_raw_ptr->mField;
686 auto fe_ptr = static_cast<Fe *>(fe_raw_ptr);
687 auto fe_handle = fe_ptr->getFEEntityHandle();
688
689 auto set_base_quadrature = [&]() {
691 int rule = funRule(order_data);
692 if (rule < QUAD_2D_TABLE_SIZE) {
693 if (QUAD_2D_TABLE[rule]->dim != 2) {
694 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
695 }
696 if (QUAD_2D_TABLE[rule]->order < rule) {
697 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
698 "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
699 }
700 const size_t nb_gauss_pts = QUAD_2D_TABLE[rule]->npoints;
701 fe_ptr->gaussPts.resize(3, nb_gauss_pts, false);
702 cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[1], 3,
703 &fe_ptr->gaussPts(0, 0), 1);
704 cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[2], 3,
705 &fe_ptr->gaussPts(1, 0), 1);
706 cblas_dcopy(nb_gauss_pts, QUAD_2D_TABLE[rule]->weights, 1,
707 &fe_ptr->gaussPts(2, 0), 1);
708 auto &data = fe_ptr->dataOnElement[H1];
709 data->dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 3,
710 false);
711 double *shape_ptr =
712 &*data->dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
713 cblas_dcopy(3 * nb_gauss_pts, QUAD_2D_TABLE[rule]->points, 1, shape_ptr,
714 1);
715 data->dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(3, 2, false);
716 std::copy(
718 data->dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).data().begin());
719
720 } else {
721 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
722 "rule > quadrature order %d < %d", rule, QUAD_3D_TABLE_SIZE);
723 }
725 };
726
727 CHKERR set_base_quadrature();
728
730
731 auto get_singular_nodes = [&]() {
732 int num_nodes;
733 const EntityHandle *conn;
734 CHKERR m_field.get_moab().get_connectivity(fe_handle, conn, num_nodes,
735 true);
736 std::bitset<numNodes> singular_nodes;
737 for (auto nn = 0; nn != numNodes; ++nn) {
738 if (frontNodes->find(conn[nn]) != frontNodes->end()) {
739 singular_nodes.set(nn);
740 } else {
741 singular_nodes.reset(nn);
742 }
743 }
744 return singular_nodes;
745 };
746
747 auto get_singular_edges = [&]() {
748 std::bitset<numEdges> singular_edges;
749 for (int ee = 0; ee != numEdges; ee++) {
750 EntityHandle edge;
751 CHKERR m_field.get_moab().side_element(fe_handle, 1, ee, edge);
752 if (frontEdges->find(edge) != frontEdges->end()) {
753 singular_edges.set(ee);
754 } else {
755 singular_edges.reset(ee);
756 }
757 }
758 return singular_edges;
759 };
760
761 auto set_gauss_pts = [&](auto &ref_gauss_pts) {
763 fe_ptr->gaussPts.swap(ref_gauss_pts);
764 const size_t nb_gauss_pts = fe_ptr->gaussPts.size2();
765 auto &data = fe_ptr->dataOnElement[H1];
766 data->dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 4);
767 double *shape_ptr =
768 &*data->dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
769 CHKERR ShapeMBTRI(shape_ptr, &fe_ptr->gaussPts(0, 0),
770 &fe_ptr->gaussPts(1, 0), nb_gauss_pts);
772 };
773
774 auto singular_nodes = get_singular_nodes();
775 if (singular_nodes.count()) {
776 auto it_map_ref_coords = mapRefCoords.find(singular_nodes.to_ulong());
777 if (it_map_ref_coords != mapRefCoords.end()) {
778 CHKERR set_gauss_pts(it_map_ref_coords->second);
780 } else {
781
782 auto refine_quadrature = [&]() {
784
785 const int max_level = refinementLevels;
786
787 moab::Core moab_ref;
788 double base_coords[] = {0, 0, 0, 1, 0, 0, 0, 1, 0};
789 EntityHandle nodes[numNodes];
790 for (int nn = 0; nn != numNodes; nn++)
791 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
792 EntityHandle tri;
793 CHKERR moab_ref.create_element(MBTRI, nodes, numNodes, tri);
794 MoFEM::CoreTmp<-1> mofem_ref_core(moab_ref, PETSC_COMM_SELF, -2);
795 MoFEM::Interface &m_field_ref = mofem_ref_core;
796 {
797 Range tris(tri, tri);
798 Range edges;
799 CHKERR m_field_ref.get_moab().get_adjacencies(
800 tris, 1, true, edges, moab::Interface::UNION);
801 CHKERR m_field_ref.getInterface<BitRefManager>()->setBitRefLevel(
802 tris, BitRefLevel().set(0), false, VERBOSE);
803 }
804
805 Range nodes_at_front;
806 for (int nn = 0; nn != numNodes; nn++) {
807 if (singular_nodes[nn]) {
808 EntityHandle ent;
809 CHKERR moab_ref.side_element(tri, 0, nn, ent);
810 nodes_at_front.insert(ent);
811 }
812 }
813
814 auto singular_edges = get_singular_edges();
815
816 EntityHandle meshset;
817 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
818 for (int ee = 0; ee != numEdges; ee++) {
819 if (singular_edges[ee]) {
820 EntityHandle ent;
821 CHKERR moab_ref.side_element(tri, 1, ee, ent);
822 CHKERR moab_ref.add_entities(meshset, &ent, 1);
823 }
824 }
825
826 // refine mesh
827 auto *m_ref = m_field_ref.getInterface<MeshRefinement>();
828 for (int ll = 0; ll != max_level; ll++) {
829 Range edges;
830 CHKERR m_field_ref.getInterface<BitRefManager>()
831 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
832 BitRefLevel().set(), MBEDGE,
833 edges);
834 Range ref_edges;
835 CHKERR moab_ref.get_adjacencies(
836 nodes_at_front, 1, true, ref_edges, moab::Interface::UNION);
837 ref_edges = intersect(ref_edges, edges);
838 Range ents;
839 CHKERR moab_ref.get_entities_by_type(meshset, MBEDGE, ents, true);
840 ref_edges = intersect(ref_edges, ents);
841 Range tris;
842 CHKERR m_field_ref.getInterface<BitRefManager>()
843 ->getEntitiesByTypeAndRefLevel(
844 BitRefLevel().set(ll), BitRefLevel().set(), MBTRI, tris);
845 CHKERR m_ref->addVerticesInTheMiddleOfEdges(
846 ref_edges, BitRefLevel().set(ll + 1));
847 CHKERR m_ref->refineTris(tris, BitRefLevel().set(ll + 1));
848 CHKERR m_field_ref.getInterface<BitRefManager>()
849 ->updateMeshsetByEntitiesChildren(meshset,
850 BitRefLevel().set(ll + 1),
851 meshset, MBEDGE, true);
852 }
853
854 // get ref coords
855 Range tris;
856 CHKERR m_field_ref.getInterface<BitRefManager>()
857 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(max_level),
858 BitRefLevel().set(), MBTRI,
859 tris);
860
861 if (debug) {
862 CHKERR save_range(moab_ref, "ref_tris.vtk", tris);
863 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "debug");
864 }
865
866 MatrixDouble ref_coords(tris.size(), 9, false);
867 int tt = 0;
868 for (Range::iterator tit = tris.begin(); tit != tris.end();
869 tit++, tt++) {
870 int num_nodes;
871 const EntityHandle *conn;
872 CHKERR moab_ref.get_connectivity(*tit, conn, num_nodes, false);
873 CHKERR moab_ref.get_coords(conn, num_nodes, &ref_coords(tt, 0));
874 }
875
876 auto &data = fe_ptr->dataOnElement[H1];
877 const size_t nb_gauss_pts = fe_ptr->gaussPts.size2();
878 MatrixDouble ref_gauss_pts(3, nb_gauss_pts * ref_coords.size1());
879 MatrixDouble &shape_n =
880 data->dataOnEntities[MBVERTEX][0].getN(NOBASE);
881 int gg = 0;
882 for (size_t tt = 0; tt != ref_coords.size1(); tt++) {
883 double *tri_coords = &ref_coords(tt, 0);
885 CHKERR Tools::getTriNormal(tri_coords, &t_normal(0));
886 auto det = t_normal.l2();
887 for (size_t ggg = 0; ggg != nb_gauss_pts; ++ggg, ++gg) {
888 for (int dd = 0; dd != 2; dd++) {
889 ref_gauss_pts(dd, gg) =
890 shape_n(ggg, 0) * tri_coords[3 * 0 + dd] +
891 shape_n(ggg, 1) * tri_coords[3 * 1 + dd] +
892 shape_n(ggg, 2) * tri_coords[3 * 2 + dd];
893 }
894 ref_gauss_pts(2, gg) = fe_ptr->gaussPts(2, ggg) * det;
895 }
896 }
897
898 mapRefCoords[singular_nodes.to_ulong()].swap(ref_gauss_pts);
899 CHKERR set_gauss_pts(mapRefCoords[singular_nodes.to_ulong()]);
900
902 };
903
904 CHKERR refine_quadrature();
905 }
906 }
907 }
908
910 }
@ VERBOSE
@ NOBASE
Definition definitions.h:59
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ H1
continuous field
Definition definitions.h:85
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ 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.
constexpr int order
static const bool debug
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition fem_tools.c:182
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition ddTensor0.hpp:33
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
#define QUAD_2D_TABLE_SIZE
Definition quad.h:174
#define QUAD_3D_TABLE_SIZE
Definition quad.h:186
static QUAD *const QUAD_2D_TABLE[]
Definition quad.h:175
static PetscBool setSingularity
static std::map< long int, MatrixDouble > mapRefCoords
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
Deprecated interface functions.
Mesh refinement interface.
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition Tools.hpp:104
static MoFEMErrorCode getTriNormal(const double *coords, double *normal, double *d_normal=nullptr)
Get the Tri Normal objectGet triangle normal.
Definition Tools.cpp:353
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
int npoints
Definition quad.h:29
auto save_range

Member Data Documentation

◆ frontEdges

boost::shared_ptr<Range> EshelbianPlasticity::SetIntegrationAtFrontFace::frontEdges
private

◆ frontNodes

boost::shared_ptr<Range> EshelbianPlasticity::SetIntegrationAtFrontFace::frontNodes
private

◆ funRule

FunRule EshelbianPlasticity::SetIntegrationAtFrontFace::funRule = [](int p) { return 2 * (p + 1) + 1; }

◆ mapRefCoords

std::map<long int, MatrixDouble> EshelbianPlasticity::SetIntegrationAtFrontFace::mapRefCoords
inlinestaticprivate

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