1 // Copyright(C) 1999-2021 National Technology & Engineering Solutions
2 // of Sandia, LLC (NTESS).  Under the terms of Contract DE-NA0003525 with
3 // NTESS, the U.S. Government retains certain rights in this software.
4 //
5 // See packages/seacas/LICENSE for details
6 
7 #include <Ioss_Assembly.h>
8 #include <Ioss_Blob.h>
9 #include <Ioss_CommSet.h>
10 #include <Ioss_CoordinateFrame.h>
11 #include <Ioss_DBUsage.h>
12 #include <Ioss_DatabaseIO.h>
13 #include <Ioss_EdgeBlock.h>
14 #include <Ioss_EdgeSet.h>
15 #include <Ioss_ElementBlock.h>
16 #include <Ioss_ElementSet.h>
17 #include <Ioss_ElementTopology.h>
18 #include <Ioss_EntityBlock.h>
19 #include <Ioss_EntityType.h>
20 #include <Ioss_FaceBlock.h>
21 #include <Ioss_FaceSet.h>
22 #include <Ioss_Field.h>
23 #include <Ioss_GroupingEntity.h>
24 #include <Ioss_NodeBlock.h>
25 #include <Ioss_NodeSet.h>
26 #include <Ioss_Property.h>
27 #include <Ioss_PropertyManager.h>
28 #include <Ioss_Region.h>
29 #include <Ioss_SideBlock.h>
30 #include <Ioss_SideSet.h>
31 #include <Ioss_SmartAssert.h>
32 #include <Ioss_Sort.h>
33 #include <Ioss_State.h>
34 #include <Ioss_StructuredBlock.h>
35 
36 #include <algorithm>
37 #include <cctype>
38 #include <climits>
39 #include <cstddef>
40 #include <fmt/ostream.h>
41 #include <iomanip>
42 #include <iostream>
43 #include <map>
44 #include <string>
45 #ifndef _MSC_VER
46 #include <unistd.h>
47 #endif
48 #include <utility>
49 #include <vector>
50 
51 namespace {
id_str()52   std::string id_str() { return std::string("id"); }
db_name_str()53   std::string db_name_str() { return std::string("db_name"); }
orig_topo_str()54   std::string orig_topo_str() { return std::string("original_topology_type"); }
orig_block_order()55   std::string orig_block_order() { return std::string("original_block_order"); }
56 
57   template <typename T>
get_entity_internal(int64_t id,const std::vector<T> & entities)58   Ioss::GroupingEntity *get_entity_internal(int64_t id, const std::vector<T> &entities)
59   {
60     for (auto ent : entities) {
61       if (ent->property_exists(id_str())) {
62         if (id == ent->get_property(id_str()).get_int()) {
63           return ent;
64         }
65       }
66     }
67     return nullptr;
68   }
69 
70   template <typename T>
internal_get_variable_count(const std::vector<T> & entities,Ioss::Field::RoleType role)71   size_t internal_get_variable_count(const std::vector<T> &entities, Ioss::Field::RoleType role)
72   {
73     Ioss::NameList names;
74     for (auto ent : entities) {
75       ent->field_describe(role, &names);
76     }
77     Ioss::Utils::uniquify(names);
78     return names.size();
79   }
80 
get_variable_count(const std::vector<T> & entities)81   template <typename T> size_t get_variable_count(const std::vector<T> &entities)
82   {
83     return internal_get_variable_count(entities, Ioss::Field::TRANSIENT);
84   }
85 
get_reduction_variable_count(const std::vector<T> & entities)86   template <typename T> size_t get_reduction_variable_count(const std::vector<T> &entities)
87   {
88     return internal_get_variable_count(entities, Ioss::Field::REDUCTION);
89   }
90 
get_entity_count(const std::vector<T> & entities)91   template <typename T> int64_t get_entity_count(const std::vector<T> &entities)
92   {
93     int64_t count = 0;
94     for (auto ent : entities) {
95       count += ent->entity_count();
96     }
97     return count;
98   }
99 
update_database(const Ioss::Region * region,Ioss::GroupingEntity * entity)100   void update_database(const Ioss::Region *region, Ioss::GroupingEntity *entity)
101   {
102     entity->reset_database(region->get_database());
103   }
104 
update_database(const Ioss::Region * region,Ioss::SideSet * sset)105   void update_database(const Ioss::Region *region, Ioss::SideSet *sset)
106   {
107     sset->reset_database(region->get_database());
108     const auto &blocks = sset->get_side_blocks();
109     for (const auto &block : blocks) {
110       block->reset_database(region->get_database());
111     }
112   }
113 
check_for_duplicate_names(const Ioss::Region * region,const Ioss::GroupingEntity * entity)114   void check_for_duplicate_names(const Ioss::Region *region, const Ioss::GroupingEntity *entity)
115   {
116     const std::string &name = entity->name();
117 
118     // See if any alias with this name...
119     std::string alias = region->get_alias__(name);
120 
121     if (!alias.empty()) {
122       // There is an entity with this name...
123       const Ioss::GroupingEntity *old_ge = region->get_entity(name);
124 
125       if (old_ge != nullptr &&
126           !(old_ge->type() == Ioss::SIDEBLOCK || old_ge->type() == Ioss::SIDESET)) {
127         std::string        filename = region->get_database()->get_filename();
128         int64_t            id1      = entity->get_optional_property(id_str(), 0);
129         int64_t            id2      = old_ge->get_optional_property(id_str(), 0);
130         std::ostringstream errmsg;
131         fmt::print(errmsg,
132                    "ERROR: There are multiple blocks or sets with the same name defined in the "
133                    "database file '{}'.\n"
134                    "\tBoth {} {} and {} {} are named '{}'.  All names must be unique.",
135                    filename, entity->type_string(), id1, old_ge->type_string(), id2, name);
136         IOSS_ERROR(errmsg);
137       }
138     }
139   }
140 
numberOfBits(unsigned x)141   constexpr unsigned numberOfBits(unsigned x) { return x < 2 ? x : 1 + numberOfBits(x >> 1); }
142 
compute_hash(Ioss::GroupingEntity * entity,size_t which)143   size_t compute_hash(Ioss::GroupingEntity *entity, size_t which)
144   {
145     // Can add more properties and or fields later.  For now just do
146     // name and optional id.
147     size_t hash = entity->hash();
148     if (entity->property_exists(id_str())) {
149       hash += which * entity->get_property(id_str()).get_int();
150     }
151     return hash;
152   }
153 
154   template <typename T>
compute_hashes(const std::vector<T> & entities,std::array<size_t,Ioss::entityTypeCount> & hashes,Ioss::EntityType type)155   void compute_hashes(const std::vector<T> &                     entities,
156                       std::array<size_t, Ioss::entityTypeCount> &hashes, Ioss::EntityType type)
157   {
158     auto index = numberOfBits(type) - 1;
159     SMART_ASSERT(index < hashes.size())(type)(index)(hashes.size());
160 
161     size_t which = 1;
162     for (const auto &entity : entities) {
163       hashes[index] += compute_hash(entity, which++);
164     }
165   }
166 
check_hashes(const std::vector<size_t> & min_hash,const std::vector<size_t> & max_hash,Ioss::EntityType type)167   bool check_hashes(const std::vector<size_t> &min_hash, const std::vector<size_t> &max_hash,
168                     Ioss::EntityType type)
169   {
170     auto index = numberOfBits(type) - 1;
171     SMART_ASSERT(index < min_hash.size())(type)(index)(min_hash.size());
172     return (min_hash[index] == max_hash[index]);
173   }
174 
175   template <typename T>
report_inconsistency(const std::vector<T> & entities,Ioss::ParallelUtils & util)176   void report_inconsistency(const std::vector<T> &entities, Ioss::ParallelUtils &util)
177   {
178     // Know that there is some mismatch in name or (optional)id.  Let user know where...
179     std::vector<size_t> hashes;
180 
181     size_t which = 1;
182     hashes.reserve(entities.size());
183     for (const auto &entity : entities) {
184       hashes.push_back(compute_hash(entity, which++));
185     }
186 
187     std::ostringstream errmsg;
188     fmt::print(errmsg, "IOSS: ERROR: Parallel Consistency Error.\n\t\t");
189 
190     auto min_hash = hashes;
191     auto max_hash = hashes;
192     // Now find mismatched location...
193     util.global_array_minmax(min_hash, Ioss::ParallelUtils::DO_MIN);
194     util.global_array_minmax(max_hash, Ioss::ParallelUtils::DO_MAX);
195 
196     if (util.parallel_rank() == 0) {
197       int count = 0;
198       for (size_t i = 0; i < hashes.size(); i++) {
199         if (min_hash[i] != max_hash[i]) {
200           auto ge = entities[i];
201           if (count == 0) {
202             fmt::print(errmsg, "{}(s) ", ge->type_string());
203           }
204           else {
205             fmt::print(errmsg, ", ");
206           }
207           fmt::print(errmsg, "'{}'", ge->name());
208           count++;
209         }
210       }
211       fmt::print(errmsg,
212                  " {} not consistently defined on all processors.\n\t\t"
213                  "Check that name and id matches across processors.\n",
214                  (count == 1 ? "is" : "are"));
215       IOSS_ERROR(errmsg);
216     }
217   }
218 
check_parallel_consistency(const Ioss::Region & region)219   bool check_parallel_consistency(const Ioss::Region &region)
220   {
221     if (!region.get_database()->is_parallel()) {
222       return true;
223     }
224 
225     // Want a good approximate test that the grouping entity lists on
226     // all processor contain the same entities in the same order.
227     // We will say an entity is the same if the name and optional id match.
228     //
229     // Hash the name and multiply it by position in list and add id+1.
230     // Do this for each type separately...  Then verify that they
231     // match on all processors...
232     std::array<size_t, Ioss::entityTypeCount> hashes{};
233 
234     compute_hashes(region.get_node_blocks(), hashes, Ioss::NODEBLOCK);
235     compute_hashes(region.get_edge_blocks(), hashes, Ioss::EDGEBLOCK);
236     compute_hashes(region.get_face_blocks(), hashes, Ioss::FACEBLOCK);
237     compute_hashes(region.get_element_blocks(), hashes, Ioss::ELEMENTBLOCK);
238     compute_hashes(region.get_nodesets(), hashes, Ioss::NODESET);
239     compute_hashes(region.get_edgesets(), hashes, Ioss::EDGESET);
240     compute_hashes(region.get_facesets(), hashes, Ioss::FACESET);
241     compute_hashes(region.get_elementsets(), hashes, Ioss::ELEMENTSET);
242     compute_hashes(region.get_sidesets(), hashes, Ioss::SIDESET);
243     compute_hashes(region.get_commsets(), hashes, Ioss::COMMSET);
244     compute_hashes(region.get_structured_blocks(), hashes, Ioss::STRUCTUREDBLOCK);
245     compute_hashes(region.get_assemblies(), hashes, Ioss::ASSEMBLY);
246     compute_hashes(region.get_blobs(), hashes, Ioss::BLOB);
247 
248     auto                util = region.get_database()->util();
249     std::vector<size_t> min_hash(hashes.begin(), hashes.end());
250     std::vector<size_t> max_hash(hashes.begin(), hashes.end());
251     util.global_array_minmax(min_hash, Ioss::ParallelUtils::DO_MIN);
252     util.global_array_minmax(max_hash, Ioss::ParallelUtils::DO_MAX);
253 
254     bool differ = false;
255     if (!check_hashes(min_hash, max_hash, Ioss::NODEBLOCK)) {
256       report_inconsistency(region.get_node_blocks(), util);
257       differ = true;
258     }
259     if (!check_hashes(min_hash, max_hash, Ioss::EDGEBLOCK)) {
260       report_inconsistency(region.get_edge_blocks(), util);
261       differ = true;
262     }
263     if (!check_hashes(min_hash, max_hash, Ioss::FACEBLOCK)) {
264       report_inconsistency(region.get_face_blocks(), util);
265       differ = true;
266     }
267     if (!check_hashes(min_hash, max_hash, Ioss::ELEMENTBLOCK)) {
268       report_inconsistency(region.get_element_blocks(), util);
269       differ = true;
270     }
271     if (!check_hashes(min_hash, max_hash, Ioss::NODESET)) {
272       report_inconsistency(region.get_nodesets(), util);
273       differ = true;
274     }
275     if (!check_hashes(min_hash, max_hash, Ioss::EDGESET)) {
276       report_inconsistency(region.get_edgesets(), util);
277       differ = true;
278     }
279     if (!check_hashes(min_hash, max_hash, Ioss::FACESET)) {
280       report_inconsistency(region.get_facesets(), util);
281       differ = true;
282     }
283     if (!check_hashes(min_hash, max_hash, Ioss::ELEMENTSET)) {
284       report_inconsistency(region.get_elementsets(), util);
285       differ = true;
286     }
287     if (!check_hashes(min_hash, max_hash, Ioss::SIDESET)) {
288       report_inconsistency(region.get_sidesets(), util);
289       differ = true;
290     }
291     if (!check_hashes(min_hash, max_hash, Ioss::COMMSET)) {
292       report_inconsistency(region.get_commsets(), util);
293       differ = true;
294     }
295     if (!check_hashes(min_hash, max_hash, Ioss::STRUCTUREDBLOCK)) {
296       report_inconsistency(region.get_structured_blocks(), util);
297       differ = true;
298     }
299     if (!check_hashes(min_hash, max_hash, Ioss::ASSEMBLY)) {
300       report_inconsistency(region.get_assemblies(), util);
301       differ = true;
302     }
303     if (!check_hashes(min_hash, max_hash, Ioss::BLOB)) {
304       report_inconsistency(region.get_blobs(), util);
305       differ = true;
306     }
307     return !differ;
308   }
309 
is_input_or_appending_output(const Ioss::DatabaseIO * iodatabase)310   bool is_input_or_appending_output(const Ioss::DatabaseIO *iodatabase)
311   {
312     return iodatabase->is_input() || iodatabase->open_create_behavior() == Ioss::DB_APPEND ||
313            iodatabase->open_create_behavior() == Ioss::DB_MODIFY;
314   }
315 } // namespace
316 
317 namespace Ioss {
318 
319   /** \brief Constructor reads in all metadata from disk.
320    *
321    *  This constructor connects this region to the database, opens the
322    *  underlying file, reads all the metadata in the file into Region
323    *  and its subentities, and closes the underlying file. Region properties,
324    *  such as spatial_dimension, element_block_count, element_count, etc, are
325    *  also added to the Region's property manager.
326    *
327    *  \param[in] iodatabase The name of the database associated with the Region.
328    *  \param[in] my_name The name of the Region.
329    *
330    */
Region(DatabaseIO * iodatabase,const std::string & my_name)331   Region::Region(DatabaseIO *iodatabase, const std::string &my_name)
332       : GroupingEntity(iodatabase, my_name, 1)
333   {
334     SMART_ASSERT(iodatabase != nullptr);
335     iodatabase->set_region(this);
336 
337     if (iodatabase->usage() != Ioss::WRITE_HEARTBEAT &&
338         (is_input_or_appending_output(iodatabase))) {
339       Region::begin_mode(STATE_DEFINE_MODEL);
340       iodatabase->read_meta_data();
341       Region::end_mode(STATE_DEFINE_MODEL);
342       if (iodatabase->open_create_behavior() != Ioss::DB_APPEND &&
343           iodatabase->open_create_behavior() != Ioss::DB_MODIFY) {
344         modelDefined     = true;
345         transientDefined = true;
346         Region::begin_mode(STATE_READONLY);
347       }
348     }
349 
350     properties.add(Property(this, "spatial_dimension", Property::INTEGER));
351     properties.add(Property(this, "node_block_count", Property::INTEGER));
352     properties.add(Property(this, "edge_block_count", Property::INTEGER));
353     properties.add(Property(this, "face_block_count", Property::INTEGER));
354     properties.add(Property(this, "element_block_count", Property::INTEGER));
355     properties.add(Property(this, "structured_block_count", Property::INTEGER));
356     properties.add(Property(this, "assembly_count", Property::INTEGER));
357     properties.add(Property(this, "blob_count", Property::INTEGER));
358     properties.add(Property(this, "side_set_count", Property::INTEGER));
359     properties.add(Property(this, "node_set_count", Property::INTEGER));
360     properties.add(Property(this, "edge_set_count", Property::INTEGER));
361     properties.add(Property(this, "face_set_count", Property::INTEGER));
362     properties.add(Property(this, "element_set_count", Property::INTEGER));
363     properties.add(Property(this, "comm_set_count", Property::INTEGER));
364     properties.add(Property(this, "node_count", Property::INTEGER));
365     properties.add(Property(this, "edge_count", Property::INTEGER));
366     properties.add(Property(this, "face_count", Property::INTEGER));
367     properties.add(Property(this, "element_count", Property::INTEGER));
368     properties.add(Property(this, "coordinate_frame_count", Property::INTEGER));
369     properties.add(Property(this, "state_count", Property::INTEGER));
370     properties.add(Property(this, "current_state", Property::INTEGER));
371     properties.add(Property(this, "database_name", Property::STRING));
372   }
373 
~Region()374   Region::~Region()
375   {
376     // Do anything to the database to make it consistent prior to closing and destructing...
377     get_database()->finalize_database();
378 
379     // Region owns all sub-grouping entities it contains...
380     try {
381       IOSS_FUNC_ENTER(m_);
382       for (auto nb : nodeBlocks) {
383         delete (nb);
384       }
385 
386       for (auto eb : edgeBlocks) {
387         delete (eb);
388       }
389 
390       for (auto fb : faceBlocks) {
391         delete (fb);
392       }
393 
394       for (auto eb : elementBlocks) {
395         delete (eb);
396       }
397 
398       for (auto sb : structuredBlocks) {
399         delete (sb);
400       }
401 
402       for (auto ss : sideSets) {
403         delete (ss);
404       }
405 
406       for (auto ns : nodeSets) {
407         delete (ns);
408       }
409 
410       for (auto es : edgeSets) {
411         delete (es);
412       }
413 
414       for (auto fs : faceSets) {
415         delete (fs);
416       }
417 
418       for (auto es : elementSets) {
419         delete (es);
420       }
421 
422       for (auto cs : commSets) {
423         delete (cs);
424       }
425 
426       for (auto as : assemblies) {
427         delete (as);
428       }
429 
430       for (auto bl : blobs) {
431         delete (bl);
432       }
433 
434       // Region owns the database pointer even though other entities use it.
435       GroupingEntity::really_delete_database();
436     }
437     catch (...) {
438     }
439   }
440 
delete_database()441   void Region::delete_database() { GroupingEntity::really_delete_database(); }
442 
node_major()443   bool Region::node_major() const { return get_database()->node_major(); }
444 
mesh_type()445   MeshType Region::mesh_type() const
446   {
447     if (elementBlocks.empty() && structuredBlocks.empty()) {
448       return MeshType::UNSTRUCTURED;
449     }
450     if (!elementBlocks.empty() && !structuredBlocks.empty()) {
451       return MeshType::HYBRID;
452     }
453     if (!structuredBlocks.empty()) {
454       return MeshType::STRUCTURED;
455     }
456     SMART_ASSERT(!elementBlocks.empty());
457     return MeshType::UNSTRUCTURED;
458   }
459 
mesh_type_string()460   const std::string Region::mesh_type_string() const
461   {
462     switch (mesh_type()) {
463     case MeshType::UNKNOWN: return "Unknown";
464     case MeshType::HYBRID: return "Hybrid";
465     case MeshType::STRUCTURED: return "Structured";
466     case MeshType::UNSTRUCTURED: return "Unstructured";
467     }
468     SMART_ASSERT(1 == 0 && "Program Error");
469     return "Invalid";
470   }
471 
472   /** \brief Print a summary of entities in the region.
473    *
474    *  \param[in,out] strm The output stream to use for printing.
475    *  \param[in]     do_transient deprecated and ignored
476    */
output_summary(std::ostream & strm,bool)477   void Region::output_summary(std::ostream &strm, bool /* do_transient */) const
478   {
479     int64_t total_cells       = get_entity_count(get_structured_blocks());
480     int64_t total_fs_faces    = get_entity_count(get_facesets());
481     int64_t total_ns_nodes    = get_entity_count(get_nodesets());
482     int64_t total_es_edges    = get_entity_count(get_edgesets());
483     int64_t total_es_elements = get_entity_count(get_elementsets());
484 
485     int64_t                       total_sides = 0;
486     const Ioss::SideSetContainer &sss         = get_sidesets();
487     for (auto fs : sss) {
488       total_sides += get_entity_count(fs->get_side_blocks());
489     }
490 
491     int64_t total_nodes    = get_property("node_count").get_int();
492     int64_t total_elements = get_property("element_count").get_int();
493     auto    max_entity = std::max({total_sides, total_es_elements, total_fs_faces, total_es_edges,
494                                 total_ns_nodes, total_cells, total_nodes, total_elements});
495 
496     int64_t num_ts = get_property("state_count").get_int();
497     auto    max_sb = std::max(
498         {get_property("spatial_dimension").get_int(), get_property("node_block_count").get_int(),
499          get_property("edge_block_count").get_int(), get_property("face_block_count").get_int(),
500          get_property("element_block_count").get_int(),
501          get_property("structured_block_count").get_int(), get_property("node_set_count").get_int(),
502          get_property("edge_set_count").get_int(), get_property("face_set_count").get_int(),
503          get_property("element_set_count").get_int(), get_property("side_set_count").get_int(),
504          get_property("assembly_count").get_int(), get_property("blob_count").get_int(), num_ts});
505 
506     // Global variables transitioning from TRANSIENT to REDUCTION..
507     size_t num_glo_vars  = field_count(Ioss::Field::TRANSIENT);
508     size_t num_nod_vars  = get_variable_count(get_node_blocks());
509     size_t num_edg_vars  = get_variable_count(get_edge_blocks());
510     size_t num_fac_vars  = get_variable_count(get_face_blocks());
511     size_t num_ele_vars  = get_variable_count(get_element_blocks());
512     size_t num_str_vars  = get_variable_count(get_structured_blocks());
513     size_t num_ns_vars   = get_variable_count(get_nodesets());
514     size_t num_es_vars   = get_variable_count(get_edgesets());
515     size_t num_fs_vars   = get_variable_count(get_facesets());
516     size_t num_els_vars  = get_variable_count(get_elementsets());
517     size_t num_asm_vars  = get_variable_count(get_assemblies());
518     size_t num_blob_vars = get_variable_count(get_blobs());
519 
520     size_t num_glo_red_vars  = field_count(Ioss::Field::REDUCTION);
521     size_t num_nod_red_vars  = get_reduction_variable_count(get_node_blocks());
522     size_t num_edg_red_vars  = get_reduction_variable_count(get_edge_blocks());
523     size_t num_fac_red_vars  = get_reduction_variable_count(get_face_blocks());
524     size_t num_ele_red_vars  = get_reduction_variable_count(get_element_blocks());
525     size_t num_str_red_vars  = get_reduction_variable_count(get_structured_blocks());
526     size_t num_ns_red_vars   = get_reduction_variable_count(get_nodesets());
527     size_t num_es_red_vars   = get_reduction_variable_count(get_edgesets());
528     size_t num_fs_red_vars   = get_reduction_variable_count(get_facesets());
529     size_t num_els_red_vars  = get_reduction_variable_count(get_elementsets());
530     size_t num_asm_red_vars  = get_reduction_variable_count(get_assemblies());
531     size_t num_blob_red_vars = get_reduction_variable_count(get_blobs());
532 
533     size_t                       num_ss_vars = 0;
534     const Ioss::SideSetContainer fss         = get_sidesets();
535     for (auto fs : fss) {
536       num_ss_vars += get_variable_count(fs->get_side_blocks());
537     }
538 
539     auto max_vr    = std::max({num_glo_vars,     num_nod_vars,     num_ele_vars,     num_str_vars,
540                             num_ns_vars,      num_ss_vars,      num_edg_vars,     num_fac_vars,
541                             num_es_vars,      num_fs_vars,      num_els_vars,     num_blob_vars,
542                             num_asm_vars,     num_glo_red_vars, num_nod_red_vars, num_edg_red_vars,
543                             num_fac_red_vars, num_ele_red_vars, num_str_red_vars, num_ns_red_vars,
544                             num_es_red_vars,  num_fs_red_vars,  num_els_red_vars, num_asm_red_vars,
545                             num_blob_red_vars});
546     int  vr_width  = Ioss::Utils::number_width(max_vr, true) + 2;
547     int  num_width = Ioss::Utils::number_width(max_entity, true) + 2;
548     int  sb_width  = Ioss::Utils::number_width(max_sb, true) + 2;
549 
550     // clang-format off
551     fmt::print(
552         strm,
553         "\n Database: {0}\n"
554         " Mesh Type = {1}, {39}\n"
555         "                      {38:{24}s}\t                 {38:{23}s}\t Variables : Transient / Reduction\n"
556         " Spatial dimensions = {2:{24}L}\t                 {38:{23}s}\t Global     = {26:{25}L}\t{44:{25}L}\n"
557         " Node blocks        = {7:{24}L}\t Nodes         = {3:{23}L}\t Nodal      = {27:{25}L}\t{45:{25}L}\n"
558         " Edge blocks        = {8:{24}L}\t Edges         = {4:{23}L}\t Edge       = {33:{25}L}\t{46:{25}L}\n"
559         " Face blocks        = {9:{24}L}\t Faces         = {5:{23}L}\t Face       = {34:{25}L}\t{47:{25}L}\n"
560         " Element blocks     = {10:{24}L}\t Elements      = {6:{23}L}\t Element    = {28:{25}L}\t{48:{25}L}\n"
561         " Structured blocks  = {11:{24}L}\t Cells         = {17:{23}L}\t Structured = {29:{25}L}\t{49:{25}L}\n"
562         " Node sets          = {12:{24}L}\t Node list     = {18:{23}L}\t Nodeset    = {30:{25}L}\t{50:{25}L}\n"
563         " Edge sets          = {13:{24}L}\t Edge list     = {19:{23}L}\t Edgeset    = {35:{25}L}\t{51:{25}L}\n"
564         " Face sets          = {14:{24}L}\t Face list     = {20:{23}L}\t Faceset    = {36:{25}L}\t{52:{25}L}\n"
565         " Element sets       = {15:{24}L}\t Element list  = {21:{23}L}\t Elementset = {37:{25}L}\t{53:{25}L}\n"
566         " Element side sets  = {16:{24}L}\t Element sides = {22:{23}L}\t Sideset    = {31:{25}L}\n"
567         " Assemblies         = {40:{24}L}\t                 {38:{23}s}\t Assembly   = {41:{25}L}\t{54:{25}L}\n"
568         " Blobs              = {42:{24}L}\t                 {38:{23}s}\t Blob       = {43:{25}L}\t{55:{25}L}\n\n"
569         " Time steps         = {32:{24}L}\n",
570         get_database()->get_filename(), mesh_type_string(),
571         get_property("spatial_dimension").get_int(), get_property("node_count").get_int(),
572         get_property("edge_count").get_int(), get_property("face_count").get_int(),
573         get_property("element_count").get_int(), get_property("node_block_count").get_int(),
574         get_property("edge_block_count").get_int(), get_property("face_block_count").get_int(),
575         get_property("element_block_count").get_int(),
576         get_property("structured_block_count").get_int(), get_property("node_set_count").get_int(),
577         get_property("edge_set_count").get_int(), get_property("face_set_count").get_int(),
578         get_property("element_set_count").get_int(), get_property("side_set_count").get_int(),
579         total_cells, total_ns_nodes, total_es_edges, total_fs_faces, total_es_elements, total_sides,
580         num_width, sb_width, vr_width, num_glo_vars, num_nod_vars, num_ele_vars, num_str_vars,
581         num_ns_vars, num_ss_vars, num_ts, num_edg_vars, num_fac_vars, num_es_vars, num_fs_vars,
582         num_els_vars, " ", get_database()->get_format(), get_property("assembly_count").get_int(),
583         num_asm_vars, get_property("blob_count").get_int(), num_blob_vars, num_glo_red_vars,
584         num_nod_red_vars, num_edg_red_vars, num_fac_red_vars, num_ele_red_vars, num_str_red_vars,
585         num_ns_red_vars, num_es_red_vars, num_fs_red_vars, num_els_red_vars, num_asm_red_vars,
586         num_blob_red_vars);
587     // clang-format on
588   }
589 
590   /** \brief Set the Region and the associated DatabaseIO to the given State.
591    *
592    *  All transitions must begin from the 'STATE_CLOSED' state or be to
593    *  the 'STATE_CLOSED' state (There are no nested begin/end pairs at
594    *  this time.)
595    *
596    *  \param[in] new_state The new State to which the Region and DatabaseIO should be set.
597    *  \returns True if successful.
598    *
599    */
begin_mode(State new_state)600   bool Region::begin_mode(State new_state)
601   {
602     bool success = false;
603     {
604       IOSS_FUNC_ENTER(m_);
605       success = begin_mode__(new_state);
606     }
607     // Pass the 'begin state' message on to the database so it can do any
608     // cleanup/data checking/manipulations it needs to do.
609     if (success) {
610       DatabaseIO *db = get_database();
611 
612       if (new_state == STATE_DEFINE_TRANSIENT && db->usage() == Ioss::WRITE_HISTORY &&
613           !(is_input_or_appending_output(db))) {
614         set_state(STATE_CLOSED);
615         Ioss::Utils::generate_history_mesh(this);
616         set_state(new_state);
617       }
618       success = db->begin(new_state);
619     }
620     return success;
621   }
622 
begin_mode__(State new_state)623   bool Region::begin_mode__(State new_state)
624   {
625     bool success = false;
626     if (new_state == STATE_CLOSED) {
627       success = set_state(new_state);
628     }
629     else {
630       switch (get_state()) {
631       case STATE_CLOSED:
632         // Make sure we can go to the specified state.
633         switch (new_state) {
634         default: success = set_state(new_state);
635         }
636         break;
637 
638       // For the invalid transitions; provide a more meaningful
639       // message in certain cases...
640       case STATE_READONLY: {
641         std::ostringstream errmsg;
642         fmt::print(errmsg, "Cannot change state of an input (readonly) database in {}",
643                    get_database()->get_filename());
644         IOSS_ERROR(errmsg);
645       }
646 
647       break;
648       default: {
649         std::ostringstream errmsg;
650         fmt::print(errmsg, "Invalid nesting of begin/end pairs in {}",
651                    get_database()->get_filename());
652         IOSS_ERROR(errmsg);
653       }
654       }
655     }
656     return success;
657   }
658 
659   /** \brief Return the Region and the associated DatabaseIO to STATE_CLOSED.
660    *
661    *  \param[in] current_state The State to end.
662    *  \returns True if successful.
663    *
664    */
end_mode(State current_state)665   bool Region::end_mode(State current_state)
666   {
667     bool success = true;
668     {
669       IOSS_FUNC_ENTER(m_);
670       success = end_mode__(current_state);
671     }
672 
673     // Pass the 'end state' message on to the database so it can do any
674     // cleanup/data checking/manipulations it needs to do.
675     success = get_database()->end(current_state);
676     begin_mode(STATE_CLOSED);
677     return success;
678   }
679 
end_mode__(State current_state)680   bool Region::end_mode__(State current_state)
681   {
682     bool success = true;
683     // Check that 'current_state' matches the current state of the
684     // Region (that is, we are leaving the state we are in).
685     if (get_state() != current_state) {
686       std::ostringstream errmsg;
687       fmt::print(errmsg,
688                  "ERROR: Specified end state does not match currently open state\n"
689                  "       [{}]\n",
690                  get_database()->get_filename());
691       IOSS_ERROR(errmsg);
692     }
693 
694     if (current_state == STATE_DEFINE_MODEL) {
695       if (is_input_or_appending_output(get_database())) {
696         auto sortName = [](const Ioss::EntityBlock *b1, const Ioss::EntityBlock *b2) {
697           return (b1->name() < b2->name());
698         };
699         Ioss::sort(structuredBlocks.begin(), structuredBlocks.end(), sortName);
700       }
701       else {
702         // Sort the element blocks based on the idOffset field, followed by
703         // name...
704         auto lessOffset = [](const Ioss::EntityBlock *b1, const Ioss::EntityBlock *b2) {
705           SMART_ASSERT(b1->property_exists(orig_block_order()));
706           SMART_ASSERT(b2->property_exists(orig_block_order()));
707           int64_t b1_orderInt = b1->get_property(orig_block_order()).get_int();
708           int64_t b2_orderInt = b2->get_property(orig_block_order()).get_int();
709           return ((b1_orderInt == b2_orderInt) ? (b1->name() < b2->name())
710                                                : (b1_orderInt < b2_orderInt));
711         };
712 
713         Ioss::sort(elementBlocks.begin(), elementBlocks.end(), lessOffset);
714         Ioss::sort(faceBlocks.begin(), faceBlocks.end(), lessOffset);
715         Ioss::sort(edgeBlocks.begin(), edgeBlocks.end(), lessOffset);
716 
717         // Now update the block offsets based on this new order...
718         {
719           int64_t offset = 0;
720           for (auto eb : elementBlocks) {
721             eb->set_offset(offset);
722             offset += eb->entity_count();
723           }
724         }
725         {
726           int64_t offset = 0;
727           for (auto fb : faceBlocks) {
728             fb->set_offset(offset);
729             offset += fb->entity_count();
730           }
731         }
732         {
733           int64_t offset = 0;
734           for (auto eb : edgeBlocks) {
735             eb->set_offset(offset);
736             offset += eb->entity_count();
737           }
738         }
739       }
740 
741       // GroupingEntity consistency check:
742       // -- debug and parallel     -- default to true; can disable via environment variable
743       // -- non-debug and parallel -- default to false; can enable via environment variable
744 #ifndef NDEBUG
745       bool check_consistency = true;
746 #else
747       bool check_consistency = false;
748 #endif
749       Ioss::Utils::check_set_bool_property(get_database()->get_property_manager(),
750                                            "CHECK_PARALLEL_CONSISTENCY", check_consistency);
751       if (check_consistency) {
752         bool ok = check_parallel_consistency(*this);
753         if (!ok) {
754           std::ostringstream errmsg;
755           fmt::print(errmsg, "ERROR: Parallel Consistency Failure for {} database '{}'.",
756                      (get_database()->is_input() ? "input" : "output"),
757                      get_database()->get_filename());
758           IOSS_ERROR(errmsg);
759         }
760       }
761 
762       modelDefined = true;
763     }
764     else if (current_state == STATE_DEFINE_TRANSIENT) {
765       transientDefined = true;
766     }
767 
768     return success;
769   }
770 
771   /** \brief Add a state for a specified time.
772    *
773    *  The states in the region will be 1-based.
774    *
775    *  \param[in] time The time at the new state.
776    *  \returns The state index (1-based).
777    */
add_state__(double time)778   int Region::add_state__(double time)
779   {
780     static bool warning_output = false;
781 
782     // NOTE:  For restart input databases, it is possible that the time
783     //        is not monotonically increasing...
784     if (!get_database()->is_input() && !stateTimes.empty() && time <= stateTimes.back()) {
785       // Check that time is increasing...
786       if (!warning_output) {
787         fmt::print(Ioss::WARNING(),
788                    "Current time {} is not greater than previous time {} in\n\t{}.\n"
789                    "This may cause problems in applications that assume monotonically increasing "
790                    "time values.\n",
791                    time, stateTimes.back(), get_database()->get_filename());
792         warning_output = true;
793       }
794     }
795 
796     if (get_database()->is_input() || get_database()->usage() == WRITE_RESULTS ||
797         get_database()->usage() == WRITE_RESTART) {
798       stateTimes.push_back(time);
799       SMART_ASSERT((int)stateTimes.size() == stateCount + 1)(stateTimes.size())(stateCount);
800     }
801     else {
802 
803       // Keep only the last time in the vector... This is to avoid
804       // memory growth for output databases that write lots of steps
805       // (heartbeat, history).  There is no need to keep a list of
806       // times that have been written since they are just streamed out
807       // and never read We do sometimes need the list of times written
808       // to restart or results files though...
809       if (stateTimes.empty()) {
810         stateTimes.push_back(time);
811       }
812       else {
813         stateTimes[0] = time;
814       }
815     }
816     return ++stateCount;
817     ;
818   }
819 
820   /** \brief Get the time corresponding to the specified state or the currently active state.
821    *
822    *  \param[in] state The state index (1-based) or -1 for the currently active state.
823    *  \returns The time at the specified state or the currently active state.
824    */
get_state_time(int state)825   double Region::get_state_time(int state) const
826   {
827     IOSS_FUNC_ENTER(m_);
828     double time = 0.0;
829     if (state == -1) {
830       if (get_database()->is_input() || get_database()->usage() == WRITE_RESULTS ||
831           get_database()->usage() == WRITE_RESTART) {
832         if (currentState == -1) {
833           std::ostringstream errmsg;
834           fmt::print(errmsg, "ERROR: No currently active state.\n       [{}]\n",
835                      get_database()->get_filename());
836           IOSS_ERROR(errmsg);
837         }
838         else {
839           SMART_ASSERT((int)stateTimes.size() >= currentState)(stateTimes.size())(currentState);
840           time = stateTimes[currentState - 1];
841         }
842       }
843       else {
844         SMART_ASSERT(!stateTimes.empty());
845         time = stateTimes[0];
846       }
847     }
848     else if (state <= 0 || state > stateCount) {
849       std::ostringstream errmsg;
850       fmt::print(errmsg,
851                  "ERROR: Requested state ({}) is invalid. State must be between 1 and {}.\n"
852                  "       [{}]\n",
853                  state, stateCount, get_database()->get_filename());
854       IOSS_ERROR(errmsg);
855     }
856     else {
857       if (get_database()->is_input() || get_database()->usage() == WRITE_RESULTS ||
858           get_database()->usage() == WRITE_RESTART) {
859         SMART_ASSERT((int)stateTimes.size() >= state)(stateTimes.size())(state);
860         time = stateTimes[state - 1];
861       }
862       else {
863         SMART_ASSERT(!stateTimes.empty());
864         time = stateTimes[0];
865       }
866     }
867     return time;
868   }
869 
870   /** \brief Get the maximum time step index (1-based) and time for the region.
871    *
872    *  \returns A pair consisting of the step (1-based) corresponding to
873    *           the maximum time on the database and the corresponding maximum
874    *           time value. Note that this may not necessarily be the last step
875    *           on the database if cycle and overlay are being used.
876    */
get_max_time()877   std::pair<int, double> Region::get_max_time() const
878   {
879     IOSS_FUNC_ENTER(m_);
880     if (!get_database()->is_input() && get_database()->usage() != WRITE_RESULTS &&
881         get_database()->usage() != WRITE_RESTART) {
882       return std::make_pair(currentState, stateTimes[0]);
883     }
884     // Cleanout the stateTimes vector and reload with current data in
885     // case the database is being read and written at the same time.
886     // This is rare, but is a supported use case.
887     stateCount = 0;
888     Ioss::Utils::clear(stateTimes);
889     DatabaseIO *db = get_database();
890     db->get_step_times();
891 
892     int    step     = -1;
893     double max_time = -1.0;
894     for (int i = 0; i < static_cast<int>(stateTimes.size()); i++) {
895       if (stateTimes[i] > max_time) {
896         step     = i;
897         max_time = stateTimes[i];
898       }
899     }
900     return std::make_pair(step + 1, max_time);
901   }
902 
903   /** \brief Get the minimum time step index (1-based) and time for the region.
904    *
905    *  \returns A pair consisting of the step (1-based) corresponding to
906    *           the minimum time on the database and the corresponding minimum
907    *           time value. Note that this may not necessarily be the first step
908    *           on the database if cycle and overlay are being used.
909    */
get_min_time()910   std::pair<int, double> Region::get_min_time() const
911   {
912     IOSS_FUNC_ENTER(m_);
913     if (!get_database()->is_input() && get_database()->usage() != WRITE_RESULTS &&
914         get_database()->usage() != WRITE_RESTART) {
915       return std::make_pair(currentState, stateTimes[0]);
916     }
917     // Cleanout the stateTimes vector and reload with current data in
918     // case the database is being read and written at the same time.
919     // This is rare, but is a supported use case.
920     stateCount = 0;
921     Ioss::Utils::clear(stateTimes);
922     DatabaseIO *db = get_database();
923     db->get_step_times();
924 
925     int    step     = 0;
926     double min_time = stateTimes[0];
927     for (int i = 1; i < static_cast<int>(stateTimes.size()); i++) {
928       if (stateTimes[i] < min_time) {
929         step     = i;
930         min_time = stateTimes[i];
931       }
932     }
933     return std::make_pair(step + 1, min_time);
934   }
935 
936   /** \brief Begin a state (moment in time).
937    *
938    *  \param[in] state The state index (1-based).
939    *  \returns The time of this state.
940    */
begin_state(int state)941   double Region::begin_state(int state)
942   {
943     double time = 0.0;
944     if (get_database()->is_input() && stateCount == 0) {
945       std::ostringstream errmsg;
946       fmt::print(errmsg,
947                  "ERROR: There are no states (time steps) on the input database.\n"
948                  "       [{}]\n",
949                  get_database()->get_filename());
950       IOSS_ERROR(errmsg);
951     }
952     if (state <= 0 || state > stateCount) {
953       std::ostringstream errmsg;
954       fmt::print(errmsg,
955                  "ERROR: Requested state ({}) is invalid. State must be between 1 and {}.\n"
956                  "       [{}]\n",
957                  state, stateCount, get_database()->get_filename());
958       IOSS_ERROR(errmsg);
959     }
960     else if (currentState != -1 && !get_database()->is_input()) {
961       std::ostringstream errmsg;
962       fmt::print(errmsg, "ERROR: State {} was not ended. Can not begin new state.\n       [{}]\n",
963                  currentState, get_database()->get_filename());
964       IOSS_ERROR(errmsg);
965     }
966     else {
967       {
968         IOSS_FUNC_ENTER(m_);
969 
970         SMART_ASSERT(state <= stateCount)(state)(stateCount);
971         if (get_database()->is_input() || get_database()->usage() == WRITE_RESULTS ||
972             get_database()->usage() == WRITE_RESTART) {
973           SMART_ASSERT((int)stateTimes.size() >= state)(stateTimes.size())(state);
974           time = stateTimes[state - 1];
975         }
976         else {
977           SMART_ASSERT(!stateTimes.empty());
978           time = stateTimes[0];
979         }
980         currentState = state;
981       }
982       DatabaseIO *db = get_database();
983       db->begin_state(state, time);
984     }
985     return time;
986   }
987 
988   /** \brief End a state (moment in time).
989    *
990    *  \param[in] state The state index (1-based).
991    *  \returns The time of this state.
992    */
end_state(int state)993   double Region::end_state(int state)
994   {
995     if (state != currentState) {
996       std::ostringstream errmsg;
997       fmt::print(errmsg,
998                  "ERROR: The current database state ({}) does not match the ending state ({}).\n"
999                  "       [{}]\n",
1000                  currentState, state, get_database()->get_filename());
1001       IOSS_ERROR(errmsg);
1002     }
1003     DatabaseIO *db   = get_database();
1004     double      time = 0.0;
1005     {
1006       IOSS_FUNC_ENTER(m_);
1007       if (get_database()->is_input() || get_database()->usage() == WRITE_RESULTS ||
1008           get_database()->usage() == WRITE_RESTART) {
1009         SMART_ASSERT((int)stateTimes.size() >= state)(stateTimes.size())(state);
1010         time = stateTimes[state - 1];
1011       }
1012       else {
1013         SMART_ASSERT(!stateTimes.empty());
1014         time = stateTimes[0];
1015       }
1016     }
1017     db->end_state(state, time);
1018     currentState = -1;
1019     return time;
1020   }
1021 
1022   /** \brief Add a structured block to the region.
1023    *
1024    *  \param[in] structured_block The structured block to add
1025    *  \returns True if successful.
1026    */
add(StructuredBlock * structured_block)1027   bool Region::add(StructuredBlock *structured_block)
1028   {
1029     check_for_duplicate_names(this, structured_block);
1030     update_database(this, structured_block);
1031     IOSS_FUNC_ENTER(m_);
1032 
1033     // Check that region is in correct state for adding entities
1034     if (get_state() == STATE_DEFINE_MODEL) {
1035       // Add node and cell offsets based on the node_count and
1036       // cell_count of the previous block.  Only add if there is more
1037       // than one block; use default for first block (or user-defined
1038       // values)
1039       if (!structuredBlocks.empty()) {
1040         auto   prev_block = structuredBlocks.back();
1041         size_t num_node   = prev_block->get_property("node_count").get_int();
1042         size_t num_cell   = prev_block->get_property("cell_count").get_int();
1043         num_node += prev_block->get_node_offset();
1044         num_cell += prev_block->get_cell_offset();
1045 
1046         structured_block->set_node_offset(num_node);
1047         structured_block->set_cell_offset(num_cell);
1048 
1049         size_t global_num_node = prev_block->get_property("global_node_count").get_int();
1050         size_t global_num_cell = prev_block->get_property("global_cell_count").get_int();
1051         global_num_node += prev_block->get_node_global_offset();
1052         global_num_cell += prev_block->get_cell_global_offset();
1053 
1054         structured_block->set_node_global_offset(global_num_node);
1055         structured_block->set_cell_global_offset(global_num_cell);
1056       }
1057 
1058       structured_block->property_add(
1059           Ioss::Property(orig_block_order(), (int)structuredBlocks.size()));
1060       structuredBlocks.push_back(structured_block);
1061 
1062       // This will possibly be overwritten at a later time when the block is output
1063       // to the cgns file
1064       structured_block->property_add(Ioss::Property("zone", (int)structuredBlocks.size()));
1065       structured_block->property_add(Ioss::Property("base", 1));
1066       // Add name as alias to itself to simplify later uses...
1067       add_alias__(structured_block);
1068       return true;
1069     }
1070     return false;
1071   }
1072 
1073   /** \brief Add a node block to the region.
1074    *
1075    *  \param[in] node_block The node block to add
1076    *  \returns True if successful.
1077    */
add(NodeBlock * node_block)1078   bool Region::add(NodeBlock *node_block)
1079   {
1080     check_for_duplicate_names(this, node_block);
1081     update_database(this, node_block);
1082     IOSS_FUNC_ENTER(m_);
1083 
1084     // Check that region is in correct state for adding entities
1085     if (get_state() == STATE_DEFINE_MODEL) {
1086       nodeBlocks.push_back(node_block);
1087       // Add name as alias to itself to simplify later uses...
1088       add_alias__(node_block);
1089 
1090       return true;
1091     }
1092     return false;
1093   }
1094 
1095   /** \brief Remove an assembly to the region.
1096    *
1097    *  \param[in] removal The assembly to remove
1098    *  \returns True if successful.
1099    *
1100    *  Checks other assemblies for uses of this assembly and removes
1101    * if from their member lists.
1102    */
remove(Assembly * removal)1103   bool Region::remove(Assembly *removal)
1104   {
1105     IOSS_FUNC_ENTER(m_);
1106 
1107     bool changed = false;
1108     // Check that region is in correct state for modifying entities
1109     if (get_state() == STATE_DEFINE_MODEL) {
1110       // Check all existing assemblies to see if any contain this assembly:
1111       for (auto *assembly : assemblies) {
1112         bool removed = assembly->remove(removal);
1113         if (removed) {
1114           changed = true;
1115         }
1116       }
1117 
1118       // Now remove this assembly...
1119       for (size_t i = 0; i < assemblies.size(); i++) {
1120         if (assemblies[i] == removal) {
1121           assemblies.erase(assemblies.begin() + i);
1122           changed = true;
1123         }
1124       }
1125     }
1126     return changed;
1127   }
1128 
1129   /** \brief Add an assembly to the region.
1130    *
1131    *  \param[in] assembly The assembly to add
1132    *  \returns True if successful.
1133    */
add(Assembly * assembly)1134   bool Region::add(Assembly *assembly)
1135   {
1136     check_for_duplicate_names(this, assembly);
1137     update_database(this, assembly);
1138     IOSS_FUNC_ENTER(m_);
1139 
1140     // Check that region is in correct state for adding entities
1141     if (get_state() == STATE_DEFINE_MODEL) {
1142       assemblies.push_back(assembly);
1143       // Add name as alias to itself to simplify later uses...
1144       add_alias__(assembly);
1145 
1146       return true;
1147     }
1148     return false;
1149   }
1150 
1151   /** \brief Add an blob to the region.
1152    *
1153    *  \param[in] blob The blob to add
1154    *  \returns True if successful.
1155    */
add(Blob * blob)1156   bool Region::add(Blob *blob)
1157   {
1158     check_for_duplicate_names(this, blob);
1159     update_database(this, blob);
1160     IOSS_FUNC_ENTER(m_);
1161 
1162     // Check that region is in correct state for adding entities
1163     if (get_state() == STATE_DEFINE_MODEL) {
1164       blobs.push_back(blob);
1165       // Add name as alias to itself to simplify later uses...
1166       add_alias__(blob);
1167 
1168       return true;
1169     }
1170     return false;
1171   }
1172 
1173   /** \brief Add a coordinate frame to the region.
1174    *
1175    *  \param[in] frame The coordinate frame to add
1176    *  \returns True if successful.
1177    */
add(const CoordinateFrame & frame)1178   bool Region::add(const CoordinateFrame &frame)
1179   {
1180     IOSS_FUNC_ENTER(m_);
1181     // Check that region is in correct state for adding entities
1182     if (get_state() == STATE_DEFINE_MODEL) {
1183       coordinateFrames.push_back(frame);
1184       return true;
1185     }
1186     return false;
1187   }
1188 
1189   /** \brief Add an element block to the region.
1190    *
1191    *  \param[in] element_block The element block to add
1192    *  \returns True if successful.
1193    */
add(ElementBlock * element_block)1194   bool Region::add(ElementBlock *element_block)
1195   {
1196     check_for_duplicate_names(this, element_block);
1197     update_database(this, element_block);
1198     IOSS_FUNC_ENTER(m_);
1199 
1200     // Check that region is in correct state for adding entities
1201     if (get_state() == STATE_DEFINE_MODEL) {
1202       // Add name as alias to itself to simplify later uses...
1203       add_alias__(element_block);
1204 
1205       // An input database defines these in the order matching the order
1206       // on the "file".  For output, we need to order based on the
1207       // "original_block_order" property and calculate the offset at that
1208       // point.  This is done in "end".
1209 
1210       if (is_input_or_appending_output(get_database())) {
1211         size_t  nblocks = elementBlocks.size();
1212         int64_t offset  = 0;
1213         if (nblocks > 0) {
1214           offset =
1215               elementBlocks[nblocks - 1]->get_offset() + elementBlocks[nblocks - 1]->entity_count();
1216         }
1217         SMART_ASSERT(offset >= 0)(offset);
1218         element_block->set_offset(offset);
1219       }
1220 #if 0
1221       // Would like to use this, but gives issue in legacy contact...
1222       // If this is enabled, then remove all settings of
1223       // "orig_block_order()" from individual DatabaseIO classes.
1224       element_block->property_add(Ioss::Property(orig_block_order(), (int)elementBlocks.size()));
1225 #else
1226       else {
1227         // Check whether the "original_block_order" property exists on
1228         // this element block. If it isn't there, then add it with a
1229         // large value. If this is an element block read from the
1230         // input mesh, then the value will be updated during the
1231         // 'synchronize_id_and_name' function; if it is a block
1232         // created by the application during execution, then this
1233         // value will persist.  Add the property with a very large
1234         // number such that it will later be sorted after all
1235         // "original" blocks.  Note that it doesn't matter if two of
1236         // the "new" blocks have the same value since there is no
1237         // ordering of new blocks that must be preserved. (Use
1238         // int_MAX/2 just to avoid some paranoia about strange issue
1239         // that might arise from int_MAX)
1240         if (!element_block->property_exists(orig_block_order())) {
1241           element_block->property_add(Property(orig_block_order(), INT_MAX / 2));
1242         }
1243       }
1244 #endif
1245       elementBlocks.push_back(element_block);
1246       return true;
1247     }
1248     return false;
1249   }
1250 
1251   /** \brief Add a face block to the region.
1252    *
1253    *  \param[in] face_block The face block to add
1254    *  \returns True if successful.
1255    */
add(FaceBlock * face_block)1256   bool Region::add(FaceBlock *face_block)
1257   {
1258     check_for_duplicate_names(this, face_block);
1259     update_database(this, face_block);
1260     IOSS_FUNC_ENTER(m_);
1261 
1262     // Check that region is in correct state for adding entities
1263     if (get_state() == STATE_DEFINE_MODEL) {
1264       // Add name as alias to itself to simplify later uses...
1265       add_alias__(face_block);
1266 
1267       // An input database defines these in the order matching the order
1268       // on the "file".  For output, we need to order based on the
1269       // "original_block_order" property and calculate the offset at that
1270       // point.  This is done in "end".
1271 
1272       if (is_input_or_appending_output(get_database())) {
1273         size_t  nblocks = faceBlocks.size();
1274         int64_t offset  = 0;
1275         if (nblocks > 0) {
1276           offset = faceBlocks[nblocks - 1]->get_offset() + faceBlocks[nblocks - 1]->entity_count();
1277         }
1278         face_block->set_offset(offset);
1279       }
1280       face_block->property_add(Ioss::Property(orig_block_order(), (int)faceBlocks.size()));
1281       faceBlocks.push_back(face_block);
1282       return true;
1283     }
1284     return false;
1285   }
1286 
1287   /** \brief Add an edge block to the region.
1288    *
1289    *  \param[in] edge_block The edge block to add
1290    *  \returns True if successful.
1291    */
add(EdgeBlock * edge_block)1292   bool Region::add(EdgeBlock *edge_block)
1293   {
1294     check_for_duplicate_names(this, edge_block);
1295     update_database(this, edge_block);
1296     IOSS_FUNC_ENTER(m_);
1297 
1298     // Check that region is in correct state for adding entities
1299     if (get_state() == STATE_DEFINE_MODEL) {
1300       // Add name as alias to itself to simplify later uses...
1301       add_alias__(edge_block);
1302 
1303       // An input database defines these in the order matching the order
1304       // on the "file".  For output, we need to order based on the
1305       // "original_block_order" property and calculate the offset at that
1306       // point.  This is done in "end".
1307 
1308       if (is_input_or_appending_output(get_database())) {
1309         size_t  nblocks = edgeBlocks.size();
1310         int64_t offset  = 0;
1311         if (nblocks > 0) {
1312           offset = edgeBlocks[nblocks - 1]->get_offset() + edgeBlocks[nblocks - 1]->entity_count();
1313         }
1314         edge_block->set_offset(offset);
1315       }
1316       edge_block->property_add(Ioss::Property(orig_block_order(), (int)edgeBlocks.size()));
1317       edgeBlocks.push_back(edge_block);
1318       return true;
1319     }
1320     return false;
1321   }
1322 
1323   /** \brief Add a side set to the region.
1324    *
1325    *  \param[in] sideset The side set to add
1326    *  \returns True if successful.
1327    */
add(SideSet * sideset)1328   bool Region::add(SideSet *sideset)
1329   {
1330     check_for_duplicate_names(this, sideset);
1331     update_database(this, sideset);
1332     IOSS_FUNC_ENTER(m_);
1333     // Check that region is in correct state for adding entities
1334     if (get_state() == STATE_DEFINE_MODEL) {
1335       // Add name as alias to itself to simplify later uses...
1336       add_alias__(sideset);
1337       sideSets.push_back(sideset);
1338       return true;
1339     }
1340     return false;
1341   }
1342 
1343   /** \brief Add a node set to the region.
1344    *
1345    *  \param[in] nodeset The node set to add
1346    *  \returns True if successful.
1347    */
add(NodeSet * nodeset)1348   bool Region::add(NodeSet *nodeset)
1349   {
1350     check_for_duplicate_names(this, nodeset);
1351     update_database(this, nodeset);
1352     IOSS_FUNC_ENTER(m_);
1353     // Check that region is in correct state for adding entities
1354     if (get_state() == STATE_DEFINE_MODEL) {
1355       // Add name as alias to itself to simplify later uses...
1356       add_alias__(nodeset);
1357       nodeSets.push_back(nodeset);
1358       return true;
1359     }
1360     return false;
1361   }
1362 
1363   /** \brief Add an edge set to the region.
1364    *
1365    *  \param[in] edgeset The edge set to add
1366    *  \returns True if successful.
1367    */
add(EdgeSet * edgeset)1368   bool Region::add(EdgeSet *edgeset)
1369   {
1370     check_for_duplicate_names(this, edgeset);
1371     update_database(this, edgeset);
1372     IOSS_FUNC_ENTER(m_);
1373     // Check that region is in correct state for adding entities
1374     if (get_state() == STATE_DEFINE_MODEL) {
1375       // Add name as alias to itself to simplify later uses...
1376       add_alias__(edgeset);
1377       edgeSets.push_back(edgeset);
1378       return true;
1379     }
1380     return false;
1381   }
1382 
1383   /** \brief Add a face set to the region.
1384    *
1385    *  \param[in] faceset The face set to add
1386    *  \returns True if successful.
1387    */
add(FaceSet * faceset)1388   bool Region::add(FaceSet *faceset)
1389   {
1390     check_for_duplicate_names(this, faceset);
1391     update_database(this, faceset);
1392     IOSS_FUNC_ENTER(m_);
1393     // Check that region is in correct state for adding entities
1394     if (get_state() == STATE_DEFINE_MODEL) {
1395       // Add name as alias to itself to simplify later uses...
1396       add_alias__(faceset);
1397       faceSets.push_back(faceset);
1398       return true;
1399     }
1400     return false;
1401   }
1402 
1403   /** \brief Add an element set to the region.
1404    *
1405    *  \param[in] elementset The element set to add
1406    *  \returns True if successful.
1407    */
add(ElementSet * elementset)1408   bool Region::add(ElementSet *elementset)
1409   {
1410     check_for_duplicate_names(this, elementset);
1411     update_database(this, elementset);
1412     IOSS_FUNC_ENTER(m_);
1413     // Check that region is in correct state for adding entities
1414     if (get_state() == STATE_DEFINE_MODEL) {
1415       // Add name as alias to itself to simplify later uses...
1416       add_alias__(elementset);
1417       elementSets.push_back(elementset);
1418       return true;
1419     }
1420     return false;
1421   }
1422 
1423   /** \brief Add a comm set to the region.
1424    *
1425    *  \param[in] commset The comm set to add
1426    *  \returns True if successful.
1427    */
add(CommSet * commset)1428   bool Region::add(CommSet *commset)
1429   {
1430     check_for_duplicate_names(this, commset);
1431     update_database(this, commset);
1432     IOSS_FUNC_ENTER(m_);
1433     // Check that region is in correct state for adding entities
1434     if (get_state() == STATE_DEFINE_MODEL) {
1435       // Add name as alias to itself to simplify later uses...
1436       add_alias__(commset);
1437       commSets.push_back(commset);
1438       return true;
1439     }
1440     return false;
1441   }
1442 
1443   /** \brief Get all the region's Assembly objects.
1444    *
1445    *  \returns A vector of all the region's Assembly objects.
1446    */
get_assemblies()1447   const AssemblyContainer &Region::get_assemblies() const { return assemblies; }
1448 
1449   /** \brief Get all the region's Blob objects.
1450    *
1451    *  \returns A vector of all the region's Blob objects.
1452    */
get_blobs()1453   const BlobContainer &Region::get_blobs() const { return blobs; }
1454 
1455   /** \brief Get all the region's NodeBlock objects.
1456    *
1457    *  \returns A vector of all the region's NodeBlock objects.
1458    */
get_node_blocks()1459   const NodeBlockContainer &Region::get_node_blocks() const { return nodeBlocks; }
1460 
1461   /** \brief Get all the region's EdgeBlock objects.
1462    *
1463    *  \returns A vector of all the region's EdgeBlock objects.
1464    */
get_edge_blocks()1465   const EdgeBlockContainer &Region::get_edge_blocks() const { return edgeBlocks; }
1466 
1467   /** \brief Get all the region's FaceBlock objects.
1468    *
1469    *  \returns A vector of all the region's FaceBlock objects.
1470    */
get_face_blocks()1471   const FaceBlockContainer &Region::get_face_blocks() const { return faceBlocks; }
1472 
1473   /** \brief Get all the region's ElementBlock objects.
1474    *
1475    *  \returns A vector of all the region's ElementBlock objects.
1476    */
get_element_blocks()1477   const ElementBlockContainer &Region::get_element_blocks() const { return elementBlocks; }
1478 
1479   /** \brief Get all the region's StructuredBlock objects.
1480    *
1481    *  \returns A vector of all the region's StructuredBlock objects.
1482    */
get_structured_blocks()1483   const StructuredBlockContainer &Region::get_structured_blocks() const { return structuredBlocks; }
1484 
1485   /** \brief Get all the region's SideSet objects.
1486    *
1487    *  \returns A vector of all the region's SideSet objects.
1488    */
get_sidesets()1489   const SideSetContainer &Region::get_sidesets() const { return sideSets; }
1490 
1491   /** \brief Get all the region's NodeSet objects.
1492    *
1493    *  \returns A vector of all the region's NodeSet objects.
1494    */
get_nodesets()1495   const NodeSetContainer &Region::get_nodesets() const { return nodeSets; }
1496 
1497   /** \brief Get all the region's EdgeSet objects.
1498    *
1499    *  \returns A vector of all the region's EdgeSet objects.
1500    */
get_edgesets()1501   const EdgeSetContainer &Region::get_edgesets() const { return edgeSets; }
1502 
1503   /** \brief Get all the region's FaceSet objects.
1504    *
1505    *  \returns A vector of all the region's FaceSet objects.
1506    */
get_facesets()1507   const FaceSetContainer &Region::get_facesets() const { return faceSets; }
1508 
1509   /** \brief Get all the region's ElementSet objects.
1510    *
1511    *  \returns A vector of all the region's ElementSet objects.
1512    */
get_elementsets()1513   const ElementSetContainer &Region::get_elementsets() const { return elementSets; }
1514 
1515   /** \brief Get all the region's CommSet objects.
1516    *
1517    *  \returns A vector of all the region's CommSet objects.
1518    */
get_commsets()1519   const CommSetContainer &Region::get_commsets() const { return commSets; }
1520 
1521   /** \brief Get all the region's CoordinateFrame objects.
1522    *
1523    *  \returns A vector of all the region's CoordinateFrame objects.
1524    */
get_coordinate_frames()1525   const CoordinateFrameContainer &Region::get_coordinate_frames() const { return coordinateFrames; }
1526 
1527   /** \brief Add a grouping entity's name as an alias for itself.
1528    *
1529    *  \param[in] ge The grouping entity.
1530    *  \returns True if successful
1531    */
add_alias(const GroupingEntity * ge)1532   bool Region::add_alias(const GroupingEntity *ge)
1533   {
1534     IOSS_FUNC_ENTER(m_);
1535     return add_alias__(ge);
1536   }
1537 
add_alias__(const GroupingEntity * ge)1538   bool Region::add_alias__(const GroupingEntity *ge)
1539   {
1540     // See if an entity with this name already exists...
1541     const std::string &db_name = ge->name();
1542     std::string        alias   = get_alias__(db_name);
1543 
1544     if (!alias.empty()) {
1545       const GroupingEntity *old_ge = get_entity(db_name);
1546       if (old_ge != nullptr && ge != old_ge) {
1547         if (!((old_ge->type() == SIDEBLOCK && ge->type() == SIDESET) ||
1548               (ge->type() == SIDEBLOCK && old_ge->type() == SIDESET))) {
1549           ssize_t            old_id = old_ge->get_optional_property(id_str(), -1);
1550           ssize_t            new_id = ge->get_optional_property(id_str(), -1);
1551           std::ostringstream errmsg;
1552           fmt::print(errmsg,
1553                      "\n\nERROR: Duplicate names detected.\n"
1554                      "       The name '{}' was found for both {} {} and {} {}.\n"
1555                      "       Names must be unique over all types in a finite element model.\n\n",
1556                      db_name, old_ge->type_string(), old_id, ge->type_string(), new_id);
1557           IOSS_ERROR(errmsg);
1558         }
1559       }
1560     }
1561     bool success = add_alias__(db_name, db_name);
1562 
1563     // "db_name" property is used with the canonical name setting.
1564     if (success && ge->property_exists("db_name")) {
1565       std::string canon_name = ge->get_property("db_name").get_string();
1566       if (canon_name != db_name) {
1567         success = add_alias__(db_name, canon_name);
1568       }
1569     }
1570 
1571     return success;
1572   }
1573 
1574   /** \brief Add an alias for a name in a region.
1575    *
1576    *  For use with the USTRING type in Sierra, create an uppercase
1577    *  version of all aliases...
1578    *
1579    *  \param[in] db_name The original name.
1580    *  \param[in] alias the alias
1581    *  \returns True if successful
1582    */
add_alias(const std::string & db_name,const std::string & alias)1583   bool Region::add_alias(const std::string &db_name, const std::string &alias)
1584   {
1585     IOSS_FUNC_ENTER(m_);
1586     return add_alias__(db_name, alias);
1587   }
1588 
add_alias__(const std::string & db_name,const std::string & alias)1589   bool Region::add_alias__(const std::string &db_name, const std::string &alias)
1590   {
1591     // Possible that 'db_name' is itself an alias, resolve down to "canonical"
1592     // name...
1593     std::string canon = db_name;
1594     if (db_name != alias) {
1595       canon = get_alias__(db_name);
1596     }
1597 
1598     if (!canon.empty()) {
1599       std::string uname = Ioss::Utils::uppercase(alias);
1600       if (uname != alias) {
1601         aliases_.insert(std::make_pair(uname, canon));
1602       }
1603 
1604       bool result;
1605       std::tie(std::ignore, result) = aliases_.insert(std::make_pair(alias, canon));
1606       return result;
1607     }
1608     std::ostringstream errmsg;
1609     fmt::print(errmsg,
1610                "\n\nERROR: The entity named '{}' which is being aliased to '{}' does not exist in "
1611                "region '{}'.\n",
1612                db_name, alias, name());
1613     IOSS_ERROR(errmsg);
1614   }
1615 
1616   /** \brief Get the original name for an alias.
1617    *
1618    *  \param[in] alias The alias name.
1619    *  \returns The original name.
1620    */
get_alias(const std::string & alias)1621   std::string Region::get_alias(const std::string &alias) const
1622   {
1623     IOSS_FUNC_ENTER(m_);
1624     return get_alias__(alias);
1625   }
1626 
get_alias__(const std::string & alias)1627   std::string Region::get_alias__(const std::string &alias) const
1628   {
1629     std::string ci_alias = Ioss::Utils::uppercase(alias);
1630     auto        I        = aliases_.find(ci_alias);
1631     if (I == aliases_.end()) {
1632       return "";
1633     }
1634     return (*I).second;
1635   }
1636 
1637   /** \brief Get all aliases for a name in the region.
1638    *
1639    *  \param[in] my_name The original name.
1640    *  \param[in,out] aliases On input, any vector of strings.
1641    *                         On output, all aliases for my_name are appended.
1642    *  \returns The number of aliases that were appended.
1643    *
1644    */
get_aliases(const std::string & my_name,std::vector<std::string> & aliases)1645   int Region::get_aliases(const std::string &my_name, std::vector<std::string> &aliases) const
1646   {
1647     IOSS_FUNC_ENTER(m_);
1648     size_t size = aliases.size();
1649     for (const auto &alias_pair : aliases_) {
1650       std::string alias = alias_pair.first;
1651       std::string base  = alias_pair.second;
1652       if (base == my_name) {
1653         aliases.push_back(alias);
1654       }
1655     }
1656     return static_cast<int>(aliases.size() - size);
1657   }
1658 
1659   /** \brief Get all original name / alias pairs for the region.
1660    *
1661    *  \returns All original name / alias pairs for the region.
1662    */
get_alias_map()1663   const AliasMap &Region::get_alias_map() const { return aliases_; }
1664 
1665   /** \brief Get an entity of a known EntityType
1666    *
1667    *  \param[in] my_name The name of the entity to get
1668    *  \param[in] io_type The known type of the entity.
1669    *  \returns The entity with the given name of the given type, or nullptr if not found.
1670    */
get_entity(const std::string & my_name,EntityType io_type)1671   GroupingEntity *Region::get_entity(const std::string &my_name, EntityType io_type) const
1672   {
1673     if (io_type == NODEBLOCK) {
1674       return get_node_block(my_name);
1675     }
1676     if (io_type == ELEMENTBLOCK) {
1677       return get_element_block(my_name);
1678     }
1679     if (io_type == STRUCTUREDBLOCK) {
1680       return get_structured_block(my_name);
1681     }
1682     if (io_type == FACEBLOCK) {
1683       return get_face_block(my_name);
1684     }
1685     if (io_type == EDGEBLOCK) {
1686       return get_edge_block(my_name);
1687     }
1688     if (io_type == SIDESET) {
1689       return get_sideset(my_name);
1690     }
1691     if (io_type == NODESET) {
1692       return get_nodeset(my_name);
1693     }
1694     else if (io_type == EDGESET) {
1695       return get_edgeset(my_name);
1696     }
1697     else if (io_type == FACESET) {
1698       return get_faceset(my_name);
1699     }
1700     else if (io_type == ELEMENTSET) {
1701       return get_elementset(my_name);
1702     }
1703     else if (io_type == COMMSET) {
1704       return get_commset(my_name);
1705     }
1706     else if (io_type == SIDEBLOCK) {
1707       return get_sideblock(my_name);
1708     }
1709     else if (io_type == ASSEMBLY) {
1710       return get_assembly(my_name);
1711     }
1712     else if (io_type == BLOB) {
1713       return get_blob(my_name);
1714     }
1715     return nullptr;
1716   }
1717 
1718   /** \brief Get an entity of a unknown EntityType
1719    *
1720    *  Searches for an entity with the given name in a fixed order of
1721    *  entity types. First NODEBLOCK entities are searched. Then ELEMENTBLOCK
1722    *  entities, etc.
1723    *
1724    *  \param[in] my_name The name of the entity to get
1725    *  \returns The entity with the given name, or nullptr if not found.
1726    */
get_entity(const std::string & my_name)1727   GroupingEntity *Region::get_entity(const std::string &my_name) const
1728   {
1729     GroupingEntity *entity = get_node_block(my_name);
1730     if (entity != nullptr) {
1731       return entity;
1732     }
1733     entity = get_element_block(my_name);
1734     if (entity != nullptr) {
1735       return entity;
1736     }
1737     entity = get_structured_block(my_name);
1738     if (entity != nullptr) {
1739       return entity;
1740     }
1741     entity = get_face_block(my_name);
1742     if (entity != nullptr) {
1743       return entity;
1744     }
1745     entity = get_edge_block(my_name);
1746     if (entity != nullptr) {
1747       return entity;
1748     }
1749     entity = get_sideset(my_name);
1750     if (entity != nullptr) {
1751       return entity;
1752     }
1753     entity = get_nodeset(my_name);
1754     if (entity != nullptr) {
1755       return entity;
1756     }
1757     entity = get_edgeset(my_name);
1758     if (entity != nullptr) {
1759       return entity;
1760     }
1761     entity = get_faceset(my_name);
1762     if (entity != nullptr) {
1763       return entity;
1764     }
1765     entity = get_elementset(my_name);
1766     if (entity != nullptr) {
1767       return entity;
1768     }
1769     entity = get_commset(my_name);
1770     if (entity != nullptr) {
1771       return entity;
1772     }
1773     entity = get_sideblock(my_name);
1774     if (entity != nullptr) {
1775       return entity;
1776     }
1777     entity = get_assembly(my_name);
1778     if (entity != nullptr) {
1779       return entity;
1780     }
1781     entity = get_blob(my_name);
1782     if (entity != nullptr) {
1783       return entity;
1784     }
1785 
1786     return entity;
1787   }
1788 
1789   /** \brief Get an entity of a known EntityType and specified id
1790    *
1791    *  \param[in] id The id of the entity to get
1792    *  \param[in] io_type The known type of the entity.
1793    *  \returns The entity with the given id of the given type, or nullptr if not found or if ids not
1794    * supported.
1795    */
get_entity(int64_t id,EntityType io_type)1796   GroupingEntity *Region::get_entity(int64_t id, EntityType io_type) const
1797   {
1798     if (io_type == NODEBLOCK) {
1799       return get_entity_internal(id, get_node_blocks());
1800     }
1801     if (io_type == ELEMENTBLOCK) {
1802       return get_entity_internal(id, get_element_blocks());
1803     }
1804     if (io_type == STRUCTUREDBLOCK) {
1805       return get_entity_internal(id, get_structured_blocks());
1806     }
1807     if (io_type == FACEBLOCK) {
1808       return get_entity_internal(id, get_face_blocks());
1809     }
1810     if (io_type == EDGEBLOCK) {
1811       return get_entity_internal(id, get_edge_blocks());
1812     }
1813     if (io_type == SIDESET) {
1814       return get_entity_internal(id, get_sidesets());
1815     }
1816     if (io_type == NODESET) {
1817       return get_entity_internal(id, get_nodesets());
1818     }
1819     else if (io_type == EDGESET) {
1820       return get_entity_internal(id, get_edgesets());
1821     }
1822     else if (io_type == FACESET) {
1823       return get_entity_internal(id, get_facesets());
1824     }
1825     else if (io_type == ELEMENTSET) {
1826       return get_entity_internal(id, get_elementsets());
1827     }
1828     else if (io_type == COMMSET) {
1829       return get_entity_internal(id, get_commsets());
1830     }
1831     else if (io_type == ASSEMBLY) {
1832       return get_entity_internal(id, get_assemblies());
1833     }
1834     else if (io_type == BLOB) {
1835       return get_entity_internal(id, get_blobs());
1836     }
1837     return nullptr;
1838   }
1839 
1840   /** \brief Get the assembly with the given name.
1841    *
1842    *  \param[in] my_name The name of the assembly to get.
1843    *  \returns The assembly, or nullptr if not found.
1844    */
get_assembly(const std::string & my_name)1845   Assembly *Region::get_assembly(const std::string &my_name) const
1846   {
1847     IOSS_FUNC_ENTER(m_);
1848     const std::string db_name = get_alias__(my_name);
1849     unsigned int      db_hash = Ioss::Utils::hash(db_name);
1850 
1851     Assembly *ge = nullptr;
1852     for (auto as : assemblies) {
1853       if (db_hash == as->hash() && as->name() == db_name) {
1854         ge = as;
1855         break;
1856       }
1857     }
1858     return ge;
1859   }
1860 
1861   /** \brief Get the blob with the given name.
1862    *
1863    *  \param[in] my_name The name of the blob to get.
1864    *  \returns The blob, or nullptr if not found.
1865    */
get_blob(const std::string & my_name)1866   Blob *Region::get_blob(const std::string &my_name) const
1867   {
1868     IOSS_FUNC_ENTER(m_);
1869     const std::string db_name = get_alias__(my_name);
1870     unsigned int      db_hash = Ioss::Utils::hash(db_name);
1871 
1872     Blob *ge = nullptr;
1873     for (auto bl : blobs) {
1874       if (db_hash == bl->hash() && bl->name() == db_name) {
1875         ge = bl;
1876         break;
1877       }
1878     }
1879     return ge;
1880   }
1881 
1882   /** \brief Get the node block with the given name.
1883    *
1884    *  \param[in] my_name The name of the node block to get.
1885    *  \returns The node block, or nullptr if not found.
1886    */
get_node_block(const std::string & my_name)1887   NodeBlock *Region::get_node_block(const std::string &my_name) const
1888   {
1889     IOSS_FUNC_ENTER(m_);
1890     const std::string db_name = get_alias__(my_name);
1891     unsigned int      db_hash = Ioss::Utils::hash(db_name);
1892 
1893     NodeBlock *ge = nullptr;
1894     for (auto nb : nodeBlocks) {
1895       if (db_hash == nb->hash() && nb->name() == db_name) {
1896         ge = nb;
1897         break;
1898       }
1899     }
1900     return ge;
1901   }
1902 
1903   /** \brief Get the edge block with the given name.
1904    *
1905    *  \param[in] my_name The name of the edge block to get.
1906    *  \returns The edge block, or nullptr if not found.
1907    */
get_edge_block(const std::string & my_name)1908   EdgeBlock *Region::get_edge_block(const std::string &my_name) const
1909   {
1910     IOSS_FUNC_ENTER(m_);
1911     const std::string db_name = get_alias__(my_name);
1912     unsigned int      db_hash = Ioss::Utils::hash(db_name);
1913 
1914     EdgeBlock *ge = nullptr;
1915     for (auto eb : edgeBlocks) {
1916       if (db_hash == eb->hash() && eb->name() == db_name) {
1917         ge = eb;
1918         break;
1919       }
1920     }
1921     return ge;
1922   }
1923 
1924   /** \brief Get the face block with the given name.
1925    *
1926    *  \param[in] my_name The name of the face block to get.
1927    *  \returns The face block, or nullptr if not found.
1928    */
get_face_block(const std::string & my_name)1929   FaceBlock *Region::get_face_block(const std::string &my_name) const
1930   {
1931     IOSS_FUNC_ENTER(m_);
1932     const std::string db_name = get_alias__(my_name);
1933     unsigned int      db_hash = Ioss::Utils::hash(db_name);
1934 
1935     FaceBlock *ge = nullptr;
1936     for (auto fb : faceBlocks) {
1937       if (db_hash == fb->hash() && fb->name() == db_name) {
1938         ge = fb;
1939         break;
1940       }
1941     }
1942     return ge;
1943   }
1944 
1945   /** \brief Get the element block with the given name.
1946    *
1947    *  \param[in] my_name The name of the element block to get.
1948    *  \returns The element block, or nullptr if not found.
1949    */
get_element_block(const std::string & my_name)1950   ElementBlock *Region::get_element_block(const std::string &my_name) const
1951   {
1952     IOSS_FUNC_ENTER(m_);
1953     const std::string db_name = get_alias__(my_name);
1954     unsigned int      db_hash = Ioss::Utils::hash(db_name);
1955 
1956     ElementBlock *ge = nullptr;
1957     for (auto eb : elementBlocks) {
1958       if (db_hash == eb->hash() && eb->name() == db_name) {
1959         ge = eb;
1960         break;
1961       }
1962     }
1963     return ge;
1964   }
1965 
1966   /** \brief Get the structured block with the given name.
1967    *
1968    *  \param[in] my_name The name of the structured block to get.
1969    *  \returns The structured block, or nullptr if not found.
1970    */
get_structured_block(const std::string & my_name)1971   StructuredBlock *Region::get_structured_block(const std::string &my_name) const
1972   {
1973     IOSS_FUNC_ENTER(m_);
1974     const std::string db_name = get_alias__(my_name);
1975     unsigned int      db_hash = Ioss::Utils::hash(db_name);
1976 
1977     StructuredBlock *ge = nullptr;
1978     for (auto sb : structuredBlocks) {
1979       if (db_hash == sb->hash() && sb->name() == db_name) {
1980         ge = sb;
1981         break;
1982       }
1983     }
1984     return ge;
1985   }
1986 
1987   /** \brief Get the side set with the given name.
1988    *
1989    *  \param[in] my_name The name of the side set to get.
1990    *  \returns The side set, or nullptr if not found.
1991    */
get_sideset(const std::string & my_name)1992   SideSet *Region::get_sideset(const std::string &my_name) const
1993   {
1994     IOSS_FUNC_ENTER(m_);
1995     const std::string db_name = get_alias__(my_name);
1996     unsigned int      db_hash = Ioss::Utils::hash(db_name);
1997 
1998     SideSet *ge = nullptr;
1999     for (auto ss : sideSets) {
2000       if (db_hash == ss->hash() && ss->name() == db_name) {
2001         ge = ss;
2002         break;
2003       }
2004     }
2005     return ge;
2006   }
2007 
2008   /** \brief Get the side block with the given name.
2009    *
2010    *  \param[in] my_name The name of the side block to get.
2011    *  \returns The side block, or nullptr if not found.
2012    */
get_sideblock(const std::string & my_name)2013   SideBlock *Region::get_sideblock(const std::string &my_name) const
2014   {
2015     IOSS_FUNC_ENTER(m_);
2016     SideBlock *ge = nullptr;
2017     for (auto ss : sideSets) {
2018       ge = ss->get_side_block(my_name);
2019       if (ge != nullptr) {
2020         break;
2021       }
2022     }
2023     return ge;
2024   }
2025 
2026   /** \brief Get the node set with the given name.
2027    *
2028    *  \param[in] my_name The name of the node set to get.
2029    *  \returns The node set, or nullptr if not found.
2030    */
get_nodeset(const std::string & my_name)2031   NodeSet *Region::get_nodeset(const std::string &my_name) const
2032   {
2033     IOSS_FUNC_ENTER(m_);
2034     const std::string db_name = get_alias__(my_name);
2035     unsigned int      db_hash = Ioss::Utils::hash(db_name);
2036 
2037     NodeSet *ge = nullptr;
2038     for (auto ns : nodeSets) {
2039       if (db_hash == ns->hash() && ns->name() == db_name) {
2040         ge = ns;
2041         break;
2042       }
2043     }
2044     return ge;
2045   }
2046 
2047   /** \brief Get the edge set with the given name.
2048    *
2049    *  \param[in] my_name The name of the edge set to get.
2050    *  \returns The edge set, or nullptr if not found.
2051    */
get_edgeset(const std::string & my_name)2052   EdgeSet *Region::get_edgeset(const std::string &my_name) const
2053   {
2054     IOSS_FUNC_ENTER(m_);
2055     const std::string db_name = get_alias__(my_name);
2056     unsigned int      db_hash = Ioss::Utils::hash(db_name);
2057 
2058     EdgeSet *ge = nullptr;
2059     for (auto es : edgeSets) {
2060       if (db_hash == es->hash() && es->name() == db_name) {
2061         ge = es;
2062         break;
2063       }
2064     }
2065     return ge;
2066   }
2067 
2068   /** \brief Get the face set with the given name.
2069    *
2070    *  \param[in] my_name The name of the face set to get.
2071    *  \returns The face set, or nullptr if not found.
2072    */
get_faceset(const std::string & my_name)2073   FaceSet *Region::get_faceset(const std::string &my_name) const
2074   {
2075     IOSS_FUNC_ENTER(m_);
2076     const std::string db_name = get_alias__(my_name);
2077     unsigned int      db_hash = Ioss::Utils::hash(db_name);
2078 
2079     FaceSet *ge = nullptr;
2080     for (auto fs : faceSets) {
2081       if (db_hash == fs->hash() && fs->name() == db_name) {
2082         ge = fs;
2083         break;
2084       }
2085     }
2086     return ge;
2087   }
2088 
2089   /** \brief Get the element set with the given name.
2090    *
2091    *  \param[in] my_name The name of the element set to get.
2092    *  \returns The element set, or nullptr if not found.
2093    */
get_elementset(const std::string & my_name)2094   ElementSet *Region::get_elementset(const std::string &my_name) const
2095   {
2096     IOSS_FUNC_ENTER(m_);
2097     const std::string db_name = get_alias__(my_name);
2098     unsigned int      db_hash = Ioss::Utils::hash(db_name);
2099 
2100     ElementSet *ge = nullptr;
2101     for (auto es : elementSets) {
2102       if (db_hash == es->hash() && es->name() == db_name) {
2103         ge = es;
2104         break;
2105       }
2106     }
2107     return ge;
2108   }
2109 
2110   /** \brief Get the comm set with the given name.
2111    *
2112    *  \param[in] my_name The name of the comm set to get.
2113    *  \returns The comm set, or nullptr if not found.
2114    */
get_commset(const std::string & my_name)2115   CommSet *Region::get_commset(const std::string &my_name) const
2116   {
2117     IOSS_FUNC_ENTER(m_);
2118     const std::string db_name = get_alias__(my_name);
2119     unsigned int      db_hash = Ioss::Utils::hash(db_name);
2120 
2121     CommSet *ge = nullptr;
2122     for (auto cs : commSets) {
2123       if (db_hash == cs->hash() && cs->name() == db_name) {
2124         ge = cs;
2125         break;
2126       }
2127     }
2128     return ge;
2129   }
2130 
2131   /** \brief Get the coordinate frame with the given id
2132    *
2133    *  \param[in] id The id of the coordinate frame to get.
2134    *  \returns The coordinate frame, or nullptr if not found.
2135    */
get_coordinate_frame(int64_t id)2136   const CoordinateFrame &Region::get_coordinate_frame(int64_t id) const
2137   {
2138     IOSS_FUNC_ENTER(m_);
2139     for (auto &coor_frame : coordinateFrames) {
2140       if (coor_frame.id() == id) {
2141         return coor_frame;
2142       }
2143     }
2144     std::ostringstream errmsg;
2145     fmt::print(errmsg, "Error: Invalid id {} specified for coordinate frame.", id);
2146     IOSS_ERROR(errmsg);
2147   }
2148 
2149   /** \brief Determine whether the entity with the given name and type exists.
2150    *
2151    *  \param[in] my_name The name of the entity to search for.
2152    *  \param[in] io_type The type of the entity.
2153    *  \param[out] my_type A string representing the type if the entity is found.
2154    *                      "INVALID" if the entity is not found or the type is invalid.
2155    *  \returns True if the type is valid and the entity is found.
2156    */
is_valid_io_entity(const std::string & my_name,unsigned int io_type,std::string * my_type)2157   bool Region::is_valid_io_entity(const std::string &my_name, unsigned int io_type,
2158                                   std::string *my_type) const
2159   {
2160     // Search all entities defined on this region for the name 'my_name'.
2161     // If found, then set 'type' (if non-nullptr) to the type of the entity
2162     // (the 'type' values are from client code that was developed prior
2163     // to this function, so they are somewhat exodusII specific...).
2164     if (((io_type & NODEBLOCK) != 0u) && get_node_block(my_name) != nullptr) {
2165       if (my_type != nullptr) {
2166         *my_type = "NODE_BLOCK";
2167       }
2168       return true;
2169     }
2170     if (((io_type & ASSEMBLY) != 0u) && get_assembly(my_name) != nullptr) {
2171       if (my_type != nullptr) {
2172         *my_type = "ASSEMBLY";
2173       }
2174       return true;
2175     }
2176     if (((io_type & BLOB) != 0u) && get_blob(my_name) != nullptr) {
2177       if (my_type != nullptr) {
2178         *my_type = "BLOB";
2179       }
2180       return true;
2181     }
2182     if (((io_type & EDGEBLOCK) != 0u) && get_edge_block(my_name) != nullptr) {
2183       if (my_type != nullptr) {
2184         *my_type = "EDGE_BLOCK";
2185       }
2186       return true;
2187     }
2188     if (((io_type & FACEBLOCK) != 0u) && get_face_block(my_name) != nullptr) {
2189       if (my_type != nullptr) {
2190         *my_type = "FACE_BLOCK";
2191       }
2192       return true;
2193     }
2194     if (((io_type & ELEMENTBLOCK) != 0u) && get_element_block(my_name) != nullptr) {
2195       if (my_type != nullptr) {
2196         *my_type = "ELEMENT_BLOCK";
2197       }
2198       return true;
2199     }
2200     if (((io_type & STRUCTUREDBLOCK) != 0u) && get_structured_block(my_name) != nullptr) {
2201       if (my_type != nullptr) {
2202         *my_type = "STRUCTURED_BLOCK";
2203       }
2204       return true;
2205     }
2206     if (((io_type & SIDESET) != 0u) && get_sideset(my_name) != nullptr) {
2207       if (my_type != nullptr) {
2208         *my_type = "SURFACE";
2209       }
2210       return true;
2211     }
2212     else if (((io_type & NODESET) != 0u) && get_nodeset(my_name) != nullptr) {
2213       if (my_type != nullptr) {
2214         *my_type = "NODESET";
2215       }
2216       return true;
2217     }
2218     else if (((io_type & EDGESET) != 0u) && get_edgeset(my_name) != nullptr) {
2219       if (my_type != nullptr) {
2220         *my_type = "EDGESET";
2221       }
2222       return true;
2223     }
2224     else if (((io_type & FACESET) != 0u) && get_faceset(my_name) != nullptr) {
2225       if (my_type != nullptr) {
2226         *my_type = "FACESET";
2227       }
2228       return true;
2229     }
2230     else if (((io_type & ELEMENTSET) != 0u) && get_elementset(my_name) != nullptr) {
2231       if (my_type != nullptr) {
2232         *my_type = "ELEMENTSET";
2233       }
2234       return true;
2235     }
2236     else if (((io_type & COMMSET) != 0u) && get_commset(my_name) != nullptr) {
2237       if (my_type != nullptr) {
2238         *my_type = "COMMSET";
2239       }
2240       return true;
2241     }
2242     if (my_type != nullptr) {
2243       *my_type = "INVALID";
2244     }
2245     return false;
2246   }
2247 
2248   /** \brief Get the element block containing a specified element.
2249    *
2250    *  \param[in] local_id The local database id (1-based), not the global id.
2251    *  \returns The element block, or nullptr if no element block contains this
2252    *           element (local_id <= 0 or greater than number of elements in database)
2253    */
get_element_block(size_t local_id)2254   ElementBlock *Region::get_element_block(size_t local_id) const
2255   {
2256     IOSS_FUNC_ENTER(m_);
2257     for (auto eb : elementBlocks) {
2258       if (eb->contains(local_id)) {
2259         return eb;
2260       }
2261     }
2262     // Should not reach this point...
2263     std::ostringstream errmsg;
2264     fmt::print(errmsg,
2265                "ERROR: In Ioss::Region::get_element_block, an invalid local_id of {} is specified. "
2266                " The valid range is 1 to {}",
2267                local_id, get_implicit_property("element_count").get_int());
2268     IOSS_ERROR(errmsg);
2269   }
2270 
2271   /** \brief Get the structured block containing a specified global-offset-node.
2272    *
2273    *  \param[in] global_offset The offset of cell-nodes for all blocks; 0-based.
2274    *  \returns The structured block, or nullptr if no structured block contains this
2275    *           node (local_id <= 0 or greater than number of cell-nodes in database)
2276    */
get_structured_block(size_t global_offset)2277   StructuredBlock *Region::get_structured_block(size_t global_offset) const
2278   {
2279     IOSS_FUNC_ENTER(m_);
2280     for (auto sb : structuredBlocks) {
2281       if (sb->contains(global_offset)) {
2282         return sb;
2283       }
2284     }
2285     // Should not reach this point...
2286     std::ostringstream errmsg;
2287     fmt::print(errmsg,
2288                "ERROR: In Ioss::Region::get_structured_block, an invalid global_offset of {} is "
2289                "specified.",
2290                global_offset);
2291     IOSS_ERROR(errmsg);
2292   }
2293 
2294   /** \brief Get an implicit property -- These are calcuated from data stored
2295    *         in the grouping entity instead of having an explicit value assigned.
2296    *
2297    *  An example would be 'element_block_count' for a region.
2298    *
2299    *  \param[in] my_name The property name.
2300    *  \returns The property.
2301    */
get_implicit_property(const std::string & my_name)2302   Property Region::get_implicit_property(const std::string &my_name) const
2303   {
2304     if (my_name == "spatial_dimension") {
2305       if (!nodeBlocks.empty()) {
2306         return nodeBlocks[0]->get_property("component_degree");
2307       }
2308 
2309       return Property(my_name, 0);
2310     }
2311 
2312     if (my_name == "node_block_count") {
2313       return Property(my_name, static_cast<int>(nodeBlocks.size()));
2314     }
2315 
2316     if (my_name == "edge_block_count") {
2317       return Property(my_name, static_cast<int>(edgeBlocks.size()));
2318     }
2319 
2320     if (my_name == "face_block_count") {
2321       return Property(my_name, static_cast<int>(faceBlocks.size()));
2322     }
2323 
2324     if (my_name == "element_block_count") {
2325       return Property(my_name, static_cast<int>(elementBlocks.size()));
2326     }
2327 
2328     if (my_name == "structured_block_count") {
2329       return Property(my_name, static_cast<int>(structuredBlocks.size()));
2330     }
2331 
2332     if (my_name == "assembly_count") {
2333       return Property(my_name, static_cast<int>(assemblies.size()));
2334     }
2335 
2336     if (my_name == "blob_count") {
2337       return Property(my_name, static_cast<int>(blobs.size()));
2338     }
2339 
2340     if (my_name == "side_set_count") {
2341       return Property(my_name, static_cast<int>(sideSets.size()));
2342     }
2343 
2344     if (my_name == "node_set_count") {
2345       return Property(my_name, static_cast<int>(nodeSets.size()));
2346     }
2347 
2348     if (my_name == "edge_set_count") {
2349       return Property(my_name, static_cast<int>(edgeSets.size()));
2350     }
2351 
2352     if (my_name == "face_set_count") {
2353       return Property(my_name, static_cast<int>(faceSets.size()));
2354     }
2355 
2356     if (my_name == "element_set_count") {
2357       return Property(my_name, static_cast<int>(elementSets.size()));
2358     }
2359 
2360     if (my_name == "comm_set_count") {
2361       return Property(my_name, static_cast<int>(commSets.size()));
2362     }
2363 
2364     if (my_name == "coordinate_frame_count") {
2365       return Property(my_name, static_cast<int>(coordinateFrames.size()));
2366     }
2367 
2368     if (my_name == "state_count") {
2369       return Property(my_name, stateCount);
2370     }
2371 
2372     if (my_name == "current_state") {
2373       return Property(my_name, currentState);
2374     }
2375 
2376     if (my_name == "element_count") {
2377       int64_t count = 0;
2378       for (auto eb : elementBlocks) {
2379         count += eb->entity_count();
2380       }
2381       return Property(my_name, count);
2382     }
2383 
2384     if (my_name == "cell_count") {
2385       int64_t count = 0;
2386       for (auto eb : structuredBlocks) {
2387         count += eb->get_property("cell_count").get_int();
2388       }
2389       return Property(my_name, count);
2390     }
2391 
2392     if (my_name == "face_count") {
2393       int64_t count = 0;
2394       for (auto fb : faceBlocks) {
2395         count += fb->entity_count();
2396       }
2397       return Property(my_name, count);
2398     }
2399 
2400     if (my_name == "edge_count") {
2401       int64_t count = 0;
2402       for (auto eb : edgeBlocks) {
2403         count += eb->entity_count();
2404       }
2405       return Property(my_name, count);
2406     }
2407 
2408     if (my_name == "node_count") {
2409       int64_t count = 0;
2410       for (auto nb : nodeBlocks) {
2411         count += nb->entity_count();
2412       }
2413       return Property(my_name, count);
2414     }
2415 
2416     if (my_name == "database_name") {
2417       std::string filename = get_database()->get_filename();
2418       return Property(my_name, filename);
2419     }
2420 
2421     {
2422       return GroupingEntity::get_implicit_property(my_name);
2423     }
2424   }
2425 
internal_get_field_data(const Field & field,void * data,size_t data_size)2426   int64_t Region::internal_get_field_data(const Field &field, void *data, size_t data_size) const
2427   {
2428     return get_database()->get_field(this, field, data, data_size);
2429   }
2430 
internal_put_field_data(const Field & field,void * data,size_t data_size)2431   int64_t Region::internal_put_field_data(const Field &field, void *data, size_t data_size) const
2432   {
2433     return get_database()->put_field(this, field, data, data_size);
2434   }
2435 
2436   /** \brief Transfer all relevant aliases from this region to another region
2437    *
2438    *  \param[in] to The region to which the aliases are to be transferred.
2439    */
transfer_mesh_aliases(Region * to)2440   void Region::transfer_mesh_aliases(Region *to) const
2441   {
2442     IOSS_FUNC_ENTER(m_);
2443     // Iterate through list, [ returns <alias, base_entity_name> ], if
2444     // 'base_entity_name' is defined on the restart file, add 'alias' as
2445     // an alias for it...
2446     for (const auto &alias_pair : aliases_) {
2447       std::string alias = alias_pair.first;
2448       std::string base  = alias_pair.second;
2449       if (alias != base && to->get_entity(base) != nullptr) {
2450         to->add_alias__(base, alias);
2451       }
2452     }
2453   }
2454 
2455   /** \brief Ensure that the restart and results files have the same ids.
2456    *
2457    *  There is very little connection between an input (mesh) database
2458    *  and an output (results/restart) database.  Basically, the entity
2459    *  names are the same between the two files.  This works fine in the
2460    *  case that an input database has 'generated' entity names of the
2461    *  form 'block_10' or 'surface_32' since then the output database
2462    *  can de-generate or decode the name and infer that the block
2463    *  should have an id of 10 and the surface an id of 32.
2464    *
2465    *  However, if alias or other renaming happens, then the output
2466    *  block may have a name of the form 'fireset' and the underlying
2467    *  database cannot infer that the id of the block should be 10.
2468    *  Instead, it sets the id to an arbitrary number (1,2,...).  This
2469    *  is annoying in the case of the results file since there is no
2470    *  correspondence between the mesh numbering and the results
2471    *  numbering. In the case of the restart output file, it can be
2472    *  disastrous since when the file is used to restart the analysis,
2473    *  there is no match between the mesh blocks and those found on the
2474    *  restart file and the restart fails.
2475    *
2476    *  So... We need to somehow ensure that the restart (and results)
2477    *  files have the same ids.  To do this, we do the following:
2478    *
2479    *  1. The mesh database will set the property 'id' on input.
2480    *
2481    *  2. The results/restart files will have a 'name' based on either
2482    *     the true name or an alias name.  Whichever, that alias will
2483    *     appear on the mesh database also, so we can query the mesh
2484    *     database aliases to get the entity.
2485    *
2486    *  3. Once we have the entity, we can query the 'id' property and
2487    *     add the same 'id' property to the results/restart database.
2488    *  4. Also set the 'name' property to the base 'name' on the output file.
2489    *
2490    *  5. Note that a property may already exist and must be removed
2491    *    before the 'correct' value is set.
2492    */
synchronize_id_and_name(const Region * from,bool sync_attribute_field_names)2493   void Region::synchronize_id_and_name(const Region *from, bool sync_attribute_field_names)
2494   {
2495     for (const auto &alias_pair : aliases_) {
2496       std::string alias = alias_pair.first;
2497       std::string base  = alias_pair.second;
2498 
2499       if (alias == base) {
2500 
2501         // Query the 'from' database to get the entity (if any) referred
2502         // to by the 'alias'
2503         GroupingEntity *ge = from->get_entity(base);
2504 
2505         if (ge != nullptr) {
2506           // Get the entity from this region... Must be non-nullptr
2507           GroupingEntity *this_ge = get_entity(base);
2508           if (this_ge == nullptr) {
2509             std::ostringstream errmsg;
2510             fmt::print(errmsg,
2511                        "INTERNAL ERROR: Could not find entity '{}' in synchronize_id_and_name() "
2512                        "                [{}]\n",
2513                        base, get_database()->get_filename());
2514             IOSS_ERROR(errmsg);
2515           }
2516 
2517           // See if there is an 'id' property...
2518           if (ge->property_exists(id_str())) {
2519             int64_t id = ge->get_property(id_str()).get_int();
2520             this_ge->property_update(id_str(), id);
2521           }
2522           else {
2523             // No id, make sure the base name matches in both databases...
2524             // There is always a 'name' property on an entity
2525             if (this_ge->name() != base) {
2526               this_ge->set_name(base);
2527             }
2528           }
2529 
2530           // See if there is an 'db_name' property...
2531           if (ge->property_exists(db_name_str())) {
2532             std::string db_name = ge->get_property(db_name_str()).get_string();
2533             // Set the new property
2534             this_ge->property_update(db_name_str(), db_name);
2535           }
2536 
2537           // See if there is a 'original_topology_type' property...
2538           if (ge->property_exists(orig_topo_str())) {
2539             std::string oes = ge->get_property(orig_topo_str()).get_string();
2540             this_ge->property_update(orig_topo_str(), oes);
2541           }
2542 
2543           // Specific to entity blocks. Transfer the "original_block_order"
2544           // property.
2545           if (ge->property_exists(orig_block_order())) {
2546             int64_t offset = ge->get_property(orig_block_order()).get_int();
2547             this_ge->property_update(orig_block_order(), offset);
2548           }
2549 
2550           if (sync_attribute_field_names) {
2551             // If there are any attribute fields, then copy those over
2552             // to the new entity in order to maintain the same order
2553             // since some codes access attributes by implicit order and
2554             // not name... (typically, element blocks only)
2555             size_t count = this_ge->entity_count();
2556 
2557             Ioss::NameList attr_fields;
2558             ge->field_describe(Ioss::Field::ATTRIBUTE, &attr_fields);
2559             for (auto &field_name : attr_fields) {
2560               const Ioss::Field &field = ge->get_fieldref(field_name);
2561               if (this_ge->field_exists(field_name)) {
2562                 // If the field is already defined on the entity, make
2563                 // sure that the attribute index matches...
2564                 size_t             index      = field.get_index();
2565                 const Ioss::Field &this_field = this_ge->get_fieldref(field_name);
2566                 this_field.set_index(index);
2567               }
2568               else {
2569                 // If the field does not already exist, add it to the
2570                 // output node block
2571                 if (field.raw_count() != count) {
2572                   Ioss::Field new_field(field);
2573                   new_field.reset_count(count);
2574                   this_ge->field_add(new_field);
2575                 }
2576                 else {
2577                   this_ge->field_add(field);
2578                 }
2579               }
2580             }
2581           }
2582         }
2583       }
2584     }
2585 
2586     for (const auto &alias_pair : aliases_) {
2587       std::string alias = alias_pair.first;
2588       std::string base  = alias_pair.second;
2589 
2590       if (alias != base) {
2591         GroupingEntity *ge = get_entity(base);
2592         if (ge != nullptr) {
2593           add_alias__(base, alias);
2594         }
2595       }
2596     }
2597   }
2598 } // namespace Ioss
2599