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"
11using namespace std;
12#include <iostream>
13#include <vector>
14#include <algorithm>
15#include <iterator>
16#include "mpi.h"
17#include "IO/ParmParse.H"
18
19namespace IC
20{
21class PSRead : public IC<Set::Scalar>
22{
23public:
24 static constexpr const char* name = "psread";
25
26 enum Type
27 {
29 Values
30 };
31
32 PSRead(amrex::Vector<amrex::Geometry>& _geom) : IC(_geom) {}
33
34 PSRead(amrex::Vector<amrex::Geometry>& _geom, IO::ParmParse& pp, std::string name) : IC(_geom)
35 {
36 pp_queryclass(name, *this);
37 }
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 pp_query("eps", value.eps); // Diffuseness of the sphere boundary
111 pp_query_file("filename", filename); // Location of .xyzr file
112 pp_query("verbose", verbose); // Verbosity (used in parser only)
113 pp_query("mult",value.mult); // Coordinate multiplier
114 pp.query("invert",value.invert); // Coordinate multiplier
115 if (pp.contains("x0"))
116 pp_queryarr("x0",value.x0); // Coordinate offset
117 std::ifstream datafile(filename);
118 std::string line;
119 if (datafile.is_open())
120 {
121 value.X.clear();
122 value.R.clear();
123
124 while (getline(datafile, line))
125 {
126 std::istringstream in(line);
127
128 std::string strx, stry, strz, strR;
129 in >> strx >> stry >> strz >> strR;
130 value.X.push_back(value.x0 + value.mult*Set::Vector(AMREX_D_DECL(std::stod(strx), std::stod(stry), std::stod(strz))));
131 value.R.push_back(value.mult*std::stod(strR));
132 if (verbose > 0)
133 Util::Message(INFO, "x=", value.X.back().transpose(), " r=", value.R.back());
134 }
135 datafile.close();
136 }
137 else
138 {
139 Util::Abort(INFO, "Unable to open file ", filename);
140 }
141 }
142
143private:
144 std::vector<Set::Vector> X;
145 std::vector<Set::Scalar> R;
148 Set::Vector x0 = Set::Vector::Zero();
149 bool invert = false;
150};
151}
152#endif
#define pp_queryarr(...)
Definition ParmParse.H:103
#define pp_query(...)
Definition ParmParse.H:106
#define pp_query_file(...)
Definition ParmParse.H:102
#define pp_queryclass(...)
Definition ParmParse.H:107
#define INFO
Definition Util.H:20
amrex::Vector< amrex::Geometry > & geom
Definition IC.H:61
PSRead(amrex::Vector< amrex::Geometry > &_geom, IO::ParmParse &pp, std::string name)
Definition PSRead.H:34
static void Parse(PSRead &value, IO::ParmParse &pp)
Definition PSRead.H:106
std::vector< Set::Vector > X
Definition PSRead.H:144
void Define()
Definition PSRead.H:39
bool invert
Definition PSRead.H:149
@ Partition
Definition PSRead.H:28
PSRead(amrex::Vector< amrex::Geometry > &_geom)
Definition PSRead.H:32
Set::Scalar eps
Definition PSRead.H:146
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:24
Set::Scalar mult
Definition PSRead.H:147
Set::Vector x0
Definition PSRead.H:148
std::vector< Set::Scalar > R
Definition PSRead.H:145
bool contains(std::string name)
Definition ParmParse.H:154
Initialize a spherical inclusion.
Definition BMP.H:19
amrex::Real Scalar
Definition Base.H:19
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:20
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:121
AMREX_FORCE_INLINE Vector Size(const amrex::Geometry &geom)
Definition Base.H:161
void Abort(const char *msg)
Definition Util.cpp:170
void Message(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:141