Line data Source code
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 :
10 : #include "PolycrystalVoronoiVoidIC.h"
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 :
19 : InputParameters
20 651 : PolycrystalVoronoiVoidIC::actionParameters()
21 : {
22 651 : InputParameters params = MultiSmoothCircleIC::validParams();
23 :
24 1302 : params.addRequiredParam<unsigned int>("op_num", "Number of order parameters");
25 :
26 1302 : params.addParam<bool>(
27 1302 : "columnar_3D", false, "3D microstructure will be columnar in the z-direction?");
28 :
29 651 : return params;
30 0 : }
31 :
32 : registerMooseObject("PhaseFieldApp", PolycrystalVoronoiVoidIC);
33 :
34 : InputParameters
35 612 : PolycrystalVoronoiVoidIC::validParams()
36 : {
37 612 : InputParameters params = PolycrystalVoronoiVoidIC::actionParameters();
38 1224 : MooseEnum structure_options("grains voids");
39 1224 : params.addRequiredParam<MooseEnum>("structure_type",
40 : structure_options,
41 : "Which structure type is being initialized, grains or voids");
42 1224 : params.addParam<unsigned int>("op_index",
43 1224 : 0,
44 : "The index for the current order parameter, "
45 : "not needed if structure_type = voids");
46 1224 : params.addRequiredParam<UserObjectName>(
47 : "polycrystal_ic_uo", "UserObject for obtaining the polycrystal grain structure.");
48 1224 : 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 612 : return params;
53 612 : }
54 :
55 332 : PolycrystalVoronoiVoidIC::PolycrystalVoronoiVoidIC(const InputParameters & parameters)
56 : : MultiSmoothCircleIC(parameters),
57 332 : _structure_type(getParam<MooseEnum>("structure_type")),
58 664 : _op_num(getParam<unsigned int>("op_num")),
59 664 : _op_index(getParam<unsigned int>("op_index")),
60 664 : _columnar_3D(getParam<bool>("columnar_3D")),
61 332 : _poly_ic_uo(getUserObject<PolycrystalVoronoi>("polycrystal_ic_uo")),
62 996 : _file_name(getParam<FileName>("file_name"))
63 : {
64 332 : if (_invalue < _outvalue)
65 0 : mooseWarning("Detected invalue < outvalue in PolycrystalVoronoiVoidIC. Please make sure that's "
66 : "the intended usage for representing voids.");
67 332 : if (_numbub == 0)
68 0 : 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 332 : }
73 :
74 : void
75 322 : PolycrystalVoronoiVoidIC::initialSetup()
76 : {
77 322 : if (_op_num <= _op_index)
78 0 : mooseError("op_index is too large in CircleGrainVoidIC");
79 :
80 : // Obtain total number and centerpoints of the grains
81 322 : _grain_num = _poly_ic_uo.getNumGrains();
82 322 : _centerpoints = _poly_ic_uo.getGrainCenters();
83 :
84 : // Call initial setup from MultiSmoothCircleIC to create _centers and _radii
85 : // for voids
86 322 : MultiSmoothCircleIC::initialSetup();
87 322 : }
88 :
89 : void
90 322 : PolycrystalVoronoiVoidIC::computeCircleCenters()
91 : {
92 322 : _centers.resize(_numbub);
93 :
94 : // This Code will place void center points on grain boundaries
95 3776 : 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 8924 : num_tries++;
104 :
105 8924 : if (num_tries > _max_num_tries)
106 0 : mooseError("Too many tries of assigning void centers in "
107 : "PolycrystalVoronoiVoidIC");
108 :
109 : Point rand_point;
110 :
111 35696 : for (const auto i : make_range(Moose::dim))
112 26772 : 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 8924 : std::vector<PolycrystalVoronoiVoidIC::DistancePoint> diff(_grain_num);
117 :
118 113524 : for (unsigned int gr = 0; gr < _grain_num; ++gr)
119 : {
120 104600 : diff[gr].d = _mesh.minPeriodicDistance(_var, rand_point, _centerpoints[gr]);
121 104600 : diff[gr].gr = gr;
122 : }
123 :
124 8924 : std::sort(diff.begin(), diff.end(), _customLess);
125 :
126 8924 : Point closest_point = _centerpoints[diff[0].gr];
127 8924 : 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 8924 : Point pa = rand_point + _mesh.minPeriodicVector(_var, rand_point, closest_point);
132 8924 : Point pb = rand_point + _mesh.minPeriodicVector(_var, rand_point, next_closest_point);
133 : Point diff_centerpoints = pb - pa;
134 :
135 8924 : Point diff_rand_center = _mesh.minPeriodicVector(_var, closest_point, rand_point);
136 8924 : Point normal_vector = diff_centerpoints.cross(diff_rand_center);
137 8924 : Point slope = normal_vector.cross(diff_centerpoints);
138 :
139 : // Midpoint position vector between two center points
140 : Point midpoint = closest_point + (0.5 * diff_centerpoints);
141 :
142 : // Solve for the scalar multiplier solution on the line
143 : Real lambda = 0;
144 8924 : Point mid_rand_vector = _mesh.minPeriodicVector(_var, midpoint, rand_point);
145 :
146 : Real slope_dot = slope * slope;
147 : mooseAssert(slope_dot > 0, "The dot product of slope with itself is zero");
148 35696 : for (const auto i : make_range(Moose::dim))
149 26772 : lambda += (mid_rand_vector(i) * slope(i)) / slope_dot;
150 :
151 : // Assigning points to vector
152 8924 : _centers[vp] = slope * lambda + midpoint;
153 :
154 : // Checking to see if points are in the domain ONLY WORKS FOR PERIODIC
155 35696 : for (const auto i : make_range(Moose::dim))
156 26772 : if ((_centers[vp](i) > _top_right(i)) || (_centers[vp](i) < _bottom_left(i)))
157 : try_again = true;
158 :
159 98822 : for (unsigned int i = 0; i < vp; ++i)
160 : {
161 89898 : Real dist = _mesh.minPeriodicDistance(_var, _centers[vp], _centers[i]);
162 :
163 89898 : if (dist < _bubspac)
164 : try_again = true;
165 : }
166 :
167 : // Two algorithms are available for screening bubbles falling in grain
168 : // interior. They produce
169 : // nearly identical results.
170 : // Here only one is listed. The other one is available upon request.
171 :
172 : // Use circle center for checking whether voids are at GBs
173 8924 : if (try_again == false)
174 : {
175 : Real min_rij_1, min_rij_2, rij, rij_diff_tol;
176 :
177 3532 : min_rij_1 = _range.norm();
178 3532 : min_rij_2 = _range.norm();
179 :
180 3532 : rij_diff_tol = 0.1 * _radius;
181 :
182 63758 : for (unsigned int gr = 0; gr < _grain_num; ++gr)
183 : {
184 60226 : rij = _mesh.minPeriodicDistance(_var, _centers[vp], _centerpoints[gr]);
185 :
186 60226 : if (rij < min_rij_1)
187 : {
188 : min_rij_2 = min_rij_1;
189 : min_rij_1 = rij;
190 : }
191 48458 : else if (rij < min_rij_2)
192 : min_rij_2 = rij;
193 : }
194 :
195 3532 : if (std::abs(min_rij_1 - min_rij_2) > rij_diff_tol)
196 : try_again = true;
197 : }
198 :
199 8924 : } while (try_again == true);
200 : }
201 322 : }
202 :
203 : Real
204 3264560 : PolycrystalVoronoiVoidIC::value(const Point & p)
205 : {
206 : Real value = 0.0;
207 :
208 : // Determine value for voids
209 3264560 : Real void_value = MultiSmoothCircleIC::value(p);
210 :
211 : // Determine value for grains
212 3264560 : Real grain_value = _poly_ic_uo.getVariableValue(_op_index, p);
213 :
214 3264560 : switch (_structure_type)
215 : {
216 2792160 : case 0: // assigning values for grains (order parameters)
217 2792160 : if (grain_value == 0) // Not in this grain
218 : value = grain_value;
219 : else // in this grain, but might be in a void
220 496040 : if (void_value == _outvalue) // Not in a void
221 : value = grain_value;
222 58792 : else if (void_value > _outvalue && void_value < _invalue) // On void interface
223 52872 : value = grain_value * (_invalue - void_value) / (_invalue - _outvalue);
224 : else if (void_value == _invalue) // In a void, so op = 0
225 : value = 0.0;
226 : break;
227 :
228 : case 1: // assigning values for voids (concentration)
229 : value = void_value;
230 : break;
231 : }
232 :
233 3264560 : return value;
234 : }
235 :
236 : RealGradient
237 0 : PolycrystalVoronoiVoidIC::gradient(const Point & p)
238 : {
239 : RealGradient gradient;
240 0 : RealGradient void_gradient = MultiSmoothCircleIC::gradient(p);
241 :
242 : // Order parameter assignment assumes zero gradient (sharp interface)
243 0 : switch (_structure_type)
244 : {
245 0 : case 1: // assigning gradient for voids
246 0 : gradient = void_gradient;
247 0 : break;
248 : }
249 :
250 0 : return gradient;
251 : }
|