libMesh
nemesis_io_helper.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 // C++ headers
20 #include <iomanip>
21 #include <set>
22 #include <sstream>
23 
24 // Libmesh headers
25 #include "libmesh/dof_map.h"
26 #include "libmesh/equation_systems.h"
27 #include "libmesh/nemesis_io_helper.h"
28 #include "libmesh/node.h"
29 #include "libmesh/elem.h"
30 #include "libmesh/boundary_info.h"
31 #include "libmesh/parallel.h"
32 #include "libmesh/mesh_base.h"
33 #include "libmesh/numeric_vector.h"
34 #include "libmesh/int_range.h"
35 #include "libmesh/utility.h"
36 
37 #if defined(LIBMESH_HAVE_NEMESIS_API) && defined(LIBMESH_HAVE_EXODUS_API)
38 
39 namespace libMesh
40 {
41 
42 
43 // Initialize the various integer members to zero. We can check
44 // these later to see if they've been properly initialized...
45 // The parent ExodusII_IO_Helper is created with the run_only_on_proc0
46 // flag set to false, so that we can make use of its functionality
47 // on multiple processors.
49  bool verbose_in, bool single_precision) :
50  ExodusII_IO_Helper(parent, verbose_in, /*run_only_on_proc0=*/false, /*single_precision=*/single_precision),
51  nemesis_err_flag(0),
52  num_nodes_global(0),
53  num_elems_global(0),
54  num_elem_blks_global(0),
55  num_node_sets_global(0),
56  num_side_sets_global(0),
57  num_proc(0),
58  num_proc_in_file(0),
59  ftype('\0'),
60  num_internal_nodes(0),
61  num_border_nodes(0),
62  num_external_nodes(0),
63  num_internal_elems(0),
64  num_border_elems(0),
65  num_node_cmaps(0),
66  num_elem_cmaps(0)
67 {
68  // Warn about using untested code!
69  libmesh_experimental();
70 }
71 
72 
74 {
75  // Our destructor is called from Nemesis_IO. We close the Exodus file here since we have
76  // responsibility for managing the file's lifetime. Only call ex_update() if the file was
77  // opened for writing!
78  if (this->opened_for_writing)
79  {
80  this->ex_err = exII::ex_update(this->ex_id);
81  EX_EXCEPTIONLESS_CHECK_ERR(ex_err, "Error flushing buffers to file.");
82  }
83  this->close();
84 }
85 
86 
87 
89 {
91  Nemesis::ne_get_init_global(ex_id,
97  EX_CHECK_ERR(nemesis_err_flag, "Error reading initial global data!");
98 
99  if (verbose)
100  {
101  libMesh::out << "[" << this->processor_id() << "] " << "num_nodes_global=" << num_nodes_global << std::endl;
102  libMesh::out << "[" << this->processor_id() << "] " << "num_elems_global=" << num_elems_global << std::endl;
103  libMesh::out << "[" << this->processor_id() << "] " << "num_elem_blks_global=" << num_elem_blks_global << std::endl;
104  libMesh::out << "[" << this->processor_id() << "] " << "num_node_sets_global=" << num_node_sets_global << std::endl;
105  libMesh::out << "[" << this->processor_id() << "] " << "num_side_sets_global=" << num_side_sets_global << std::endl;
106  }
107 }
108 
109 
110 
112 {
113  if (num_side_sets_global > 0)
114  {
117 
118  // df = "distribution factor", not really sure what that is. I don't yet have a file
119  // which has distribution factors so I guess we'll worry about it later...
121 
123  Nemesis::ne_get_ss_param_global(ex_id,
124  global_sideset_ids.data(),
125  num_global_side_counts.data(),
127  EX_CHECK_ERR(nemesis_err_flag, "Error reading global sideset parameters!");
128 
129  if (verbose)
130  {
131  libMesh::out << "[" << this->processor_id() << "] " << "Global Sideset IDs, Side Counts, and DF counts:" << std::endl;
132  for (auto bn : index_range(global_sideset_ids))
133  {
134  libMesh::out << " [" << this->processor_id() << "] "
135  << "global_sideset_ids["<<bn<<"]=" << global_sideset_ids[bn]
136  << ", num_global_side_counts["<<bn<<"]=" << num_global_side_counts[bn]
137  << ", num_global_side_df_counts["<<bn<<"]=" << num_global_side_df_counts[bn]
138  << std::endl;
139  }
140  }
141  }
142 }
143 
144 
145 
146 
148 {
149  if (num_node_sets_global > 0)
150  {
154 
156  Nemesis::ne_get_ns_param_global(ex_id,
157  global_nodeset_ids.data(),
158  num_global_node_counts.data(),
160  EX_CHECK_ERR(nemesis_err_flag, "Error reading global nodeset parameters!");
161 
162  if (verbose)
163  {
164  libMesh::out << "[" << this->processor_id() << "] " << "Global Nodeset IDs, Node Counts, and DF counts:" << std::endl;
165  for (auto bn : index_range(global_nodeset_ids))
166  {
167  libMesh::out << " [" << this->processor_id() << "] "
168  << "global_nodeset_ids["<<bn<<"]=" << global_nodeset_ids[bn]
169  << ", num_global_node_counts["<<bn<<"]=" << num_global_node_counts[bn]
170  << ", num_global_node_df_counts["<<bn<<"]=" << num_global_node_df_counts[bn]
171  << std::endl;
172  }
173  }
174  }
175 }
176 
177 
178 
180 {
183 
184  if (num_elem_blks_global > 0)
185  {
187  Nemesis::ne_get_eb_info_global(ex_id,
188  global_elem_blk_ids.data(),
189  global_elem_blk_cnts.data());
190  EX_CHECK_ERR(nemesis_err_flag, "Error reading global element block info!");
191  }
192 
193  if (verbose)
194  {
195  libMesh::out << "[" << this->processor_id() << "] " << "Global Element Block IDs and Counts:" << std::endl;
196  for (auto bn : index_range(global_elem_blk_ids))
197  {
198  libMesh::out << " [" << this->processor_id() << "] "
199  << "global_elem_blk_ids[" << bn << "]=" << global_elem_blk_ids[bn]
200  << ", global_elem_blk_cnts[" << bn << "]=" << global_elem_blk_cnts[bn]
201  << std::endl;
202  }
203  }
204 }
205 
206 
207 
209 {
211  Nemesis::ne_get_init_info(ex_id,
212  &num_proc,
214  &ftype);
215  EX_CHECK_ERR(nemesis_err_flag, "Error reading initial info!");
216 
217  if (verbose)
218  {
219  libMesh::out << "[" << this->processor_id() << "] " << "num_proc=" << num_proc << std::endl;
220  libMesh::out << "[" << this->processor_id() << "] " << "num_proc_in_file=" << num_proc_in_file << std::endl;
221  libMesh::out << "[" << this->processor_id() << "] " << "ftype=" << ftype << std::endl;
222  }
223 }
224 
225 
226 
228 {
230  Nemesis::ne_get_loadbal_param(ex_id,
238  this->processor_id() // The ID of the processor for which info is to be read
239  );
240  EX_CHECK_ERR(nemesis_err_flag, "Error reading load balance parameters!");
241 
242 
243  if (verbose)
244  {
245  libMesh::out << "[" << this->processor_id() << "] " << "num_internal_nodes=" << num_internal_nodes << std::endl;
246  libMesh::out << "[" << this->processor_id() << "] " << "num_border_nodes=" << num_border_nodes << std::endl;
247  libMesh::out << "[" << this->processor_id() << "] " << "num_external_nodes=" << num_external_nodes << std::endl;
248  libMesh::out << "[" << this->processor_id() << "] " << "num_internal_elems=" << num_internal_elems << std::endl;
249  libMesh::out << "[" << this->processor_id() << "] " << "num_border_elems=" << num_border_elems << std::endl;
250  libMesh::out << "[" << this->processor_id() << "] " << "num_node_cmaps=" << num_node_cmaps << std::endl;
251  libMesh::out << "[" << this->processor_id() << "] " << "num_elem_cmaps=" << num_elem_cmaps << std::endl;
252  }
253 }
254 
255 
256 
258 {
260  elem_mapb.resize(num_border_elems);
261 
263  Nemesis::ne_get_elem_map(ex_id,
264  elem_mapi.empty() ? nullptr : elem_mapi.data(),
265  elem_mapb.empty() ? nullptr : elem_mapb.data(),
266  this->processor_id()
267  );
268  EX_CHECK_ERR(nemesis_err_flag, "Error reading element maps!");
269 
270 
271  if (verbose)
272  {
273  libMesh::out << "[" << this->processor_id() << "] elem_mapi[i] = ";
274  for (auto i : IntRange<int>(0, num_internal_elems-1))
275  libMesh::out << elem_mapi[i] << ", ";
276  libMesh::out << "... " << elem_mapi.back() << std::endl;
277 
278  libMesh::out << "[" << this->processor_id() << "] elem_mapb[i] = ";
279  for (auto i : IntRange<int>(0, std::min(10, num_border_elems-1)))
280  libMesh::out << elem_mapb[i] << ", ";
281  libMesh::out << "... " << elem_mapb.back() << std::endl;
282  }
283 }
284 
285 
286 
287 
289 {
291  node_mapb.resize(num_border_nodes);
293 
295  Nemesis::ne_get_node_map(ex_id,
296  node_mapi.empty() ? nullptr : node_mapi.data(),
297  node_mapb.empty() ? nullptr : node_mapb.data(),
298  node_mape.empty() ? nullptr : node_mape.data(),
299  this->processor_id()
300  );
301  EX_CHECK_ERR(nemesis_err_flag, "Error reading node maps!");
302 
303  if (verbose)
304  {
305  // Remark: The Exodus/Nemesis node numbering is always (?) 1-based! This means the first interior node id will
306  // always be == 1.
307  libMesh::out << "[" << this->processor_id() << "] " << "first interior node id=" << node_mapi[0] << std::endl;
308  libMesh::out << "[" << this->processor_id() << "] " << "last interior node id=" << node_mapi.back() << std::endl;
309 
310  libMesh::out << "[" << this->processor_id() << "] " << "first boundary node id=" << node_mapb[0] << std::endl;
311  libMesh::out << "[" << this->processor_id() << "] " << "last boundary node id=" << node_mapb.back() << std::endl;
312 
313  // The number of external nodes is sometimes zero, don't try to access
314  // node_mape.back() in this case!
315  if (num_external_nodes > 0)
316  {
317  libMesh::out << "[" << this->processor_id() << "] " << "first external node id=" << node_mape[0] << std::endl;
318  libMesh::out << "[" << this->processor_id() << "] " << "last external node id=" << node_mape.back() << std::endl;
319  }
320  }
321 }
322 
323 
324 
326 {
331 
333  Nemesis::ne_get_cmap_params(ex_id,
334  node_cmap_ids.empty() ? nullptr : node_cmap_ids.data(),
335  node_cmap_node_cnts.empty() ? nullptr : node_cmap_node_cnts.data(),
336  elem_cmap_ids.empty() ? nullptr : elem_cmap_ids.data(),
337  elem_cmap_elem_cnts.empty() ? nullptr : elem_cmap_elem_cnts.data(),
338  this->processor_id());
339  EX_CHECK_ERR(nemesis_err_flag, "Error reading cmap parameters!");
340 
341 
342  if (verbose)
343  {
344  libMesh::out << "[" << this->processor_id() << "] ";
345  for (auto i : index_range(node_cmap_ids))
346  libMesh::out << "node_cmap_ids[" << i << "]=" << node_cmap_ids[i] << " ";
347  libMesh::out << std::endl;
348 
349  libMesh::out << "[" << this->processor_id() << "] ";
350  for (auto i : index_range(node_cmap_node_cnts))
351  libMesh::out << "node_cmap_node_cnts[" << i << "]=" << node_cmap_node_cnts[i] << " ";
352  libMesh::out << std::endl;
353 
354  libMesh::out << "[" << this->processor_id() << "] ";
355  for (auto i : index_range(elem_cmap_ids))
356  libMesh::out << "elem_cmap_ids[" << i << "]=" << elem_cmap_ids[i] << " ";
357  libMesh::out << std::endl;
358 
359  libMesh::out << "[" << this->processor_id() << "] ";
360  for (auto i : index_range(elem_cmap_elem_cnts))
361  libMesh::out << "elem_cmap_elem_cnts[" << i << "]=" << elem_cmap_elem_cnts[i] << " ";
362  libMesh::out << std::endl;
363  }
364 }
365 
366 
367 
369 {
372 
373  for (auto i : index_range(node_cmap_node_ids))
374  {
377 
378  // Don't call ne_get_node_cmap() if there is nothing there to
379  // get, Nemesis throws an error in this case.
380  if (node_cmap_node_cnts[i] > 0)
381  {
383  Nemesis::ne_get_node_cmap(ex_id,
384  node_cmap_ids[i],
387  this->processor_id());
388  EX_CHECK_ERR(nemesis_err_flag, "Error reading node cmap node and processor ids!");
389  }
390 
391  if (verbose)
392  {
393  libMesh::out << "[" << this->processor_id() << "] node_cmap_node_ids[" << i << "]=";
394  for (const auto & dof : node_cmap_node_ids[i])
395  libMesh::out << dof << " ";
396  libMesh::out << std::endl;
397 
398  // This is basically a vector, all entries of which are = node_cmap_ids[i]
399  // Not sure if it's always guaranteed to be that or what...
400  libMesh::out << "[" << this->processor_id() << "] node_cmap_proc_ids[" << i << "]=";
401  for (const auto & dof : node_cmap_proc_ids[i])
402  libMesh::out << dof << " ";
403  libMesh::out << std::endl;
404  }
405  }
406 }
407 
408 
409 
411 {
415 
416  for (auto i : index_range(elem_cmap_elem_ids))
417  {
421 
422  if (elem_cmap_elem_cnts[i] > 0)
423  {
425  Nemesis::ne_get_elem_cmap(ex_id,
426  elem_cmap_ids[i],
430  this->processor_id());
431  EX_CHECK_ERR(nemesis_err_flag, "Error reading elem cmap elem, side, and processor ids!");
432  }
433 
434  if (verbose)
435  {
436  libMesh::out << "[" << this->processor_id() << "] elem_cmap_elem_ids[" << i << "]=";
437  for (const auto & dof : elem_cmap_elem_ids[i])
438  libMesh::out << dof << " ";
439  libMesh::out << std::endl;
440 
441  // These must be the (local) side IDs (in the ExodusII face numbering scheme)
442  // of the sides shared across processors.
443  libMesh::out << "[" << this->processor_id() << "] elem_cmap_side_ids[" << i << "]=";
444  for (const auto & dof : elem_cmap_side_ids[i])
445  libMesh::out << dof << " ";
446  libMesh::out << std::endl;
447 
448  // This is basically a vector, all entries of which are = elem_cmap_ids[i]
449  // Not sure if it's always guaranteed to be that or what...
450  libMesh::out << "[" << this->processor_id() << "] elem_cmap_proc_ids[" << i << "]=";
451  for (const auto & dof : elem_cmap_proc_ids[i])
452  libMesh::out << dof << " ";
453  libMesh::out << std::endl;
454  }
455  }
456 }
457 
458 
459 
460 
461 void Nemesis_IO_Helper::put_init_info(unsigned num_proc_in,
462  unsigned num_proc_in_file_in,
463  const char * ftype_in)
464 {
466  Nemesis::ne_put_init_info(ex_id,
467  num_proc_in,
468  num_proc_in_file_in,
469  const_cast<char *>(ftype_in));
470 
471  EX_CHECK_ERR(nemesis_err_flag, "Error writing initial information!");
472 }
473 
474 
475 
476 
478  dof_id_type num_elems_global_in,
479  unsigned num_elem_blks_global_in,
480  unsigned num_node_sets_global_in,
481  unsigned num_side_sets_global_in)
482 {
484  Nemesis::ne_put_init_global(ex_id,
485  num_nodes_global_in,
486  num_elems_global_in,
487  num_elem_blks_global_in,
488  num_node_sets_global_in,
489  num_side_sets_global_in);
490 
491  EX_CHECK_ERR(nemesis_err_flag, "Error writing initial global data!");
492 }
493 
494 
495 
496 void Nemesis_IO_Helper::put_eb_info_global(std::vector<int> & global_elem_blk_ids_in,
497  std::vector<int> & global_elem_blk_cnts_in)
498 {
500  Nemesis::ne_put_eb_info_global(ex_id,
501  global_elem_blk_ids_in.data(),
502  global_elem_blk_cnts_in.data());
503 
504  EX_CHECK_ERR(nemesis_err_flag, "Error writing global element block information!");
505 }
506 
507 
508 
509 
510 void Nemesis_IO_Helper::put_ns_param_global(std::vector<int> & global_nodeset_ids_in,
511  std::vector<int> & num_global_node_counts_in,
512  std::vector<int> & num_global_node_df_counts_in)
513 {
514  // Only add nodesets if there are some
515  if (global_nodeset_ids.size())
516  {
518  Nemesis::ne_put_ns_param_global(ex_id,
519  global_nodeset_ids_in.data(),
520  num_global_node_counts_in.data(),
521  num_global_node_df_counts_in.data());
522  }
523 
524  EX_CHECK_ERR(nemesis_err_flag, "Error writing global nodeset parameters!");
525 }
526 
527 
528 
529 
530 void Nemesis_IO_Helper::put_ss_param_global(std::vector<int> & global_sideset_ids_in,
531  std::vector<int> & num_global_side_counts_in,
532  std::vector<int> & num_global_side_df_counts_in)
533 {
534  // Only add sidesets if there are some
535  if (global_sideset_ids.size())
536  {
538  Nemesis::ne_put_ss_param_global(ex_id,
539  global_sideset_ids_in.data(),
540  num_global_side_counts_in.data(),
541  num_global_side_df_counts_in.data());
542  }
543 
544  EX_CHECK_ERR(nemesis_err_flag, "Error writing global sideset parameters!");
545 }
546 
547 
548 
549 
550 void Nemesis_IO_Helper::put_loadbal_param(unsigned num_internal_nodes_in,
551  unsigned num_border_nodes_in,
552  unsigned num_external_nodes_in,
553  unsigned num_internal_elems_in,
554  unsigned num_border_elems_in,
555  unsigned num_node_cmaps_in,
556  unsigned num_elem_cmaps_in)
557 {
559  Nemesis::ne_put_loadbal_param(ex_id,
560  num_internal_nodes_in,
561  num_border_nodes_in,
562  num_external_nodes_in,
563  num_internal_elems_in,
564  num_border_elems_in,
565  num_node_cmaps_in,
566  num_elem_cmaps_in,
567  this->processor_id());
568 
569  EX_CHECK_ERR(nemesis_err_flag, "Error writing loadbal parameters!");
570 }
571 
572 
573 
574 
575 
576 void Nemesis_IO_Helper::put_cmap_params(std::vector<int> & node_cmap_ids_in,
577  std::vector<int> & node_cmap_node_cnts_in,
578  std::vector<int> & elem_cmap_ids_in,
579  std::vector<int> & elem_cmap_elem_cnts_in)
580 {
581  libmesh_assert(!node_cmap_ids_in.empty());
582  libmesh_assert(!node_cmap_node_cnts_in.empty());
583  libmesh_assert(!elem_cmap_ids_in.empty());
584  libmesh_assert(!elem_cmap_elem_cnts_in.empty());
585 
587  Nemesis::ne_put_cmap_params(ex_id,
588  node_cmap_ids_in.data(),
589  node_cmap_node_cnts_in.data(),
590  elem_cmap_ids_in.data(),
591  elem_cmap_elem_cnts_in.data(),
592  this->processor_id());
593 
594  EX_CHECK_ERR(nemesis_err_flag, "Error writing cmap parameters!");
595 }
596 
597 
598 
599 
600 void Nemesis_IO_Helper::put_node_cmap(std::vector<std::vector<int>> & node_cmap_node_ids_in,
601  std::vector<std::vector<int>> & node_cmap_proc_ids_in)
602 {
603 
604  // Print to screen what we are about to print to Nemesis file
605  if (verbose)
606  {
607  for (auto i : index_range(node_cmap_node_ids_in))
608  {
609  libMesh::out << "[" << this->processor_id() << "] put_node_cmap() : nodes communicated to proc "
610  << this->node_cmap_ids[i]
611  << " = ";
612  for (const auto & node_id : node_cmap_node_ids_in[i])
613  libMesh::out << node_id << " ";
614  libMesh::out << std::endl;
615  }
616 
617  for (auto i : index_range(node_cmap_node_ids_in))
618  {
619  libMesh::out << "[" << this->processor_id() << "] put_node_cmap() : processor IDs = ";
620  for (const auto & proc_id : node_cmap_proc_ids_in[i])
621  libMesh::out << proc_id << " ";
622  libMesh::out << std::endl;
623  }
624  }
625 
626  for (auto i : index_range(node_cmap_node_ids_in))
627  {
628  int * node_ids_ptr = node_cmap_node_ids_in[i].empty() ?
629  nullptr : node_cmap_node_ids_in[i].data();
630  int * proc_ids_ptr = node_cmap_proc_ids_in[i].empty() ?
631  nullptr : node_cmap_proc_ids_in[i].data();
632 
634  Nemesis::ne_put_node_cmap(ex_id, this->node_cmap_ids[i],
635  node_ids_ptr, proc_ids_ptr,
636  this->processor_id());
637 
638  EX_CHECK_ERR(nemesis_err_flag, "Error writing node communication map to file!");
639  }
640 }
641 
642 
643 
644 
645 void Nemesis_IO_Helper::put_node_map(std::vector<int> & node_mapi_in,
646  std::vector<int> & node_mapb_in,
647  std::vector<int> & node_mape_in)
648 {
650  Nemesis::ne_put_node_map(ex_id,
651  node_mapi_in.empty() ? nullptr : node_mapi_in.data(),
652  node_mapb_in.empty() ? nullptr : node_mapb_in.data(),
653  node_mape_in.empty() ? nullptr : node_mape_in.data(),
654  this->processor_id());
655 
656  EX_CHECK_ERR(nemesis_err_flag, "Error writing Nemesis internal and border node maps to file!");
657 }
658 
659 
660 
661 
662 void Nemesis_IO_Helper::put_elem_cmap(std::vector<std::vector<int>> & elem_cmap_elem_ids_in,
663  std::vector<std::vector<int>> & elem_cmap_side_ids_in,
664  std::vector<std::vector<int>> & elem_cmap_proc_ids_in)
665 {
666  for (auto i : index_range(elem_cmap_ids))
667  {
669  Nemesis::ne_put_elem_cmap(ex_id,
670  this->elem_cmap_ids[i],
671  elem_cmap_elem_ids_in[i].data(),
672  elem_cmap_side_ids_in[i].data(),
673  elem_cmap_proc_ids_in[i].data(),
674  this->processor_id());
675 
676  EX_CHECK_ERR(nemesis_err_flag, "Error writing elem communication map to file!");
677  }
678 }
679 
680 
681 
682 
683 void Nemesis_IO_Helper::put_elem_map(std::vector<int> & elem_mapi_in,
684  std::vector<int> & elem_mapb_in)
685 {
687  Nemesis::ne_put_elem_map(ex_id,
688  elem_mapi_in.empty() ? nullptr : elem_mapi_in.data(),
689  elem_mapb_in.empty() ? nullptr : elem_mapb_in.data(),
690  this->processor_id());
691 
692  EX_CHECK_ERR(nemesis_err_flag, "Error writing Nemesis internal and border element maps to file!");
693 }
694 
695 
696 
697 
698 
699 
700 void Nemesis_IO_Helper::put_n_coord(unsigned start_node_num,
701  unsigned num_nodes_in,
702  std::vector<Real> & x_coor,
703  std::vector<Real> & y_coor,
704  std::vector<Real> & z_coor)
705 {
706  if (num_nodes_in)
707  {
709  Nemesis::ne_put_n_coord(ex_id,
710  start_node_num,
711  num_nodes_in,
712  x_coor.empty() ? nullptr : x_coor.data(),
713  y_coor.empty() ? nullptr : y_coor.data(),
714  z_coor.empty() ? nullptr : z_coor.data());
715 
716  EX_CHECK_ERR(nemesis_err_flag, "Error writing coords to file!");
717  }
718 }
719 
720 
721 
722 
723 
724 
725 
726 
727 // Note: we can't reuse the ExodusII_IO_Helper code directly, since it only runs
728 // on processor 0. TODO: We could have the body of this function as a separate
729 // function and then ExodusII_IO_Helper would only call it if on processor 0...
730 void Nemesis_IO_Helper::create(std::string filename)
731 {
732  // Fall back on double precision when necessary since ExodusII
733  // doesn't seem to support long double
734  int
735  comp_ws = 0,
736  io_ws = 0;
737 
738  if (_single_precision)
739  {
740  comp_ws = sizeof(float);
741  io_ws = sizeof(float);
742  }
743  else
744  {
745  comp_ws = cast_int<int>(std::min(sizeof(Real), sizeof(double)));
746  io_ws = cast_int<int>(std::min(sizeof(Real), sizeof(double)));
747  }
748 
749  this->ex_id = exII::ex_create(filename.c_str(), EX_CLOBBER, &comp_ws, &io_ws);
750 
751  EX_CHECK_ERR(ex_id, "Error creating Nemesis mesh file.");
752 
753  if (verbose)
754  libMesh::out << "File created successfully." << std::endl;
755 
756  this->opened_for_writing = true;
757 }
758 
759 
760 
761 
762 
763 
764 
765 
766 void Nemesis_IO_Helper::initialize(std::string title_in, const MeshBase & mesh, bool /*use_discontinuous*/)
767 {
768  // Make sure that the reference passed in is really a DistributedMesh
769  // const DistributedMesh & pmesh = cast_ref<const DistributedMesh &>(mesh);
770  const MeshBase & pmesh = mesh;
771 
772  // According to Nemesis documentation, first call when writing should be to
773  // ne_put_init_info(). Our reader doesn't actually call this, but we should
774  // strive to be as close to a normal nemesis file as possible...
775  this->put_init_info(this->n_processors(), 1, "p");
776 
777 
778  // Gather global "initial" information for Nemesis. This consists of
779  // three parts labeled I, II, and III below...
780 
781  //
782  // I.) Need to compute the number of global element blocks. To be consistent with
783  // Exodus, we also incorrectly associate the number of element blocks with the
784  // number of libmesh subdomains...
785  //
786  this->compute_num_global_elem_blocks(pmesh);
787 
788  //
789  // II.) Determine the global number of nodesets by communication.
790  // This code relies on BoundaryInfo storing side and node
791  // boundary IDs separately at the time they are added to the
792  // BoundaryInfo object.
793  //
794  this->compute_num_global_nodesets(pmesh);
795 
796  //
797  // III.) Need to compute the global number of sidesets by communication:
798  // This code relies on BoundaryInfo storing side and node
799  // boundary IDs separately at the time they are added to the
800  // BoundaryInfo object.
801  //
802  this->compute_num_global_sidesets(pmesh);
803 
804  // Now write the global data obtained in steps I, II, and III to the Nemesis file
805  this->put_init_global(pmesh.parallel_n_nodes(),
806  pmesh.parallel_n_elem(),
807  this->num_elem_blks_global, /* I. */
808  this->num_node_sets_global, /* II. */
809  this->num_side_sets_global /* III. */
810  );
811 
812  // Next, we'll write global element block information to the file. This was already
813  // gathered in step I. above
815  this->global_elem_blk_cnts);
816 
817 
818  // Next, write global nodeset information to the file. This was already gathered in
819  // step II. above.
820  this->num_global_node_df_counts.clear();
821  this->num_global_node_df_counts.resize(this->global_nodeset_ids.size()); // distribution factors all zero...
825 
826 
827  // Next, write global sideset information to the file. This was already gathered in
828  // step III. above.
829  this->num_global_side_df_counts.clear();
830  this->num_global_side_df_counts.resize(this->global_sideset_ids.size()); // distribution factors all zero...
834 
835 
836  // Before we go any further we need to derive consistent node and
837  // element numbering schemes for all local elems and nodes connected
838  // to local elements.
839  //
840  // Must be called *after* the local_subdomain_counts map has been constructed
841  // by the compute_num_global_elem_blocks() function!
842  this->build_element_and_node_maps(pmesh);
843 
844  // Next step is to write "load balance" parameters. Several things need to
845  // be computed first though...
846 
847  // First we'll collect IDs of border nodes.
848  this->compute_border_node_ids(pmesh);
849 
850  // Next we'll collect numbers of internal and border elements, and internal nodes.
851  // Note: "A border node does not a border element make...", that is, just because one
852  // of an element's nodes has been identified as a border node, the element is not
853  // necessarily a border element. It must have a side on the boundary between processors,
854  // i.e. have a face neighbor with a different processor id...
856 
857  // Finally we are ready to write the loadbal information to the file
859  this->num_border_nodes,
860  this->num_external_nodes,
861  this->num_internal_elems,
862  this->num_border_elems,
863  this->num_node_cmaps,
864  this->num_elem_cmaps);
865 
866 
867  // Now we need to compute the "communication map" parameters. These are basically
868  // lists of nodes and elements which need to be communicated between different processors
869  // when the mesh file is read back in.
871 
872  // Do we have communication maps to write? Note that
873  // ne_put_cmap_params expects us to have either *both* node and elem
874  // cmaps or *neither*
875  if (!this->node_cmap_ids.empty() &&
876  !this->node_cmap_node_cnts.empty() &&
877  !this->elem_cmap_ids.empty() &&
878  !this->elem_cmap_elem_cnts.empty())
879  {
880  // Write communication map parameters to file.
881  this->put_cmap_params(this->node_cmap_ids,
882  this->node_cmap_node_cnts,
883  this->elem_cmap_ids,
884  this->elem_cmap_elem_cnts);
885 
886  // Ready the node communication maps. The node IDs which
887  // are communicated are the ones currently stored in
888  // proc_nodes_touched_intersections.
890 
891  // Write the packed node communication vectors to file.
892  this->put_node_cmap(this->node_cmap_node_ids,
893  this->node_cmap_proc_ids);
894 
895  // Ready the node maps. These have nothing to do with communication, they map
896  // the nodes to internal, border, and external nodes in the file.
897  this->compute_node_maps();
898 
899  // Call the Nemesis API to write the node maps to file.
900  this->put_node_map(this->node_mapi,
901  this->node_mapb,
902  this->node_mape);
903 
904  // Ready the element communication maps. This includes border
905  // element IDs, sides which are on the border, and the processors to which
906  // they are to be communicated...
908 
909  // Call the Nemesis API to write the packed element communication maps vectors to file
910  this->put_elem_cmap(this->elem_cmap_elem_ids,
911  this->elem_cmap_side_ids,
912  this->elem_cmap_proc_ids);
913  }
914 
915 
916  // Ready the Nemesis element maps (internal and border) for writing to file.
917  this->compute_element_maps();
918 
919  // Call the Nemesis API to write the internal and border element IDs.
920  this->put_elem_map(this->elem_mapi,
921  this->elem_mapb);
922 
923  // Now write Exodus-specific initialization information, some of which is
924  // different when you are using Nemesis.
925  this->write_exodus_initialization_info(pmesh, title_in);
926 } // end initialize()
927 
928 
929 
930 
931 
932 
934  const std::string & title_in)
935 {
936  // This follows the convention of Exodus: we always write out the mesh as LIBMESH_DIM-dimensional,
937  // even if it is 2D...
938  this->num_dim = LIBMESH_DIM;
939 
940  this->num_elem = static_cast<unsigned int>(std::distance (pmesh.active_local_elements_begin(),
941  pmesh.active_local_elements_end()));
942 
943  // Exodus will also use *global* number of side and node sets,
944  // though it will not write out entries for all of them...
945  this->num_side_sets =
946  cast_int<int>(this->global_sideset_ids.size());
947  this->num_node_sets =
948  cast_int<int>(this->global_nodeset_ids.size());
949 
950  // We need to write the global number of blocks, even though this processor might not have
951  // elements in some of them!
952  this->num_elem_blk = this->num_elem_blks_global;
953 
954  ex_err = exII::ex_put_init(ex_id,
955  title_in.c_str(),
956  this->num_dim,
957  this->num_nodes,
958  this->num_elem,
959  this->num_elem_blk,
960  this->num_node_sets,
961  this->num_side_sets);
962 
963  EX_CHECK_ERR(ex_err, "Error initializing new Nemesis file.");
964 }
965 
966 
967 
968 
969 
971 {
972  // Make sure we don't have any leftover info
973  this->elem_mapi.clear();
974  this->elem_mapb.clear();
975 
976  // Copy set contents into vectors
977  this->elem_mapi.resize(this->internal_elem_ids.size());
978  this->elem_mapb.resize(this->border_elem_ids.size());
979 
980  {
981  unsigned cnt = 0;
982  for (const auto & id : this->internal_elem_ids)
983  this->elem_mapi[cnt++] = libmesh_map_find(libmesh_elem_num_to_exodus, id);
984  }
985 
986  {
987  unsigned cnt = 0;
988  for (const auto & id : this->border_elem_ids)
989  this->elem_mapb[cnt++] = libmesh_map_find(libmesh_elem_num_to_exodus, id);
990  }
991 }
992 
993 
994 
996 {
997  // Make sure there is no leftover information
998  this->elem_cmap_elem_ids.clear();
999  this->elem_cmap_side_ids.clear();
1000  this->elem_cmap_proc_ids.clear();
1001 
1002  // Allocate enough space for all our element maps
1003  this->elem_cmap_elem_ids.resize(this->num_elem_cmaps);
1004  this->elem_cmap_side_ids.resize(this->num_elem_cmaps);
1005  this->elem_cmap_proc_ids.resize(this->num_elem_cmaps);
1006  {
1007  unsigned cnt=0; // Index into vectors
1009  it = this->proc_border_elem_sets.begin(),
1010  end = this->proc_border_elem_sets.end();
1011 
1012  for (; it != end; ++it)
1013  {
1014  // Make sure the current elem_cmap_id matches the index in our map of node intersections
1015  libmesh_assert_equal_to (static_cast<unsigned>(this->elem_cmap_ids[cnt]), it->first);
1016 
1017  // Get reference to the set of IDs to be packed into the vector
1018  std::set<std::pair<unsigned,unsigned>> & elem_set = it->second;
1019 
1020  // Resize the vectors to receive their payload
1021  this->elem_cmap_elem_ids[cnt].resize(elem_set.size());
1022  this->elem_cmap_side_ids[cnt].resize(elem_set.size());
1023  this->elem_cmap_proc_ids[cnt].resize(elem_set.size());
1024 
1025  std::set<std::pair<unsigned,unsigned>>::iterator elem_set_iter = elem_set.begin();
1026 
1027  // Pack the vectors with elem IDs, side IDs, and processor IDs.
1028  for (std::size_t j=0, eceis=this->elem_cmap_elem_ids[cnt].size(); j<eceis; ++j, ++elem_set_iter)
1029  {
1030  this->elem_cmap_elem_ids[cnt][j] =
1031  libmesh_map_find(libmesh_elem_num_to_exodus, elem_set_iter->first);
1032  this->elem_cmap_side_ids[cnt][j] = elem_set_iter->second; // Side ID, this has already been converted above
1033  this->elem_cmap_proc_ids[cnt][j] = it->first; // All have the same processor ID
1034  }
1035 
1036  // increment vector index to go to next processor
1037  cnt++;
1038  }
1039  } // end scope for packing
1040 }
1041 
1042 
1043 
1044 
1045 
1047 {
1048  // Make sure we don't have any leftover information
1049  this->node_mapi.clear();
1050  this->node_mapb.clear();
1051  this->node_mape.clear();
1052 
1053  // Make sure there's enough space to hold all our node IDs
1054  this->node_mapi.resize(this->internal_node_ids.size());
1055  this->node_mapb.resize(this->border_node_ids.size());
1056 
1057  // Copy set contents into vectors
1058  {
1059  unsigned cnt = 0;
1060  for (const auto & id : this->internal_node_ids)
1061  this->node_mapi[cnt++] = libmesh_map_find(libmesh_node_num_to_exodus, id);
1062  }
1063 
1064  {
1065  unsigned cnt=0;
1066  for (const auto & id : this->border_node_ids)
1067  this->node_mapb[cnt++] = libmesh_map_find(libmesh_node_num_to_exodus, id);
1068  }
1069 }
1070 
1071 
1072 
1073 
1074 
1076 {
1077  // Make sure there's no left-over information
1078  this->node_cmap_node_ids.clear();
1079  this->node_cmap_proc_ids.clear();
1080 
1081  libmesh_assert_less_equal
1082  (this->proc_nodes_touched_intersections.size(),
1083  std::size_t(this->num_node_cmaps));
1084 
1085  // Allocate enough space for all our node maps
1086  this->node_cmap_node_ids.resize(this->num_node_cmaps);
1087  this->node_cmap_proc_ids.resize(this->num_node_cmaps);
1088  {
1089  unsigned cnt=0; // Index into vectors
1091  it = this->proc_nodes_touched_intersections.begin(),
1092  end = this->proc_nodes_touched_intersections.end();
1093 
1094  for (; it != end; ++it)
1095  {
1096  // Make sure the current node_cmap_id matches the index in our map of node intersections
1097  libmesh_assert_equal_to (static_cast<unsigned>(this->node_cmap_ids[cnt]), it->first);
1098 
1099  // Get reference to the set of IDs to be packed into the vector.
1100  std::set<unsigned> & node_set = it->second;
1101 
1102  // Resize the vectors to receive their payload
1103  this->node_cmap_node_ids[cnt].resize(node_set.size());
1104  this->node_cmap_proc_ids[cnt].resize(node_set.size());
1105 
1106  std::set<unsigned>::iterator node_set_iter = node_set.begin();
1107 
1108  // Pack the vectors with node IDs and processor IDs.
1109  for (std::size_t j=0, nceis=this->node_cmap_node_ids[cnt].size(); j<nceis; ++j, ++node_set_iter)
1110  {
1111  this->node_cmap_node_ids[cnt][j] =
1112  libmesh_map_find(libmesh_node_num_to_exodus, *node_set_iter);
1113  this->node_cmap_proc_ids[cnt][j] = it->first;
1114  }
1115 
1116  // increment vector index to go to next processor
1117  cnt++;
1118  }
1119  } // end scope for packing
1120 
1121  // Print out the vectors we just packed
1122  if (verbose)
1123  {
1124  for (auto i : index_range(this->node_cmap_node_ids))
1125  {
1126  libMesh::out << "[" << this->processor_id() << "] nodes communicated to proc "
1127  << this->node_cmap_ids[i]
1128  << " = ";
1129  for (const auto & node_id : this->node_cmap_node_ids[i])
1130  libMesh::out << node_id << " ";
1131  libMesh::out << std::endl;
1132  }
1133 
1134  for (const auto & id_vec : this->node_cmap_node_ids)
1135  {
1136  libMesh::out << "[" << this->processor_id() << "] processor ID node communicated to = ";
1137  for (const auto & proc_id : id_vec)
1138  libMesh::out << proc_id << " ";
1139  libMesh::out << std::endl;
1140  }
1141  }
1142 }
1143 
1144 
1145 
1146 
1148 {
1149  // For the nodes, these are the number of entries in the sets in proc_nodes_touched_intersections
1150  // map computed above. Note: this map does not contain self-intersections so we can loop over it
1151  // directly.
1152  this->node_cmap_node_cnts.clear(); // Make sure we don't have any leftover information...
1153  this->node_cmap_ids.clear(); // Make sure we don't have any leftover information...
1154  this->node_cmap_node_cnts.resize(this->num_node_cmaps);
1155  this->node_cmap_ids.resize(this->num_node_cmaps);
1156 
1157  {
1158  unsigned cnt=0; // Index into the vector
1160  it = this->proc_nodes_touched_intersections.begin(),
1161  end = this->proc_nodes_touched_intersections.end();
1162 
1163  for (; it != end; ++it)
1164  {
1165  this->node_cmap_ids[cnt] = it->first; // The ID of the proc we communicate with
1166  this->node_cmap_node_cnts[cnt] = cast_int<int>(it->second.size()); // The number of nodes we communicate
1167  cnt++; // increment vector index!
1168  }
1169  }
1170 
1171  // Print the packed vectors we just filled
1172  if (verbose)
1173  {
1174  libMesh::out << "[" << this->processor_id() << "] node_cmap_node_cnts = ";
1175  for (const auto & node_cnt : node_cmap_node_cnts)
1176  libMesh::out << node_cnt << ", ";
1177  libMesh::out << std::endl;
1178 
1179  libMesh::out << "[" << this->processor_id() << "] node_cmap_ids = ";
1180  for (const auto & node_id : node_cmap_ids)
1181  libMesh::out << node_id << ", ";
1182  libMesh::out << std::endl;
1183  }
1184 
1185  // For the elements, we have not yet computed all this information..
1186  this->elem_cmap_elem_cnts.clear(); // Make sure we don't have any leftover information...
1187  this->elem_cmap_ids.clear(); // Make sure we don't have any leftover information...
1188  this->elem_cmap_elem_cnts.resize(this->num_elem_cmaps);
1189  this->elem_cmap_ids.resize(this->num_elem_cmaps);
1190 
1191  // Pack the elem_cmap_ids and elem_cmap_elem_cnts vectors
1192  {
1193  unsigned cnt=0; // Index into the vectors we're filling
1195  it = this->proc_border_elem_sets.begin(),
1196  end = this->proc_border_elem_sets.end();
1197 
1198  for (; it != end; ++it)
1199  {
1200  this->elem_cmap_ids[cnt] = it->first; // The ID of the proc we communicate with
1201  this->elem_cmap_elem_cnts[cnt] = cast_int<int>(it->second.size()); // The number of elems we communicate to/from that proc
1202  cnt++; // increment vector index!
1203  }
1204  }
1205 
1206  // Print the packed vectors we just filled
1207  if (verbose)
1208  {
1209  libMesh::out << "[" << this->processor_id() << "] elem_cmap_elem_cnts = ";
1210  for (const auto & elem_cnt : elem_cmap_elem_cnts)
1211  libMesh::out << elem_cnt << ", ";
1212  libMesh::out << std::endl;
1213 
1214  libMesh::out << "[" << this->processor_id() << "] elem_cmap_ids = ";
1215  for (const auto & elem_id : elem_cmap_ids)
1216  libMesh::out << elem_id << ", ";
1217  libMesh::out << std::endl;
1218  }
1219 }
1220 
1221 
1222 
1223 
1224 void
1226 {
1227  // Set of all local, active element IDs. After we have identified border element
1228  // IDs, the set_difference between this set and the border_elem_ids set will give us
1229  // the set of internal_elem_ids.
1230  std::set<unsigned> all_elem_ids;
1231 
1232  // A set of processor IDs which elements on this processor have as
1233  // neighbors. The size of this set will determine the number of
1234  // element communication maps in Exodus.
1235  std::set<unsigned> neighboring_processor_ids;
1236 
1237  for (const auto & elem : pmesh.active_local_element_ptr_range())
1238  {
1239  // Add this Elem's ID to all_elem_ids, later we will take the difference
1240  // between this set and the set of border_elem_ids, to get the set of
1241  // internal_elem_ids.
1242  all_elem_ids.insert(elem->id());
1243 
1244  // Will be set to true if element is determined to be a border element
1245  bool is_border_elem = false;
1246 
1247  // Construct a conversion object for this Element. This will help us map
1248  // Libmesh numberings into Nemesis numberings for sides.
1249  const auto & conv = get_conversion(elem->type());
1250 
1251  // Add all this element's node IDs to the set of all node IDs.
1252  // The set of internal_node_ids will be the set difference between
1253  // the set of all nodes and the set of border nodes.
1254  //
1255  // In addition, if any node of a local node is listed in the
1256  // border nodes list, then this element goes into the proc_border_elem_sets.
1257  // Note that there is not a 1:1 correspondence between
1258  // border_elem_ids and the entries which go into proc_border_elem_sets.
1259  // The latter is for communication purposes, ie determining which elements
1260  // should be shared between processors.
1261  for (auto node : elem->node_index_range())
1262  this->nodes_attached_to_local_elems.insert(elem->node_id(node));
1263 
1264  // Loop over element's neighbors, see if it has a neighbor which is off-processor
1265  for (auto n : elem->side_index_range())
1266  {
1267  if (elem->neighbor_ptr(n) != nullptr)
1268  {
1269  unsigned neighbor_proc_id = elem->neighbor_ptr(n)->processor_id();
1270 
1271  // If my neighbor has a different processor ID, I must be a border element.
1272  // Also track the neighboring processor ID if it is are different from our processor ID
1273  if (neighbor_proc_id != this->processor_id())
1274  {
1275  is_border_elem = true;
1276  neighboring_processor_ids.insert(neighbor_proc_id);
1277 
1278  // Convert libmesh side(n) of this element into a side ID for Nemesis
1279  unsigned nemesis_side_id = conv.get_inverse_side_map(n);
1280 
1281  if (verbose)
1282  libMesh::out << "[" << this->processor_id() << "] LibMesh side "
1283  << n
1284  << " mapped to (1-based) Exodus side "
1285  << nemesis_side_id
1286  << std::endl;
1287 
1288  // Add this element's ID and the ID of the side which is on the boundary
1289  // to the set of border elements for this processor.
1290  // Note: if the set does not already exist, this creates it.
1291  this->proc_border_elem_sets[ neighbor_proc_id ].insert( std::make_pair(elem->id(), nemesis_side_id) );
1292  }
1293  }
1294  } // end for loop over neighbors
1295 
1296  // If we're on a border element, add it to the set
1297  if (is_border_elem)
1298  this->border_elem_ids.insert( elem->id() );
1299 
1300  } // end for loop over active local elements
1301 
1302  // Take the set_difference between all elements and border elements to get internal
1303  // element IDs
1304  std::set_difference(all_elem_ids.begin(), all_elem_ids.end(),
1305  this->border_elem_ids.begin(), this->border_elem_ids.end(),
1306  std::inserter(this->internal_elem_ids, this->internal_elem_ids.end()));
1307 
1308  // Take the set_difference between all nodes and border nodes to get internal nodes
1309  std::set_difference(this->nodes_attached_to_local_elems.begin(), this->nodes_attached_to_local_elems.end(),
1310  this->border_node_ids.begin(), this->border_node_ids.end(),
1311  std::inserter(this->internal_node_ids, this->internal_node_ids.end()));
1312 
1313  if (verbose)
1314  {
1315  libMesh::out << "[" << this->processor_id() << "] neighboring_processor_ids = ";
1316  for (const auto & id : neighboring_processor_ids)
1317  libMesh::out << id << " ";
1318  libMesh::out << std::endl;
1319  }
1320 
1321  // The size of the neighboring_processor_ids set should be the number of element communication maps
1322  this->num_elem_cmaps =
1323  cast_int<int>(neighboring_processor_ids.size());
1324 
1325  if (verbose)
1326  libMesh::out << "[" << this->processor_id() << "] "
1327  << "Number of neighboring processor IDs="
1328  << this->num_elem_cmaps
1329  << std::endl;
1330 
1331  if (verbose)
1332  {
1333  // Print out counts of border elements for each processor
1334  for (const auto & pr : proc_border_elem_sets)
1335  {
1336  libMesh::out << "[" << this->processor_id() << "] "
1337  << "Proc "
1338  << pr.first << " communicates "
1339  << pr.second.size() << " elements." << std::endl;
1340  }
1341  }
1342 
1343  // Store the number of internal and border elements, and the number of internal nodes,
1344  // to be written to the Nemesis file.
1345  this->num_internal_elems =
1346  cast_int<int>(this->internal_elem_ids.size());
1347  this->num_border_elems =
1348  cast_int<int>(this->border_elem_ids.size());
1349  this->num_internal_nodes =
1350  cast_int<int>(this->internal_node_ids.size());
1351 
1352  if (verbose)
1353  {
1354  libMesh::out << "[" << this->processor_id() << "] num_internal_nodes=" << this->num_internal_nodes << std::endl;
1355  libMesh::out << "[" << this->processor_id() << "] num_border_nodes=" << this->num_border_nodes << std::endl;
1356  libMesh::out << "[" << this->processor_id() << "] num_border_elems=" << this->num_border_elems << std::endl;
1357  libMesh::out << "[" << this->processor_id() << "] num_internal_elems=" << this->num_internal_elems << std::endl;
1358  }
1359 }
1360 
1361 
1362 
1364 {
1365  // 1.) Get reference to the set of side boundary IDs
1366  std::set<boundary_id_type> global_side_boundary_ids
1367  (pmesh.get_boundary_info().get_side_boundary_ids().begin(),
1368  pmesh.get_boundary_info().get_side_boundary_ids().end());
1369 
1370  // 2.) Gather boundary side IDs from other processors
1371  this->comm().set_union(global_side_boundary_ids);
1372 
1373  // 3.) Now global_side_boundary_ids actually contains a global list of all side boundary IDs
1374  this->num_side_sets_global =
1375  cast_int<int>(global_side_boundary_ids.size());
1376 
1377  // 4.) Pack these sidesets into a vector so they can be written by Nemesis
1378  this->global_sideset_ids.clear(); // Make sure there is no leftover information
1379  this->global_sideset_ids.insert(this->global_sideset_ids.end(),
1380  global_side_boundary_ids.begin(),
1381  global_side_boundary_ids.end());
1382 
1383  if (verbose)
1384  {
1385  libMesh::out << "[" << this->processor_id() << "] global_sideset_ids = ";
1386  for (const auto & id : this->global_sideset_ids)
1387  libMesh::out << id << ", ";
1388  libMesh::out << std::endl;
1389  }
1390 
1391  // We also need global counts of sides in each of the sidesets.
1392  // Build a list of (elem, side, bc) tuples.
1393  typedef std::tuple<dof_id_type, unsigned short int, boundary_id_type> Tuple;
1394  std::vector<Tuple> bc_triples = pmesh.get_boundary_info().build_side_list();
1395 
1396  // Iterators to the beginning and end of the current range.
1397  std::vector<Tuple>::iterator
1398  it = bc_triples.begin(),
1399  new_end = bc_triples.end();
1400 
1401  while (it != new_end)
1402  {
1403  if (pmesh.elem_ref(std::get<0>(*it)).processor_id() != this->processor_id())
1404  {
1405  // Back up the new end iterators to prepare for swap
1406  --new_end;
1407 
1408  // Swap places, the non-local elem will now be "past-the-end"
1409  std::swap (*it, *new_end);
1410  }
1411  else // elem is local, go to next
1412  ++it;
1413  }
1414 
1415  // Erase from "new" end to old.
1416  bc_triples.erase(new_end, bc_triples.end());
1417 
1418  this->num_global_side_counts.clear(); // Make sure we don't have any leftover information
1419  this->num_global_side_counts.resize(this->global_sideset_ids.size());
1420 
1421  // Get the count for each global sideset ID
1422  for (auto i : index_range(global_sideset_ids))
1423  {
1424  int id = global_sideset_ids[i];
1425  this->num_global_side_counts[i] =
1426  cast_int<int>(std::count_if(bc_triples.begin(),
1427  bc_triples.end(),
1428  [id](const Tuple & t)->bool { return std::get<2>(t) == id; }));
1429  }
1430 
1431  if (verbose)
1432  {
1433  libMesh::out << "[" << this->processor_id() << "] num_global_side_counts = ";
1434  for (const auto & cnt : this->num_global_side_counts)
1435  libMesh::out << cnt << ", ";
1436  libMesh::out << std::endl;
1437  }
1438 
1439  // Finally sum up the result
1440  this->comm().sum(this->num_global_side_counts);
1441 
1442  if (verbose)
1443  {
1444  libMesh::out << "[" << this->processor_id() << "] num_global_side_counts = ";
1445  for (const auto & cnt : this->num_global_side_counts)
1446  libMesh::out << cnt << ", ";
1447  libMesh::out << std::endl;
1448  }
1449 }
1450 
1451 
1452 
1453 
1454 
1455 
1457 {
1458  std::set<boundary_id_type> local_node_boundary_ids;
1459 
1460  // 1.) Get reference to the set of node boundary IDs *for this processor*
1461  std::set<boundary_id_type> global_node_boundary_ids
1462  (pmesh.get_boundary_info().get_node_boundary_ids().begin(),
1463  pmesh.get_boundary_info().get_node_boundary_ids().end());
1464 
1465  // Save a copy of the local_node_boundary_ids...
1466  local_node_boundary_ids = global_node_boundary_ids;
1467 
1468  // 2.) Gather boundary node IDs from other processors
1469  this->comm().set_union(global_node_boundary_ids);
1470 
1471  // 3.) Now global_node_boundary_ids actually contains a global list of all node boundary IDs
1472  this->num_node_sets_global =
1473  cast_int<int>(global_node_boundary_ids.size());
1474 
1475  // 4.) Create a vector<int> from the global_node_boundary_ids set
1476  this->global_nodeset_ids.clear();
1477  this->global_nodeset_ids.insert(this->global_nodeset_ids.end(),
1478  global_node_boundary_ids.begin(),
1479  global_node_boundary_ids.end());
1480 
1481  if (verbose)
1482  {
1483  libMesh::out << "[" << this->processor_id() << "] global_nodeset_ids = ";
1484  for (const auto & id : global_nodeset_ids)
1485  libMesh::out << id << ", ";
1486  libMesh::out << std::endl;
1487 
1488  libMesh::out << "[" << this->processor_id() << "] local_node_boundary_ids = ";
1489  for (const auto & id : local_node_boundary_ids)
1490  libMesh::out << id << ", ";
1491  libMesh::out << std::endl;
1492  }
1493 
1494  // 7.) We also need to know the number of nodes which is in each of the nodesets, globally.
1495 
1496  // Build list of (node-id, bc-id) tuples.
1497  typedef std::tuple<dof_id_type, boundary_id_type> Tuple;
1498  std::vector<Tuple> bc_tuples = pmesh.get_boundary_info().build_node_list();
1499 
1500  if (verbose)
1501  {
1502  libMesh::out << "[" << this->processor_id() << "] boundary_node_list.size()="
1503  << bc_tuples.size() << std::endl;
1504  libMesh::out << "[" << this->processor_id() << "] (boundary_node_id, boundary_id) = ";
1505  for (const auto & t : bc_tuples)
1506  libMesh::out << "(" << std::get<0>(t) << ", " << std::get<1>(t) << ") ";
1507  libMesh::out << std::endl;
1508  }
1509 
1510  // Now get the global information. In this case, we only want to count boundary
1511  // information for nodes *owned* by this processor, so there are no duplicates.
1512 
1513  // Make sure we don't have any left over information
1514  this->num_global_node_counts.clear();
1515  this->num_global_node_counts.resize(this->global_nodeset_ids.size());
1516 
1517  // Unfortunately, we can't just count up all occurrences of a given id,
1518  // that would give us duplicate entries when we do the parallel summation.
1519  // So instead, only count entries for nodes owned by this processor.
1520  // Start by getting rid of all non-local node entries from the vectors.
1521  std::vector<Tuple>::iterator
1522  it = bc_tuples.begin(),
1523  new_end = bc_tuples.end();
1524 
1525  while (it != new_end)
1526  {
1527  if (pmesh.node_ptr(std::get<0>(*it))->processor_id() != this->processor_id())
1528  {
1529  // Back up the new end iterators to prepare for swap
1530  --new_end;
1531 
1532  // Swap places, the non-local node will now be "past-the-end"
1533  std::swap(*it, *new_end);
1534  }
1535  else // node is local, go to next
1536  ++it;
1537  }
1538 
1539  // Erase from "new" end to old end.
1540  bc_tuples.erase(new_end, bc_tuples.end());
1541 
1542  // Now we can do the local count for each ID...
1543  for (auto i : index_range(global_nodeset_ids))
1544  {
1545  int id = this->global_nodeset_ids[i];
1546  this->num_global_node_counts[i] =
1547  cast_int<int>(std::count_if(bc_tuples.begin(),
1548  bc_tuples.end(),
1549  [id](const Tuple & t)->bool { return std::get<1>(t) == id; }));
1550  }
1551 
1552  // And finally we can sum them up
1553  this->comm().sum(this->num_global_node_counts);
1554 
1555  if (verbose)
1556  {
1557  libMesh::out << "[" << this->processor_id() << "] num_global_node_counts = ";
1558  for (const auto & cnt : num_global_node_counts)
1559  libMesh::out << cnt << ", ";
1560  libMesh::out << std::endl;
1561  }
1562 }
1563 
1564 
1565 
1566 
1568 {
1569  // 1.) Loop over active local elements, build up set of subdomain IDs.
1570  std::set<subdomain_id_type> global_subdomain_ids;
1571 
1572  // This map keeps track of the number of elements in each subdomain over all processors
1573  std::map<subdomain_id_type, unsigned> global_subdomain_counts;
1574 
1575  for (const auto & elem : pmesh.active_local_element_ptr_range())
1576  {
1577  subdomain_id_type cur_subdomain = elem->subdomain_id();
1578 
1579  /*
1580  // We can't have a zero subdomain ID in Exodus (for some reason?)
1581  // so map zero subdomains to a max value...
1582  if (cur_subdomain == 0)
1583  cur_subdomain = std::numeric_limits<subdomain_id_type>::max();
1584  */
1585 
1586  global_subdomain_ids.insert(cur_subdomain);
1587 
1588  // Increment the count of elements in this subdomain
1589  global_subdomain_counts[cur_subdomain]++;
1590  }
1591 
1592  // We're next going to this->comm().sum the subdomain counts, so save the local counts
1593  this->local_subdomain_counts = global_subdomain_counts;
1594 
1595  {
1596  // 2.) Copy local subdomain IDs into a vector for communication
1597  std::vector<subdomain_id_type> global_subdomain_ids_vector(global_subdomain_ids.begin(),
1598  global_subdomain_ids.end());
1599 
1600  // 3.) Gather them into an enlarged vector
1601  this->comm().allgather(global_subdomain_ids_vector);
1602 
1603  // 4.) Insert any new IDs into the set (any duplicates will be dropped)
1604  global_subdomain_ids.insert(global_subdomain_ids_vector.begin(),
1605  global_subdomain_ids_vector.end());
1606  }
1607 
1608  // 5.) Now global_subdomain_ids actually contains a global list of all subdomain IDs
1609  this->num_elem_blks_global =
1610  cast_int<int>(global_subdomain_ids.size());
1611 
1612  // Print the number of elements found locally in each subdomain
1613  if (verbose)
1614  {
1615  libMesh::out << "[" << this->processor_id() << "] ";
1616  for (const auto & pr : global_subdomain_counts)
1617  {
1618  libMesh::out << "ID: "
1619  << static_cast<unsigned>(pr.first)
1620  << ", Count: " << pr.second << ", ";
1621  }
1622  libMesh::out << std::endl;
1623  }
1624 
1625  // 6.) this->comm().sum up the number of elements in each block. We know the global
1626  // subdomain IDs, so pack them into a vector one by one. Use a vector of int since
1627  // that is what Nemesis wants
1628  this->global_elem_blk_cnts.resize(global_subdomain_ids.size());
1629 
1630  unsigned cnt=0;
1631  // Find the entry in the local map, note: if not found, will be created with 0 default value, which is OK...
1632  for (const auto & id : global_subdomain_ids)
1633  this->global_elem_blk_cnts[cnt++] = global_subdomain_counts[id];
1634 
1635  // Sum up subdomain counts from all processors
1636  this->comm().sum(this->global_elem_blk_cnts);
1637 
1638  if (verbose)
1639  {
1640  libMesh::out << "[" << this->processor_id() << "] global_elem_blk_cnts = ";
1641  for (const auto & bc : this->global_elem_blk_cnts)
1642  libMesh::out << bc << ", ";
1643  libMesh::out << std::endl;
1644  }
1645 
1646  // 7.) Create a vector<int> from the global_subdomain_ids set, for passing to Nemesis
1647  this->global_elem_blk_ids.clear();
1648  this->global_elem_blk_ids.insert(this->global_elem_blk_ids.end(), // pos
1649  global_subdomain_ids.begin(),
1650  global_subdomain_ids.end());
1651 
1652  if (verbose)
1653  {
1654  libMesh::out << "[" << this->processor_id() << "] global_elem_blk_ids = ";
1655  for (const auto & id : this->global_elem_blk_ids)
1656  libMesh::out << id << ", ";
1657  libMesh::out << std::endl;
1658  }
1659 
1660 
1661  // 8.) We will call put_eb_info_global later, it must be called after this->put_init_global().
1662 }
1663 
1664 
1665 
1666 
1668 {
1669  // If we don't have any local subdomains, it had better be because
1670  // we don't have any local elements
1671 #ifdef DEBUG
1672  if (local_subdomain_counts.empty())
1673  {
1675  pmesh.active_local_elements_end());
1677  }
1678 #endif
1679 
1680  // Elements have to be numbered contiguously based on what block
1681  // number they are in. Therefore we have to do a bit of work to get
1682  // the block (ie subdomain) numbers first and store them off as
1683  // block_ids.
1684 
1685  // Make sure there is no leftover information in the subdomain_map, and reserve
1686  // enough space to store the elements we need.
1687  this->subdomain_map.clear();
1688  for (const auto & pr : local_subdomain_counts)
1689  {
1690  if (verbose)
1691  {
1692  libMesh::out << "[" << this->processor_id() << "] "
1693  << "local_subdomain_counts [" << static_cast<unsigned>(pr.first) << "]= "
1694  << pr.second
1695  << std::endl;
1696  }
1697 
1698  // pr.first is the subdomain ID, pr.second is the number of elements it contains
1699  this->subdomain_map[pr.first].reserve(pr.second);
1700  }
1701 
1702 
1703  // First loop over the elements to figure out which elements are in which subdomain
1704  for (const auto & elem : pmesh.active_local_element_ptr_range())
1705  {
1706  // Grab the nodes while we're here.
1707  for (auto n : elem->node_index_range())
1708  this->nodes_attached_to_local_elems.insert( elem->node_id(n) );
1709 
1710  subdomain_id_type cur_subdomain = elem->subdomain_id();
1711 
1712  this->subdomain_map[cur_subdomain].push_back(elem->id());
1713  }
1714 
1715  // Set num_nodes which is used by exodusII_io_helper
1716  this->num_nodes =
1717  cast_int<int>(this->nodes_attached_to_local_elems.size());
1718 
1719  // Now come up with a 1-based numbering for these nodes
1720  this->exodus_node_num_to_libmesh.clear(); // Make sure it's empty
1721  this->exodus_node_num_to_libmesh.reserve(this->nodes_attached_to_local_elems.size());
1722 
1723  // Also make sure there's no leftover information in the map which goes the
1724  // other direction.
1725  this->libmesh_node_num_to_exodus.clear();
1726 
1727  // Set the map for nodes
1728  for (const auto & id : nodes_attached_to_local_elems)
1729  {
1730  // I.e. given exodus_node_id,
1731  // exodus_node_num_to_libmesh[ exodus_node_id ] returns the libmesh ID for that node.
1732  // Note that even though most of Exodus is 1-based, this code will map an Exodus ID of
1733  // zero to some libmesh node ID. Is that a problem?
1734  this->exodus_node_num_to_libmesh.push_back(id);
1735 
1736  // Likewise, given libmesh_node_id,
1737  // libmesh_node_num_to_exodus[ libmesh_node_id ] returns the *Exodus* ID for that node.
1738  // Unlike the exodus_node_num_to_libmesh vector above, this one is a std::map
1739  this->libmesh_node_num_to_exodus[id] =
1740  cast_int<int>(this->exodus_node_num_to_libmesh.size()); // should never be zero...
1741  }
1742 
1743  // Now we're going to loop over the subdomain map and build a few things right
1744  // now that we'll use later.
1745 
1746  // First make sure our data structures don't have any leftover data...
1747  this->exodus_elem_num_to_libmesh.clear();
1748  this->block_ids.clear();
1749  this->libmesh_elem_num_to_exodus.clear();
1750 
1751  // Now loop over each subdomain and get a unique numbering for the elements
1752  for (auto & pr : subdomain_map)
1753  {
1754  block_ids.push_back(pr.first);
1755 
1756  // Vector of element IDs for this subdomain
1757  std::vector<dof_id_type> & elem_ids_this_subdomain = pr.second;
1758 
1759  // The code below assumes this subdomain block is not empty, make sure that's the case!
1760  if (elem_ids_this_subdomain.size() == 0)
1761  libmesh_error_msg("Error, no element IDs found in subdomain " << pr.first);
1762 
1763  // Use the first element in this block to get representative information.
1764  // Note that Exodus assumes all elements in a block are of the same type!
1765  // We are using that same assumption here!
1766  const auto & conv = get_conversion
1767  (pmesh.elem_ref(elem_ids_this_subdomain[0]).type());
1768  this->num_nodes_per_elem =
1769  pmesh.elem_ref(elem_ids_this_subdomain[0]).n_nodes();
1770 
1771  // Get a reference to the connectivity vector for this subdomain. This vector
1772  // is most likely empty, we are going to fill it up now.
1773  std::vector<int> & current_block_connectivity = this->block_id_to_elem_connectivity[pr.first];
1774 
1775  // Just in case it's not already empty...
1776  current_block_connectivity.clear();
1777  current_block_connectivity.resize(elem_ids_this_subdomain.size() * this->num_nodes_per_elem);
1778 
1779  for (auto i : index_range(elem_ids_this_subdomain))
1780  {
1781  auto elem_id = elem_ids_this_subdomain[i];
1782 
1783  // Set the number map for elements
1784  this->exodus_elem_num_to_libmesh.push_back(elem_id);
1785  this->libmesh_elem_num_to_exodus[elem_id] =
1786  cast_int<int>(this->exodus_elem_num_to_libmesh.size());
1787 
1788  const Elem & elem = pmesh.elem_ref(elem_id);
1789 
1790  // Exodus/Nemesis want every block to have the same element type
1791  // libmesh_assert_equal_to (elem->type(), conv.libmesh_elem_type());
1792 
1793  // But we can get away with writing e.g. HEX8 and INFHEX8 in
1794  // the same block...
1795  libmesh_assert_equal_to (elem.n_nodes(), Elem::build(conv.libmesh_elem_type(), nullptr)->n_nodes());
1796 
1797  for (auto j : IntRange<int>(0, this->num_nodes_per_elem))
1798  {
1799  const unsigned int connect_index = (i*this->num_nodes_per_elem)+j;
1800  const unsigned int elem_node_index = conv.get_node_map(j);
1801 
1802  current_block_connectivity[connect_index] =
1803  libmesh_map_find(libmesh_node_num_to_exodus,
1804  elem.node_id(elem_node_index));
1805  }
1806  } // End loop over elems in this subdomain
1807  } // end loop over subdomain_map
1808 }
1809 
1810 
1811 
1812 
1813 
1815 {
1816  // The set which will eventually contain the IDs of "border nodes". These are nodes
1817  // that lie on the boundary between one or more processors.
1818  //std::set<unsigned> border_node_ids;
1819 
1820  // map from processor ID to set of nodes which elements from this processor "touch",
1821  // that is,
1822  // proc_nodes_touched[p] = (set all node IDs found in elements owned by processor p)
1823  std::map<unsigned, std::set<unsigned>> proc_nodes_touched;
1824 
1825 
1826  // We are going to create a lot of intermediate data structures here, so make sure
1827  // as many as possible all cleaned up by creating scope!
1828  {
1829  // Loop over active (not just active local) elements, make sets of node IDs for each
1830  // processor which has an element that "touches" a node.
1831  for (const auto & elem : pmesh.active_element_ptr_range())
1832  {
1833  // Get reference to the set for this processor. If it does not exist
1834  // it will be created.
1835  std::set<unsigned> & set_p = proc_nodes_touched[ elem->processor_id() ];
1836 
1837  // Insert all nodes touched by this element into the set
1838  for (auto node : elem->node_index_range())
1839  set_p.insert(elem->node_id(node));
1840  }
1841 
1842  // The number of node communication maps is the number of other processors
1843  // with which we share nodes. (I think.) This is just the size of the map we just
1844  // created, minus 1 unless this processor has no nodes of its own.
1845  this->num_node_cmaps =
1846  cast_int<int>(proc_nodes_touched.size() -
1847  proc_nodes_touched.count(this->processor_id()));
1848 
1849  // We can't be connecting to more processors than exist outside
1850  // ourselves
1851  libmesh_assert_less (this->num_node_cmaps, this->n_processors());
1852 
1853  if (verbose)
1854  {
1855  libMesh::out << "[" << this->processor_id()
1856  << "] proc_nodes_touched contains "
1857  << proc_nodes_touched.size()
1858  << " sets of nodes."
1859  << std::endl;
1860 
1861  for (const auto & pr : proc_nodes_touched)
1862  libMesh::out << "[" << this->processor_id()
1863  << "] proc_nodes_touched[" << pr.first << "] has "
1864  << pr.second.size()
1865  << " entries."
1866  << std::endl;
1867  }
1868 
1869 
1870  // Loop over all the sets we just created and compute intersections with the
1871  // this processor's set. Obviously, don't intersect with ourself.
1872  for (auto & pr : proc_nodes_touched)
1873  {
1874  // Don't compute intersections with ourself
1875  if (pr.first == this->processor_id())
1876  continue;
1877 
1878  // Otherwise, compute intersection with other processor and ourself
1879  std::set<unsigned> & my_set = proc_nodes_touched[this->processor_id()];
1880  std::set<unsigned> & other_set = pr.second;
1881  std::set<unsigned> & result_set = this->proc_nodes_touched_intersections[pr.first]; // created if does not exist
1882 
1883  std::set_intersection(my_set.begin(), my_set.end(),
1884  other_set.begin(), other_set.end(),
1885  std::inserter(result_set, result_set.end()));
1886  }
1887 
1888  if (verbose)
1889  {
1890  for (const auto & pr : proc_nodes_touched_intersections)
1891  libMesh::out << "[" << this->processor_id()
1892  << "] this->proc_nodes_touched_intersections[" << pr.first << "] has "
1893  << pr.second.size()
1894  << " entries."
1895  << std::endl;
1896  }
1897 
1898  // Compute the set_union of all the preceding intersections. This will be the set of
1899  // border node IDs for this processor.
1900  for (auto & pr : proc_nodes_touched_intersections)
1901  {
1902  std::set<unsigned> & other_set = pr.second;
1903  std::set<unsigned> intermediate_result; // Don't think we can insert into one of the sets we're unioning...
1904 
1905  std::set_union(this->border_node_ids.begin(), this->border_node_ids.end(),
1906  other_set.begin(), other_set.end(),
1907  std::inserter(intermediate_result, intermediate_result.end()));
1908 
1909  // Swap our intermediate result into the final set
1910  this->border_node_ids.swap(intermediate_result);
1911  }
1912 
1913  libmesh_assert_less_equal
1914  (this->proc_nodes_touched_intersections.size(),
1915  std::size_t(this->num_node_cmaps));
1916 
1917  if (verbose)
1918  {
1919  libMesh::out << "[" << this->processor_id()
1920  << "] border_node_ids.size()=" << this->border_node_ids.size()
1921  << std::endl;
1922  }
1923  } // end scope for border node ID creation
1924 
1925  // Store the number of border node IDs to be written to Nemesis file
1926  this->num_border_nodes = cast_int<int>(this->border_node_ids.size());
1927 }
1928 
1929 
1930 
1931 
1932 
1934 {
1935  // Write the nodesets. In Nemesis, the idea is to "create space" for the global
1936  // set of boundary nodesets, but to only write node IDs which are local to the current
1937  // processor. This is what is done in Nemesis files created by the "loadbal" script.
1938 
1939  // Store a map of vectors for boundary node IDs on this processor.
1940  // Use a vector of int here so it can be passed directly to Exodus.
1941  std::map<boundary_id_type, std::vector<int>> local_node_boundary_id_lists;
1942 
1943  // FIXME: We should build this list only one time!! We already built it above, but we
1944  // did not have the libmesh to exodus node mapping at that time... for now we'll just
1945  // build it here again, hopefully it's small relative to the size of the entire mesh.
1946 
1947  // Build list of (node-id, bc-id) tuples.
1948  typedef std::tuple<dof_id_type, boundary_id_type> Tuple;
1949  std::vector<Tuple> bc_tuples = mesh.get_boundary_info().build_node_list();
1950 
1951  if (verbose)
1952  {
1953  libMesh::out << "[" << this->processor_id() << "] boundary_node_list.size()="
1954  << bc_tuples.size() << std::endl;
1955  libMesh::out << "[" << this->processor_id() << "] (boundary_node_id, boundary_id) = ";
1956  for (const auto & t : bc_tuples)
1957  libMesh::out << "(" << std::get<0>(t) << ", " << std::get<1>(t) << ") ";
1958  libMesh::out << std::endl;
1959  }
1960 
1961  // For each node in the node list, add it to the vector of node IDs for that
1962  // set for the local processor. This will be used later when writing Exodus
1963  // nodesets.
1964  for (const auto & t : bc_tuples)
1965  {
1966  // Don't try to grab a reference to the vector unless the current node is attached
1967  // to a local element. Otherwise, another processor will be responsible for writing it in its nodeset.
1968  std::map<int, int>::iterator it = this->libmesh_node_num_to_exodus.find(std::get<0>(t));
1969 
1970  if (it != this->libmesh_node_num_to_exodus.end())
1971  {
1972  // Get reference to the vector where this node ID will be inserted. If it
1973  // doesn't yet exist, this will create it.
1974  std::vector<int> & current_id_set = local_node_boundary_id_lists[std::get<1>(t)];
1975 
1976  // Push back Exodus-mapped node ID for this set
1977  // TODO: reserve space in these vectors somehow.
1978  current_id_set.push_back( it->second );
1979  }
1980  }
1981 
1982  // See what we got
1983  if (verbose)
1984  {
1985  for (const auto & pr : local_node_boundary_id_lists)
1986  {
1987  libMesh::out << "[" << this->processor_id() << "] ID: " << pr.first << ", ";
1988 
1989  // Libmesh node ID (Exodus Node ID)
1990  for (const auto & id : pr.second)
1991  libMesh::out << id << ", ";
1992  libMesh::out << std::endl;
1993  }
1994  }
1995 
1996  // Loop over *global* nodeset IDs, call the Exodus API. Note that some nodesets may be empty
1997  // for a given processor.
1998  if (global_nodeset_ids.size() > 0)
1999  {
2000  NamesData names_table(global_nodeset_ids.size(), MAX_STR_LENGTH);
2001 
2002  for (const auto & nodeset_id : this->global_nodeset_ids)
2003  {
2004  const std::string & current_ns_name =
2005  mesh.get_boundary_info().get_nodeset_name
2006  (cast_int<boundary_id_type>(nodeset_id));
2007 
2008  // Store this name in a data structure that will be used to
2009  // write sideset names to file.
2010  names_table.push_back_entry(current_ns_name);
2011 
2012  if (verbose)
2013  {
2014  libMesh::out << "[" << this->processor_id()
2015  << "] Writing out Exodus nodeset info for ID: " << nodeset_id
2016  << ", Name: " << current_ns_name
2017  << std::endl;
2018  }
2019 
2020  // Convert current global_nodeset_id into an exodus ID, which can't be zero...
2021  int exodus_id = nodeset_id;
2022 
2023  /*
2024  // Exodus can't handle zero nodeset IDs (?) Use max short here since
2025  // when libmesh reads it back in, it will want to store it as a short...
2026  if (exodus_id==0)
2027  exodus_id = std::numeric_limits<short>::max();
2028  */
2029 
2030  // Try to find this boundary ID in the local list we created
2031  auto it =
2032  local_node_boundary_id_lists.find (cast_int<boundary_id_type>(nodeset_id));
2033 
2034  // No nodes found for this boundary ID on this processor
2035  if (it == local_node_boundary_id_lists.end())
2036  {
2037  if (verbose)
2038  libMesh::out << "[" << this->processor_id()
2039  << "] No nodeset data for ID: " << nodeset_id
2040  << " on this processor." << std::endl;
2041 
2042  // Call the Exodus interface to write the parameters of this node set
2043  this->ex_err = exII::ex_put_node_set_param(this->ex_id,
2044  exodus_id,
2045  0, /* No nodes for this ID */
2046  0 /* No distribution factors */);
2047  EX_CHECK_ERR(this->ex_err, "Error writing nodeset parameters in Nemesis");
2048 
2049  }
2050  else // Boundary ID *was* found in list
2051  {
2052  // Get reference to the vector of node IDs
2053  std::vector<int> & current_nodeset_ids = it->second;
2054 
2055  // Call the Exodus interface to write the parameters of this node set
2056  this->ex_err = exII::ex_put_node_set_param(this->ex_id,
2057  exodus_id,
2058  current_nodeset_ids.size(),
2059  0 /* No distribution factors */);
2060 
2061  EX_CHECK_ERR(this->ex_err, "Error writing nodeset parameters in Nemesis");
2062 
2063  // Call Exodus interface to write the actual node IDs for this boundary ID
2064  this->ex_err = exII::ex_put_node_set(this->ex_id,
2065  exodus_id,
2066  current_nodeset_ids.data());
2067 
2068  EX_CHECK_ERR(this->ex_err, "Error writing nodesets in Nemesis");
2069 
2070  }
2071  } // end loop over global nodeset IDs
2072 
2073  // Write out the nodeset names
2074  ex_err = exII::ex_put_names(ex_id,
2075  exII::EX_NODE_SET,
2076  names_table.get_char_star_star());
2077  EX_CHECK_ERR(ex_err, "Error writing nodeset names");
2078  } // end for loop over global nodeset IDs
2079 }
2080 
2081 
2082 
2083 
2085 {
2086  // Write the sidesets. In Nemesis, the idea is to "create space" for the global
2087  // set of boundary sidesets, but to only write sideset IDs which are local to the current
2088  // processor. This is what is done in Nemesis files created by the "loadbal" script.
2089  // See also: ExodusII_IO_Helper::write_sidesets()...
2090 
2091 
2092  // Store a map of vectors for boundary side IDs on this processor.
2093  // Use a vector of int here so it can be passed directly to Exodus.
2094  std::map<boundary_id_type, std::vector<int>> local_elem_boundary_id_lists;
2095  std::map<boundary_id_type, std::vector<int>> local_elem_boundary_id_side_lists;
2096 
2097  // FIXME: We already built this list once, we should reuse that information!
2098  std::vector<std::tuple<dof_id_type, unsigned short int, boundary_id_type>> bndry_elem_side_id_list =
2099  mesh.get_boundary_info().build_side_list();
2100 
2101  // Integer looping, skipping non-local elements
2102  for (const auto & t : bndry_elem_side_id_list)
2103  {
2104  // Get pointer to current Elem
2105  const Elem * elem = mesh.elem_ptr(std::get<0>(t));
2106 
2107  std::vector<const Elem *> family;
2108 #ifdef LIBMESH_ENABLE_AMR
2109  // We need to build up active elements if AMR is enabled and add
2110  // them to the exodus sidesets instead of the potentially inactive "parent" elements
2111  // Technically we don't need to "reset" the tree since the vector was just created.
2112  elem->active_family_tree_by_side(family, std::get<1>(t), /*reset tree=*/false);
2113 #else
2114  // If AMR is not even enabled, just push back the element itself
2115  family.push_back( elem );
2116 #endif
2117 
2118  // Loop over all the elements in the family tree, store their converted IDs
2119  // and side IDs to the map's vectors. TODO: Somehow reserve enough space for these
2120  // push_back's...
2121  for (const auto & tree_elem : family)
2122  {
2123  const dof_id_type f_id = tree_elem->id();
2124  const Elem & f = mesh.elem_ref(f_id);
2125 
2126  // If element is local, process it
2127  if (f.processor_id() == this->processor_id())
2128  {
2129  const auto & conv = get_conversion(f.type());
2130 
2131  // Use the libmesh to exodus data structure map to get the proper sideset IDs
2132  // The data structure contains the "collapsed" contiguous ids.
2133  //
2134  // We know the parent element is local, but let's be absolutely sure that all the children have been
2135  // actually mapped to Exodus IDs before we blindly try to add them...
2136  local_elem_boundary_id_lists[ std::get<2>(t) ].push_back( libmesh_map_find(libmesh_elem_num_to_exodus, f_id) );
2137  local_elem_boundary_id_side_lists[ std::get<2>(t) ].push_back(conv.get_inverse_side_map( std::get<1>(t) ));
2138  }
2139  }
2140  }
2141 
2142 
2143  // Loop over *global* sideset IDs, call the Exodus API. Note that some sidesets may be empty
2144  // for a given processor.
2145  if (global_sideset_ids.size() > 0)
2146  {
2147  NamesData names_table(global_sideset_ids.size(), MAX_STR_LENGTH);
2148 
2149  for (const auto & exodus_id : this->global_sideset_ids)
2150  {
2151  const std::string & current_ss_name =
2152  mesh.get_boundary_info().get_sideset_name
2153  (cast_int<boundary_id_type>(exodus_id));
2154 
2155  // Store this name in a data structure that will be used to
2156  // write sideset names to file.
2157  names_table.push_back_entry(current_ss_name);
2158 
2159  if (verbose)
2160  {
2161  libMesh::out << "[" << this->processor_id()
2162  << "] Writing out Exodus sideset info for ID: " << exodus_id
2163  << ", Name: " << current_ss_name
2164  << std::endl;
2165  }
2166 
2167  // Try to find this boundary ID in the local list we created
2168  auto it =
2169  local_elem_boundary_id_lists.find (cast_int<boundary_id_type>(exodus_id));
2170 
2171  // No sides found for this boundary ID on this processor
2172  if (it == local_elem_boundary_id_lists.end())
2173  {
2174  if (verbose)
2175  libMesh::out << "[" << this->processor_id()
2176  << "] No sideset data for ID: " << exodus_id
2177  << " on this processor." << std::endl;
2178 
2179  // Call the Exodus interface to write the parameters of this side set
2180  this->ex_err = exII::ex_put_side_set_param(this->ex_id,
2181  exodus_id,
2182  0, /* No sides for this ID */
2183  0 /* No distribution factors */);
2184  EX_CHECK_ERR(this->ex_err, "Error writing sideset parameters in Nemesis");
2185 
2186  }
2187  else // Boundary ID *was* found in list
2188  {
2189  // Get reference to the vector of elem IDs
2190  std::vector<int> & current_sideset_elem_ids = it->second;
2191 
2192  // Get reference to the vector of side IDs
2193  std::vector<int> & current_sideset_side_ids =
2194  libmesh_map_find(local_elem_boundary_id_side_lists,
2195  cast_int<boundary_id_type>(exodus_id));
2196 
2197  // Call the Exodus interface to write the parameters of this side set
2198  this->ex_err = exII::ex_put_side_set_param(this->ex_id,
2199  exodus_id,
2200  current_sideset_elem_ids.size(),
2201  0 /* No distribution factors */);
2202 
2203  EX_CHECK_ERR(this->ex_err, "Error writing sideset parameters in Nemesis");
2204 
2205  // Call Exodus interface to write the actual side IDs for this boundary ID
2206  this->ex_err = exII::ex_put_side_set(this->ex_id,
2207  exodus_id,
2208  current_sideset_elem_ids.data(),
2209  current_sideset_side_ids.data());
2210 
2211  EX_CHECK_ERR(this->ex_err, "Error writing sidesets in Nemesis");
2212  }
2213  } // end for loop over global sideset IDs
2214 
2215  // Write sideset names to file. Some of these may be blank strings
2216  // if the current processor didn't have all the sideset names for
2217  // any reason...
2218  ex_err = exII::ex_put_names(this->ex_id,
2219  exII::EX_SIDE_SET,
2220  names_table.get_char_star_star());
2221  EX_CHECK_ERR(ex_err, "Error writing sideset names");
2222 
2223  } // end if (global_sideset_ids.size() > 0)
2224 }
2225 
2226 
2227 
2228 void Nemesis_IO_Helper::write_nodal_coordinates(const MeshBase & mesh, bool /*use_discontinuous*/)
2229 {
2230  auto local_num_nodes = this->exodus_node_num_to_libmesh.size();
2231 
2232  x.resize(local_num_nodes);
2233  y.resize(local_num_nodes);
2234  z.resize(local_num_nodes);
2235 
2236  // Just loop over our list outputting the nodes the way we built the map
2237  for (auto i : IntRange<int>(0, local_num_nodes))
2238  {
2239  const Point & pt = mesh.point(this->exodus_node_num_to_libmesh[i]);
2240  x[i]=pt(0);
2241  y[i]=pt(1);
2242  z[i]=pt(2);
2243  }
2244 
2245  if (local_num_nodes)
2246  {
2247  if (_single_precision)
2248  {
2249  std::vector<float>
2250  x_single(x.begin(), x.end()),
2251  y_single(y.begin(), y.end()),
2252  z_single(z.begin(), z.end());
2253 
2254  ex_err = exII::ex_put_coord(ex_id,
2255  x_single.data(),
2256  y_single.data(),
2257  z_single.data());
2258  }
2259  else
2260  {
2261  // Call Exodus API to write nodal coordinates...
2262  ex_err = exII::ex_put_coord(ex_id,
2263  x.data(),
2264  y.data(),
2265  z.data());
2266  }
2267  EX_CHECK_ERR(ex_err, "Error writing node coordinates");
2268 
2269  // And write the nodal map we created for them
2270  ex_err = exII::ex_put_node_num_map(ex_id, this->exodus_node_num_to_libmesh.data());
2271  EX_CHECK_ERR(ex_err, "Error writing node num map");
2272  }
2273  else // Does the Exodus API want us to write empty nodal coordinates?
2274  {
2275  ex_err = exII::ex_put_coord(ex_id, nullptr, nullptr, nullptr);
2276  EX_CHECK_ERR(ex_err, "Error writing empty node coordinates");
2277 
2278  ex_err = exII::ex_put_node_num_map(ex_id, nullptr);
2279  EX_CHECK_ERR(ex_err, "Error writing empty node num map");
2280  }
2281 }
2282 
2283 
2284 
2285 
2286 
2287 void Nemesis_IO_Helper::write_elements(const MeshBase & mesh, bool /*use_discontinuous*/)
2288 {
2289  // Only write elements if there are elements blocks available.
2290  if (this->num_elem_blks_global > 0)
2291  {
2292  // Data structure to store element block names that will be used to
2293  // write the element block names to file.
2294  NamesData names_table(this->num_elem_blks_global, MAX_STR_LENGTH);
2295 
2296  // Loop over all blocks, even if we don't have elements in each block.
2297  // If we don't have elements we need to write out a 0 for that block...
2298  for (auto i : IntRange<int>(0, this->num_elem_blks_global))
2299  {
2300  // Even if there are no elements for this block on the current
2301  // processor, we still want to write its name to file, if
2302  // possible. MeshBase::subdomain_name() will just return an
2303  // empty string if there is no name associated with the current
2304  // block.
2305  names_table.push_back_entry
2306  (mesh.subdomain_name(cast_int<subdomain_id_type>(this->global_elem_blk_ids[i])));
2307 
2308  // Search for the current global block ID in the map
2309  std::map<int, std::vector<int>>::iterator it =
2310  this->block_id_to_elem_connectivity.find( this->global_elem_blk_ids[i] );
2311 
2312  // If not found, write a zero to file....
2313  if (it == this->block_id_to_elem_connectivity.end())
2314  {
2315  this->ex_err = exII::ex_put_elem_block(this->ex_id,
2316  this->global_elem_blk_ids[i],
2317  "Empty",
2318  0, /* n. elements in this block */
2319  0, /* n. nodes per element */
2320  0); /* number of attributes per element */
2321 
2322  EX_CHECK_ERR(this->ex_err, "Error writing element block from Nemesis.");
2323  }
2324 
2325  // Otherwise, write the actual block information and connectivity to file
2326  else
2327  {
2328  subdomain_id_type block =
2329  cast_int<subdomain_id_type>(it->first);
2330  std::vector<int> & this_block_connectivity = it->second;
2331  std::vector<dof_id_type> & elements_in_this_block = subdomain_map[block];
2332 
2333  // Use the first element in this block to get representative information.
2334  // Note that Exodus assumes all elements in a block are of the same type!
2335  // We are using that same assumption here!
2336  const auto & conv =
2337  get_conversion(mesh.elem_ref(elements_in_this_block[0]).type());
2338 
2339  this->num_nodes_per_elem =
2340  mesh.elem_ref(elements_in_this_block[0]).n_nodes();
2341 
2342  ex_err = exII::ex_put_elem_block(ex_id,
2343  block,
2344  conv.exodus_elem_type().c_str(),
2345  elements_in_this_block.size(),
2347  0);
2348  EX_CHECK_ERR(ex_err, "Error writing element block from Nemesis.");
2349 
2350  ex_err = exII::ex_put_elem_conn(ex_id,
2351  block,
2352  this_block_connectivity.data());
2353  EX_CHECK_ERR(ex_err, "Error writing element connectivities from Nemesis.");
2354  }
2355  } // end loop over global block IDs
2356 
2357  // Only call this once, not in the loop above!
2358  ex_err = exII::ex_put_elem_num_map(ex_id,
2359  exodus_elem_num_to_libmesh.empty() ? nullptr : exodus_elem_num_to_libmesh.data());
2360  EX_CHECK_ERR(ex_err, "Error writing element map");
2361 
2362  // Write the element block names to file.
2363  ex_err = exII::ex_put_names(ex_id, exII::EX_ELEM_BLOCK, names_table.get_char_star_star());
2364  EX_CHECK_ERR(ex_err, "Error writing element block names");
2365  } // end if (this->num_elem_blks_global > 0)
2366 }
2367 
2368 
2369 
2370 
2371 
2372 void Nemesis_IO_Helper::write_nodal_solution(const std::vector<Number> & values,
2373  const std::vector<std::string> & names,
2374  int timestep)
2375 {
2376  int num_vars = cast_int<int>(names.size());
2377  //int num_values = values.size(); // Not used?
2378 
2379  for (int c=0; c<num_vars; c++)
2380  {
2381 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
2382  std::vector<Real> real_parts(num_nodes);
2383  std::vector<Real> imag_parts(num_nodes);
2384  std::vector<Real> magnitudes(num_nodes);
2385 
2386  for (int i=0; i<num_nodes; ++i)
2387  {
2388  Number value = values[this->exodus_node_num_to_libmesh[i]*num_vars + c];
2389  real_parts[i] = value.real();
2390  imag_parts[i] = value.imag();
2391  magnitudes[i] = std::abs(value);
2392  }
2393  write_nodal_values(3*c+1,real_parts,timestep);
2394  write_nodal_values(3*c+2,imag_parts,timestep);
2395  write_nodal_values(3*c+3,magnitudes,timestep);
2396 #else
2397  std::vector<Number> cur_soln(num_nodes);
2398 
2399  // Copy out this variable's solution
2400  for (int i=0; i<num_nodes; i++)
2401  cur_soln[i] = values[this->exodus_node_num_to_libmesh[i]*num_vars + c];
2402 
2403  write_nodal_values(c+1,cur_soln,timestep);
2404 #endif
2405  }
2406 }
2407 
2408 
2409 
2411  const std::vector<std::string> & names,
2412  int timestep,
2413  const std::vector<std::string> & output_names)
2414 {
2415  int num_vars = cast_int<int>(names.size());
2416 
2417  for (int c=0; c<num_vars; c++)
2418  {
2419  // Find the position of names[c] in the output_names vector, if it exists.
2420  auto pos = std::find(output_names.begin(), output_names.end(), names[c]);
2421 
2422  // Skip names[c] if it's not supposed to be output.
2423  if (pos == output_names.end())
2424  continue;
2425 
2426  // Compute the (zero-based) index which determines which
2427  // variable this will be as far as Nemesis is concerned. This
2428  // will be used below in the write_nodal_values() call.
2429  int variable_name_position =
2430  cast_int<int>(std::distance(output_names.begin(), pos));
2431 
2432  // Fill up a std::vector with the dofs for the current variable
2433  std::vector<numeric_index_type> required_indices(num_nodes);
2434 
2435  for (int i=0; i<num_nodes; i++)
2436  required_indices[i] = static_cast<dof_id_type>(this->exodus_node_num_to_libmesh[i]) * num_vars + c;
2437 
2438  // Get the dof values required to write just our local part of
2439  // the solution vector.
2440  std::vector<Number> local_soln;
2441  parallel_soln.localize(local_soln, required_indices);
2442 
2443 #ifndef LIBMESH_USE_COMPLEX_NUMBERS
2444  // Call the ExodusII_IO_Helper function to write the data.
2445  write_nodal_values(variable_name_position + 1, local_soln, timestep);
2446 #else
2447  // We have the local (complex) values. Now extract the real,
2448  // imaginary, and magnitude values from them.
2449  std::vector<Real> real_parts(num_nodes);
2450  std::vector<Real> imag_parts(num_nodes);
2451  std::vector<Real> magnitudes(num_nodes);
2452 
2453  for (int i=0; i<num_nodes; ++i)
2454  {
2455  real_parts[i] = local_soln[i].real();
2456  imag_parts[i] = local_soln[i].imag();
2457  magnitudes[i] = std::abs(local_soln[i]);
2458  }
2459 
2460  // Write the real, imaginary, and magnitude values to file.
2461  write_nodal_values(3 * variable_name_position + 1, real_parts, timestep);
2462  write_nodal_values(3 * variable_name_position + 2, imag_parts, timestep);
2463  write_nodal_values(3 * variable_name_position + 3, magnitudes, timestep);
2464 #endif
2465  }
2466 }
2467 
2468 
2469 
2471  const std::vector<std::pair<unsigned int, unsigned int>> & var_nums,
2472  int timestep,
2473  const std::vector<std::string> & output_names)
2474 {
2475  const MeshBase & mesh = es.get_mesh();
2476 
2477  for (auto & var_num : var_nums)
2478  {
2479  const unsigned int sys_num = var_num.first;
2480  const unsigned int var = var_num.second;
2481 
2482  const System & sys = es.get_system(sys_num);
2483  const std::string & name = sys.variable_name(var);
2484 
2485  auto pos = std::find(output_names.begin(), output_names.end(), name);
2486 
2487  // Skip this name if it's not supposed to be output.
2488  if (pos == output_names.end())
2489  continue;
2490 
2491  // Compute the (zero-based) index which determines which
2492  // variable this will be as far as Nemesis is concerned. This
2493  // will be used below in the write_nodal_values() call.
2494  int variable_name_position =
2495  cast_int<int>(std::distance(output_names.begin(), pos));
2496 
2497  // Fill up a std::vector with the dofs for the current variable
2498  std::vector<numeric_index_type> required_indices(num_nodes);
2499 
2500  const FEType type = sys.variable_type(var);
2501  if (type.family == SCALAR)
2502  {
2503  std::vector<numeric_index_type> scalar_indices;
2504  sys.get_dof_map().SCALAR_dof_indices(scalar_indices, var);
2505  for (int i=0; i<num_nodes; i++)
2506  required_indices[i] = scalar_indices[0];
2507  }
2508  else
2509  for (int i=0; i<num_nodes; i++)
2510  {
2511  const Node & node = mesh.node_ref(this->exodus_node_num_to_libmesh[i]);
2512  required_indices[i] = node.dof_number(sys_num, var, 0);
2513  }
2514 
2515  // Get the dof values required to write just our local part of
2516  // the solution vector.
2517  std::vector<Number> local_soln;
2518  sys.current_local_solution->get(required_indices, local_soln);
2519 
2520 #ifndef LIBMESH_USE_COMPLEX_NUMBERS
2521  // Call the ExodusII_IO_Helper function to write the data.
2522  write_nodal_values(variable_name_position + 1, local_soln, timestep);
2523 #else
2524  // We have the local (complex) values. Now extract the real,
2525  // imaginary, and magnitude values from them.
2526  std::vector<Real> real_parts(num_nodes);
2527  std::vector<Real> imag_parts(num_nodes);
2528  std::vector<Real> magnitudes(num_nodes);
2529 
2530  for (int i=0; i<num_nodes; ++i)
2531  {
2532  real_parts[i] = local_soln[i].real();
2533  imag_parts[i] = local_soln[i].imag();
2534  magnitudes[i] = std::abs(local_soln[i]);
2535  }
2536 
2537  // Write the real, imaginary, and magnitude values to file.
2538  write_nodal_values(3 * variable_name_position + 1, real_parts, timestep);
2539  write_nodal_values(3 * variable_name_position + 2, imag_parts, timestep);
2540  write_nodal_values(3 * variable_name_position + 3, magnitudes, timestep);
2541 #endif
2542  }
2543 }
2544 
2545 
2546 
2547 void
2549  const std::vector<std::set<subdomain_id_type>> & vars_active_subdomains)
2550 {
2551  // Quick return if there are no element variables to write
2552  if (names.size() == 0)
2553  return;
2554 
2555  // Quick return if we have already called this function
2557  return;
2558 
2559  // Be sure that variables in the file match what we are asking for
2560  if (num_elem_vars > 0)
2561  {
2562  this->check_existing_vars(ELEMENTAL, names, this->elem_var_names);
2563  return;
2564  }
2565 
2566  // Set the flag so we can skip this stuff on subsequent calls to
2567  // initialize_element_variables()
2568  _elem_vars_initialized = true;
2569 
2570  this->write_var_names(ELEMENTAL, names);
2571 
2572  // Create a truth table from global_elem_blk_ids and the information
2573  // in vars_active_subdomains. Create a truth table of
2574  // size global_elem_blk_ids.size() * names.size().
2575  std::vector<int> truth_tab(global_elem_blk_ids.size() * names.size());
2576  for (auto blk : index_range(global_elem_blk_ids))
2577  for (auto var : index_range(names))
2578  if (vars_active_subdomains[var].empty() ||
2579  vars_active_subdomains[var].count(cast_int<subdomain_id_type>(global_elem_blk_ids[blk])))
2580  truth_tab[names.size() * blk + var] = 1;
2581 
2582  // Write truth table to file.
2583  if (truth_tab.size())
2584  {
2585  ex_err = exII::ex_put_elem_var_tab(ex_id,
2586  cast_int<int>(global_elem_blk_ids.size()),
2587  cast_int<int>(names.size()),
2588  truth_tab.data());
2589  EX_CHECK_ERR(ex_err, "Error writing element truth table.");
2590  }
2591 }
2592 
2593 
2594 
2595 void
2597  const EquationSystems & es,
2598  const std::vector<std::pair<unsigned int, unsigned int>> &var_nums,
2599  int timestep,
2600  const std::vector<std::set<subdomain_id_type>> & vars_active_subdomains)
2601 {
2602  // For each variable in names,
2603  // For each subdomain in subdomain_map,
2604  // If this (subdomain, variable) combination is active
2605  // Extract our values into local_soln (localize is a collective)
2606  // Write local_soln to file
2607  for (auto v : index_range(var_nums))
2608  {
2609  // Get list of active subdomains for variable v
2610  const auto & active_subdomains = vars_active_subdomains[v];
2611 
2612  for (const int sbd_id_int : global_elem_blk_ids)
2613  {
2614  const subdomain_id_type sbd_id =
2615  cast_int<subdomain_id_type>(sbd_id_int);
2616  auto it = subdomain_map.find(sbd_id);
2617  const std::vector<dof_id_type> empty_vec;
2618  const std::vector<dof_id_type> & elem_ids =
2619  (it == subdomain_map.end()) ? empty_vec : it->second;
2620 
2621  // Possibly skip this (variable, subdomain) combination
2622  if (active_subdomains.empty() || active_subdomains.count(sbd_id))
2623  {
2624  std::vector<numeric_index_type> required_indices;
2625  required_indices.reserve(elem_ids.size());
2626 
2627  const unsigned int sys_num = var_nums[v].first;
2628 
2629  for (const auto & id : elem_ids)
2630  required_indices.push_back
2631  (mesh.elem_ref(id).dof_number
2632  (sys_num, var_nums[v].second, 0));
2633 
2634  std::vector<Number> local_soln;
2635  es.get_system(sys_num).current_local_solution->get
2636  (required_indices, local_soln);
2637 
2638  // It's possible that there's nothing for us to write:
2639  // we may not be responsible for any elements on the
2640  // current subdomain. We did still have to participate
2641  // in the localize() call above, however, since it is a
2642  // collective.
2643  if (local_soln.size())
2644  {
2645  ex_err = exII::ex_put_elem_var(ex_id,
2646  timestep,
2647  static_cast<int>(v+1),
2648  static_cast<int>(sbd_id),
2649  static_cast<int>(local_soln.size()),
2650  local_soln.data());
2651  EX_CHECK_ERR(ex_err, "Error writing element values.");
2652  }
2653  }
2654  }
2655  } // end loop over vars
2656 
2657  ex_err = exII::ex_update(ex_id);
2658  EX_CHECK_ERR(ex_err, "Error flushing buffers to file.");
2659 }
2660 
2661 
2662 
2663 std::string Nemesis_IO_Helper::construct_nemesis_filename(const std::string & base_filename)
2664 {
2665  // Build a filename for this processor. This code is cut-n-pasted from the read function
2666  // and should probably be put into a separate function...
2667  std::ostringstream file_oss;
2668 
2669  // We have to be a little careful here: Nemesis left pads its file
2670  // numbers based on the number of processors, so for example on 10
2671  // processors, we'd have:
2672  // mesh.e.10.00
2673  // mesh.e.10.01
2674  // mesh.e.10.02
2675  // ...
2676  // mesh.e.10.09
2677 
2678  // And on 100 you'd have:
2679  // mesh.e.100.000
2680  // mesh.e.100.001
2681  // ...
2682  // mesh.e.128.099
2683 
2684  // Find the length of the highest processor ID
2685  file_oss << (this->n_processors());
2686  unsigned int field_width = cast_int<unsigned int>(file_oss.str().size());
2687 
2688  if (verbose)
2689  libMesh::out << "field_width=" << field_width << std::endl;
2690 
2691  file_oss.str(""); // reset the string stream
2692  file_oss << base_filename
2693  << '.' << this->n_processors()
2694  << '.' << std::setfill('0') << std::setw(field_width) << this->processor_id();
2695 
2696  // Return the resulting string
2697  return file_oss.str();
2698 }
2699 
2700 } // namespace libMesh
2701 
2702 #endif // #if defined(LIBMESH_HAVE_NEMESIS_API) && defined(LIBMESH_HAVE_EXODUS_API)
libMesh::ExodusII_IO_Helper::opened_for_writing
bool opened_for_writing
Definition: exodusII_io_helper.h:667
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::Nemesis_IO_Helper::num_node_sets_global
int num_node_sets_global
Definition: nemesis_io_helper.h:374
libMesh::Nemesis_IO_Helper::get_init_info
void get_init_info()
Definition: nemesis_io_helper.C:208
libMesh::Elem::active_family_tree_by_side
void active_family_tree_by_side(std::vector< const Elem * > &family, unsigned int side, bool reset=true) const
Same as the active_family_tree() member, but only adds elements which are next to side.
Definition: elem.C:1507
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::Nemesis_IO_Helper::num_internal_nodes
int num_internal_nodes
To be used with the Nemesis::ne_get_loadbal_param() routine.
Definition: nemesis_io_helper.h:453
libMesh::ExodusII_IO_Helper::get_conversion
const ExodusII_IO_Helper::Conversion & get_conversion(const ElemType type) const
Definition: exodusII_io_helper.C:383
libMesh::Nemesis_IO_Helper::get_cmap_params
void get_cmap_params()
Definition: nemesis_io_helper.C:325
libMesh::ExodusII_IO_Helper::_single_precision
bool _single_precision
Definition: exodusII_io_helper.h:740
libMesh::Nemesis_IO_Helper::num_side_sets_global
int num_side_sets_global
Definition: nemesis_io_helper.h:375
libMesh::ExodusII_IO_Helper::libmesh_node_num_to_exodus
std::map< int, int > libmesh_node_num_to_exodus
Definition: exodusII_io_helper.h:623
libMesh::Nemesis_IO_Helper::put_cmap_params
void put_cmap_params(std::vector< int > &node_cmap_ids, std::vector< int > &node_cmap_node_cnts, std::vector< int > &elem_cmap_ids, std::vector< int > &elem_cmap_elem_cnts)
Outputs initial information for communication maps.
Definition: nemesis_io_helper.C:576
libMesh::Nemesis_IO_Helper::node_cmap_node_cnts
std::vector< int > node_cmap_node_cnts
Definition: nemesis_io_helper.h:543
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::Nemesis_IO_Helper::elem_cmap_proc_ids
std::vector< std::vector< int > > elem_cmap_proc_ids
Definition: nemesis_io_helper.h:567
libMesh::ExodusII_IO_Helper::check_existing_vars
void check_existing_vars(ExodusVarType type, std::vector< std::string > &names, std::vector< std::string > &names_from_file)
When appending: during initialization, check that variable names in the file match those you attempt ...
Definition: exodusII_io_helper.C:2522
libMesh::Nemesis_IO_Helper::write_elements
virtual void write_elements(const MeshBase &mesh, bool use_discontinuous=false) override
This function is specialized to write the connectivity.
Definition: nemesis_io_helper.C:2287
libMesh::FEType::family
FEFamily family
The type of finite element.
Definition: fe_type.h:203
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
libMesh::Nemesis_IO_Helper::num_proc
int num_proc
The number of processors for which the NEMESIS I file was created.
Definition: nemesis_io_helper.h:381
libMesh::Nemesis_IO_Helper::compute_num_global_sidesets
void compute_num_global_sidesets(const MeshBase &pmesh)
This function uses global communication routines to determine the number of sidesets across the entir...
Definition: nemesis_io_helper.C:1363
libMesh::ExodusII_IO_Helper::y
std::vector< Real > y
Definition: exodusII_io_helper.h:605
libMesh::ExodusII_IO_Helper::write_var_names
void write_var_names(ExodusVarType type, const std::vector< std::string > &names)
Wraps calls to exII::ex_put_var_names() and exII::ex_put_var_param().
Definition: exodusII_io_helper.C:1408
libMesh::Nemesis_IO_Helper::elem_cmap_elem_cnts
std::vector< int > elem_cmap_elem_cnts
Definition: nemesis_io_helper.h:545
libMesh::Nemesis_IO_Helper::elem_cmap_side_ids
std::vector< std::vector< int > > elem_cmap_side_ids
Definition: nemesis_io_helper.h:566
libMesh::Nemesis_IO_Helper::~Nemesis_IO_Helper
virtual ~Nemesis_IO_Helper()
Destructor.
Definition: nemesis_io_helper.C:73
libMesh::Elem::n_nodes
virtual unsigned int n_nodes() const =0
libMesh::Nemesis_IO_Helper::compute_element_maps
void compute_element_maps()
This function computes element maps (really just packs vectors) which map the elements to internal an...
Definition: nemesis_io_helper.C:970
libMesh::Nemesis_IO_Helper::initialize_element_variables
virtual void initialize_element_variables(std::vector< std::string > names, const std::vector< std::set< subdomain_id_type >> &vars_active_subdomains) override
Override the Exodus Helper's implementation of this function so that it works correctly in parallel.
Definition: nemesis_io_helper.C:2548
libMesh::Nemesis_IO_Helper::num_elems_global
int num_elems_global
Definition: nemesis_io_helper.h:372
libMesh::Nemesis_IO_Helper::Nemesis_IO_Helper
Nemesis_IO_Helper(const ParallelObject &parent, bool verbose=false, bool single_precision=false)
Constructor.
Definition: nemesis_io_helper.C:48
libMesh::Nemesis_IO_Helper::global_sideset_ids
std::vector< int > global_sideset_ids
Containers for reading global sideset (boundary conditions) information.
Definition: nemesis_io_helper.h:404
libMesh::ExodusII_IO_Helper::block_ids
std::vector< int > block_ids
Definition: exodusII_io_helper.h:541
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::Nemesis_IO_Helper::put_elem_cmap
void put_elem_cmap(std::vector< std::vector< int >> &elem_cmap_elem_ids, std::vector< std::vector< int >> &elem_cmap_side_ids, std::vector< std::vector< int >> &elem_cmap_proc_ids)
Writes information about elemental communication map.
Definition: nemesis_io_helper.C:662
libMesh::MeshBase::active_local_elements_begin
virtual element_iterator active_local_elements_begin()=0
libMesh::Nemesis_IO_Helper::proc_nodes_touched_iterator
std::map< unsigned, std::set< unsigned > >::iterator proc_nodes_touched_iterator
Typedef for an iterator into the data structure above.
Definition: nemesis_io_helper.h:593
libMesh::Nemesis_IO_Helper::num_external_nodes
int num_external_nodes
The number of FEM nodes that reside on another processor but whose element partially resides on the c...
Definition: nemesis_io_helper.h:467
libMesh::Nemesis_IO_Helper::compute_border_node_ids
void compute_border_node_ids(const MeshBase &pmesh)
This function constructs the set of border node IDs present on the current mesh.
Definition: nemesis_io_helper.C:1814
libMesh::ExodusII_IO_Helper::libmesh_elem_num_to_exodus
std::map< int, int > libmesh_elem_num_to_exodus
Definition: exodusII_io_helper.h:618
libMesh::ExodusII_IO_Helper
This is the ExodusII_IO_Helper class.
Definition: exodusII_io_helper.h:80
libMesh::MeshBase::active_element_ptr_range
virtual SimpleRange< element_iterator > active_element_ptr_range()=0
libMesh::MeshBase::elem_ref
virtual const Elem & elem_ref(const dof_id_type i) const
Definition: mesh_base.h:521
libMesh::Nemesis_IO_Helper::create
virtual void create(std::string filename) override
This function is specialized from ExodusII_IO_Helper to create the nodal coordinates stored on the lo...
Definition: nemesis_io_helper.C:730
libMesh::System::variable_name
const std::string & variable_name(const unsigned int i) const
Definition: system.h:2203
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Nemesis_IO_Helper::construct_nemesis_filename
std::string construct_nemesis_filename(const std::string &base_filename)
Given base_filename, foo.e, constructs the Nemesis filename foo.e.X.Y, where X=n.
Definition: nemesis_io_helper.C:2663
libMesh::Nemesis_IO_Helper::compute_num_global_nodesets
void compute_num_global_nodesets(const MeshBase &pmesh)
This function uses global communication routines to determine the number of nodesets across the entir...
Definition: nemesis_io_helper.C:1456
libMesh::Nemesis_IO_Helper::put_elem_map
void put_elem_map(std::vector< int > &elem_mapi, std::vector< int > &elem_mapb)
Outputs IDs of internal and border elements.
Definition: nemesis_io_helper.C:683
libMesh::ExodusII_IO_Helper::num_nodes_per_elem
int num_nodes_per_elem
Definition: exodusII_io_helper.h:532
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::Nemesis_IO_Helper::get_ss_param_global
void get_ss_param_global()
Fills: global_sideset_ids, num_global_side_counts, num_global_side_df_counts Call after: get_init_glo...
Definition: nemesis_io_helper.C:111
end
IterBase * end
Also have a polymorphic pointer to the end object, this prevents iterating past the end.
Definition: variant_filter_iterator.h:343
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::Nemesis_IO_Helper::border_elem_ids
std::set< unsigned > border_elem_ids
A set of border elem IDs for this processor.
Definition: nemesis_io_helper.h:618
libMesh::BoundaryInfo::get_side_boundary_ids
const std::set< boundary_id_type > & get_side_boundary_ids() const
Definition: boundary_info.h:794
libMesh::Nemesis_IO_Helper::elem_mapb
std::vector< int > elem_mapb
Vector which stores border element IDs.
Definition: nemesis_io_helper.h:510
libMesh::Nemesis_IO_Helper::node_mapi
std::vector< int > node_mapi
Vector which stores internal node IDs.
Definition: nemesis_io_helper.h:519
libMesh::ExodusII_IO_Helper::num_nodes
int num_nodes
Definition: exodusII_io_helper.h:507
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::node_ptr
virtual const Node * node_ptr(const dof_id_type i) const =0
libMesh::Nemesis_IO_Helper::global_elem_blk_ids
std::vector< int > global_elem_blk_ids
Read the global element block IDs and counts.
Definition: nemesis_io_helper.h:426
libMesh::Nemesis_IO_Helper::num_border_nodes
int num_border_nodes
The number of FEM nodes local to a processor but residing in an element which also has FEM nodes on o...
Definition: nemesis_io_helper.h:460
libMesh::MeshBase::parallel_n_elem
virtual dof_id_type parallel_n_elem() const =0
libMesh::if
if(subdm)
Definition: petsc_dm_wrapper.C:77
libMesh::ExodusII_IO_Helper::num_side_sets
int num_side_sets
Definition: exodusII_io_helper.h:526
libMesh::Nemesis_IO_Helper::put_node_map
void put_node_map(std::vector< int > &node_mapi, std::vector< int > &node_mapb, std::vector< int > &node_mape)
Outputs IDs of internal, border, and external nodes.
Definition: nemesis_io_helper.C:645
libMesh::Nemesis_IO_Helper::put_ns_param_global
void put_ns_param_global(std::vector< int > &global_nodeset_ids, std::vector< int > &num_global_node_counts, std::vector< int > &num_global_node_df_counts)
This function writes information about global node sets.
Definition: nemesis_io_helper.C:510
libMesh::Nemesis_IO_Helper::compute_communication_map_parameters
void compute_communication_map_parameters()
This function determines the communication map parameters which will eventually be written to file.
Definition: nemesis_io_helper.C:1147
libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes
void compute_internal_and_border_elems_and_internal_nodes(const MeshBase &pmesh)
This function constructs the set of border and internal element IDs and internal node IDs present on ...
Definition: nemesis_io_helper.C:1225
libMesh::Nemesis_IO_Helper::num_proc_in_file
int num_proc_in_file
The number of processors for which the NEMESIS I file stores information.
Definition: nemesis_io_helper.h:388
libMesh::DofObject::processor_id
processor_id_type processor_id() const
Definition: dof_object.h:829
libMesh::Nemesis_IO_Helper::write_exodus_initialization_info
void write_exodus_initialization_info(const MeshBase &pmesh, const std::string &title)
This function writes exodus-specific initialization information.
Definition: nemesis_io_helper.C:933
libMesh::Nemesis_IO_Helper::num_global_node_counts
std::vector< int > num_global_node_counts
Definition: nemesis_io_helper.h:417
libMesh::BoundaryInfo::build_side_list
void build_side_list(std::vector< dof_id_type > &element_id_list, std::vector< unsigned short int > &side_list, std::vector< boundary_id_type > &bc_id_list) const
Creates a list of element numbers, sides, and ids for those sides.
Definition: boundary_info.C:1976
libMesh::Nemesis_IO_Helper::compute_elem_communication_maps
void compute_elem_communication_maps()
This function computes element communication maps (really just packs vectors) in preparation for writ...
Definition: nemesis_io_helper.C:995
libMesh::NumericVector< Number >
libMesh::ExodusII_IO_Helper::num_dim
int num_dim
Definition: exodusII_io_helper.h:498
libMesh::DofObject::dof_number
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:956
libMesh::Nemesis_IO_Helper::ftype
char ftype
The type of file to be written.
Definition: nemesis_io_helper.h:395
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::Nemesis_IO_Helper::global_nodeset_ids
std::vector< int > global_nodeset_ids
Containers for reading global nodeset information.
Definition: nemesis_io_helper.h:416
libMesh::Nemesis_IO_Helper::elem_cmap_elem_ids
std::vector< std::vector< int > > elem_cmap_elem_ids
3 vectors of vectors for storing element communication IDs for this processor.
Definition: nemesis_io_helper.h:565
libMesh::Nemesis_IO_Helper::get_elem_map
void get_elem_map()
Definition: nemesis_io_helper.C:257
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::NumericVector::localize
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.
libMesh::Nemesis_IO_Helper::global_elem_blk_cnts
std::vector< int > global_elem_blk_cnts
Definition: nemesis_io_helper.h:427
libMesh::Nemesis_IO_Helper::initialize
virtual void initialize(std::string title, const MeshBase &mesh, bool use_discontinuous=false) override
Specialization of the initialize function from ExodusII_IO_Helper that also writes global initial dat...
Definition: nemesis_io_helper.C:766
libMesh::Nemesis_IO_Helper::num_elem_cmaps
int num_elem_cmaps
The number of elemental communication maps for this processor.
Definition: nemesis_io_helper.h:494
libMesh::Nemesis_IO_Helper::node_mape
std::vector< int > node_mape
Vector which stores external node IDs.
Definition: nemesis_io_helper.h:533
libMesh::ExodusII_IO_Helper::elem_var_names
std::vector< std::string > elem_var_names
Definition: exodusII_io_helper.h:645
libMesh::ParallelObject::n_processors
processor_id_type n_processors() const
Definition: parallel_object.h:100
libMesh::Nemesis_IO_Helper::write_element_values
void write_element_values(const MeshBase &mesh, const EquationSystems &es, const std::vector< std::pair< unsigned int, unsigned int >> &var_nums, int timestep, const std::vector< std::set< subdomain_id_type >> &vars_active_subdomains)
Writes the vector of elemental variable values, one variable and one subdomain at a time.
Definition: nemesis_io_helper.C:2596
libMesh::Nemesis_IO_Helper::node_cmap_proc_ids
std::vector< std::vector< int > > node_cmap_proc_ids
Definition: nemesis_io_helper.h:557
libMesh::Nemesis_IO_Helper::proc_nodes_touched_intersections
std::map< unsigned, std::set< unsigned > > proc_nodes_touched_intersections
Another map to store sets of intersections with each other processor (other than ourself,...
Definition: nemesis_io_helper.h:588
libMesh::Nemesis_IO_Helper::internal_elem_ids
std::set< unsigned > internal_elem_ids
A set of internal elem IDs for this processor.
Definition: nemesis_io_helper.h:613
libMesh::Nemesis_IO_Helper::put_node_cmap
void put_node_cmap(std::vector< std::vector< int >> &node_cmap_node_ids, std::vector< std::vector< int >> &node_cmap_proc_ids)
Outputs all of the nodal communication maps for this processor.
Definition: nemesis_io_helper.C:600
node_set
std::unordered_set< const Node * > & node_set
Definition: mesh_tools.C:2135
libMesh::Nemesis_IO_Helper::elem_mapi
std::vector< int > elem_mapi
Vector which stores internal element IDs.
Definition: nemesis_io_helper.h:503
libMesh::Nemesis_IO_Helper::put_ss_param_global
void put_ss_param_global(std::vector< int > &global_sideset_ids, std::vector< int > &num_global_side_counts, std::vector< int > &num_global_side_df_counts)
This function writes information about global side sets.
Definition: nemesis_io_helper.C:530
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::ExodusII_IO_Helper::_elem_vars_initialized
bool _elem_vars_initialized
Definition: exodusII_io_helper.h:719
libMesh::Nemesis_IO_Helper::get_node_map
void get_node_map()
Definition: nemesis_io_helper.C:288
libMesh::System::current_local_solution
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
Definition: system.h:1551
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::DofMap::SCALAR_dof_indices
void SCALAR_dof_indices(std::vector< dof_id_type > &di, const unsigned int vn, const bool old_dofs=false) const
Fills the vector di with the global degree of freedom indices corresponding to the SCALAR variable vn...
Definition: dof_map.C:2488
libMesh::MeshBase::active_local_elements_end
virtual element_iterator active_local_elements_end()=0
libMesh::BoundaryInfo::get_node_boundary_ids
const std::set< boundary_id_type > & get_node_boundary_ids() const
Definition: boundary_info.h:819
libMesh::Nemesis_IO_Helper::write_sidesets
virtual void write_sidesets(const MeshBase &mesh) override
Writes the sidesets for this processor.
Definition: nemesis_io_helper.C:2084
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::Nemesis_IO_Helper::node_cmap_ids
std::vector< int > node_cmap_ids
Vectors for storing the communication map parameters.
Definition: nemesis_io_helper.h:542
libMesh::ExodusII_IO_Helper::num_node_sets
int num_node_sets
Definition: exodusII_io_helper.h:523
libMesh::MeshBase::parallel_n_nodes
virtual dof_id_type parallel_n_nodes() const =0
libMesh::ExodusII_IO_Helper::num_elem
int num_elem
Definition: exodusII_io_helper.h:510
libMesh::ExodusII_IO_Helper::exodus_elem_num_to_libmesh
std::vector< int > exodus_elem_num_to_libmesh
Definition: exodusII_io_helper.h:619
libMesh::Nemesis_IO_Helper::compute_node_maps
void compute_node_maps()
Compute the node maps (really just pack vectors) which map the nodes to internal, border,...
Definition: nemesis_io_helper.C:1046
libMesh::Nemesis_IO_Helper::elem_cmap_ids
std::vector< int > elem_cmap_ids
Definition: nemesis_io_helper.h:544
libMesh::Nemesis_IO_Helper::put_n_coord
void put_n_coord(unsigned start_node_num, unsigned num_nodes, std::vector< Real > &x_coor, std::vector< Real > &y_coor, std::vector< Real > &z_coor)
Writes the specified number of coordinate values starting at the specified index.
Definition: nemesis_io_helper.C:700
libMesh::Nemesis_IO_Helper::num_nodes_global
int num_nodes_global
Global initial information.
Definition: nemesis_io_helper.h:371
libMesh::Nemesis_IO_Helper::get_eb_info_global
void get_eb_info_global()
Definition: nemesis_io_helper.C:179
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
libMesh::Nemesis_IO_Helper::num_elem_blks_global
int num_elem_blks_global
Definition: nemesis_io_helper.h:373
libMesh::Nemesis_IO_Helper::local_subdomain_counts
std::map< subdomain_id_type, unsigned > local_subdomain_counts
This map keeps track of the number of elements in each subdomain (block) for this processor.
Definition: nemesis_io_helper.h:575
libMesh::Nemesis_IO_Helper::node_mapb
std::vector< int > node_mapb
Vector which stores border node IDs.
Definition: nemesis_io_helper.h:526
libMesh::Nemesis_IO_Helper::nemesis_err_flag
int nemesis_err_flag
Member data.
Definition: nemesis_io_helper.h:365
libMesh::ExodusII_IO_Helper::write_nodal_values
void write_nodal_values(int var_id, const std::vector< Real > &values, int timestep)
Writes the vector of values to a nodal variable.
Definition: exodusII_io_helper.C:3070
libMesh::ExodusII_IO_Helper::close
void close()
Closes the ExodusII mesh file.
Definition: exodusII_io_helper.C:1228
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::Nemesis_IO_Helper::proc_border_elem_sets
std::map< unsigned, std::set< std::pair< unsigned, unsigned > > > proc_border_elem_sets
Map between processor ID and (element,side) pairs bordering that processor ID.
Definition: nemesis_io_helper.h:598
libMesh::System::variable_type
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2233
libMesh::Nemesis_IO_Helper::nodes_attached_to_local_elems
std::set< int > nodes_attached_to_local_elems
libMesh numbered node ids attached to local elems.
Definition: nemesis_io_helper.h:432
libMesh::Nemesis_IO_Helper::get_ns_param_global
void get_ns_param_global()
Definition: nemesis_io_helper.C:147
libMesh::Nemesis_IO_Helper::border_node_ids
std::set< unsigned > border_node_ids
The set which will eventually contain the IDs of "border nodes".
Definition: nemesis_io_helper.h:581
libMesh::ExodusII_IO_Helper::ELEMENTAL
Definition: exodusII_io_helper.h:692
libMesh::ExodusII_IO_Helper::exodus_node_num_to_libmesh
std::vector< int > exodus_node_num_to_libmesh
Definition: exodusII_io_helper.h:624
distance
Real distance(const Point &p)
Definition: subdomains_ex3.C:50
libMesh::ExodusII_IO_Helper::z
std::vector< Real > z
Definition: exodusII_io_helper.h:608
libMesh::ExodusII_IO_Helper::NamesData::get_char_star_star
char ** get_char_star_star()
Provide access to the underlying C data table.
Definition: exodusII_io_helper.C:3392
libMesh::Nemesis_IO_Helper::put_init_info
void put_init_info(unsigned num_proc, unsigned num_proc_in_file, const char *ftype)
Writing functions.
Definition: nemesis_io_helper.C:461
libMesh::Nemesis_IO_Helper::proc_border_elem_sets_iterator
std::map< unsigned, std::set< std::pair< unsigned, unsigned > > >::iterator proc_border_elem_sets_iterator
Typedef for an iterator into the data structure above.
Definition: nemesis_io_helper.h:603
swap
void swap(Iterator &lhs, Iterator &rhs)
swap, used to implement op=
Definition: variant_filter_iterator.h:478
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::ExodusII_IO_Helper::num_elem_blk
int num_elem_blk
Definition: exodusII_io_helper.h:513
libMesh::Nemesis_IO_Helper::num_global_side_counts
std::vector< int > num_global_side_counts
Definition: nemesis_io_helper.h:405
value
static const bool value
Definition: xdr_io.C:56
libMesh::Nemesis_IO_Helper::num_node_cmaps
int num_node_cmaps
The number of nodal communication maps for this processor.
Definition: nemesis_io_helper.h:487
libMesh::BoundaryInfo::build_node_list
void build_node_list(std::vector< dof_id_type > &node_id_list, std::vector< boundary_id_type > &bc_id_list) const
Creates a list of nodes and ids for those nodes.
Definition: boundary_info.C:1704
libMesh::Nemesis_IO_Helper::write_nodal_solution
void write_nodal_solution(const NumericVector< Number > &parallel_soln, const std::vector< std::string > &names, int timestep, const std::vector< std::string > &output_names)
Takes a parallel solution vector containing the node-major solution vector for all variables and outp...
Definition: nemesis_io_helper.C:2410
libMesh::Nemesis_IO_Helper::internal_node_ids
std::set< unsigned > internal_node_ids
A set of internal node IDs for this processor.
Definition: nemesis_io_helper.h:608
libMesh::ExodusII_IO_Helper::verbose
bool verbose
Definition: exodusII_io_helper.h:663
libMesh::Nemesis_IO_Helper::write_nodal_coordinates
virtual void write_nodal_coordinates(const MeshBase &mesh, bool use_discontinuous=false) override
This function is specialized from ExodusII_IO_Helper to write only the nodal coordinates stored on th...
Definition: nemesis_io_helper.C:2228
libMesh::Nemesis_IO_Helper::write_nodesets
virtual void write_nodesets(const MeshBase &mesh) override
Writes the nodesets for this processor.
Definition: nemesis_io_helper.C:1933
libMesh::Nemesis_IO_Helper::get_elem_cmap
void get_elem_cmap()
Definition: nemesis_io_helper.C:410
libMesh::Nemesis_IO_Helper::put_loadbal_param
void put_loadbal_param(unsigned num_internal_nodes, unsigned num_border_nodes, unsigned num_external_nodes, unsigned num_internal_elems, unsigned num_border_elems, unsigned num_node_cmaps, unsigned num_elem_cmaps)
Writes load balance parameters, some of which are described below: .) num_internal_nodes - nodes "who...
Definition: nemesis_io_helper.C:550
libMesh::ExodusII_IO_Helper::num_elem_vars
int num_elem_vars
Definition: exodusII_io_helper.h:642
libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks
void compute_num_global_elem_blocks(const MeshBase &pmesh)
This function uses global communication routines to determine the number of element blocks across the...
Definition: nemesis_io_helper.C:1567
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::Nemesis_IO_Helper::put_init_global
void put_init_global(dof_id_type num_nodes_global, dof_id_type num_elems_global, unsigned num_elem_blks_global, unsigned num_node_sets_global, unsigned num_side_sets_global)
Writes global information including: .) global number of nodes .) global number of elems ....
Definition: nemesis_io_helper.C:477
libMesh::Nemesis_IO_Helper::get_init_global
void get_init_global()
Reading functions.
Definition: nemesis_io_helper.C:88
data
IterBase * data
Ideally this private member data should have protected access.
Definition: variant_filter_iterator.h:337
libMesh::Nemesis_IO_Helper::num_internal_elems
int num_internal_elems
The number of internal FEM elements.
Definition: nemesis_io_helper.h:473
libMesh::ExodusII_IO_Helper::NamesData
This class is useful for managing anything that requires a char ** input/output in ExodusII file.
Definition: exodusII_io_helper.h:984
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::ExodusII_IO_Helper::NamesData::push_back_entry
void push_back_entry(const std::string &name)
Adds another name to the current data table.
Definition: exodusII_io_helper.C:3376
libMesh::SCALAR
Definition: enum_fe_family.h:58
libMesh::Nemesis_IO_Helper::compute_node_communication_maps
void compute_node_communication_maps()
Compute the node communication maps (really just pack vectors) in preparation for writing them to fil...
Definition: nemesis_io_helper.C:1075
libMesh::TestClass
Definition: id_types.h:33
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::Nemesis_IO_Helper::num_global_side_df_counts
std::vector< int > num_global_side_df_counts
Definition: nemesis_io_helper.h:406
libMesh::ParallelObject
An object whose state is distributed along a set of processors.
Definition: parallel_object.h:55
libMesh::Nemesis_IO_Helper::subdomain_map
std::map< subdomain_id_type, std::vector< dof_id_type > > subdomain_map
Map of subdomains to element numbers.
Definition: nemesis_io_helper.h:437
libMesh::Elem::build
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:246
libMesh::out
OStreamProxy out
libMesh::Nemesis_IO_Helper::node_cmap_node_ids
std::vector< std::vector< int > > node_cmap_node_ids
2 vectors of vectors for storing the node communication IDs for this processor.
Definition: nemesis_io_helper.h:556
libMesh::ExodusII_IO_Helper::x
std::vector< Real > x
Definition: exodusII_io_helper.h:602
libMesh::Nemesis_IO_Helper::num_global_node_df_counts
std::vector< int > num_global_node_df_counts
Definition: nemesis_io_helper.h:418
libMesh::Nemesis_IO_Helper::block_id_to_elem_connectivity
std::map< int, std::vector< int > > block_id_to_elem_connectivity
This is the block connectivity, i.e.
Definition: nemesis_io_helper.h:443
libMesh::ExodusII_IO_Helper::ex_err
int ex_err
Definition: exodusII_io_helper.h:495
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::Nemesis_IO_Helper::get_node_cmap
void get_node_cmap()
Definition: nemesis_io_helper.C:368
libMesh::Nemesis_IO_Helper::put_eb_info_global
void put_eb_info_global(std::vector< int > &global_elem_blk_ids, std::vector< int > &global_elem_blk_cnts)
Writes global block information to the file .) global_elem_blk_ids - list of block IDs for all blocks...
Definition: nemesis_io_helper.C:496
libMesh::Nemesis_IO_Helper::build_element_and_node_maps
void build_element_and_node_maps(const MeshBase &pmesh)
This function builds the libmesh -> exodus and exodus -> libmesh node and element maps.
Definition: nemesis_io_helper.C:1667
libMesh::Quality::name
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
libMesh::ExodusII_IO_Helper::ex_id
int ex_id
Definition: exodusII_io_helper.h:492
libMesh::Nemesis_IO_Helper::get_loadbal_param
void get_loadbal_param()
Definition: nemesis_io_helper.C:227
libMesh::Nemesis_IO_Helper::num_border_elems
int num_border_elems
The number of border FEM elements.
Definition: nemesis_io_helper.h:480