v0.14.0
Loading...
Searching...
No Matches
remove_mofem_meshsets.cpp
Go to the documentation of this file.
1/** \file field_to_vertices.cpp
2 \brief Field to vertices
3 \example field_to_vertices.cpp
4
5*/
6
7
8
9#include <MoFEM.hpp>
10
11using namespace MoFEM;
12
13static char help[] = "...\n\n";
14
15int main(int argc, char *argv[]) {
16 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
17
18 try {
19
20 // global variables
21 char mesh_file_name[255];
22 PetscBool flg_file = PETSC_FALSE;
23 char mesh_out_file[255] = "out.h5m";
24
25 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "Field to vertices options",
26 "none");
27 CHKERR PetscOptionsString("-file_name", "mesh file name", "", "mesh.h5m",
28 mesh_file_name, 255, &flg_file);
29 CHKERR PetscOptionsString("-output_file", "output mesh file name", "",
30 "out.h5m", mesh_out_file, 255, PETSC_NULL);
31
32 ierr = PetscOptionsEnd(); CHKERRG(ierr);
33
34 moab::Core mb_instance;
35 moab::Interface &moab = mb_instance;
36 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
37 if (pcomm == NULL)
38 pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
39 const char *option;
40 option = "";
41 CHKERR moab.load_file(mesh_file_name, 0, option);
42
43 // Create MoFEM database
44 MoFEM::Core core(moab);
45 MoFEM::Interface &m_field = core;
46
47 // Add logging channel for example
48 auto core_log = logging::core::get();
49 core_log->add_sink(
51 LogManager::setLog("REMOVER");
52 MOFEM_LOG_TAG("REMOVER", "remover");
53
54 auto prb_ptr = m_field.get_problems();
55 std::vector<std::string> prb_list;
56 for(auto &it : *prb_ptr)
57 prb_list.push_back(it.getName());
58
59 for (auto &it : prb_list) {
60 MOFEM_LOG("REMOVER", Sev::inform) << "Delete problem " << it;
61 CHKERR m_field.delete_problem(it);
62 }
63
64 auto fe_ptr = m_field.get_finite_elements();
65 std::vector<std::string> fe_list;
66 for (auto &it : *fe_ptr)
67 fe_list.push_back(it->getName());
68
69 for (auto &it : fe_list) {
70 MOFEM_LOG("REMOVER", Sev::inform)
71 << "Delete finite element " << it;
72 CHKERR m_field.delete_finite_element(it);
73 }
74
75 auto field_ptr = m_field.get_fields();
76 std::vector<std::string> field_list;
77 for (auto &it : *field_ptr)
78 field_list.push_back(it->getName());
79
80 for (auto &it : field_list) {
81 MOFEM_LOG("REMOVER", Sev::inform) << "Delete field " << it;
82 CHKERR m_field.delete_field(it);
83 }
84
85 CHKERR moab.write_file(mesh_out_file);
86
87 }
89
91
92 return 0;
93}
int main()
#define CATCH_ERRORS
Catch errors.
#define MYPCOMM_INDEX
default communicator number PCOMM
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
#define CHKERR
Inline error check.
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
virtual const Problem_multiIndex * get_problems() const =0
Get the problems object.
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
virtual MoFEMErrorCode delete_problem(const std::string name)=0
Delete problem.
char mesh_file_name[255]
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
static char help[]
virtual MoFEMErrorCode delete_finite_element(const std::string name, int verb=DEFAULT_VERBOSITY)=0
delete finite element from mofem database
virtual MoFEMErrorCode delete_field(const std::string name, int verb=DEFAULT_VERBOSITY)=0
Delete field.
Core (interface) class.
Definition Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.