Line data Source code
1 : // The libMesh Finite Element Library.
2 : // Copyright (C) 2002-2025 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 : // Local includes
21 : #include "libmesh/libmesh_config.h"
22 : #include "libmesh/libmesh_logging.h"
23 : #include "libmesh/tecplot_io.h"
24 : #include "libmesh/mesh_base.h"
25 : #include "libmesh/elem.h"
26 : #include "libmesh/enum_io_package.h"
27 : #include "libmesh/int_range.h"
28 :
29 : #ifdef LIBMESH_HAVE_TECPLOT_API
30 : extern "C" {
31 : # include <TECIO.h>
32 : }
33 : #endif
34 :
35 : // C++ includes
36 : #include <fstream>
37 : #include <iomanip>
38 : #include <sstream>
39 :
40 : namespace libMesh
41 : {
42 :
43 :
44 : //--------------------------------------------------------
45 : // Macros for handling Tecplot API data
46 :
47 : #ifdef LIBMESH_HAVE_TECPLOT_API
48 :
49 : namespace
50 : {
51 : class TecplotMacros
52 : {
53 : public:
54 : TecplotMacros(const dof_id_type n_nodes,
55 : const unsigned int n_vars,
56 : const dof_id_type n_cells,
57 : const unsigned int n_vert);
58 : float & nd(const std::size_t i, const std::size_t j);
59 : int & cd(const std::size_t i, const std::size_t j);
60 : std::vector<float> nodalData;
61 : std::vector<int> connData;
62 : //float* nodalData;
63 : //int * connData;
64 :
65 : void set_n_cells (const dof_id_type nc);
66 :
67 : const dof_id_type n_nodes;
68 : const unsigned int n_vars;
69 : dof_id_type n_cells;
70 : const unsigned int n_vert;
71 : };
72 :
73 :
74 :
75 : inline
76 : TecplotMacros::TecplotMacros(const dof_id_type nn,
77 : const unsigned int nvar,
78 : const dof_id_type nc,
79 : const unsigned int nvrt) :
80 : n_nodes(nn),
81 : n_vars(nvar),
82 : n_cells(nc),
83 : n_vert(nvrt)
84 : {
85 : nodalData.resize(n_nodes*n_vars);
86 : connData.resize(n_cells*n_vert);
87 : }
88 :
89 :
90 :
91 : inline
92 : float & TecplotMacros::nd(const std::size_t i, const std::size_t j)
93 : {
94 : return nodalData[(i)*(n_nodes) + (j)];
95 : }
96 :
97 :
98 :
99 : inline
100 : int & TecplotMacros::cd(const std::size_t i, const std::size_t j)
101 : {
102 : return connData[(i) + (j)*(n_vert)];
103 : }
104 :
105 :
106 : inline
107 : void TecplotMacros::set_n_cells (const dof_id_type nc)
108 : {
109 : n_cells = nc;
110 : connData.resize(n_cells*n_vert);
111 : }
112 : }
113 :
114 : #endif
115 : //--------------------------------------------------------
116 :
117 :
118 :
119 : // ------------------------------------------------------------
120 : // TecplotIO members
121 1293 : TecplotIO::TecplotIO (const MeshBase & mesh_in,
122 : const bool binary_in,
123 : const double time_in,
124 1293 : const int strand_offset_in) :
125 : MeshOutput<MeshBase> (mesh_in),
126 201 : _binary (binary_in),
127 201 : _time (time_in),
128 201 : _strand_offset (strand_offset_in),
129 402 : _zone_title ("zone"),
130 1293 : _ascii_append(false)
131 : {
132 : // Gather a list of subdomain ids in the mesh.
133 : // We must do this now, while we have every
134 : // processor's attention
135 : // (some of the write methods only execute on processor 0).
136 1293 : mesh_in.subdomain_ids (_subdomain_ids);
137 1293 : }
138 :
139 :
140 :
141 1116 : bool & TecplotIO::binary ()
142 : {
143 1116 : return _binary;
144 : }
145 :
146 :
147 :
148 0 : double & TecplotIO::time ()
149 : {
150 0 : return _time;
151 : }
152 :
153 :
154 :
155 0 : int & TecplotIO::strand_offset ()
156 : {
157 0 : return _strand_offset;
158 : }
159 :
160 :
161 :
162 0 : std::string & TecplotIO::zone_title ()
163 : {
164 0 : return _zone_title;
165 : }
166 :
167 :
168 0 : bool & TecplotIO::ascii_append ()
169 : {
170 0 : return _ascii_append;
171 : }
172 :
173 :
174 1080 : void TecplotIO::write (const std::string & fname)
175 : {
176 1620 : if (this->mesh().processor_id() == 0)
177 : {
178 1080 : if (this->binary())
179 0 : this->write_binary (fname);
180 : else
181 1080 : this->write_ascii (fname);
182 : }
183 1080 : }
184 :
185 :
186 :
187 213 : void TecplotIO::write_nodal_data (const std::string & fname,
188 : const std::vector<Number> & soln,
189 : const std::vector<std::string> & names)
190 : {
191 12 : LOG_SCOPE("write_nodal_data()", "TecplotIO");
192 :
193 219 : if (this->mesh().processor_id() == 0)
194 : {
195 36 : if (this->binary())
196 0 : this->write_binary (fname, &soln, &names);
197 : else
198 36 : this->write_ascii (fname, &soln, &names);
199 : }
200 213 : }
201 :
202 :
203 :
204 1116 : unsigned TecplotIO::elem_dimension()
205 : {
206 : // Get a constant reference to the mesh.
207 1086 : const MeshBase & the_mesh = MeshOutput<MeshBase>::mesh();
208 :
209 1659 : std::vector<unsigned> elem_dims(3);
210 :
211 : // Loop over all the elements and mark the proper dimension entry in
212 : // the elem_dims vector.
213 323019 : for (const auto & elem : the_mesh.active_element_ptr_range())
214 176880 : elem_dims[elem->dim() - 1] = 1;
215 :
216 : // Detect and disallow (for now) the writing of mixed dimension meshes.
217 1116 : libmesh_error_msg_if(std::count(elem_dims.begin(), elem_dims.end(), 1) > 1,
218 : "Error, cannot write Mesh with mixed element dimensions to Tecplot file!");
219 :
220 1116 : if (elem_dims[0])
221 0 : return 1;
222 1116 : else if (elem_dims[1])
223 3 : return 2;
224 1080 : else if (elem_dims[2])
225 540 : return 3;
226 : else
227 0 : libmesh_error_msg("No 1, 2, or 3D elements detected!");
228 : }
229 :
230 :
231 :
232 1116 : void TecplotIO::write_ascii (const std::string & fname,
233 : const std::vector<Number> * v,
234 : const std::vector<std::string> * solution_names)
235 : {
236 : // Should only do this on processor 0!
237 543 : libmesh_assert_equal_to (this->mesh().processor_id(), 0);
238 :
239 : // Create an output stream, possibly in append mode.
240 2775 : std::ofstream out_stream(fname.c_str(), _ascii_append ? std::ofstream::app : std::ofstream::out);
241 :
242 : // Make sure it opened correctly
243 1116 : if (!out_stream.good())
244 0 : libmesh_file_error(fname.c_str());
245 :
246 : // Get a constant reference to the mesh.
247 1086 : const MeshBase & the_mesh = MeshOutput<MeshBase>::mesh();
248 :
249 : // Write header to stream
250 : {
251 : {
252 : // TODO: We used to print out the SVN revision here when we did keyword expansions...
253 : out_stream << "# For a description of the Tecplot format see the Tecplot User's guide.\n"
254 1116 : << "#\n";
255 : }
256 :
257 1116 : out_stream << "Variables=x,y,z";
258 :
259 1116 : if (solution_names != nullptr)
260 72 : for (const auto & val : *solution_names)
261 : {
262 : #ifdef LIBMESH_USE_REAL_NUMBERS
263 :
264 : // Write variable names for real variables
265 30 : out_stream << "," << val;
266 :
267 : #else
268 :
269 : // Write variable names for complex variables
270 : out_stream << "," << "r_" << val
271 : << "," << "i_" << val
272 9 : << "," << "a_" << val;
273 :
274 : #endif
275 : }
276 :
277 1116 : out_stream << '\n';
278 :
279 1689 : out_stream << "Zone f=fepoint, n=" << the_mesh.n_nodes() << ", e=" << the_mesh.n_active_sub_elem();
280 :
281 : // We cannot choose the element type simply based on the mesh
282 : // dimension... there might be 1D elements living in a 3D mesh.
283 : // So look at the elements which are actually in the Mesh, and
284 : // choose either "lineseg", "quadrilateral", or "brick" depending
285 : // on if the elements are 1, 2, or 3D.
286 :
287 : // Write the element type we've determined to the header.
288 1116 : out_stream << ", et=";
289 :
290 1116 : switch (this->elem_dimension())
291 : {
292 0 : case 1:
293 0 : out_stream << "lineseg";
294 0 : break;
295 36 : case 2:
296 36 : out_stream << "quadrilateral";
297 3 : break;
298 1080 : case 3:
299 1080 : out_stream << "brick";
300 540 : break;
301 0 : default:
302 0 : libmesh_error_msg("Unsupported element dimension: " << this->elem_dimension());
303 : }
304 :
305 : // Output the time in the header
306 1116 : out_stream << ", t=\"T " << _time << "\"";
307 :
308 : // Use default mesh color = black
309 1116 : out_stream << ", c=black\n";
310 :
311 : } // finished writing header
312 :
313 115715 : for (auto i : make_range(the_mesh.n_nodes()))
314 : {
315 : // Print the point without a newline
316 114572 : the_mesh.point(i).write_unformatted(out_stream, false);
317 :
318 114572 : if ((v != nullptr) && (solution_names != nullptr))
319 : {
320 18054 : const std::size_t n_vars = solution_names->size();
321 :
322 :
323 216648 : for (std::size_t c=0; c<n_vars; c++)
324 : {
325 : #ifdef LIBMESH_USE_REAL_NUMBERS
326 : // Write real data
327 99297 : out_stream << std::setprecision(this->ascii_precision())
328 108324 : << (*v)[i*n_vars + c] << " ";
329 :
330 : #else
331 : // Write complex data
332 9027 : out_stream << std::setprecision(this->ascii_precision())
333 9027 : << (*v)[i*n_vars + c].real() << " "
334 9027 : << (*v)[i*n_vars + c].imag() << " "
335 9027 : << std::abs((*v)[i*n_vars + c]) << " ";
336 :
337 : #endif
338 : }
339 : }
340 :
341 : // Write a new line after the data for this node
342 114572 : out_stream << '\n';
343 : }
344 :
345 323019 : for (const auto & elem : the_mesh.active_element_ptr_range())
346 176880 : elem->write_connectivity(out_stream, TECPLOT);
347 1116 : }
348 :
349 :
350 :
351 0 : void TecplotIO::write_binary (const std::string & fname,
352 : const std::vector<Number> * vec,
353 : const std::vector<std::string> * solution_names)
354 : {
355 : //-----------------------------------------------------------
356 : // Call the ASCII output function if configure did not detect
357 : // the Tecplot binary API
358 : #ifndef LIBMESH_HAVE_TECPLOT_API
359 :
360 0 : libMesh::err << "WARNING: Tecplot Binary files require the Tecplot API." << std::endl
361 0 : << "Continuing with ASCII output."
362 0 : << std::endl;
363 :
364 0 : if (this->mesh().processor_id() == 0)
365 0 : this->write_ascii (fname, vec, solution_names);
366 0 : return;
367 :
368 :
369 :
370 : //------------------------------------------------------------
371 : // New binary formats, time aware and whatnot
372 : #elif defined(LIBMESH_HAVE_TECPLOT_API_112)
373 :
374 : // Get a constant reference to the mesh.
375 : const MeshBase & the_mesh = MeshOutput<MeshBase>::mesh();
376 :
377 : // Required variables
378 : std::string tecplot_variable_names;
379 : int
380 : ierr = 0,
381 : file_type = 0, // full
382 : is_double = 0,
383 : #ifdef DEBUG
384 : tec_debug = 1,
385 : #else
386 : tec_debug = 0,
387 : #endif
388 : cell_type = -1,
389 : nn_per_elem = -1;
390 :
391 : switch (this->elem_dimension())
392 : {
393 : case 1:
394 : cell_type = 1; // FELINESEG
395 : nn_per_elem = 2;
396 : break;
397 :
398 : case 2:
399 : cell_type = 3; // FEQUADRILATERAL
400 : nn_per_elem = 4;
401 : break;
402 :
403 : case 3:
404 : cell_type = 5; // FEBRICK
405 : nn_per_elem = 8;
406 : break;
407 :
408 : default:
409 : libmesh_error_msg("Unsupported element dimension: " << this->elem_dimension());
410 : }
411 :
412 : // Build a string containing all the variable names to pass to Tecplot
413 : {
414 : tecplot_variable_names += "x, y, z";
415 :
416 : if (solution_names != nullptr)
417 : {
418 : for (const auto & val : *solution_names)
419 : {
420 : #ifdef LIBMESH_USE_REAL_NUMBERS
421 :
422 : tecplot_variable_names += ", ";
423 : tecplot_variable_names += val;
424 :
425 : #else
426 :
427 : tecplot_variable_names += ", ";
428 : tecplot_variable_names += "r_";
429 : tecplot_variable_names += val;
430 : tecplot_variable_names += ", ";
431 : tecplot_variable_names += "i_";
432 : tecplot_variable_names += val;
433 : tecplot_variable_names += ", ";
434 : tecplot_variable_names += "a_";
435 : tecplot_variable_names += val;
436 :
437 : #endif
438 : }
439 : }
440 : }
441 :
442 : // Instantiate a TecplotMacros interface. In 2D the most nodes per
443 : // face should be 4, in 3D it's 8.
444 :
445 :
446 : TecplotMacros tm(the_mesh.n_nodes(),
447 : #ifdef LIBMESH_USE_REAL_NUMBERS
448 : (3 + ((solution_names == nullptr) ? 0 :
449 : cast_int<unsigned int>(solution_names->size()))),
450 : #else
451 : (3 + 3*((solution_names == nullptr) ? 0 :
452 : cast_int<unsigned int>(solution_names->size()))),
453 : #endif
454 : the_mesh.n_active_sub_elem(),
455 : nn_per_elem
456 : );
457 :
458 :
459 : // Copy the nodes and data to the TecplotMacros class. Note that we store
460 : // everything as a float here since the eye doesn't require a double to
461 : // understand what is going on
462 : for (auto v : make_range(the_mesh.n_nodes()))
463 : {
464 : tm.nd(0,v) = static_cast<float>(the_mesh.point(v)(0));
465 : tm.nd(1,v) = static_cast<float>(the_mesh.point(v)(1));
466 : tm.nd(2,v) = static_cast<float>(the_mesh.point(v)(2));
467 :
468 : if ((vec != nullptr) &&
469 : (solution_names != nullptr))
470 : {
471 : const std::size_t n_vars = solution_names->size();
472 :
473 : for (std::size_t c=0; c<n_vars; c++)
474 : {
475 : #ifdef LIBMESH_USE_REAL_NUMBERS
476 :
477 : tm.nd((3+c),v) = static_cast<float>((*vec)[v*n_vars + c]);
478 : #else
479 : tm.nd((3+3*c),v) = static_cast<float>((*vec)[v*n_vars + c].real());
480 : tm.nd((3+3*c+1),v) = static_cast<float>((*vec)[v*n_vars + c].imag());
481 : tm.nd((3+3*c+2),v) = static_cast<float>(std::abs((*vec)[v*n_vars + c]));
482 : #endif
483 : }
484 : }
485 : }
486 :
487 :
488 : // Initialize the file
489 : ierr = TECINI112 (nullptr,
490 : const_cast<char *>(tecplot_variable_names.c_str()),
491 : const_cast<char *>(fname.c_str()),
492 : const_cast<char *>("."),
493 : &file_type,
494 : &tec_debug,
495 : &is_double);
496 :
497 : if (ierr)
498 : libmesh_file_error(fname);
499 :
500 : // A zone for each subdomain
501 : bool firstzone=true;
502 : for (const auto & sbd_id : _subdomain_ids)
503 : {
504 : // Copy the connectivity for this subdomain
505 : {
506 : unsigned int n_subcells_in_subdomain=0;
507 :
508 : for (const auto & elem :
509 : as_range(the_mesh.active_subdomain_elements_begin(sbd_id),
510 : the_mesh.active_subdomain_elements_end(sbd_id)))
511 : n_subcells_in_subdomain += elem->n_sub_elem();
512 :
513 : // update the connectivity array to include only the elements in this subdomain
514 : tm.set_n_cells (n_subcells_in_subdomain);
515 :
516 : unsigned int te = 0;
517 :
518 : for (const auto & elem :
519 : as_range(the_mesh.active_subdomain_elements_begin(sbd_id),
520 : the_mesh.active_subdomain_elements_end(sbd_id)))
521 : {
522 : std::vector<dof_id_type> conn;
523 : for (auto se : make_range(elem->n_sub_elem()))
524 : {
525 : elem->connectivity(se, TECPLOT, conn);
526 :
527 : for (auto node : index_range(conn))
528 : tm.cd(node,te) = conn[node];
529 :
530 : te++;
531 : }
532 : }
533 : }
534 :
535 :
536 : // Ready to call the Tecplot API for this subdomain
537 : {
538 : int
539 : num_nodes = static_cast<int>(the_mesh.n_nodes()),
540 : num_cells = static_cast<int>(tm.n_cells),
541 : num_faces = 0,
542 : i_cell_max = 0,
543 : j_cell_max = 0,
544 : k_cell_max = 0,
545 : strand_id = std::max(sbd_id, static_cast<subdomain_id_type>(1)) + this->strand_offset(),
546 : parent_zone = 0,
547 : is_block = 1,
548 : num_face_connect = 0,
549 : face_neighbor_mode = 0,
550 : tot_num_face_nodes = 0,
551 : num_connect_boundary_faces = 0,
552 : tot_num_boundary_connect = 0,
553 : share_connect_from_zone=0;
554 :
555 : std::vector<int>
556 : passive_var_list (tm.n_vars, 0),
557 : share_var_from_zone (tm.n_vars, 1); // We only write data for the first zone, all other
558 : // zones will share from this one.
559 :
560 : // get the subdomain name from libMesh, if there is one.
561 : std::string subdomain_name = the_mesh.subdomain_name(sbd_id);
562 : std::ostringstream zone_name;
563 : zone_name << this->zone_title();
564 :
565 : // We will title this
566 : // "{zone_title()}_{subdomain_name}", or
567 : // "{zone_title()}_{subdomain_id}", or
568 : // "{zone_title()}"
569 : if (subdomain_name.size())
570 : {
571 : zone_name << "_";
572 : zone_name << subdomain_name;
573 : }
574 : else if (_subdomain_ids.size() > 1)
575 : {
576 : zone_name << "_";
577 : zone_name << sbd_id;
578 : }
579 :
580 : ierr = TECZNE112 (const_cast<char *>(zone_name.str().c_str()),
581 : &cell_type,
582 : &num_nodes,
583 : &num_cells,
584 : &num_faces,
585 : &i_cell_max,
586 : &j_cell_max,
587 : &k_cell_max,
588 : &_time,
589 : &strand_id,
590 : &parent_zone,
591 : &is_block,
592 : &num_face_connect,
593 : &face_neighbor_mode,
594 : &tot_num_face_nodes,
595 : &num_connect_boundary_faces,
596 : &tot_num_boundary_connect,
597 : passive_var_list.data(),
598 : nullptr, // = all are node centered
599 : (firstzone) ? nullptr : share_var_from_zone.data(),
600 : &share_connect_from_zone);
601 :
602 : if (ierr)
603 : libmesh_file_error(fname);
604 :
605 : // Write *all* the data for the first zone, then share it with the others
606 : if (firstzone)
607 : {
608 : int total = cast_int<int>
609 : #ifdef LIBMESH_USE_REAL_NUMBERS
610 : ((3 + ((solution_names == nullptr) ? 0 : solution_names->size()))*num_nodes);
611 : #else
612 : ((3 + 3*((solution_names == nullptr) ? 0 : solution_names->size()))*num_nodes);
613 : #endif
614 :
615 :
616 : ierr = TECDAT112 (&total,
617 : tm.nodalData.data(),
618 : &is_double);
619 :
620 : if (ierr)
621 : libmesh_file_error(fname);
622 : }
623 :
624 : // Write the connectivity
625 : ierr = TECNOD112 (tm.connData.data());
626 :
627 : if (ierr)
628 : libmesh_file_error(fname);
629 : }
630 :
631 : firstzone = false;
632 : }
633 :
634 : // Done, close the file.
635 : ierr = TECEND112 ();
636 :
637 : if (ierr)
638 : libmesh_file_error(fname);
639 :
640 :
641 :
642 :
643 : //------------------------------------------------------------
644 : // Legacy binary format
645 : #else
646 :
647 : // Get a constant reference to the mesh.
648 : const MeshBase & the_mesh = MeshOutput<MeshBase>::mesh();
649 :
650 : // Tecplot binary output only good for dim=2,3
651 : if (the_mesh.mesh_dimension() == 1)
652 : {
653 : this->write_ascii (fname, vec, solution_names);
654 :
655 : return;
656 : }
657 :
658 : // Required variables
659 : std::string tecplot_variable_names;
660 : int is_double = 0,
661 : tec_debug = 0,
662 : cell_type = ((the_mesh.mesh_dimension()==2) ? (1) : (3));
663 :
664 : // Build a string containing all the variable names to pass to Tecplot
665 : {
666 : tecplot_variable_names += "x, y, z";
667 :
668 : if (solution_names != nullptr)
669 : {
670 : for (const auto & val : *solution_names)
671 : {
672 : #ifdef LIBMESH_USE_REAL_NUMBERS
673 :
674 : tecplot_variable_names += ", ";
675 : tecplot_variable_names += val;
676 :
677 : #else
678 :
679 : tecplot_variable_names += ", ";
680 : tecplot_variable_names += "r_";
681 : tecplot_variable_names += val;
682 : tecplot_variable_names += ", ";
683 : tecplot_variable_names += "i_";
684 : tecplot_variable_names += val;
685 : tecplot_variable_names += ", ";
686 : tecplot_variable_names += "a_";
687 : tecplot_variable_names += val;
688 :
689 : #endif
690 : }
691 : }
692 : }
693 :
694 : // Instantiate a TecplotMacros interface. In 2D the most nodes per
695 : // face should be 4, in 3D it's 8.
696 :
697 :
698 : TecplotMacros tm(cast_int<unsigned int>(the_mesh.n_nodes()),
699 : cast_int<unsigned int>
700 : #ifdef LIBMESH_USE_REAL_NUMBERS
701 : (3 + ((solution_names == nullptr) ? 0 : solution_names->size())),
702 : #else
703 : (3 + 3*((solution_names == nullptr) ? 0 : solution_names->size())),
704 : #endif
705 : cast_int<unsigned int>
706 : (the_mesh.n_active_sub_elem()),
707 : ((the_mesh.mesh_dimension() == 2) ? 4 : 8)
708 : );
709 :
710 :
711 : // Copy the nodes and data to the TecplotMacros class. Note that we store
712 : // everything as a float here since the eye doesn't require a double to
713 : // understand what is going on
714 : for (auto v : make_range(the_mesh.n_nodes()))
715 : {
716 : tm.nd(0,v) = static_cast<float>(the_mesh.point(v)(0));
717 : tm.nd(1,v) = static_cast<float>(the_mesh.point(v)(1));
718 : tm.nd(2,v) = static_cast<float>(the_mesh.point(v)(2));
719 :
720 : if ((vec != nullptr) &&
721 : (solution_names != nullptr))
722 : {
723 : const std::size_t n_vars = solution_names->size();
724 :
725 : for (std::size_t c=0; c<n_vars; c++)
726 : {
727 : #ifdef LIBMESH_USE_REAL_NUMBERS
728 :
729 : tm.nd((3+c),v) = static_cast<float>((*vec)[v*n_vars + c]);
730 : #else
731 : tm.nd((3+3*c),v) = static_cast<float>((*vec)[v*n_vars + c].real());
732 : tm.nd((3+3*c+1),v) = static_cast<float>((*vec)[v*n_vars + c].imag());
733 : tm.nd((3+3*c+2),v) = static_cast<float>(std::abs((*vec)[v*n_vars + c]));
734 : #endif
735 : }
736 : }
737 : }
738 :
739 :
740 : // Copy the connectivity
741 : {
742 : unsigned int te = 0;
743 :
744 : for (const auto & elem : the_mesh.active_element_ptr_range())
745 : {
746 : std::vector<dof_id_type> conn;
747 : for (auto se : make_range(elem->n_sub_elem()))
748 : {
749 : elem->connectivity(se, TECPLOT, conn);
750 :
751 : for (auto node : index_range(conn))
752 : tm.cd(node,te) = conn[node];
753 :
754 : te++;
755 : }
756 : }
757 : }
758 :
759 :
760 : // Ready to call the Tecplot API
761 : {
762 : int ierr = 0,
763 : num_nodes = static_cast<int>(the_mesh.n_nodes()),
764 : num_cells = static_cast<int>(the_mesh.n_active_sub_elem());
765 :
766 :
767 : ierr = TECINI (nullptr,
768 : (char *) tecplot_variable_names.c_str(),
769 : (char *) fname.c_str(),
770 : (char *) ".",
771 : &tec_debug,
772 : &is_double);
773 :
774 : if (ierr)
775 : libmesh_file_error(fname);
776 :
777 :
778 : ierr = TECZNE (nullptr,
779 : &num_nodes,
780 : &num_cells,
781 : &cell_type,
782 : (char *) "FEBLOCK",
783 : nullptr);
784 :
785 : if (ierr)
786 : libmesh_file_error(fname);
787 :
788 :
789 : int total =
790 : #ifdef LIBMESH_USE_REAL_NUMBERS
791 : ((3 + ((solution_names == nullptr) ? 0 : solution_names->size()))*num_nodes);
792 : #else
793 : ((3 + 3*((solution_names == nullptr) ? 0 : solution_names->size()))*num_nodes);
794 : #endif
795 :
796 :
797 : ierr = TECDAT (&total,
798 : tm.nodalData.data(),
799 : &is_double);
800 :
801 : if (ierr)
802 : libmesh_file_error(fname);
803 :
804 : ierr = TECNOD (tm.connData.data());
805 :
806 : if (ierr)
807 : libmesh_file_error(fname);
808 :
809 : ierr = TECEND ();
810 :
811 : if (ierr)
812 : libmesh_file_error(fname);
813 : }
814 :
815 : #endif
816 : }
817 :
818 : } // namespace libMesh
|