GetFEM++  5.3
getfem_export.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2001-2017 Yves Renard, Julien Pommier
5 
6  This file is a part of GetFEM++
7 
8  GetFEM++ is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program; if not, write to the Free Software Foundation,
19  Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
20 
21  As a special exception, you may use this file as it is a part of a free
22  software library without restriction. Specifically, if other files
23  instantiate templates or use macros or inline functions from this file,
24  or you compile this file and link it with other files to produce an
25  executable, this file does not by itself cause the resulting executable
26  to be covered by the GNU Lesser General Public License. This exception
27  does not however invalidate any other reasons why the executable file
28  might be covered by the GNU Lesser General Public License.
29 
30 ===========================================================================*/
31 
32 /**@file getfem_export.h
33  @author Yves Renard <Yves.Renard@insa-lyon.fr>,
34  @author Julien Pommier <Julien.Pommier@insa-toulouse.fr>
35  @date October 15, 2001.
36  @brief Export solutions to various formats.
37 */
38 #ifndef GETFEM_EXPORT_H__
39 #define GETFEM_EXPORT_H__
40 
41 #include "getfem_interpolation.h"
42 #include "getfem_mesh_slice.h"
43 #include <list>
44 
45 namespace getfem {
46 
47  /* ********************************************************************* */
48  /* */
49  /* Save a solution in a file with a Pk interpolation. */
50  /* */
51  /* ********************************************************************* */
52 
53  inline std::string remove_spaces(const std::string &s) {
54  std::string s2(s);
55  for (unsigned i=0; i < s.size(); ++i)
56  if (s2[i] <= ' ') s2[i] = '_';
57  return s2;
58  }
59 
60  /** @brief VTK export.
61 
62  export class to VTK ( http://www.kitware.com/vtk.html ) file format
63  (not the XML format, but the old format)
64 
65  A vtk_export can store multiple scalar/vector fields.
66  */
67  class vtk_export {
68  protected:
69  std::ostream &os;
70  char header[256]; // hard limit in vtk
71  bool ascii;
72  const stored_mesh_slice *psl;
73  std::unique_ptr<mesh_fem> pmf;
74  dal::bit_vector pmf_dof_used;
75  std::vector<unsigned> pmf_mapping_type;
76  std::ofstream real_os;
77  dim_type dim_;
78  bool reverse_endian;
79  enum { EMPTY, HEADER_WRITTEN, STRUCTURE_WRITTEN, IN_CELL_DATA,
80  IN_POINT_DATA } state;
81 
82  template<class T> void write_val(T v);
83  template<class V> void write_vec(V p, size_type qdim);
84  template<class IT> void write_3x3tensor(IT p);
85  void write_separ();
86 
87  public:
88  typedef enum { VTK_VERTEX = 1,
89  VTK_LINE = 3,
90  VTK_TRIANGLE = 5,
91  VTK_PIXEL = 8,
92  VTK_QUAD = 9,
93  VTK_TETRA = 10,
94  VTK_VOXEL = 11,
95  VTK_HEXAHEDRON = 12,
96  VTK_WEDGE = 13,
97  VTK_PYRAMID = 14,
98  VTK_QUADRATIC_EDGE = 21,
99  VTK_QUADRATIC_TRIANGLE = 22,
100  VTK_QUADRATIC_QUAD = 23,
101  VTK_QUADRATIC_TETRA = 24,
102  VTK_QUADRATIC_HEXAHEDRON = 25,
103  VTK_QUADRATIC_WEDGE = 26,
104  VTK_QUADRATIC_PYRAMID = 27,
105  VTK_BIQUADRATIC_QUAD = 28,
106  VTK_TRIQUADRATIC_HEXAHEDRON = 29,
107  VTK_BIQUADRATIC_QUADRATIC_WEDGE = 32 } vtk_cell_type;
108  vtk_export(const std::string& fname, bool ascii_ = false);
109  vtk_export(std::ostream &os_, bool ascii_ = false);
110 
111  /** should be called before write_*_data */
112  void exporting(const mesh& m);
113  void exporting(const mesh_fem& mf);
114  void exporting(const stored_mesh_slice& sl);
115 
116  /** the header is the second line of text in the exported file,
117  you can put whatever you want -- call this before any write_dataset
118  or write_mesh */
119  void set_header(const std::string& s);
120  void write_mesh();
121 
122  /** append a new scalar or vector field defined on mf to the .vtk file. If
123  you are exporting a slice, or if mf != get_exported_mesh_fem(), U will
124  be interpolated on the slice, or on get_exported_mesh_fem().
125 
126  Note that vectors should be written AFTER scalars, and tensors
127  after vectors
128 
129  NO SPACE ALLOWED in 'name' */
130  template<class VECT> void write_point_data(const getfem::mesh_fem &mf,
131  const VECT& U0,
132  const std::string& name);
133 
134  /** append a new scalar or vector field to .vtk file. The Uslice vector is
135  the field interpolated on the exported mesh_slice This function should
136  not be used if you are not exporting a slice! NO SPACE ALLOWED in
137  'name' */
138  template<class VECT> void write_sliced_point_data(const VECT& Uslice,
139  const std::string& name,
140  size_type qdim=1);
141  /** export data which is constant over each element. You should not use
142  this function if you are exporting a slice. U should have
143  convex_index().card() elements. */
144 
145  template<class VECT> void write_cell_data(const VECT& U,
146  const std::string& name,
147  size_type qdim = 1);
148  /** export a data_set correspounding to measures of quality for each convex
149  of the supplied mesh (which should have the same number of convex than
150  the one used in the vtk_export)
151 
152  If a slice is being exported, the convex quality is written as
153  point_data (TO IMPROVE ONEDAY), if a mesh/mesh_fem is being exported,
154  it is written as cell_data
155  */
156  void write_mesh_quality(const mesh &m);
157  void write_normals();
158  const stored_mesh_slice& get_exported_slice() const;
159  const mesh_fem& get_exported_mesh_fem() const;
160  private:
161  void init();
162  void check_header();
163  void write_mesh_structure_from_slice();
164  void write_mesh_structure_from_mesh_fem();
165  void switch_to_cell_data();
166  void switch_to_point_data();
167  template<class VECT> void write_dataset_(const VECT& U,
168  const std::string& name,
169  size_type qdim,
170  bool cell_data=false);
171  };
172 
173  template<class T> void vtk_export::write_val(T v) {
174  if (ascii) os << " " << v;
175  else {
176  char *p = (char*)&v;
177  if (reverse_endian)
178  for (size_type i=0; i < sizeof(v)/2; ++i)
179  std::swap(p[i], p[sizeof(v)-i-1]);
180  os.write(p, sizeof(T));
181  }
182  }
183 
184  template<class IT> void vtk_export::write_vec(IT p, size_type qdim) {
185  float v[3];
186  for (size_type i=0; i < qdim; ++i) {
187  v[i] = float(p[i]);
188  }
189  for (size_type i=qdim; i < 3; ++i) v[i] = 0.0f;
190  write_val(v[0]);write_val(v[1]);write_val(v[2]);
191  }
192 
193  template<class IT> void vtk_export::write_3x3tensor(IT p) {
194  float v[3][3];
195  memset(v, 0, sizeof v);
196  for (size_type i=0; i < dim_; ++i) {
197  for (size_type j=0; j < dim_; ++j)
198  v[i][j] = float(p[i + j*dim_]);
199  }
200  for (size_type i=0; i < 3; ++i) {
201  for (size_type j=0; j < 3; ++j) {
202  write_val(v[i][j]);
203  }
204  if (ascii) os << "\n";
205  }
206  }
207 
208  template<class VECT>
209  void vtk_export::write_point_data(const getfem::mesh_fem &mf, const VECT& U,
210  const std::string& name) {
211  size_type Q = (gmm::vect_size(U) / mf.nb_dof()) * mf.get_qdim();
212  size_type qdim = mf.get_qdim();
213  if (psl) {
214  std::vector<scalar_type> Uslice(Q*psl->nb_points());
215  psl->interpolate(mf, U, Uslice);
216  write_dataset_(Uslice, name, qdim);
217  } else {
218  std::vector<scalar_type> V(pmf->nb_dof() * Q);
219  if (&mf != &(*pmf)) {
220  interpolation(mf, *pmf, U, V);
221  } else gmm::copy(U,V);
222  size_type cnt = 0;
223  for (dal::bv_visitor d(pmf_dof_used); !d.finished(); ++d, ++cnt) {
224  if (cnt != d)
225  for (size_type q=0; q < Q; ++q) {
226  V[cnt*Q + q] = V[d*Q + q];
227  }
228  }
229  V.resize(Q*pmf_dof_used.card());
230  write_dataset_(V, name, qdim);
231  }
232  }
233 
234  template<class VECT>
235  void vtk_export::write_cell_data(const VECT& U, const std::string& name,
236  size_type qdim) {
237  write_dataset_(U, name, qdim, true);
238  }
239 
240  template<class VECT>
242  const std::string& name,
243  size_type qdim) {
244  write_dataset_(U, name, qdim, false);
245  }
246 
247  template<class VECT>
248  void vtk_export::write_dataset_(const VECT& U, const std::string& name,
249  size_type qdim, bool cell_data) {
250  write_mesh();
251  size_type nb_val = 0;
252  if (cell_data) {
253  switch_to_cell_data();
254  nb_val = psl ? psl->linked_mesh().convex_index().card()
255  : pmf->linked_mesh().convex_index().card();
256  } else {
257  switch_to_point_data();
258  nb_val = psl ? psl->nb_points() : pmf_dof_used.card();
259  }
260  size_type Q = qdim;
261  if (Q == 1) Q = gmm::vect_size(U) / nb_val;
262  GMM_ASSERT1(gmm::vect_size(U) == nb_val*Q,
263  "inconsistency in the size of the dataset: "
264  << gmm::vect_size(U) << " != " << nb_val << "*" << Q);
265  write_separ();
266  if (Q == 1) {
267  os << "SCALARS " << remove_spaces(name) << " float 1\n";
268  os << "LOOKUP_TABLE default\n";
269  for (size_type i=0; i < nb_val; ++i) {
270  write_val(float(U[i]));
271  }
272  } else if (Q <= 3) {
273  os << "VECTORS " << remove_spaces(name) << " float\n";
274  for (size_type i=0; i < nb_val; ++i) {
275  write_vec(U.begin() + i*Q, Q);
276  }
277  } else if (Q == gmm::sqr(dim_)) {
278  /* tensors : coef are supposed to be stored in FORTRAN order
279  in the VTK file, they are written with C (row major) order
280  */
281  os << "TENSORS " << remove_spaces(name) << " float\n";
282  for (size_type i=0; i < nb_val; ++i) {
283  write_3x3tensor(U.begin() + i*Q);
284  }
285  } else GMM_ASSERT1(false, "vtk does not accept vectors of dimension > 3");
286  write_separ();
287  }
288 
289 
290  /** @brief A (quite large) class for exportation of data to IBM OpenDX.
291 
292  http://www.opendx.org/
293 
294  This class is more capable than the VTK export, as it is
295  possible to export many different meshes/slices, with their
296  edges, datasets, and create series of dataset for animations
297  etc, in a single '.dx' file.
298 
299  Moreover, it is able to reopen a '.dx' file and append new data
300  into it. Hence it is possible, if many time-steps are to be
301  saved, to view intermediate results in OpenDX during the
302  computation.
303  */
304  class dx_export {
305  std::ostream &os;
306  char header[256];
307  bool ascii;
308  const stored_mesh_slice *psl;
309  bool psl_use_merged; /* flag enabled if we merge the points of
310  psl before export */
311  std::unique_ptr<mesh_fem> pmf;
312  dal::bit_vector pmf_dof_used;
313  std::vector<unsigned> pmf_cell_type;
314  std::fstream real_os;
315  dim_type dim_, connections_dim;
316  struct dxSeries {
317  std::string name;
318  std::list<std::string> members;
319  };
320  struct dxObject {
321  std::string name;
322  std::string mesh;
323  };
324  struct dxMesh {
325  unsigned flags;
326  typedef enum { NONE=0, WITH_EDGES=1, STRUCTURE_WRITTEN=2 } flags_t;
327  std::string name;
328  dxMesh() : flags(NONE) {}
329  };
330  std::list<dxObject> objects;
331  std::list<dxSeries> series;
332  std::list<dxMesh> meshes;
333  bool header_written;
334  public:
335  dx_export(const std::string& fname, bool ascii_ = false,
336  bool append_ = false);
337  dx_export(std::ostream &os_, bool ascii_ = false);
338  ~dx_export(); /* the file is not complete until the destructor
339  has been executed */
340  void exporting(const mesh& m, std::string name = std::string());
341  void exporting(const mesh_fem& mf, std::string name = std::string());
342  void exporting(const stored_mesh_slice& sl, bool merge_points = true,
343  std::string name = std::string());
344  /** append edges information (if you want to draw the mesh and are
345  using a refined slice. Should be called just after exporting(..) */
346  void exporting_mesh_edges(bool with_slice_edge = true);
347 
348  /** the header is the second line of text in the exported file,
349  you can put whatever you want -- call this before any write_dataset
350  or write_mesh */
351  void set_header(const std::string& s);
352  void write_mesh();
353  /** add an object (typically the name of a data field) to a
354  'series', i.e. an aggregate of consecutive objects. Using
355  'series' is useful for animations in opendx
356 
357  If 'field_name' corresponds to a data_set whose mesh edges have
358  been exported, a second series called serie_name + '_edges'
359  will be filled, which will allow you to view the mesh edges.
360  */
361  void serie_add_object(const std::string& serie_name,
362  const std::string& object_name);
363  void serie_add_object(const std::string& serie_name)
364  { serie_add_object(serie_name, current_data_name()); }
365  /** return the name of current mesh (use exporting(...) to change
366  the current mesh) */
367  std::string current_mesh_name() { return current_mesh().name; }
368  /** return the name of last written data_set */
369  std::string current_data_name() { return current_data().name; }
370  template<class VECT> void
372  const VECT& U0, std::string name = std::string());
373  template<class VECT> void
374  write_sliced_point_data(const VECT& Uslice,
375  std::string name = std::string());
376  /* TOBEDONE !!!!!!!!!!!
377  template<class VECT> void
378  write_cell_data(const VECT& U, std::string name = std::string());
379  void write_mesh_quality(const mesh &m);*/
380  void write_normals();
381  const stored_mesh_slice& get_exported_slice() const;
382  const mesh_fem& get_exported_mesh_fem() const;
383 
384  private:
385  void init();
386  void reread_metadata();
387  void update_metadata(std::ios::pos_type);
388  void write_series();
389  void serie_add_object_(const std::string &serie_name,
390  const std::string &object_name);
391  void write_separ();
392  std::string default_name(std::string s, int count,
393  const char *default_prefix) {
394  if (s.size() == 0) {
395  std::stringstream ss; ss << default_prefix << count; return ss.str();
396  } else return s;
397  }
398  template<class T> void write_val(T v) {
399  if (ascii) os << " " << v;
400  else os.write((char*)&v, sizeof(T));
401  }
402  static const char* endianness() {
403  static int i=0x12345678;
404  char *p = (char*)&i;
405  if (*p == 0x12) return "msb";
406  else if (*p == 0x78) return "lsb";
407  else return "this is very strange..";
408  }
409  bool new_mesh(std::string &name);
410  std::list<dxMesh>::iterator get_mesh(const std::string& name,
411  bool raise_error = true);
412  std::list<dxObject>::iterator get_object(const std::string& name,
413  bool raise_error = true);
414  dxMesh &current_mesh() {
415  if (meshes.size()) return meshes.back();
416  else GMM_ASSERT1(false, "no mesh!");
417  }
418  dxObject &current_data() {
419  if (objects.size()) return objects.back();
420  else GMM_ASSERT1(false, "no data!");
421  }
422 
423  std::string name_of_pts_array(const std::string &meshname)
424  { return meshname + std::string("_pts"); }
425  std::string name_of_conn_array(const std::string &meshname)
426  { return meshname + std::string("_conn"); }
427  std::string name_of_edges_array(const std::string &meshname)
428  { return meshname + std::string("_edges"); }
429  void check_header();
430  const char *dxname_of_convex_structure(bgeot::pconvex_structure cvs);
431  void write_convex_attributes(bgeot::pconvex_structure cvs);
432  void write_mesh_structure_from_slice();
433  void write_mesh_structure_from_mesh_fem();
434  void write_mesh_edges_from_slice(bool with_slice_edge);
435  void write_mesh_edges_from_mesh();
436  template <class VECT>
437  void smooth_field(const VECT& U, base_vector &sU);
438  template<class VECT>
439  void write_dataset_(const VECT& U, std::string name, bool cell_data=false);
440  };
441 
442  template <class VECT>
443  void dx_export::smooth_field(const VECT& U, base_vector &sU) {
444  size_type Q = gmm::vect_size(U) / psl->nb_points();
445  sU.clear(); sU.resize(Q*psl->nb_merged_nodes());
446  for (size_type i=0; i < psl->nb_merged_nodes(); ++i) {
447  for (size_type j=0; j < psl->merged_point_cnt(i); ++j)
448  for (size_type q=0; q < Q; ++q)
449  sU[i*Q+q] += U[psl->merged_point_nodes(i)[j].pos*Q+q];
450  for (size_type q=0; q < Q; ++q)
451  sU[i*Q+q] /= double(psl->merged_point_cnt(i));
452  }
453  }
454 
455  template<class VECT>
456  void dx_export::write_point_data(const getfem::mesh_fem &mf, const VECT& U,
457  std::string name) {
458  size_type Q = (gmm::vect_size(U) / mf.nb_dof())*mf.get_qdim();
459  if (psl) {
460  std::vector<scalar_type> Uslice(Q*psl->nb_points());
461  psl->interpolate(mf, U, Uslice);
462  write_sliced_point_data(Uslice,name);
463  } else {
464  std::vector<scalar_type> V(pmf->nb_dof() * Q);
465  if (&mf != &(*pmf)) {
466  interpolation(mf, *pmf, U, V);
467  } else gmm::copy(U,V);
468  size_type cnt = 0;
469  for (dal::bv_visitor d(pmf_dof_used); !d.finished(); ++d, ++cnt) {
470  if (cnt != d)
471  for (size_type q=0; q < Q; ++q) {
472  V[cnt*Q + q] = V[d*Q + q];
473  }
474  }
475  V.resize(Q*pmf_dof_used.card());
476  write_dataset_(V, name);
477  }
478  }
479 
480  template<class VECT> void
481  dx_export::write_sliced_point_data(const VECT& Uslice, std::string name) {
482  if (!psl_use_merged)
483  write_dataset_(Uslice, name, false);
484  else {
485  base_vector Umerged; smooth_field(Uslice,Umerged);
486  write_dataset_(Umerged, name, false);
487  }
488  }
489 
490  template<class VECT> void
491  dx_export::write_dataset_(const VECT& U, std::string name, bool cell_data) {
492  write_mesh();
493  objects.push_back(dxObject());
494  name = default_name(name, int(objects.size()), "gf_field");
495  objects.back().name = name;
496  objects.back().mesh = current_mesh_name();
497  size_type nb_val = 0;
498  if (cell_data) {
499  nb_val = psl ? psl->linked_mesh().convex_index().card()
500  : pmf->linked_mesh().convex_index().card();
501  } else {
502  nb_val = psl ? (psl_use_merged ? psl->nb_merged_nodes() : psl->nb_points())
503  : pmf_dof_used.card();
504  }
505  size_type Q = gmm::vect_size(U) / nb_val;
506  GMM_ASSERT1(gmm::vect_size(U) == nb_val*Q,
507  "inconsistency in the size of the dataset: "
508  << gmm::vect_size(U) << " != " << nb_val << "*" << Q);
509 
510  os << "\nobject \"" << name << "_data\" class array type float rank ";
511  if (Q == 1) { os << "0"; } /* scalar data */
512  else if (Q == 4) { os << "2 shape 2 2"; } /* or 2x2 tensor data */
513  else if (Q == 9) { os << "2 shape 3 3"; } /* or 2x2 tensor data */
514  else { os << "1 shape " << Q; } /* fallback: vector data */
515  os << " items " << nb_val;
516  if (!ascii) os << " " << endianness() << " binary";
517  os << " data follows" << endl;
518  for (size_type i=0; i < nb_val*Q; ++i) {
519  write_val(float(U[i]));
520  if (((i+1) % (Q > 1 ? Q : 10)) == 0) write_separ();
521  }
522  write_separ();
523 
524  if (!cell_data)
525  os << "\n attribute \"dep\" string \"positions\"\n";
526  else os << "\n attribute \"dep\" string \"connections\"\n";
527  os << "\n";
528 
529  if (current_mesh().flags & dxMesh::WITH_EDGES) {
530  os << "\nobject \"" << name << "_edges\" class field\n"
531  << " component \"positions\" value \""
532  << name_of_pts_array(current_mesh_name()) << "\"\n"
533  << " component \"connections\" value \""
534  << name_of_conn_array(name_of_edges_array(current_mesh_name()))
535  << "\"\n"
536  << " component \"data\" value \"" << name << "_data\"\n";
537  }
538 
539  /* write footer */
540  os << "\nobject \"" << name << "\" class field\n"
541  << " component \"positions\" value \""
542  << name_of_pts_array(current_mesh_name()) << "\"\n"
543  << " component \"connections\" value \""
544  << name_of_conn_array(current_mesh_name()) << "\"\n"
545  << " component \"data\" value \"" << name << "_data\"\n";
546  }
547 
548  /** @brief POS export.
549 
550  export class to Gmsh post-processing file format.
551 
552  ( http://geuz.org/gmsh )
553 
554  A pos_export can store multiple scalar/vector/tensor fields.
555  */
556 
557  class pos_export {
558  protected:
559  std::ostream& os;
560  char header[256];
561 
562  std::vector<std::vector<float> > pos_pts;
563  std::vector<unsigned> pos_cell_type;
564  std::vector<std::vector<unsigned> > pos_cell_dof;
565 
566  std::unique_ptr<mesh_fem> pmf;
567  const stored_mesh_slice *psl;
568 
569  size_type view;
570  dim_type dim;
571  enum { EMPTY, HEADER_WRITTEN, STRUCTURE_WRITTEN, IN_CELL_DATA} state;
572  std::ofstream real_os;
573 
574  public:
575  typedef enum {
576  POS_PT = 0, // point
577  POS_LN = 1, // line
578  POS_TR = 2, // triangles
579  POS_QU = 3, // quadrangles
580  POS_SI = 4, // tetrahedra
581  POS_HE = 5, // hexahedra
582  POS_PR = 6, // prisms
583  POS_PY = 7 // pyramids
584  } pos_cell_types;
585 
586  pos_export(const std::string& fname);
587  pos_export(std::ostream& osname);
588 
589  void set_header(const std::string& s);
590 
591  void exporting(const mesh& m);
592  void exporting(const mesh_fem& mf);
593  void exporting(const stored_mesh_slice& sl);
594 
595  void write(const mesh& m, const std::string& name="");
596  void write(const mesh_fem& mf, const std::string& name="");
597  void write(const stored_mesh_slice& sl, const std::string& name="");
598 
599  template <class VECT>
600  void write(const mesh_fem& mf,const VECT& U, const std::string& name);
601  template <class VECT>
602  void write(const stored_mesh_slice& sl,const VECT& U, const std::string& name);
603 
604  private:
605  void init();
606  void check_header();
607 
608  template <class VECT>
609  void write(const VECT& V, const size_type qdim_v);
610 
611  template <class VECT>
612  void write_cell(const int& t, const std::vector<unsigned>& dof,
613  const VECT& val);
614  };
615 
616  template <class VECT>
617  void pos_export::write(const mesh_fem& mf,const VECT& U,
618  const std::string& name){
619  check_header();
620  exporting(mf);
621 
622  os << "View \"" << name.c_str() <<"\" {\n";
623 
624  size_type nb_points = mf.nb_dof()/mf.get_qdim();
625  size_type qdim_u = gmm::vect_size(U)/nb_points;
626  if (psl){
627  std::vector<scalar_type> Uslice(psl->nb_points()*qdim_u);
628  psl->interpolate(mf, U, Uslice);
629  qdim_u = gmm::vect_size(Uslice)/psl->nb_points();
630  write(Uslice, qdim_u);
631  }else {
632  std::vector<scalar_type> V(pmf->nb_dof()*qdim_u);
633  if (&mf != &(*pmf)) {
634  interpolation(mf, *pmf, U, V);
635  } else gmm::copy(U,V);
636  /*for (dal::bv_visitor d(pmf_dof_used); !d.finished(); ++d, ++cnt) {
637  if (cnt != d)
638  for (size_type q=0; q < Q; ++q) {
639  V[cnt*Q + q] = V[d*Q + q];
640  }
641  }
642  V.resize(Q*pmf_dof_used.card());*/
643  nb_points = pmf->nb_dof()/pmf->get_qdim();
644  qdim_u = gmm::vect_size(V)/nb_points;
645  write(V, qdim_u);
646  }
647 
648  os << "};\n";
649  os << "View[" << view << "].ShowScale = 1;\n";
650  os << "View[" << view << "].ShowElement = 0;\n";
651  os << "View[" << view << "].DrawScalars = 1;\n";
652  os << "View[" << view << "].DrawVectors = 1;\n";
653  os << "View[" << view++ << "].DrawTensors = 1;\n";
654  }
655 
656  template <class VECT>
657  void pos_export::write(const stored_mesh_slice& sl,const VECT& V,
658  const std::string& name){
659  check_header();
660  exporting(sl);
661 
662  os << "View \"" << name.c_str() <<"\" {\n";
663 
664  size_type qdim_v = gmm::vect_size(V)/psl->nb_points();
665  write(V, qdim_v);
666 
667  os << "};\n";
668  os << "View[" << view << "].ShowScale = 1;\n";
669  os << "View[" << view << "].ShowElement = 0;\n";
670  os << "View[" << view << "].DrawScalars = 1;\n";
671  os << "View[" << view << "].DrawVectors = 1;\n";
672  os << "View[" << view++ << "].DrawTensors = 1;\n";
673  }
674 
675  template <class VECT>
676  void pos_export::write(const VECT& V, const size_type qdim_v) {
677  int t;
678  std::vector<unsigned> cell_dof;
679  std::vector<scalar_type> cell_dof_val;
680  for (size_type cell = 0; cell < pos_cell_type.size(); ++cell) {
681  t = pos_cell_type[cell];
682  cell_dof = pos_cell_dof[cell];
683  cell_dof_val.resize(cell_dof.size()*qdim_v, scalar_type(0));
684  for (size_type i=0; i< cell_dof.size(); ++i)
685  for (size_type j=0; j< qdim_v; ++j)
686  cell_dof_val[i*qdim_v+j] = scalar_type(V[cell_dof[i]*qdim_v+j]);
687  write_cell(t,cell_dof,cell_dof_val);
688  }
689  }
690 
691  template <class VECT>
692  void pos_export::write_cell(const int& t, const std::vector<unsigned>& dof,
693  const VECT& val) {
694  size_type qdim_cell = val.size()/dof.size();
695  size_type dim3D = size_type(-1);
696  if (1==qdim_cell){
697  dim3D = size_type(1);
698  os << "S";
699  } else if (2==qdim_cell || 3==qdim_cell){
700  dim3D = size_type(3);
701  os << "V";
702  } else if (4<=qdim_cell && qdim_cell<=9){
703  dim3D = size_type(9);
704  os << "T";
705  }
706  switch (t){
707  case POS_PT: os << "P("; break; // point
708  case POS_LN: os << "L("; break; // line
709  case POS_TR: os << "T("; break; // triangle
710  case POS_QU: os << "Q("; break; // quadrangle
711  case POS_SI: os << "S("; break; // tetrahedra (simplex)
712  case POS_HE: os << "H("; break; // hexahedra
713  case POS_PR: os << "I("; break; // prism
714  case POS_PY: os << "Y("; break; // pyramid
715  }
716  for (size_type i=0; i<dof.size(); ++i){
717  for(size_type j=0; j<dim; ++j){
718  if(0!=i || 0!=j) os << ",";
719  os << pos_pts[dof[i]][j];
720  }
721  for (size_type j=dim; j<3; ++j){
722  os << ",0.00";
723  }
724  }
725 
726  os << "){";
727  for (size_type i=0; i<dof.size(); ++i){
728  for(size_type j=0; j<qdim_cell; ++j){
729  if(0!=i || 0!=j) os << ",";
730  os << val[i*qdim_cell+j];
731  }
732  for (size_type j=qdim_cell; j< dim3D; ++j){
733  os << ",0.00";
734  }
735  }
736  os << "};\n";
737  }
738 } /* end of namespace getfem. */
739 
740 #endif /* GETFEM_EXPORT_H__ */
const mesh & linked_mesh() const
return a pointer to the original mesh
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...
size_type nb_points() const
Return the number of nodes in the slice.
void interpolation(const mesh_fem &mf_source, const mesh_fem &mf_target, const VECTU &U, VECTV &V, int extrapolation=0, double EPS=1E-10, mesh_region rg_source=mesh_region::all_convexes(), mesh_region rg_target=mesh_region::all_convexes())
interpolation/extrapolation of (mf_source, U) on mf_target.
std::string current_mesh_name()
return the name of current mesh (use exporting(...) to change the current mesh)
std::string current_data_name()
return the name of last written data_set
Describe a mesh (collection of convexes (elements) and points).
Definition: getfem_mesh.h:95
void exporting(const mesh &m)
should be called before write_*_data
void write_sliced_point_data(const VECT &Uslice, const std::string &name, size_type qdim=1)
append a new scalar or vector field to .vtk file.
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.
Interpolation of fields from a mesh_fem onto another.
virtual dim_type get_qdim() const
Return the Q dimension.
void interpolate(const getfem::mesh_fem &mf, const V1 &UU, V2 &V) const
Interpolation of a mesh_fem on a slice.
The output of a getfem::mesh_slicer which has been recorded.
A (quite large) class for exportation of data to IBM OpenDX.
GEneric Tool for Finite Element Methods.
Describe a finite element method linked to a mesh.
Define the class getfem::stored_mesh_slice.
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...
VTK export.
Definition: getfem_export.h:67
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.