GetFEM++  5.3
getfem_export.cc
1 /*===========================================================================
2 
3  Copyright (C) 2001-2017 Yves Renard
4 
5  This file is a part of GetFEM++
6 
7  GetFEM++ is free software; you can redistribute it and/or modify it
8  under the terms of the GNU Lesser General Public License as published
9  by the Free Software Foundation; either version 3 of the License, or
10  (at your option) any later version along with the GCC Runtime Library
11  Exception either version 3.1 or (at your option) any later version.
12  This program is distributed in the hope that it will be useful, but
13  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
14  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
15  License and GCC Runtime Library Exception for more details.
16  You should have received a copy of the GNU Lesser General Public License
17  along with this program; if not, write to the Free Software Foundation,
18  Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
19 
20 ===========================================================================*/
21 
22 #include <iomanip>
23 #include "getfem/dal_singleton.h"
25 #include "getfem/getfem_export.h"
26 
27 namespace getfem
28 {
29  /* -------------------------------------------------------------
30  * VTK export
31  * ------------------------------------------------------------- */
32 
33  struct gf2vtk_dof_mapping : public std::vector<std::vector<unsigned> > {};
34  struct gf2vtk_vtk_type : public std::vector<int> {};
35 
36  typedef enum { NO_VTK_MAPPING,
37  N1_TO_VTK_VERTEX,
38  N2_TO_VTK_LINE,
39  N3_TO_VTK_TRIANGLE,
40  N4_TO_VTK_PIXEL,
41  N4_TO_VTK_QUAD,
42  N4_TO_VTK_TETRA,
43  N8_TO_VTK_VOXEL,
44  N8_TO_VTK_HEXAHEDRON,
45  N6_TO_VTK_WEDGE,
46  N5_TO_VTK_PYRAMID,
47  N3_TO_VTK_QUADRATIC_EDGE,
48  N6_TO_VTK_QUADRATIC_TRIANGLE,
49  N8_TO_VTK_QUADRATIC_QUAD,
50  N10_TO_VTK_QUADRATIC_TETRA,
51  N20_TO_VTK_QUADRATIC_HEXAHEDRON,
52  N15_TO_VTK_QUADRATIC_WEDGE,
53  N13_TO_VTK_QUADRATIC_PYRAMID,
54  N14_TO_VTK_QUADRATIC_PYRAMID,
55  N9_TO_VTK_BIQUADRATIC_QUAD,
56  N27_TO_VTK_TRIQUADRATIC_HEXAHEDRON,
57  N18_TO_VTK_BIQUADRATIC_QUADRATIC_WEDGE } vtk_mapping_type;
58 
59 
60  void init_gf2vtk() {
61  // see https://www.kitware.com/products/books/VTKUsersGuide.pdf
62  // and https://www.kitware.com/products/books/VTKTextbook.pdf
63  // (there are some conflicts between the two)
64  gf2vtk_dof_mapping &vtkmaps = dal::singleton<gf2vtk_dof_mapping>::instance();
65  gf2vtk_vtk_type &vtktypes = dal::singleton<gf2vtk_vtk_type>::instance();
66  vtkmaps.resize(25);
67  vtktypes.resize(25);
68 
69  vtktypes[N1_TO_VTK_VERTEX] = vtk_export::VTK_VERTEX;
70  vtkmaps [N1_TO_VTK_VERTEX] = {0};
71  vtktypes[N2_TO_VTK_LINE] = vtk_export::VTK_LINE;
72  vtkmaps [N2_TO_VTK_LINE] = {0, 1};
73  vtktypes[N3_TO_VTK_TRIANGLE] = vtk_export::VTK_TRIANGLE;
74  vtkmaps [N3_TO_VTK_TRIANGLE] = {0, 1, 2};
75  vtktypes[N4_TO_VTK_PIXEL] = vtk_export::VTK_PIXEL;
76  vtkmaps [N4_TO_VTK_PIXEL] = {0, 1, 2, 3};
77  vtktypes[N4_TO_VTK_QUAD] = vtk_export::VTK_QUAD;
78  vtkmaps [N4_TO_VTK_QUAD] = {0, 1, 3, 2};
79  vtktypes[N4_TO_VTK_TETRA] = vtk_export::VTK_TETRA;
80  vtkmaps [N4_TO_VTK_TETRA] = {0, 1, 2, 3};
81  vtktypes[N8_TO_VTK_VOXEL] = vtk_export::VTK_VOXEL;
82  vtkmaps [N8_TO_VTK_VOXEL] = {0, 1, 2, 3, 4, 5, 6, 7};
83  vtktypes[N8_TO_VTK_HEXAHEDRON] = vtk_export::VTK_HEXAHEDRON;
84  vtkmaps [N8_TO_VTK_HEXAHEDRON] = {0, 1, 3, 2, 4, 5, 7, 6};
85  vtktypes[N6_TO_VTK_WEDGE] = vtk_export::VTK_WEDGE;
86  vtkmaps [N6_TO_VTK_WEDGE] = {0, 1, 2, 3, 4, 5};
87  vtktypes[N5_TO_VTK_PYRAMID] = vtk_export::VTK_PYRAMID;
88  vtkmaps [N5_TO_VTK_PYRAMID] = {0, 1, 3, 2, 4};
89  vtktypes[N3_TO_VTK_QUADRATIC_EDGE] = vtk_export::VTK_QUADRATIC_EDGE;
90  vtkmaps [N3_TO_VTK_QUADRATIC_EDGE] = {0, 2, 1};
91  vtktypes[N6_TO_VTK_QUADRATIC_TRIANGLE] = vtk_export::VTK_QUADRATIC_TRIANGLE;
92  vtkmaps [N6_TO_VTK_QUADRATIC_TRIANGLE] = {0, 2, 5, 1, 4, 3};
93  vtktypes[N8_TO_VTK_QUADRATIC_QUAD] = vtk_export::VTK_QUADRATIC_QUAD;
94  vtkmaps [N8_TO_VTK_QUADRATIC_QUAD] = {0, 2, 7, 5, 1, 4, 6, 3};
95  vtktypes[N10_TO_VTK_QUADRATIC_TETRA] = vtk_export::VTK_QUADRATIC_TETRA;
96  vtkmaps [N10_TO_VTK_QUADRATIC_TETRA] = {0, 2, 5, 9, 1, 4, 3, 6, 7, 8};
97  vtktypes[N20_TO_VTK_QUADRATIC_HEXAHEDRON] = vtk_export::VTK_QUADRATIC_HEXAHEDRON;
98  vtkmaps [N20_TO_VTK_QUADRATIC_HEXAHEDRON] = {0, 2, 7, 5, 12, 14, 19, 17, 1, 4, 6, 3, 13, 16, 18, 15, 8, 9, 11, 10};
99  vtktypes[N15_TO_VTK_QUADRATIC_WEDGE] = vtk_export::VTK_QUADRATIC_WEDGE;
100  vtkmaps [N15_TO_VTK_QUADRATIC_WEDGE] = {0, 2, 5, 9, 11, 14, 1, 4, 3, 10, 13, 12, 6, 7, 8};
101  // = {0, 5, 2, 9, 14, 11, 3, 4, 1, 12, 13, 10, 6, 8, 7};
102  vtktypes[N13_TO_VTK_QUADRATIC_PYRAMID] = vtk_export::VTK_QUADRATIC_PYRAMID;
103  vtkmaps [N13_TO_VTK_QUADRATIC_PYRAMID] = {0, 2, 7, 5, 12, 1, 4, 6, 3, 8, 9, 11, 10};
104  vtktypes[N14_TO_VTK_QUADRATIC_PYRAMID] = vtk_export::VTK_QUADRATIC_PYRAMID;
105  vtkmaps [N14_TO_VTK_QUADRATIC_PYRAMID] = {0, 2, 8, 6, 13, 1, 5, 7, 3, 9, 10, 12, 11};
106  vtktypes[N9_TO_VTK_BIQUADRATIC_QUAD] = vtk_export::VTK_BIQUADRATIC_QUAD;
107  vtkmaps [N9_TO_VTK_BIQUADRATIC_QUAD] = {0, 2, 8, 6, 1, 5, 7, 3, 4};
108  vtktypes[N27_TO_VTK_TRIQUADRATIC_HEXAHEDRON] = vtk_export::VTK_TRIQUADRATIC_HEXAHEDRON;
109  vtkmaps [N27_TO_VTK_TRIQUADRATIC_HEXAHEDRON] = {0, 2, 8, 6, 18, 20, 26, 24, 1, 5, 7, 3, 19, 23, 25, 21, 9, 11, 17, 15, 12, 14, 10, 16, 4, 22};
110  vtktypes[N18_TO_VTK_BIQUADRATIC_QUADRATIC_WEDGE] = vtk_export::VTK_BIQUADRATIC_QUADRATIC_WEDGE;
111  vtkmaps [N18_TO_VTK_BIQUADRATIC_QUADRATIC_WEDGE] = {0, 2, 5, 12, 14, 17, 1, 4, 3, 13, 16, 15, 6, 8, 11, 7, 10, 9};
112  }
113 
114  static const std::vector<unsigned> &
115  select_vtk_dof_mapping(unsigned t) {
116  gf2vtk_dof_mapping &vtkmaps = dal::singleton<gf2vtk_dof_mapping>::instance();
117  if (vtkmaps.size() == 0) init_gf2vtk();
118  return vtkmaps[t];
119  }
120 
121  int select_vtk_type(unsigned t) {
122  gf2vtk_vtk_type &vtktypes = dal::singleton<gf2vtk_vtk_type>::instance();
123  if (vtktypes.size() == 0) init_gf2vtk();
124  return vtktypes[t];
125  }
126 
127  vtk_export::vtk_export(std::ostream &os_, bool ascii_)
128  : os(os_), ascii(ascii_) { init(); }
129 
130  vtk_export::vtk_export(const std::string& fname, bool ascii_)
131  : os(real_os), ascii(ascii_),
132  real_os(fname.c_str(), !ascii ? std::ios_base::binary | std::ios_base::out :
133  std::ios_base::out) {
134  GMM_ASSERT1(real_os, "impossible to write to vtk file '" << fname << "'");
135  init();
136  }
137 
138  void vtk_export::init() {
139  static int test_endian = 0x01234567;
140  strcpy(header, "Exported by getfem++");
141  psl = 0; dim_ = dim_type(-1);
142  if (*((char*)&test_endian) == 0x67)
143  reverse_endian = true;
144  else reverse_endian = false;
145  state = EMPTY;
146  }
147 
148  void vtk_export::switch_to_cell_data() {
149  if (state != IN_CELL_DATA) {
150  state = IN_CELL_DATA;
151  write_separ();
152  if (psl) {
153  os << "CELL_DATA " << psl->nb_simplexes(0) + psl->nb_simplexes(1)
154  + psl->nb_simplexes(2) + psl->nb_simplexes(3) << "\n";
155  } else {
156  os << "CELL_DATA " << pmf->convex_index().card() << "\n";
157  }
158  write_separ();
159  }
160  }
161 
162  void vtk_export::switch_to_point_data() {
163  if (state != IN_POINT_DATA) {
164  state = IN_POINT_DATA;
165  write_separ();
166  if (psl) {
167  write_separ(); os << "POINT_DATA " << psl->nb_points() << "\n";
168  } else {
169  os << "POINT_DATA " << pmf_dof_used.card() << "\n";
170  }
171  write_separ();
172  }
173  }
174 
175 
176  /* try to check if a quad or hexahedric cell is "rectangular" and oriented
177  along the axes */
178  template<typename C> static bool check_voxel(const C& c) {
179  scalar_type h[3];
180  unsigned N = c[0].size();
181  if (c.size() != (1U << N)) return false;
182  const base_node P0 = c[0];
183  h[0] = c[1][0] - P0[0];
184  h[1] = c[2][0] - P0[0];
185  if (c.size() != 4) h[2] = c[4][0] - P0[0];
186  for (unsigned i=1; i < c.size(); ++i) {
187  const base_node d = c[i] - P0;
188  for (unsigned j=0; j < N; ++j)
189  if (gmm::abs(d[j]) > 1e-7*h[j] && gmm::abs(d[j] - h[j]) > 1e-7*h[j])
190  return false;
191  }
192  return true;
193  }
194 
195 
196  void vtk_export::exporting(const stored_mesh_slice& sl) {
197  psl = &sl; dim_ = dim_type(sl.dim());
198  GMM_ASSERT1(psl->dim() <= 3, "attempt to export a " << int(dim_)
199  << "D slice (not supported)");
200  }
201 
202  void vtk_export::exporting(const mesh& m) {
203  dim_ = m.dim();
204  GMM_ASSERT1(dim_ <= 3, "attempt to export a " << int(dim_)
205  << "D slice (not supported)");
206  pmf = std::make_unique<mesh_fem>(const_cast<mesh&>(m), dim_type(1));
207  for (dal::bv_visitor cv(m.convex_index()); !cv.finished(); ++cv) {
208  bgeot::pgeometric_trans pgt = m.trans_of_convex(cv);
209  pfem pf = getfem::classical_fem(pgt, pgt->complexity() > 1 ? 2 : 1);
210  pmf->set_finite_element(cv, pf);
211  }
212  exporting(*pmf);
213  }
214 
215  void vtk_export::exporting(const mesh_fem& mf) {
216  dim_ = mf.linked_mesh().dim();
217  GMM_ASSERT1(dim_ <= 3, "attempt to export a " << int(dim_)
218  << "D slice (not supported)");
219  if (&mf != pmf.get())
220  pmf = std::make_unique<mesh_fem>(mf.linked_mesh());
221  /* initialize pmf with finite elements suitable for VTK (which only knows
222  isoparametric FEMs of order 1 and 2) */
223  for (dal::bv_visitor cv(mf.convex_index()); !cv.finished(); ++cv) {
224  bgeot::pgeometric_trans pgt = mf.linked_mesh().trans_of_convex(cv);
225  pfem pf = mf.fem_of_element(cv);
226 
227  if (pf == fem_descriptor("FEM_Q2_INCOMPLETE(2)") ||
228  pf == fem_descriptor("FEM_Q2_INCOMPLETE(3)") ||
229  pf == fem_descriptor("FEM_PYRAMID_Q2_INCOMPLETE") ||
230  pf == fem_descriptor("FEM_PYRAMID_Q2_INCOMPLETE_DISCONTINUOUS") ||
231  pf == fem_descriptor("FEM_PRISM_INCOMPLETE_P2") ||
232  pf == fem_descriptor("FEM_PRISM_INCOMPLETE_P2_DISCONTINUOUS"))
233  pmf->set_finite_element(cv, pf);
234  else {
235  bool discontinuous = false;
236  for (unsigned i=0; i < pf->nb_dof(cv); ++i) {
237  /* could be a better test for discontinuity .. */
238  if (!dof_linkable(pf->dof_types()[i])) { discontinuous = true; break; }
239  }
240 
241  pfem classical_pf1 = discontinuous ? classical_discontinuous_fem(pgt, 1)
242  : classical_fem(pgt, 1);
243 
244  short_type degree = 1;
245  if ((pf != classical_pf1 && pf->estimated_degree() > 1) ||
246  pgt->structure() != pgt->basic_structure())
247  degree = 2;
248 
249  pmf->set_finite_element(cv, discontinuous ?
250  classical_discontinuous_fem(pgt, degree, 0, true) :
251  classical_fem(pgt, degree, true));
252  }
253  }
254  /* find out which dof will be exported to VTK */
255 
256  const mesh &m = pmf->linked_mesh();
257  pmf_mapping_type.resize(pmf->convex_index().last_true() + 1, unsigned(-1));
258  pmf_dof_used.sup(0, pmf->nb_basic_dof());
259  for (dal::bv_visitor cv(pmf->convex_index()); !cv.finished(); ++cv) {
260  vtk_mapping_type t = NO_VTK_MAPPING;
261  size_type nbd = pmf->fem_of_element(cv)->nb_dof(cv);
262  switch (pmf->fem_of_element(cv)->dim()) {
263  case 0: t = N1_TO_VTK_VERTEX; break;
264  case 1:
265  if (nbd == 2) t = N2_TO_VTK_LINE;
266  else if (nbd == 3) t = N3_TO_VTK_QUADRATIC_EDGE;
267  break;
268  case 2:
269  if (nbd == 3) t = N3_TO_VTK_TRIANGLE;
270  else if (nbd == 4)
271  t = check_voxel(m.points_of_convex(cv)) ? N4_TO_VTK_PIXEL
272  : N4_TO_VTK_QUAD;
273  else if (nbd == 6) t = N6_TO_VTK_QUADRATIC_TRIANGLE;
274  else if (nbd == 8) t = N8_TO_VTK_QUADRATIC_QUAD;
275  else if (nbd == 9) t = N9_TO_VTK_BIQUADRATIC_QUAD;
276  break;
277  case 3:
278  if (nbd == 4) t = N4_TO_VTK_TETRA;
279  else if (nbd == 10) t = N10_TO_VTK_QUADRATIC_TETRA;
280  else if (nbd == 8)
281  t = check_voxel(m.points_of_convex(cv)) ? N8_TO_VTK_VOXEL
282  : N8_TO_VTK_HEXAHEDRON;
283  else if (nbd == 20) t = N20_TO_VTK_QUADRATIC_HEXAHEDRON;
284  else if (nbd == 27) t = N27_TO_VTK_TRIQUADRATIC_HEXAHEDRON;
285  else if (nbd == 5) t = N5_TO_VTK_PYRAMID;
286  else if (nbd == 13) t = N13_TO_VTK_QUADRATIC_PYRAMID;
287  else if (nbd == 14) t = N14_TO_VTK_QUADRATIC_PYRAMID;
288  else if (nbd == 6) t = N6_TO_VTK_WEDGE;
289  else if (nbd == 15) t = N15_TO_VTK_QUADRATIC_WEDGE;
290  else if (nbd == 18) t = N18_TO_VTK_BIQUADRATIC_QUADRATIC_WEDGE;
291  break;
292  }
293  GMM_ASSERT1(t != -1, "semi internal error. Could not map " <<
294  name_of_fem(pmf->fem_of_element(cv))
295  << " to a VTK cell type");
296  pmf_mapping_type[cv] = t;
297 
298  const std::vector<unsigned> &dmap = select_vtk_dof_mapping(t);
299  //cout << "nbd = " << nbd << ", t = " << t << ", dmap = "<<dmap << "\n";
300  GMM_ASSERT1(dmap.size() <= pmf->nb_basic_dof_of_element(cv),
301  "inconsistency in vtk_dof_mapping");
302  for (unsigned i=0; i < dmap.size(); ++i)
303  pmf_dof_used.add(pmf->ind_basic_dof_of_element(cv)[dmap[i]]);
304  }
305  // cout << "mf.nb_dof = " << mf.nb_dof() << ", pmf->nb_dof="
306  // << pmf->nb_dof() << ", dof_used = " << pmf_dof_used.card() << "\n";
307  }
308 
309 
310  const stored_mesh_slice& vtk_export::get_exported_slice() const {
311  GMM_ASSERT1(psl, "no slice!");
312  return *psl;
313  }
314 
315  const mesh_fem& vtk_export::get_exported_mesh_fem() const {
316  GMM_ASSERT1(pmf.get(), "no mesh_fem!");
317  return *pmf;
318  }
319 
320  void vtk_export::set_header(const std::string& s)
321  {
322  strncpy(header, s.c_str(), 256);
323  header[255] = 0;
324  }
325 
326  void vtk_export::check_header() {
327  if (state >= HEADER_WRITTEN) return;
328  os << "# vtk DataFile Version 2.0\n";
329  os << header << "\n";
330  if (ascii) os << "ASCII\n"; else os << "BINARY\n";
331  state = HEADER_WRITTEN;
332  }
333 
334  void vtk_export::write_separ()
335  { if (ascii) os << "\n"; }
336 
337  void vtk_export::write_mesh() {
338  if (psl) write_mesh_structure_from_slice();
339  else write_mesh_structure_from_mesh_fem();
340  }
341 
342  /* export the slice data as an unstructured mesh composed of simplexes */
343  void vtk_export::write_mesh_structure_from_slice() {
344  /* element type code for (linear) simplexes of dimensions 0,1,2,3 in VTK */
345  static int vtk_simplex_code[4] = { VTK_VERTEX, VTK_LINE, VTK_TRIANGLE, VTK_TETRA };
346  if (state >= STRUCTURE_WRITTEN) return;
347  check_header();
348  /* possible improvement: detect structured grids */
349  os << "DATASET UNSTRUCTURED_GRID\n";
350  os << "POINTS " << psl->nb_points() << " float\n";
351  /*
352  points are not merge, vtk is mostly fine with that (except for
353  transparency and normals at vertices of 3D elements: all simplex faces
354  are considered to be "boundary" faces by vtk)
355  */
356  for (size_type ic=0; ic < psl->nb_convex(); ++ic) {
357  for (size_type i=0; i < psl->nodes(ic).size(); ++i)
358  write_vec(psl->nodes(ic)[i].pt.begin(),psl->nodes(ic)[i].pt.size());
359  write_separ();
360  }
361  /* count total number of simplexes, and total number of entries
362  * in the CELLS section */
363  size_type cells_cnt = 0, splx_cnt = 0;
364  for (size_type ic=0; ic < psl->nb_convex(); ++ic) {
365  for (size_type i=0; i < psl->simplexes(ic).size(); ++i)
366  cells_cnt += psl->simplexes(ic)[i].dim() + 2;
367  splx_cnt += psl->simplexes(ic).size();
368  }
369  size_type nodes_cnt = 0;
370  write_separ(); os << "CELLS " << splx_cnt << " " << cells_cnt << "\n";
371  for (size_type ic=0; ic < psl->nb_convex(); ++ic) {
372  const getfem::mesh_slicer::cs_simplexes_ct& s = psl->simplexes(ic);
373  for (size_type i=0; i < s.size(); ++i) {
374  write_val(int(s[i].dim()+1));
375  for (size_type j=0; j < s[i].dim()+1; ++j)
376  write_val(int(s[i].inodes[j] + nodes_cnt));
377  write_separ();
378  }
379  nodes_cnt += psl->nodes(ic).size();
380  }
381  assert(nodes_cnt == psl->nb_points()); // sanity check
382  write_separ(); os << "CELL_TYPES " << splx_cnt << "\n";
383  for (size_type ic=0; ic < psl->nb_convex(); ++ic) {
384  const getfem::mesh_slicer::cs_simplexes_ct& s = psl->simplexes(ic);
385  for (size_type i=0; i < s.size(); ++i) {
386  write_val(int(vtk_simplex_code[s[i].dim()]));
387  }
388  write_separ();
389  splx_cnt -= s.size();
390  }
391  assert(splx_cnt == 0); // sanity check
392  state = STRUCTURE_WRITTEN;
393  }
394 
395 
396  void vtk_export::write_mesh_structure_from_mesh_fem() {
397  if (state >= STRUCTURE_WRITTEN) return;
398  check_header();
399  /* possible improvement: detect structured grids */
400  os << "DATASET UNSTRUCTURED_GRID\n";
401  os << "POINTS " << pmf_dof_used.card() << " float\n";
402  std::vector<int> dofmap(pmf->nb_dof());
403  int cnt = 0;
404  for (dal::bv_visitor d(pmf_dof_used); !d.finished(); ++d) {
405  dofmap[d] = cnt++;
406  base_node P = pmf->point_of_basic_dof(d);
407  write_vec(P.const_begin(),P.size());
408  write_separ();
409  }
410 
411  size_type nb_cell_values = 0;
412  for (dal::bv_visitor cv(pmf->convex_index()); !cv.finished(); ++cv)
413  nb_cell_values += (1 + select_vtk_dof_mapping(pmf_mapping_type[cv]).size());
414 
415  write_separ(); os << "CELLS " << pmf->convex_index().card() << " " << nb_cell_values << "\n";
416 
417  for (dal::bv_visitor cv(pmf->convex_index()); !cv.finished(); ++cv) {
418  const std::vector<unsigned> &dmap = select_vtk_dof_mapping(pmf_mapping_type[cv]);
419  write_val(int(dmap.size()));
420  for (size_type i=0; i < dmap.size(); ++i)
421  write_val(int(dofmap[pmf->ind_basic_dof_of_element(cv)[dmap[i]]]));
422  write_separ();
423  }
424 
425  write_separ(); os << "CELL_TYPES " << pmf->convex_index().card() << "\n";
426  for (dal::bv_visitor cv(pmf->convex_index()); !cv.finished(); ++cv) {
427  write_val(select_vtk_type(pmf_mapping_type[cv]));
428  write_separ();
429  }
430 
431  state = STRUCTURE_WRITTEN;
432  }
433 
435  if (psl) {
436  mesh_fem mf(const_cast<mesh&>(m),1);
438  std::vector<scalar_type> q(mf.nb_dof());
439  for (size_type d=0; d < mf.nb_dof(); ++d) {
441  }
442  write_point_data(mf, q, "convex_quality");
443  } else {
444  std::vector<scalar_type> q(pmf->convex_index().card());
445  for (dal::bv_visitor cv(pmf->convex_index()); !cv.finished(); ++cv) {
446  q[cv] = m.convex_quality_estimate(cv);
447  }
448  write_cell_data(q, "convex_quality");
449  }
450  }
451 
452 
453  /* -------------------------------------------------------------
454  * OPENDX export
455  * ------------------------------------------------------------- */
456 
457  dx_export::dx_export(std::ostream &os_, bool ascii_)
458  : os(os_), ascii(ascii_) { init(); }
459 
460  dx_export::dx_export(const std::string& fname, bool ascii_, bool append_)
461  : os(real_os), ascii(ascii_) {
462  real_os.open(fname.c_str(),
463  std::ios_base::openmode(std::ios_base::in |
464  std::ios_base::out |
465  (append_ ? std::ios_base::ate : std::ios_base::trunc)));
466  GMM_ASSERT1(real_os.good(), "impossible to write to dx file '"
467  << fname << "'");
468  init();
469  if (append_) { reread_metadata(); header_written = true; }
470  }
471 
472  dx_export::~dx_export() {
473  std::ios::pos_type p = os.tellp();
474  write_series();
475  os << "\n# --end of getfem export\nend\n";
476  update_metadata(p);
477  }
478 
479  void dx_export::init() {
480  strcpy(header, "Exported by getfem++");
481  psl = 0; dim_ = dim_type(-1); connections_dim = dim_type(-1);
482  psl_use_merged = false;
483  header_written = false;
484  }
485 
486  void dx_export::write_separ()
487  { if (ascii) os << "\n"; }
488 
489  template<typename T> static typename std::list<T>::iterator
490  get_from_name(std::list<T> &c,
491  const std::string& name, bool raise_error) {
492  for (typename std::list<T>::iterator it = c.begin();
493  it != c.end(); ++it) {
494  if (it->name == name) return it;
495  }
496  GMM_ASSERT1(!raise_error, "object not found in dx file: " << name);
497  return c.end();
498  }
499 
500  std::list<dx_export::dxMesh>::iterator
501  dx_export::get_mesh(const std::string& name, bool raise_error) {
502  return get_from_name(meshes,name,raise_error);
503  }
504  std::list<dx_export::dxObject>::iterator
505  dx_export::get_object(const std::string& name, bool raise_error) {
506  return get_from_name(objects,name,raise_error);
507  }
508 
509 
510  bool dx_export::new_mesh(std::string &name) {
511  name = default_name(name, int(meshes.size()), "mesh");
512  std::list<dxMesh>::iterator it = get_mesh(name, false);
513  if (it != meshes.end()) {
514  if (&(*it) != &current_mesh())
515  std::swap(current_mesh(),*it);
516  return false;
517  } else {
518  meshes.push_back(dxMesh()); meshes.back().name = name;
519  return true;
520  }
521  }
522 
523  void dx_export::exporting(const stored_mesh_slice& sl, bool merge_points,
524  std::string name) {
525  if (!new_mesh(name)) return;
526  psl_use_merged = merge_points;
527  if (merge_points) sl.merge_nodes();
528  psl = &sl; dim_ = dim_type(sl.dim());
529  GMM_ASSERT1(psl->dim() <= 3, "4D slices and more are not supported");
530  for (dim_type d = 0; d <= psl->dim(); ++d) {
531  if (psl->nb_simplexes(d)) {
532  if (connections_dim == dim_type(-1)) connections_dim = d;
533  else GMM_ASSERT1(false, "Cannot export a slice containing "
534  "simplexes of different dimensions");
535  }
536  }
537  GMM_ASSERT1(connections_dim != dim_type(-1), "empty slice!");
538  }
539 
540 
541  void dx_export::exporting(const mesh_fem& mf, std::string name) {
542  name = default_name(name, int(meshes.size()), "mesh");
543  if (!new_mesh(name)) return;
544  const mesh &m = mf.linked_mesh();
545  GMM_ASSERT1(mf.linked_mesh().convex_index().card() != 0,
546  "won't export an empty mesh");
547 
548  dim_ = m.dim();
549  GMM_ASSERT1(dim_ <= 3, "4D meshes and more are not supported");
550  if (&mf != pmf.get())
551  pmf = std::make_unique<mesh_fem>(const_cast<mesh&>(m), dim_type(1));
552  bgeot::pgeometric_trans pgt = m.trans_of_convex(m.convex_index().first_true());
553  GMM_ASSERT1(dxname_of_convex_structure
554  (basic_structure(pgt->structure())) != 0,
555  "DX Cannot handle " <<
556  bgeot::name_of_geometric_trans(pgt) << ", use slices");
557  /* initialize pmf with finite elements suitable for OpenDX */
558  for (dal::bv_visitor cv(mf.convex_index()); !cv.finished(); ++cv) {
559  bgeot::pgeometric_trans pgt2 = mf.linked_mesh().trans_of_convex(cv);
560  GMM_ASSERT1(basic_structure(pgt->structure()) ==
561  basic_structure(pgt2->structure()),
562  "Cannot export this mesh to opendx, it contains "
563  "different convex types. Slice it first.");
564  pfem pf = mf.fem_of_element(cv);
565  bool discontinuous = false;
566  for (unsigned i=0; i < pf->nb_dof(cv); ++i) {
567  /* could be a better test for discontinuity .. */
568  if (!dof_linkable(pf->dof_types()[i])) { discontinuous = true; break; }
569  }
570  pfem classical_pf1 = discontinuous ? classical_discontinuous_fem(pgt, 1)
571  : classical_fem(pgt, 1);
572  pmf->set_finite_element(cv, classical_pf1);
573  }
574  pmf_dof_used.add(0, pmf->nb_basic_dof());
575  connections_dim = dim_type(pmf->nb_basic_dof_of_element(m.convex_index().first_true()));
576  }
577 
578  void dx_export::exporting(const mesh& m, std::string name) {
579  dim_ = m.dim();
580  GMM_ASSERT1(dim_ <= 3, "4D meshes and more are not supported");
581  pmf = std::make_unique<mesh_fem>(const_cast<mesh&>(m), dim_type(1));
582  pmf->set_classical_finite_element(1);
583  exporting(*pmf, name);
584  }
585 
586  void dx_export::write_series() {
587  for (std::list<dxSeries>::const_iterator it = series.begin();
588  it != series.end(); ++it) {
589  if (it->members.size() == 0) continue;
590  size_type count = 0;
591  os << "\nobject \"" << it->name << "\" class series\n";
592  for (std::list<std::string>::const_iterator ito = it->members.begin();
593  ito != it->members.end(); ++ito, ++count) {
594  os << " member " << count << " \"" << (*ito) << "\"\n";
595  }
596  }
597  }
598 
599  void dx_export::serie_add_object_(const std::string &serie_name,
600  const std::string &object_name) {
601  std::list<dxSeries>::iterator it = series.begin();
602  while (it != series.end() && it->name != serie_name) ++it;
603  if (it == series.end()) {
604  series.push_back(dxSeries()); it = series.end(); --it;
605  it->name = serie_name;
606  }
607  it->members.push_back(object_name);
608  }
609 
610  void dx_export::serie_add_object(const std::string &serie_name,
611  const std::string &object_name) {
612  /* create a series for edge data if possible (the cost is null
613  and it may be useful) */
614  std::list<dxObject>::iterator ito = get_object(object_name, false);
615  if (ito != objects.end()) {
616  std::list<dxMesh>::iterator itm = get_mesh(ito->mesh);
617  if (itm != meshes.end() && (itm->flags & dxMesh::WITH_EDGES)) {
618  serie_add_object_(serie_name + "_edges",
619  object_name + "_edges");
620  }
621  }
622  /* fill the real serie */
623  serie_add_object_(serie_name, object_name);
624  }
625 
626  void dx_export::set_header(const std::string& s)
627  { strncpy(header, s.c_str(), 256); header[255] = 0; }
628 
629  void dx_export::check_header() {
630  if (header_written) return;
631  header_written = true;
632  os << "# data file for IBM OpenDX, generated by GetFem++ v "
633  << GETFEM_VERSION << "\n";
634  os << "# " << header << "\n";
635  }
636 
637  void dx_export::update_metadata(std::ios::pos_type pos_series) {
638  os.seekp(0,std::ios::end);
639  os << "# This file contains the following objects\n";
640  std::ios::pos_type pos_end = os.tellp();
641  for (std::list<dxSeries>::const_iterator it = series.begin();
642  it != series.end(); ++it) {
643  os << "#S \"" << it->name << "\" which contains:\n";
644  for (std::list<std::string>::const_iterator its = it->members.begin();
645  its != it->members.end(); ++its)
646  os << "#+ \"" << *its << "\"\n";
647  }
648  for (std::list<dxObject>::const_iterator it = objects.begin();
649  it != objects.end(); ++it) {
650  os << "#O \"" << it->name << "\" \"" << it->mesh << "\"\n";
651  }
652  for (std::list<dxMesh>::const_iterator it = meshes.begin();
653  it != meshes.end(); ++it) {
654  os << "#M \"" << it->name << "\" " << it->flags << "\n";
655  }
656  os << "#E \"THE_END\" " << std::setw(20) << pos_series << std::setw(20) << pos_end << "\n";
657  }
658 
659  void dx_export::reread_metadata() {
660  char line[512];
661  real_os.seekg(0, std::ios::end);
662  int count=0; char c;
663  unsigned long lu_end, lu_series;
664  do {
665  real_os.seekg(-1, std::ios::cur);
666  c = char(real_os.peek());
667  } while (++count < 512 && c != '#');
668  real_os.getline(line, sizeof line);
669  if (sscanf(line, "#E \"THE_END\" %lu %lu", &lu_series, &lu_end) != 2)
670  GMM_ASSERT1(false, "this file was not generated by getfem, "
671  "cannot append data to it!\n");
672  real_os.seekg(lu_end, std::ios::beg);
673  do {
674  char name[512]; unsigned n;
675  int pos;
676  real_os.getline(line, sizeof line);
677  if (sscanf(line, "#%c \"%512[^\"]\"%n", &c, name, &pos) < 1)
678  GMM_ASSERT1(false, "corrupted file! your .dx file is broken\n");
679  if (c == 'S') {
680  series.push_back(dxSeries()); series.back().name = name;
681  } else if (c == '+') {
682  series.back().members.push_back(name);
683  } else if (c == 'O') {
684  objects.push_back(dxObject()); objects.back().name = name;
685  sscanf(line+pos, " \"%512[^\"]\"", name); objects.back().mesh = name;
686  } else if (c == 'M') {
687  meshes.push_back(dxMesh()); meshes.back().name = name;
688  sscanf(line+pos, "%u", &n); meshes.back().flags = n;
689  } else if (c == 'E') {
690  break;
691  } else GMM_ASSERT1(false, "corrupted file! your .dx file is broken\n");
692  } while (1);
693  real_os.seekp(lu_series, std::ios::beg);
694  }
695 
696  void dx_export::write_convex_attributes(bgeot::pconvex_structure cvs) {
697  const char *s_elem_type = dxname_of_convex_structure(cvs);
698  if (!s_elem_type)
699  GMM_WARNING1("OpenDX won't handle this kind of convexes");
700  os << "\n attribute \"element type\" string \"" << s_elem_type << "\"\n"
701  << " attribute \"ref\" string \"positions\"\n\n";
702  }
703 
704  const char *dx_export::dxname_of_convex_structure(bgeot::pconvex_structure cvs) {
705  const char *s_elem_type = 0;
706  switch (cvs->dim()) {
707  /* TODO: do something for point data */
708  case 1: s_elem_type = "lines"; break;
709  case 2:
710  if (cvs->nb_points() == 3)
711  s_elem_type = "triangles";
712  else if (cvs->nb_points() == 4)
713  s_elem_type = "quads";
714  break;
715  case 3:
716  if (cvs->nb_points() == 4)
717  s_elem_type = "tetrahedra";
718  else if (cvs->nb_points() == 8)
719  s_elem_type = "cubes";
720  break;
721  }
722  return s_elem_type;
723  }
724 
725  void dx_export::write_mesh() {
726  check_header();
727  if (current_mesh().flags & dxMesh::STRUCTURE_WRITTEN) return;
728  if (psl) write_mesh_structure_from_slice();
729  else write_mesh_structure_from_mesh_fem();
730 
731  os << "\nobject \"" << current_mesh_name() << "\" class field\n"
732  << " component \"positions\" value \""
733  << name_of_pts_array(current_mesh_name()) << "\"\n"
734  << " component \"connections\" value \""
735  << name_of_conn_array(current_mesh_name()) << "\"\n";
736  current_mesh().flags |= dxMesh::STRUCTURE_WRITTEN;
737  }
738 
739  /* export the slice data as an unstructured mesh composed of simplexes */
740  void dx_export::write_mesh_structure_from_slice() {
741  os << "\nobject \"" << name_of_pts_array(current_mesh_name())
742  << "\" class array type float rank 1 shape "
743  << int(psl->dim())
744  << " items " << (psl_use_merged ? psl->nb_merged_nodes() : psl->nb_points());
745  if (!ascii) os << " " << endianness() << " binary";
746  os << " data follows\n";
747  if (psl_use_merged) {
748  for (size_type i=0; i < psl->nb_merged_nodes(); ++i) {
749  for (size_type k=0; k < psl->dim(); ++k)
750  write_val(float(psl->merged_point(i)[k]));
751  write_separ();
752  }
753  } else {
754  for (size_type ic=0; ic < psl->nb_convex(); ++ic) {
755  for (size_type i=0; i < psl->nodes(ic).size(); ++i)
756  for (size_type k=0; k < psl->dim(); ++k)
757  write_val(float(psl->nodes(ic)[i].pt[k]));
758  write_separ();
759  }
760  }
761 
762  os << "\nobject \"" << name_of_conn_array(current_mesh_name())
763  << "\" class array type int rank 1 shape "
764  << int(connections_dim+1)
765  << " items " << psl->nb_simplexes(connections_dim);
766  if (!ascii) os << " " << endianness() << " binary";
767  os << " data follows\n";
768 
769  size_type nodes_cnt = 0; /* <- a virer , global_index le remplace */
770  for (size_type ic=0; ic < psl->nb_convex(); ++ic) {
771  const getfem::mesh_slicer::cs_simplexes_ct& s = psl->simplexes(ic);
772  for (size_type i=0; i < s.size(); ++i) {
773  if (s[i].dim() == connections_dim) {
774  for (size_type j=0; j < s[i].dim()+1; ++j) {
775  size_type k;
776  if (psl_use_merged)
777  k = psl->merged_index(ic, s[i].inodes[j]);
778  else k = psl->global_index(ic, s[i].inodes[j]);
779  write_val(int(k));
780  }
781  write_separ();
782  }
783  }
784  nodes_cnt += psl->nodes(ic).size();
785  }
786 
787  write_convex_attributes(bgeot::simplex_structure(connections_dim));
788  assert(nodes_cnt == psl->nb_points()); // sanity check
789  }
790 
791  void dx_export::write_mesh_structure_from_mesh_fem() {
792  os << "\nobject \"" << name_of_pts_array(current_mesh_name())
793  << "\" class array type float rank 1 shape "
794  << int(pmf->linked_mesh().dim())
795  << " items " << pmf->nb_dof();
796  if (!ascii) os << " " << endianness() << " binary";
797  os << " data follows\n";
798 
799  /* possible improvement: detect structured grids */
800  for (size_type d = 0; d < pmf->nb_basic_dof(); ++d) {
801  const base_node P = pmf->point_of_basic_dof(d);
802  for (size_type k=0; k < dim_; ++k)
803  write_val(float(P[k]));
804  write_separ();
805  }
806 
807  os << "\nobject \"" << name_of_conn_array(current_mesh_name())
808  << "\" class array type int rank 1 shape "
809  << int(connections_dim)
810  << " items " << pmf->convex_index().card();
811  if (!ascii) os << " " << endianness() << " binary";
812  os << " data follows\n";
813 
814  for (dal::bv_visitor cv(pmf->convex_index()); !cv.finished(); ++cv) {
815  for (size_type i=0; i < connections_dim; ++i)
816  write_val(int(pmf->ind_basic_dof_of_element(cv)[i]));
817  write_separ();
818  }
819  write_convex_attributes(basic_structure(pmf->linked_mesh().structure_of_convex(pmf->convex_index().first_true())));
820  }
821 
822  void dx_export::exporting_mesh_edges(bool with_slice_edges) {
823  write_mesh();
824  if (current_mesh().flags & dxMesh::WITH_EDGES) return;
825  if (psl) write_mesh_edges_from_slice(with_slice_edges);
826  else write_mesh_edges_from_mesh();
827  current_mesh().flags |= dxMesh::WITH_EDGES;
828  os << "\nobject \"" << name_of_edges_array(current_mesh_name())
829  << "\" class field\n"
830  << " component \"positions\" value \""
831  << name_of_pts_array(current_mesh_name()) << "\"\n"
832  << " component \"connections\" value \""
833  << name_of_conn_array(name_of_edges_array(current_mesh_name()))
834  << "\"\n";
835  }
836 
837  void dx_export::write_mesh_edges_from_slice(bool with_slice_edges) {
838  std::vector<size_type> edges;
839  dal::bit_vector slice_edges;
840  psl->get_edges(edges, slice_edges, psl_use_merged);
841  if (with_slice_edges) slice_edges.clear();
842  os << "\nobject \""
843  << name_of_conn_array(name_of_edges_array(current_mesh_name()))
844  << "\" class array type int rank 1 shape 2"
845  << " items " << edges.size()/2 - slice_edges.card();
846  if (!ascii) os << " " << endianness() << " binary";
847  os << " data follows\n";
848  for (size_type i=0; i < edges.size()/2; ++i) {
849  if (!slice_edges.is_in(i)) {
850  write_val(int(edges[2*i]));
851  write_val(int(edges[2*i+1]));
852  }
853  if ((i+1)%10 == 0) write_separ();
854  }
855  write_separ();
856  write_convex_attributes(bgeot::simplex_structure(1));
857  }
858 
859  void dx_export::write_mesh_edges_from_mesh() {
860  bgeot::mesh_structure ms(pmf->linked_mesh()); ms.to_edges();
861  os << "\nobject \""
862  << name_of_conn_array(name_of_edges_array(current_mesh_name()))
863  << "\" class array type int rank 1 shape 2"
864  << " items " << ms.convex_index().card();
865  if (!ascii) os << " " << endianness() << " binary";
866  os << " data follows\n";
867  for (dal::bv_visitor cv(ms.convex_index()); !cv.finished(); ++cv) {
868  write_val(int(ms.ind_points_of_convex(cv)[0]));
869  write_val(int(ms.ind_points_of_convex(cv)[1]));
870  if ((cv+1)%20 == 0) write_separ();
871  }
872  write_separ();
873  write_convex_attributes(bgeot::simplex_structure(1));
874  }
875 
876 
877  /* -------------------------------------------------------------
878  * POS export (Gmsh post-processing format)
879  * ------------------------------------------------------------- */
880  struct gf2pos_dof_mapping : public std::vector<std::vector<unsigned> > {};
881 
882  static const std::vector<unsigned>& getfem_to_pos_dof_mapping(int t) {
883  gf2pos_dof_mapping &dm = dal::singleton<gf2pos_dof_mapping>::instance();
884  if (dm.size() == 0) {
885  dm.resize(8);
886  dm[pos_export::POS_PT] = {0};
887  dm[pos_export::POS_LN] = {0, 1};
888  dm[pos_export::POS_TR] = {0, 1, 2};
889  dm[pos_export::POS_QU] = {0, 1, 3, 2};
890  dm[pos_export::POS_SI] = {0, 1, 2, 3};
891  dm[pos_export::POS_HE] = {0, 1, 3, 2, 4, 5, 7, 6};
892  dm[pos_export::POS_PR] = {0, 1, 2, 3, 4, 5};
893  dm[pos_export::POS_PY] = {0, 1, 3, 2, 4};
894  }
895  return dm[t];
896  }
897 
898  pos_export::pos_export(std::ostream& osname): os(osname) {
899  init();
900  }
901 
902  pos_export::pos_export(const std::string& fname)
903  : os(real_os), real_os(fname.c_str()) {
904  GMM_ASSERT1(real_os, "impossible to write to pos file '" << fname << "'");
905  init();
906  }
907 
908  void pos_export::init() {
909  strcpy(header, "Exported by GetFEM++");
910  state = EMPTY;
911  view = 0;
912  }
913 
914  void pos_export::set_header(const std::string& s){
915  strncpy(header, s.c_str(), 256);
916  header[255] = 0;
917  }
918 
919  void pos_export::check_header() {
920  if (state >= HEADER_WRITTEN) return;
921  os << "/* " << header << " */\n";
922  os << "General.FastRedraw = 0;\n";
923  os << "General.ColorScheme = 1;\n";
924  state = HEADER_WRITTEN;
925  }
926 
927  void pos_export::exporting(const mesh& m) {
928  if (state >= STRUCTURE_WRITTEN) return;
929  dim = dim_type(m.dim());
930  GMM_ASSERT1(int(dim) <= 3, "attempt to export a "
931  << int(dim) << "D mesh (not supported)");
932  pmf = std::make_unique<mesh_fem>(const_cast<mesh&>(m), dim_type(1));
933  pmf->set_classical_finite_element(1);
934  exporting(*pmf);
935  state = STRUCTURE_WRITTEN;
936  }
937 
938  void pos_export::exporting(const mesh_fem& mf) {
939  if (state >= STRUCTURE_WRITTEN) return;
940  dim = dim_type(mf.linked_mesh().dim());
941  GMM_ASSERT1(int(dim) <= 3, "attempt to export a "
942  << int(dim) << "D mesh_fem (not supported)");
943  if (&mf != pmf.get())
944  pmf = std::make_unique<mesh_fem>(mf.linked_mesh(), dim_type(1));
945 
946  /* initialize pmf with finite elements suitable for Gmsh */
947  for (dal::bv_visitor cv(mf.convex_index()); !cv.finished(); ++cv) {
948  bgeot::pgeometric_trans pgt = mf.linked_mesh().trans_of_convex(cv);
949  pfem pf = mf.fem_of_element(cv);
950 
951  bool discontinuous = false;
952  for (unsigned i=0; i < pf->nb_dof(cv); ++i) {
953  // could be a better test for discontinuity ...
954  if (!dof_linkable(pf->dof_types()[i])) { discontinuous = true; break; }
955  }
956  pfem classical_pf1 = discontinuous ?
957  classical_discontinuous_fem(pgt, 1) : classical_fem(pgt, 1);
958  pmf->set_finite_element(cv, classical_pf1);
959  }
960  psl = NULL;
961 
962  /* find out which dof will be exported to Gmsh */
963  for (dal::bv_visitor cv(pmf->convex_index()); !cv.finished(); ++cv) {
964  int t = -1;
965  switch (pmf->fem_of_element(cv)->nb_dof(cv)){
966  case 1: t = POS_PT; break;
967  case 2: t = POS_LN; break;
968  case 3: t = POS_TR; break;
969  case 4:
970  if ( 2 == pmf->fem_of_element(cv)->dim() ) t = POS_QU;
971  else if (3 == pmf->fem_of_element(cv)->dim()) t = POS_SI;
972  break;
973  case 6: t = POS_PR; break;
974  case 8: t = POS_HE; break;
975  case 5: t = POS_PY; break;
976  }
977  GMM_ASSERT1(t != -1, "semi internal error. Could not map "
978  << name_of_fem(pmf->fem_of_element(cv))
979  << " to a POS primitive type");
980  pos_cell_type.push_back(unsigned(t));
981 
982  const std::vector<unsigned>& dmap = getfem_to_pos_dof_mapping(t);
983  GMM_ASSERT1(dmap.size() <= pmf->nb_basic_dof_of_element(cv),
984  "inconsistency in pos_dof_mapping");
985  std::vector<unsigned> cell_dof;
986  cell_dof.resize(dmap.size(),unsigned(-1));
987  for (size_type i=0; i < dmap.size(); ++i)
988  cell_dof[i] = unsigned(pmf->ind_basic_dof_of_element(cv)[dmap[i]]);
989  pos_cell_dof.push_back(cell_dof);
990  }
991  for (size_type i=0; i< pmf->nb_basic_dof(); ++i){
992  std::vector<float> pt;
993  pt.resize(dim,float(0));
994  for (size_type j=0; j<dim; ++j)
995  pt[j] = float(pmf->point_of_basic_dof(i)[j]);
996  pos_pts.push_back(pt);
997  }
998  state = STRUCTURE_WRITTEN;
999  }
1000 
1001  void pos_export::exporting(const stored_mesh_slice& sl) {
1002  if (state >= STRUCTURE_WRITTEN) return;
1003  psl = &sl;
1004  dim = dim_type(sl.dim());
1005  GMM_ASSERT1(int(dim) <= 3, "attempt to export a "
1006  << int(dim) << "D slice (not supported)");
1007 
1008  for (size_type ic=0, pcnt=0; ic < psl->nb_convex(); ++ic) {
1009  for (getfem::mesh_slicer::cs_simplexes_ct::const_iterator it=psl->simplexes(ic).begin();
1010  it != psl->simplexes(ic).end(); ++it) {
1011  int t = -1;
1012  switch ((*it).dim()){
1013  case 0: t = POS_PT; break;
1014  case 1: t = POS_LN; break;
1015  case 2: t = POS_TR; break;
1016  case 3: t = POS_SI; break;
1017  }
1018  GMM_ASSERT1(t != -1, "semi internal error.. could not map to a POS primitive type");
1019  pos_cell_type.push_back(unsigned(t));
1020 
1021  const std::vector<unsigned>& dmap = getfem_to_pos_dof_mapping(t);
1022  GMM_ASSERT1(dmap.size() <= size_type(t+1), "inconsistency in pos_dof_mapping");
1023 
1024  std::vector<unsigned> cell_dof;
1025  cell_dof.resize(dmap.size(),unsigned(-1));
1026  for (size_type i=0; i < dmap.size(); ++i)
1027  cell_dof[i] = unsigned(it->inodes[dmap[i]] + pcnt);
1028  pos_cell_dof.push_back(cell_dof);
1029  }
1030  for(getfem::mesh_slicer::cs_nodes_ct::const_iterator it=psl->nodes(ic).begin();
1031  it != psl->nodes(ic).end(); ++it) {
1032  std::vector<float> pt;
1033  pt.resize(dim,float(0));
1034  for (size_type i=0; i<dim; ++i)
1035  pt[i] = float(it->pt[i]);
1036  pos_pts.push_back(pt);
1037  }
1038  pcnt += psl->nodes(ic).size();
1039  }
1040  state = STRUCTURE_WRITTEN;
1041  }
1042 
1043  void pos_export::write(const mesh& m, const std::string &name){
1044  if (state >= IN_CELL_DATA) return;
1045  GMM_ASSERT1(int(m.dim()) <= 3, "attempt to export a "
1046  << int(m.dim()) << "D mesh (not supported)");
1047  pmf = std::make_unique<mesh_fem>(const_cast<mesh&>(m), dim_type(1));
1048  pmf->set_classical_finite_element(1);
1049  write(*pmf,name);
1050  state = IN_CELL_DATA;
1051  }
1052 
1053  void pos_export::write(const mesh_fem& mf, const std::string &name){
1054  if (state >= IN_CELL_DATA) return;
1055  check_header();
1056  exporting(mf);
1057 
1058  if (""==name) os << "View \"mesh " << view <<"\" {\n";
1059  else os << "View \"" << name <<"\" {\n";
1060 
1061  int t;
1062  std::vector<unsigned> cell_dof;
1063  std::vector<float> cell_dof_val;
1064  for (size_type cell = 0; cell < pos_cell_type.size(); ++cell) {
1065  t = pos_cell_type[cell];
1066  cell_dof = pos_cell_dof[cell];
1067  cell_dof_val.resize(cell_dof.size(),float(0));
1068  write_cell(t,cell_dof,cell_dof_val);
1069  }
1070 
1071  os << "};\n";
1072  os << "View[" << view << "].ShowScale = 0;\n";
1073  os << "View[" << view << "].ShowElement = 1;\n";
1074  os << "View[" << view << "].DrawScalars = 0;\n";
1075  os << "View[" << view << "].DrawVectors = 0;\n";
1076  os << "View[" << view++ << "].DrawTensors = 0;\n";
1077  state = IN_CELL_DATA;
1078  }
1079 
1080  void pos_export::write(const stored_mesh_slice& sl, const std::string &name){
1081  if (state >= IN_CELL_DATA) return;
1082  check_header();
1083  exporting(sl);
1084 
1085  if (""==name) os << "View \"mesh " << view <<"\" {\n";
1086  else os << "View \"" << name <<"\" {\n";
1087 
1088  int t;
1089  std::vector<unsigned> cell_dof;
1090  std::vector<float> cell_dof_val;
1091  for (size_type cell = 0; cell < pos_cell_type.size(); ++cell) {
1092  t = pos_cell_type[cell];
1093  cell_dof = pos_cell_dof[cell];
1094  cell_dof_val.resize(cell_dof.size(),float(0));
1095  write_cell(t,cell_dof,cell_dof_val);
1096  }
1097 
1098  os << "};\n";
1099  os << "View[" << view << "].ShowScale = 0;\n";
1100  os << "View[" << view << "].ShowElement = 1;\n";
1101  os << "View[" << view << "].DrawScalars = 0;\n";
1102  os << "View[" << view << "].DrawVectors = 0;\n";
1103  os << "View[" << view++ << "].DrawTensors = 0;\n";
1104  state = IN_CELL_DATA;
1105  }
1106 } /* end of namespace getfem. */
virtual size_type nb_dof() const
Return the total number of degrees of freedom.
void set_header(const std::string &s)
the header is the second line of text in the exported file, you can put whatever you want – call thi...
void set_header(const std::string &s)
the header is the second line of text in the exported file, you can put whatever you want – call thi...
const mesh_slicer::cs_nodes_ct & nodes(size_type ic) const
Return the list of nodes for the &#39;ic&#39;th convex of the slice.
Export solutions to various formats.
void to_edges()
build a new mesh, such that its convexes are the edges of the convexes of the previous one ...
size_type nb_points() const
Return the number of nodes in the slice.
size_type nb_merged_nodes() const
Return the number of merged nodes in slice.
pfem fem_descriptor(const std::string &name)
get a fem descriptor from its string name.
Definition: getfem_fem.cc:3950
void exporting_mesh_edges(bool with_slice_edge=true)
append edges information (if you want to draw the mesh and are using a refined slice.
virtual pfem fem_of_element(size_type cv) const
Return the basic fem associated with an element (if no fem is associated, the function will crash! us...
Describe a mesh (collection of convexes (elements) and points).
Definition: getfem_mesh.h:95
pfem classical_discontinuous_fem(bgeot::pgeometric_trans pg, short_type k, scalar_type alpha=0, bool complete=false)
Give a pointer on the structures describing the classical polynomial discontinuous fem of degree k on...
Definition: getfem_fem.cc:3872
void serie_add_object(const std::string &serie_name, const std::string &object_name)
add an object (typically the name of a data field) to a &#39;series&#39;, i.e.
void exporting(const mesh &m)
should be called before write_*_data
size_type dim() const
return the slice dimension
static T & instance()
Instance from the current thread.
const base_node merged_point(size_type i_merged) const
Return the physical position of the merged node.
void nb_simplexes(std::vector< size_type > &c) const
return the simplex count, in an array.
std::shared_ptr< const convex_structure > pconvex_structure
Pointer on a convex structure description.
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:49
void write_point_data(const getfem::mesh_fem &mf, const VECT &U0, const std::string &name)
append a new scalar or vector field defined on mf to the .vtk file.
scalar_type convex_quality_estimate(size_type ic) const
Return an estimate of the convex quality (0 <= Q <= 1).
Definition: getfem_mesh.cc:409
std::string name_of_fem(pfem p)
get the string name of a fem descriptor.
Definition: getfem_fem.cc:3959
size_type nb_convex() const
return the number of convexes of the original mesh referenced in the slice
A simple singleton implementation.
void get_edges(std::vector< size_type > &edges, dal::bit_vector &slice_edges, bool from_merged_nodes) const
Extract the list of mesh edges.
void merge_nodes() const
build a list of merged nodes.
void set_classical_finite_element(size_type cv, dim_type fem_degree, bool complete=false)
Set a classical (i.e.
The output of a getfem::mesh_slicer which has been recorded.
GEneric Tool for Finite Element Methods.
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
Definition: getfem_fem.h:239
bool dof_linkable(pdof_description)
Says if the dof is linkable.
Definition: getfem_fem.cc:594
Describe a finite element method linked to a mesh.
gmm::uint16_type short_type
used as the common short type integer in the library
Definition: bgeot_config.h:79
std::string name_of_geometric_trans(pgeometric_trans p)
Get the string name of a geometric transformation.
virtual size_type first_convex_of_basic_dof(size_type d) const
Shortcut for convex_to_dof(d)[0].
convenient initialization of vectors via overload of "operator,".
pconvex_structure simplex_structure(dim_type nc)
Give a pointer on the structures of a simplex of dimension d.
const mesh_slicer::cs_simplexes_ct & simplexes(size_type ic) const
Return the list of simplexes for the &#39;ic&#39;th convex of the slice.
pconvex_structure basic_structure(pconvex_structure cv)
Original structure (if concerned)
Mesh structure definition.
void write_mesh_quality(const mesh &m)
export a data_set correspounding to measures of quality for each convex of the supplied mesh (which s...
pfem classical_fem(bgeot::pgeometric_trans pgt, short_type k, bool complete=false)
Give a pointer on the structures describing the classical polynomial fem of degree k on a given conve...
Definition: getfem_fem.cc:3867
const mesh & linked_mesh() const
Return a reference to the underlying mesh.
const dal::bit_vector & convex_index() const
Get the set of convexes where a finite element has been assigned.
void write_cell_data(const VECT &U, const std::string &name, size_type qdim=1)
export data which is constant over each element.
const dal::bit_vector & convex_index() const
Return the list of valid convex IDs.
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation