$extrastylesheet
tecplot_io.C
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 
00020 // C++ includes
00021 #include <fstream>
00022 #include <iomanip>
00023 #include <sstream>
00024 
00025 // Local includes
00026 #include "libmesh/libmesh_config.h"
00027 #include "libmesh/libmesh_logging.h"
00028 #include "libmesh/tecplot_io.h"
00029 #include "libmesh/mesh_base.h"
00030 #include "libmesh/elem.h"
00031 #include "libmesh/parallel.h"
00032 
00033 #ifdef LIBMESH_HAVE_TECPLOT_API
00034 extern "C" {
00035 # include <TECIO.h>
00036 }
00037 #endif
00038 
00039 
00040 namespace libMesh
00041 {
00042 
00043 
00044 //--------------------------------------------------------
00045 // Macros for handling Tecplot API data
00046 
00047 #ifdef LIBMESH_HAVE_TECPLOT_API
00048 
00049 namespace
00050 {
00051 class TecplotMacros
00052 {
00053 public:
00054   TecplotMacros(const dof_id_type n_nodes,
00055                 const unsigned int n_vars,
00056                 const dof_id_type n_cells,
00057                 const unsigned int n_vert);
00058   float & nd(const std::size_t i, const std::size_t j);
00059   int   & cd(const std::size_t i, const std::size_t j);
00060   std::vector<float> nodalData;
00061   std::vector<int>   connData;
00062   //float* nodalData;
00063   //int*   connData;
00064 
00065   void set_n_cells (const dof_id_type nc);
00066 
00067   const dof_id_type n_nodes;
00068   const unsigned int n_vars;
00069   dof_id_type n_cells;
00070   const unsigned int n_vert;
00071 };
00072 }
00073 
00074 
00075 
00076 inline
00077 TecplotMacros::TecplotMacros(const dof_id_type nn,
00078                              const unsigned int nvar,
00079                              const dof_id_type nc,
00080                              const unsigned int nvrt) :
00081   n_nodes(nn),
00082   n_vars(nvar),
00083   n_cells(nc),
00084   n_vert(nvrt)
00085 {
00086   nodalData.resize(n_nodes*n_vars);
00087   connData.resize(n_cells*n_vert);
00088 }
00089 
00090 
00091 
00092 inline
00093 float & TecplotMacros::nd(const std::size_t i, const std::size_t j)
00094 {
00095   return nodalData[(i)*(n_nodes) + (j)];
00096 }
00097 
00098 
00099 
00100 inline
00101 int & TecplotMacros::cd(const std::size_t i, const std::size_t j)
00102 {
00103   return connData[(i) + (j)*(n_vert)];
00104 }
00105 
00106 
00107 inline
00108 void TecplotMacros::set_n_cells (const dof_id_type nc)
00109 {
00110   n_cells = nc;
00111   connData.resize(n_cells*n_vert);
00112 }
00113 #endif
00114 //--------------------------------------------------------
00115 
00116 
00117 
00118 // ------------------------------------------------------------
00119 // TecplotIO  members
00120 TecplotIO::TecplotIO (const MeshBase& mesh_in,
00121                       const bool binary_in,
00122                       const double time_in,
00123                       const int strand_offset_in) :
00124   MeshOutput<MeshBase> (mesh_in),
00125   _binary (binary_in),
00126   _time (time_in),
00127   _strand_offset (strand_offset_in),
00128   _zone_title ("zone"),
00129   _ascii_append(false)
00130 {
00131   // Gather a list of subdomain ids in the mesh.
00132   // We must do this now, while we have every
00133   // processor's attention
00134   // (some of the write methods only execute on processor 0).
00135   mesh_in.subdomain_ids (_subdomain_ids);
00136 }
00137 
00138 
00139 
00140 bool & TecplotIO::binary ()
00141 {
00142   return _binary;
00143 }
00144 
00145 
00146 
00147 double & TecplotIO::time ()
00148 {
00149   return _time;
00150 }
00151 
00152 
00153 
00154 int & TecplotIO::strand_offset ()
00155 {
00156   return _strand_offset;
00157 }
00158 
00159 
00160 
00161 std::string & TecplotIO::zone_title ()
00162 {
00163   return _zone_title;
00164 }
00165 
00166 
00167 bool & TecplotIO::ascii_append ()
00168 {
00169   return _ascii_append;
00170 }
00171 
00172 
00173 void TecplotIO::write (const std::string& fname)
00174 {
00175   if (this->mesh().processor_id() == 0)
00176     {
00177       if (this->binary())
00178         this->write_binary (fname);
00179       else
00180         this->write_ascii  (fname);
00181     }
00182 }
00183 
00184 
00185 
00186 void TecplotIO::write_nodal_data (const std::string& fname,
00187                                   const std::vector<Number>& soln,
00188                                   const std::vector<std::string>& names)
00189 {
00190   START_LOG("write_nodal_data()", "TecplotIO");
00191 
00192   if (this->mesh().processor_id() == 0)
00193     {
00194       if (this->binary())
00195         this->write_binary (fname, &soln, &names);
00196       else
00197         this->write_ascii  (fname, &soln, &names);
00198     }
00199 
00200   STOP_LOG("write_nodal_data()", "TecplotIO");
00201 }
00202 
00203 
00204 
00205 unsigned TecplotIO::elem_dimension()
00206 {
00207   // Get a constant reference to the mesh.
00208   const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
00209 
00210   std::vector<unsigned> elem_dims(3);
00211 
00212   // Loop over all the elements and mark the proper dimension entry in
00213   // the elem_dims vector.
00214   MeshBase::const_element_iterator       it  = the_mesh.active_elements_begin();
00215   const MeshBase::const_element_iterator end = the_mesh.active_elements_end();
00216   for ( ; it != end; ++it)
00217     elem_dims[(*it)->dim() - 1] = 1;
00218 
00219   // Detect and disallow (for now) the writing of mixed dimension meshes.
00220   if (std::count(elem_dims.begin(), elem_dims.end(), 1) > 1)
00221     libmesh_error_msg("Error, cannot write Mesh with mixed element dimensions to Tecplot file!");
00222 
00223   if (elem_dims[0])
00224     return 1;
00225   else if (elem_dims[1])
00226     return 2;
00227   else if (elem_dims[2])
00228     return 3;
00229   else
00230     libmesh_error_msg("No 1, 2, or 3D elements detected!");
00231 }
00232 
00233 
00234 
00235 void TecplotIO::write_ascii (const std::string& fname,
00236                              const std::vector<Number>* v,
00237                              const std::vector<std::string>* solution_names)
00238 {
00239   // Should only do this on processor 0!
00240   libmesh_assert_equal_to (this->mesh().processor_id(), 0);
00241 
00242   // Create an output stream, possibly in append mode.
00243   std::ofstream out_stream(fname.c_str(), _ascii_append ? std::ofstream::app : std::ofstream::out);
00244 
00245   // Make sure it opened correctly
00246   if (!out_stream.good())
00247     libmesh_file_error(fname.c_str());
00248 
00249   // Get a constant reference to the mesh.
00250   const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
00251 
00252   // Write header to stream
00253   {
00254     {
00255       // TODO: We used to print out the SVN revision here when we did keyword expansions...
00256       out_stream << "# For a description of the Tecplot format see the Tecplot User's guide.\n"
00257                  << "#\n";
00258     }
00259 
00260     out_stream << "Variables=x,y,z";
00261 
00262     if (solution_names != NULL)
00263       for (unsigned int n=0; n<solution_names->size(); n++)
00264         {
00265 #ifdef LIBMESH_USE_REAL_NUMBERS
00266 
00267           // Write variable names for real variables
00268           out_stream << "," << (*solution_names)[n];
00269 
00270 #else
00271 
00272           // Write variable names for complex variables
00273           out_stream << "," << "r_"   << (*solution_names)[n]
00274                      << "," << "i_"   << (*solution_names)[n]
00275                      << "," << "a_"   << (*solution_names)[n];
00276 
00277 #endif
00278         }
00279 
00280     out_stream << '\n';
00281 
00282     out_stream << "Zone f=fepoint, n=" << the_mesh.n_nodes() << ", e=" << the_mesh.n_active_sub_elem();
00283 
00284     // We cannot choose the element type simply based on the mesh
00285     // dimension... there might be 1D elements living in a 3D mesh.
00286     // So look at the elements which are actually in the Mesh, and
00287     // choose either "lineseg", "quadrilateral", or "brick" depending
00288     // on if the elements are 1, 2, or 3D.
00289 
00290     // Write the element type we've determined to the header.
00291     out_stream << ", et=";
00292 
00293     switch (this->elem_dimension())
00294       {
00295       case 1:
00296         out_stream << "lineseg";
00297         break;
00298       case 2:
00299         out_stream << "quadrilateral";
00300         break;
00301       case 3:
00302         out_stream << "brick";
00303         break;
00304       default:
00305         libmesh_error_msg("Unsupported element dimension: " << this->elem_dimension());
00306       }
00307 
00308     // Output the time in the header
00309     out_stream << ", t=\"T " << _time << "\"";
00310 
00311     // Use default mesh color = black
00312     out_stream << ", c=black\n";
00313 
00314   } // finished writing header
00315 
00316   for (unsigned int i=0; i<the_mesh.n_nodes(); i++)
00317     {
00318       // Print the point without a newline
00319       the_mesh.point(i).write_unformatted(out_stream, false);
00320 
00321       if ((v != NULL) && (solution_names != NULL))
00322         {
00323           const std::size_t n_vars = solution_names->size();
00324 
00325 
00326           for (std::size_t c=0; c<n_vars; c++)
00327             {
00328 #ifdef LIBMESH_USE_REAL_NUMBERS
00329               // Write real data
00330               out_stream << std::setprecision(this->ascii_precision())
00331                          << (*v)[i*n_vars + c] << " ";
00332 
00333 #else
00334               // Write complex data
00335               out_stream << std::setprecision(this->ascii_precision())
00336                          << (*v)[i*n_vars + c].real() << " "
00337                          << (*v)[i*n_vars + c].imag() << " "
00338                          << std::abs((*v)[i*n_vars + c]) << " ";
00339 
00340 #endif
00341             }
00342         }
00343 
00344       // Write a new line after the data for this node
00345       out_stream << '\n';
00346     }
00347 
00348   MeshBase::const_element_iterator       it  = the_mesh.active_elements_begin();
00349   const MeshBase::const_element_iterator end = the_mesh.active_elements_end();
00350 
00351   for ( ; it != end; ++it)
00352     (*it)->write_connectivity(out_stream, TECPLOT);
00353 }
00354 
00355 
00356 
00357 void TecplotIO::write_binary (const std::string& fname,
00358                               const std::vector<Number>* vec,
00359                               const std::vector<std::string>* solution_names)
00360 {
00361   //-----------------------------------------------------------
00362   // Call the ASCII output function if configure did not detect
00363   // the Tecplot binary API
00364 #ifndef LIBMESH_HAVE_TECPLOT_API
00365 
00366   libMesh::err << "WARNING: Tecplot Binary files require the Tecplot API." << std::endl
00367                << "Continuing with ASCII output."
00368                << std::endl;
00369 
00370   if (this->mesh().processor_id() == 0)
00371     this->write_ascii (fname, vec, solution_names);
00372   return;
00373 
00374 
00375 
00376   //------------------------------------------------------------
00377   // New binary formats, time aware and whatnot
00378 #elif defined(LIBMESH_HAVE_TECPLOT_API_112)
00379 
00380   // Get a constant reference to the mesh.
00381   const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
00382 
00383   // Required variables
00384   std::string tecplot_variable_names;
00385   int
00386     ierr      =  0,
00387     file_type =  0, // full
00388     is_double =  0,
00389 #ifdef DEBUG
00390     tec_debug =  1,
00391 #else
00392     tec_debug =  0,
00393 #endif
00394     cell_type   = -1,
00395     nn_per_elem = -1;
00396 
00397   switch (this->elem_dimension())
00398     {
00399     case 1:
00400       cell_type   = 1;  // FELINESEG
00401       nn_per_elem = 2;
00402       break;
00403 
00404     case 2:
00405       cell_type   = 3; // FEQUADRILATERAL
00406       nn_per_elem = 4;
00407       break;
00408 
00409     case 3:
00410       cell_type   = 5; // FEBRICK
00411       nn_per_elem = 8;
00412       break;
00413 
00414     default:
00415       libmesh_error_msg("Unsupported element dimension: " << this->elem_dimension());
00416     }
00417 
00418   // Build a string containing all the variable names to pass to Tecplot
00419   {
00420     tecplot_variable_names += "x, y, z";
00421 
00422     if (solution_names != NULL)
00423       {
00424         for (unsigned int name=0; name<solution_names->size(); name++)
00425           {
00426 #ifdef LIBMESH_USE_REAL_NUMBERS
00427 
00428             tecplot_variable_names += ", ";
00429             tecplot_variable_names += (*solution_names)[name];
00430 
00431 #else
00432 
00433             tecplot_variable_names += ", ";
00434             tecplot_variable_names += "r_";
00435             tecplot_variable_names += (*solution_names)[name];
00436             tecplot_variable_names += ", ";
00437             tecplot_variable_names += "i_";
00438             tecplot_variable_names += (*solution_names)[name];
00439             tecplot_variable_names += ", ";
00440             tecplot_variable_names += "a_";
00441             tecplot_variable_names += (*solution_names)[name];
00442 
00443 #endif
00444           }
00445       }
00446   }
00447 
00448   // Instantiate a TecplotMacros interface.  In 2D the most nodes per
00449   // face should be 4, in 3D it's 8.
00450 
00451 
00452   TecplotMacros tm(the_mesh.n_nodes(),
00453 #ifdef LIBMESH_USE_REAL_NUMBERS
00454                    (3 + ((solution_names == NULL) ? 0 :
00455                          cast_int<unsigned int>(solution_names->size()))),
00456 #else
00457                    (3 + 3*((solution_names == NULL) ? 0 :
00458                            cast_int<unsigned int>(solution_names->size()))),
00459 #endif
00460                    the_mesh.n_active_sub_elem(),
00461                    nn_per_elem
00462                    );
00463 
00464 
00465   // Copy the nodes and data to the TecplotMacros class. Note that we store
00466   // everything as a float here since the eye doesn't require a double to
00467   // understand what is going on
00468   for (unsigned int v=0; v<the_mesh.n_nodes(); v++)
00469     {
00470       tm.nd(0,v) = static_cast<float>(the_mesh.point(v)(0));
00471       tm.nd(1,v) = static_cast<float>(the_mesh.point(v)(1));
00472       tm.nd(2,v) = static_cast<float>(the_mesh.point(v)(2));
00473 
00474       if ((vec != NULL) &&
00475           (solution_names != NULL))
00476         {
00477           const std::size_t n_vars = solution_names->size();
00478 
00479           for (std::size_t c=0; c<n_vars; c++)
00480             {
00481 #ifdef LIBMESH_USE_REAL_NUMBERS
00482 
00483               tm.nd((3+c),v)     = static_cast<float>((*vec)[v*n_vars + c]);
00484 #else
00485               tm.nd((3+3*c),v)   = static_cast<float>((*vec)[v*n_vars + c].real());
00486               tm.nd((3+3*c+1),v) = static_cast<float>((*vec)[v*n_vars + c].imag());
00487               tm.nd((3+3*c+2),v) = static_cast<float>(std::abs((*vec)[v*n_vars + c]));
00488 #endif
00489             }
00490         }
00491     }
00492 
00493 
00494   // Initialize the file
00495   ierr = TECINI112 (NULL,
00496                     const_cast<char*>(tecplot_variable_names.c_str()),
00497                     const_cast<char*>(fname.c_str()),
00498                     const_cast<char*>("."),
00499                     &file_type,
00500                     &tec_debug,
00501                     &is_double);
00502 
00503   if (ierr)
00504     libmesh_file_error(fname);
00505 
00506   // A zone for each subdomain
00507   bool firstzone=true;
00508   for (std::set<subdomain_id_type>::const_iterator sbd_it=_subdomain_ids.begin();
00509        sbd_it!=_subdomain_ids.end(); ++sbd_it)
00510     {
00511       // Copy the connectivity for this subdomain
00512       {
00513         MeshBase::const_element_iterator       it  = the_mesh.active_subdomain_elements_begin (*sbd_it);
00514         const MeshBase::const_element_iterator end = the_mesh.active_subdomain_elements_end   (*sbd_it);
00515 
00516         unsigned int n_subcells_in_subdomain=0;
00517 
00518         for (; it != end; ++it)
00519           n_subcells_in_subdomain += (*it)->n_sub_elem();
00520 
00521         // update the connectivty array to include only the elements in this subdomain
00522         tm.set_n_cells (n_subcells_in_subdomain);
00523 
00524         unsigned int te = 0;
00525 
00526         for (it  = the_mesh.active_subdomain_elements_begin (*sbd_it);
00527              it != end; ++it)
00528           {
00529             std::vector<dof_id_type> conn;
00530             for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
00531               {
00532                 (*it)->connectivity(se, TECPLOT, conn);
00533 
00534                 for (unsigned int node=0; node<conn.size(); node++)
00535                   tm.cd(node,te) = conn[node];
00536 
00537                 te++;
00538               }
00539           }
00540       }
00541 
00542 
00543       // Ready to call the Tecplot API for this subdomain
00544       {
00545         int
00546           num_nodes   = static_cast<int>(the_mesh.n_nodes()),
00547           num_cells   = static_cast<int>(tm.n_cells),
00548           num_faces   = 0,
00549           i_cell_max  = 0,
00550           j_cell_max  = 0,
00551           k_cell_max  = 0,
00552           strand_id   = std::max(*sbd_it,static_cast<subdomain_id_type>(1)) + this->strand_offset(),
00553           parent_zone = 0,
00554           is_block    = 1,
00555           num_face_connect   = 0,
00556           face_neighbor_mode = 0,
00557           tot_num_face_nodes = 0,
00558           num_connect_boundary_faces = 0,
00559           tot_num_boundary_connect   = 0,
00560           share_connect_from_zone=0;
00561 
00562         std::vector<int>
00563           passive_var_list    (tm.n_vars, 0),
00564           share_var_from_zone (tm.n_vars, 1); // We only write data for the first zone, all other
00565         // zones will share from this one.
00566 
00567         // get the subdomain name from libMesh, if there is one.
00568         std::string subdomain_name = the_mesh.subdomain_name(*sbd_it);
00569         std::ostringstream zone_name;
00570         zone_name << this->zone_title();
00571 
00572         // We will title this
00573         // "{zone_title()}_{subdomain_name}", or
00574         // "{zone_title()}_{subdomain_id}", or
00575         // "{zone_title()}"
00576         if (subdomain_name.size())
00577           {
00578             zone_name << "_";
00579             zone_name << subdomain_name;
00580           }
00581         else if (_subdomain_ids.size() > 1)
00582           {
00583             zone_name << "_";
00584             zone_name << *sbd_it;
00585           }
00586 
00587         ierr = TECZNE112 (const_cast<char*>(zone_name.str().c_str()),
00588                           &cell_type,
00589                           &num_nodes,
00590                           &num_cells,
00591                           &num_faces,
00592                           &i_cell_max,
00593                           &j_cell_max,
00594                           &k_cell_max,
00595                           &_time,
00596                           &strand_id,
00597                           &parent_zone,
00598                           &is_block,
00599                           &num_face_connect,
00600                           &face_neighbor_mode,
00601                           &tot_num_face_nodes,
00602                           &num_connect_boundary_faces,
00603                           &tot_num_boundary_connect,
00604                           &passive_var_list[0],
00605                           NULL, // = all are node centered
00606                           (firstzone) ? NULL : &share_var_from_zone[0],
00607                           &share_connect_from_zone);
00608 
00609         if (ierr)
00610           libmesh_file_error(fname);
00611 
00612         // Write *all* the data for the first zone, then share it with the others
00613         if (firstzone)
00614           {
00615             int total = cast_int<int>
00616 #ifdef LIBMESH_USE_REAL_NUMBERS
00617               ((3 + ((solution_names == NULL) ? 0 : solution_names->size()))*num_nodes);
00618 #else
00619             ((3 + 3*((solution_names == NULL) ? 0 : solution_names->size()))*num_nodes);
00620 #endif
00621 
00622 
00623             ierr = TECDAT112 (&total,
00624                               &tm.nodalData[0],
00625                               &is_double);
00626 
00627             if (ierr)
00628               libmesh_file_error(fname);
00629           }
00630 
00631         // Write the connectivity
00632         ierr = TECNOD112 (&tm.connData[0]);
00633 
00634         if (ierr)
00635           libmesh_file_error(fname);
00636       }
00637 
00638       firstzone = false;
00639     }
00640 
00641   // Done, close the file.
00642   ierr = TECEND112 ();
00643 
00644   if (ierr)
00645     libmesh_file_error(fname);
00646 
00647 
00648 
00649 
00650   //------------------------------------------------------------
00651   // Legacy binary format
00652 #else
00653 
00654   // Get a constant reference to the mesh.
00655   const MeshBase& the_mesh = MeshOutput<MeshBase>::mesh();
00656 
00657   // Tecplot binary output only good for dim=2,3
00658   if (the_mesh.mesh_dimension() == 1)
00659     {
00660       this->write_ascii (fname, vec, solution_names);
00661 
00662       return;
00663     }
00664 
00665   // Required variables
00666   std::string tecplot_variable_names;
00667   int is_double =  0,
00668     tec_debug =  0,
00669     cell_type = ((the_mesh.mesh_dimension()==2) ? (1) : (3));
00670 
00671   // Build a string containing all the variable names to pass to Tecplot
00672   {
00673     tecplot_variable_names += "x, y, z";
00674 
00675     if (solution_names != NULL)
00676       {
00677         for (unsigned int name=0; name<solution_names->size(); name++)
00678           {
00679 #ifdef LIBMESH_USE_REAL_NUMBERS
00680 
00681             tecplot_variable_names += ", ";
00682             tecplot_variable_names += (*solution_names)[name];
00683 
00684 #else
00685 
00686             tecplot_variable_names += ", ";
00687             tecplot_variable_names += "r_";
00688             tecplot_variable_names += (*solution_names)[name];
00689             tecplot_variable_names += ", ";
00690             tecplot_variable_names += "i_";
00691             tecplot_variable_names += (*solution_names)[name];
00692             tecplot_variable_names += ", ";
00693             tecplot_variable_names += "a_";
00694             tecplot_variable_names += (*solution_names)[name];
00695 
00696 #endif
00697           }
00698       }
00699   }
00700 
00701   // Instantiate a TecplotMacros interface.  In 2D the most nodes per
00702   // face should be 4, in 3D it's 8.
00703 
00704 
00705   TecplotMacros tm(cast_int<unsigned int>(the_mesh.n_nodes()),
00706                    cast_int<unsigned int>
00707 #ifdef LIBMESH_USE_REAL_NUMBERS
00708                    (3 + ((solution_names == NULL) ? 0 : solution_names->size())),
00709 #else
00710                    (3 + 3*((solution_names == NULL) ? 0 : solution_names->size())),
00711 #endif
00712                    cast_int<unsigned int>
00713                    (the_mesh.n_active_sub_elem()),
00714                    ((the_mesh.mesh_dimension() == 2) ? 4 : 8)
00715                    );
00716 
00717 
00718   // Copy the nodes and data to the TecplotMacros class. Note that we store
00719   // everything as a float here since the eye doesn't require a double to
00720   // understand what is going on
00721   for (unsigned int v=0; v<the_mesh.n_nodes(); v++)
00722     {
00723       tm.nd(0,v) = static_cast<float>(the_mesh.point(v)(0));
00724       tm.nd(1,v) = static_cast<float>(the_mesh.point(v)(1));
00725       tm.nd(2,v) = static_cast<float>(the_mesh.point(v)(2));
00726 
00727       if ((vec != NULL) &&
00728           (solution_names != NULL))
00729         {
00730           const unsigned int n_vars = solution_names->size();
00731 
00732           for (unsigned int c=0; c<n_vars; c++)
00733             {
00734 #ifdef LIBMESH_USE_REAL_NUMBERS
00735 
00736               tm.nd((3+c),v)     = static_cast<float>((*vec)[v*n_vars + c]);
00737 #else
00738               tm.nd((3+3*c),v)   = static_cast<float>((*vec)[v*n_vars + c].real());
00739               tm.nd((3+3*c+1),v) = static_cast<float>((*vec)[v*n_vars + c].imag());
00740               tm.nd((3+3*c+2),v) = static_cast<float>(std::abs((*vec)[v*n_vars + c]));
00741 #endif
00742             }
00743         }
00744     }
00745 
00746 
00747   // Copy the connectivity
00748   {
00749     unsigned int te = 0;
00750 
00751     MeshBase::const_element_iterator       it  = the_mesh.active_elements_begin();
00752     const MeshBase::const_element_iterator end = the_mesh.active_elements_end();
00753 
00754     for ( ; it != end; ++it)
00755       {
00756         std::vector<dof_id_type> conn;
00757         for (unsigned int se=0; se<(*it)->n_sub_elem(); se++)
00758           {
00759             (*it)->connectivity(se, TECPLOT, conn);
00760 
00761             for (unsigned int node=0; node<conn.size(); node++)
00762               tm.cd(node,te) = conn[node];
00763 
00764             te++;
00765           }
00766       }
00767   }
00768 
00769 
00770   // Ready to call the Tecplot API
00771   {
00772     int ierr = 0,
00773       num_nodes = static_cast<int>(the_mesh.n_nodes()),
00774       num_cells = static_cast<int>(the_mesh.n_active_sub_elem());
00775 
00776 
00777     ierr = TECINI (NULL,
00778                    (char*) tecplot_variable_names.c_str(),
00779                    (char*) fname.c_str(),
00780                    (char*) ".",
00781                    &tec_debug,
00782                    &is_double);
00783 
00784     if (ierr)
00785       libmesh_file_error(fname);
00786 
00787 
00788     ierr = TECZNE (NULL,
00789                    &num_nodes,
00790                    &num_cells,
00791                    &cell_type,
00792                    (char*) "FEBLOCK",
00793                    NULL);
00794 
00795     if (ierr)
00796       libmesh_file_error(fname);
00797 
00798 
00799     int total =
00800 #ifdef LIBMESH_USE_REAL_NUMBERS
00801       ((3 + ((solution_names == NULL) ? 0 : solution_names->size()))*num_nodes);
00802 #else
00803     ((3 + 3*((solution_names == NULL) ? 0 : solution_names->size()))*num_nodes);
00804 #endif
00805 
00806 
00807     ierr = TECDAT (&total,
00808                    &tm.nodalData[0],
00809                    &is_double);
00810 
00811     if (ierr)
00812       libmesh_file_error(fname);
00813 
00814     ierr = TECNOD (&tm.connData[0]);
00815 
00816     if (ierr)
00817       libmesh_file_error(fname);
00818 
00819     ierr = TECEND ();
00820 
00821     if (ierr)
00822       libmesh_file_error(fname);
00823   }
00824 
00825 #endif
00826 }
00827 
00828 } // namespace libMesh