www.mooseframework.org
EBSDReader.h
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 #pragma once
11 
12 #include "EulerAngleProvider.h"
13 #include "EBSDAccessFunctors.h"
14 
15 class EBSDReader;
16 
17 template <>
18 InputParameters validParams<EBSDReader>();
19 
37 {
38 public:
39  EBSDReader(const InputParameters & params);
40  virtual ~EBSDReader();
41 
42  virtual void readFile();
43 
44  virtual void initialize() {}
45  virtual void execute() {}
46  virtual void finalize() {}
47 
51  const EBSDPointData & getData(const Point & p) const;
52 
56  const EBSDAvgData & getAvgData(unsigned int i) const;
57 
61  const EBSDAvgData & getAvgData(unsigned int phase, unsigned int local_id) const;
62 
66  virtual const EulerAngles & getEulerAngles(unsigned int) const;
67 
71  virtual unsigned int getGrainNum() const;
72 
76  virtual unsigned int getPhaseNum() const { return _global_id.size(); }
77 
81  virtual unsigned int getGrainNum(unsigned int phase) const;
82 
84  unsigned int getFeatureID(unsigned int phase, unsigned int local_id) const
85  {
86  return _avg_data[_global_id[phase][local_id]]._feature_id;
87  }
89  unsigned int getFeatureID(unsigned int global_id) const
90  {
91  return _avg_data[global_id]._feature_id;
92  }
93 
95  virtual unsigned int getGlobalID(unsigned int phase, unsigned int local_id) const
96  {
97  return _global_id[phase][local_id];
98  }
100  virtual unsigned int getGlobalID(unsigned int feature_id) const;
101 
103  MooseSharedPointer<EBSDPointDataFunctor>
104  getPointDataAccessFunctor(const MooseEnum & field_name) const;
106  MooseSharedPointer<EBSDAvgDataFunctor>
107  getAvgDataAccessFunctor(const MooseEnum & field_name) const;
108 
113  const std::map<dof_id_type, std::vector<Real>> & getNodeToGrainWeightMap() const;
114 
119  const std::map<dof_id_type, std::vector<Real>> & getNodeToPhaseWeightMap() const;
120 
122  void meshChanged();
123 
124 protected:
126  MooseMesh & _mesh;
127  NonlinearSystemBase & _nl;
129 
131  unsigned int _grain_num;
133  Point _top_right;
134  Point _range;
136 
138  unsigned int _custom_columns;
139 
141  std::vector<EBSDPointData> _data;
142 
144  std::vector<EBSDAvgData> _avg_data;
145 
147  std::vector<EulerAngles> _avg_angles;
148 
150  std::map<unsigned int, unsigned int> _global_id_map;
151 
153  std::vector<std::vector<unsigned int>> _global_id;
154 
156  std::map<dof_id_type, std::vector<Real>> _node_to_grain_weight_map;
157 
159  std::map<dof_id_type, std::vector<Real>> _node_to_phase_weight_map;
160 
162  const int & _time_step;
163 
165  unsigned int _mesh_dimension;
166 
168  unsigned int _bins;
169 
171  unsigned int _L_norm;
172 
174  unsigned _nx, _ny, _nz;
175 
177  Real _dx, _dy, _dz;
178 
180  Real _minx, _miny, _minz;
181 
183  Real _maxx, _maxy, _maxz;
184 
186  unsigned indexFromPoint(const Point & p) const;
187 
189  unsigned indexFromIndex(unsigned int var) const;
190 
192  void buildNodeWeightMaps();
193 };
EBSDReader::getFeatureID
unsigned int getFeatureID(unsigned int global_id) const
Return the EBSD feature id for a given (global) grain number.
Definition: EBSDReader.h:89
EBSDReader::_bins
unsigned int _bins
number of bins for each quaternion component
Definition: EBSDReader.h:168
EBSDReader::meshChanged
void meshChanged()
Maps need to be updated when the mesh changes.
Definition: EBSDReader.C:453
EBSDReader::_maxz
Real _maxz
Definition: EBSDReader.h:183
EBSDReader::_data
std::vector< EBSDPointData > _data
Logically three-dimensional data indexed by geometric points in a 1D vector.
Definition: EBSDReader.h:141
EulerAngleProvider.h
EBSDReader::~EBSDReader
virtual ~EBSDReader()
Definition: EBSDReader.C:342
EBSDReader::getGlobalID
virtual unsigned int getGlobalID(unsigned int phase, unsigned int local_id) const
Return the (global) grain id for a given phase and (local) grain number.
Definition: EBSDReader.h:95
EBSDReader::finalize
virtual void finalize()
Definition: EBSDReader.h:46
EBSDReader::_mesh
MooseMesh & _mesh
MooseMesh Variables.
Definition: EBSDReader.h:126
EBSDReader::_top_right
Point _top_right
Definition: EBSDReader.h:133
EBSDReader::indexFromIndex
unsigned indexFromIndex(unsigned int var) const
Transfer the index into the _avg_data array from given index.
Definition: EBSDReader.C:415
EBSDReader::_time_step
const int & _time_step
current timestep. Maps are only rebuild on mesh change during time step zero
Definition: EBSDReader.h:162
EBSDReader::_dz
Real _dz
Definition: EBSDReader.h:177
EBSDReader::getNodeToGrainWeightMap
const std::map< dof_id_type, std::vector< Real > > & getNodeToGrainWeightMap() const
Returns a map consisting of the node index followd by a vector of all grain weights for that node.
Definition: EBSDReader.C:432
EBSDReader::getNodeToPhaseWeightMap
const std::map< dof_id_type, std::vector< Real > > & getNodeToPhaseWeightMap() const
Returns a map consisting of the node index followd by a vector of all phase weights for that node.
Definition: EBSDReader.C:438
EBSDReader::_node_to_grain_weight_map
std::map< dof_id_type, std::vector< Real > > _node_to_grain_weight_map
Map of grain weights per node.
Definition: EBSDReader.h:156
EBSDReader::_ny
unsigned _ny
Definition: EBSDReader.h:174
EBSDReader::_dy
Real _dy
Definition: EBSDReader.h:177
EBSDReader::_custom_columns
unsigned int _custom_columns
number of additional custom data columns
Definition: EBSDReader.h:138
EBSDReader::_mesh_dimension
unsigned int _mesh_dimension
Dimension of the problem domain.
Definition: EBSDReader.h:165
EBSDReader::_nl
NonlinearSystemBase & _nl
Definition: EBSDReader.h:127
EBSDReader::getGrainNum
virtual unsigned int getGrainNum() const
Return the total number of grains.
Definition: EBSDReader.C:369
EBSDReader::_range
Point _range
Definition: EBSDReader.h:134
EBSDAccessFunctors
Mix-in class that adds so called access functors to select a field from an EBSDPointData or EBSDPoint...
Definition: EBSDAccessFunctors.h:23
EBSDReader::_avg_angles
std::vector< EulerAngles > _avg_angles
Euler Angles by (global) grain ID.
Definition: EBSDReader.h:147
EulerAngleProvider
Abstract base class for user objects that implement the Euler Angle provider interface.
Definition: EulerAngleProvider.h:24
EBSDReader
A GeneralUserObject that reads an EBSD file and stores the centroid data in a data structure which in...
Definition: EBSDReader.h:36
EBSDReader::_global_id
std::vector< std::vector< unsigned int > > _global_id
global ID for given phases and grains
Definition: EBSDReader.h:153
EBSDReader::_maxx
Real _maxx
Maximum grid extent.
Definition: EBSDReader.h:183
EBSDReader::_grain_num
unsigned int _grain_num
Variables needed to determine reduced order parameter values.
Definition: EBSDReader.h:131
EBSDReader::_minx
Real _minx
Grid origin.
Definition: EBSDReader.h:180
EBSDReader::readFile
virtual void readFile()
Definition: EBSDReader.C:59
EBSDReader::execute
virtual void execute()
Definition: EBSDReader.h:45
EBSDReader::_maxy
Real _maxy
Definition: EBSDReader.h:183
EBSDReader::buildNodeWeightMaps
void buildNodeWeightMaps()
Build grain and phase weight maps.
Definition: EBSDReader.C:461
EBSDReader::getPhaseNum
virtual unsigned int getPhaseNum() const
Return the total number of phases.
Definition: EBSDReader.h:76
EBSDReader::_nz
unsigned _nz
Definition: EBSDReader.h:174
EBSDReader::_minz
Real _minz
Definition: EBSDReader.h:180
EBSDReader::_miny
Real _miny
Definition: EBSDReader.h:180
EBSDReader::initialize
virtual void initialize()
Definition: EBSDReader.h:44
EBSDReader::getEulerAngles
virtual const EulerAngles & getEulerAngles(unsigned int) const
EulerAngleProvider interface implementation to fetch a triplet of Euler angles.
Definition: EBSDReader.C:357
EBSDReader::getAvgData
const EBSDAvgData & getAvgData(unsigned int i) const
Get the requested type of average data for (global) grain number i.
Definition: EBSDReader.C:351
EBSDReader::getAvgDataAccessFunctor
MooseSharedPointer< EBSDAvgDataFunctor > getAvgDataAccessFunctor(const MooseEnum &field_name) const
Factory function to return a average functor specified by name.
Definition: EBSDReader.C:556
EBSDReader::getData
const EBSDPointData & getData(const Point &p) const
Get the requested type of data at the point p.
Definition: EBSDReader.C:345
EBSDReader::getPointDataAccessFunctor
MooseSharedPointer< EBSDPointDataFunctor > getPointDataAccessFunctor(const MooseEnum &field_name) const
Factory function to return a point functor specified by name.
Definition: EBSDReader.C:510
EBSDReader::_bottom_left
Point _bottom_left
Definition: EBSDReader.h:132
EBSDReader::getFeatureID
unsigned int getFeatureID(unsigned int phase, unsigned int local_id) const
Return the EBSD feature id for a given phase and phase (local) grain number.
Definition: EBSDReader.h:84
validParams< EBSDReader >
InputParameters validParams< EBSDReader >()
Definition: EBSDReader.C:23
EBSDReader::_dx
Real _dx
The spacing of the values in x, y and z directions.
Definition: EBSDReader.h:177
EBSDReader::_nx
unsigned _nx
The number of values in the x, y and z directions.
Definition: EBSDReader.h:174
EBSDReader::_node_to_phase_weight_map
std::map< dof_id_type, std::vector< Real > > _node_to_phase_weight_map
Map of phase weights per node.
Definition: EBSDReader.h:159
EulerAngles
Euler angle triplet.
Definition: EulerAngles.h:22
EBSDReader::indexFromPoint
unsigned indexFromPoint(const Point &p) const
Computes a global index in the _data array given an input centroid point.
Definition: EBSDReader.C:381
EBSDReader::EBSDReader
EBSDReader(const InputParameters &params)
Definition: EBSDReader.C:35
EBSDAccessFunctors.h
EBSDReader::_L_norm
unsigned int _L_norm
L_norm value for averaging.
Definition: EBSDReader.h:171
EBSDReader::_avg_data
std::vector< EBSDAvgData > _avg_data
Averages by (global) grain ID.
Definition: EBSDReader.h:144
EBSDReader::_global_id_map
std::map< unsigned int, unsigned int > _global_id_map
map from feature_id to global_id
Definition: EBSDReader.h:150