Alamo
PSRead.H
Go to the documentation of this file.
1//
2// Fill a domain (region where field=0) with packed spheres (regions where field=1).
3// Sphere locations and radii are determined from an xyzr file.
4//
5
6#ifndef IC_PSREAD_H_
7#define IC_PSREAD_H_
8
9#include "Set/Set.H"
10#include "IC/IC.H"
11#include "Util/Util.H"
12using namespace std;
13#include <iostream>
14#include <vector>
15#include <algorithm>
16#include <iterator>
17#include "mpi.h"
18#include "IO/ParmParse.H"
19
20namespace IC
21{
22class PSRead : public IC<Set::Scalar>
23{
24public:
25 static constexpr const char* name = "psread";
26
27 enum Type
28 {
30 Values
31 };
32
33 PSRead(amrex::Vector<amrex::Geometry>& _geom) : IC(_geom) {}
34
35 PSRead(amrex::Vector<amrex::Geometry>& _geom, IO::ParmParse& pp, std::string name) : IC(_geom)
36 {
37 pp.queryclass(name, *this);
38 }
39 void Define() {
40
41 };
42
43 void Add(const int& lev, Set::Field<Set::Scalar>& a_phi, Set::Scalar)
44 {
45 Set::Vector size = Set::Size(geom[lev]);
46 int ncomp = a_phi[lev]->nComp();
47
48 for (amrex::MFIter mfi(*a_phi[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
49 {
50 amrex::Box bx;
51 amrex::IndexType type = a_phi[lev]->ixType();
52 if (type == amrex::IndexType::TheCellType()) bx = mfi.growntilebox();
53 else if (type == amrex::IndexType::TheNodeType()) bx = mfi.grownnodaltilebox();
54 else Util::Abort(INFO, "Unkonwn index type");
55
56 amrex::Array4<Set::Scalar> const& phi = a_phi[lev]->array(mfi);
57 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
58 {
59 Set::Vector x = Set::Position(i, j, k, geom[lev], type);
60 Set::Scalar min_grain_id = 0;
61
62 for (unsigned int n = 0; n < X.size(); n++)
63 {
65 Set::Scalar d1 = std::numeric_limits<Set::Scalar>::infinity();
66
67 d = (x - X[n]).lpNorm<2>();
68
69 if (geom[0].isPeriodic(0))
70 {
71 d1 = std::min((x - X[n] + size(0) * Set::Vector::Unit(0)).lpNorm<2>(),
72 (x - X[n] - size(0) * Set::Vector::Unit(0)).lpNorm<2>());
73 }
74#if AMREX_SPACEDIM > 1
75 if (geom[0].isPeriodic(1))
76 {
77 d1 = std::min((x - X[n] + size(1) * Set::Vector::Unit(1)).lpNorm<2>(),
78 (x - X[n] - size(1) * Set::Vector::Unit(1)).lpNorm<2>());
79 }
80#endif
81#if AMREX_SPACEDIM > 2
82 if (geom[0].isPeriodic(2))
83 {
84 d1 = std::min((x - X[n] + size(2) * Set::Vector::Unit(2)).lpNorm<2>(),
85 (x - X[n] - size(2) * Set::Vector::Unit(2)).lpNorm<2>());
86 }
87#endif
88
89 d = std::min(d, d1);
90
91 if (d <= (R[n] + eps))
92
93 {
94 Set::Scalar m = 0.5 * (1 + erf((-d + R[n]) / (eps)));
95 min_grain_id = min_grain_id + m * (1. - min_grain_id);
96 }
97 }
98
99 phi(i, j, k, 0) = min_grain_id;
100 if (invert) phi(i,j,k,0) = 1.0 - min_grain_id;
101 if (ncomp > 1) phi(i, j, k, 1) = 1.0 - min_grain_id;
102 });
103 }
104 }
105
106 static void Parse(PSRead& value, IO::ParmParse& pp)
107 {
108 std::string filename;
109 int verbose = 0;
110 // Diffuseness of the sphere boundary
111 pp.query_default("eps", value.eps, "0.0", Unit::Length());
112
113 pp.forbid("filename", "use file.name instead");
114 // Location of .xyzr file
115 pp.query_file("file.name", filename);
116
117 // Verbosity (used in parser only)
118 pp.query_default("verbose", verbose, 0);
119
120 // Unitless coordinate multiplier
121 pp.query_default("mult",value.mult, "1.0", Unit::Less());
122
123 // Whether to invert (1 outside instead of inside spheres)
124 pp.query_default("invert",value.invert,false);
125
126 // X-offset
127 pp.queryarr_default("x0",value.x0,"0 0 0", Unit::Length());
128
129 if (pp.contains("x0"))
130 pp_queryarr("x0",value.x0); // Coordinate offset
131
132 Unit unit = Unit::Length();
133 // Units of length in the file
134 pp.queryunit("file.unit",unit);
136 "Unit must be of type length but got unit ",unit.normalized_unitstring());
137
138 std::ifstream datafile(filename);
139 std::string line;
140 if (datafile.is_open())
141 {
142 value.X.clear();
143 value.R.clear();
144
145 while (getline(datafile, line))
146 {
147 std::istringstream in(line);
148
149 std::string strx, stry, strz, strR;
150 in >> strx >> stry >> strz >> strR;
151
152 Set::Scalar x = (std::stod(strx) * unit).normalized_value();
153 Set::Scalar y = (std::stod(stry) * unit).normalized_value();
154 #if AMREX_SPACEDIM > 2
155 Set::Scalar z = (std::stod(strz) * unit).normalized_value();
156 #endif
157 Set::Scalar r = (std::stod(strR) * unit).normalized_value();
158
159 Set::Vector X(AMREX_D_DECL(x,y,z));
160 X = value.x0 + value.mult*X;
161
162 value.X.push_back(X);
163 value.R.push_back(r);
164 if (verbose > 0)
165 Util::Message(INFO, "x=", value.X.back().transpose(), " r=", value.R.back());
166 }
167 datafile.close();
168 }
169 else
170 {
171 Util::Abort(INFO, "Unable to open file ", filename);
172 }
173 }
174
175private:
176 std::vector<Set::Vector> X;
177 std::vector<Set::Scalar> R;
180 Set::Vector x0 = Set::Vector::Zero();
181 bool invert = false;
182};
183}
184#endif
#define pp_queryarr(...)
Definition ParmParse.H:105
#define TEST(x)
Definition Util.H:22
#define INFO
Definition Util.H:21
amrex::Vector< amrex::Geometry > & geom
Definition IC.H:61
PSRead(amrex::Vector< amrex::Geometry > &_geom, IO::ParmParse &pp, std::string name)
Definition PSRead.H:35
static void Parse(PSRead &value, IO::ParmParse &pp)
Definition PSRead.H:106
std::vector< Set::Vector > X
Definition PSRead.H:176
void Define()
Definition PSRead.H:39
bool invert
Definition PSRead.H:181
@ Partition
Definition PSRead.H:29
PSRead(amrex::Vector< amrex::Geometry > &_geom)
Definition PSRead.H:33
Set::Scalar eps
Definition PSRead.H:178
void Add(const int &lev, Set::Field< Set::Scalar > &a_phi, Set::Scalar)
Definition PSRead.H:43
static constexpr const char * name
Definition PSRead.H:25
Set::Scalar mult
Definition PSRead.H:179
Set::Vector x0
Definition PSRead.H:180
std::vector< Set::Scalar > R
Definition PSRead.H:177
bool contains(std::string name)
Definition ParmParse.H:173
int queryunit(std::string name, Unit &value)
Definition ParmParse.H:192
int queryarr_default(std::string name, std::vector< std::string > &value, std::vector< std::string > defaultvalue)
Definition ParmParse.H:660
void forbid(std::string name, std::string explanation)
Definition ParmParse.H:160
int query_file(std::string name, std::string &value, bool copyfile, bool checkfile)
Definition ParmParse.H:486
void queryclass(std::string name, T *value)
Definition ParmParse.H:960
int query_default(std::string name, T &value, T defaultvalue)
Definition ParmParse.H:293
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
AMREX_FORCE_INLINE Vector Position(const int &i, const int &j, const int &k, const amrex::Geometry &geom, const amrex::IndexType &ixType)
Definition Base.H:120
AMREX_FORCE_INLINE Vector Size(const amrex::Geometry &geom)
Definition Base.H:159
AMREX_FORCE_INLINE void AssertException(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition Util.H:251
void Abort(const char *msg)
Definition Util.cpp:243
void Message(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:126
Definition Unit.H:21
bool isType(const Unit &test) const
Definition Unit.H:425
std::string normalized_unitstring() const
Definition Unit.H:510
static Unit Length()
Definition Unit.H:198
static Unit Less()
Definition Unit.H:197