Alamo
Voronoi.H
Go to the documentation of this file.
1// Initialize an N component field with a Voronoi tessellation.
2// (Mostly used by :code:`Integrator::PhaseFieldMicrostructure`).
3
4#ifndef IC_VORONOI_H_
5#define IC_VORONOI_H_
6
7#include "Set/Set.H"
8#include "IC/IC.H"
9#include "IO/ParmParse.H"
10
11namespace IC
12{
13class Voronoi : public IC<Set::Scalar>
14{
15public:
16 static constexpr const char* name = "voronoi";
17
19
20 Voronoi (amrex::Vector<amrex::Geometry> &_geom) : IC(_geom) {}
21 Voronoi (amrex::Vector<amrex::Geometry> &_geom,IO::ParmParse &pp, std::string name) : IC(_geom)
22 {pp_queryclass(name,*this);}
23 Voronoi (amrex::Vector<amrex::Geometry> &_geom, int _number_of_grains, Set::Scalar a_alpha) : IC(_geom)
24 {
25 Define(_number_of_grains,a_alpha);
26 }
27 Voronoi (amrex::Vector<amrex::Geometry> &a_geom, int a_number_of_grains) : IC(a_geom)
28 {
29 Define(a_number_of_grains,1.0);
30 }
31
32 void Define (int a_number_of_grains, Set::Scalar a_alpha)
33 {
34 Define(a_number_of_grains, std::vector<Set::Scalar>(a_number_of_grains,a_alpha),Type::Partition);
35 }
36
37 void Define (int a_number_of_grains,
38 std::vector<Set::Scalar> a_alpha,
39 Type a_type = Type::Values)
40 {
41 number_of_grains = a_number_of_grains;
42 alpha = a_alpha;
43 type = a_type;
44
46 srand(seed);
47 for (int n = 0; n<number_of_grains; n++)
48 {
49 AMREX_D_TERM(voronoi[n](0) = geom[0].ProbLo(0) + (geom[0].ProbHi(0)-geom[0].ProbLo(0))*Util::Random();,
50 voronoi[n](1) = geom[0].ProbLo(1) + (geom[0].ProbHi(1)-geom[0].ProbLo(1))*Util::Random();,
51 voronoi[n](2) = geom[0].ProbLo(2) + (geom[0].ProbHi(2)-geom[0].ProbLo(2))*Util::Random(););
52 }
53 };
54
55 void Add(const int &lev, Set::Field<Set::Scalar> &a_field, Set::Scalar)
56 {
57 Set::Vector size;
58 AMREX_D_TERM(size(0) = geom[0].ProbHi()[0] - geom[0].ProbLo()[0];,
59 size(1) = geom[0].ProbHi()[1] - geom[0].ProbLo()[1];,
60 size(2) = geom[0].ProbHi()[2] - geom[0].ProbLo()[2];)
61
62 for (amrex::MFIter mfi(*a_field[lev],amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
63 {
64 amrex::Box bx;
65 amrex::IndexType ixtype = a_field[lev]->ixType();
66 if (ixtype == amrex::IndexType::TheCellType()) bx = mfi.growntilebox();
67 else if (ixtype == amrex::IndexType::TheNodeType()) bx = mfi.grownnodaltilebox();
68
69 int ncomp = a_field[lev]->nComp();
70 amrex::Array4<Set::Scalar> const& field = a_field[lev]->array(mfi);
71 amrex::ParallelFor (bx,[=] AMREX_GPU_DEVICE(int i, int j, int k) {
72
74 AMREX_D_TERM(x(0) = geom[lev].ProbLo()[0] + ((amrex::Real)(i) + 0.5) * geom[lev].CellSize()[0];,
75 x(1) = geom[lev].ProbLo()[1] + ((amrex::Real)(j) + 0.5) * geom[lev].CellSize()[1];,
76 x(2) = geom[lev].ProbLo()[2] + ((amrex::Real)(k) + 0.5) * geom[lev].CellSize()[2];);
77
78 amrex::Real min_distance = std::numeric_limits<amrex::Real>::infinity();
79 int min_grain_id = -1;
80
81 for (int n = 0; n<number_of_grains; n++)
82 {
83 Set::Scalar d = (x - voronoi[n]).lpNorm<2>();
84
85 if (geom[0].isPeriodic(0))
86 {
87 d = std::min(d,
88 std::min( (x-voronoi[n] + size(0)*Set::Vector::Unit(0)).lpNorm<2>(),
89 (x-voronoi[n] - size(0)*Set::Vector::Unit(0)).lpNorm<2>()));
90 }
91#if AMREX_SPACEDIM>1
92 if (geom[0].isPeriodic(1))
93 {
94 d = std::min(d,
95 std::min( (x-voronoi[n] + size(0)*Set::Vector::Unit(1)).lpNorm<2>(),
96 (x-voronoi[n] - size(0)*Set::Vector::Unit(1)).lpNorm<2>()));
97 }
98#endif
99#if AMREX_SPACEDIM>2
100 if (geom[0].isPeriodic(2))
101 {
102 d = std::min(d,
103 std::min( (x-voronoi[n] + size(0)*Set::Vector::Unit(2)).lpNorm<2>(),
104 (x-voronoi[n] - size(0)*Set::Vector::Unit(2)).lpNorm<2>()));
105 }
106#endif
107 if (d<min_distance)
108 {
109 min_distance = d;
110 min_grain_id = n;
111 }
112 }
113
114 if (type == Type::Values) field(i,j,k) = alpha[min_grain_id];
115 else if (type == Type::Partition) field(i,j,k,min_grain_id % ncomp) = alpha[min_grain_id];
116 });
117 }
118 }
119
120 static void Parse(Voronoi &value, IO::ParmParse &pp)
121 {
122 pp_query("number_of_grains",value.number_of_grains); // Number of grains
123 if (pp.contains("alpha"))
124 {
125 pp_queryarr("alpha",value.alpha); // Value to take in the region [1.0]
126 value.Define(value.number_of_grains,value.alpha);
127 }
128 else
129 {
130 value.Define(value.number_of_grains,1.0);
131 }
132 pp_query("seed",value.seed); // Random seed to use
133 }
134
135private:
137 int seed = 1;
138 std::vector<Set::Scalar> alpha;
139 std::vector<Set::Vector> voronoi;
141};
142}
143#endif
#define pp_queryarr(...)
Definition ParmParse.H:105
#define pp_query(...)
Definition ParmParse.H:108
#define pp_queryclass(...)
Definition ParmParse.H:109
amrex::Vector< amrex::Geometry > & geom
Definition IC.H:61
void Add(const int &lev, Set::Field< Set::Scalar > &a_field, Set::Scalar)
Definition Voronoi.H:55
static void Parse(Voronoi &value, IO::ParmParse &pp)
Definition Voronoi.H:120
std::vector< Set::Scalar > alpha
Definition Voronoi.H:138
int number_of_grains
Definition Voronoi.H:136
static constexpr const char * name
Definition Voronoi.H:16
Voronoi(amrex::Vector< amrex::Geometry > &_geom, int _number_of_grains, Set::Scalar a_alpha)
Definition Voronoi.H:23
void Define(int a_number_of_grains, std::vector< Set::Scalar > a_alpha, Type a_type=Type::Values)
Definition Voronoi.H:37
Type type
Definition Voronoi.H:140
Voronoi(amrex::Vector< amrex::Geometry > &_geom)
Definition Voronoi.H:20
Voronoi(amrex::Vector< amrex::Geometry > &a_geom, int a_number_of_grains)
Definition Voronoi.H:27
void Define(int a_number_of_grains, Set::Scalar a_alpha)
Definition Voronoi.H:32
std::vector< Set::Vector > voronoi
Definition Voronoi.H:139
Voronoi(amrex::Vector< amrex::Geometry > &_geom, IO::ParmParse &pp, std::string name)
Definition Voronoi.H:21
bool contains(std::string name)
Definition ParmParse.H:173
Initialize a spherical inclusion.
Definition BMP.H:20
amrex::Real Scalar
Definition Base.H:18
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:19
Set::Scalar Random()
Definition Set.cpp:34