1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
17 
18 
19 
20 #ifndef LIBMESH_GMV_IO_H
21 #define LIBMESH_GMV_IO_H
22 
23 // Local includes
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/mesh_output.h"
26 #include "libmesh/mesh_input.h"
27 
28 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS
29 namespace libMesh
30 {
31 enum ElemType : int;
32 }
33 #else
34 #include "libmesh/enum_elem_type.h"
35 #endif
36 
37 // C++ includes
38 #include <map>
39 
40 namespace libMesh
41 {
42 
43 // Forward declarations
44 class MeshBase;
45 
46 /**
47  * This class implements writing meshes in the GMV format.
48  * For a full description of the GMV format and to obtain the
49  * GMV software see http://www.generalmeshviewer.com
50  *
51  * \author Benjamin S. Kirk
52  * \date 2004
53  */
54 class GMVIO : public MeshInput<MeshBase>,
55               public MeshOutput<MeshBase>
56 {
57 public:
58 
59   /**
60    * Constructor.  Takes a reference to a constant mesh object.
61    * This constructor will only allow us to write the mesh.
62    */
63   explicit
64   GMVIO (const MeshBase &);
65 
66   /**
67    * Constructor.  Takes a writable reference to a mesh object.
68    * This constructor is required to let us read in a mesh.
69    */
70   explicit
71   GMVIO (MeshBase &);
72 
73   /**
74    * This method implements writing a mesh to a specified file.
75    */
76   virtual void write (const std::string &) override;
77 
78   /**
79    * This method implements reading a mesh from a specified file.
80    */
81   virtual void read (const std::string & mesh_file) override;
82 
83   /**
84    * Bring in base class functionality for name resolution and to
85    * avoid warnings about hidden overloaded virtual functions.
86    */
87   using MeshOutput<MeshBase>::write_nodal_data;
88 
89   /**
90    * This method implements writing a mesh with nodal data to a
91    * specified file where the nodal data and variable names are provided.
92    */
93   virtual void write_nodal_data (const std::string &,
94                                  const std::vector<Number> &,
95                                  const std::vector<std::string> &) override;
96 
97   /**
98    * Flag indicating whether or not to write a binary file.  While
99    * binary files may end up being smaller than equivalent ASCII
100    * files, they are harder to debug if anything goes wrong, since
101    * they are not human-readable.
102    */
binary()103   bool & binary () { return _binary; }
104 
105   /**
106    * Flag indicating whether or not to write the mesh
107    * as discontinuous cell patches
108    */
discontinuous()109   bool & discontinuous() { return _discontinuous; }
110 
111   /**
112    * Flag indicating whether or not to write the partitioning
113    * information for the mesh.
114    */
partitioning()115   bool & partitioning() { return _partitioning; }
116 
117   /**
118    * Flag to write element subdomain_id's as GMV "materials" instead
119    * of element processor_id's.  Allows you to generate exploded views
120    * on user-defined subdomains, potentially creating a pretty picture.
121    */
write_subdomain_id_as_material()122   bool & write_subdomain_id_as_material() { return _write_subdomain_id_as_material; }
123 
124   /**
125    * Flag indicating whether or not to subdivide second order
126    * elements
127    */
subdivide_second_order()128   bool & subdivide_second_order() { return _subdivide_second_order; }
129 
130   /**
131    * Flag indicating whether or not to write p level
132    * information for p refined meshes
133    */
p_levels()134   bool & p_levels() { return _p_levels; }
135 
136   /**
137    * Writes a GMV file with discontinuous data
138    */
139   void write_discontinuous_gmv (const std::string & name,
140                                 const EquationSystems & es,
141                                 const bool write_partitioning,
142                                 const std::set<std::string> * system_names=nullptr) const;
143 
144 
145   /**
146    * This method implements writing a mesh with nodal data to a
147    * specified file where the nodal data and variable names are optionally
148    * provided.  This will write an ASCII file.   This is the new implementation
149    * (without subcells).
150    */
151   void write_ascii_new_impl (const std::string &,
152                              const std::vector<Number> * = nullptr,
153                              const std::vector<std::string> * = nullptr);
154 
155   /**
156    * Takes a vector of cell-centered data to be plotted.
157    * You must ensure that for every active element e,
158    * v[e->id()] is a valid number.  You can add an arbitrary
159    * number of different cell-centered data sets by calling
160    * this function multiple times.
161    *
162    * .) GMV does not like spaces in the cell_centered_data_name
163    * .) No matter what order you add cell-centered data, it will be
164    *    output alphabetically.
165    */
166   void add_cell_centered_data (const std::string &       cell_centered_data_name,
167                                const std::vector<Real> * cell_centered_data_vals);
168 
169   /**
170    * If we read in a nodal solution while reading in a mesh, we can attempt
171    * to copy that nodal solution into an EquationSystems object.
172    */
173   void copy_nodal_solution(EquationSystems & es);
174 
175 private:
176 
177   /**
178    * This method implements writing a mesh with nodal data to a
179    * specified file where the nodal data and variable names are optionally
180    * provided.  This will write an ASCII file.  This is the old implementation
181    * (using subcells) which was the default in libMesh-0.4.3-rc2.
182    */
183   void write_ascii_old_impl (const std::string &,
184                              const std::vector<Number> * = nullptr,
185                              const std::vector<std::string> * = nullptr);
186 
187   /**
188    * This method implements writing a mesh with nodal data to a
189    * specified file where the nodal data and variable names are optionally
190    * provided.
191    */
192   void write_binary (const std::string &,
193                      const std::vector<Number> * = nullptr,
194                      const std::vector<std::string> * = nullptr);
195 
196   /**
197    * Flag to write binary data.
198    */
199   bool _binary;
200 
201   /**
202    * Flag to write the mesh as discontinuous patches.
203    */
204   bool _discontinuous;
205 
206   /**
207    * Flag to write the mesh partitioning.
208    */
209   bool _partitioning;
210 
211   /**
212    * Flag to write element subdomain_id's as GMV "materials" instead
213    * of element processor_id's.
214    */
215   bool _write_subdomain_id_as_material;
216 
217   /**
218    * Flag to subdivide second order elements
219    */
220   bool _subdivide_second_order;
221 
222   /**
223    * Flag to write the mesh p refinement levels.
224    */
225   bool _p_levels;
226 
227   /**
228    * Storage for arbitrary cell-centered data.  Ex: You can use this
229    * to plot the effectivity index for a given cell.  The map is
230    * between the string representing the variable name and a pointer
231    * to a vector containing the data.
232    */
233   std::map<std::string, const std::vector<Real> * > _cell_centered_data;
234 
235   /**
236    * Helper functions for reading nodes/cells from a GMV file
237    */
238   void _read_nodes();
239   unsigned int _next_elem_id;
240   void _read_one_cell();
241   ElemType gmv_elem_to_libmesh_elem(std::string elemname);
242   void _read_materials();
243   void _read_var();
244   std::map<std::string, std::vector<Number>> _nodal_data;
245 
246   /**
247    * Static map from string -> ElementType for
248    * use during reading.
249    */
250   static std::map<std::string, ElemType> _reading_element_map;
251 
252   /**
253    * Static function used to build the _reading_element_map.
254    */
255   static std::map<std::string, ElemType> build_reading_element_map();
256 };
257 
258 } // namespace libMesh
259 
260 
261 #endif // LIBMESH_GMV_IO_H
262