v0.13.2
Loading...
Searching...
No Matches
Functions | Variables
record_series_atom.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Functions

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

Variables

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

Function Documentation

◆ main()

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

Definition at line 13 of file record_series_atom.cpp.

13 {
14
15 MoFEM::Core::Initialize(&argc, &argv, PETSC_NULL, help);
16
17 try {
18
19 moab::Core mb_instance;
20 moab::Interface &moab = mb_instance;
21 int rank;
22 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
23
24 // Reade parameters from line command
25 PetscBool flg = PETSC_TRUE;
26 char mesh_file_name[255];
27#if PETSC_VERSION_GE(3, 6, 4)
28 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
29 255, &flg);
30#else
31 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
32 mesh_file_name, 255, &flg);
33#endif
34 if (flg != PETSC_TRUE) {
35 SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
36 }
37 PetscInt order = 3;
38#if PETSC_VERSION_GE(3, 6, 4)
39 CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-my_order", &order, &flg);
40#else
41 CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
42 &flg);
43#endif
44
45 // Read mesh to MOAB
46 const char *option;
47 option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
48 CHKERR moab.load_file(mesh_file_name, 0, option);
49
50 // Create MoFEM
51 // Starting from one MoFEM::CoreTmp<1> for testing purposes
52 MoFEM::CoreTmp<0> core(moab);
53 MoFEM::Interface &m_field = core;
54
55 // stl::bitset see for more details
56 auto bit_level0 = BitRefLevel().set(0);
57 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
58 0, 3, BitRefLevel().set(0));
59
60 /***/
61 // Define problem
62
63 // Fields
64 CHKERR m_field.add_field("FIELD_A", H1, AINSWORTH_LEGENDRE_BASE, 3);
65 CHKERR m_field.add_field("FIELD_B", H1, AINSWORTH_LEGENDRE_BASE, 3);
66
67 CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "FIELD_A");
68 CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "FIELD_B");
69
70 CHKERR m_field.set_field_order(0, MBTET, "FIELD_A", order);
71 CHKERR m_field.set_field_order(0, MBTRI, "FIELD_A", order);
72 CHKERR m_field.set_field_order(0, MBEDGE, "FIELD_A", order);
73 CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD_A", 1);
74
75 CHKERR m_field.set_field_order(0, MBTET, "FIELD_B", order);
76 CHKERR m_field.set_field_order(0, MBTRI, "FIELD_B", order);
77 CHKERR m_field.set_field_order(0, MBEDGE, "FIELD_B", order);
78 CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD_B", 1);
79
80 // build field
81 CHKERR m_field.build_fields();
82
83 CHKERR m_field.getInterface<FieldBlas>()->setField(0, MBVERTEX, "FIELD_B");
84 CHKERR m_field.getInterface<FieldBlas>()->setField(1, MBVERTEX, "FIELD_A");
85
86 SeriesRecorder *recorder_ptr;
87 CHKERR m_field.getInterface(recorder_ptr);
88
89 CHKERR recorder_ptr->add_series_recorder("TEST_SERIES1");
90 CHKERR recorder_ptr->add_series_recorder("TEST_SERIES2");
91
92 // initialize
93 CHKERR recorder_ptr->initialize_series_recorder("TEST_SERIES1");
94
95 CHKERR recorder_ptr->record_begin("TEST_SERIES1");
96 CHKERR recorder_ptr->record_field("TEST_SERIES1", "FIELD_B", bit_level0,
97 bit_level0);
98 CHKERR recorder_ptr->record_end("TEST_SERIES1", 1);
99
100 CHKERR m_field.getInterface<FieldBlas>()->fieldAxpy(1., "FIELD_A",
101 "FIELD_B");
102 CHKERR recorder_ptr->record_begin("TEST_SERIES1");
103 CHKERR recorder_ptr->record_field("TEST_SERIES1", "FIELD_B", bit_level0,
104 bit_level0);
105
106 CHKERR recorder_ptr->initialize_series_recorder("TEST_SERIES2");
107 CHKERR recorder_ptr->record_begin("TEST_SERIES2");
108 CHKERR recorder_ptr->record_field("TEST_SERIES2", "FIELD_A", bit_level0,
109 bit_level0);
110 CHKERR recorder_ptr->record_field("TEST_SERIES2", "FIELD_B", bit_level0,
111 bit_level0);
112 CHKERR recorder_ptr->record_end("TEST_SERIES2", 1);
113 CHKERR recorder_ptr->finalize_series_recorder("TEST_SERIES2");
114
115 CHKERR recorder_ptr->record_end("TEST_SERIES1", 2);
116
117 // finalize
118 CHKERR recorder_ptr->finalize_series_recorder("TEST_SERIES1");
119 CHKERR recorder_ptr->print_series_steps();
120
121 CHKERR m_field.getInterface<FieldBlas>()->fieldScale(2, "FIELD_A");
122
123 // Setting MoFEM::CoreTmp<-1> for testing purposes
124 MoFEM::CoreTmp<-1> core2(moab);
125 MoFEM::Interface &m_field2 = core2;
126
127 // build field
128 CHKERR m_field2.build_fields();
129
130 typedef tee_device<std::ostream, std::ofstream> TeeDevice;
131 typedef stream<TeeDevice> TeeStream;
132 std::ofstream ofs("record_series_atom.txt");
133 TeeDevice my_tee(std::cout, ofs);
134 TeeStream my_split(my_tee);
135
136 SeriesRecorder *recorder2_ptr;
137 CHKERR m_field2.getInterface(recorder2_ptr);
138 CHKERR recorder2_ptr->print_series_steps();
139
140 auto dofs_ptr = m_field.get_dofs();;
141
142 my_split << "TEST_SERIES1" << std::endl;
143 for (_IT_SERIES_STEPS_BY_NAME_FOR_LOOP_(recorder2_ptr, "TEST_SERIES1",
144 sit)) {
145
146 CHKERR recorder2_ptr->load_series_data("TEST_SERIES1",
147 sit->get_step_number());
148
149 my_split << "next step:\n";
150 my_split << *sit << std::endl;
151
152 {
154 for (_IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_(m_field2, "FIELD_B", dof)) {
155 dofs_view.insert(*dof);
156 }
157 for (DofEntity_multiIndex_uid_view::iterator dit = dofs_view.begin();
158 dit != dofs_view.end(); dit++) {
159 my_split << **dit << endl;
160 }
161 }
162 }
163
164 my_split << "TEST_SERIES2" << std::endl;
165 for (_IT_SERIES_STEPS_BY_NAME_FOR_LOOP_(recorder2_ptr, "TEST_SERIES2",
166 sit)) {
167
168 CHKERR recorder2_ptr->load_series_data("TEST_SERIES2",
169 sit->get_step_number());
170
171 my_split << "next step:\n";
172 {
174 for (_IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_(m_field2, "FIELD_A", dof)) {
175 dofs_view.insert(*dof);
176 }
177 for (DofEntity_multiIndex_uid_view::iterator dit = dofs_view.begin();
178 dit != dofs_view.end(); dit++) {
179 my_split << **dit << endl;
180 }
181 }
182 {
184 for (_IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_(m_field2, "FIELD_B", dof)) {
185 dofs_view.insert(*dof);
186 }
187 for (DofEntity_multiIndex_uid_view::iterator dit = dofs_view.begin();
188 dit != dofs_view.end(); dit++) {
189 my_split << **dit << endl;
190 }
191 }
192 }
193 }
195
197 return 0;
198}
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
@ H1
continuous field
Definition: definitions.h:85
#define CHKERR
Inline error check.
Definition: definitions.h:535
tee_device< std::ostream, std::ofstream > TeeDevice
multi_index_container< boost::shared_ptr< DofEntity >, indexed_by< ordered_unique< const_mem_fun< DofEntity, UId, &DofEntity::getLocalUniqueId > > > > DofEntity_multiIndex_uid_view
multi-index view on DofEntity by uid
virtual const DofEntity_multiIndex * get_dofs() const =0
Get the dofs object.
#define _IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_(MFIELD, NAME, IT)
Definition: Interface.hpp:1843
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
virtual MoFEMErrorCode load_series_data(const std::string &serie_name, const int step_number)
virtual MoFEMErrorCode print_series_steps()
virtual MoFEMErrorCode finalize_series_recorder(const std::string &serie_name)
virtual MoFEMErrorCode record_field(const std::string &serie_name, const std::string &field_name, const BitRefLevel &bit, const BitRefLevel &mask)
virtual MoFEMErrorCode initialize_series_recorder(const std::string &serie_name)
#define _IT_SERIES_STEPS_BY_NAME_FOR_LOOP_(RECORDER, NAME, IT)
loop over recorded series step
virtual MoFEMErrorCode add_series_recorder(const std::string &series_name)
virtual MoFEMErrorCode record_begin(const std::string &serie_name)
virtual MoFEMErrorCode record_end(const std::string &serie_name, double time=0)
char mesh_file_name[255]
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
static char help[]
Managing BitRefLevels.
virtual MoFEMErrorCode add_field(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_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
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.
Basic algebra on fields.
Definition: FieldBlas.hpp:21
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

Variable Documentation

◆ help

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

Definition at line 11 of file record_series_atom.cpp.