v0.15.0
Loading...
Searching...
No Matches
Classes | Public Member Functions | Private Attributes | Static Private Attributes | List of all members
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 363 of file EshelbianPlasticity.cpp.

Constructor & Destructor Documentation

◆ SetIntegrationAtFrontVolume()

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

Definition at line 365 of file EshelbianPlasticity.cpp.

367 : 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 369 of file EshelbianPlasticity.cpp.

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

◆ frontEdges

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

Definition at line 614 of file EshelbianPlasticity.cpp.

◆ frontNodes

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

Definition at line 613 of file EshelbianPlasticity.cpp.

◆ mapRefCoords

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

Definition at line 616 of file EshelbianPlasticity.cpp.


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