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::SetIntegrationAtFrontVolume Struct Reference
Collaboration diagram for EshelbianPlasticity::SetIntegrationAtFrontVolume:
[legend]

Classes

struct  Fe
 

Public Types

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

Public Member Functions

 SetIntegrationAtFrontVolume (boost::shared_ptr< Range > front_nodes, boost::shared_ptr< Range > front_edges, boost::shared_ptr< CGGUserPolynomialBase::CachePhi > cache_phi=nullptr)
 
 SetIntegrationAtFrontVolume (boost::shared_ptr< Range > front_nodes, boost::shared_ptr< Range > front_edges, FunRule fun_rule, boost::shared_ptr< CGGUserPolynomialBase::CachePhi > cache_phi=nullptr)
 
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; }
 

Private Attributes

boost::shared_ptr< RangefrontNodes
 
boost::shared_ptr< RangefrontEdges
 
boost::shared_ptr< CGGUserPolynomialBase::CachePhicachePhi
 

Static Private Attributes

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

Detailed Description

Definition at line 363 of file EshelbianPlasticity.cpp.

Member Typedef Documentation

◆ FunRule

Constructor & Destructor Documentation

◆ SetIntegrationAtFrontVolume() [1/2]

EshelbianPlasticity::SetIntegrationAtFrontVolume::SetIntegrationAtFrontVolume ( boost::shared_ptr< Range front_nodes,
boost::shared_ptr< Range front_edges,
boost::shared_ptr< CGGUserPolynomialBase::CachePhi cache_phi = nullptr 
)
inline

Definition at line 368 of file EshelbianPlasticity.cpp.

372 : frontNodes(front_nodes), frontEdges(front_edges), cachePhi(cache_phi){};
boost::shared_ptr< CGGUserPolynomialBase::CachePhi > cachePhi

◆ SetIntegrationAtFrontVolume() [2/2]

EshelbianPlasticity::SetIntegrationAtFrontVolume::SetIntegrationAtFrontVolume ( boost::shared_ptr< Range front_nodes,
boost::shared_ptr< Range front_edges,
FunRule  fun_rule,
boost::shared_ptr< CGGUserPolynomialBase::CachePhi cache_phi = nullptr 
)
inline

Definition at line 374 of file EshelbianPlasticity.cpp.

378 : frontNodes(front_nodes), frontEdges(front_edges), funRule(fun_rule),
379 cachePhi(cache_phi){};

Member Function Documentation

◆ operator()()

MoFEMErrorCode EshelbianPlasticity::SetIntegrationAtFrontVolume::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 381 of file EshelbianPlasticity.cpp.

