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 "DiscreteNucleationMap.h" 11 : #include "MooseMesh.h" 12 : #include "Executioner.h" 13 : 14 : #include "libmesh/quadrature.h" 15 : 16 : registerMooseObject("PhaseFieldApp", DiscreteNucleationMap); 17 : 18 : InputParameters 19 230 : DiscreteNucleationMap::validParams() 20 : { 21 230 : InputParameters params = ElementUserObject::validParams(); 22 230 : params.addClassDescription("Generates a spatial smoothed map of all nucleation sites with the " 23 : "data of the DiscreteNucleationInserter for use by the " 24 : "DiscreteNucleation material."); 25 460 : params.addParam<Real>("int_width", 0.0, "Nucleus interface width for smooth nuclei"); 26 460 : params.addRequiredParam<UserObjectName>("inserter", "DiscreteNucleationInserter user object"); 27 460 : params.addCoupledVar("periodic", 28 : "Use the periodicity settings of this variable to populate the grain map"); 29 : // the mapping needs to run at timestep begin, which is after the adaptivity 30 : // run of the previous timestep. 31 230 : params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_BEGIN; 32 230 : return params; 33 0 : } 34 : 35 124 : DiscreteNucleationMap::DiscreteNucleationMap(const InputParameters & parameters) 36 : : ElementUserObject(parameters), 37 124 : _mesh_changed(false), 38 124 : _inserter(getUserObject<DiscreteNucleationInserterBase>("inserter")), 39 226 : _periodic_var(isCoupled("periodic") ? getFieldVar("periodic", 0) : nullptr), 40 248 : _int_width(getParam<Real>("int_width")), 41 248 : _nucleus_list(_inserter.getNucleusList()) 42 : { 43 124 : _zero_map.assign(_fe_problem.getMaxQps(), 0.0); 44 124 : } 45 : 46 : void 47 923 : DiscreteNucleationMap::initialize() 48 : { 49 923 : if (_inserter.isMapUpdateRequired() || _mesh_changed || _force_rebuild_map) 50 : { 51 : // If the last solve didn't converge, postpone the rebuild until the next timestep 52 419 : if (_app.getExecutioner()->lastSolveConverged()) 53 : { 54 410 : _rebuild_map = true; 55 : _nucleus_map.clear(); 56 410 : _force_rebuild_map = false; 57 : } 58 : else 59 : { 60 9 : _rebuild_map = false; 61 9 : _force_rebuild_map = true; 62 : } 63 : } 64 : else 65 504 : _rebuild_map = false; 66 : 67 923 : _mesh_changed = false; 68 923 : } 69 : 70 : void 71 1375858 : DiscreteNucleationMap::execute() 72 : { 73 1375858 : if (_rebuild_map) 74 : { 75 : // reserve space for each quadrature point in the element 76 717459 : _elem_map.assign(_qrule->n_points(), 0); 77 : 78 : // store a random number for each quadrature point 79 : unsigned int active_nuclei = 0; 80 3550095 : for (unsigned int qp = 0; qp < _qrule->n_points(); ++qp) 81 : { 82 : Real r = std::numeric_limits<Real>::max(); 83 : 84 : // find the distance to the closest nucleus 85 : Real local_radius = 0.0; 86 11108296 : for (unsigned i = 0; i < _nucleus_list.size(); ++i) 87 : { 88 : // use a non-periodic or periodic distance 89 8275660 : r = _periodic_var 90 8275660 : ? _mesh.minPeriodicDistance(*_periodic_var, _q_point[qp], _nucleus_list[i].center) 91 1327600 : : (_q_point[qp] - _nucleus_list[i].center).norm(); 92 : 93 : // grab the radius of the nucleus that this qp is closest to 94 8275660 : local_radius = _nucleus_list[i].radius; 95 : 96 : // compute intensity value with smooth interface 97 : Real value = 0.0; 98 8275660 : if (r <= local_radius - _int_width / 2.0) // Inside circle 99 : { 100 931124 : active_nuclei++; 101 : value = 1.0; 102 : } 103 7344536 : else if (r < local_radius + _int_width / 2.0) // Smooth interface 104 : { 105 41869 : Real int_pos = (r - local_radius + _int_width / 2.0) / _int_width; 106 41869 : active_nuclei++; 107 41869 : value = (1.0 + std::cos(int_pos * libMesh::pi)) / 2.0; 108 : } 109 8275660 : if (value > _elem_map[qp]) 110 773289 : _elem_map[qp] = value; 111 : } 112 : } 113 : 114 : // if the map is not empty insert it 115 717459 : if (active_nuclei > 0) 116 : _nucleus_map.insert( 117 401122 : std::pair<dof_id_type, std::vector<Real>>(_current_elem->id(), _elem_map)); 118 : } 119 1375858 : } 120 : 121 : void 122 142 : DiscreteNucleationMap::threadJoin(const UserObject & y) 123 : { 124 : // if the map needs to be updated we merge the maps from all threads 125 142 : if (_rebuild_map) 126 : { 127 : const auto & uo = static_cast<const DiscreteNucleationMap &>(y); 128 64 : _nucleus_map.insert(uo._nucleus_map.begin(), uo._nucleus_map.end()); 129 : } 130 142 : } 131 : 132 : void 133 21 : DiscreteNucleationMap::meshChanged() 134 : { 135 21 : _mesh_changed = true; 136 21 : } 137 : 138 : const std::vector<Real> & 139 1910974 : DiscreteNucleationMap::nuclei(const Elem * elem) const 140 : { 141 1910974 : NucleusMap::const_iterator i = _nucleus_map.find(elem->id()); 142 : 143 : // if no entry in the map was found the element contains no nucleus 144 1910974 : if (i == _nucleus_map.end()) 145 1524386 : return _zero_map; 146 : 147 386588 : return i->second; 148 : }