38 #ifndef GETFEM_EXPORT_H__ 
   39 #define GETFEM_EXPORT_H__ 
   53   inline std::string remove_spaces(
const std::string &s) {
 
   55     for (
unsigned i=0; i < s.size(); ++i)
 
   56       if (s2[i] <= 
' ') s2[i] = 
'_';
 
   75     std::unique_ptr<mesh_fem> pmf;
 
   76     dal::bit_vector pmf_dof_used;
 
   77     std::vector<unsigned> pmf_mapping_type;
 
   78     std::ofstream real_os;
 
   81     std::vector<unsigned char> vals;
 
   82     enum { EMPTY, HEADER_WRITTEN, STRUCTURE_WRITTEN, IN_CELL_DATA,
 
   83            IN_POINT_DATA } state;
 
   85     template<
class T> 
void write_val(T v);
 
   86     template<
class V> 
void write_vec(V p, 
size_type qdim);
 
   87     template<
class IT> 
void write_3x3tensor(IT p);
 
   93     vtk_export(
const std::string& fname, 
bool ascii_ = 
false, 
bool vtk_= 
true);
 
   94     vtk_export(std::ostream &os_, 
bool ascii_ = 
false, 
bool vtk_ = 
true);
 
  117                                                const std::string& name);
 
  124                                                       const std::string& name,
 
  131                                               const std::string& name,
 
  142     void write_normals();
 
  144     const mesh_fem& get_exported_mesh_fem() 
const;
 
  148     void write_mesh_structure_from_slice();
 
  149     void write_mesh_structure_from_mesh_fem();
 
  150     void switch_to_cell_data();
 
  151     void switch_to_point_data();
 
  152     template<
class VECT> 
void write_dataset_(
const VECT& U,
 
  153                                              const std::string& name,
 
  155                                              bool cell_data=
false);
 
  158   template<
