v0.14.0
Classes | Functions | Variables
loop_entities.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Classes

struct  TestEntityMethod
 

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)
Examples
loop_entities.cpp.

Definition at line 47 of file loop_entities.cpp.

47  {
48  // initialize petsc
49  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
50 
51  try {
52 
53  // read mesh and create moab and mofem data structures
54  moab::Core mb_instance;
55  moab::Interface &moab = mb_instance;
56 
57  MoFEM::Core core(moab, PETSC_COMM_WORLD);
58  MoFEM::Interface &m_field = core;
59 
60  // Register DM Manager
61  DMType dm_name = "DMMOFEM";
62  CHKERR DMRegister_MoFEM(dm_name);
63 
64  EntityHandle root_set = moab.get_root_set();
65 
66  // Simple interface
67  Simple *simple_interface;
68  CHKERR m_field.getInterface(simple_interface);
69  {
70  // get options from command line
71  CHKERR simple_interface->getOptions();
72  // load mesh file
73  CHKERR simple_interface->loadFile();
74  // add fields
75  CHKERR simple_interface->addDomainField("FIELD1", H1,
77  // set fields order
78  CHKERR simple_interface->setFieldOrder("FIELD1", 4);
79  // setup problem
80  CHKERR simple_interface->setUp();
81 
82  Range all_ents;
83  CHKERR m_field.get_moab().get_entities_by_handle(root_set, all_ents);
84  all_ents = subtract(all_ents, all_ents.subset_by_type(MBENTITYSET));
85 
86  auto testingEntitiesInDatabase = [&]() {
88  TestEntityMethod method(all_ents);
89  CHKERR m_field.loop_entities("FIELD1", method);
91  };
92 
93  CHKERR testingEntitiesInDatabase();
94 
95  auto testingEntitiesInDatabaseSubRange = [&]() {
97  Range edges = all_ents.subset_by_type(MBEDGE);
98  TestEntityMethod method(edges);
99  CHKERR m_field.loop_entities("FIELD1", method, &edges);
101  };
102 
103  CHKERR testingEntitiesInDatabaseSubRange();
104 
105  auto dm = simple_interface->getDM();
106 
107  auto testingEntitiesInDM = [&]() {
109  TestEntityMethod method(all_ents);
110  const MoFEM::Problem *prb_ptr;
111  CHKERR DMMoFEMGetProblemPtr(dm, &prb_ptr);
112  CHKERR m_field.loop_entities(prb_ptr, "FIELD1", ROW, method, 0,
113  m_field.get_comm_size());
115  };
116 
117  CHKERR testingEntitiesInDM();
118 
119  }
120  }
121  CATCH_ERRORS;
122 
123  // finish work cleaning memory, getting statistics, ect.
125 
126  return 0;
127 }

Variable Documentation

◆ help

char help[] = "...\n\n"
static
Examples
loop_entities.cpp.

Definition at line 13 of file loop_entities.cpp.

MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
help
static char help[]
Definition: loop_entities.cpp:13
EntityHandle
MoFEM::Simple::loadFile
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
ROW
@ ROW
Definition: definitions.h:136
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MoFEM::Simple::getOptions
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEM::Simple::getDM
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:747
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
TestEntityMethod
Definition: loop_entities.cpp:15
MoFEM::Simple::addDomainField
MoFEMErrorCode addDomainField(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_ZERO, int verb=-1)
Add field on domain.
Definition: Simple.cpp:264
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::CoreInterface::get_comm_size
virtual int get_comm_size() const =0
MoFEM::Simple::setFieldOrder
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:545
Range
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::DMMoFEMGetProblemPtr
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition: DMMoFEM.cpp:426
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
MoFEM::Simple::setUp
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:683
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::CoreInterface::loop_entities
virtual MoFEMErrorCode loop_entities(const std::string field_name, EntityMethod &method, Range const *const ents=nullptr, int verb=DEFAULT_VERBOSITY)=0
Loop over field entities.