Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ComponentComsol.cc
Go to the documentation of this file.
1// Copied and modified ComponentAnsys123.cc
2
3#include <math.h>
4#include <stdlib.h>
5#include <fstream>
6#include <iostream>
7#include <map>
8#include <sstream>
9
11
12namespace Garfield {
13
15 m_className = "ComponentComsol";
16}
17
18ComponentComsol::ComponentComsol(std::string mesh, std::string mplist,
19 std::string field)
21 m_className = "ComponentComsol";
22 Initialise(mesh, mplist, field);
23}
24
25bool ends_with(std::string s, std::string t) {
26 return s.size() >= t.size() && s.substr(s.size() - t.size(), t.size()) == t;
27}
28
29int readInt(std::string s) {
30 std::istringstream iss(s);
31 int ret;
32 iss >> ret;
33 return ret;
34}
35
36bool ComponentComsol::Initialise(std::string mesh, std::string mplist,
37 std::string field) {
38 m_ready = false;
39 m_warning = false;
40 m_nWarnings = 0;
41
42 double unit = 100.0; // m
43
44 std::string line;
45
46 // Open the materials file.
47 materials.clear();
48 std::ifstream fmplist;
49 fmplist.open(mplist.c_str(), std::ios::in);
50 if (fmplist.fail()) {
51 std::cerr << m_className << "::Initialise:\n";
52 std::cerr << " Could not open result file " << mplist
53 << " for reading.\n";
54 return false;
55 }
56 fmplist >> m_nMaterials;
57 for (unsigned int i = 0; i < m_nMaterials; ++i) {
58 Material newMaterial;
59 newMaterial.driftmedium = true;
60 newMaterial.medium = nullptr;
61 newMaterial.ohm = -1;
62 fmplist >> newMaterial.eps;
63 materials.push_back(newMaterial);
64 }
65 {
66 // add default material
67 Material newMaterial;
68 newMaterial.driftmedium = false;
69 newMaterial.medium = nullptr;
70 newMaterial.eps = newMaterial.ohm = -1;
71 materials.push_back(newMaterial);
73 }
74 std::map<int, int> domain2material;
75 int d2msize;
76 fmplist >> d2msize;
77 for (int i = 0; i < d2msize; ++i) {
78 int domain;
79 fmplist >> domain;
80 fmplist >> domain2material[domain];
81 }
82 fmplist.close();
83
84 nodes.clear();
85 std::ifstream fmesh;
86 fmesh.open(mesh.c_str(), std::ios::in);
87 if (fmesh.fail()) {
88 std::cerr << m_className << "::Initialise:\n";
89 std::cerr << " Could not open nodes file " << mesh << " for reading.\n";
90 return false;
91 }
92
93 do {
94 std::getline(fmesh, line);
95 } while (!ends_with(line, "# number of mesh points"));
96 nNodes = readInt(line);
97
98 std::cout << m_className << "::Initialise:\n";
99 std::cout << " Read " << nNodes << " nodes from file " << mesh << ".\n";
100 do {
101 std::getline(fmesh, line);
102 } while (line != "# Mesh point coordinates");
103 double minx = 1e100, miny = 1e100, minz = 1e100, maxx = -1e100, maxy = -1e100,
104 maxz = -1e100;
105 for (int i = 0; i < nNodes; ++i) {
106 Node newNode;
107 fmesh >> newNode.x >> newNode.y >> newNode.z;
108 newNode.x *= unit;
109 newNode.y *= unit;
110 newNode.z *= unit;
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);
118 }
119 std::cout << minx << " < x < " << maxx << "\n";
120 std::cout << miny << " < y < " << maxy << "\n";
121 std::cout << minz << " < z < " << maxz << "\n";
122
123 do {
124 std::getline(fmesh, line);
125 } while (line != "4 tet2 # type name");
126 do {
127 std::getline(fmesh, line);
128 } while (!ends_with(line, "# number of elements"));
129 nElements = readInt(line);
130 elements.clear();
131 std::cout << m_className << "::Initialise:\n";
132 std::cout << " Read " << nElements << " elements from file " << mesh
133 << ".\n";
134 std::getline(fmesh, line);
135 // elements 6 & 7 are swapped due to differences in COMSOL and ANSYS
136 // representation
137 int perm[10] = {0, 1, 2, 3, 4, 5, 7, 6, 8, 9};
138 for (int i = 0; i < nElements; ++i) {
139 Element newElement;
140 newElement.degenerate = false;
141 for (int j = 0; j < 10; ++j) {
142 fmesh >> newElement.emap[perm[j]];
143 }
144 elements.push_back(newElement);
145 }
146
147 do {
148 std::getline(fmesh, line);
149 } while (line != "# Geometric entity indices");
150 for (int i = 0; i < nElements; ++i) {
151 int domain;
152 fmesh >> domain;
153 elements[i].matmap = domain2material.count(domain) ? domain2material[domain]
154 : m_nMaterials - 1;
155 }
156 fmesh.close();
157
158 std::map<Node, std::vector<int>, nodeCmp> nodeIdx;
159 for (int i = 0; i < nNodes; ++i) {
160 nodeIdx[nodes[i]].push_back(i);
161 }
162 std::cout << "Map size: " << nodeIdx.size() << std::endl;
163
164 std::ifstream ffield;
165 ffield.open(field.c_str(), std::ios::in);
166 if (ffield.fail()) {
167 std::cerr << m_className << "::Initialise:\n";
168 std::cerr << " Could not open field potentials file " << field
169 << " for reading.\n";
170 return false;
171 }
172 do {
173 std::getline(ffield, line);
174 } while (line.substr(0, 81) !=
175 "% x y z "
176 " V (V)");
177 {
178 std::istringstream sline(line);
179 std::string token;
180 sline >> token; // %
181 sline >> token; // x
182 sline >> token; // y
183 sline >> token; // z
184 sline >> token; // V
185 sline >> token; // (V)
186 while (sline >> token) {
187 std::cout << m_className << "::Initialise:\n";
188 std::cout << " Reading data for weighting field " << token << ".\n";
190 wfields.push_back(token);
191 wfieldsOk.push_back(true);
192 sline >> token; // (V)
193 }
194 }
195 for (int i = 0; i < nNodes; ++i) {
196 Node tmp;
197 ffield >> tmp.x >> tmp.y >> tmp.z >> tmp.v;
198 tmp.x *= unit;
199 tmp.y *= unit;
200 tmp.z *= unit;
201 for (int j = 0; j < nWeightingFields; ++j) {
202 double w;
203 ffield >> w;
204 tmp.w.push_back(w);
205 }
206 int closest = -1;
207 double closestDist = 1;
208 const unsigned int nIdx = nodeIdx[tmp].size();
209 // for (int j : nodeIdx[tmp]) {
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) +
213 (tmp.y - nodes[j].y) * (tmp.y - nodes[j].y) +
214 (tmp.z - nodes[j].z) * (tmp.z - nodes[j].z);
215 if (dist < closestDist) {
216 closestDist = dist;
217 closest = j;
218 }
219 }
220 if (closest == -1) {
221 std::cerr << m_className << "::Initialise:\n";
222 std::cerr << " Could not match the node from field potentials file: "
223 << tmp.x << " " << tmp.y << " " << tmp.z << "\n.";
224 return false;
225 }
226 nodes[closest].v = tmp.v;
227 nodes[closest].w = tmp.w;
228 }
229
230 m_ready = true;
231
232 // for (int i = 0; i < nNodes; ++i) {
233 // double ex, ey, ez, v;
234 // Medium* m;
235 // int status;
236 // ElectricField(nodes[i].x, nodes[i].y, nodes[i].z, ex, ey, ez, v, m,
237 // status);
238 // std::cout << "Field at " << nodes[i].x << " " << nodes[i].y << " " <<
239 // nodes[i].z << ": " << ex << " " << ey << " " << ez << " " << v << "\n";
240 // }
241
242 // Establish the ranges.
243 SetRange();
245 return true;
246}
247
248bool ComponentComsol::SetWeightingField(std::string field, std::string label) {
249 double unit = 100.0; // m;
250
251 if (!m_ready) {
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";
255 return false;
256 }
257
258 // Open the voltage list.
259 std::ifstream ffield;
260 ffield.open(field.c_str(), std::ios::in);
261 if (ffield.fail()) {
262 std::cerr << m_className << "::Initialise:\n";
263 std::cerr << " Could not open field potentials file " << field
264 << " for reading.\n";
265 return false;
266 }
267
268 // Check if a weighting field with the same label alm_ready exists.
269 int iw = nWeightingFields;
270 for (int i = nWeightingFields; i--;) {
271 if (wfields[i] == label) {
272 iw = i;
273 break;
274 }
275 }
276 if (iw == nWeightingFields) {
280 for (int j = 0; j < nNodes; ++j) {
281 nodes[j].w.resize(nWeightingFields);
282 }
283 } else {
284 std::cout << m_className << "::SetWeightingField:\n";
285 std::cout << " Replacing existing weighting field " << label << ".\n";
286 }
287 wfields[iw] = label;
288 wfieldsOk[iw] = false;
289 std::map<Node, std::vector<int>, nodeCmp> nodeIdx;
290 for (int i = 0; i < nNodes; ++i) {
291 nodeIdx[nodes[i]].push_back(i);
292 }
293 std::cout << "Map size: " << nodeIdx.size() << std::endl;
294
295 std::string line;
296 do {
297 std::getline(ffield, line);
298 } while (line !=
299 "% x y z "
300 " V (V)");
301 for (int i = 0; i < nNodes; ++i) {
302 Node tmp;
303 ffield >> tmp.x >> tmp.y >> tmp.z >> tmp.v;
304 tmp.x *= unit;
305 tmp.y *= unit;
306 tmp.z *= unit;
307 int closest = -1;
308 double closestDist = 1;
309 const unsigned int nIdx = nodeIdx[tmp].size();
310 // for (int j : nodeIdx[tmp]) {
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) +
314 (tmp.y - nodes[j].y) * (tmp.y - nodes[j].y) +
315 (tmp.z - nodes[j].z) * (tmp.z - nodes[j].z);
316 if (dist < closestDist) {
317 closestDist = dist;
318 closest = j;
319 }
320 }
321 if (closest == -1) {
322 std::cerr << m_className << "::Initialise:\n";
323 std::cerr << " Could not match the node from field potentials file: "
324 << tmp.x << " " << tmp.y << " " << tmp.z << "\n.";
325 return false;
326 }
327 nodes[closest].w[iw] = tmp.v;
328 }
329
330 return true;
331}
332
333void ComponentComsol::ElectricField(const double x, const double y,
334 const double z, double& ex, double& ey,
335 double& ez, Medium*& m, int& status) {
336 double v = 0.;
337 ElectricField(x, y, z, ex, ey, ez, v, m, status);
338}
339
340void ComponentComsol::ElectricField(const double xin, const double yin,
341 const double zin, double& ex, double& ey,
342 double& ez, double& volt, Medium*& m,
343 int& status) {
344 // Copy the coordinates
345 double x = xin, y = yin, z = zin;
346
347 // Map the coordinates onto field map coordinates
348 bool xmirr, ymirr, zmirr;
349 double rcoordinate, rotation;
350 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
351
352 // Initial values
353 ex = ey = ez = volt = 0.;
354 status = 0;
355 m = NULL;
356
357 // Do not proceed if not properly initialised.
358 if (!m_ready) {
359 status = -10;
360 PrintNotReady("ElectricField");
361 return;
362 }
363
364 if (m_warning) PrintWarning("ElectricField");
365
366 // Find the element that contains this point
367 double t1, t2, t3, t4, jac[4][4], det;
368 const int imap = FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
369 if (imap < 0) {
370 if (m_debug) {
371 std::cout << m_className << "::ElectricField:\n";
372 std::cout << " Point (" << x << ", " << y << ", " << z
373 << " not in the mesh.\n";
374 }
375 status = -6;
376 return;
377 }
378
379 const Element& element = elements[imap];
380 if (m_debug) {
381 PrintElement("ElectricField", x, y, z, t1, t2, t3, t4, element, 10);
382 }
383 const Node& n0 = nodes[element.emap[0]];
384 const Node& n1 = nodes[element.emap[1]];
385 const Node& n2 = nodes[element.emap[2]];
386 const Node& n3 = nodes[element.emap[3]];
387 const Node& n4 = nodes[element.emap[4]];
388 const Node& n5 = nodes[element.emap[5]];
389 const Node& n6 = nodes[element.emap[6]];
390 const Node& n7 = nodes[element.emap[7]];
391 const Node& n8 = nodes[element.emap[8]];
392 const Node& n9 = nodes[element.emap[9]];
393 // Tetrahedral field
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])) /
406 det;
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])) /
415 det;
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])) /
424 det;
425
426 // Transform field to global coordinates
427 UnmapFields(ex, ey, ez, x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
428 // std::cout << "ef @(" << xin << ", " << yin << ", " << zin << ") = " <<
429 // volt << "\n";
430
431 // Drift medium?
432 if (m_debug) {
433 std::cout << m_className << "::ElectricField:\n";
434 std::cout << " Material " << element.matmap << ", drift flag "
435 << materials[element.matmap].driftmedium << "\n";
436 }
437 m = materials[element.matmap].medium;
438 status = -5;
439 if (materials[element.matmap].driftmedium) {
440 if (m && m->IsDriftable()) status = 0;
441 }
442}
443
444void ComponentComsol::WeightingField(const double xin, const double yin,
445 const double zin, double& wx, double& wy,
446 double& wz, const std::string& label) {
447 // Initial values
448 wx = wy = wz = 0;
449
450 // Do not proceed if not properly initialised.
451 if (!m_ready) return;
452
453 // Look for the label.
454 int iw = 0;
455 bool found = false;
456 for (int i = nWeightingFields; i--;) {
457 if (wfields[i] == label) {
458 iw = i;
459 found = true;
460 break;
461 }
462 }
463
464 // Do not proceed if the requested weighting field does not exist.
465 if (!found) return;
466 // Check if the weighting field is properly initialised.
467 if (!wfieldsOk[iw]) return;
468
469 // Copy the coordinates.
470 double x = xin, y = yin, z = zin;
471
472 // Map the coordinates onto field map coordinates
473 bool xmirr, ymirr, zmirr;
474 double rcoordinate, rotation;
475 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
476
477 if (m_warning) PrintWarning("WeightingField");
478
479 // Find the element that contains this point.
480 double t1, t2, t3, t4, jac[4][4], det;
481 const int imap = FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
482 // Check if the point is in the mesh.
483 if (imap < 0) return;
484
485 const Element& element = elements[imap];
486 if (m_debug) {
487 PrintElement("WeightingField", x, y, z, t1, t2, t3, t4, element, 10, iw);
488 }
489 const Node& n0 = nodes[element.emap[0]];
490 const Node& n1 = nodes[element.emap[1]];
491 const Node& n2 = nodes[element.emap[2]];
492 const Node& n3 = nodes[element.emap[3]];
493 const Node& n4 = nodes[element.emap[4]];
494 const Node& n5 = nodes[element.emap[5]];
495 const Node& n6 = nodes[element.emap[6]];
496 const Node& n7 = nodes[element.emap[7]];
497 const Node& n8 = nodes[element.emap[8]];
498 const Node& n9 = nodes[element.emap[9]];
499 // Tetrahedral field
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])) /
510 det;
511
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])) /
522 det;
523
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])) /
534 det;
535
536 // Transform field to global coordinates
537 UnmapFields(wx, wy, wz, x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
538}
539
540double ComponentComsol::WeightingPotential(const double xin, const double yin,
541 const double zin,
542 const std::string& label) {
543 // Do not proceed if not properly initialised.
544 if (!m_ready) return 0.;
545
546 // Look for the label.
547 int iw = 0;
548 bool found = false;
549 for (int i = nWeightingFields; i--;) {
550 if (wfields[i] == label) {
551 iw = i;
552 found = true;
553 break;
554 }
555 }
556
557 // Do not proceed if the requested weighting field does not exist.
558 if (!found) return 0.;
559 // Check if the weighting field is properly initialised.
560 if (!wfieldsOk[iw]) return 0.;
561
562 // Copy the coordinates.
563 double x = xin, y = yin, z = zin;
564
565 // Map the coordinates onto field map coordinates.
566 bool xmirr, ymirr, zmirr;
567 double rcoordinate, rotation;
568 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
569
570 if (m_warning) PrintWarning("WeightingPotential");
571
572 // Find the element that contains this point.
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.;
576
577 const Element& element = elements[imap];
578 if (m_debug) {
579 PrintElement("WeightingPotential", x, y, z, t1, t2, t3, t4, element, 10,
580 iw);
581 }
582 const Node& n0 = nodes[element.emap[0]];
583 const Node& n1 = nodes[element.emap[1]];
584 const Node& n2 = nodes[element.emap[2]];
585 const Node& n3 = nodes[element.emap[3]];
586 const Node& n4 = nodes[element.emap[4]];
587 const Node& n5 = nodes[element.emap[5]];
588 const Node& n6 = nodes[element.emap[6]];
589 const Node& n7 = nodes[element.emap[7]];
590 const Node& n8 = nodes[element.emap[8]];
591 const Node& n9 = nodes[element.emap[9]];
592 // Tetrahedral field
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;
598}
599
600Medium* ComponentComsol::GetMedium(const double xin, const double yin,
601 const double zin) {
602 // Copy the coordinates
603 double x = xin, y = yin, z = zin;
604
605 // Map the coordinates onto field map coordinates
606 bool xmirr, ymirr, zmirr;
607 double rcoordinate, rotation;
608 MapCoordinates(x, y, z, xmirr, ymirr, zmirr, rcoordinate, rotation);
609
610 // Do not proceed if not properly initialised.
611 if (!m_ready) {
612 PrintNotReady("GetMedium");
613 return nullptr;
614 }
615 if (m_warning) PrintWarning("GetMedium");
616
617 // Find the element that contains this point
618 double t1, t2, t3, t4, jac[4][4], det;
619 const int imap = FindElement13(x, y, z, t1, t2, t3, t4, jac, det);
620 if (imap < 0) {
621 if (m_debug) {
622 std::cout << m_className << "::GetMedium:\n";
623 std::cout << " Point (" << x << ", " << y << ", " << z
624 << ") not in the mesh.\n";
625 }
626 return nullptr;
627 }
628 const Element& element = elements[imap];
629 if (element.matmap >= m_nMaterials) {
630 if (m_debug) {
631 std::cerr << m_className << "::GetMedium:\n";
632 std::cerr << " Point (" << x << ", " << y
633 << ") has out of range material number " << imap << ".\n";
634 }
635 return nullptr;
636 }
637
638 if (m_debug) {
639 PrintElement("GetMedium", x, y, z, t1, t2, t3, t4, element, 10);
640 }
641
642 return materials[element.matmap].medium;
643}
644
645double ComponentComsol::GetElementVolume(const unsigned int i) {
646 if (i >= elements.size()) return 0.;
647 const Element& element = elements[i];
648 const Node& n0 = nodes[element.emap[0]];
649 const Node& n1 = nodes[element.emap[1]];
650 const Node& n2 = nodes[element.emap[2]];
651 const Node& n3 = nodes[element.emap[3]];
652
653 // Uses formula V = |a (dot) b x c|/6
654 // with a => "3", b => "1", c => "2" and origin "0"
655 const double vol =
656 fabs((n3.x - n0.x) *
657 ((n1.y - n0.y) * (n2.z - n0.z) - (n2.y - n0.y) * (n1.z - n0.z)) +
658 (n3.y - n0.y) *
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))) /
662 6.;
663 return vol;
664}
665
666void ComponentComsol::GetAspectRatio(const unsigned int i, double& dmin,
667 double& dmax) {
668 if (i >= elements.size()) {
669 dmin = dmax = 0.;
670 return;
671 }
672
673 const Element& element = elements[i];
674 const int np = 4;
675 // Loop over all pairs of vertices.
676 for (int j = 0; j < np - 1; ++j) {
677 const Node& nj = nodes[element.emap[j]];
678 for (int k = j + 1; k < np; ++k) {
679 const Node& nk = nodes[element.emap[k]];
680 // Compute distance.
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);
685 if (k == 1) {
686 dmin = dmax = dist;
687 } else {
688 if (dist < dmin) dmin = dist;
689 if (dist > dmax) dmax = dist;
690 }
691 }
692 }
693}
694
695} // namespace Garfield
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.
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)
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.
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.
Definition: Medium.hh:13
bool IsDriftable() const
Is charge carrier transport enabled in this medium?
Definition: Medium.hh:74
int readInt(std::string s)
bool ends_with(std::string s, std::string t)