www.mooseframework.org
GrainTracker.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "GrainTracker.h"
11 
12 // MOOSE includes
14 #include "GeneratedMesh.h"
15 #include "MooseMesh.h"
16 #include "MooseVariable.h"
17 #include "NonlinearSystem.h"
18 
19 #include "libmesh/periodic_boundary_base.h"
20 
21 // C++ includes
22 #include <algorithm>
23 #include <limits>
24 #include <numeric>
25 
26 template <>
27 void
28 dataStore(std::ostream & stream, GrainTracker::PartialFeatureData & feature, void * context)
29 {
30  storeHelper(stream, feature.boundary_intersection, context);
31  storeHelper(stream, feature.id, context);
32  storeHelper(stream, feature.centroid, context);
33  storeHelper(stream, feature.status, context);
34 }
35 
36 template <>
37 void
38 dataLoad(std::istream & stream, GrainTracker::PartialFeatureData & feature, void * context)
39 {
40  loadHelper(stream, feature.boundary_intersection, context);
41  loadHelper(stream, feature.id, context);
42  loadHelper(stream, feature.centroid, context);
43  loadHelper(stream, feature.status, context);
44 }
45 
46 registerMooseObject("PhaseFieldApp", GrainTracker);
47 
48 template <>
49 InputParameters
51 {
52  InputParameters params = validParams<FeatureFloodCount>();
54 
55  // FeatureFloodCount adds a relationship manager, but we need to extend that for GrainTracker
56  params.clearRelationshipManagers();
57 
58  params.addRelationshipManager(
59  "ElementSideNeighborLayers",
60  Moose::RelationshipManagerType::GEOMETRIC,
61 
62  [](const InputParameters & obj_params, InputParameters & rm_params) {
63  rm_params.set<unsigned short>("layers") = obj_params.get<unsigned short>("halo_level");
64  }
65 
66  );
67 
68  params.addRelationshipManager("ElementSideNeighborLayers",
69  Moose::RelationshipManagerType::ALGEBRAIC);
70 
71  // The GrainTracker requires non-volatile storage for tracking grains across invocations.
72  params.set<bool>("restartable_required") = true;
73 
74  params.addClassDescription("Grain Tracker object for running reduced order parameter simulations "
75  "without grain coalescence.");
76 
77  return params;
78 }
79 
80 GrainTracker::GrainTracker(const InputParameters & parameters)
81  : FeatureFloodCount(parameters),
83  _tracking_step(getParam<int>("tracking_step")),
84  _halo_level(getParam<unsigned short>("halo_level")),
85  _max_remap_recursion_depth(getParam<unsigned short>("max_remap_recursion_depth")),
86  _n_reserve_ops(getParam<unsigned short>("reserve_op")),
87  _reserve_op_index(_n_reserve_ops <= _n_vars ? _n_vars - _n_reserve_ops : 0),
88  _reserve_op_threshold(getParam<Real>("reserve_op_threshold")),
89  _remap(getParam<bool>("remap_grains")),
90  _tolerate_failure(getParam<bool>("tolerate_failure")),
91  _nl(_fe_problem.getNonlinearSystemBase()),
92  _poly_ic_uo(parameters.isParamValid("polycrystal_ic_uo")
93  ? &getUserObject<PolycrystalUserObjectBase>("polycrystal_ic_uo")
94  : nullptr),
95  _verbosity_level(getParam<short>("verbosity_level")),
96  _first_time(declareRestartableData<bool>("first_time", true)),
97  _error_on_grain_creation(getParam<bool>("error_on_grain_creation")),
98  _reserve_grain_first_index(0),
99  _old_max_grain_id(0),
100  _max_curr_grain_id(declareRestartableData<unsigned int>("max_curr_grain_id", invalid_id)),
101  _is_transient(_subproblem.isTransient()),
102  _finalize_timer(registerTimedSection("finalize", 1)),
103  _remap_timer(registerTimedSection("remapGrains", 2)),
104  _track_grains(registerTimedSection("trackGrains", 2)),
105  _broadcast_update(registerTimedSection("broadCastUpdate", 2)),
106  _update_field_info(registerTimedSection("updateFieldInfo", 2))
107 {
108  if (_tolerate_failure)
109  paramInfo("tolerate_failure",
110  "Tolerate failure has been set to true. Non-physical simulation results "
111  "are possible, you will be notified in the event of a failed remapping operation.");
112 
113  if (_tracking_step > 0 && _poly_ic_uo)
114  mooseError("Can't start tracking after the initial condition when using a polycrystal_ic_uo");
115 }
116 
118 
119 Real
120 GrainTracker::getEntityValue(dof_id_type entity_id,
121  FieldType field_type,
122  std::size_t var_index) const
123 {
124  if (_t_step < _tracking_step)
125  return 0;
126 
127  return FeatureFloodCount::getEntityValue(entity_id, field_type, var_index);
128 }
129 
130 const std::vector<unsigned int> &
131 GrainTracker::getVarToFeatureVector(dof_id_type elem_id) const
132 {
134 }
135 
136 unsigned int
137 GrainTracker::getFeatureVar(unsigned int feature_id) const
138 {
139  return FeatureFloodCount::getFeatureVar(feature_id);
140 }
141 
142 std::size_t
144 {
145  // Note: This value is parallel consistent, see FeatureFloodCount::communicateAndMerge()
146  return _feature_count;
147 }
148 
149 std::size_t
151 {
152  // Note: This value is parallel consistent, see assignGrains()/trackGrains()
154 }
155 
156 Point
157 GrainTracker::getGrainCentroid(unsigned int grain_id) const
158 {
159  mooseAssert(grain_id < _feature_id_to_local_index.size(), "Grain ID out of bounds");
160  auto grain_index = _feature_id_to_local_index[grain_id];
161 
162  if (grain_index != invalid_size_t)
163  {
164  mooseAssert(_feature_id_to_local_index[grain_id] < _feature_sets.size(),
165  "Grain index out of bounds");
166  // Note: This value is parallel consistent, see GrainTracker::broadcastAndUpdateGrainData()
167  return _feature_sets[_feature_id_to_local_index[grain_id]]._centroid;
168  }
169 
170  // Inactive grain
171  return Point();
172 }
173 
174 bool
175 GrainTracker::doesFeatureIntersectBoundary(unsigned int feature_id) const
176 {
177  // TODO: This data structure may need to be turned into a Multimap
178  mooseAssert(feature_id < _feature_id_to_local_index.size(), "Grain ID out of bounds");
179 
180  auto feature_index = _feature_id_to_local_index[feature_id];
181  if (feature_index != invalid_size_t)
182  {
183  mooseAssert(feature_index < _feature_sets.size(), "Grain index out of bounds");
184  return _feature_sets[feature_index]._boundary_intersection != BoundaryIntersection::NONE;
185  }
186 
187  return false;
188 }
189 
190 bool
192 {
193  // TODO: This data structure may need to be turned into a Multimap
194  mooseAssert(feature_id < _feature_id_to_local_index.size(), "Grain ID out of bounds");
195 
196  auto feature_index = _feature_id_to_local_index[feature_id];
197  if (feature_index != invalid_size_t)
198  {
199  mooseAssert(feature_index < _feature_sets.size(), "Grain index out of bounds");
200  return ((_feature_sets[feature_index]._boundary_intersection &
202  }
203 
204  return false;
205 }
206 
207 bool
208 GrainTracker::isFeaturePercolated(unsigned int feature_id) const
209 {
210  // TODO: This data structure may need to be turned into a Multimap
211  mooseAssert(feature_id < _feature_id_to_local_index.size(), "Grain ID out of bounds");
212 
213  auto feature_index = _feature_id_to_local_index[feature_id];
214  if (feature_index != invalid_size_t)
215  {
216  mooseAssert(feature_index < _feature_sets.size(), "Grain index out of bounds");
217  bool primary = ((_feature_sets[feature_index]._boundary_intersection &
220  bool secondary = ((_feature_sets[feature_index]._boundary_intersection &
223  return (primary && secondary);
224  }
225 
226  return false;
227 }
228 
229 void
231 {
232  // Don't track grains if the current simulation step is before the specified tracking step
233  if (_t_step < _tracking_step)
234  return;
235 
241  if (!_first_time)
243 
245 }
246 
247 void
249 {
250  // Update the element ID ranges for use when computing halo maps
251  if (_compute_halo_maps && _mesh.isDistributedMesh())
252  {
253  _all_ranges.clear();
254 
255  auto range = std::make_pair(std::numeric_limits<dof_id_type>::max(),
256  std::numeric_limits<dof_id_type>::min());
257  for (const auto & current_elem : _mesh.getMesh().active_local_element_ptr_range())
258  {
259  auto id = current_elem->id();
260  if (id < range.first)
261  range.first = id;
262  else if (id > range.second)
263  range.second = id;
264  }
265 
266  _communicator.gather(0, range, _all_ranges);
267  }
268 
270 }
271 
272 void
274 {
275  // Don't track grains if the current simulation step is before the specified tracking step
276  if (_t_step < _tracking_step)
277  return;
278 
279  if (_poly_ic_uo && _first_time)
280  return;
281 
283 }
284 
285 Real
286 GrainTracker::getThreshold(std::size_t var_index) const
287 {
288  // If we are inspecting a reserve op parameter, we need to make sure
289  // that there is an entity above the reserve_op threshold before
290  // starting the flood of the feature.
291  if (var_index >= _reserve_op_index)
292  return _reserve_op_threshold;
293  else
294  return _step_threshold;
295 }
296 
297 void
299 {
300  mooseAssert(_first_time, "This method should only be called on the first invocation");
301 
302  _feature_sets.clear();
303 
309  if (_is_master)
310  {
311  const auto & features = ffc_object.getFeatures();
312  for (auto & feature : features)
313  _feature_sets.emplace_back(feature.duplicate());
314 
315  _feature_count = _feature_sets.size();
316  }
317  else
318  {
319  const auto & features = ffc_object.getFeatures();
320  _partial_feature_sets[0].clear();
321  for (auto & feature : features)
322  _partial_feature_sets[0].emplace_back(feature.duplicate());
323  }
324 
325  // Make sure that feature count is communicated to all ranks
326  _communicator.broadcast(_feature_count);
327 }
328 
329 void
331 {
332  // Don't track grains if the current simulation step is before the specified tracking step
333  if (_t_step < _tracking_step)
334  return;
335 
336  TIME_SECTION(_finalize_timer);
337 
338  // Expand the depth of the halos around all grains
339  auto num_halo_layers = _halo_level >= 1
340  ? _halo_level - 1
341  : 0; // The first level of halos already exists so subtract one
342 
343  if (_poly_ic_uo && _first_time)
345  else
346  {
347  expandEdgeHalos(num_halo_layers);
348 
349  // Build up the grain map on the root processor
351  }
352 
356  if (_first_time)
357  assignGrains();
358  else
359  trackGrains();
360 
361  if (_verbosity_level > 1)
362  _console << "Finished inside of trackGrains" << std::endl;
363 
368 
372  if (_remap)
373  remapGrains();
374 
375  updateFieldInfo();
376  if (_verbosity_level > 1)
377  _console << "Finished inside of updateFieldInfo\n";
378 
379  // Set the first time flag false here (after all methods of finalize() have completed)
380  _first_time = false;
381 
382  // TODO: Release non essential memory
383  if (_verbosity_level > 0)
384  _console << "Finished inside of GrainTracker\n" << std::endl;
385 }
386 
387 void
389 {
390  TIME_SECTION(_broadcast_update);
391 
392  std::vector<PartialFeatureData> root_feature_data;
393  std::vector<std::string> send_buffer(1), recv_buffer;
394 
395  if (_is_master)
396  {
397  root_feature_data.reserve(_feature_sets.size());
398 
399  // Populate a subset of the information in a small data structure
400  std::transform(_feature_sets.begin(),
401  _feature_sets.end(),
402  std::back_inserter(root_feature_data),
403  [](FeatureData & feature) {
404  PartialFeatureData partial_feature;
405  partial_feature.boundary_intersection = feature._boundary_intersection;
406  partial_feature.id = feature._id;
407  partial_feature.centroid = feature._centroid;
408  partial_feature.status = feature._status;
409  return partial_feature;
410  });
411 
412  std::ostringstream oss;
413  dataStore(oss, root_feature_data, this);
414  send_buffer[0].assign(oss.str());
415  }
416 
417  // Broadcast the data to all ranks
418  _communicator.broadcast_packed_range((void *)(nullptr),
419  send_buffer.begin(),
420  send_buffer.end(),
421  (void *)(nullptr),
422  std::back_inserter(recv_buffer));
423 
424  // Unpack and update
425  if (!_is_master)
426  {
427  std::istringstream iss;
428  iss.str(recv_buffer[0]);
429  iss.clear();
430 
431  dataLoad(iss, root_feature_data, this);
432 
433  for (const auto & partial_data : root_feature_data)
434  {
435  // See if this processor has a record of this grain
436  if (partial_data.id < _feature_id_to_local_index.size() &&
437  _feature_id_to_local_index[partial_data.id] != invalid_size_t)
438  {
439  auto & grain = _feature_sets[_feature_id_to_local_index[partial_data.id]];
440  grain._boundary_intersection = partial_data.boundary_intersection;
441  grain._centroid = partial_data.centroid;
442  if (partial_data.status == Status::INACTIVE)
443  grain._status = Status::INACTIVE;
444  }
445  }
446  }
447 }
448 
449 void
451 {
452  mooseAssert(_first_time, "assignGrains may only be called on the first tracking step");
453 
459  if (_is_master)
460  {
461  // Find the largest grain ID, this requires sorting if the ID is not already set
462  sortAndLabel();
463 
464  if (_feature_sets.empty())
465  {
468  }
469  else
470  {
471  _max_curr_grain_id = _feature_sets.back()._id;
473  }
474 
475  for (auto & grain : _feature_sets)
476  grain._status = Status::MARKED; // Mark the grain
477 
478  } // is_master
479 
480  /*************************************************************
481  ****************** COLLECTIVE WORK SECTION ******************
482  *************************************************************/
483 
484  // Make IDs on all non-master ranks consistent
486 
487  // Build up an id to index map
488  _communicator.broadcast(_max_curr_grain_id);
490 
491  // Now trigger the newGrainCreated() callback on all ranks
493  for (unsigned int new_id = 0; new_id <= _max_curr_grain_id; ++new_id)
494  newGrainCreated(new_id);
495 }
496 
497 void
499 {
500  TIME_SECTION(_track_grains);
501 
502  mooseAssert(!_first_time, "Track grains may only be called when _tracking_step > _t_step");
503 
504  // Used to track indices for which to trigger the new grain callback on (used on all ranks)
506 
511  if (_is_master)
512  {
513  // Reset Status on active unique grains
514  std::vector<unsigned int> map_sizes(_maps_size);
515  for (auto & grain : _feature_sets_old)
516  {
517  if (grain._status != Status::INACTIVE)
518  {
519  grain._status = Status::CLEAR;
520  map_sizes[grain._var_index]++;
521  }
522  }
523 
524  // Print out stats on overall tracking changes per var_index
525  if (_verbosity_level > 0)
526  {
527  _console << "\nGrain Tracker Status:";
528  for (MooseIndex(_maps_size) map_num = 0; map_num < _maps_size; ++map_num)
529  {
530  _console << "\nGrains active index " << map_num << ": " << map_sizes[map_num] << " -> "
531  << _feature_counts_per_map[map_num];
532  if (map_sizes[map_num] > _feature_counts_per_map[map_num])
533  _console << "--";
534  else if (map_sizes[map_num] < _feature_counts_per_map[map_num])
535  _console << "++";
536  }
537  _console << '\n' << std::endl;
538  }
539 
540  // Before we track grains, lets sort them so that we get parallel consistent answers
541  std::sort(_feature_sets.begin(), _feature_sets.end());
542 
549  std::vector<std::size_t> new_grain_index_to_existing_grain_index(_feature_sets.size(),
551 
552  for (MooseIndex(_feature_sets_old) old_grain_index = 0;
553  old_grain_index < _feature_sets_old.size();
554  ++old_grain_index)
555  {
556  auto & old_grain = _feature_sets_old[old_grain_index];
557 
558  if (old_grain._status == Status::INACTIVE) // Don't try to find matches for inactive grains
559  continue;
560 
561  std::size_t closest_match_index = invalid_size_t;
562  Real min_centroid_diff = std::numeric_limits<Real>::max();
563 
569  // clang-format off
570  auto start_it =
571  std::lower_bound(_feature_sets.begin(), _feature_sets.end(), old_grain._var_index,
572  [](const FeatureData & item, std::size_t var_index)
573  {
574  return item._var_index < var_index;
575  });
576  // clang-format on
577 
578  // We only need to examine grains that have matching variable indices
579  bool any_boxes_intersect = false;
580  for (MooseIndex(_feature_sets)
581  new_grain_index = std::distance(_feature_sets.begin(), start_it);
582  new_grain_index < _feature_sets.size() &&
583  _feature_sets[new_grain_index]._var_index == old_grain._var_index;
584  ++new_grain_index)
585  {
586  auto & new_grain = _feature_sets[new_grain_index];
587 
592  if (new_grain.boundingBoxesIntersect(old_grain))
593  {
594  any_boxes_intersect = true;
595  Real curr_centroid_diff = centroidRegionDistance(old_grain._bboxes, new_grain._bboxes);
596  if (curr_centroid_diff <= min_centroid_diff)
597  {
598  closest_match_index = new_grain_index;
599  min_centroid_diff = curr_centroid_diff;
600  }
601  }
602  }
603 
604  if (_verbosity_level > 2 && !any_boxes_intersect)
605  _console << "\nNo intersecting bounding boxes found while trying to match grain "
606  << old_grain;
607 
608  // found a match
609  if (closest_match_index != invalid_size_t)
610  {
617  auto curr_index = new_grain_index_to_existing_grain_index[closest_match_index];
618  if (curr_index != invalid_size_t)
619  {
620  // The new feature being competed for
621  auto & new_grain = _feature_sets[closest_match_index];
622 
623  // The other old grain competing to match up to the same new grain
624  auto & other_old_grain = _feature_sets_old[curr_index];
625 
626  auto centroid_diff1 = centroidRegionDistance(new_grain._bboxes, old_grain._bboxes);
627  auto centroid_diff2 = centroidRegionDistance(new_grain._bboxes, other_old_grain._bboxes);
628 
629  auto & inactive_grain = (centroid_diff1 < centroid_diff2) ? other_old_grain : old_grain;
630 
631  inactive_grain._status = Status::INACTIVE;
632  if (_verbosity_level > 0)
633  {
634  _console << COLOR_GREEN << "Marking Grain " << inactive_grain._id
635  << " as INACTIVE (variable index: " << inactive_grain._var_index << ")\n"
636  << COLOR_DEFAULT;
637  if (_verbosity_level > 1)
638  _console << inactive_grain;
639  }
640 
646  if (&inactive_grain == &other_old_grain)
647  new_grain_index_to_existing_grain_index[closest_match_index] = old_grain_index;
648  }
649  else
650  new_grain_index_to_existing_grain_index[closest_match_index] = old_grain_index;
651  }
652  }
653 
654  // Mark all resolved grain matches
655  for (MooseIndex(new_grain_index_to_existing_grain_index) new_index = 0;
656  new_index < new_grain_index_to_existing_grain_index.size();
657  ++new_index)
658  {
659  auto curr_index = new_grain_index_to_existing_grain_index[new_index];
660 
661  // This may be a new grain, we'll handle that case below
662  if (curr_index == invalid_size_t)
663  continue;
664 
665  mooseAssert(_feature_sets_old[curr_index]._id != invalid_id,
666  "Invalid ID in old grain structure");
667 
668  _feature_sets[new_index]._id = _feature_sets_old[curr_index]._id; // Transfer ID
669  _feature_sets[new_index]._status = Status::MARKED; // Mark the status in the new set
670  _feature_sets_old[curr_index]._status = Status::MARKED; // Mark the status in the old set
671  }
672 
687  // Case 1 (new grains in _feature_sets):
688  for (MooseIndex(_feature_sets) grain_num = 0; grain_num < _feature_sets.size(); ++grain_num)
689  {
690  auto & grain = _feature_sets[grain_num];
691 
692  // New Grain
693  if (grain._status == Status::CLEAR)
694  {
713  // clang-format off
714  auto start_it =
715  std::lower_bound(_feature_sets.begin(), _feature_sets.end(), grain._var_index,
716  [](const FeatureData & item, std::size_t var_index)
717  {
718  return item._var_index < var_index;
719  });
720  // clang-format on
721 
722  // Loop over matching variable indices
723  for (MooseIndex(_feature_sets)
724  new_grain_index = std::distance(_feature_sets.begin(), start_it);
725  new_grain_index < _feature_sets.size() &&
726  _feature_sets[new_grain_index]._var_index == grain._var_index;
727  ++new_grain_index)
728  {
729  auto & other_grain = _feature_sets[new_grain_index];
730 
731  // Splitting grain?
732  if (grain_num != new_grain_index && // Make sure indices aren't pointing at the same grain
733  other_grain._status == Status::MARKED && // and that the other grain is indeed marked
734  other_grain.boundingBoxesIntersect(grain) && // and the bboxes intersect
735  other_grain.halosIntersect(grain)) // and the halos also intersect
736  // TODO: Inspect combined volume and see if it's "close" to the expected value
737  {
738  grain._id = other_grain._id; // Set the duplicate ID
739  grain._status = Status::MARKED; // Mark it
740 
741  if (_verbosity_level > 0)
742  _console << COLOR_YELLOW << "Split Grain Detected #" << grain._id
743  << " (variable index: " << grain._var_index << ")\n"
744  << COLOR_DEFAULT;
745  if (_verbosity_level > 1)
746  _console << grain << other_grain;
747  }
748  }
749 
750  if (grain._var_index < _reserve_op_index)
751  {
770  if (_verbosity_level > 1)
771  _console << COLOR_YELLOW
772  << "Trying harder to detect a split grain while examining grain on variable "
773  "index "
774  << grain._var_index << '\n'
775  << COLOR_DEFAULT;
776 
777  std::vector<std::size_t> old_grain_indices;
778  for (MooseIndex(_feature_sets_old) old_grain_index = 0;
779  old_grain_index < _feature_sets_old.size();
780  ++old_grain_index)
781  {
782  auto & old_grain = _feature_sets_old[old_grain_index];
783 
784  if (old_grain._status == Status::INACTIVE)
785  continue;
786 
792  if (grain._var_index == old_grain._var_index &&
793  grain.boundingBoxesIntersect(old_grain) && grain.halosIntersect(old_grain))
794  old_grain_indices.push_back(old_grain_index);
795  }
796 
797  if (old_grain_indices.size() == 1)
798  {
799  grain._id = _feature_sets_old[old_grain_indices[0]]._id;
800  grain._status = Status::MARKED;
801 
802  if (_verbosity_level > 0)
803  _console << COLOR_YELLOW << "Split Grain Detected #" << grain._id
804  << " (variable index: " << grain._var_index << ")\n"
805  << COLOR_DEFAULT;
806  }
807  else if (old_grain_indices.size() > 1)
808  _console
809  << COLOR_RED << "Split Grain Likely Detected #" << grain._id
810  << " Need more information to find correct candidate - contact a developer!\n\n"
811  << COLOR_DEFAULT;
812  }
813 
814  // Must be a nucleating grain (status is still not set)
815  if (grain._status == Status::CLEAR)
816  {
817  auto new_index = getNextUniqueID();
818  grain._id = new_index; // Set the ID
819  grain._status = Status::MARKED; // Mark it
820 
821  if (_verbosity_level > 0)
822  _console << COLOR_YELLOW << "Nucleating Grain Detected "
823  << " (variable index: " << grain._var_index << ")\n"
824  << COLOR_DEFAULT;
825  if (_verbosity_level > 1)
826  _console << grain;
827  }
828  }
829  }
830 
831  // Case 2 (inactive grains in _feature_sets_old)
832  for (auto & grain : _feature_sets_old)
833  {
834  if (grain._status == Status::CLEAR)
835  {
836  grain._status = Status::INACTIVE;
837  if (_verbosity_level > 0)
838  {
839  _console << COLOR_GREEN << "Marking Grain " << grain._id
840  << " as INACTIVE (variable index: " << grain._var_index << ")\n"
841  << COLOR_DEFAULT;
842  if (_verbosity_level > 1)
843  _console << grain;
844  }
845  }
846  }
847  } // is_master
848 
849  /*************************************************************
850  ****************** COLLECTIVE WORK SECTION ******************
851  *************************************************************/
852 
853  // Make IDs on all non-master ranks consistent
855 
856  // Build up an id to index map
857  _communicator.broadcast(_max_curr_grain_id);
859 
864  {
865  for (auto new_id = _old_max_grain_id + 1; new_id <= _max_curr_grain_id; ++new_id)
866  {
867  // Don't trigger the callback on the reserve IDs
869  {
870  // See if we've been instructed to terminate with an error
872  mooseError(
873  "Error: New grain detected and \"error_on_new_grain_creation\" is set to true");
874  else
875  newGrainCreated(new_id);
876  }
877  }
878  }
879 }
880 
881 void
882 GrainTracker::newGrainCreated(unsigned int new_grain_id)
883 {
884  if (!_first_time && _is_master)
885  {
886  mooseAssert(new_grain_id < _feature_id_to_local_index.size(), "new_grain_id is out of bounds");
887  auto grain_index = _feature_id_to_local_index[new_grain_id];
888  mooseAssert(grain_index != invalid_size_t && grain_index < _feature_sets.size(),
889  "new_grain_id appears to be invalid");
890 
891  const auto & grain = _feature_sets[grain_index];
892  _console << COLOR_YELLOW
893  << "\n*****************************************************************************"
894  << "\nCouldn't find a matching grain while working on variable index: "
895  << grain._var_index << "\nCreating new unique grain: " << new_grain_id << '\n'
896  << grain
897  << "\n*****************************************************************************\n"
898  << COLOR_DEFAULT;
899  }
900 }
901 
902 std::vector<unsigned int>
904 {
905  std::vector<unsigned int> new_ids(_max_curr_grain_id - _old_max_grain_id);
906  auto new_id = _old_max_grain_id + 1;
907 
908  // Generate the new ids
909  std::iota(new_ids.begin(), new_ids.end(), new_id);
910 
911  return new_ids;
912 }
913 
914 void
916 {
917  // Don't remap grains if the current simulation step is before the specified tracking step
918  if (_t_step < _tracking_step)
919  return;
920 
921  TIME_SECTION(_remap_timer);
922 
923  if (_verbosity_level > 1)
924  _console << "Running remap Grains\n" << std::endl;
925 
932  std::map<unsigned int, std::size_t> grain_id_to_new_var;
933 
934  // Items are added to this list when split grains are found
935  std::list<std::pair<std::size_t, std::size_t>> split_pairs;
936 
945  if (_is_master)
946  {
947  // Build the map to detect difference in _var_index mappings after the remap operation
948  std::map<unsigned int, std::size_t> grain_id_to_existing_var_index;
949  for (auto & grain : _feature_sets)
950  {
951  // Unmark the grain so it can be used in the remap loop
952  grain._status = Status::CLEAR;
953 
954  grain_id_to_existing_var_index[grain._id] = grain._var_index;
955  }
956 
957  // Make sure that all split pieces of any grain are on the same OP
958  for (MooseIndex(_feature_sets) i = 0; i < _feature_sets.size(); ++i)
959  {
960  auto & grain1 = _feature_sets[i];
961 
962  for (MooseIndex(_feature_sets) j = 0; j < _feature_sets.size(); ++j)
963  {
964  auto & grain2 = _feature_sets[j];
965  if (i == j)
966  continue;
967 
968  // The first condition below is there to prevent symmetric checks (duplicate values)
969  if (i < j && grain1._id == grain2._id)
970  {
971  split_pairs.push_front(std::make_pair(i, j));
972  if (grain1._var_index != grain2._var_index)
973  {
974  if (_verbosity_level > 0)
975  _console << COLOR_YELLOW << "Split Grain (#" << grain1._id
976  << ") detected on unmatched OPs (" << grain1._var_index << ", "
977  << grain2._var_index << ") attempting to remap to " << grain1._var_index
978  << ".\n"
979  << COLOR_DEFAULT;
980 
985  grain1._var_index = grain2._var_index;
986  grain1._status |= Status::DIRTY;
987  }
988  }
989  }
990  }
991 
995  bool any_grains_remapped = false;
996  bool grains_remapped;
997 
998  std::set<unsigned int> notify_ids;
999  do
1000  {
1001  grains_remapped = false;
1002  notify_ids.clear();
1003 
1004  for (auto & grain1 : _feature_sets)
1005  {
1006  // We need to remap any grains represented on any variable index above the cuttoff
1007  if (grain1._var_index >= _reserve_op_index)
1008  {
1009  if (_verbosity_level > 0)
1010  _console << COLOR_YELLOW << "\nGrain #" << grain1._id
1011  << " detected on a reserved order parameter #" << grain1._var_index
1012  << ", remapping to another variable\n"
1013  << COLOR_DEFAULT;
1014 
1015  for (MooseIndex(_max_remap_recursion_depth) max = 0; max <= _max_remap_recursion_depth;
1016  ++max)
1017  if (max < _max_remap_recursion_depth)
1018  {
1019  if (attemptGrainRenumber(grain1, 0, max))
1020  break;
1021  }
1022  else if (!attemptGrainRenumber(grain1, 0, max))
1023  {
1024  _console << std::flush;
1025  std::stringstream oss;
1026  oss << "Unable to find any suitable order parameters for remapping while working "
1027  << "with Grain #" << grain1._id << ", which is on a reserve order parameter.\n"
1028  << "\n\nPossible Resolutions:\n"
1029  << "\t- Add more order parameters to your simulation (8 for 2D, 28 for 3D)\n"
1030  << "\t- Increase adaptivity or reduce your grain boundary widths\n"
1031  << "\t- Make sure you are not starting with too many grains for the mesh size\n";
1032  mooseError(oss.str());
1033  }
1034 
1035  grains_remapped = true;
1036  }
1037 
1038  for (auto & grain2 : _feature_sets)
1039  {
1040  // Don't compare a grain with itself and don't try to remap inactive grains
1041  if (&grain1 == &grain2)
1042  continue;
1043 
1044  if (grain1._var_index == grain2._var_index && // grains represented by same variable?
1045  grain1._id != grain2._id && // are they part of different grains?
1046  grain1.boundingBoxesIntersect(grain2) && // do bboxes intersect (coarse level)?
1047  grain1.halosIntersect(grain2)) // do they actually overlap (fine level)?
1048  {
1049  if (_verbosity_level > 0)
1050  _console << COLOR_YELLOW << "Grain #" << grain1._id << " intersects Grain #"
1051  << grain2._id << " (variable index: " << grain1._var_index << ")\n"
1052  << COLOR_DEFAULT;
1053 
1054  for (MooseIndex(_max_remap_recursion_depth) max = 0; max <= _max_remap_recursion_depth;
1055  ++max)
1056  {
1057  if (max < _max_remap_recursion_depth)
1058  {
1059  if (attemptGrainRenumber(grain1, 0, max))
1060  {
1061  grains_remapped = true;
1062  break;
1063  }
1064  }
1065  else if (!attemptGrainRenumber(grain1, 0, max) &&
1066  !attemptGrainRenumber(grain2, 0, max))
1067  {
1068  notify_ids.insert(grain1._id);
1069  notify_ids.insert(grain2._id);
1070  }
1071  }
1072  }
1073  }
1074  }
1075  any_grains_remapped |= grains_remapped;
1076  } while (grains_remapped);
1077 
1078  if (!notify_ids.empty())
1079  {
1080  _console << std::flush;
1081  std::stringstream oss;
1082  oss << "Unable to find any suitable order parameters for remapping while working "
1083  << "with the following grain IDs:\n"
1084  << Moose::stringify(notify_ids, ", ", "", true) << "\n\nPossible Resolutions:\n"
1085  << "\t- Add more order parameters to your simulation (8 for 2D, 28 for 3D)\n"
1086  << "\t- Increase adaptivity or reduce your grain boundary widths\n"
1087  << "\t- Make sure you are not starting with too many grains for the mesh size\n";
1088 
1089  if (_tolerate_failure)
1090  mooseWarning(oss.str());
1091  else
1092  mooseError(oss.str());
1093  }
1094 
1095  // Verify that split grains are still intact
1096  for (auto & split_pair : split_pairs)
1097  if (_feature_sets[split_pair.first]._var_index != _feature_sets[split_pair.first]._var_index)
1098  mooseError("Split grain remapped - This case is currently not handled");
1099 
1105  for (auto & grain : _feature_sets)
1106  {
1107  mooseAssert(grain_id_to_existing_var_index.find(grain._id) !=
1108  grain_id_to_existing_var_index.end(),
1109  "Missing unique ID");
1110 
1111  auto old_var_index = grain_id_to_existing_var_index[grain._id];
1112 
1113  if (old_var_index != grain._var_index)
1114  {
1115  mooseAssert(static_cast<bool>(grain._status & Status::DIRTY), "grain status is incorrect");
1116 
1117  grain_id_to_new_var.emplace_hint(
1118  grain_id_to_new_var.end(),
1119  std::pair<unsigned int, std::size_t>(grain._id, grain._var_index));
1120 
1129  grain._var_index = old_var_index;
1130  // Clear the DIRTY status as well for consistency
1131  grain._status &= ~Status::DIRTY;
1132  }
1133  }
1134 
1135  if (!grain_id_to_new_var.empty())
1136  {
1137  if (_verbosity_level > 1)
1138  {
1139  _console << "Final remapping tally:\n";
1140  for (const auto & remap_pair : grain_id_to_new_var)
1141  _console << "Grain #" << remap_pair.first << " var_index "
1142  << grain_id_to_existing_var_index[remap_pair.first] << " -> "
1143  << remap_pair.second << '\n';
1144  _console << "Communicating swaps with remaining processors..." << std::endl;
1145  }
1146  }
1147  } // root processor
1148 
1149  // Communicate the std::map to all ranks
1150  _communicator.broadcast(grain_id_to_new_var);
1151 
1152  // Perform swaps if any occurred
1153  if (!grain_id_to_new_var.empty())
1154  {
1155  // Cache for holding values during swaps
1156  std::vector<std::map<Node *, CacheValues>> cache(_n_vars);
1157 
1158  // Perform the actual swaps on all processors
1159  for (auto & grain : _feature_sets)
1160  {
1161  // See if this grain was remapped
1162  auto new_var_it = grain_id_to_new_var.find(grain._id);
1163  if (new_var_it != grain_id_to_new_var.end())
1164  swapSolutionValues(grain, new_var_it->second, cache, RemapCacheMode::FILL);
1165  }
1166 
1167  for (auto & grain : _feature_sets)
1168  {
1169  // See if this grain was remapped
1170  auto new_var_it = grain_id_to_new_var.find(grain._id);
1171  if (new_var_it != grain_id_to_new_var.end())
1172  swapSolutionValues(grain, new_var_it->second, cache, RemapCacheMode::USE);
1173  }
1174 
1175  _nl.solution().close();
1176  _nl.solutionOld().close();
1177  _nl.solutionOlder().close();
1178 
1179  _fe_problem.getNonlinearSystemBase().system().update();
1180 
1181  if (_verbosity_level > 1)
1182  _console << "Swaps complete" << std::endl;
1183  }
1184 }
1186 void
1188  std::vector<std::list<GrainDistance>> & min_distances)
1189 {
1213  for (MooseIndex(_feature_sets) i = 0; i < _feature_sets.size(); ++i)
1214  {
1215  auto & other_grain = _feature_sets[i];
1216 
1217  if (other_grain._var_index == grain._var_index || other_grain._var_index >= _reserve_op_index)
1218  continue;
1219 
1220  auto target_var_index = other_grain._var_index;
1221  auto target_grain_index = i;
1222  auto target_grain_id = other_grain._id;
1223 
1224  Real curr_bbox_diff = boundingRegionDistance(grain._bboxes, other_grain._bboxes);
1225 
1226  GrainDistance grain_distance_obj(
1227  curr_bbox_diff, target_var_index, target_grain_index, target_grain_id);
1228 
1229  // To handle touching halos we penalize the top pick each time we see another
1230  if (curr_bbox_diff == -1.0 && !min_distances[target_var_index].empty())
1231  {
1232  Real last_distance = min_distances[target_var_index].begin()->_distance;
1233  if (last_distance < 0)
1234  grain_distance_obj._distance += last_distance;
1235  }
1236 
1237  // Insertion sort into a list
1238  auto insert_it = min_distances[target_var_index].begin();
1239  while (insert_it != min_distances[target_var_index].end() && !(grain_distance_obj < *insert_it))
1240  ++insert_it;
1241  min_distances[target_var_index].insert(insert_it, grain_distance_obj);
1242  }
1243 
1249  for (MooseIndex(_vars) var_index = 0; var_index < _reserve_op_index; ++var_index)
1250  {
1251  // Don't put an entry in for matching variable indices (i.e. we can't remap to ourselves)
1252  if (grain._var_index == var_index)
1253  continue;
1254 
1255  if (min_distances[var_index].empty())
1256  min_distances[var_index].emplace_front(std::numeric_limits<Real>::max(), var_index);
1257  }
1258 }
1260 bool
1261 GrainTracker::attemptGrainRenumber(FeatureData & grain, unsigned int depth, unsigned int max_depth)
1262 {
1263  // End the recursion of our breadth first search
1264  if (depth > max_depth)
1265  return false;
1266 
1267  std::size_t curr_var_index = grain._var_index;
1268 
1269  std::vector<std::map<Node *, CacheValues>> cache;
1270 
1271  std::vector<std::list<GrainDistance>> min_distances(_vars.size());
1272 
1278  computeMinDistancesFromGrain(grain, min_distances);
1279 
1290  // clang-format off
1291  std::sort(min_distances.begin(), min_distances.end(),
1292  [](const std::list<GrainDistance> & lhs, const std::list<GrainDistance> & rhs)
1293  {
1294  // Sort lists in reverse order (largest distance first)
1295  // These empty cases are here to make this comparison stable
1296  if (lhs.empty())
1297  return false;
1298  else if (rhs.empty())
1299  return true;
1300  else
1301  return lhs.begin()->_distance > rhs.begin()->_distance;
1302  });
1303  // clang-format on
1304 
1305  for (auto & list_ref : min_distances)
1306  {
1307  const auto target_it = list_ref.begin();
1308  if (target_it == list_ref.end())
1309  continue;
1310 
1311  // If the distance is positive we can just remap and be done
1312  if (target_it->_distance > 0)
1313  {
1314  if (_verbosity_level > 0)
1315  {
1316  _console << COLOR_GREEN << "- Depth " << depth << ": Remapping grain #" << grain._id
1317  << " from variable index " << curr_var_index << " to " << target_it->_var_index;
1318  if (target_it->_distance == std::numeric_limits<Real>::max())
1319  _console << " which currently contains zero grains.\n\n" << COLOR_DEFAULT;
1320  else
1321  _console << " whose closest grain (#" << target_it->_grain_id << ") is at a distance of "
1322  << std::sqrt(target_it->_distance) << "\n\n"
1323  << COLOR_DEFAULT;
1324  }
1325 
1326  grain._status |= Status::DIRTY;
1327  grain._var_index = target_it->_var_index;
1328  return true;
1329  }
1330 
1331  // If the distance isn't positive we just need to make sure that none of the grains represented
1332  // by the target variable index would intersect this one if we were to remap
1333  {
1334  auto next_target_it = target_it;
1335  bool intersection_hit = false;
1336  unsigned short num_close_targets = 0;
1337  std::ostringstream oss;
1338  while (!intersection_hit && next_target_it != list_ref.end())
1339  {
1340  if (next_target_it->_distance > 0)
1341  break;
1342 
1343  mooseAssert(next_target_it->_grain_index < _feature_sets.size(),
1344  "Error in indexing target grain in attemptGrainRenumber");
1345  FeatureData & next_target_grain = _feature_sets[next_target_it->_grain_index];
1346 
1347  // If any grains touch we're done here
1348  if (grain.halosIntersect(next_target_grain))
1349  intersection_hit = true;
1350  else
1351  {
1352  if (num_close_targets > 0)
1353  oss << ", "; // delimiter
1354  oss << "#" << next_target_it->_grain_id;
1355  }
1356 
1357  ++next_target_it;
1358  ++num_close_targets;
1359  }
1360 
1361  if (!intersection_hit)
1362  {
1363  if (_verbosity_level > 0)
1364  {
1365  _console << COLOR_GREEN << "- Depth " << depth << ": Remapping grain #" << grain._id
1366  << " from variable index " << curr_var_index << " to " << target_it->_var_index;
1367 
1368  if (num_close_targets == 1)
1369  _console << " whose closest grain (" << oss.str()
1370  << ") is inside our bounding box but whose halo is not touching.\n\n"
1371  << COLOR_DEFAULT;
1372  else
1373  _console << " whose closest grains (" << oss.str()
1374  << ") are inside our bounding box but whose halos are not touching.\n\n"
1375  << COLOR_DEFAULT;
1376  }
1377 
1378  grain._status |= Status::DIRTY;
1379  grain._var_index = target_it->_var_index;
1380  return true;
1381  }
1382  }
1383 
1384  // If we reach this part of the loop, there is no simple renumbering that can be done.
1385  mooseAssert(target_it->_grain_index < _feature_sets.size(),
1386  "Error in indexing target grain in attemptGrainRenumber");
1387  FeatureData & target_grain = _feature_sets[target_it->_grain_index];
1388 
1396  if (target_it->_distance < -1)
1397  return false;
1398 
1399  // Make sure this grain isn't marked. If it is, we can't recurse here
1400  if ((target_grain._status & Status::MARKED) == Status::MARKED)
1401  return false;
1402 
1408  grain._var_index = target_it->_var_index;
1409  grain._status |= Status::MARKED;
1410  if (attemptGrainRenumber(target_grain, depth + 1, max_depth))
1411  {
1412  // SUCCESS!
1413  if (_verbosity_level > 0)
1414  _console << COLOR_GREEN << "- Depth " << depth << ": Remapping grain #" << grain._id
1415  << " from variable index " << curr_var_index << " to " << target_it->_var_index
1416  << "\n\n"
1417  << COLOR_DEFAULT;
1418 
1419  // Now we need to mark the grain as DIRTY since the recursion succeeded.
1420  grain._status |= Status::DIRTY;
1421  return true;
1422  }
1423  else
1424  // FAILURE, We need to set our var index back after failed recursive step
1425  grain._var_index = curr_var_index;
1426 
1427  // ALWAYS "unmark" (or clear the MARKED status) after recursion so it can be used by other remap
1428  // operations
1429  grain._status &= ~Status::MARKED;
1430  }
1431 
1432  return false;
1433 }
1435 void
1437  std::size_t new_var_index,
1438  std::vector<std::map<Node *, CacheValues>> & cache,
1439  RemapCacheMode cache_mode)
1440 {
1441  MeshBase & mesh = _mesh.getMesh();
1442 
1443  // Remap the grain
1444  std::set<Node *> updated_nodes_tmp; // Used only in the elemental case
1445  for (auto entity : grain._local_ids)
1446  {
1447  if (_is_elemental)
1448  {
1449  Elem * elem = mesh.query_elem_ptr(entity);
1450  if (!elem)
1451  continue;
1452 
1453  for (unsigned int i = 0; i < elem->n_nodes(); ++i)
1454  {
1455  Node * curr_node = elem->node_ptr(i);
1456  if (updated_nodes_tmp.find(curr_node) == updated_nodes_tmp.end())
1457  {
1458  // cache this node so we don't attempt to remap it again within this loop
1459  updated_nodes_tmp.insert(curr_node);
1460  swapSolutionValuesHelper(curr_node, grain._var_index, new_var_index, cache, cache_mode);
1461  }
1462  }
1463  }
1464  else
1466  mesh.query_node_ptr(entity), grain._var_index, new_var_index, cache, cache_mode);
1467  }
1468 
1469  // Update the variable index in the unique grain datastructure after swaps are complete
1470  if (cache_mode == RemapCacheMode::USE || cache_mode == RemapCacheMode::BYPASS)
1471  grain._var_index = new_var_index;
1472 }
1474 void
1476  std::size_t curr_var_index,
1477  std::size_t new_var_index,
1478  std::vector<std::map<Node *, CacheValues>> & cache,
1479  RemapCacheMode cache_mode)
1480 {
1481  if (curr_node && curr_node->processor_id() == processor_id())
1482  {
1483  // Reinit the node so we can get and set values of the solution here
1484  _subproblem.reinitNode(curr_node, 0);
1485 
1486  // Local variables to hold values being transferred
1487  Real current, old = 0, older = 0;
1488  // Retrieve the value either from the old variable or cache
1489  if (cache_mode == RemapCacheMode::FILL || cache_mode == RemapCacheMode::BYPASS)
1490  {
1491  current = _vars[curr_var_index]->dofValues()[0];
1492  if (_is_transient)
1493  {
1494  old = _vars[curr_var_index]->dofValuesOld()[0];
1495  older = _vars[curr_var_index]->dofValuesOlder()[0];
1496  }
1497  }
1498  else // USE
1499  {
1500  const auto cache_it = cache[curr_var_index].find(curr_node);
1501  mooseAssert(cache_it != cache[curr_var_index].end(), "Error in cache");
1502  current = cache_it->second.current;
1503  old = cache_it->second.old;
1504  older = cache_it->second.older;
1505  }
1506 
1507  // Cache the value or use it!
1508  if (cache_mode == RemapCacheMode::FILL)
1509  {
1510  cache[curr_var_index][curr_node].current = current;
1511  cache[curr_var_index][curr_node].old = old;
1512  cache[curr_var_index][curr_node].older = older;
1513  }
1514  else // USE or BYPASS
1515  {
1516  const auto & dof_index = _vars[new_var_index]->nodalDofIndex();
1517 
1518  // Transfer this solution from the old to the current
1519  _nl.solution().set(dof_index, current);
1520  if (_is_transient)
1521  {
1522  _nl.solutionOld().set(dof_index, old);
1523  _nl.solutionOlder().set(dof_index, older);
1524  }
1525  }
1526 
1539  if (cache_mode == RemapCacheMode::FILL || cache_mode == RemapCacheMode::BYPASS)
1540  {
1541  const auto & dof_index = _vars[curr_var_index]->nodalDofIndex();
1542 
1543  // Set the DOF for the current variable to zero
1544  _nl.solution().set(dof_index, 0.0);
1545  if (_is_transient)
1546  {
1547  _nl.solutionOld().set(dof_index, 0.0);
1548  _nl.solutionOlder().set(dof_index, 0.0);
1549  }
1550  }
1551  }
1552 }
1554 void
1556 {
1557  TIME_SECTION(_update_field_info);
1558 
1559  for (MooseIndex(_maps_size) map_num = 0; map_num < _maps_size; ++map_num)
1560  _feature_maps[map_num].clear();
1561 
1562  std::map<dof_id_type, Real> tmp_map;
1563 
1564  for (const auto & grain : _feature_sets)
1565  {
1566  std::size_t curr_var = grain._var_index;
1567  std::size_t map_index = (_single_map_mode || _condense_map_info) ? 0 : curr_var;
1568 
1569  for (auto entity : grain._local_ids)
1570  {
1571  // Highest variable value at this entity wins
1572  Real entity_value = std::numeric_limits<Real>::lowest();
1573  if (_is_elemental)
1574  {
1575  const Elem * elem = _mesh.elemPtr(entity);
1576  std::vector<Point> centroid(1, elem->centroid());
1577  if (_poly_ic_uo && _first_time)
1578  {
1579  entity_value = _poly_ic_uo->getVariableValue(grain._var_index, centroid[0]);
1580  }
1581  else
1582  {
1583  _fe_problem.reinitElemPhys(elem, centroid, 0, /* suppress_displaced_init = */ true);
1584  entity_value = _vars[curr_var]->sln()[0];
1585  }
1586  }
1587  else
1588  {
1589  auto node_ptr = _mesh.nodePtr(entity);
1590  entity_value = _vars[curr_var]->getNodalValue(*node_ptr);
1591  }
1592 
1593  if (entity_value != std::numeric_limits<Real>::lowest() &&
1594  (tmp_map.find(entity) == tmp_map.end() || entity_value > tmp_map[entity]))
1595  {
1596  mooseAssert(grain._id != invalid_id, "Missing Grain ID");
1597  _feature_maps[map_index][entity] = grain._id;
1598 
1599  if (_var_index_mode)
1600  _var_index_maps[map_index][entity] = grain._var_index;
1601 
1602  tmp_map[entity] = entity_value;
1603  }
1604 
1606  {
1607  auto insert_pair = moose_try_emplace(
1608  _entity_var_to_features, entity, std::vector<unsigned int>(_n_vars, invalid_id));
1609  auto & vec_ref = insert_pair.first->second;
1610 
1611  if (insert_pair.second)
1612  {
1613  // insert the reserve op numbers (if appropriate)
1614  for (MooseIndex(_n_reserve_ops) reserve_index = 0; reserve_index < _n_reserve_ops;
1615  ++reserve_index)
1616  vec_ref[reserve_index] = _reserve_grain_first_index + reserve_index;
1617  }
1618  vec_ref[grain._var_index] = grain._id;
1619  }
1620  }
1621 
1622  if (_compute_halo_maps)
1623  for (auto entity : grain._halo_ids)
1624  _halo_ids[grain._var_index][entity] = grain._var_index;
1625 
1626  for (auto entity : grain._ghosted_ids)
1627  _ghosted_entity_ids[entity] = 1;
1628  }
1629 
1631 }
1633 void
1635 {
1636  if (_compute_halo_maps)
1637  {
1638  // rank var_index entity_id
1639  std::vector<std::pair<std::size_t, dof_id_type>> halo_ids_all;
1640 
1641  std::vector<int> counts;
1642  std::vector<std::pair<std::size_t, dof_id_type>> local_halo_ids;
1643  std::size_t counter = 0;
1644 
1645  const bool isDistributedMesh = _mesh.isDistributedMesh();
1646 
1647  if (_is_master)
1648  {
1649  std::vector<std::vector<std::pair<std::size_t, dof_id_type>>> root_halo_ids(_n_procs);
1650  counts.resize(_n_procs);
1651 
1652  // Loop over the _halo_ids "field" and build minimal lists for all of the other ranks
1653  for (MooseIndex(_halo_ids) var_index = 0; var_index < _halo_ids.size(); ++var_index)
1654  {
1655  for (const auto & entity_pair : _halo_ids[var_index])
1656  {
1657  auto entity_id = entity_pair.first;
1658  if (isDistributedMesh)
1659  {
1660  // Check to see which contiguous range this entity ID falls into
1661  auto range_it =
1662  std::lower_bound(_all_ranges.begin(),
1663  _all_ranges.end(),
1664  entity_id,
1665  [](const std::pair<dof_id_type, dof_id_type> range,
1666  dof_id_type entity_id) { return range.second < entity_id; });
1667 
1668  mooseAssert(range_it != _all_ranges.end(), "No range round?");
1669 
1670  // Recover the index from the iterator
1671  auto proc_id = std::distance(_all_ranges.begin(), range_it);
1672 
1673  // Now add this halo entity to the map for the corresponding proc to scatter latter
1674  root_halo_ids[proc_id].push_back(std::make_pair(var_index, entity_id));
1675  }
1676  else
1677  {
1678  DofObject * halo_entity;
1679  if (_is_elemental)
1680  halo_entity = _mesh.queryElemPtr(entity_id);
1681  else
1682  halo_entity = _mesh.queryNodePtr(entity_id);
1683 
1684  if (halo_entity)
1685  root_halo_ids[halo_entity->processor_id()].push_back(
1686  std::make_pair(var_index, entity_id));
1687  }
1688  }
1689  }
1690 
1691  // Build up the counts vector for MPI scatter
1692  std::size_t global_count = 0;
1693  for (const auto & vector_ref : root_halo_ids)
1694  {
1695  std::copy(vector_ref.begin(), vector_ref.end(), std::back_inserter(halo_ids_all));
1696  counts[counter] = vector_ref.size();
1697  global_count += counts[counter++];
1698  }
1699  }
1700 
1701  _communicator.scatter(halo_ids_all, counts, local_halo_ids);
1702 
1703  // Now add the contributions from the root process to the processor local maps
1704  for (const auto & halo_pair : local_halo_ids)
1705  _halo_ids[halo_pair.first].emplace(std::make_pair(halo_pair.second, halo_pair.first));
1706 
1713  for (const auto & grain : _feature_sets)
1714  for (auto local_id : grain._local_ids)
1715  _halo_ids[grain._var_index].erase(local_id);
1716  }
1717 }
1719 Real
1720 GrainTracker::centroidRegionDistance(std::vector<BoundingBox> & bboxes1,
1721  std::vector<BoundingBox> & bboxes2) const
1722 {
1726  auto min_distance = std::numeric_limits<Real>::max();
1727  for (const auto & bbox1 : bboxes1)
1728  {
1729  const auto centroid_point1 = (bbox1.max() + bbox1.min()) / 2.0;
1730 
1731  for (const auto & bbox2 : bboxes2)
1732  {
1733  const auto centroid_point2 = (bbox2.max() + bbox2.min()) / 2.0;
1734 
1735  // Here we'll calculate a distance between the centroids
1736  auto curr_distance = _mesh.minPeriodicDistance(_var_number, centroid_point1, centroid_point2);
1737 
1738  if (curr_distance < min_distance)
1739  min_distance = curr_distance;
1740  }
1741  }
1742 
1743  return min_distance;
1744 }
1746 Real
1747 GrainTracker::boundingRegionDistance(std::vector<BoundingBox> & bboxes1,
1748  std::vector<BoundingBox> & bboxes2) const
1749 {
1756  auto min_distance = std::numeric_limits<Real>::max();
1757  for (const auto & bbox1 : bboxes1)
1758  {
1759  for (const auto & bbox2 : bboxes2)
1760  {
1761  // AABB squared distance
1762  Real curr_distance = 0.0;
1763  bool boxes_overlap = true;
1764  for (unsigned int dim = 0; dim < LIBMESH_DIM; ++dim)
1765  {
1766  const auto & min1 = bbox1.min()(dim);
1767  const auto & max1 = bbox1.max()(dim);
1768  const auto & min2 = bbox2.min()(dim);
1769  const auto & max2 = bbox2.max()(dim);
1770 
1771  if (min1 > max2)
1772  {
1773  const auto delta = max2 - min1;
1774  curr_distance += delta * delta;
1775  boxes_overlap = false;
1776  }
1777  else if (min2 > max1)
1778  {
1779  const auto delta = max1 - min2;
1780  curr_distance += delta * delta;
1781  boxes_overlap = false;
1782  }
1783  }
1784 
1785  if (boxes_overlap)
1786  return -1.0; /* all overlaps are treated the same */
1787 
1788  if (curr_distance < min_distance)
1789  min_distance = curr_distance;
1790  }
1791  }
1792 
1793  return min_distance;
1794 }
1796 unsigned int
1798 {
1807  _reserve_grain_first_index + _n_reserve_ops /* no +1 here!*/);
1808 
1809  return _max_curr_grain_id;
1810 }
1811 
1812 /*************************************************
1813  ************** Helper Structures ****************
1814  ************************************************/
1815 GrainDistance::GrainDistance(Real distance, std::size_t var_index)
1816  : GrainDistance(distance,
1817  var_index,
1818  std::numeric_limits<std::size_t>::max(),
1819  std::numeric_limits<unsigned int>::max())
1821 }
1822 
1823 GrainDistance::GrainDistance(Real distance,
1824  std::size_t var_index,
1825  std::size_t grain_index,
1826  unsigned int grain_id)
1827  : _distance(distance), _var_index(var_index), _grain_index(grain_index), _grain_id(grain_id)
1828 {
1830 
1831 bool
1832 GrainDistance::operator<(const GrainDistance & rhs) const
1833 {
1834  return _distance < rhs._distance;
1835 }
FeatureFloodCount::expandEdgeHalos
void expandEdgeHalos(unsigned int num_layers_to_expand)
This method expands the existing halo set by some width determined by the passed in value.
Definition: FeatureFloodCount.C:1527
FeatureFloodCount::_var_index_mode
const bool _var_index_mode
This variable is used to indicate whether the maps will contain unique region information or just the...
Definition: FeatureFloodCount.h:601
GrainTracker::_reserve_grain_first_index
unsigned int _reserve_grain_first_index
Holds the first unique grain index when using _reserve_op (all the remaining indices are sequential)
Definition: GrainTracker.h:236
FeatureFloodCount::_condense_map_info
const bool _condense_map_info
Definition: FeatureFloodCount.h:593
GrainTracker::getEntityValue
virtual Real getEntityValue(dof_id_type node_id, FieldType field_type, std::size_t var_index=0) const override
Definition: GrainTracker.C:120
GrainTrackerInterface
This class defines the interface for the GrainTracking objects.
Definition: GrainTrackerInterface.h:24
validParams< FeatureFloodCount >
InputParameters validParams< FeatureFloodCount >()
Definition: FeatureFloodCount.C:104
GrainTracker::_halo_level
const unsigned short _halo_level
The thickness of the halo surrounding each grain.
Definition: GrainTracker.h:183
GrainTracker::assignGrains
void assignGrains()
When the tracking phase starts (_t_step == _tracking_step) it assigns a unique id to every FeatureDat...
Definition: GrainTracker.C:450
GrainTracker::_error_on_grain_creation
const bool _error_on_grain_creation
Boolean to terminate with an error if a new grain is created during the simulation.
Definition: GrainTracker.h:232
FeatureFloodCount::invalid_size_t
static const std::size_t invalid_size_t
Definition: FeatureFloodCount.h:93
GrainTracker::PartialFeatureData::id
unsigned int id
Definition: GrainTracker.h:41
FeatureFloodCount::scatterAndUpdateRanks
void scatterAndUpdateRanks()
Calls buildLocalToGlobalIndices to build the individual local to global indicies for each rank and sc...
Definition: FeatureFloodCount.C:717
FeatureFloodCount::_mesh
MooseMesh & _mesh
A reference to the mesh.
Definition: FeatureFloodCount.h:581
FeatureFloodCount::FeatureData::_id
unsigned int _id
An ID for this feature.
Definition: FeatureFloodCount.h:287
FeatureFloodCount::getEntityValue
virtual Real getEntityValue(dof_id_type entity_id, FieldType field_type, std::size_t var_index=0) const
Definition: FeatureFloodCount.C:918
FeatureFloodCount::BoundaryIntersection::SECONDARY_PERCOLATION_BOUNDARY
FeatureFloodCount::getVarToFeatureVector
virtual const std::vector< unsigned int > & getVarToFeatureVector(dof_id_type elem_id) const
Returns a list of active unique feature ids for a particular element.
Definition: FeatureFloodCount.C:701
GrainTracker::_feature_sets_old
std::vector< FeatureData > _feature_sets_old
This data structure holds the map of unique grains from the previous time step.
Definition: GrainTracker.h:211
GrainTracker::_old_max_grain_id
unsigned int _old_max_grain_id
The previous max grain id (needed to figure out which ids are new in a given step)
Definition: GrainTracker.h:239
GrainTracker::getTotalFeatureCount
virtual std::size_t getTotalFeatureCount() const override
Returns the total feature count (active and inactive ids, useful for sizing vectors)
Definition: GrainTracker.C:150
FeatureFloodCount::Status
Status
This enumeration is used to indicate status of the grains in the _unique_grains data structure.
Definition: FeatureFloodCount.h:120
FeatureFloodCount::communicateAndMerge
void communicateAndMerge()
This routine handles all of the serialization, communication and deserialization of the data structur...
Definition: FeatureFloodCount.C:412
GrainTracker::newGrainCreated
virtual void newGrainCreated(unsigned int new_grain_id)
This method is called when a new grain is detected.
Definition: GrainTracker.C:880
GrainDistance::_distance
Real _distance
Definition: GrainTracker.h:282
FeatureFloodCount::FieldType
FieldType
Definition: FeatureFloodCount.h:103
GrainTracker::_n_reserve_ops
const unsigned short _n_reserve_ops
The number of reserved order parameters.
Definition: GrainTracker.h:189
GrainTracker::boundingRegionDistance
Real boundingRegionDistance(std::vector< BoundingBox > &bboxes1, std::vector< BoundingBox > &bboxes2) const
This method returns the minimum periodic distance between two vectors of bounding boxes.
Definition: GrainTracker.C:1745
FeatureFloodCount::_is_master
const bool _is_master
Convenience variable for testing master rank.
Definition: FeatureFloodCount.h:732
GrainTracker::updateFieldInfo
virtual void updateFieldInfo() override
This method is used to populate any of the data structures used for storing field data (nodal or elem...
Definition: GrainTracker.C:1553
GrainTracker::doesFeatureIntersectBoundary
virtual bool doesFeatureIntersectBoundary(unsigned int feature_id) const override
Returns a Boolean indicating whether this feature intersects any boundary.
Definition: GrainTracker.C:175
GrainTracker::_broadcast_update
const PerfID _broadcast_update
Definition: GrainTracker.h:254
GrainDistance
This struct is used to hold distance information to other grains in the simulation.
Definition: GrainTracker.h:263
FeatureFloodCount
This object will mark nodes or elements of continuous regions all with a unique number for the purpos...
Definition: FeatureFloodCount.h:44
FeatureFloodCount::initialize
virtual void initialize() override
Definition: FeatureFloodCount.C:298
GrainTracker
Definition: GrainTracker.h:24
FeatureFloodCount::Status::CLEAR
GrainTracker::RemapCacheMode::USE
FeatureFloodCount::BoundaryIntersection::SPECIFIED_BOUNDARY
counter
static unsigned int counter
Definition: ContactPenetrationAuxAction.C:17
validParams< GrainTrackerInterface >
InputParameters validParams< GrainTrackerInterface >()
Definition: GrainTrackerInterface.C:15
FeatureFloodCount::sortAndLabel
void sortAndLabel()
Sort and assign ids to features based on their position in the container after sorting.
Definition: FeatureFloodCount.C:573
GrainTracker::PartialFeatureData
Definition: GrainTracker.h:38
GrainTracker::getNewGrainIDs
virtual std::vector< unsigned int > getNewGrainIDs() const override
This method returns all of the new ids generated in an invocation of the GrainTracker.
Definition: GrainTracker.C:901
GrainTracker::broadcastAndUpdateGrainData
void broadcastAndUpdateGrainData()
Broadcast essential Grain information to all processors.
Definition: GrainTracker.C:388
GrainTracker::PartialFeatureData::status
Status status
Definition: GrainTracker.h:43
GrainTracker::_poly_ic_uo
const PolycrystalUserObjectBase * _poly_ic_uo
An optional IC UserObject which can provide initial data structures to this object.
Definition: GrainTracker.h:214
PolycrystalUserObjectBase
This object provides the base capability for creating proper polycrystal ICs.
Definition: PolycrystalUserObjectBase.h:27
GrainTracker::prepopulateState
void prepopulateState(const FeatureFloodCount &ffc_object)
This method extracts the necessary state from the passed in object necessary to continue tracking gra...
Definition: GrainTracker.C:298
GrainTracker::execute
virtual void execute() override
Definition: GrainTracker.C:273
GrainTracker::computeMinDistancesFromGrain
void computeMinDistancesFromGrain(FeatureData &grain, std::vector< std::list< GrainDistance >> &min_distances)
Populates and sorts a min_distances vector with the minimum distances to all grains in the simulation...
Definition: GrainTracker.C:1185
validParams< GrainTracker >
InputParameters validParams< GrainTracker >()
Definition: GrainTracker.C:50
GrainTracker::swapSolutionValues
void swapSolutionValues(FeatureData &grain, std::size_t new_var_index, std::vector< std::map< Node *, CacheValues >> &cache, RemapCacheMode cache_mode)
A routine for moving all of the solution values from a given grain to a new variable number.
Definition: GrainTracker.C:1434
FeatureFloodCount::getFeatures
const std::vector< FeatureData > & getFeatures() const
Return a constant reference to the vector of all discovered features.
Definition: FeatureFloodCount.h:340
FeatureFloodCount::_var_number
unsigned long _var_number
This variable is used to build the periodic node map.
Definition: FeatureFloodCount.h:588
GrainDistance::GrainDistance
GrainDistance(Real distance, std::size_t var_index)
Definition: GrainTracker.C:1812
FeatureFloodCount::_feature_sets
std::vector< FeatureData > & _feature_sets
The data structure used to hold the globally unique features.
Definition: FeatureFloodCount.h:662
GrainTracker::RemapCacheMode::BYPASS
FeatureFloodCount::_entity_var_to_features
std::map< dof_id_type, std::vector< unsigned int > > _entity_var_to_features
Definition: FeatureFloodCount.h:713
GrainTracker::swapSolutionValuesHelper
void swapSolutionValuesHelper(Node *curr_node, std::size_t curr_var_index, std::size_t new_var_index, std::vector< std::map< Node *, CacheValues >> &cache, RemapCacheMode cache_mode)
Helper method for actually performing the swaps.
Definition: GrainTracker.C:1473
GrainTracker::trackGrains
void trackGrains()
On subsequent time_steps, incoming FeatureData objects are compared to previous time_step information...
Definition: GrainTracker.C:497
GrainTracker::finalize
virtual void finalize() override
Definition: GrainTracker.C:330
FeatureFloodCount::getFeatureVar
virtual unsigned int getFeatureVar(unsigned int feature_id) const
Returns the variable representing the passed in feature.
Definition: FeatureFloodCount.C:809
GrainTracker::_max_curr_grain_id
unsigned int & _max_curr_grain_id
Holds the next "regular" grain ID (a grain found or remapped to the standard op vars)
Definition: GrainTracker.h:242
FeatureFloodCount::_single_map_mode
const bool _single_map_mode
This variable is used to indicate whether or not multiple maps are used during flooding.
Definition: FeatureFloodCount.h:591
GrainTracker::_finalize_timer
const PerfID _finalize_timer
Timers.
Definition: GrainTracker.h:251
GrainTracker::_all_ranges
std::vector< std::pair< dof_id_type, dof_id_type > > _all_ranges
Data structure to hold element ID ranges when using Distributed Mesh (populated on rank 0 only)
Definition: GrainTracker.h:248
PolycrystalUserObjectBase.h
GrainTracker::_remap
const bool _remap
Inidicates whether remapping should be done or not (remapping is independent of tracking)
Definition: GrainTracker.h:199
FeatureFloodCount::_ghosted_entity_ids
std::map< dof_id_type, int > _ghosted_entity_ids
The map for holding reconstructed ghosted element information.
Definition: FeatureFloodCount.h:695
FeatureFloodCount::_maps_size
const std::size_t _maps_size
Convenience variable holding the size of all the datastructures size by the number of maps.
Definition: FeatureFloodCount.h:620
GrainTracker::RemapCacheMode::FILL
FeatureFloodCount::_partial_feature_sets
std::vector< std::list< FeatureData > > _partial_feature_sets
The data structure used to hold partial and communicated feature data, during the discovery and mergi...
Definition: FeatureFloodCount.h:655
GrainTracker::_nl
NonlinearSystemBase & _nl
A reference to the nonlinear system (used for retrieving solution vectors)
Definition: GrainTracker.h:205
GrainTracker::centroidRegionDistance
Real centroidRegionDistance(std::vector< BoundingBox > &bboxes1, std::vector< BoundingBox > &bboxes2) const
This method returns the minimum periodic distance between the centroids of two vectors of bounding bo...
Definition: GrainTracker.C:1718
dataStore
void dataStore(std::ostream &stream, GrainTracker::PartialFeatureData &feature, void *context)
Definition: GrainTracker.C:28
FeatureFloodCount::_is_elemental
const bool _is_elemental
Determines if the flood counter is elements or not (nodes)
Definition: FeatureFloodCount.h:723
GrainTracker::~GrainTracker
virtual ~GrainTracker()
Definition: GrainTracker.C:117
GrainTracker::PartialFeatureData::centroid
Point centroid
Definition: GrainTracker.h:42
GrainTracker::attemptGrainRenumber
bool attemptGrainRenumber(FeatureData &grain, unsigned int depth, unsigned int max_depth)
This is the recursive part of the remapping algorithm.
Definition: GrainTracker.C:1259
GrainTracker::PartialFeatureData::boundary_intersection
BoundaryIntersection boundary_intersection
Definition: GrainTracker.h:40
FeatureFloodCount::FeatureData
Definition: FeatureFloodCount.h:138
FeatureFloodCount::_vars
std::vector< MooseVariable * > _vars
The vector of coupled in variables cast to MooseVariable.
Definition: FeatureFloodCount.h:566
FeatureFloodCount::_halo_ids
std::vector< std::map< dof_id_type, int > > _halo_ids
The data structure for looking up halos around features.
Definition: FeatureFloodCount.h:701
GrainTracker::RemapCacheMode
RemapCacheMode
Definition: GrainTracker.h:53
GrainTracker::_update_field_info
const PerfID _update_field_info
Definition: GrainTracker.h:255
FeatureFloodCount::execute
virtual void execute() override
Definition: FeatureFloodCount.C:358
FeatureFloodCount::invalid_id
static const unsigned int invalid_id
Definition: FeatureFloodCount.h:94
GrainTracker::getVarToFeatureVector
virtual const std::vector< unsigned int > & getVarToFeatureVector(dof_id_type elem_id) const override
Returns a list of active unique feature ids for a particular element.
Definition: GrainTracker.C:131
GrainTracker::_is_transient
const bool _is_transient
Boolean to indicate whether this is a Steady or Transient solve.
Definition: GrainTracker.h:245
GrainDistance::operator<
bool operator<(const GrainDistance &rhs) const
Definition: GrainTracker.C:1829
GrainTracker::communicateHaloMap
void communicateHaloMap()
Definition: GrainTracker.C:1632
GrainTracker::GrainTracker
GrainTracker(const InputParameters &parameters)
Definition: GrainTracker.C:80
GrainTracker::_tracking_step
const int _tracking_step
The timestep to begin tracking grains.
Definition: GrainTracker.h:180
FeatureFloodCount::Status::INACTIVE
GrainTracker::_reserve_op_index
const std::size_t _reserve_op_index
The cutoff index where if variable index >= this number, no remapping TO that variable will occur.
Definition: GrainTracker.h:193
FeatureFloodCount::_feature_maps
std::vector< std::map< dof_id_type, int > > _feature_maps
The feature maps contain the raw flooded node information and eventually the unique grain numbers.
Definition: FeatureFloodCount.h:678
registerMooseObject
registerMooseObject("PhaseFieldApp", GrainTracker)
GrainTracker::_max_remap_recursion_depth
const unsigned short _max_remap_recursion_depth
Depth of renumbering recursion (a depth of zero means no recursion)
Definition: GrainTracker.h:186
GrainTracker::getGrainCentroid
virtual Point getGrainCentroid(unsigned int grain_id) const override
Returns the centroid for the given grain number.
Definition: GrainTracker.C:157
FeatureFloodCount::FeatureData::_var_index
std::size_t _var_index
The Moose variable where this feature was found (often the "order parameter")
Definition: FeatureFloodCount.h:284
GrainTracker::initialize
virtual void initialize() override
Definition: GrainTracker.C:230
FeatureFloodCount::_feature_id_to_local_index
std::vector< std::size_t > _feature_id_to_local_index
The vector recording the grain_id to local index (several indices will contain invalid_size_t)
Definition: FeatureFloodCount.h:684
GrainTracker::getThreshold
virtual Real getThreshold(std::size_t current_index) const override
Return the starting comparison threshold to use when inspecting an entity during the flood stage.
Definition: GrainTracker.C:286
FeatureFloodCount::_compute_halo_maps
const bool _compute_halo_maps
Indicates whether or not to communicate halo map information with all ranks.
Definition: FeatureFloodCount.h:604
GrainTracker::_reserve_op_threshold
const Real _reserve_op_threshold
The threshold above (or below) where a grain may be found on a reserve op field.
Definition: GrainTracker.h:196
dataLoad
void dataLoad(std::istream &stream, GrainTracker::PartialFeatureData &feature, void *context)
Definition: GrainTracker.C:38
GrainTracker::getFeatureVar
virtual unsigned int getFeatureVar(unsigned int feature_id) const override
Returns the variable representing the passed in feature.
Definition: GrainTracker.C:137
PolycrystalUserObjectBase::getVariableValue
virtual Real getVariableValue(unsigned int op_index, const Point &p) const =0
Returns the variable value for a given op_index and mesh point.
GrainTracker::_track_grains
const PerfID _track_grains
Definition: GrainTracker.h:253
FeatureFloodCount::FeatureData::halosIntersect
bool halosIntersect(const FeatureData &rhs) const
Determine if one of this FeaturesData's member sets intersects the other FeatureData's corresponding ...
Definition: FeatureFloodCount.C:1930
GrainTracker.h
GrainTracker::_verbosity_level
const short _verbosity_level
Verbosity level controlling the amount of information printed to the console.
Definition: GrainTracker.h:219
FeatureFloodCount::_feature_count
unsigned int _feature_count
The number of features seen by this object (same as summing _feature_counts_per_map)
Definition: FeatureFloodCount.h:648
FeatureFloodCount::FeatureData::_local_ids
container_type _local_ids
Holds the local ids in the interior of a feature.
Definition: FeatureFloodCount.h:272
GrainTracker::remapGrains
void remapGrains()
This method is called after trackGrains to remap grains that are too close to each other.
Definition: GrainTracker.C:913
GrainTracker::_tolerate_failure
const bool _tolerate_failure
Indicates whether we should continue after a remap failure (will result in non-physical results)
Definition: GrainTracker.h:202
GrainTracker::getNextUniqueID
unsigned int getNextUniqueID()
Retrieve the next unique grain number if a new grain is detected during trackGrains.
Definition: GrainTracker.C:1795
GrainTracker::getNumberActiveGrains
virtual std::size_t getNumberActiveGrains() const override
Returns the number of active grains current stored in the GrainTracker.
Definition: GrainTracker.C:143
FeatureFloodCount::Status::DIRTY
FeatureFloodCount::buildFeatureIdToLocalIndices
void buildFeatureIdToLocalIndices(unsigned int max_id)
This method builds a lookup map for retrieving the right local feature (by index) given a global inde...
Definition: FeatureFloodCount.C:665
GrainTracker::_remap_timer
const PerfID _remap_timer
Definition: GrainTracker.h:252
FeatureFloodCount::BoundaryIntersection::NONE
FeatureFloodCount::_compute_var_to_feature_map
const bool _compute_var_to_feature_map
Indicates whether or not the var to feature map is populated.
Definition: FeatureFloodCount.h:607
FeatureFloodCount::FeatureData::_bboxes
std::vector< BoundingBox > _bboxes
The vector of bounding boxes completely enclosing this feature (multiple used with periodic constrain...
Definition: FeatureFloodCount.h:291
FeatureFloodCount::_step_threshold
Real _step_threshold
Definition: FeatureFloodCount.h:573
FeatureFloodCount::_n_procs
const processor_id_type _n_procs
Convenience variable holding the number of processors in this simulation.
Definition: FeatureFloodCount.h:623
FeatureFloodCount::_var_index_maps
std::vector< std::map< dof_id_type, int > > _var_index_maps
This map keeps track of which variables own which nodes.
Definition: FeatureFloodCount.h:639
GrainTracker::isFeaturePercolated
virtual bool isFeaturePercolated(unsigned int feature_id) const override
Returns a Boolean indicating whether this feature is percolated (e.g.
Definition: GrainTracker.C:208
FeatureFloodCount::_n_vars
const std::size_t _n_vars
Definition: FeatureFloodCount.h:617
GrainTracker::doesFeatureIntersectSpecifiedBoundary
virtual bool doesFeatureIntersectSpecifiedBoundary(unsigned int feature_id) const override
Returns a Boolean indicating whether this feature intersects boundaries in a user-supplied list.
Definition: GrainTracker.C:191
GrainTracker::meshChanged
virtual void meshChanged() override
Definition: GrainTracker.C:248
FeatureFloodCount::meshChanged
virtual void meshChanged() override
Definition: FeatureFloodCount.C:335
FeatureFloodCount::_feature_counts_per_map
std::vector< unsigned int > _feature_counts_per_map
The number of features seen by this object per map.
Definition: FeatureFloodCount.h:645
FeatureFloodCount::BoundaryIntersection::PRIMARY_PERCOLATION_BOUNDARY
FeatureFloodCount::FeatureData::_status
Status _status
The status of a feature (used mostly in derived classes like the GrainTracker)
Definition: FeatureFloodCount.h:307
FeatureFloodCount::Status::MARKED
GrainTracker::_first_time
bool & _first_time
Boolean to indicate the first time this object executes.
Definition: GrainTracker.h:225