class T> 
void vtk_export::write_val(T v) {
 
  159     if (ascii) os << 
" " << v;
 
  164           for (
size_type i=0; i < 
sizeof(v)/2; ++i)
 
  165             std::swap(p[i], p[
sizeof(v)-i-1]);
 
  166         os.write(p, 
sizeof(T));
 
  168         union { T value; 
unsigned char bytes[
sizeof(T)]; } UNION;
 
  171           vals.push_back(UNION.bytes[i]);
 
  176   template<
class IT> 
void vtk_export::write_vec(IT p, 
size_type qdim) {
 
  181     for (
size_type i=qdim; i < 3; ++i) v[i] = 0.0f;
 
  182     write_val(v[0]);write_val(v[1]);write_val(v[2]);
 
  185   template<
class IT> 
void vtk_export::write_3x3tensor(IT p) {
 
  187     memset(v, 0, 
sizeof v);
 
  190         v[i][j] = 
float(p[i + j*dim_]);
 
  196       if (ascii) os << 
"\n";
 
  202                                     const std::string& name) {
 
  206       std::vector<scalar_type> Uslice(Q*psl->
nb_points());
 
  208       write_dataset_(Uslice, name, qdim);
 
  210       std::vector<scalar_type> V(pmf->nb_dof() * Q);
 
  211       if (&mf != &(*pmf)) {
 
  213       } 
else gmm::copy(U,V);
 
  215       for (dal::bv_visitor d(pmf_dof_used); !d.finished(); ++d, ++cnt) {
 
  218             V[cnt*Q + q] = V[d*Q + q];
 
  221       V.resize(Q*pmf_dof_used.card());
 
  222       write_dataset_(V, name, qdim);
 
  229     write_dataset_(U, name, qdim, 
true);
 
  234                                            const std::string& name,
 
  236     write_dataset_(U, name, qdim, 
false);
 
  240   void vtk_export::write_dataset_(
const VECT& U, 
const std::string& name,
 
  245       switch_to_cell_data();
 
  247                    : pmf->linked_mesh().convex_index().card();
 
  249       switch_to_point_data();
 
  250       nb_val = psl ? psl->
nb_points() : pmf_dof_used.card();
 
  253     if (Q == 1) Q = gmm::vect_size(U) / nb_val;
 
  254     GMM_ASSERT1(gmm::vect_size(U) == nb_val*Q,
 
  255                 "inconsistency in the size of the dataset: " 
  256                 << gmm::vect_size(U) << 
" != " << nb_val << 
"*" << Q);
 
  257     if (vtk) write_separ();
 
  258     if (!vtk && !ascii) write_val(
float(gmm::vect_size(U)));
 
  261         os << 
"SCALARS " << remove_spaces(name) << 
" float 1\n" 
  262            << 
"LOOKUP_TABLE default\n";
 
  264         os << 
"<DataArray type=\"Float32\" Name=\"" << remove_spaces(name) << 
"\" " 
  265            << (ascii ? 
"format=\"ascii\">\n" : 
"format=\"binary\">\n");
 
  267         write_val(
float(U[i]));
 
  270         os << 
"VECTORS " << remove_spaces(name) << 
" float\n";
 
  272         os << 
"<DataArray type=\"Float32\" Name=\"" << remove_spaces(name) << 
"\" " 
  273            << 
"NumberOfComponents=\"3\" " 
  274            << (ascii ? 
"format=\"ascii\">\n" : 
"format=\"binary\">\n");
 
  276         write_vec(U.begin() + i*Q, Q);
 
  277     } 
else if (Q == gmm::sqr(dim_)) {
 
  282         os << 
"TENSORS " << remove_spaces(name) << 
" float\n";
 
  284         os << 
"<DataArray type=\"Float32\" Name=\"" << remove_spaces(name)
 
  285            << 
"\" NumberOfComponents=\"9\" " 
  286            << (ascii ? 
"format=\"ascii\">\n" : 
"format=\"binary\">\n");
 
  288         write_3x3tensor(U.begin() + i*Q);
 
  290       GMM_ASSERT1(
false, std::string(vtk ? 
"vtk" : 
"vtu")
 
  291                          + 
" does not accept vectors of dimension > 3");
 
  293     if (vtk) write_separ();
 
  294     if (!vtk) os << 
"\n" << 
"</DataArray>\n";
 
  298   class vtu_export : 
public vtk_export {
 
  300     vtu_export(
const std::string& fname, 
bool ascii_ = 
false) : vtk_export(fname, ascii_, false) {}
 
  301     vtu_export(std::ostream &os_, 
bool ascii_ = 
false) : vtk_export(os_, ascii_, false) {}
 
  325     std::unique_ptr<mesh_fem> pmf;
 
  326     dal::bit_vector pmf_dof_used;
 
  327     std::vector<unsigned> pmf_cell_type;
 
  328     std::fstream real_os;
 
  329     dim_type dim_, connections_dim;
 
  332       std::list<std::string> members;
 
  340       typedef enum { NONE=0, WITH_EDGES=1, STRUCTURE_WRITTEN=2 } flags_t;
 
  342       dxMesh() : flags(NONE) {}
 
  344     std::list<dxObject> objects;
 
  345     std::list<dxSeries> series;
 
  346     std::list<dxMesh> meshes;
 
  349     dx_export(
const std::string& fname, 
bool ascii_ = 
false,
 
  350               bool append_ = 
false);
 
  351     dx_export(std::ostream &os_, 
bool ascii_ = 
false);
 
  354     void exporting(
const mesh& m, std::string name = std::string());
 
  355     void exporting(
const mesh_fem& mf, std::string name = std::string());
 
  357                    std::string name = std::string());
 
  376                           const std::string& object_name);
 
  384     template<
class VECT> 
void 
  386                      const VECT& U0, std::string name = std::string());
 
  387     template<
class VECT> 
void 
  388     write_sliced_point_data(
const VECT& Uslice,
 
  389                             std::string name = std::string());
 
  394     void write_normals();
 
  396     const mesh_fem& get_exported_mesh_fem() 
const;
 
  400     void reread_metadata();
 
  401     void update_metadata(std::ios::pos_type);
 
  403     void serie_add_object_(
const std::string &serie_name,
 
  404                            const std::string &object_name);
 
  406     std::string default_name(std::string s, 
int count,
 
  407                              const char *default_prefix) {
 
  409         std::stringstream ss; ss << default_prefix << count; 
return ss.str();
 
  412     template<
class T> 
void write_val(T v) {
 
  413       if (ascii) os << 
" " << v;
 
  414       else os.write((
char*)&v, 
sizeof(T));
 
  416     static const char* endianness() {
 
  417       static int i=0x12345678;
 
  419       if (*p == 0x12) 
return "msb";
 
  420       else if (*p == 0x78) 
return "lsb";
 
  421       else return "this is very strange..";
 
  423     bool new_mesh(std::string &name);
 
  424     std::list<dxMesh>::iterator get_mesh(
const std::string& name,
 
  425                                          bool raise_error = 
true);
 
  426     std::list<dxObject>::iterator get_object(
const std::string& name,
 
  427                                              bool raise_error = 
true);
 
  428     dxMesh ¤t_mesh() {
 
  429       if (meshes.size()) 
return meshes.back();
 
  430       else GMM_ASSERT1(
false, 
"no mesh!");
 
  432     dxObject ¤t_data() {
 
  433       if (objects.size()) 
return objects.back();
 
  434       else GMM_ASSERT1(
false, 
"no data!");
 
  437     std::string name_of_pts_array(
const std::string &meshname)
 
  438     { 
return meshname + std::string(
"_pts"); }
 
  439     std::string name_of_conn_array(
const std::string &meshname)
 
  440     { 
return meshname + std::string(
"_conn"); }
 
  441     std::string name_of_edges_array(
const std::string &meshname)
 
  442     { 
return meshname + std::string(
"_edges"); }
 
  446     void write_mesh_structure_from_slice();
 
  447     void write_mesh_structure_from_mesh_fem();
 
  448     void write_mesh_edges_from_slice(
bool with_slice_edge);
 
  449     void write_mesh_edges_from_mesh();
 
  450     template <
class VECT>
 
  451     void smooth_field(
const VECT& U, base_vector &sU);
 
  453     void write_dataset_(
const VECT& U, std::string name, 
bool cell_data=
false);
 
  456   template <
class VECT>
 
  457   void dx_export::smooth_field(
const VECT& U, base_vector &sU) {
 
  463           sU[i*Q+q] += U[psl->merged_point_nodes(i)[j].pos*Q+q];
 
  474       std::vector<scalar_type> Uslice(Q*psl->
nb_points());
 
  476       write_sliced_point_data(Uslice,name);
 
  478       std::vector<scalar_type> V(pmf->nb_dof() * Q);
 
  479       if (&mf != &(*pmf)) {
 
  483       for (dal::bv_visitor d(pmf_dof_used); !d.finished(); ++d, ++cnt) {
 
  486             V[cnt*Q + q] = V[d*Q + q];
 
  489       V.resize(Q*pmf_dof_used.card());
 
  490       write_dataset_(V, name);
 
  494   template<
class VECT> 
void 
  495   dx_export::write_sliced_point_data(
const VECT& Uslice, std::string name) {
 
  497       write_dataset_(Uslice, name, 
false);
 
  499       base_vector Umerged; smooth_field(Uslice,Umerged);
 
  500       write_dataset_(Umerged, name, 
false);
 
  504   template<
class VECT> 
void 
  505   dx_export::write_dataset_(
const VECT& U, std::string name, 
bool cell_data) {
 
  507     objects.push_back(dxObject());
 
  508     name = default_name(name, 
int(objects.size()), 
"gf_field");
 
  509     objects.back().name = name;
 
  514                    : pmf->linked_mesh().convex_index().card();
 
  517                    : pmf_dof_used.card();
 
  519     size_type Q = gmm::vect_size(U) / nb_val;
 
  520     GMM_ASSERT1(gmm::vect_size(U) == nb_val*Q,
 
  521                 "inconsistency in the size of the dataset: " 
  522                 << gmm::vect_size(U) << 
" != " << nb_val << 
"*" << Q);
 
  524     os << 
"\nobject \"" << name << 
"_data\" class array type float rank ";
 
  525     if (Q == 1) { os << 
"0"; } 
 
  526     else if (Q == 4) { os << 
"2 shape 2 2"; } 
 
  527     else if (Q == 9) { os << 
"2 shape 3 3"; } 
 
  528     else { os << 
"1 shape " << Q; } 
 
  529     os << 
" items " << nb_val;
 
  530     if (!ascii) os << 
" " << endianness() << 
" binary";
 
  531     os << 
" data follows" << endl;
 
  533       write_val(
float(U[i]));
 
  534       if (((i+1) % (Q > 1 ? Q : 10)) == 0) write_separ();
 
  539       os << 
"\n  attribute \"dep\" string \"positions\"\n";
 
  540     else os << 
"\n  attribute \"dep\" string \"connections\"\n";
 
  543     if (current_mesh().flags & dxMesh::WITH_EDGES) {
 
  544       os << 
"\nobject \"" << name << 
"_edges\" class field\n" 
  545          << 
"  component \"positions\" value \"" 
  547          << 
"  component \"connections\" value \"" 
  550          << 
"  component \"data\" value \"" << name << 
"_data\"\n";
 
  554     os << 
"\nobject \"" << name << 
"\" class field\n" 
  555        << 
"  component \"positions\" value \"" 
  557        << 
"  component \"connections\" value \"" 
  559        << 
"  component \"data\" value \"" << name << 
"_data\"\n";
 
  576     std::vector<std::vector<float> > pos_pts;
 
  577     std::vector<unsigned> pos_cell_type;
 
  578     std::vector<std::vector<unsigned> > pos_cell_dof;
 
  580     std::unique_ptr<mesh_fem> pmf;
 
  585     enum { EMPTY, HEADER_WRITTEN, STRUCTURE_WRITTEN, IN_CELL_DATA } state;
 
  586     std::ofstream real_os;
 
  603     void set_header(
const std::string& s);
 
  605     void exporting(
const mesh& m);
 
  609     void write(
const mesh& m, 
const std::string& name=
"");
 
  610     void write(
const mesh_fem& mf, 
const std::string& name=
"");
 
  613     template <
class VECT>
 
  614     void write(
const mesh_fem& mf,
const VECT& U, 
const std::string& name);
 
  615     template <
class VECT>
 
  622     template <
class VECT>
 
  623     void write(
const VECT& V, 
const size_type qdim_v);
 
  625     template <
class VECT>
 
  626     void write_cell(
const int& t, 
const std::vector<unsigned>& dof,
 
  630   template <
class VECT>
 
  631   void pos_export::write(
const mesh_fem& mf,
const VECT& U,
 
  632                          const std::string& name){
 
  636     os << 
"View \"" << name.c_str() <<
"\" {\n";
 
  639     size_type qdim_u = gmm::vect_size(U)/nb_points;
 
  641       std::vector<scalar_type> Uslice(psl->
nb_points()*qdim_u);
 
  643       qdim_u = gmm::vect_size(Uslice)/psl->
nb_points();
 
  644       write(Uslice, qdim_u);
 
  646       std::vector<scalar_type> V(pmf->nb_dof()*qdim_u);
 
  647       if (&mf != &(*pmf)) {
 
  657       nb_points = pmf->nb_dof()/pmf->get_qdim();
 
  658       qdim_u = gmm::vect_size(V)/nb_points;
 
  663     os << 
"View[" << view << 
"].ShowScale = 1;\n";
 
  664     os << 
"View[" << view << 
"].ShowElement = 0;\n";
 
  665     os << 
"View[" << view << 
"].DrawScalars = 1;\n";
 
  666     os << 
"View[" << view << 
"].DrawVectors = 1;\n";
 
  667     os << 
"View[" << view++ << 
"].DrawTensors = 1;\n";
 
  670   template <
class VECT>
 
  671   void pos_export::write(
const stored_mesh_slice& sl,
const VECT& V,
 
  672                          const std::string& name){
 
  676     os << 
"View \"" << name.c_str() <<
"\" {\n";
 
  682     os << 
"View[" << view << 
"].ShowScale = 1;\n";
 
  683     os << 
"View[" << view << 
"].ShowElement = 0;\n";
 
  684     os << 
"View[" << view << 
"].DrawScalars = 1;\n";
 
  685     os << 
"View[" << view << 
"].DrawVectors = 1;\n";
 
  686     os << 
"View[" << view++ << 
"].DrawTensors = 1;\n";
 
  689   template <
class VECT>
 
  690   void pos_export::write(
const VECT& V, 
const size_type qdim_v) {
 
  692     std::vector<unsigned> cell_dof;
 
  693     std::vector<scalar_type> cell_dof_val;
 
  694     for (
size_type cell = 0; cell < pos_cell_type.size(); ++cell) {
 
  695       t = pos_cell_type[cell];
 
  696       cell_dof = pos_cell_dof[cell];
 
  697       cell_dof_val.resize(cell_dof.size()*qdim_v, scalar_type(0));
 
  698       for (
size_type i=0; i< cell_dof.size(); ++i)
 
  700           cell_dof_val[i*qdim_v+j] = scalar_type(V[cell_dof[i]*qdim_v+j]);
 
  701       write_cell(t,cell_dof,cell_dof_val);
 
  705   template <
class VECT>
 
  706   void pos_export::write_cell(
const int& t, 
const std::vector<unsigned>& dof,
 
  708     size_type qdim_cell = val.size()/dof.size();
 
  713     } 
else if (2==qdim_cell || 3==qdim_cell){
 
  716     } 
else if (4<=qdim_cell && qdim_cell<=9){
 
  721       case POS_PT: os << 
"P("; 
break; 
 
  722       case POS_LN: os << 
"L("; 
break; 
 
  723       case POS_TR: os << 
"T("; 
break; 
 
  724       case POS_QU: os << 
"Q("; 
break; 
 
  725       case POS_SI: os << 
"S("; 
break; 
 
  726       case POS_HE: os << 
"H("; 
break; 
 
  727       case POS_PR: os << 
"I("; 
break; 
 
  728       case POS_PY: os << 
"Y("; 
break; 
 
  732         if(0!=i || 0!=j) os << 
",";
 
  733         os << pos_pts[dof[i]][j];
 
  743         if(0!=i || 0!=j) os << 
",";
 
  744         os << val[i*qdim_cell+j];
 
  746       for (
size_type j=qdim_cell; j< dim3D; ++j){
 
const dal::bit_vector & convex_index() const
Return the list of valid convex IDs.
 
A (quite large) class for exportation of data to IBM OpenDX.
 
std::string current_mesh_name()
return the name of current mesh (use exporting(...) to change the current mesh)
 
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 'series', i.e.
 
std::string current_data_name()
return the name of last written data_set
 
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 this...
 
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.
 
Describe a finite element method linked to a mesh.
 
virtual dim_type get_qdim() const
Return the Q dimension.
 
virtual size_type nb_dof() const
Return the total number of degrees of freedom.
 
Describe a mesh (collection of convexes (elements) and points).
 
The output of a getfem::mesh_slicer which has been recorded.
 
size_type merged_point_cnt(size_type i_merged) const
Return the number of nodes that were merged to form the merged one.
 
const mesh & linked_mesh() const
return a pointer to the original mesh
 
size_type nb_points() const
Return the number of nodes in the slice.
 
void interpolate(const getfem::mesh_fem &mf, const V1 &UU, V2 &V) const
Interpolation of a mesh_fem on a slice.
 
size_type nb_merged_nodes() const
Return the number of merged nodes in 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...
 
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 this...
 
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/.vtu file.
 
void write_cell_data(const VECT &U, const std::string &name, size_type qdim=1)
export data which is constant over each element.
 
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.
 
Interpolation of fields from a mesh_fem onto another.
 
Define the class getfem::stored_mesh_slice.
 
void copy(const L1 &l1, L2 &l2)
*/
 
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
 
GEneric Tool for Finite Element Methods.
 
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.