https://mooseframework.inl.gov
PolycrystalVoronoiVoidIC.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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 
11 
12 // MOOSE includes
13 #include "MooseMesh.h"
14 #include "MooseVariable.h"
15 #include "DelimitedFileReader.h"
16 #include "GrainTrackerInterface.h"
17 #include "PolycrystalVoronoi.h"
18 
21 {
23 
24  params.addRequiredParam<unsigned int>("op_num", "Number of order parameters");
25 
26  params.addParam<bool>(
27  "columnar_3D", false, "3D microstructure will be columnar in the z-direction?");
28 
29  return params;
30 }
31 
33 
36 {
38  MooseEnum structure_options("grains voids");
39  params.addRequiredParam<MooseEnum>("structure_type",
40  structure_options,
41  "Which structure type is being initialized, grains or voids");
42  params.addParam<unsigned int>("op_index",
43  0,
44  "The index for the current order parameter, "
45  "not needed if structure_type = voids");
46  params.addRequiredParam<UserObjectName>(
47  "polycrystal_ic_uo", "UserObject for obtaining the polycrystal grain structure.");
48  params.addParam<FileName>("file_name",
49  "",
50  "File containing grain centroids, if file_name is provided, "
51  "the centroids from the file will be used.");
52  return params;
53 }
54 
56  : MultiSmoothCircleIC(parameters),
57  _structure_type(getParam<MooseEnum>("structure_type")),
58  _op_num(getParam<unsigned int>("op_num")),
59  _op_index(getParam<unsigned int>("op_index")),
60  _columnar_3D(getParam<bool>("columnar_3D")),
61  _poly_ic_uo(getUserObject<PolycrystalVoronoi>("polycrystal_ic_uo")),
62  _file_name(getParam<FileName>("file_name"))
63 {
64  if (_invalue < _outvalue)
65  mooseWarning("Detected invalue < outvalue in PolycrystalVoronoiVoidIC. Please make sure that's "
66  "the intended usage for representing voids.");
67  if (_numbub == 0)
68  mooseError("PolycrystalVoronoiVoidIC requires numbub > 0. If you want no voids to "
69  "be represented, use invalue = outvalue. In general, you should use "
70  "PolycrystalVoronoi to represent Voronoi grain structures without "
71  "voids.");
72 }
73 
74 void
76 {
77  if (_op_num <= _op_index)
78  mooseError("op_index is too large in CircleGrainVoidIC");
79 
80  // Obtain total number and centerpoints of the grains
83 
84  // Call initial setup from MultiSmoothCircleIC to create _centers and _radii
85  // for voids
87 }
88 
89 void
91 {
92  _centers.resize(_numbub);
93 
94  // This Code will place void center points on grain boundaries
95  for (unsigned int vp = 0; vp < _numbub; ++vp)
96  {
97  bool try_again;
98  unsigned int num_tries = 0;
99 
100  do
101  {
102  try_again = false;
103  num_tries++;
104 
105  if (num_tries > _max_num_tries)
106  mooseError("Too many tries of assigning void centers in "
107  "PolycrystalVoronoiVoidIC");
108 
109  Point rand_point;
110 
111  for (const auto i : make_range(Moose::dim))
112  rand_point(i) = _bottom_left(i) + _range(i) * _random.rand(_tid);
113 
114  // Allow the vectors to be sorted based on their distance from the
115  // rand_point
116  std::vector<PolycrystalVoronoiVoidIC::DistancePoint> diff(_grain_num);
117 
118  for (unsigned int gr = 0; gr < _grain_num; ++gr)
119  {
120  diff[gr].d = _mesh.minPeriodicDistance(_var.number(), rand_point, _centerpoints[gr]);
121  diff[gr].gr = gr;
122  }
123 
124  std::sort(diff.begin(), diff.end(), _customLess);
125 
126  Point closest_point = _centerpoints[diff[0].gr];
127  Point next_closest_point = _centerpoints[diff[1].gr];
128 
129  // Find Slope of Line in the plane orthogonal to the diff_centerpoint
130  // vector
131  Point pa = rand_point + _mesh.minPeriodicVector(_var.number(), rand_point, closest_point);
132  Point pb =
133  rand_point + _mesh.minPeriodicVector(_var.number(), rand_point, next_closest_point);
134  Point diff_centerpoints = pb - pa;
135 
136  Point diff_rand_center = _mesh.minPeriodicVector(_var.number(), closest_point, rand_point);
137  Point normal_vector = diff_centerpoints.cross(diff_rand_center);
138  Point slope = normal_vector.cross(diff_centerpoints);
139 
140  // Midpoint position vector between two center points
141  Point midpoint = closest_point + (0.5 * diff_centerpoints);
142 
143  // Solve for the scalar multiplier solution on the line
144  Real lambda = 0;
145  Point mid_rand_vector = _mesh.minPeriodicVector(_var.number(), midpoint, rand_point);
146 
147  Real slope_dot = slope * slope;
148  mooseAssert(slope_dot > 0, "The dot product of slope with itself is zero");
149  for (const auto i : make_range(Moose::dim))
150  lambda += (mid_rand_vector(i) * slope(i)) / slope_dot;
151 
152  // Assigning points to vector
153  _centers[vp] = slope * lambda + midpoint;
154 
155  // Checking to see if points are in the domain ONLY WORKS FOR PERIODIC
156  for (const auto i : make_range(Moose::dim))
157  if ((_centers[vp](i) > _top_right(i)) || (_centers[vp](i) < _bottom_left(i)))
158  try_again = true;
159 
160  for (unsigned int i = 0; i < vp; ++i)
161  {
163 
164  if (dist < _bubspac)
165  try_again = true;
166  }
167 
168  // Two algorithms are available for screening bubbles falling in grain
169  // interior. They produce
170  // nearly identical results.
171  // Here only one is listed. The other one is available upon request.
172 
173  // Use circle center for checking whether voids are at GBs
174  if (try_again == false)
175  {
176  Real min_rij_1, min_rij_2, rij, rij_diff_tol;
177 
178  min_rij_1 = _range.norm();
179  min_rij_2 = _range.norm();
180 
181  rij_diff_tol = 0.1 * _radius;
182 
183  for (unsigned int gr = 0; gr < _grain_num; ++gr)
184  {
186 
187  if (rij < min_rij_1)
188  {
189  min_rij_2 = min_rij_1;
190  min_rij_1 = rij;
191  }
192  else if (rij < min_rij_2)
193  min_rij_2 = rij;
194  }
195 
196  if (std::abs(min_rij_1 - min_rij_2) > rij_diff_tol)
197  try_again = true;
198  }
199 
200  } while (try_again == true);
201  }
202 }
203 
204 Real
206 {
207  Real value = 0.0;
208 
209  // Determine value for voids
210  Real void_value = MultiSmoothCircleIC::value(p);
211 
212  // Determine value for grains
213  Real grain_value = _poly_ic_uo.getVariableValue(_op_index, p);
214 
215  switch (_structure_type)
216  {
217  case 0: // assigning values for grains (order parameters)
218  if (grain_value == 0) // Not in this grain
219  value = grain_value;
220  else // in this grain, but might be in a void
221  if (void_value == _outvalue) // Not in a void
222  value = grain_value;
223  else if (void_value > _outvalue && void_value < _invalue) // On void interface
224  value = grain_value * (_invalue - void_value) / (_invalue - _outvalue);
225  else if (void_value == _invalue) // In a void, so op = 0
226  value = 0.0;
227  break;
228 
229  case 1: // assigning values for voids (concentration)
230  value = void_value;
231  break;
232  }
233 
234  return value;
235 }
236 
239 {
241  RealGradient void_gradient = MultiSmoothCircleIC::gradient(p);
242 
243  // Order parameter assignment assumes zero gradient (sharp interface)
244  switch (_structure_type)
245  {
246  case 1: // assigning gradient for voids
247  gradient = void_gradient;
248  break;
249  }
250 
251  return gradient;
252 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
unsigned int number() const
struct PolycrystalVoronoiVoidIC::DistancePointComparator _customLess
std::vector< Point > _centers
const PolycrystalVoronoi & _poly_ic_uo
RealVectorValue minPeriodicVector(unsigned int nonlinear_var_num, Point p, Point q) const
const unsigned int _max_num_tries
virtual std::vector< Point > getGrainCenters() const
virtual Real value(const Point &p)
virtual Real getVariableValue(unsigned int op_index, const Point &p) const override
Returns the variable value for a given op_index and mesh point.
static constexpr std::size_t dim
std::vector< Point > _centerpoints
PolycrystalVoronoiVoidIC initializes either grain or void values for a voronoi tesselation with voids...
MooseVariableField< T > & _var
void mooseWarning(Args &&... args) const
void addRequiredParam(const std::string &name, const std::string &doc_string)
virtual Real value(const Point &p) override
static InputParameters validParams()
const unsigned int _numbub
TypeVector< typename CompareTypes< Real, T2 >::supertype > cross(const TypeVector< T2 > &v) const
PolycrystalVoronoiVoidIC(const InputParameters &parameters)
Real minPeriodicDistance(unsigned int nonlinear_var_num, Point p, Point q) const
virtual RealGradient gradient(const Point &p)
registerMooseObject("PhaseFieldApp", PolycrystalVoronoiVoidIC)
static InputParameters validParams()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
IntRange< T > make_range(T beg, T end)
void mooseError(Args &&... args) const
virtual void initialSetup() override
virtual unsigned int getNumGrains() const override
Must be overridden by the deriving class to provide the number of grains in the polycrystal structure...
virtual RealGradient gradient(const Point &p) override
Real rand(std::size_t i)
static InputParameters actionParameters()
void ErrorVector unsigned int
virtual void computeCircleCenters() override
MultismoothCircleIC creates multiple SmoothCircles (number = numbub) that are randomly positioned aro...
virtual void initialSetup() override