v0.15.0
Loading...
Searching...
No Matches
EshelbianPlasticity::SetIntegrationAtFrontVolume Struct Reference
Collaboration diagram for EshelbianPlasticity::SetIntegrationAtFrontVolume:
[legend]

Classes

struct  Fe
 

Public Member Functions

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

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 360 of file EshelbianPlasticity.cpp.

Constructor & Destructor Documentation

◆ SetIntegrationAtFrontVolume()

EshelbianPlasticity::SetIntegrationAtFrontVolume::SetIntegrationAtFrontVolume ( boost::shared_ptr< Range > front_nodes,
boost::shared_ptr< Range > front_edges )
inline
Examples
EshelbianPlasticity.cpp.

Definition at line 362 of file EshelbianPlasticity.cpp.

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

Member Function Documentation

◆ operator()()

MoFEMErrorCode EshelbianPlasticity::SetIntegrationAtFrontVolume::operator() ( ForcesAndSourcesCore * fe_raw_ptr,
int order_row,
int order_col,
int order_data )
inline
Examples
EshelbianPlasticity.cpp.

Definition at line 366 of file EshelbianPlasticity.cpp.

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

◆ frontEdges

boost::shared_ptr<Range> EshelbianPlasticity::SetIntegrationAtFrontVolume::frontEdges
private
Examples
EshelbianPlasticity.cpp.

Definition at line 611 of file EshelbianPlasticity.cpp.

◆ frontNodes

boost::shared_ptr<Range> EshelbianPlasticity::SetIntegrationAtFrontVolume::frontNodes
private
Examples
EshelbianPlasticity.cpp.

Definition at line 610 of file EshelbianPlasticity.cpp.

◆ mapRefCoords

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

Definition at line 613 of file EshelbianPlasticity.cpp.


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