382 {
384
385 constexpr bool debug = false;
386
387 constexpr int numNodes = 4;
388 constexpr int numEdges = 6;
389 constexpr int refinementLevels = 6;
390
391 auto &m_field = fe_raw_ptr->mField;
392 auto fe_ptr = static_cast<Fe *>(fe_raw_ptr);
393 auto fe_handle = fe_ptr->getFEEntityHandle();
394
395 auto set_base_quadrature = [&]() {
397 int rule = funRule(order_data);
398 if (rule < QUAD_3D_TABLE_SIZE) {
399 if (QUAD_3D_TABLE[rule]->dim != 3) {
400 SETERRQ(m_field.get_comm(), MOFEM_DATA_INCONSISTENCY,
401 "wrong dimension");
402 }
403 if (QUAD_3D_TABLE[rule]->order < rule) {
404 SETERRQ(m_field.get_comm(), MOFEM_DATA_INCONSISTENCY,
405 "wrong order %d != %d", QUAD_3D_TABLE[rule]->order, rule);
406 }
407 const size_t nb_gauss_pts = QUAD_3D_TABLE[rule]->npoints;
408 auto &gauss_pts = fe_ptr->gaussPts;
409 gauss_pts.resize(4, nb_gauss_pts, false);
410 cblas_dcopy(nb_gauss_pts, &QUAD_3D_TABLE[rule]->points[1], 4,
411 &gauss_pts(0, 0), 1);
412 cblas_dcopy(nb_gauss_pts, &QUAD_3D_TABLE[rule]->points[2], 4,
413 &gauss_pts(1, 0), 1);
414 cblas_dcopy(nb_gauss_pts, &QUAD_3D_TABLE[rule]->points[3], 4,
415 &gauss_pts(2, 0), 1);
416 cblas_dcopy(nb_gauss_pts, QUAD_3D_TABLE[rule]->weights, 1,
417 &gauss_pts(3, 0), 1);
418 auto &data = fe_ptr->dataOnElement[H1];
419 data->dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 4,
420 false);
421 double *shape_ptr =
422 &*data->dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
423 cblas_dcopy(4 * nb_gauss_pts, QUAD_3D_TABLE[rule]->points, 1, shape_ptr,
424 1);
425 } else {
426 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
427 "rule > quadrature order %d < %d", rule, QUAD_3D_TABLE_SIZE);
428 }
430 };
431
432 CHKERR set_base_quadrature();
433
435
436 auto get_singular_nodes = [&]() {
437 int num_nodes;
438 const EntityHandle *conn;
439 CHKERR m_field.get_moab().get_connectivity(fe_handle, conn, num_nodes,
440 true);
441 std::bitset<numNodes> singular_nodes;
442 for (auto nn = 0; nn != numNodes; ++nn) {
443 if (frontNodes->find(conn[nn]) != frontNodes->end()) {
444 singular_nodes.set(nn);
445 } else {
446 singular_nodes.reset(nn);
447 }
448 }
449 return singular_nodes;
450 };
451
452 auto get_singular_edges = [&]() {
453 std::bitset<numEdges> singular_edges;
454 for (int ee = 0; ee != numEdges; ee++) {
455 EntityHandle edge;
456 CHKERR m_field.get_moab().side_element(fe_handle, 1, ee, edge);
457 if (frontEdges->find(edge) != frontEdges->end()) {
458 singular_edges.set(ee);
459 } else {
460 singular_edges.reset(ee);
461 }
462 }
463 return singular_edges;
464 };
465
466 auto set_gauss_pts = [&](auto &ref_gauss_pts) {
468 fe_ptr->gaussPts.swap(ref_gauss_pts);
469 const size_t nb_gauss_pts = fe_ptr->gaussPts.size2();
470 auto &data = fe_ptr->dataOnElement[H1];
471 data->dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 4);
472 double *shape_ptr =
473 &*data->dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
474 CHKERR ShapeMBTET(shape_ptr, &fe_ptr->gaussPts(0, 0),
475 &fe_ptr->gaussPts(1, 0), &fe_ptr->gaussPts(2, 0),
476 nb_gauss_pts);
478 };
479
480 auto singular_nodes = get_singular_nodes();
481 if (singular_nodes.count()) {
482 auto it_map_ref_coords = mapRefCoords.find(singular_nodes.to_ulong());
483 if (it_map_ref_coords != mapRefCoords.end()) {
484 CHKERR set_gauss_pts(it_map_ref_coords->second);
486 } else {
487
488 auto refine_quadrature = [&]() {
490
491 const int max_level = refinementLevels;
492 EntityHandle tet;
493
494 moab::Core moab_ref;
495 double base_coords[] = {0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1};
496 EntityHandle nodes[4];
497 for (int nn = 0; nn != 4; nn++)
498 CHKERR moab_ref.create_vertex(&base_coords[3 * nn], nodes[nn]);
499 CHKERR moab_ref.create_element(MBTET, nodes, 4, tet);
500 MoFEM::CoreTmp<-1> mofem_ref_core(moab_ref, PETSC_COMM_SELF, -2);
501 MoFEM::Interface &m_field_ref = mofem_ref_core;
502 {
503 Range tets(tet, tet);
504 Range edges;
505 CHKERR m_field_ref.get_moab().get_adjacencies(
506 tets, 1, true, edges, moab::Interface::UNION);
507 CHKERR m_field_ref.getInterface<BitRefManager>()->setBitRefLevel(
508 tets, BitRefLevel().set(0), false, VERBOSE);
509 }
510
511 Range nodes_at_front;
512 for (int nn = 0; nn != numNodes; nn++) {
513 if (singular_nodes[nn]) {
514 EntityHandle ent;
515 CHKERR moab_ref.side_element(tet, 0, nn, ent);
516 nodes_at_front.insert(ent);
517 }
518 }
519
520 auto singular_edges = get_singular_edges();
521
522 EntityHandle meshset;
523 CHKERR moab_ref.create_meshset(MESHSET_SET, meshset);
524 for (int ee = 0; ee != numEdges; ee++) {
525 if (singular_edges[ee]) {
526 EntityHandle ent;
527 CHKERR moab_ref.side_element(tet, 1, ee, ent);
528 CHKERR moab_ref.add_entities(meshset, &ent, 1);
529 }
530 }
531
532 // refine mesh
533 auto *m_ref = m_field_ref.getInterface<MeshRefinement>();
534 for (int ll = 0; ll != max_level; ll++) {
535 Range edges;
536 CHKERR m_field_ref.getInterface<BitRefManager>()
537 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(ll),
538 BitRefLevel().set(), MBEDGE,
539 edges);
540 Range ref_edges;
541 CHKERR moab_ref.get_adjacencies(
542 nodes_at_front, 1, true, ref_edges, moab::Interface::UNION);
543 ref_edges = intersect(ref_edges, edges);
544 Range ents;
545 CHKERR moab_ref.get_entities_by_type(meshset, MBEDGE, ents, true);
546 ref_edges = intersect(ref_edges, ents);
547 Range tets;
548 CHKERR m_field_ref.getInterface<BitRefManager>()
549 ->getEntitiesByTypeAndRefLevel(
550 BitRefLevel().set(ll), BitRefLevel().set(), MBTET, tets);
551 CHKERR m_ref->addVerticesInTheMiddleOfEdges(
552 ref_edges, BitRefLevel().set(ll + 1));
553 CHKERR m_ref->refineTets(tets, BitRefLevel().set(ll + 1));
554 CHKERR m_field_ref.getInterface<BitRefManager>()
555 ->updateMeshsetByEntitiesChildren(meshset,
556 BitRefLevel().set(ll + 1),
557 meshset, MBEDGE, true);
558 }
559
560 // get ref coords
561 Range tets;
562 CHKERR m_field_ref.getInterface<BitRefManager>()
563 ->getEntitiesByTypeAndRefLevel(BitRefLevel().set(max_level),
564 BitRefLevel().set(), MBTET,
565 tets);
566
567 if (debug) {
568 CHKERR save_range(moab_ref, "ref_tets.vtk", tets);
569 }
570
571 MatrixDouble ref_coords(tets.size(), 12, false);
572 int tt = 0;
573 for (Range::iterator tit = tets.begin(); tit != tets.end();
574 tit++, tt++) {
575 int num_nodes;
576 const EntityHandle *conn;
577 CHKERR moab_ref.get_connectivity(*tit, conn, num_nodes, false);
578 CHKERR moab_ref.get_coords(conn, num_nodes, &ref_coords(tt, 0));
579 }
580
581 auto &data = fe_ptr->dataOnElement[H1];
582 const size_t nb_gauss_pts = fe_ptr->gaussPts.size2();
583 MatrixDouble ref_gauss_pts(4, nb_gauss_pts * ref_coords.size1());
584 MatrixDouble &shape_n =
585 data->dataOnEntities[MBVERTEX][0].getN(NOBASE);
586 int gg = 0;
587 for (size_t tt = 0; tt != ref_coords.size1(); tt++) {
588 double *tet_coords = &ref_coords(tt, 0);
589 double det = Tools::tetVolume(tet_coords);
590 det *= 6;
591 for (size_t ggg = 0; ggg != nb_gauss_pts; ++ggg, ++gg) {
592 for (int dd = 0; dd != 3; dd++) {
593 ref_gauss_pts(dd, gg) =
594 shape_n(ggg, 0) * tet_coords[3 * 0 + dd] +
595 shape_n(ggg, 1) * tet_coords[3 * 1 + dd] +
596 shape_n(ggg, 2) * tet_coords[3 * 2 + dd] +
597 shape_n(ggg, 3) * tet_coords[3 * 3 + dd];
598 }
599 ref_gauss_pts(3, gg) = fe_ptr->gaussPts(3, ggg) * det;
600 }
601 }
602
603 mapRefCoords[singular_nodes.to_ulong()].swap(ref_gauss_pts);
604 CHKERR set_gauss_pts(mapRefCoords[singular_nodes.to_ulong()]);
605
606 // clear cache bubble
607 cachePhi->get<0>() = 0;
608 cachePhi->get<1>() = 0;
609 // tet base cache
610 TetPolynomialBase::switchCacheBaseOff<HDIV>({fe_raw_ptr});
611 TetPolynomialBase::switchCacheBaseOn<HDIV>({fe_raw_ptr});
612
614 };
615
616 CHKERR refine_quadrature();
617 }
618 }
619 }
620
622 }
@ 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 ShapeMBTET(double *N, const double *G_X, const double *G_Y, const double *G_Z, int DIM)
calculate shape functions
Definition fem_tools.c:306
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_3D_TABLE_SIZE
Definition quad.h:186
static QUAD *const QUAD_3D_TABLE[]
Definition quad.h:187
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 double tetVolume(const double *coords)
Calculate volume of tetrahedron.
Definition Tools.cpp:30
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
int order
Definition quad.h:28
int npoints
Definition quad.h:29
auto save_range

Member Data Documentation

◆ cachePhi

boost::shared_ptr<CGGUserPolynomialBase::CachePhi> EshelbianPlasticity::SetIntegrationAtFrontVolume::cachePhi
private

◆ frontEdges

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

◆ frontNodes

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

◆ funRule

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

◆ mapRefCoords

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

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