1 #include "../../include/matrixFreePDE.h" 3 #include "../../include/FloodFiller.h" 4 #include "../../include/SimplifiedGrainRepresentation.h" 5 #include "../../include/OrderParameterRemapper.h" 9 template <
int dim,
int degree>
13 computing_timer.enter_section(
"matrixFreePDE: reassignGrains");
15 pcout <<
"Reassigning grains..." << std::endl;
18 unsigned int scalar_field_index = 0;
19 for (
unsigned int var=0; var<userInputs.number_of_variables; var++){
20 if (userInputs.var_type.at(var) ==
SCALAR){
21 scalar_field_index = var;
27 QGaussLobatto<dim> quadrature2 (degree+1);
30 std::vector<GrainSet<dim>> grain_sets;
32 unsigned int op_list_index = 0;
33 for(
unsigned int fieldIndex=0; fieldIndex<fields.size(); fieldIndex++){
34 if (op_list_index < userInputs.variables_for_remapping.size()){
35 if (fieldIndex == userInputs.variables_for_remapping.at(op_list_index)){
38 std::vector<GrainSet<dim>> single_OP_grain_sets;
39 flood_filler.
calcGrainSets(*FESet.at(scalar_field_index), *dofHandlersSet_nonconst.at(scalar_field_index), solutionSet.at(fieldIndex), userInputs.order_parameter_threshold, 1.0+userInputs.order_parameter_threshold, fieldIndex, single_OP_grain_sets);
41 grain_sets.insert(grain_sets.end(), single_OP_grain_sets.begin(), single_OP_grain_sets.end());
47 for (
unsigned int g=0; g<grain_sets.size(); g++){
48 grain_sets.at(g).setGrainIndex(g);
51 std::vector<SimplifiedGrainRepresentation<dim>> old_grain_representations = simplified_grain_representations;
52 simplified_grain_representations.clear();
53 for (
unsigned int g=0; g<grain_sets.size(); g++){
56 pcout <<
"Grain: " << simplified_grain_representation.
getGrainId() <<
" " << simplified_grain_representation.getOrderParameterId() <<
" Center: " << simplified_grain_representation.getCenter()(0) <<
" " << simplified_grain_representation.getCenter()(1) << std::endl;
58 simplified_grain_representations.push_back(simplified_grain_representation);
63 if (currentIncrement > 0 || userInputs.load_grain_structure){
64 simplified_grain_manipulator.
transferGrainIds(old_grain_representations, simplified_grain_representations);
67 simplified_grain_manipulator.
reassignGrains(simplified_grain_representations, userInputs.buffer_between_grains, userInputs.variables_for_remapping);
69 for (
unsigned int g=0; g<this->simplified_grain_representations.size(); g++){
70 pcout <<
"Grain: " << simplified_grain_representations[g].getGrainId() <<
" " << simplified_grain_representations[g].getOrderParameterId() <<
" Center: " << simplified_grain_representations[g].getCenter()(0) <<
" " << simplified_grain_representations[g].getCenter()(1) << std::endl;
75 order_parameter_remapper.
remap(simplified_grain_representations, solutionSet, *dofHandlersSet_nonconst.at(scalar_field_index), FESet.at(scalar_field_index)->dofs_per_cell, userInputs.buffer_between_grains);
77 pcout <<
"Reassigning grains completed." << std::endl << std::endl;
80 computing_timer.exit_section(
"matrixFreePDE: reassignGrains");
85 #include "../../include/matrixFreePDE_template_instantiations.h" void remap(std::vector< SimplifiedGrainRepresentation< dim >> &grain_representations, std::vector< vectorType *> &solution_fields, dealii::DoFHandler< dim > &dof_handler, unsigned int dofs_per_cell, double buffer)
void calcGrainSets(dealii::FESystem< dim > &fe, dealii::DoFHandler< dim > &dof_handler, vectorType *solution_field, double threshold_lower, double threshold_upper, unsigned int order_parameter_index, std::vector< GrainSet< dim >> &grain_sets)
void reassignGrains(std::vector< SimplifiedGrainRepresentation< dim >> &grain_representations, double buffer_distance, std::vector< unsigned int > order_parameter_id_list)
void transferGrainIds(const std::vector< SimplifiedGrainRepresentation< dim >> &old_grain_representations, std::vector< SimplifiedGrainRepresentation< dim >> &new_grain_representations) const
unsigned int getGrainId() const