www.mooseframework.org
EBSDReader.h
Go to the documentation of this file.
1 /****************************************************************/
2 /* MOOSE - Multiphysics Object Oriented Simulation Environment */
3 /* */
4 /* All contents are licensed under LGPL V2.1 */
5 /* See LICENSE for full restrictions */
6 /****************************************************************/
7 #ifndef EBSDREADER_H
8 #define EBSDREADER_H
9 
10 #include "EulerAngleProvider.h"
11 #include "EBSDAccessFunctors.h"
12 
13 class EBSDReader;
14 
15 template <>
16 InputParameters validParams<EBSDReader>();
17 
35 {
36 public:
37  EBSDReader(const InputParameters & params);
38  virtual ~EBSDReader();
39 
40  virtual void readFile();
41 
42  virtual void initialize() {}
43  virtual void execute() {}
44  virtual void finalize() {}
45 
49  const EBSDPointData & getData(const Point & p) const;
50 
54  const EBSDAvgData & getAvgData(unsigned int i) const;
55 
59  const EBSDAvgData & getAvgData(unsigned int phase, unsigned int local_id) const;
60 
64  virtual const EulerAngles & getEulerAngles(unsigned int) const;
65 
69  virtual unsigned int getGrainNum() const;
70 
74  virtual unsigned int getPhaseNum() const { return _global_id.size(); }
75 
79  unsigned int getGrainNum(unsigned int phase) const;
80 
82  unsigned int getFeatureID(unsigned int phase, unsigned int local_id) const
83  {
84  return _avg_data[_global_id[phase][local_id]]._feature_id;
85  }
87  unsigned int getFeatureID(unsigned int global_id) const
88  {
89  return _avg_data[global_id]._feature_id;
90  }
91 
93  unsigned int getGlobalID(unsigned int phase, unsigned int local_id) const
94  {
95  return _global_id[phase][local_id];
96  }
98  unsigned int getGlobalID(unsigned int feature_id) const;
99 
101  MooseSharedPointer<EBSDPointDataFunctor>
102  getPointDataAccessFunctor(const MooseEnum & field_name) const;
104  MooseSharedPointer<EBSDAvgDataFunctor>
105  getAvgDataAccessFunctor(const MooseEnum & field_name) const;
106 
111  const std::map<dof_id_type, std::vector<Real>> & getNodeToGrainWeightMap() const;
112 
117  const std::map<dof_id_type, std::vector<Real>> & getNodeToPhaseWeightMap() const;
118 
120  void meshChanged();
121 
122 protected:
124  MooseMesh & _mesh;
125  NonlinearSystemBase & _nl;
127 
129  unsigned int _grain_num;
131  Point _top_right;
132  Point _range;
134 
136  unsigned int _custom_columns;
137 
139  std::vector<EBSDPointData> _data;
140 
142  std::vector<EBSDAvgData> _avg_data;
143 
145  std::vector<EulerAngles> _avg_angles;
146 
148  std::map<unsigned int, unsigned int> _global_id_map;
149 
151  std::vector<std::vector<unsigned int>> _global_id;
152 
154  std::map<dof_id_type, std::vector<Real>> _node_to_grain_weight_map;
155 
157  std::map<dof_id_type, std::vector<Real>> _node_to_phase_weight_map;
158 
160  const int & _time_step;
161 
163  unsigned int _mesh_dimension;
164 
166  unsigned _nx, _ny, _nz;
167 
169  Real _dx, _dy, _dz;
170 
172  Real _minx, _miny, _minz;
173 
175  Real _maxx, _maxy, _maxz;
176 
178  unsigned indexFromPoint(const Point & p) const;
179 
181  unsigned indexFromIndex(unsigned int var) const;
182 
184  void buildNodeWeightMaps();
185 };
186 
187 #endif // EBSDREADER_H
unsigned int getFeatureID(unsigned int global_id) const
Return the EBSD feature id for a given (global) grain number.
Definition: EBSDReader.h:87
Point _bottom_left
Definition: EBSDReader.h:130
std::map< dof_id_type, std::vector< Real > > _node_to_phase_weight_map
Map of phase weights per node.
Definition: EBSDReader.h:157
unsigned int _custom_columns
number of additional custom data columns
Definition: EBSDReader.h:136
Real _minx
Grid origin.
Definition: EBSDReader.h:172
Point _range
Definition: EBSDReader.h:132
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:93
Real _maxx
Maximum grid extent.
Definition: EBSDReader.h:175
unsigned indexFromPoint(const Point &p) const
Computes a global index in the _data array given an input centroid point.
Definition: EBSDReader.C:260
unsigned _nz
Definition: EBSDReader.h:166
virtual const EulerAngles & getEulerAngles(unsigned int) const
EulerAngleProvider interface implementation to fetch a triplet of Euler angles.
Definition: EBSDReader.C:236
Real _miny
Definition: EBSDReader.h:172
unsigned int _grain_num
Variables needed to determine reduced order parameter values.
Definition: EBSDReader.h:129
unsigned _nx
The number of values in the x, y and z directions.
Definition: EBSDReader.h:166
virtual void execute()
Definition: EBSDReader.h:43
MooseMesh & _mesh
MooseMesh Variables.
Definition: EBSDReader.h:124
Real _minz
Definition: EBSDReader.h:172
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:306
Per element EBSD data point.
Mix-in class that adds so called access functors to select a field from an EBSDPointData or EBSDPoint...
std::map< unsigned int, unsigned int > _global_id_map
map from feature_id to global_id
Definition: EBSDReader.h:148
const EBSDPointData & getData(const Point &p) const
Get the requested type of data at the point p.
Definition: EBSDReader.C:224
virtual unsigned int getPhaseNum() const
Return the total number of phases.
Definition: EBSDReader.h:74
const int & _time_step
current timestep. Maps are only rebuild on mesh change during time step zero
Definition: EBSDReader.h:160
A GeneralUserObject that reads an EBSD file and stores the centroid data in a data structure which in...
Definition: EBSDReader.h:34
unsigned _ny
Definition: EBSDReader.h:166
Real _dx
The spacing of the values in x, y and z directions.
Definition: EBSDReader.h:169
unsigned int _mesh_dimension
Dimension of the problem domain.
Definition: EBSDReader.h:163
InputParameters validParams< EBSDReader >()
Definition: EBSDReader.C:18
EBSDReader(const InputParameters &params)
Definition: EBSDReader.C:28
virtual unsigned int getGrainNum() const
Return the total number of grains.
Definition: EBSDReader.C:248
MooseSharedPointer< EBSDAvgDataFunctor > getAvgDataAccessFunctor(const MooseEnum &field_name) const
Factory function to return a average functor specified by name.
Definition: EBSDReader.C:432
std::vector< EBSDAvgData > _avg_data
Averages by (global) grain ID.
Definition: EBSDReader.h:142
const EBSDAvgData & getAvgData(unsigned int i) const
Get the requested type of average data for (global) grain number i.
Definition: EBSDReader.C:230
Point _top_right
Definition: EBSDReader.h:131
Real _maxz
Definition: EBSDReader.h:175
std::map< dof_id_type, std::vector< Real > > _node_to_grain_weight_map
Map of grain weights per node.
Definition: EBSDReader.h:154
MooseSharedPointer< EBSDPointDataFunctor > getPointDataAccessFunctor(const MooseEnum &field_name) const
Factory function to return a point functor specified by name.
Definition: EBSDReader.C:386
Euler angle triplet.
Definition: EulerAngles.h:20
virtual ~EBSDReader()
Definition: EBSDReader.C:221
NonlinearSystemBase & _nl
Definition: EBSDReader.h:125
Real _maxy
Definition: EBSDReader.h:175
Abstract base class for user objects that implement the Euler Angle provider interface.
std::vector< EulerAngles > _avg_angles
Euler Angles by (global) grain ID.
Definition: EBSDReader.h:145
void buildNodeWeightMaps()
Build grain and phase weight maps.
Definition: EBSDReader.C:335
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:82
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:312
void meshChanged()
Maps need to be updated when the mesh changes.
Definition: EBSDReader.C:327
virtual void readFile()
Definition: EBSDReader.C:47
std::vector< std::vector< unsigned int > > _global_id
global ID for given phases and grains
Definition: EBSDReader.h:151
std::vector< EBSDPointData > _data
Logically three-dimensional data indexed by geometric points in a 1D vector.
Definition: EBSDReader.h:139
unsigned indexFromIndex(unsigned int var) const
Transfer the index into the _avg_data array from given index.
Definition: EBSDReader.C:289
virtual void initialize()
Definition: EBSDReader.h:42
virtual void finalize()
Definition: EBSDReader.h:44