v0.14.0
Loading...
Searching...
No Matches
SeriesRecorder.cpp
Go to the documentation of this file.
1/** \file SeriesRecorder.cpp
2 * \brief Record and save data series in time or load steps
3 */
4
5
6#define SeriesRecorderFunctionBegin \
7 MoFEMFunctionBegin; \
8 MOFEM_LOG_CHANNEL("WORLD"); \
9 MOFEM_LOG_CHANNEL("SYNC"); \
10 MOFEM_LOG_FUNCTION(); \
11 MOFEM_LOG_TAG("SYNC", "SeriesRecorder"); \
12 MOFEM_LOG_TAG("WORLD", "SeriesRecorder")
13
14namespace MoFEM {
15
16// const static int debug = 1;
17
19SeriesRecorder::query_interface(boost::typeindex::type_index type_index,
20 UnknownInterface **iface) const {
21 *iface = const_cast<SeriesRecorder *>(this);
22 return 0;
23}
24
26 : cOre(const_cast<Core &>(core)) {}
27
29
30 Interface &m_field = cOre;
31 moab::Interface &moab = m_field.get_moab();
33 const int def_val_len = 0;
34 CHKERR moab.tag_get_handle(
35 "_SeriesName", def_val_len, MB_TYPE_OPAQUE, th_SeriesName,
36 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
38}
39
42 sEries.clear();
43 seriesSteps.clear();
45}
46
48 Interface &m_field = cOre;
49 moab::Interface &moab = m_field.get_moab();
51 Range meshsets;
52 const int def_val_len = 0;
53 CHKERR moab.tag_get_handle(
54 "_SeriesName", def_val_len, MB_TYPE_OPAQUE, th_SeriesName,
55 MB_TAG_CREAT | MB_TAG_BYTES | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
56
57 CHKERR moab.get_entities_by_type(0, MBENTITYSET, meshsets, false);
58 for (Range::iterator mit = meshsets.begin(); mit != meshsets.end(); mit++) {
59 const void *tagName;
60 int tagNameSize;
61 rval = moab.tag_get_by_ptr(th_SeriesName, &*mit, 1, (const void **)&tagName,
62 &tagNameSize);
63 if (rval == MB_SUCCESS) {
64 std::pair<Series_multiIndex::iterator, bool> p =
65 sEries.insert(FieldSeries(moab, *mit));
66 if (verb > QUIET)
67 MOFEM_LOG("SYNC", Sev::inform) << "read series " << *p.first;
68 }
69 }
70 // //build series steps
71 for (Series_multiIndex::iterator sit = sEries.begin(); sit != sEries.end();
72 sit++) {
73 int nb_steps;
74 CHKERR sit->get_nb_steps(moab, nb_steps);
75 int ss = 0;
76 for (; ss < nb_steps; ss++) {
77 std::pair<SeriesStep_multiIndex::iterator, bool> p =
78 seriesSteps.insert(FieldSeriesStep(moab, &*sit, ss));
79 if (verb > QUIET)
80 MOFEM_LOG("SYNC", Sev::inform) << "add series step " << *p.first;
81 }
82 }
83 MOFEM_LOG_CHANNEL("WORLD");
84 MOFEM_LOG_CHANNEL("SYNC");
85 MOFEM_LOG_CHANNEL("SELF");
87}
88
90SeriesRecorder::add_series_recorder(const std::string &series_name) {
91
92 Interface &m_field = cOre;
93 moab::Interface &moab = m_field.get_moab();
95 EntityHandle meshset;
96 CHKERR moab.create_meshset(MESHSET_SET, meshset);
97 void const *tag_data[] = {series_name.c_str()};
98 int tag_sizes[1];
99 tag_sizes[0] = series_name.size();
100 CHKERR moab.tag_set_by_ptr(th_SeriesName, &meshset, 1, tag_data, tag_sizes);
101 std::pair<Series_multiIndex::iterator, bool> p =
102 sEries.insert(FieldSeries(moab, meshset));
103 if (!p.second) {
104 SETERRQ1(PETSC_COMM_SELF, 1, "series recorder <%s> is already there",
105 series_name.c_str());
106 }
108}
109
111SeriesRecorder::delete_recorder_series(const std::string &series_name) {
112
113 Interface &m_field = cOre;
114 moab::Interface &moab = m_field.get_moab();
116 SeriesStep_multiIndex::index<SeriesName_mi_tag>::type::iterator ssit, hi_ssit;
117 ssit = seriesSteps.get<SeriesName_mi_tag>().lower_bound(series_name);
118 hi_ssit = seriesSteps.get<SeriesName_mi_tag>().upper_bound(series_name);
119 seriesSteps.get<SeriesName_mi_tag>().erase(ssit, hi_ssit);
120 Series_multiIndex::index<SeriesName_mi_tag>::type::iterator sit;
121 sit = sEries.get<SeriesName_mi_tag>().find(series_name);
122 if (sit == sEries.get<SeriesName_mi_tag>().end()) {
123 SETERRQ1(m_field.get_comm(), MOFEM_DATA_INCONSISTENCY,
124 "series recorder <%s> not exist and can be deleted",
125 series_name.c_str());
126 }
127 EntityHandle series_meshset = sit->getMeshset();
128 CHKERR moab.tag_delete(sit->th_SeriesTime);
129 CHKERR moab.tag_delete(sit->th_SeriesDataHandles);
130 CHKERR moab.tag_delete(sit->th_SeriesDataUIDs);
131 CHKERR moab.tag_delete(sit->th_SeriesData);
132 sEries.get<SeriesName_mi_tag>().erase(sit);
133 std::vector<EntityHandle> contained;
134 CHKERR moab.get_contained_meshsets(series_meshset, contained);
135 CHKERR moab.remove_entities(series_meshset, &contained[0], contained.size());
136 CHKERR moab.delete_entities(&contained[0], contained.size());
137 CHKERR moab.delete_entities(&series_meshset, 1);
139}
140
141MoFEMErrorCode SeriesRecorder::record_problem(const std::string &serie_name,
142 const Problem *problemPtr,
143 RowColData rc) {
145
146 Series_multiIndex::index<SeriesName_mi_tag>::type::iterator sit =
147 sEries.get<SeriesName_mi_tag>().find(serie_name);
148 if (sit == sEries.get<SeriesName_mi_tag>().end()) {
149 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
150 "series recorder <%s> not exist", serie_name.c_str());
151 }
152 switch (rc) {
153 case ROW:
154 CHKERR const_cast<FieldSeries *>(&*sit)->push_dofs(
155 problemPtr->numeredRowDofsPtr->begin(),
156 problemPtr->numeredRowDofsPtr->end());
157 break;
158 case COL:
159 CHKERR const_cast<FieldSeries *>(&*sit)->push_dofs(
160 problemPtr->numeredColDofsPtr->begin(),
161 problemPtr->numeredColDofsPtr->end());
162 break;
163 default:
164 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
165 }
167}
168
169MoFEMErrorCode SeriesRecorder::record_problem(const std::string &serie_name,
170 const std::string &problem_name,
171 RowColData rc) {
172
173 Interface &m_field = cOre;
175 const Problem *problem_ptr;
176 CHKERR m_field.get_problem(problem_name, &problem_ptr);
177 CHKERR record_problem(serie_name, problem_ptr, rc);
179}
180
181MoFEMErrorCode SeriesRecorder::record_field(const std::string &serie_name,
182 const std::string &field_name,
183 const BitRefLevel &bit,
184 const BitRefLevel &mask) {
185 Interface &m_field = cOre;
186 auto dofs_ptr = m_field.get_dofs();
188 Series_multiIndex::index<SeriesName_mi_tag>::type::iterator sit =
189 sEries.get<SeriesName_mi_tag>().find(serie_name);
190 if (sit == sEries.get<SeriesName_mi_tag>().end()) {
191 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND, "serie recorder <%s> not exist",
192 serie_name.c_str());
193 }
194
195 const auto bit_number = m_field.get_field_bit_number(field_name);
196
197 auto dit = dofs_ptr->get<Unique_mi_tag>().lower_bound(
199 if (dit == dofs_ptr->get<Unique_mi_tag>().end())
200 SETERRQ1(PETSC_COMM_SELF, MOFEM_NOT_FOUND, "field <%s> not exist", field_name.c_str());
201 auto hi_dit = dofs_ptr->get<Unique_mi_tag>().upper_bound(
203
204 for (; dit != hi_dit; dit++) {
205 const BitRefLevel &dof_bit = (*dit)->getBitRefLevel();
206 if ((dof_bit & mask) != dof_bit)
207 continue;
208 if ((dof_bit & bit).any()) {
209 EntityHandle ent = (*dit)->getEnt();
210 ShortId uid = (*dit)->getNonNonuniqueShortId();
211 FieldData val = (*dit)->getFieldData();
212 CHKERR const_cast<FieldSeries *>(&*sit)->push_dofs(ent, uid, val);
213 }
214 }
216}
217
218MoFEMErrorCode SeriesRecorder::record_begin(const std::string &serie_name) {
219
221 Series_multiIndex::index<SeriesName_mi_tag>::type::iterator sit =
222 sEries.get<SeriesName_mi_tag>().find(serie_name);
223 if (sit == sEries.get<SeriesName_mi_tag>().end()) {
224 SETERRQ1(PETSC_COMM_SELF, 1, "serie recorder <%s> not exist",
225 serie_name.c_str());
226 }
227 CHKERR const_cast<FieldSeries *>(&*sit)->begin();
229}
230
231MoFEMErrorCode SeriesRecorder::record_end(const std::string &serie_name,
232 double time) {
233
235 Series_multiIndex::index<SeriesName_mi_tag>::type::iterator sit =
236 sEries.get<SeriesName_mi_tag>().find(serie_name);
237 if (sit == sEries.get<SeriesName_mi_tag>().end()) {
238 SETERRQ1(PETSC_COMM_SELF, 1, "serie recorder <%s> not exist",
239 serie_name.c_str());
240 }
241 CHKERR const_cast<FieldSeries *>(&*sit)->end(time);
243}
244
246SeriesRecorder::initialize_series_recorder(const std::string &serie_name) {
247
248 Interface &m_field = cOre;
249 moab::Interface &moab = m_field.get_moab();
251 Series_multiIndex::index<SeriesName_mi_tag>::type::iterator sit =
252 sEries.get<SeriesName_mi_tag>().find(serie_name);
253 if (sit == sEries.get<SeriesName_mi_tag>().end()) {
254 SETERRQ1(PETSC_COMM_SELF, 1, "serie recorder <%s> not exist",
255 serie_name.c_str());
256 }
257 CHKERR const_cast<FieldSeries *>(&*sit)->read(moab);
258 const_cast<FieldSeries *>(&*sit)->record_begin = false;
259 const_cast<FieldSeries *>(&*sit)->record_end = false;
261}
262
264SeriesRecorder::finalize_series_recorder(const std::string &serie_name) {
265
266 Interface &m_field = cOre;
267 moab::Interface &moab = m_field.get_moab();
269 Series_multiIndex::index<SeriesName_mi_tag>::type::iterator sit =
270 sEries.get<SeriesName_mi_tag>().find(serie_name);
271 if (sit == sEries.get<SeriesName_mi_tag>().end()) {
272 SETERRQ1(PETSC_COMM_SELF, 1, "serie recorder <%s> not exist",
273 serie_name.c_str());
274 }
275 CHKERR sit->save(moab);
276 int nb_steps;
277 CHKERR sit->get_nb_steps(moab, nb_steps);
278 int ss = 0;
279 for (; ss < nb_steps; ss++) {
280 /*std::pair<SeriesStep_multiIndex::iterator,bool> p =*/seriesSteps.insert(
281 FieldSeriesStep(moab, &*sit, ss));
282 }
284}
285
288
289 for (auto &sit : sEries.get<SeriesName_mi_tag>())
290 MOFEM_LOG("SYNC", Sev::inform) << "series " << sit;
291
292 for (auto &ssit : seriesSteps.get<SeriesName_mi_tag>())
293 MOFEM_LOG("SYNC", Sev::inform) << "series steps " << ssit;
294
295 MOFEM_LOG_CHANNEL("WORLD");
296 MOFEM_LOG_CHANNEL("SYNC");
297 MOFEM_LOG_CHANNEL("SELF");
299}
300
301bool SeriesRecorder::check_series(const std::string &name) const {
302 Series_multiIndex::index<SeriesName_mi_tag>::type::iterator sit =
303 sEries.get<SeriesName_mi_tag>().find(name);
304 if (sit != sEries.get<SeriesName_mi_tag>().end())
305 return true;
306 return false;
307}
308
310 const int step_number) {
311
312 Interface &m_field = cOre;
313 moab::Interface &moab = m_field.get_moab();
314 auto dofs_ptr = m_field.get_dofs();
316 SeriesStep_multiIndex::index<
317 Composite_SeriesName_And_Step_mi_tag>::type::iterator sit;
319 boost::make_tuple(serie_name, step_number));
320 if (sit == seriesSteps.get<Composite_SeriesName_And_Step_mi_tag>().end()) {
321 SETERRQ2(PETSC_COMM_SELF, 1, "series <%s> and step %d not found",
322 serie_name.c_str(), step_number);
323 }
324 CHKERR sit->get(moab, *(const_cast<DofEntity_multiIndex *>(dofs_ptr)));
326}
327
328SeriesStep_multiIndex::index<SeriesName_mi_tag>::type::iterator
330 return seriesSteps.get<SeriesName_mi_tag>().lower_bound(name);
331}
332
333SeriesStep_multiIndex::index<SeriesName_mi_tag>::type::iterator
335 return seriesSteps.get<SeriesName_mi_tag>().upper_bound(name);
336}
337
338} // namespace MoFEM
static Index< 'p', 3 > p
#define SeriesRecorderFunctionBegin
@ QUIET
Definition: definitions.h:208
RowColData
RowColData.
Definition: definitions.h:123
@ COL
Definition: definitions.h:123
@ ROW
Definition: definitions.h:123
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
multi_index_container< boost::shared_ptr< DofEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< DofEntity, UId, &DofEntity::getLocalUniqueId > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< DofEntity, EntityHandle, &DofEntity::getEnt > > > > DofEntity_multiIndex
MultiIndex container keeps DofEntity.
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
virtual const DofEntity_multiIndex * get_dofs() const =0
Get the dofs object.
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
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 bool check_series(const std::string &name) const
check if series is in database
virtual MoFEMErrorCode record_problem(const std::string &serie_name, const Problem *problemPtr, RowColData rc)
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)
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)
virtual MoFEMErrorCode delete_recorder_series(const std::string &series_name)
auto bit
set bit
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
int ShortId
Unique Id in the field.
Definition: Types.hpp:32
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
constexpr auto field_name
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
Core (interface) class.
Definition: Core.hpp:82
Deprecated interface functions.
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
Structure for recording (time) series.
Structure for keeping time and step.
keeps basic data about problem
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredRowDofsPtr
store DOFs on rows for this problem
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredColDofsPtr
store DOFs on columns for this problem
virtual SeriesStep_multiIndex::index< SeriesName_mi_tag >::type::iterator get_series_steps_byName_begin(const std::string &name)
Series_multiIndex sEries
recorded series
MoFEMErrorCode clearMap()
clear multi-index container
SeriesStep_multiIndex seriesSteps
recorded series steps
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
SeriesRecorder(const MoFEM::Core &core)
Tag th_SeriesName
Recorded series name.
MoFEMErrorCode getTags(int verb=-1)
get tags handlers used on meshsets
virtual SeriesStep_multiIndex::index< SeriesName_mi_tag >::type::iterator get_series_steps_byName_end(const std::string &name)
MoFEMErrorCode initialiseDatabaseFromMesh(int verb=0)
base class for all interface classes