26 return s.size() >= t.size() && s.substr(s.size() - t.size(), t.size()) == t;
30 std::istringstream iss(s);
48 std::ifstream fmplist;
49 fmplist.open(mplist.c_str(), std::ios::in);
52 std::cerr <<
" Could not open result file " << mplist
60 newMaterial.
medium =
nullptr;
62 fmplist >> newMaterial.
eps;
69 newMaterial.
medium =
nullptr;
70 newMaterial.
eps = newMaterial.
ohm = -1;
74 std::map<int, int> domain2material;
77 for (
int i = 0; i < d2msize; ++i) {
80 fmplist >> domain2material[domain];
86 fmesh.open(mesh.c_str(), std::ios::in);
89 std::cerr <<
" Could not open nodes file " << mesh <<
" for reading.\n";
94 std::getline(fmesh, line);
95 }
while (!
ends_with(line,
"# number of mesh points"));
99 std::cout <<
" Read " <<
nNodes <<
" nodes from file " << mesh <<
".\n";
101 std::getline(fmesh, line);
102 }
while (line !=
"# Mesh point coordinates");
103 double minx = 1e100, miny = 1e100, minz = 1e100, maxx = -1e100, maxy = -1e100,
105 for (
int i = 0; i <
nNodes; ++i) {
107 fmesh >> newNode.
x >> newNode.
y >> newNode.
z;
111 nodes.push_back(newNode);
112 minx = std::min(minx, newNode.
x);
113 maxx = std::max(maxx, newNode.
x);
114 miny = std::min(miny, newNode.
y);
115 maxy = std::max(maxy, newNode.
y);
116 minz = std::min(minz, newNode.
z);
117 maxz = std::max(maxz, newNode.
z);
119 std::cout << minx <<
" < x < " << maxx <<
"\n";
120 std::cout << miny <<
" < y < " << maxy <<
"\n";
121 std::cout << minz <<
" < z < " << maxz <<
"\n";
124 std::getline(fmesh, line);
125 }
while (line !=
"4 tet2 # type name");
127 std::getline(fmesh, line);
128 }
while (!
ends_with(line,
"# number of elements"));
132 std::cout <<
" Read " <<
nElements <<
" elements from file " << mesh
134 std::getline(fmesh, line);
137 int perm[10] = {0, 1, 2, 3, 4, 5, 7, 6, 8, 9};
141 for (
int j = 0; j < 10; ++j) {
142 fmesh >> newElement.
emap[perm[j]];
148 std::getline(fmesh, line);
149 }
while (line !=
"# Geometric entity indices");
153 elements[i].matmap = domain2material.count(domain) ? domain2material[domain]
158 std::map<Node, std::vector<int>,
nodeCmp> nodeIdx;
159 for (
int i = 0; i <
nNodes; ++i) {
160 nodeIdx[
nodes[i]].push_back(i);
162 std::cout <<
"Map size: " << nodeIdx.size() << std::endl;
164 std::ifstream ffield;
165 ffield.open(field.c_str(), std::ios::in);
168 std::cerr <<
" Could not open field potentials file " << field
169 <<
" for reading.\n";
173 std::getline(ffield, line);
174 }
while (line.substr(0, 81) !=
178 std::istringstream sline(line);
186 while (sline >> token) {
188 std::cout <<
" Reading data for weighting field " << token <<
".\n";
195 for (
int i = 0; i <
nNodes; ++i) {
197 ffield >> tmp.
x >> tmp.
y >> tmp.
z >> tmp.
v;
207 double closestDist = 1;
208 const unsigned int nIdx = nodeIdx[tmp].size();
210 for (
unsigned int k = 0; k < nIdx; ++k) {
211 int j = nodeIdx[tmp][k];
212 double dist = (tmp.
x -
nodes[j].x) * (tmp.
x -
nodes[j].x) +
215 if (dist < closestDist) {
222 std::cerr <<
" Could not match the node from field potentials file: "
223 << tmp.
x <<
" " << tmp.
y <<
" " << tmp.
z <<
"\n.";
252 std::cerr <<
m_className <<
"::SetWeightingField:\n";
253 std::cerr <<
" No valid field map is present.\n";
254 std::cerr <<
" Weighting field cannot be added.\n";
259 std::ifstream ffield;
260 ffield.open(field.c_str(), std::ios::in);
263 std::cerr <<
" Could not open field potentials file " << field
264 <<
" for reading.\n";
280 for (
int j = 0; j <
nNodes; ++j) {
284 std::cout <<
m_className <<
"::SetWeightingField:\n";
285 std::cout <<
" Replacing existing weighting field " << label <<
".\n";
289 std::map<Node, std::vector<int>,
nodeCmp> nodeIdx;
290 for (
int i = 0; i <
nNodes; ++i) {
291 nodeIdx[
nodes[i]].push_back(i);
293 std::cout <<
"Map size: " << nodeIdx.size() << std::endl;
297 std::getline(ffield, line);
301 for (
int i = 0; i <
nNodes; ++i) {
303 ffield >> tmp.
x >> tmp.
y >> tmp.
z >> tmp.
v;
308 double closestDist = 1;
309 const unsigned int nIdx = nodeIdx[tmp].size();
311 for (
unsigned int k = 0; k < nIdx; ++k) {
312 int j = nodeIdx[tmp][k];
313 double dist = (tmp.
x -
nodes[j].x) * (tmp.
x -
nodes[j].x) +
316 if (dist < closestDist) {
323 std::cerr <<
" Could not match the node from field potentials file: "
324 << tmp.
x <<
" " << tmp.
y <<
" " << tmp.
z <<
"\n.";
327 nodes[closest].w[iw] = tmp.
v;
334 const double z,
double& ex,
double& ey,
335 double& ez,
Medium*& m,
int& status) {
341 const double zin,
double& ex,
double& ey,
342 double& ez,
double& volt,
Medium*& m,
345 double x = xin, y = yin, z = zin;
348 bool xmirr, ymirr, zmirr;
349 double rcoordinate, rotation;
350 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
353 ex = ey = ez = volt = 0.;
367 double t1, t2, t3, t4, jac[4][4], det;
368 const int imap =
FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
372 std::cout <<
" Point (" << x <<
", " << y <<
", " << z
373 <<
" not in the mesh.\n";
381 PrintElement(
"ElectricField", x, y, z, t1, t2, t3, t4, element, 10);
394 volt = n0.
v * t1 * (2 * t1 - 1) + n1.
v * t2 * (2 * t2 - 1) +
395 n2.
v * t3 * (2 * t3 - 1) + n3.
v * t4 * (2 * t4 - 1) +
396 4 * n4.
v * t1 * t2 + 4 * n5.
v * t1 * t3 + 4 * n6.
v * t1 * t4 +
397 4 * n7.
v * t2 * t3 + 4 * n8.
v * t2 * t4 + 4 * n9.
v * t3 * t4;
398 ex = -(n0.
v * (4 * t1 - 1) * jac[0][1] + n1.
v * (4 * t2 - 1) * jac[1][1] +
399 n2.
v * (4 * t3 - 1) * jac[2][1] + n3.
v * (4 * t4 - 1) * jac[3][1] +
400 n4.
v * (4 * t2 * jac[0][1] + 4 * t1 * jac[1][1]) +
401 n5.
v * (4 * t3 * jac[0][1] + 4 * t1 * jac[2][1]) +
402 n6.
v * (4 * t4 * jac[0][1] + 4 * t1 * jac[3][1]) +
403 n7.
v * (4 * t3 * jac[1][1] + 4 * t2 * jac[2][1]) +
404 n8.
v * (4 * t4 * jac[1][1] + 4 * t2 * jac[3][1]) +
405 n9.
v * (4 * t4 * jac[2][1] + 4 * t3 * jac[3][1])) /
407 ey = -(n0.
v * (4 * t1 - 1) * jac[0][2] + n1.
v * (4 * t2 - 1) * jac[1][2] +
408 n2.
v * (4 * t3 - 1) * jac[2][2] + n3.
v * (4 * t4 - 1) * jac[3][2] +
409 n4.
v * (4 * t2 * jac[0][2] + 4 * t1 * jac[1][2]) +
410 n5.
v * (4 * t3 * jac[0][2] + 4 * t1 * jac[2][2]) +
411 n6.
v * (4 * t4 * jac[0][2] + 4 * t1 * jac[3][2]) +
412 n7.
v * (4 * t3 * jac[1][2] + 4 * t2 * jac[2][2]) +
413 n8.
v * (4 * t4 * jac[1][2] + 4 * t2 * jac[3][2]) +
414 n9.
v * (4 * t4 * jac[2][2] + 4 * t3 * jac[3][2])) /
416 ez = -(n0.
v * (4 * t1 - 1) * jac[0][3] + n1.
v * (4 * t2 - 1) * jac[1][3] +
417 n2.
v * (4 * t3 - 1) * jac[2][3] + n3.
v * (4 * t4 - 1) * jac[3][3] +
418 n4.
v * (4 * t2 * jac[0][3] + 4 * t1 * jac[1][3]) +
419 n5.
v * (4 * t3 * jac[0][3] + 4 * t1 * jac[2][3]) +
420 n6.
v * (4 * t4 * jac[0][3] + 4 * t1 * jac[3][3]) +
421 n7.
v * (4 * t3 * jac[1][3] + 4 * t2 * jac[2][3]) +
422 n8.
v * (4 * t4 * jac[1][3] + 4 * t2 * jac[3][3]) +
423 n9.
v * (4 * t4 * jac[2][3] + 4 * t3 * jac[3][3])) /
427 UnmapFields(ex, ey, ez, x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
434 std::cout <<
" Material " << element.
matmap <<
", drift flag "
445 const double zin,
double& wx,
double& wy,
446 double& wz,
const std::string& label) {
470 double x = xin, y = yin, z = zin;
473 bool xmirr, ymirr, zmirr;
474 double rcoordinate, rotation;
475 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
480 double t1, t2, t3, t4, jac[4][4], det;
481 const int imap =
FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
483 if (imap < 0)
return;
487 PrintElement(
"WeightingField", x, y, z, t1, t2, t3, t4, element, 10, iw);
500 wx = -(n0.
w[iw] * (4 * t1 - 1) * jac[0][1] +
501 n1.
w[iw] * (4 * t2 - 1) * jac[1][1] +
502 n2.
w[iw] * (4 * t3 - 1) * jac[2][1] +
503 n3.
w[iw] * (4 * t4 - 1) * jac[3][1] +
504 n4.
w[iw] * (4 * t2 * jac[0][1] + 4 * t1 * jac[1][1]) +
505 n5.
w[iw] * (4 * t3 * jac[0][1] + 4 * t1 * jac[2][1]) +
506 n6.
w[iw] * (4 * t4 * jac[0][1] + 4 * t1 * jac[3][1]) +
507 n7.
w[iw] * (4 * t3 * jac[1][1] + 4 * t2 * jac[2][1]) +
508 n8.
w[iw] * (4 * t4 * jac[1][1] + 4 * t2 * jac[3][1]) +
509 n9.
w[iw] * (4 * t4 * jac[2][1] + 4 * t3 * jac[3][1])) /
512 wy = -(n0.
w[iw] * (4 * t1 - 1) * jac[0][2] +
513 n1.
w[iw] * (4 * t2 - 1) * jac[1][2] +
514 n2.
w[iw] * (4 * t3 - 1) * jac[2][2] +
515 n3.
w[iw] * (4 * t4 - 1) * jac[3][2] +
516 n4.
w[iw] * (4 * t2 * jac[0][2] + 4 * t1 * jac[1][2]) +
517 n5.
w[iw] * (4 * t3 * jac[0][2] + 4 * t1 * jac[2][2]) +
518 n6.
w[iw] * (4 * t4 * jac[0][2] + 4 * t1 * jac[3][2]) +
519 n7.
w[iw] * (4 * t3 * jac[1][2] + 4 * t2 * jac[2][2]) +
520 n8.
w[iw] * (4 * t4 * jac[1][2] + 4 * t2 * jac[3][2]) +
521 n9.
w[iw] * (4 * t4 * jac[2][2] + 4 * t3 * jac[3][2])) /
524 wz = -(n0.
w[iw] * (4 * t1 - 1) * jac[0][3] +
525 n1.
w[iw] * (4 * t2 - 1) * jac[1][3] +
526 n2.
w[iw] * (4 * t3 - 1) * jac[2][3] +
527 n3.
w[iw] * (4 * t4 - 1) * jac[3][3] +
528 n4.
w[iw] * (4 * t2 * jac[0][3] + 4 * t1 * jac[1][3]) +
529 n5.
w[iw] * (4 * t3 * jac[0][3] + 4 * t1 * jac[2][3]) +
530 n6.
w[iw] * (4 * t4 * jac[0][3] + 4 * t1 * jac[3][3]) +
531 n7.
w[iw] * (4 * t3 * jac[1][3] + 4 * t2 * jac[2][3]) +
532 n8.
w[iw] * (4 * t4 * jac[1][3] + 4 * t2 * jac[3][3]) +
533 n9.
w[iw] * (4 * t4 * jac[2][3] + 4 * t3 * jac[3][3])) /
537 UnmapFields(wx, wy, wz, x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
542 const std::string& label) {
558 if (!found)
return 0.;
563 double x = xin, y = yin, z = zin;
566 bool xmirr, ymirr, zmirr;
567 double rcoordinate, rotation;
568 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
573 double t1, t2, t3, t4, jac[4][4], det;
574 const int imap =
FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
575 if (imap < 0)
return 0.;
579 PrintElement(
"WeightingPotential", x, y, z, t1, t2, t3, t4, element, 10,
593 return n0.
w[iw] * t1 * (2 * t1 - 1) + n1.
w[iw] * t2 * (2 * t2 - 1) +
594 n2.
w[iw] * t3 * (2 * t3 - 1) + n3.
w[iw] * t4 * (2 * t4 - 1) +
595 4 * n4.
w[iw] * t1 * t2 + 4 * n5.
w[iw] * t1 * t3 +
596 4 * n6.
w[iw] * t1 * t4 + 4 * n7.
w[iw] * t2 * t3 +
597 4 * n8.
w[iw] * t2 * t4 + 4 * n9.
w[iw] * t3 * t4;
603 double x = xin, y = yin, z = zin;
606 bool xmirr, ymirr, zmirr;
607 double rcoordinate, rotation;
608 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
618 double t1, t2, t3, t4, jac[4][4], det;
619 const int imap =
FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
623 std::cout <<
" Point (" << x <<
", " << y <<
", " << z
624 <<
") not in the mesh.\n";
632 std::cerr <<
" Point (" << x <<
", " << y
633 <<
") has out of range material number " << imap <<
".\n";
639 PrintElement(
"GetMedium", x, y, z, t1, t2, t3, t4, element, 10);
646 if (i >=
elements.size())
return 0.;
657 ((n1.
y - n0.
y) * (n2.
z - n0.
z) - (n2.
y - n0.
y) * (n1.
z - n0.
z)) +
659 ((n1.
z - n0.
z) * (n2.
x - n0.
x) - (n2.
z - n0.
z) * (n1.
x - n0.
x)) +
660 (n3.
z - n0.
z) * ((n1.
x - n0.
x) * (n2.
y - n0.
y) -
661 (n3.
x - n0.
x) * (n1.
y - n0.
y))) /
676 for (
int j = 0; j < np - 1; ++j) {
678 for (
int k = j + 1; k < np; ++k) {
681 const double dx = nj.
x - nk.
x;
682 const double dy = nj.
y - nk.
y;
683 const double dz = nj.
z - nk.
z;
684 const double dist = sqrt(dx * dx + dy * dy + dz * dz);
688 if (dist < dmin) dmin = dist;
689 if (dist > dmax) dmax = dist;
bool m_ready
Ready for use?
std::string m_className
Class name.
bool m_debug
Switch on/off debugging messages.
void UpdatePeriodicity() override
Verify periodicities.
Medium * GetMedium(const double x, const double y, const double z) override
Get the medium at a given location (x, y, z).
void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) override
double WeightingPotential(const double x, const double y, const double z, const std::string &label) override
ComponentComsol()
Default constructor.
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
double GetElementVolume(const unsigned int i) override
void GetAspectRatio(const unsigned int i, double &dmin, double &dmax) override
bool SetWeightingField(std::string file, std::string label)
bool Initialise(std::string header="mesh.mphtxt", std::string mplist="dielectrics.dat", std::string field="field.txt")
Base class for components based on finite-element field maps.
std::vector< bool > wfieldsOk
void MapCoordinates(double &xpos, double &ypos, double &zpos, bool &xmirrored, bool &ymirrored, bool &zmirrored, double &rcoordinate, double &rotation) const
Move (xpos, ypos, zpos) to field map coordinates.
void PrintWarning(const std::string &header)
unsigned int m_nMaterials
std::vector< Material > materials
int FindElement13(const double x, const double y, const double z, double &t1, double &t2, double &t3, double &t4, double jac[4][4], double &det)
Find the element for a point in curved quadratic tetrahedra.
std::vector< Node > nodes
void PrintElement(const std::string &header, const double x, const double y, const double z, const double t1, const double t2, const double t3, const double t4, const Element &element, const unsigned int n, const int iw=-1) const
virtual void SetRange()
Calculate x, y, z, V and angular ranges.
std::vector< Element > elements
std::vector< std::string > wfields
void UnmapFields(double &ex, double &ey, double &ez, double &xpos, double &ypos, double &zpos, bool &xmirrored, bool &ymirrored, bool &zmirrored, double &rcoordinate, double &rotation) const
Move (ex, ey, ez) to global coordinates.
void PrintNotReady(const std::string &header) const
Abstract base class for media.
bool IsDriftable() const
Is charge carrier transport enabled in this medium?
int readInt(std::string s)
bool ends_with(std::string s, std::string t)