12void PrintErrorReadingFile(
const std::string& hdr,
const std::string& file,
14 std::cerr << hdr <<
"\n Error reading file " << file <<
" (line " << line
18void PrintErrorOpeningFile(
const std::string& hdr,
const std::string& filetype,
19 const std::string& filename) {
20 std::cerr << hdr <<
"\n Could not open " << filetype <<
" file "
21 << filename <<
" for reading.\n";
22 std::cerr <<
" The file perhaps does not exist.\n";
33 const std::string& elist,
34 const std::string& nlist,
35 const std::string& mplist,
36 const std::string& volt,
const std::string& unit)
39 Initialise(header, elist, nlist, mplist, volt, unit);
43 const std::string& elist,
44 const std::string& nlist,
45 const std::string& mplist,
46 const std::string& volt,
47 const std::string& unit) {
48 const std::string hdr =
m_className +
"::Initialise:";
62 std::ifstream fheader;
63 fheader.open(header.c_str(), std::ios::in);
65 PrintErrorOpeningFile(hdr,
"header", header);
70 bool readerror =
false;
71 bool readstop =
false;
75 fheader.getline(line, size,
'\n');
76 token = strtok(line,
" ");
78 token = strtok(NULL,
" ");
81 <<
" elements from file " << header <<
".\n";
83 PrintErrorReadingFile(hdr, header, il);
93 fnodes.open(nlist.c_str(), std::ios::in);
95 PrintErrorOpeningFile(hdr,
"nodes", nlist);
100 if (strcmp(unit.c_str(),
"mum") == 0 || strcmp(unit.c_str(),
"micron") == 0 ||
101 strcmp(unit.c_str(),
"micrometer") == 0) {
103 }
else if (strcmp(unit.c_str(),
"mm") == 0 ||
104 strcmp(unit.c_str(),
"millimeter") == 0) {
106 }
else if (strcmp(unit.c_str(),
"cm") == 0 ||
107 strcmp(unit.c_str(),
"centimeter") == 0) {
109 }
else if (strcmp(unit.c_str(),
"m") == 0 ||
110 strcmp(unit.c_str(),
"meter") == 0) {
113 std::cerr << hdr <<
" Unknown length unit " << unit <<
".\n";
117 if (
m_debug) std::cout << hdr <<
" Unit scaling factor = " << funit <<
".\n";
120 for (il = 0; il <
nNodes; il++) {
122 fnodes.getline(line, size,
'\n');
125 token = strtok(line,
" ");
126 token = strtok(NULL,
" ");
129 token = strtok(NULL,
" ");
130 double xnode =
ReadDouble(token, -1, readerror);
131 token = strtok(NULL,
" ");
132 double ynode =
ReadDouble(token, -1, readerror);
133 token = strtok(NULL,
" ");
134 double znode =
ReadDouble(token, -1, readerror);
136 PrintErrorReadingFile(hdr, nlist, il);
144 newNode.
x = xnode * funit;
145 newNode.
y = ynode * funit;
146 newNode.
z = znode * funit;
147 nodes.push_back(std::move(newNode));
155 fvolt.open(volt.c_str(), std::ios::in);
157 PrintErrorOpeningFile(hdr,
"potentials", volt);
164 while (!readstop && fvolt.getline(line, size,
'\n')) {
165 token = strtok(line,
" ");
166 if (strcmp(token,
"Perm:") == 0) readstop =
true;
172 std::cerr << hdr <<
"\n Error reading past header of potentials file "
179 for (
int tl = 0; tl <
nNodes; tl++) {
180 fvolt.getline(line, size,
'\n');
185 for (
int tl = 0; tl <
nNodes; tl++) {
186 fvolt.getline(line, size,
'\n');
187 token = strtok(line,
" ");
190 PrintErrorReadingFile(hdr, volt, il);
202 std::ifstream fmplist;
203 fmplist.open(mplist.c_str(), std::ios::in);
204 if (fmplist.fail()) {
205 PrintErrorOpeningFile(hdr,
"materials", mplist);
209 fmplist.getline(line, size,
'\n');
210 token = strtok(line,
" ");
212 std::cerr << hdr <<
"\n Error reading number of materials from "
226 fmplist.getline(line, size,
'\n');
227 token = strtok(line,
" ");
229 token = strtok(NULL,
" ");
230 double dc =
ReadDouble(token, -1.0, readerror);
232 PrintErrorReadingFile(hdr, mplist, il);
238 std::cout << hdr <<
"\n Set material " << il - 2 <<
" of "
247 unsigned int iepsmin = 0;
248 for (
unsigned int imat = 0; imat <
m_nMaterials; ++imat) {
251 std::cerr << hdr <<
"\n Material " << imat
252 <<
" has been assigned a permittivity equal to zero in\n "
255 }
else if (epsmin < 0. || epsmin >
materials[imat].eps) {
262 std::cerr << hdr <<
"\n No material with positive permittivity found \n"
263 <<
" in material list " << mplist <<
".\n";
266 for (
unsigned int imat = 0; imat <
m_nMaterials; ++imat) {
267 materials[imat].driftmedium = imat == iepsmin ? true :
false;
272 std::ifstream felems;
273 felems.open(elist.c_str(), std::ios::in);
275 PrintErrorOpeningFile(hdr,
"elements", elist);
284 felems.getline(line, size,
'\n');
287 token = strtok(line,
" ");
295 token = strtok(NULL,
" ");
297 token = strtok(NULL,
" ");
298 token = strtok(NULL,
" ");
300 token = strtok(NULL,
" ");
302 token = strtok(NULL,
" ");
304 token = strtok(NULL,
" ");
306 token = strtok(NULL,
" ");
308 token = strtok(NULL,
" ");
310 token = strtok(NULL,
" ");
312 token = strtok(NULL,
" ");
314 token = strtok(NULL,
" ");
316 token = strtok(NULL,
" ");
320 std::cout <<
" Read nodes " << in0 <<
", " << in1 <<
", " << in2
321 <<
", " << in3 <<
", ... from element " << il + 1 <<
" of "
322 <<
nElements <<
" with mat " << imat <<
".\n";
327 PrintErrorReadingFile(hdr, elist, il);
335 std::cerr << hdr <<
"\n Out-of-range material number on file " << elist
336 <<
" (line " << il <<
").\n";
337 std::cerr <<
" Element: " << il <<
", material: " << imat <<
"\n";
338 std::cerr <<
" nodes: (" << in0 <<
", " << in1 <<
", " << in2 <<
", "
339 << in3 <<
", " << in4 <<
", " << in5 <<
", " << in6 <<
", "
340 << in7 <<
", " << in8 <<
", " << in9 <<
")\n";
344 std::cerr << hdr <<
"\n Element " << il <<
" in element list " << elist
345 <<
"\n uses material " << imat
346 <<
" which has not been assigned a positive permittivity in "
352 if (in0 < 1 || in1 < 1 || in2 < 1 || in3 < 1 || in4 < 1 || in5 < 1 ||
353 in6 < 1 || in7 < 1 || in8 < 1 || in9 < 1) {
354 std::cerr << hdr <<
"\n Found a node number < 1 on file " << elist
355 <<
" (line " << il <<
").\n Element: " << il
356 <<
", material: " << imat <<
"\n nodes: (" << in0 <<
", "
357 << in1 <<
", " << in2 <<
", " << in3 <<
", " << in4 <<
", "
358 << in5 <<
", " << in6 <<
", " << in7 <<
", " << in8 <<
", "
362 if (in0 > highestnode) highestnode = in0;
363 if (in1 > highestnode) highestnode = in1;
364 if (in2 > highestnode) highestnode = in2;
365 if (in3 > highestnode) highestnode = in3;
366 if (in4 > highestnode) highestnode = in4;
367 if (in5 > highestnode) highestnode = in5;
368 if (in6 > highestnode) highestnode = in6;
369 if (in7 > highestnode) highestnode = in7;
370 if (in8 > highestnode) highestnode = in8;
371 if (in9 > highestnode) highestnode = in9;
374 if (in0 == in1 || in0 == in2 || in0 == in3 || in0 == in4 || in0 == in5 ||
375 in0 == in6 || in0 == in7 || in0 == in8 || in0 == in9 || in1 == in2 ||
376 in1 == in3 || in1 == in4 || in1 == in5 || in1 == in6 || in1 == in7 ||
377 in1 == in8 || in1 == in9 || in2 == in3 || in2 == in4 || in2 == in5 ||
378 in2 == in6 || in2 == in7 || in2 == in8 || in2 == in9 || in3 == in4 ||
379 in3 == in5 || in3 == in6 || in3 == in7 || in3 == in8 || in3 == in9 ||
380 in4 == in5 || in4 == in6 || in4 == in7 || in4 == in8 || in4 == in9 ||
381 in5 == in6 || in5 == in7 || in5 == in8 || in5 == in9 || in6 == in7 ||
382 in6 == in8 || in6 == in9 || in7 == in8 || in7 == in9 || in8 == in9) {
383 std::cerr << hdr <<
"\n Element " << il <<
" of file " << elist
384 <<
" is degenerate,\n"
385 <<
" no such elements are allowed in this type of map.\n";
395 newElement.
emap[0] = in0 - 1;
396 newElement.
emap[1] = in1 - 1;
397 newElement.
emap[2] = in2 - 1;
398 newElement.
emap[3] = in3 - 1;
399 newElement.
emap[4] = in4 - 1;
400 newElement.
emap[7] = in5 - 1;
401 newElement.
emap[5] = in6 - 1;
402 newElement.
emap[6] = in7 - 1;
403 newElement.
emap[8] = in8 - 1;
404 newElement.
emap[9] = in9 - 1;
415 std::cerr << hdr <<
"\n Field map could not be "
416 <<
"read and cannot be interpolated.\n";
420 std::cout << hdr <<
" Finished.\n";
434 const std::string hdr =
m_className +
"::SetWeightingField:";
437 std::cerr <<
" Weighting field cannot be added.\n";
445 std::ifstream fwvolt;
446 fwvolt.open(wvolt.c_str(), std::ios::in);
448 PrintErrorOpeningFile(hdr,
"potential", wvolt);
464 for (
int j =
nNodes; j--;) {
468 std::cout << hdr <<
"\n Replacing existing weighting field " << label
475 const int size = 100;
478 bool readerror =
false;
479 bool readstop =
false;
483 while (!readstop && fwvolt.getline(line, size,
'\n')) {
484 token = strtok(line,
" ");
485 if (strcmp(token,
"Perm:") == 0) readstop =
true;
491 std::cerr << hdr <<
"\n Error reading past header of potentials file "
499 for (
int tl = 0; tl <
nNodes; tl++) {
500 fwvolt.getline(line, size,
'\n');
505 for (
int tl = 0; tl <
nNodes; tl++) {
507 fwvolt.getline(line, size,
'\n');
508 token = strtok(line,
" ");
511 PrintErrorReadingFile(hdr, wvolt, il);
522 std::cout << hdr <<
"\n Read potentials from file " << wvolt <<
".\n";
527 std::cerr << hdr <<
"\n Field map could not "
528 <<
"be read and cannot be interpolated.\n";
535 const double z,
double& ex,
double& ey,
536 double& ez,
Medium*& m,
int& status) {
542 const double zin,
double& ex,
double& ey,
543 double& ez,
double& volt,
Medium*& m,
546 double x = xin, y = yin, z = zin;
549 bool xmirr, ymirr, zmirr;
550 double rcoordinate, rotation;
551 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
554 ex = ey = ez = volt = 0.;
568 double t1, t2, t3, t4, jac[4][4], det;
569 const int imap =
FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
572 std::cout <<
m_className <<
"::ElectricField:\n Point (" << x <<
", "
573 << y <<
", " << z <<
") is not in the mesh.\n";
581 PrintElement(
"ElectricField", x, y, z, t1, t2, t3, t4, element, 10);
594 const double fourt1 = 4 * t1;
595 const double fourt2 = 4 * t2;
596 const double fourt3 = 4 * t3;
597 const double fourt4 = 4 * t4;
598 const double invdet = 1. / det;
600 volt = n0.
v * t1 * (2 * t1 - 1) + n1.
v * t2 * (2 * t2 - 1) +
601 n2.
v * t3 * (2 * t3 - 1) + n3.
v * t4 * (2 * t4 - 1) +
602 4 * n4.
v * t1 * t2 + 4 * n5.
v * t1 * t3 + 4 * n6.
v * t1 * t4 +
603 4 * n7.
v * t2 * t3 + 4 * n8.
v * t2 * t4 + 4 * n9.
v * t3 * t4;
604 ex = -(n0.
v * (fourt1 - 1) * jac[0][1] + n1.
v * (fourt2 - 1) * jac[1][1] +
605 n2.
v * (fourt3 - 1) * jac[2][1] + n3.
v * (fourt4 - 1) * jac[3][1] +
606 n4.
v * (fourt2 * jac[0][1] + fourt1 * jac[1][1]) +
607 n5.
v * (fourt3 * jac[0][1] + fourt1 * jac[2][1]) +
608 n6.
v * (fourt4 * jac[0][1] + fourt1 * jac[3][1]) +
609 n7.
v * (fourt3 * jac[1][1] + fourt2 * jac[2][1]) +
610 n8.
v * (fourt4 * jac[1][1] + fourt2 * jac[3][1]) +
611 n9.
v * (fourt4 * jac[2][1] + fourt3 * jac[3][1])) *
613 ey = -(n0.
v * (fourt1 - 1) * jac[0][2] + n1.
v * (fourt2 - 1) * jac[1][2] +
614 n2.
v * (fourt3 - 1) * jac[2][2] + n3.
v * (fourt4 - 1) * jac[3][2] +
615 n4.
v * (fourt2 * jac[0][2] + fourt1 * jac[1][2]) +
616 n5.
v * (fourt3 * jac[0][2] + fourt1 * jac[2][2]) +
617 n6.
v * (fourt4 * jac[0][2] + fourt1 * jac[3][2]) +
618 n7.
v * (fourt3 * jac[1][2] + fourt2 * jac[2][2]) +
619 n8.
v * (fourt4 * jac[1][2] + fourt2 * jac[3][2]) +
620 n9.
v * (fourt4 * jac[2][2] + fourt3 * jac[3][2])) *
622 ez = -(n0.
v * (fourt1 - 1) * jac[0][3] + n1.
v * (fourt2 - 1) * jac[1][3] +
623 n2.
v * (fourt3 - 1) * jac[2][3] + n3.
v * (fourt4 - 1) * jac[3][3] +
624 n4.
v * (fourt2 * jac[0][3] + fourt1 * jac[1][3]) +
625 n5.
v * (fourt3 * jac[0][3] + fourt1 * jac[2][3]) +
626 n6.
v * (fourt4 * jac[0][3] + fourt1 * jac[3][3]) +
627 n7.
v * (fourt3 * jac[1][3] + fourt2 * jac[2][3]) +
628 n8.
v * (fourt4 * jac[1][3] + fourt2 * jac[3][3]) +
629 n9.
v * (fourt4 * jac[2][3] + fourt3 * jac[3][3])) *
633 UnmapFields(ex, ey, ez, x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
638 std::cout <<
m_className <<
"::ElectricField:\n Material "
647 const double zin,
double& wx,
double& wy,
648 double& wz,
const std::string& label) {
672 double x = xin, y = yin, z = zin;
675 bool xmirr, ymirr, zmirr;
676 double rcoordinate, rotation;
677 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
682 double t1, t2, t3, t4, jac[4][4], det;
683 const int imap =
FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
685 if (imap < 0)
return;
689 PrintElement(
"WeightingField", x, y, z, t1, t2, t3, t4, element, 10, iw);
702 const double fourt1 = 4 * t1;
703 const double fourt2 = 4 * t2;
704 const double fourt3 = 4 * t3;
705 const double fourt4 = 4 * t4;
706 const double invdet = 1. / det;
708 wx = -(n0.
w[iw] * (fourt1 - 1) * jac[0][1] +
709 n1.
w[iw] * (fourt2 - 1) * jac[1][1] +
710 n2.
w[iw] * (fourt3 - 1) * jac[2][1] +
711 n3.
w[iw] * (fourt4 - 1) * jac[3][1] +
712 n4.
w[iw] * (fourt2 * jac[0][1] + fourt1 * jac[1][1]) +
713 n5.
w[iw] * (fourt3 * jac[0][1] + fourt1 * jac[2][1]) +
714 n6.
w[iw] * (fourt4 * jac[0][1] + fourt1 * jac[3][1]) +
715 n7.
w[iw] * (fourt3 * jac[1][1] + fourt2 * jac[2][1]) +
716 n8.
w[iw] * (fourt4 * jac[1][1] + fourt2 * jac[3][1]) +
717 n9.
w[iw] * (fourt4 * jac[2][1] + fourt3 * jac[3][1])) *
719 wy = -(n0.
w[iw] * (fourt1 - 1) * jac[0][2] +
720 n1.
w[iw] * (fourt2 - 1) * jac[1][2] +
721 n2.
w[iw] * (fourt3 - 1) * jac[2][2] +
722 n3.
w[iw] * (fourt4 - 1) * jac[3][2] +
723 n4.
w[iw] * (fourt2 * jac[0][2] + fourt1 * jac[1][2]) +
724 n5.
w[iw] * (fourt3 * jac[0][2] + fourt1 * jac[2][2]) +
725 n6.
w[iw] * (fourt4 * jac[0][2] + fourt1 * jac[3][2]) +
726 n7.
w[iw] * (fourt3 * jac[1][2] + fourt2 * jac[2][2]) +
727 n8.
w[iw] * (fourt4 * jac[1][2] + fourt2 * jac[3][2]) +
728 n9.
w[iw] * (fourt4 * jac[2][2] + fourt3 * jac[3][2])) *
730 wz = -(n0.
w[iw] * (fourt1 - 1) * jac[0][3] +
731 n1.
w[iw] * (fourt2 - 1) * jac[1][3] +
732 n2.
w[iw] * (fourt3 - 1) * jac[2][3] +
733 n3.
w[iw] * (fourt4 - 1) * jac[3][3] +
734 n4.
w[iw] * (fourt2 * jac[0][3] + fourt1 * jac[1][3]) +
735 n5.
w[iw] * (fourt3 * jac[0][3] + fourt1 * jac[2][3]) +
736 n6.
w[iw] * (fourt4 * jac[0][3] + fourt1 * jac[3][3]) +
737 n7.
w[iw] * (fourt3 * jac[1][3] + fourt2 * jac[2][3]) +
738 n8.
w[iw] * (fourt4 * jac[1][3] + fourt2 * jac[3][3]) +
739 n9.
w[iw] * (fourt4 * jac[2][3] + fourt3 * jac[3][3])) *
743 UnmapFields(wx, wy, wz, x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
748 const std::string& label) {
764 if (!found)
return 0.;
769 double x = xin, y = yin, z = zin;
772 bool xmirr, ymirr, zmirr;
773 double rcoordinate, rotation;
774 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
779 double t1, t2, t3, t4, jac[4][4], det;
780 const int imap =
FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
781 if (imap < 0)
return 0.;
785 PrintElement(
"WeightingPotential", x, y, z, t1, t2, t3, t4, element, 10,
799 return n0.
w[iw] * t1 * (2 * t1 - 1) + n1.
w[iw] * t2 * (2 * t2 - 1) +
800 n2.
w[iw] * t3 * (2 * t3 - 1) + n3.
w[iw] * t4 * (2 * t4 - 1) +
801 4 * n4.
w[iw] * t1 * t2 + 4 * n5.
w[iw] * t1 * t3 +
802 4 * n6.
w[iw] * t1 * t4 + 4 * n7.
w[iw] * t2 * t3 +
803 4 * n8.
w[iw] * t2 * t4 + 4 * n9.
w[iw] * t3 * t4;
809 double x = xin, y = yin, z = zin;
812 bool xmirr, ymirr, zmirr;
813 double rcoordinate, rotation;
814 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
824 double t1, t2, t3, t4, jac[4][4], det;
825 const int imap =
FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
828 std::cout <<
m_className <<
"::GetMedium:\n Point (" << x <<
", " << y
829 <<
", " << z <<
") is not in the mesh.\n";
836 std::cerr <<
m_className <<
"::GetMedium:\n Point (" << x <<
", " << y
837 <<
", " << z <<
") has out of range material number " << imap
849 if (i >=
elements.size())
return 0.;
860 ((n1.
y - n0.
y) * (n2.
z - n0.
z) - (n2.
y - n0.
y) * (n1.
z - n0.
z)) +
862 ((n1.
z - n0.
z) * (n2.
x - n0.
x) - (n2.
z - n0.
z) * (n1.
x - n0.
x)) +
863 (n3.
z - n0.
z) * ((n1.
x - n0.
x) * (n2.
y - n0.
y) -
864 (n3.
x - n0.
x) * (n1.
y - n0.
y))) /
879 for (
int j = 0; j < np - 1; ++j) {
881 for (
int k = j + 1; k < np; ++k) {
884 const double dx = nj.
x - nk.
x;
885 const double dy = nj.
y - nk.
y;
886 const double dz = nj.
z - nk.
z;
887 const double dist = sqrt(dx * dx + dy * dy + dz * dz);
891 if (dist < dmin) dmin = dist;
892 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 WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) override
bool Initialise(const std::string &header="mesh.header", const std::string &elist="mesh.elements", const std::string &nlist="mesh.nodes", const std::string &mplist="dielectrics.dat", const std::string &volt="out.result", const std::string &unit="cm")
void UpdatePeriodicity() override
Verify periodicities.
bool SetWeightingField(std::string prnsol, std::string label)
Import a list of voltages to be used as weighting field.
ComponentElmer()
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
double WeightingPotential(const double x, const double y, const double z, const std::string &label) override
Medium * GetMedium(const double x, const double y, const double z) override
Get the medium at a given location (x, y, z).
Base class for components based on finite-element field maps.
double ReadDouble(char *token, double def, bool &error)
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)
int ReadInteger(char *token, int def, bool &error)
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?