PRISMS-PF Manual v3.0-pre
Loading...
Searching...
No Matches
SimplifiedGrainRepresentation.h
1#ifndef INCLUDE_SIMPLIFIEDGRAINREPRESENTATION_H_
2#define INCLUDE_SIMPLIFIEDGRAINREPRESENTATION_H_
3
4#include <grains/FloodFiller.h>
5
13template <int dim>
15{
16public:
24
28 [[nodiscard]] dealii::Point<dim>
29 getCenter() const;
30
34 [[nodiscard]] double
35 getRadius() const;
36
40 [[nodiscard]] unsigned int
41 getGrainId() const;
42
46 void
47 setGrainId(unsigned int _grain_id);
48
52 [[nodiscard]] unsigned int
53 getOrderParameterId() const;
54
58 void
59 setOrderParameterId(unsigned int _order_parameter_id);
60
65 [[nodiscard]] unsigned int
67
72 void
73 setDistanceToNeighbor(double dist);
74
79 [[nodiscard]] double
81
82protected:
86 dealii::Point<dim> center;
87
91 double radius;
92
97 unsigned int grain_id;
98
102 unsigned int order_parameter_id;
103
109
116};
117
125template <int dim>
127{
128public:
133 void
134 reassignGrains(std::vector<SimplifiedGrainRepresentation<dim>> &grain_representations,
135 double buffer_distance,
136 std::vector<unsigned int> &order_parameter_id_list);
137
144 void
146 const std::vector<SimplifiedGrainRepresentation<dim>> &old_grain_representations,
147 std::vector<SimplifiedGrainRepresentation<dim>> &new_grain_representations) const;
148
149protected:
150};
151
152#endif
Definition FloodFiller.h:20
Definition SimplifiedGrainRepresentation.h:127
void reassignGrains(std::vector< SimplifiedGrainRepresentation< dim > > &grain_representations, double buffer_distance, std::vector< unsigned int > &order_parameter_id_list)
Definition SimplifiedGrainRepresentation.cc:146
void transferGrainIds(const std::vector< SimplifiedGrainRepresentation< dim > > &old_grain_representations, std::vector< SimplifiedGrainRepresentation< dim > > &new_grain_representations) const
Definition SimplifiedGrainRepresentation.cc:273
Definition SimplifiedGrainRepresentation.h:15
unsigned int order_parameter_id
Definition SimplifiedGrainRepresentation.h:102
double radius
Definition SimplifiedGrainRepresentation.h:91
double getRadius() const
Definition SimplifiedGrainRepresentation.cc:86
dealii::Point< dim > center
Definition SimplifiedGrainRepresentation.h:86
void setOrderParameterId(unsigned int _order_parameter_id)
Definition SimplifiedGrainRepresentation.cc:114
unsigned int getGrainId() const
Definition SimplifiedGrainRepresentation.cc:93
unsigned int old_order_parameter_id
Definition SimplifiedGrainRepresentation.h:108
dealii::Point< dim > getCenter() const
Definition SimplifiedGrainRepresentation.cc:79
void setDistanceToNeighbor(double dist)
Definition SimplifiedGrainRepresentation.cc:128
unsigned int getOrderParameterId() const
Definition SimplifiedGrainRepresentation.cc:107
double distance_to_neighbor_sharing_op
Definition SimplifiedGrainRepresentation.h:115
unsigned int grain_id
Definition SimplifiedGrainRepresentation.h:97
void setGrainId(unsigned int _grain_id)
Definition SimplifiedGrainRepresentation.cc:100
double getDistanceToNeighbor() const
Definition SimplifiedGrainRepresentation.cc:135
unsigned int getOldOrderParameterId() const
Definition SimplifiedGrainRepresentation.cc:121