14 : m_className(
"ViewFEMesh"),
19 m_hasExternalCanvas(false),
29 m_plotMeshBorders(false),
40 m_axes->SetStats(
false);
41 m_axes->GetXaxis()->SetTitle(
"x");
42 m_axes->GetYaxis()->SetTitle(
"y");
47 if (!m_hasExternalCanvas && m_canvas)
delete m_canvas;
53 std::cerr << m_className <<
"::SetComponent:\n";
54 std::cerr <<
" Component pointer is null.\n";
64 if (!m_hasExternalCanvas && m_canvas) {
69 m_hasExternalCanvas =
true;
75 double ymax,
double zmax) {
78 if (xmin == xmax || ymin == ymax) {
79 std::cout << m_className <<
"::SetArea:\n";
80 std::cout <<
" Null area range not permitted.\n";
83 m_xMin = std::min(xmin, xmax);
84 m_yMin = std::min(ymin, ymax);
85 m_zMin = std::min(zmin, zmax);
86 m_xMax = std::max(xmin, xmax);
87 m_yMax = std::max(ymin, ymax);
88 m_zMax = std::max(zmin, zmax);
100 std::cerr << m_className <<
"::Plot:\n";
101 std::cerr <<
" Component is not defined.\n";
105 double pmin = 0., pmax = 0.;
107 std::cerr << m_className <<
"::Plot:\n";
108 std::cerr <<
" Component is not ready.\n";
113 if (!m_hasUserArea) {
114 std::cerr << m_className <<
"::Plot:\n";
115 std::cerr <<
" Bounding box cannot be determined.\n";
116 std::cerr <<
" Call SetArea first.\n";
122 m_canvas =
new TCanvas();
123 m_canvas->SetTitle(m_label.c_str());
124 if (m_hasExternalCanvas) m_hasExternalCanvas =
false;
126 m_canvas->Range(m_xMin, m_yMin, m_xMax, m_yMax);
131 std::cout << m_className <<
"::Plot:\n";
132 std::cout <<
" The given component is a CST component.\n";
133 std::cout <<
" Method PlotCST is called now!.\n";
134 DrawCST(componentCST);
146 const double x0,
const double y0,
const double z0) {
149 double fnorm = sqrt(fx * fx + fy * fy + fz * fz);
150 double dist = fx * x0 + fy * y0 + fz * z0;
151 if (fnorm > 0 && fx * fx + fz * fz > 0) {
152 project[0][0] = fz / sqrt(fx * fx + fz * fz);
154 project[0][2] = -fx / sqrt(fx * fx + fz * fz);
155 project[1][0] = -fx * fy / (sqrt(fx * fx + fz * fz) * fnorm);
156 project[1][1] = (fx * fx + fz * fz) / (sqrt(fx * fx + fz * fz) * fnorm);
157 project[1][2] = -fy * fz / (sqrt(fx * fx + fz * fz) * fnorm);
158 project[2][0] = dist * fx / (fnorm * fnorm);
159 project[2][1] = dist * fy / (fnorm * fnorm);
160 project[2][2] = dist * fz / (fnorm * fnorm);
161 }
else if (fnorm > 0 && fy * fy + fz * fz > 0) {
162 project[0][0] = (fy * fy + fz * fz) / (sqrt(fy * fy + fz * fz) * fnorm);
163 project[0][1] = -fx * fz / (sqrt(fy * fy + fz * fz) * fnorm);
164 project[0][2] = -fy * fz / (sqrt(fy * fy + fz * fz) * fnorm);
166 project[1][1] = fz / sqrt(fy * fy + fz * fz);
167 project[1][2] = -fy / sqrt(fy * fy + fz * fz);
168 project[2][0] = dist * fx / (fnorm * fnorm);
169 project[2][1] = dist * fy / (fnorm * fnorm);
170 project[2][2] = dist * fz / (fnorm * fnorm);
172 std::cout << m_className <<
"::SetPlane:\n";
173 std::cout <<
" Normal vector has zero norm.\n";
174 std::cout <<
" No new projection set.\n";
213 m_axes->GetXaxis()->SetTitle(xtitle);
218 m_axes->GetYaxis()->SetTitle(ytitle);
225 m_xaxis =
new TGaxis(m_xMin + std::abs(m_xMax - m_xMin) * 0.1,
226 m_yMin + std::abs(m_yMax - m_yMin) * 0.1,
227 m_xMax - std::abs(m_xMax - m_xMin) * 0.1,
228 m_yMin + std::abs(m_yMax - m_yMin) * 0.1,
229 m_xMin + std::abs(m_xMax - m_xMin) * 0.1,
230 m_xMax - std::abs(m_xMax - m_xMin) * 0.1, 2405,
"x");
231 m_yaxis =
new TGaxis(m_xMin + std::abs(m_xMax - m_xMin) * 0.1,
232 m_yMin + std::abs(m_yMax - m_yMin) * 0.1,
233 m_xMin + std::abs(m_xMax - m_xMin) * 0.1,
234 m_yMax - std::abs(m_yMax - m_yMin) * 0.1,
235 m_yMin + std::abs(m_yMax - m_yMin) * 0.1,
236 m_yMax - std::abs(m_yMax - m_yMin) * 0.1, 2405,
"y");
239 m_xaxis->SetLabelSize(0.025);
240 m_yaxis->SetLabelSize(0.025);
243 m_xaxis->SetTitleSize(0.03);
244 m_xaxis->SetTitle(
"x [cm]");
245 m_yaxis->SetTitleSize(0.03);
246 m_yaxis->SetTitle(
"y [cm]");
251void ViewFEMesh::DrawElements() {
254 double mapxmax = m_component->
mapxmax;
255 double mapxmin = m_component->
mapxmin;
256 double mapymax = m_component->
mapymax;
257 double mapymin = m_component->
mapymin;
258 double mapzmax = m_component->
mapzmax;
259 double mapzmin = m_component->
mapzmin;
262 double sx = mapxmax - mapxmin;
263 double sy = mapymax - mapymin;
264 double sz = mapzmax - mapzmin;
271 m_driftLines.clear();
276 sqrt(plane[0] * plane[0] + plane[1] * plane[1] + plane[2] * plane[2]);
278 dataProj[0] = project[0][0];
279 dataProj[1] = project[1][0];
280 dataProj[2] = plane[0] / fnorm;
281 dataProj[3] = project[0][1];
282 dataProj[4] = project[1][1];
283 dataProj[5] = plane[1] / fnorm;
284 dataProj[6] = project[0][2];
285 dataProj[7] = project[1][2];
286 dataProj[8] = plane[2] / fnorm;
287 TMatrixD projMat(3, 3, dataProj.GetArray());
292 (projMat(1, 1) * projMat(2, 2) - projMat(1, 2) * projMat(2, 1)) -
294 (projMat(1, 0) * projMat(2, 2) - projMat(1, 2) * projMat(2, 0)) +
296 (projMat(1, 0) * projMat(2, 1) - projMat(1, 1) * projMat(2, 0));
303 std::cerr << m_className <<
"::DrawElements:\n";
304 std::cerr <<
" Projection matrix is not invertible.\n";
305 std::cerr <<
" Finite element mesh will not be drawn.\n";
309 double fx = plane[0];
310 double fy = plane[1];
311 double fz = plane[2];
312 double dist = plane[3];
318 int nMaxX = 0, nMinX = 0;
319 int nMaxY = 0, nMinY = 0;
320 int nMaxZ = 0, nMinZ = 0;
322 nMinX = int(m_xMin / sx) - 1;
323 nMaxX = int(m_xMax / sx) + 1;
326 nMinY = int(m_yMin / sy) - 1;
327 nMaxY = int(m_yMax / sy) + 1;
330 nMinZ = int(m_zMin / sz) - 1;
331 nMaxZ = int(m_zMax / sz) + 1;
335 for (
int elem = 0; elem < m_component->
nElements; elem++) {
340 !(m_plotMeshBorders))
343 if (m_disabledMaterial[m_component->
elements[elem].matmap])
continue;
347 double vx1, vy1, vz1;
348 double vx2, vy2, vz2;
349 double vx3, vy3, vz3;
350 double vx4, vy4, vz4;
353 int colorID = m_colorMap.count(m_component->
elements[elem].matmap);
355 colorID = m_colorMap[m_component->
elements[elem].matmap];
361 m_colorMap_fill.count(m_component->
elements[elem].matmap);
362 if (colorID_fill != 0)
363 colorID_fill = m_colorMap_fill[m_component->
elements[elem].matmap];
365 colorID_fill = colorID;
368 for (
int nx = nMinX; nx <= nMaxX; nx++) {
390 m_component->
nodes[m_component->
elements[elem].emap[0]].x + sx * nx;
392 m_component->
nodes[m_component->
elements[elem].emap[1]].x + sx * nx;
394 m_component->
nodes[m_component->
elements[elem].emap[2]].x + sx * nx;
396 m_component->
nodes[m_component->
elements[elem].emap[3]].x + sx * nx;
400 for (
int ny = nMinY; ny <= nMaxY; ny++) {
421 vy1 = m_component->
nodes[m_component->
elements[elem].emap[0]].y +
423 vy2 = m_component->
nodes[m_component->
elements[elem].emap[1]].y +
425 vy3 = m_component->
nodes[m_component->
elements[elem].emap[2]].y +
427 vy4 = m_component->
nodes[m_component->
elements[elem].emap[3]].y +
432 for (
int nz = nMinZ; nz <= nMaxZ; nz++) {
453 vz1 = m_component->
nodes[m_component->
elements[elem].emap[0]].z +
455 vz2 = m_component->
nodes[m_component->
elements[elem].emap[1]].z +
457 vz3 = m_component->
nodes[m_component->
elements[elem].emap[2]].z +
459 vz4 = m_component->
nodes[m_component->
elements[elem].emap[3]].z +
464 std::vector<double> vX;
465 std::vector<double> vY;
469 bool in1 =
false, in2 =
false, in3 =
false, in4 =
false;
472 bool planeCut =
false;
475 double pcf = std::max(
476 std::max(std::max(std::max(std::max(std::max(std::abs(vx1),
485 if (std::abs(fx * vx1 + fy * vy1 + fz * vz1 - dist) < 1.0e-4 * pcf) {
488 if (std::abs(fx * vx2 + fy * vy2 + fz * vz2 - dist) < 1.0e-4 * pcf) {
491 if (std::abs(fx * vx3 + fy * vy3 + fz * vz3 - dist) < 1.0e-4 * pcf) {
494 if (std::abs(fx * vx4 + fy * vy4 + fz * vz4 - dist) < 1.0e-4 * pcf) {
501 PlaneCoords(vx1, vy1, vz1, projMat, xMat);
502 vX.push_back(xMat(0, 0));
503 vY.push_back(xMat(1, 0));
506 PlaneCoords(vx2, vy2, vz2, projMat, xMat);
507 vX.push_back(xMat(0, 0));
508 vY.push_back(xMat(1, 0));
511 PlaneCoords(vx3, vy3, vz3, projMat, xMat);
512 vX.push_back(xMat(0, 0));
513 vY.push_back(xMat(1, 0));
516 PlaneCoords(vx4, vy4, vz4, projMat, xMat);
517 vX.push_back(xMat(0, 0));
518 vY.push_back(xMat(1, 0));
523 planeCut = PlaneCut(vx1, vy1, vz1, vx2, vy2, vz2, xMat);
525 vX.push_back(xMat(0, 0));
526 vY.push_back(xMat(1, 0));
530 planeCut = PlaneCut(vx1, vy1, vz1, vx3, vy3, vz3, xMat);
532 vX.push_back(xMat(0, 0));
533 vY.push_back(xMat(1, 0));
537 planeCut = PlaneCut(vx1, vy1, vz1, vx4, vy4, vz4, xMat);
539 vX.push_back(xMat(0, 0));
540 vY.push_back(xMat(1, 0));
544 planeCut = PlaneCut(vx2, vy2, vz2, vx3, vy3, vz3, xMat);
546 vX.push_back(xMat(0, 0));
547 vY.push_back(xMat(1, 0));
551 planeCut = PlaneCut(vx2, vy2, vz2, vx4, vy4, vz4, xMat);
553 vX.push_back(xMat(0, 0));
554 vY.push_back(xMat(1, 0));
558 planeCut = PlaneCut(vx3, vy3, vz3, vx4, vy4, vz4, xMat);
560 vX.push_back(xMat(0, 0));
561 vY.push_back(xMat(1, 0));
566 if (vX.size() >= 3) {
570 RemoveCrossings(vX, vY);
573 std::vector<double> cX;
574 std::vector<double> cY;
577 ClipToView(vX, vY, cX, cY);
583 RemoveCrossings(cX, cY);
587 TPolyLine* poly =
new TPolyLine();
588 poly->SetLineColor(colorID);
589 poly->SetFillColor(colorID_fill);
590 poly->SetLineWidth(3);
593 for (
int pt = 0; pt < (int)cX.size(); pt++) {
596 poly->SetPoint(polyPts, cX[pt], cY[pt]);
601 m_mesh.push_back(poly);
611 const int nDlines = m_viewDrift->m_driftLines.size();
612 for (
int dline = 0; dline < nDlines; dline++) {
614 const unsigned int npts = m_viewDrift->m_driftLines[dline].vect.size();
616 TPolyLine* poly =
new TPolyLine();
617 poly->SetLineColor(m_viewDrift->m_driftLines[dline].n);
619 for (
unsigned int pt = 0; pt < npts; pt++) {
621 double ptx = m_viewDrift->m_driftLines[dline].vect[pt].x;
622 double pty = m_viewDrift->m_driftLines[dline].vect[pt].y;
623 double ptz = m_viewDrift->m_driftLines[dline].vect[pt].z;
625 PlaneCoords(ptx, pty, ptz, projMat, xMat);
627 if (xMat(0, 0) >= m_xMin && xMat(0, 0) <= m_xMax &&
628 xMat(1, 0) >= m_yMin && xMat(1, 0) <= m_yMax) {
629 poly->SetPoint(polyPts, xMat(0, 0), xMat(1, 0));
635 m_driftLines.push_back(poly);
644 if (!m_xaxis && !m_yaxis && m_drawAxes) {
645 m_axes->GetXaxis()->SetLimits(m_xMin, m_xMax);
646 m_axes->GetYaxis()->SetLimits(m_yMin, m_yMax);
651 if (m_xaxis && m_drawAxes) m_xaxis->Draw();
652 if (m_yaxis && m_drawAxes) m_yaxis->Draw();
655 for (
int m = m_mesh.size(); m--;) {
656 if (m_plotMeshBorders || !m_fillMesh) m_mesh[m]->Draw(
"same");
657 if (m_fillMesh) m_mesh[m]->Draw(
"f:same");
661 for (
int m = m_driftLines.size(); m--;) {
662 m_driftLines[m]->Draw(
"same");
666void ViewFEMesh::DrawCST(ComponentCST* componentCST) {
676 double mapxmax = m_component->
mapxmax;
677 double mapxmin = m_component->
mapxmin;
678 double mapymax = m_component->
mapymax;
679 double mapymin = m_component->
mapymin;
680 double mapzmax = m_component->
mapzmax;
681 double mapzmin = m_component->
mapzmin;
684 double sx = mapxmax - mapxmin;
685 double sy = mapymax - mapymin;
686 double sz = mapzmax - mapzmin;
693 m_driftLines.clear();
698 sqrt(plane[0] * plane[0] + plane[1] * plane[1] + plane[2] * plane[2]);
700 dataProj[0] = project[0][0];
701 dataProj[1] = project[1][0];
702 dataProj[2] = plane[0] / fnorm;
703 dataProj[3] = project[0][1];
704 dataProj[4] = project[1][1];
705 dataProj[5] = plane[1] / fnorm;
706 dataProj[6] = project[0][2];
707 dataProj[7] = project[1][2];
708 dataProj[8] = plane[2] / fnorm;
709 TMatrixD projMat(3, 3, dataProj.GetArray());
714 (projMat(1, 1) * projMat(2, 2) - projMat(1, 2) * projMat(2, 1)) -
716 (projMat(1, 0) * projMat(2, 2) - projMat(1, 2) * projMat(2, 0)) +
718 (projMat(1, 0) * projMat(2, 1) - projMat(1, 1) * projMat(2, 0));
725 std::cerr << m_className <<
"::DrawCST:\n";
726 std::cerr <<
" Projection matrix is not invertible.\n";
727 std::cerr <<
" Finite element mesh will not be drawn.\n";
734 int nMaxX = 0, nMinX = 0;
735 int nMaxY = 0, nMinY = 0;
736 int nMaxZ = 0, nMinZ = 0;
738 nMinX = int(m_xMin / sx) - 1;
739 nMaxX = int(m_xMax / sx) + 1;
742 nMinY = int(m_yMin / sy) - 1;
743 nMaxY = int(m_yMax / sy) + 1;
746 nMinZ = int(m_zMin / sz) - 1;
747 nMaxZ = int(m_zMax / sz) + 1;
751 std::vector<PolygonInfo> elements;
752 int nMinU = 0, nMaxU = 0, nMinV = 0, nMaxV = 0;
753 double mapumin = 0., mapumax = 0., mapvmin = 0., mapvmax = 0.;
754 double su = 0., sv = 0.;
755 bool mirroru =
false, mirrorv =
false;
756 double uMin, vMin, uMax, vMax;
757 unsigned int n_x, n_y, n_z;
758 componentCST->GetNumberOfMeshLines(n_x, n_y, n_z);
759 double e_xmin, e_xmax, e_ymin, e_ymax, e_zmin, e_zmax;
761 if (plane[0] == 0 && plane[1] == 0 && plane[2] == 1) {
762 std::cout << m_className <<
"::DrawCST:\n";
763 std::cout <<
" Creating x-y mesh view.\n";
767 unsigned int i, j, z;
768 if (!componentCST->Coordinate2Index(0, 0, project[2][2], i, j, z)) {
769 std::cerr <<
"Could determine the position of the plane in z direction."
773 std::cout <<
" The plane position in z direction is: " << project[2][2]
792 for (
unsigned int y = 0; y < (n_y - 1); y++) {
793 for (
unsigned int x = 0; x < (n_x - 1); x++) {
794 elem = componentCST->Index2Element(x, y, z);
795 componentCST->GetElementBoundaries(elem, e_xmin, e_xmax, e_ymin, e_ymax,
797 PolygonInfo tmp_info;
798 tmp_info.element = elem;
799 tmp_info.p1[0] = e_xmin;
800 tmp_info.p2[0] = e_xmax;
801 tmp_info.p3[0] = e_xmax;
802 tmp_info.p4[0] = e_xmin;
803 tmp_info.p1[1] = e_ymin;
804 tmp_info.p2[1] = e_ymin;
805 tmp_info.p3[1] = e_ymax;
806 tmp_info.p4[1] = e_ymax;
807 tmp_info.material = componentCST->GetElementMaterial(elem);
808 elements.push_back(tmp_info);
812 }
else if (plane[0] == 0 && plane[1] == -1 && plane[2] == 0) {
813 std::cout << m_className <<
"::DrawCST:\n";
814 std::cout <<
" Creating x-z mesh view.\n";
818 unsigned int i = 0, j = 0, y = 0;
819 if (!componentCST->Coordinate2Index(0, project[2][1], 0, i, y, j)) {
820 std::cerr <<
"Could determine the position of the plane in y direction."
824 std::cout <<
" The plane position in y direction is: " << project[2][1]
844 for (
unsigned int z = 0; z < (n_z - 1); z++) {
845 for (
unsigned int x = 0; x < (n_x - 1); x++) {
846 elem = componentCST->Index2Element(x, y, z);
847 componentCST->GetElementBoundaries(elem, e_xmin, e_xmax, e_ymin, e_ymax,
849 PolygonInfo tmp_info;
850 tmp_info.element = elem;
851 tmp_info.p1[0] = e_xmin;
852 tmp_info.p2[0] = e_xmax;
853 tmp_info.p3[0] = e_xmax;
854 tmp_info.p4[0] = e_xmin;
855 tmp_info.p1[1] = e_zmin;
856 tmp_info.p2[1] = e_zmin;
857 tmp_info.p3[1] = e_zmax;
858 tmp_info.p4[1] = e_zmax;
859 tmp_info.material = componentCST->GetElementMaterial(elem);
860 elements.push_back(tmp_info);
865 }
else if (plane[0] == -1 && plane[1] == 0 && plane[2] == 0) {
866 std::cout << m_className <<
"::DrawCST:\n";
867 std::cout <<
" Creating z-y mesh view.\n";
871 unsigned int i, j, x;
872 if (!componentCST->Coordinate2Index(project[2][0], 0, 0, x, i, j)) {
873 std::cerr <<
"Could determine the position of the plane in x direction."
877 std::cout <<
" The plane position in x direction is: " << project[2][0]
896 for (
unsigned int z = 0; z < (n_z - 1); z++) {
897 for (
unsigned int y = 0; y < (n_y - 1); y++) {
898 elem = componentCST->Index2Element(x, y, z);
899 componentCST->GetElementBoundaries(elem, e_xmin, e_xmax, e_ymin, e_ymax,
901 PolygonInfo tmp_info;
902 tmp_info.element = elem;
903 tmp_info.p1[0] = e_zmin;
904 tmp_info.p2[0] = e_zmax;
905 tmp_info.p3[0] = e_zmax;
906 tmp_info.p4[0] = e_zmin;
907 tmp_info.p1[1] = e_ymin;
908 tmp_info.p2[1] = e_ymin;
909 tmp_info.p3[1] = e_ymax;
910 tmp_info.p4[1] = e_ymax;
911 tmp_info.material = componentCST->GetElementMaterial(elem);
913 elements.push_back(tmp_info);
917 std::cerr << m_className <<
"::DrawCST:\n";
918 std::cerr <<
" The given plane name is not known.\n";
919 std::cerr <<
" Please choose one of the following: xy, xz, yz.\n";
922 std::cout << m_className <<
"::PlotCST:\n";
923 std::cout <<
" Number of elements in the projection of the unit cell:"
924 << elements.size() << std::endl;
925 std::vector<PolygonInfo>::iterator it;
926 std::vector<PolygonInfo>::iterator itend = elements.end();
928 for (
int nu = nMinU; nu <= nMaxU; nu++) {
929 for (
int nv = nMinV; nv <= nMaxV; nv++) {
930 it = elements.begin();
931 while (it != itend) {
932 if (m_disabledMaterial[(*it).material]) {
937 int colorID = m_colorMap.count((*it).material);
939 colorID = m_colorMap[(*it).material];
944 int colorID_fill = m_colorMap_fill.count((*it).material);
945 if (colorID_fill != 0)
946 colorID_fill = m_colorMap_fill[(*it).material];
948 colorID_fill = colorID;
950 TPolyLine* poly =
new TPolyLine();
951 poly->SetLineColor(colorID);
952 poly->SetFillColor(colorID_fill);
953 if (m_plotMeshBorders)
954 poly->SetLineWidth(3);
956 poly->SetLineWidth(1);
958 Double_t tmp_u[4], tmp_v[4];
959 if (mirroru && nu != 2 * (nu / 2)) {
961 tmp_u[0] = mapumin + (mapumax - (*it).p1[0]) + su * nu;
962 tmp_u[1] = mapumin + (mapumax - (*it).p2[0]) + su * nu;
963 tmp_u[2] = mapumin + (mapumax - (*it).p3[0]) + su * nu;
964 tmp_u[3] = mapumin + (mapumax - (*it).p4[0]) + su * nu;
967 tmp_u[0] = (*it).p1[0] + su * nu;
968 tmp_u[1] = (*it).p2[0] + su * nu;
969 tmp_u[2] = (*it).p3[0] + su * nu;
970 tmp_u[3] = (*it).p4[0] + su * nu;
972 if (mirrorv && nv != 2 * (nv / 2)) {
973 tmp_v[0] = mapvmin + (mapvmax - (*it).p1[1]) + sv * nv;
974 tmp_v[1] = mapvmin + (mapvmax - (*it).p2[1]) + sv * nv;
975 tmp_v[2] = mapvmin + (mapvmax - (*it).p3[1]) + sv * nv;
976 tmp_v[3] = mapvmin + (mapvmax - (*it).p4[1]) + sv * nv;
978 tmp_v[0] = (*it).p1[1] + sv * nv;
979 tmp_v[1] = (*it).p2[1] + sv * nv;
980 tmp_v[2] = (*it).p3[1] + sv * nv;
981 tmp_v[3] = (*it).p4[1] + sv * nv;
983 if (tmp_u[0] < uMin || tmp_u[1] > uMax || tmp_v[0] < vMin ||
988 poly->SetPoint(0, tmp_u[0], tmp_v[0]);
989 poly->SetPoint(1, tmp_u[1], tmp_v[1]);
990 poly->SetPoint(2, tmp_u[2], tmp_v[2]);
991 poly->SetPoint(3, tmp_u[3], tmp_v[3]);
993 m_mesh.push_back(poly);
998 std::cout << m_className <<
"::PlotCST:\n";
999 std::cout <<
" Number of polygons to be drawn:" << m_mesh.size()
1003 const int nDlines = m_viewDrift->m_driftLines.size();
1004 for (
int dline = 0; dline < nDlines; dline++) {
1006 const unsigned int npts = m_viewDrift->m_driftLines[dline].vect.size();
1008 TPolyLine* poly =
new TPolyLine();
1009 poly->SetLineColor(m_viewDrift->m_driftLines[dline].n);
1011 for (
unsigned int pt = 0; pt < npts; pt++) {
1013 double ptx = m_viewDrift->m_driftLines[dline].vect[pt].x;
1014 double pty = m_viewDrift->m_driftLines[dline].vect[pt].y;
1015 double ptz = m_viewDrift->m_driftLines[dline].vect[pt].z;
1017 PlaneCoords(ptx, pty, ptz, projMat, xMat);
1019 if (xMat(0, 0) >= uMin && xMat(0, 0) <= uMax && xMat(1, 0) >= vMin &&
1020 xMat(1, 0) <= vMax) {
1021 poly->SetPoint(polyPts, xMat(0, 0), xMat(1, 0));
1026 m_driftLines.push_back(poly);
1033 if (!m_xaxis && !m_yaxis && m_drawAxes) {
1034 m_axes->GetXaxis()->SetLimits(uMin, uMax);
1035 m_axes->GetYaxis()->SetLimits(vMin, vMax);
1039 if (m_xaxis && m_drawAxes) m_xaxis->Draw(
"");
1040 if (m_yaxis && m_drawAxes) m_yaxis->Draw(
"");
1042 for (
int m = m_mesh.size(); m--;) {
1043 if (m_plotMeshBorders || !m_fillMesh) m_mesh[m]->Draw(
"same");
1044 if (m_fillMesh) m_mesh[m]->Draw(
"f:sames");
1048 for (
int m = m_driftLines.size(); m--;) {
1049 m_driftLines[m]->Draw(
"sames");
1055bool ViewFEMesh::InView(
double x,
double y) {
1056 return (x >= m_xMin && x <= m_xMax && y >= m_yMin && y <= m_yMax);
1065void ViewFEMesh::RemoveCrossings(std::vector<double>& x,
1066 std::vector<double>& y) {
1069 double xmin = x[0], xmax = x[0];
1070 double ymin = y[0], ymax = y[0];
1071 for (
int i = 1; i < (int)x.size(); i++) {
1072 if (x[i] < xmin) xmin = x[i];
1073 if (x[i] > xmax) xmax = x[i];
1074 if (y[i] < ymin) ymin = y[i];
1075 if (y[i] > ymax) ymax = y[i];
1079 double xtol = 1e-10 * std::abs(xmax - xmin);
1080 double ytol = 1e-10 * std::abs(ymax - ymin);
1081 for (
int i = 0; i < (int)x.size(); i++) {
1082 for (
int j = i + 1; j < (int)x.size(); j++) {
1083 if (std::abs(x[i] - x[j]) < xtol && std::abs(y[i] - y[j]) < ytol) {
1084 x.erase(x.begin() + j);
1085 y.erase(y.begin() + j);
1092 if (x.size() <= 3)
return;
1102 bool crossings =
true;
1103 while (crossings && (attempts < NN)) {
1108 for (
int i = 1; i <= NN; i++) {
1109 for (
int j = i + 2; j <= NN; j++) {
1112 if ((j + 1) > NN && 1 + (j % NN) >= i)
break;
1118 double xc = 0., yc = 0.;
1119 if (LinesCrossed(x[(i - 1) % NN], y[(i - 1) % NN], x[i % NN],
1120 y[i % NN], x[(j - 1) % NN], y[(j - 1) % NN],
1121 x[j % NN], y[j % NN], xc, yc)) {
1125 for (
int k = 1; k <= (j - i) / 2; k++) {
1127 double xs = x[(i + k - 1) % NN];
1128 double ys = y[(i + k - 1) % NN];
1129 x[(i + k - 1) % NN] = x[(j - k) % NN];
1130 y[(i + k - 1) % NN] = y[(j - k) % NN];
1131 x[(j - k) % NN] = xs;
1132 y[(j - k) % NN] = ys;
1147 if (attempts > NN) {
1148 std::cerr << m_className <<
"::RemoveCrossings:\n";
1150 <<
" WARNING: Maximum attempts reached - crossings not removed.\n";
1161bool ViewFEMesh::LinesCrossed(
double x1,
double y1,
double x2,
double y2,
1162 double u1,
double v1,
double u2,
double v2,
1163 double& xc,
double& yc) {
1168 std::max(std::abs(x1),
1169 std::max(std::abs(x2), std::max(std::abs(u1), std::abs(u2))));
1172 std::max(std::abs(y1),
1173 std::max(std::abs(y2), std::max(std::abs(v1), std::abs(v2))));
1174 if (xtol <= 0) xtol = 1.0e-10;
1175 if (ytol <= 0) ytol = 1.0e-10;
1178 double dy = y2 - y1;
1179 double dv = v2 - v1;
1180 double dx = x1 - x2;
1181 double du = u1 - u2;
1182 double det = dy * du - dx * dv;
1185 if (OnLine(x1, y1, x2, y2, u1, v1)) {
1189 }
else if (OnLine(x1, y1, x2, y2, u2, v2)) {
1193 }
else if (OnLine(u1, v1, u2, v2, x1, y1)) {
1197 }
else if (OnLine(u1, v1, u2, v2, x2, y2)) {
1203 else if (std::abs(det) < xtol * ytol)
1209 xc = (du * (x1 * y2 - x2 * y1) - dx * (u1 * v2 - u2 * v1)) / det;
1210 yc = ((-1 * dv) * (x1 * y2 - x2 * y1) + dy * (u1 * v2 - u2 * v1)) / det;
1213 if (OnLine(x1, y1, x2, y2, xc, yc) && OnLine(u1, v1, u2, v2, xc, yc))
1227bool ViewFEMesh::OnLine(
double x1,
double y1,
double x2,
double y2,
double u,
1232 1.0e-10 * std::max(std::abs(x1), std::max(std::abs(x2), std::abs(u)));
1234 1.0e-10 * std::max(std::abs(y1), std::max(std::abs(y2), std::abs(v)));
1235 if (xtol <= 0) xtol = 1.0e-10;
1236 if (ytol <= 0) ytol = 1.0e-10;
1239 double xc = 0, yc = 0;
1242 if ((std::abs(x1 - u) <= xtol && std::abs(y1 - v) <= ytol) ||
1243 (std::abs(x2 - u) <= xtol && std::abs(y2 - v) <= ytol)) {
1247 else if (std::abs(x1 - x2) <= xtol && std::abs(y1 - y2) <= ytol) {
1251 else if (std::abs(u - x1) + std::abs(v - y1) <
1252 std::abs(u - x2) + std::abs(v - y2)) {
1256 double dpar = ((u - x1) * (x2 - x1) + (v - y1) * (y2 - y1)) /
1257 ((x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1));
1263 }
else if (dpar > 1.0) {
1267 xc = x1 + dpar * (x2 - x1);
1268 yc = y1 + dpar * (y2 - y1);
1276 double dpar = ((u - x2) * (x1 - x2) + (v - y2) * (y1 - y2)) /
1277 ((x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1));
1283 }
else if (dpar > 1.0) {
1287 xc = x2 + dpar * (x1 - x2);
1288 yc = y2 + dpar * (y1 - y2);
1293 if (std::abs(u - xc) < xtol && std::abs(v - yc) < ytol)
return true;
1303bool ViewFEMesh::PlaneCut(
double x1,
double y1,
double z1,
double x2,
double y2,
1304 double z2, TMatrixD& xMat) {
1308 TMatrixD cutMat(3, 3);
1309 dataCut[0] = project[0][0];
1310 dataCut[1] = project[1][0];
1311 dataCut[2] = x1 - x2;
1312 dataCut[3] = project[0][1];
1313 dataCut[4] = project[1][1];
1314 dataCut[5] = y1 - y2;
1315 dataCut[6] = project[0][2];
1316 dataCut[7] = project[1][2];
1317 dataCut[8] = z1 - z2;
1318 cutMat.SetMatrixArray(dataCut.GetArray());
1323 (cutMat(1, 1) * cutMat(2, 2) - cutMat(1, 2) * cutMat(2, 1)) -
1325 (cutMat(1, 0) * cutMat(2, 2) - cutMat(1, 2) * cutMat(2, 0)) +
1327 (cutMat(1, 0) * cutMat(2, 1) - cutMat(1, 1) * cutMat(2, 0));
1330 if (cutDet == 0)
return false;
1333 TArrayD dataCoords(3);
1334 TMatrixD coordMat(3, 1);
1335 dataCoords[0] = x1 - project[2][0];
1336 dataCoords[1] = y1 - project[2][1];
1337 dataCoords[2] = z1 - project[2][2];
1338 coordMat.SetMatrixArray(dataCoords.GetArray());
1342 xMat = cutMat * coordMat;
1345 if (xMat(2, 0) < 0 || xMat(2, 0) > 1)
return false;
1353bool ViewFEMesh::PlaneCoords(
double x,
double y,
double z,
1354 const TMatrixD& projMat, TMatrixD& xMat) {
1357 TArrayD dataCoords(3);
1358 TMatrixD coordMat(3, 1);
1362 coordMat.SetMatrixArray(dataCoords.GetArray());
1363 xMat = projMat * coordMat;
1375bool ViewFEMesh::IsInPolygon(
double x,
double y, std::vector<double>& px,
1376 std::vector<double>& py,
bool& edge) {
1379 int pN = (int)px.size();
1382 if (pN < 2)
return false;
1384 if (pN == 2)
return OnLine(px[0], py[0], px[1], py[1], x, y);
1387 double px_min = px[0], py_min = py[0];
1388 double px_max = px[0], py_max = py[0];
1389 for (
int i = 0; i < pN; i++) {
1390 px_min = std::min(px_min, px[i]);
1391 py_min = std::min(py_min, py[i]);
1392 px_max = std::max(px_max, px[i]);
1393 py_max = std::max(py_max, py[i]);
1397 double xtol = 1.0e-10 * std::max(std::abs(px_min), std::abs(px_max));
1398 double ytol = 1.0e-10 * std::max(std::abs(py_min), std::abs(py_max));
1399 if (xtol <= 0) xtol = 1.0e-10;
1400 if (ytol <= 0) ytol = 1.0e-10;
1403 if (std::abs(px_max - px_min) < xtol) {
1404 edge = (y > (py_min - ytol) && y < (py_max + ytol) &&
1405 std::abs(px_max + px_min - 2 * x) < xtol);
1409 if (std::abs(py_max - py_min) < ytol) {
1410 edge = (x > (px_min - xtol) && x < (px_max + xtol) &&
1411 std::abs(py_max + py_min - 2 * y) < ytol);
1416 double xinf = px_min - std::abs(px_max - px_min);
1417 double yinf = py_min - std::abs(py_max - py_min);
1424 while (!done && niter < 100) {
1432 for (
int i = 0; (done && i < pN); i++) {
1435 if (OnLine(px[i % pN], py[i % pN], px[(i + 1) % pN], py[(i + 1) % pN], x,
1443 double xc = 0., yc = 0.;
1444 if (LinesCrossed(x, y, xinf, yinf, px[i % pN], py[i % pN],
1445 px[(i + 1) % pN], py[(i + 1) % pN], xc, yc))
1450 if (OnLine(x, y, xinf, yinf, px[i], py[i])) {
1453 xinf = px_min -
RndmUniform() * std::abs(px_max - xinf);
1454 yinf = py_min -
RndmUniform() * std::abs(py_max - yinf);
1465 std::cerr << m_className <<
"::IsInPolygon: unable to determine whether ("
1466 << x <<
", " << y <<
") is inside a polygon. Returning false.\n";
1471 return (ncross != 2 * (ncross / 2));
1480void ViewFEMesh::ClipToView(std::vector<double>& px, std::vector<double>& py,
1481 std::vector<double>& cx, std::vector<double>& cy) {
1484 int pN = (int)px.size();
1491 std::vector<double> vx;
1492 vx.push_back(m_xMin);
1493 vx.push_back(m_xMax);
1494 vx.push_back(m_xMax);
1495 vx.push_back(m_xMin);
1496 std::vector<double> vy;
1497 vy.push_back(m_yMax);
1498 vy.push_back(m_yMax);
1499 vy.push_back(m_yMin);
1500 vy.push_back(m_yMin);
1501 int vN = (int)vx.size();
1507 for (
int i = 0; i < pN; i++) {
1513 for (
int j = 0; j < vN; j++) {
1517 if (OnLine(vx[j % vN], vy[j % vN], vx[(j + 1) % vN], vy[(j + 1) % vN],
1521 cx.push_back(px[i]);
1522 cy.push_back(py[i]);
1530 if (OnLine(px[i % pN], py[i % pN], px[(i + 1) % pN], py[(i + 1) % pN],
1534 cx.push_back(vx[j]);
1535 cy.push_back(vy[j]);
1547 for (
int j = 0; j < vN; j++) {
1551 double xc = 0., yc = 0.;
1552 if (LinesCrossed(vx[j % vN], vy[j % vN], vx[(j + 1) % vN],
1553 vy[(j + 1) % vN], px[i % pN], py[i % pN],
1554 px[(i + 1) % pN], py[(i + 1) % pN], xc, yc)) {
1565 for (
int j = 0; j < vN; j++) {
1570 if (IsInPolygon(vx[j], vy[j], px, py, edge)) {
1573 cx.push_back(vx[j]);
1574 cy.push_back(vy[j]);
1579 for (
int i = 0; i < pN; i++) {
1584 if (IsInPolygon(px[i], py[i], vx, vy, edge)) {
1587 cx.push_back(px[i]);
1588 cy.push_back(py[i]);
bool m_zMirrorPeriodic
Mirror periodicity in z.
bool m_yPeriodic
Simple periodicity in y.
bool m_yMirrorPeriodic
Mirror periodicity in y.
bool m_xPeriodic
Simple periodicity in x.
bool m_zPeriodic
Simple periodicity in z.
bool m_xMirrorPeriodic
Mirror periodicity in x.
Base class for components based on finite-element field maps.
std::vector< Material > materials
std::vector< Node > nodes
virtual bool GetVoltageRange(double &vmin, double &vmax)
Calculate the voltage range [V].
std::vector< Element > elements
void SetXaxisTitle(const char *xtitle)
void CreateDefaultAxes()
Create a default set of custom-made axes.
void SetCanvas(TCanvas *c)
Set the canvas to be painted on.
void SetComponent(ComponentFieldMap *comp)
Set the component from which to retrieve the mesh and field.
void SetYaxisTitle(const char *ytitle)
bool Plot()
Plot method to be called by user.
void SetPlane(double fx, double fy, double fz, double x0, double y0, double z0)
Set the projection plane.
void SetArea()
Set area to be plotted to the default.
void SetXaxis(TGaxis *ax)
void SetYaxis(TGaxis *ay)
void SetDefaultProjection()
Reset the projection plane.
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
PlottingEngineRoot plottingEngine
DoubleAc sqrt(const DoubleAc &f)