Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
REveGeoTopNode.cxx
Go to the documentation of this file.
1
3
6
8#include <ROOT/RGeomData.hxx>
9#include <ROOT/RWebWindow.hxx>
10#include <ROOT/REveManager.hxx>
12
14
15#include <ROOT/REveUtil.hxx>
16#include <ROOT/RLogger.hxx>
17#include <ROOT/REveUtil.hxx>
18#include "TBufferJSON.h"
19#include "TMath.h"
20
21#include "TGeoCompositeShape.h"
22#include "TGeoManager.h"
23#include "TClass.h"
24#include "TGeoNode.h"
25#include "TGeoMatrix.h"
26#include "TBase64.h"
27#include "TStopwatch.h"
28
29
30#include <cassert>
31#include <iostream>
32#include <regex>
33#include <nlohmann/json.hpp>
34
35
36using namespace ROOT::Experimental;
37
38thread_local ElementId_t gSelId;
39
40#define REVEGEO_DEBUG
41#ifdef REVEGEO_DEBUG
42#define REVEGEO_DEBUG_PRINT(fmt, ...) printf(fmt, ##__VA_ARGS__)
43#else
44#define REVEGEO_DEBUG_PRINT(fmt, ...)
45#endif
46
47
49
50/*
51namespace {
52void PrintStackPath(const std::vector<int>& stack)
53{
54 printf("Path: ");
55
56 for (auto idx : stack)
57 printf("/%d", idx);
58
59 printf("\n");
60}
61}*/
62
63
64bool REveGeomDescription::ChangeEveVisibility(const std::vector<int> &stack, ERnrFlags flags, bool on)
65{
66 std::vector<RGeomNodeVisibility> &visVec = (flags == kRnrSelf) ? fVisibilitySelf : fVisibilityRec;
67
68 for (auto iter = visVec.begin(); iter != visVec.end(); iter++) {
69 if (iter->stack == stack) {
70 // AMT TODO remove path fom the vsibilirt vector if it is true
71 iter->visible = on;
72 return true;
73 }
74 }
75
76 visVec.emplace_back(stack, on);
77 return true;
78}
79
81{
82 std::vector<int> stack = fApex.GetIndexStack();
83 stack.insert(stack.end(), iStack.begin(), iStack.end());
84
85 auto isVisible = [&stack](std::vector<RGeomNodeVisibility> &visVec) -> bool {
86 for (auto &visVecEl : visVec) {
87 /*
88 printf("compare ======\n");
89 PrintStackPath(stack);
90 PrintStackPath(visVecEl.stack);
91*/
92 if (stack == visVecEl.stack)
93 return visVecEl.visible ? 1 : 0;
94 }
95 return true;
96 };
97
100
101 item.SetLogicalVisibility(visRec);
102 item.SetPhysicalVisibility(visSelf);
103
104 //return RGeoItem(node.name, node.chlds.size(), node.id, node.color, node.material,
105 // visRec, vis);
106}
107
108void REveGeomDescription::InitPath(const std::vector<std::string>& path)
109{
110 fApex.SetFromPath(path);
111 Build(fApex.GetNode()->GetVolume()); // rebuild geo-webviewer
112}
113
115{
116 // visibility self
117 for (auto &visVecEl : fVisibilitySelf) {
118 if (nodeStack == visVecEl.stack) {
119 return false;
120 }
121 }
122
123 // visibility recurse/children
124 for (auto &visVecEl : fVisibilityRec) {
125 bool inside =
126 nodeStack.size() >= visVecEl.stack.size() && std::equal(visVecEl.stack.begin(), visVecEl.stack.end(), nodeStack.begin());
127 if (inside)
128 return false;
129 }
130
131 return true;
132}
133
134////////////////////////////////////////////////////////////////////////////////
135/// Table signal handling
136
137void REveGeomHierarchy::WebWindowCallback(unsigned connid, const std::string &arg)
138{
139 using namespace std::string_literals;
141
143 if (arg.compare(0, 6, "CDTOP:") == 0)
144 {
145 std::vector<std::string> ep;
146 eveDesc.InitPath(ep);
147 fDesc.IssueSignal(this, "CdTop");
148 fWebWindow->Send(connid, "RELOAD"s);
149 }
150 else if (arg.compare(0, 5, "CDUP:") == 0)
151 {
152 std::vector<std::string> result = eveDesc.GetApexPath();
153 result.pop_back();
154 eveDesc.InitPath(result);
155 fDesc.IssueSignal(this, "CdUp");
156 fWebWindow->Send(connid, "RELOAD"s);
157 }
158 else if (arg.compare(0, 8, "SETAPEX:") == 0) {
159 auto path = TBufferJSON::FromJSON<std::vector<std::string>>(arg.substr(8));
160
161 //const std::vector<int> &sstack = fDesc.GetSelectedStack();
162 // std::vector<std::string> sspath = fDesc.MakePathByStack(sstack);
163 std::vector<std::string> result = eveDesc.GetApexPath();
164 if (path->size() > 1) {
165 result.insert(result.end(), path->begin() + 1, path->end());
166 eveDesc.InitPath(result);
167 fDesc.IssueSignal(this, "SelectTop");
168 fWebWindow->Send(connid, "RELOAD"s);
169 }
170 }
171 else if ((arg.compare(0, 7, "SETVI0:") == 0) || (arg.compare(0, 7, "SETVI1:") == 0)) {
172 {
174 bool on = (arg[5] == '1');
175 auto path = TBufferJSON::FromJSON<std::vector<std::string>>(arg.substr(7));
176 // Get integer stack from string stack
177 std::vector<int> base = eveDesc.GetIndexStack();
178 std::vector<int> stack = fDesc.MakeStackByPath(*path);
179 stack.insert(stack.begin(), base.begin(), base.end());
180
181 if (eveDesc.ChangeEveVisibility(stack, REveGeomDescription::kRnrChildren , on)) {
182 std::cout << "Set visibilty rnr CHIDLREN \n";
184 }
185 }
186 }
187 else if ((arg.compare(0, 5, "SHOW:") == 0) || (arg.compare(0, 5, "HIDE:") == 0)) {
188 {
189 auto path = TBufferJSON::FromJSON<std::vector<std::string>>(arg.substr(5));
190 bool on = (arg.compare(0, 5, "SHOW:") == 0);
191 // Get integer stack from string stack
192
193 std::vector<int> base = eveDesc.GetIndexStack();
194 std::vector<int> stack = fDesc.MakeStackByPath(*path);
195 stack.insert(stack.begin(), base.begin(), base.end());
196
197 if (path && eveDesc.ChangeEveVisibility(stack, REveGeomDescription::kRnrSelf, on)) {
198 std::cout << "Set visibilty rnr PHY \n";
201 }
202 }
203 }
204
205 else {
207 }
208}
209
210
211////////////////////////////////////////////////////////////////////////////////
212////////////////////////////////////////////////////////////////////////////////
213////////////////////////////////////////////////////////////////////////////////
214
220
222{
223 fPath = absPath;
225}
226
227TGeoNode *REveGeomDescription::Apex::LocateNodeWithPath(const std::vector<std::string> &path) const
228{
229 TGeoNode *top = REveGeomDescription::GetGeoManager()->GetTopNode();
230 // printf("Top node name from geoData name (%s)\n", top->GetName());
231 for (size_t t = 0; t < path.size(); t++) {
232 std::string s = path[t];
233 std::cout << s << std::endl;
234 TGeoNode *ntop = top->GetVolume()->FindNode(s.c_str());
235 if (!ntop)
236 throw std::runtime_error("Apex::LocateNodeWithPath(), can't locate node with path " + s);
237 top = ntop;
238 }
239 return top;
240}
241
243{
244 if (fPath.empty())
245 return "";
246
247 std::ostringstream oss;
248
249 oss << fPath[0];
250
251 for (size_t i = 1; i < fPath.size(); ++i)
252 oss << "/" << fPath[i];
253
254 return oss.str();
255}
256
258{
259 std::vector<int> indexStack;
260
261 TGeoNode* current = REveGeomDescription::GetGeoManager()->GetTopNode();
262
263 // optional: skip first if it is top itself
264 size_t start = 0;
265 std::vector<std::string> nameStack = fPath;
266 if (!nameStack.empty() && nameStack[0] == current->GetName())
267 start = 1;
268
269 for (size_t i = start; i < nameStack.size(); ++i)
270 {
271 const std::string& targetName = nameStack[i];
272
273 TGeoVolume* vol = current->GetVolume();
274
275 int nd = vol->GetNdaughters();
276
277 int foundIndex = -1;
278
279 for (int j = 0; j < nd; ++j)
280 {
281 TGeoNode* daughter = vol->GetNode(j);
282
283 if (targetName == daughter->GetName())
284 {
285 foundIndex = j;
286 current = daughter;
287 break;
288 }
289 }
290
291 if (foundIndex == -1)
292 {
293 std::cerr << "Node not found: " << targetName << std::endl;
294 return {};
295 }
296
297 indexStack.push_back(foundIndex);
298 }
299
300 // PrintStackPath(indexStack);
301 return indexStack;
302}
303
305{
307 if (!s_geoManager) {
308 throw std::runtime_error("Critical Error: Failed to import geometry file");
309 }
310}
311
312////////////////////////////////////////////////////////////////////////////////
313///
314/// Constructor.
315
317{
318 // this below will be obsolete
319 fDesc.AddSignalHandler(this, [this](const std::string &kind) { ProcessSignal(kind); });
320 fDesc.ImportFile(filename);
321
322
323 fWebHierarchy = std::make_shared<REveGeomHierarchy>(fDesc, true);
324 fWebHierarchy->SetReceiver(this);
325}
326
327void REveGeoTopNodeData::InitPath(const std::string &path)
328{
329 std::regex re(R"([/\\]+)"); // split on one or more slashes
330 std::sregex_token_iterator it(path.begin(), path.end(), re, -1);
331 std::sregex_token_iterator end;
332 std::vector<std::string> result;
333
334 for (; it != end; ++it) {
335 if (!it->str().empty()) { // skip empty parts
336 result.push_back(*it);
337 }
338 }
339
340 fDesc.InitPath(result);
341
342 for (auto &el : fNieces) {
343 REveGeoTopNodeViz *etn = dynamic_cast<REveGeoTopNodeViz *>(el);
344 etn->BuildDesc();
345 }
346}
347
349{
350
351 for (auto &el : fNieces) {
352 REveGeoTopNodeViz *etn = dynamic_cast<REveGeoTopNodeViz *>(el);
353 etn->VisibilityChanged(on, flag, path);
354 }
355}
356
357////////////////////////////////////////////////////////////////////////////////
358
359void REveGeoTopNodeData::SetChannel(unsigned connid, int chid)
360{
361 fWebHierarchy->Show({gEve->GetWebWindow(), connid, chid});
362}
363
364////////////////////////////////////////////////////////////////////////////////
365void REveGeoTopNodeData::ProcessSignal(const std::string &kind)
366{
368 if ((kind == "SelectTop") || (kind == "CdTop") || (kind == "CdUp"))
369 {
370 for (auto &el : fNieces) {
371 REveGeoTopNodeViz *etn = dynamic_cast<REveGeoTopNodeViz *>(el);
372 etn->BuildDesc();
373 }
374 }
375 else if (kind == "HighlightItem") {
376 /*
377 printf("REveGeoTopNodeData element highlighted --------------------------------"\n);
378 */
379
380 } else if (kind == "ClickItem") {
381 printf("REveGeoTopNodeData element CLICKED selected --------------------------------\n");
382 auto sstack = fDesc.GetClickedItem();
383 std::set<int> ss;
384
385 for (auto &n : fNieces) {
386 REveGeoTopNodeViz* viz = dynamic_cast<REveGeoTopNodeViz*>(n);
387 viz->GetIndicesFromBrowserStack(sstack, ss);
388 bool multi = false;
389 bool secondary = true;
390 gEve->GetSelection()->NewElementPicked(n->GetElementId(), multi, secondary, ss);
391 }
392 }
393}
394
395////////////////////////////////////////////////////////////////////////////////
396/// Fill core part of JSON representation.
397
404
405////////////////////////////////////////////////////////////////////////////////
406// REveGeoTopNodeViz
407////////////////////////////////////////////////////////////////////////////////
408////////////////////////////////////////////////////////////////////////////////
409
414
419
421{
423 {
424 if (it.GetLevel() > fGeoData->fDesc.GetVisLevel()) {
425 it.Skip();
426 return false;
427 }
428 }
429 else if (fMode == EMode::kModeLeafOnly)
430 {
431 // printf("accep mkod eleaf node ptr %p \n", (void*)it.GetNode(it.GetLevel()));
432 if (it.GetNode(it.GetLevel())->GetNdaughters())
433 return false;
434 }
435 else if (fMode == EMode::kModeMixed)
436 {
437 if (it.GetLevel() > fGeoData->fDesc.GetVisLevel()) {
438 if (skip) it.Skip();
439 return false;
440 }
441 // printf("accep mkod eleaf node ptr %p \n", (void*)it.GetNode(it.GetLevel()));
442 if (it.GetNode(it.GetLevel())->GetNdaughters())
443 return false;
444 }
445
446 return true;
447}
448
449std::string REveGeoTopNodeViz::GetHighlightTooltip(const std::set<int> & set) const
450{
452 if (set.empty()) {
453 return "";
454 } else {
455 auto it = set.begin();
456 int pos = *it;
457 //const BNode &bn = fNodes[pos];
458 std::cout << "highlight node with ID " << pos << "\n";
459
460 std::string res = "GeoNode name";
461
463 TGeoIterator git(top->GetVolume());
464 TGeoNode *node;
465 int i = 0;
466 TString path;
467 while ((node = git.Next()))
468 {
469 if (!AcceptNode(git))
470 continue;
471 if (i == pos) {
472 git.GetPath(path);
473 res = path;
474 break;
475 }
476 i++;
477 }
478 return res;
479 }
480}
481
483{
484 // locate top node
486
487 fNodes.clear();
488 fShapes.clear();
489 // shape array
490 std::set<TGeoShape *> shapes;
492 timer.Start();
493 CollectShapes(top, shapes, fShapes);
494 std::cout << "Shape size " << shapes.size() << "\n";
495
496 timer.Stop();
497
498 printf("Real time: %.3f s\n", timer.RealTime());
499 printf("CPU time: %.3f s\n", timer.CpuTime());
500
501 // node array
502 timer.Start();
504 std::cout << "Node size " << fNodes.size() << "\n";
505
506 timer.Stop();
507
508 printf("NODES Real time: %.3f s\n", timer.RealTime());
509 printf("NODES CPU time: %.3f s\n", timer.CpuTime());
510
512}
513
514void REveGeoTopNodeViz::CollectNodes(TGeoVolume *volume, std::vector<BNode> &bnl, std::vector<BShape> &browsables)
515{
516 printf("collect nodes \n");
517 TGeoIterator it(volume);
518 TGeoNode *node;
519 int nodeId = 0;
520
521 std::vector<int> apexStack = fGeoData->RefDescription().GetIndexStack();
522
523 // get top node transformation
525 {
527 for (int idx : apexStack) {
528 inode = inode->GetDaughter(idx);
529 global.Multiply(inode->GetMatrix());
530 }
531 }
532
533 while ((node = it.Next()))
534 {
535 if (!AcceptNode(it))
536 continue;
537
538 TGeoHMatrix full = global; // identity if global is identity
539 full.Multiply(it.GetCurrentMatrix());
540 const TGeoMatrix *mat = &full;
541
542 // const TGeoMatrix *mat = it.GetCurrentMatrix();
543 const Double_t *t = mat->GetTranslation(); // size 3
544 const Double_t *r = mat->GetRotationMatrix(); // size 9 (3x3)
545
546 Double_t m[16];
547 if (mat->IsScale()) {
548 const Double_t *s = mat->GetScale();
549 m[0] = r[0] * s[0];
550 m[1] = r[3] * s[0];
551 m[2] = r[6] * s[0];
552 m[3] = 0;
553 m[4] = r[1] * s[1];
554 m[5] = r[4] * s[1];
555 m[6] = r[7] * s[1];
556 m[7] = 0;
557 m[8] = r[2] * s[2];
558 m[9] = r[5] * s[2];
559 m[10] = r[8] * s[2];
560 m[11] = 0;
561 m[12] = t[0];
562 m[13] = t[1];
563 m[14] = t[2];
564 m[15] = 1;
565 } else {
566 m[0] = r[0];
567 m[1] = r[3];
568 m[2] = r[6];
569 m[3] = 0;
570 m[4] = r[1];
571 m[5] = r[4];
572 m[6] = r[7];
573 m[7] = 0;
574 m[8] = r[2];
575 m[9] = r[5];
576 m[10] = r[8];
577 m[11] = 0;
578 m[12] = t[0];
579 m[13] = t[1];
580 m[14] = t[2];
581 m[15] = 1;
582 }
583
584 BNode b;
585 b.node = node;
586 b.nodeId = nodeId;
587 b.color = node->GetVolume()->GetLineColor();
588
589 // TString path; it.GetPath(path);
590 // printf("[%d] %d %s \n", node->GetNdaughters(), it.GetLevel(), path.Data());
591
592
593 // set BNode transformation matrix
594 for (int i = 0; i < 16; ++i)
595 b.trans[i] = m[i];
596
597 // find shape
598 TGeoShape *shape = node->GetVolume()->GetShape();
599 b.shapeId = -1; // mark invalid at start
600 for (size_t i = 0; i < browsables.size(); i++) {
601 if (shape == browsables[i].shape) {
602 b.shapeId = i;
603 break;
604 }
605 }
606 assert(b.shapeId >= 0);
607
608
609 // set visibility flag
610 std::vector<int> visStack = apexStack;
611 for (int i = 1; i <= it.GetLevel(); ++i)
612 visStack.push_back(it.GetIndex(i));
613 // PrintStackPath(visStack);
615
616 // printf("Node %d shape id %d \n", (int)bnl.size(), b.shapeId);
617 bnl.push_back(b);
618 nodeId++;
619
620 if (nodeId > 300000) {
621 R__LOG_ERROR(REveLog()) << "Max number of nodes reached ... breaking the loop \n";
622 printf("num nodes locked !!! \n");
623 break;
624 }
625 }
626}
627
628void REveGeoTopNodeViz::CollectShapes(TGeoNode *tnode, std::set<TGeoShape *> &shapes, std::vector<BShape> &browsables)
629{
630 printf("collect shapes \n");
631 TGeoIterator geoit(tnode->GetVolume());
632 TGeoNode *node = nullptr;
633 while ((node = geoit.Next()))
634 {
635 if (!AcceptNode(geoit))
636 continue;
637
638 TGeoVolume *vol = node->GetVolume();
639 if (vol) {
640 TGeoShape *shape = vol->GetShape();
641 if (shape) {
642 auto it = shapes.find(shape);
643 if (it == shapes.end()) {
644 shapes.insert(shape); // use set to avoid duplicates
646 TGeoCompositeShape *compositeShape = dynamic_cast<TGeoCompositeShape *>(shape);
647 int n_seg = 60; // default value in the geo manager and poly shape
648 if (compositeShape)
649 polyShape.BuildFromComposite(compositeShape, n_seg);
650 else
651 polyShape.BuildFromShape(shape, n_seg);
652
653 // printf("[%d] Shape name %s %s \n",(int)browsables.size(), shape->GetName(), shape->ClassName());
654
655 // printf("vertices %lu: \n", polyShape.fVertices.size());
656
657 // create browser shape
659 browserShape.shape = shape;
660 browsables.push_back(browserShape);
661
662 // copy vertices transform vec double to float
663 browsables.back().vertices.reserve(polyShape.GetVertices().size());
664 for (size_t i = 0; i < polyShape.GetVertices().size(); i++)
665 browsables.back().vertices.push_back(polyShape.GetVertices()[i]);
666
667 // copy indices kip the first integer in the sequence of 4
668 for (size_t i = 0; i < polyShape.GetPolyDesc().size(); i += 4) {
669 browsables.back().indices.push_back(polyShape.GetPolyDesc()[i + 1]);
670 browsables.back().indices.push_back(polyShape.GetPolyDesc()[i + 2]);
671 browsables.back().indices.push_back(polyShape.GetPolyDesc()[i + 3]);
672 }
673 // printf("last browsable size indices size %lu \n", browsables.back().indices.size());
674 }
675 }
676 }
677 }
678}
679
681{
682 fRenderData = std::make_unique<REveRenderData>("makeGeoTopNode");
683 for (size_t i = 0; i < fNodes.size(); ++i) {
684
685 UChar_t c[4] = {1, 2, 3, 4};
687 // if (i < 400) printf("%d > %d %d %d %d \n",fNodes[i].color, c[0], c[1], c[2], c[3]);
688 uint32_t v = (c[0] << 16) + (c[1] << 8) + c[2];
689 float pc;
690 std::memcpy(&pc, &v, sizeof(pc));
691 GetRenderData()->PushV(pc);
692 }
693}
694//------------------------------------------------------------------------------
695//------------------------------------------------------------------------------
702
703//------------------------------------------------------------------------------
704
706{
709
710 if (!fGeoData) {
711 j["dataId"] = -1;
712 } else {
713 // std::string json = fGeoData->fDesc.ProduceJson();
714 // j["geomDescription"] = TBase64::Encode(json.c_str());
715 j["dataId"] = fGeoData->GetElementId();
716 }
717 j["visLevel"] = fGeoData ? fGeoData->fDesc.GetVisLevel() : 0;
718
719 // put shapes vector in json array
720 using namespace nlohmann;
721
722 json shapeVertexArr = json::array();
723 int vertexOff = 0;
724
725 json shapeIndexArr = json::array();
726 json shapePolySizeArr = json::array();
727 json shapePolyOffArr = json::array();
728
729 json nodeVisibility = json::array();
730
731 int polyOff = 0;
732
733 // need four integers for
734 for (size_t i = 0; i < fShapes.size(); ++i) {
735 // vertices
736
737 std::copy(fShapes[i].vertices.begin(), fShapes[i].vertices.end(), std::back_inserter(shapeVertexArr));
738
739 int numVertices = int(fShapes[i].vertices.size());
740 // indices
741 // write shape indices with the vertexOff
742 for (size_t p = 0; p < fShapes[i].indices.size(); ++p)
743 shapeIndexArr.push_back(fShapes[i].indices[p] + vertexOff);
744
745 int numIndices = int(fShapes[i].indices.size());
746 shapePolySizeArr.push_back(numIndices);
747 shapePolyOffArr.push_back(polyOff);
748
749 // printf("shape [%d] numIndices %d \n", i, numIndices);
750
752 vertexOff += numVertices / 3;
753 }
754
755 // write vector of shape ids for visible nodes
756 json nodeShapeIds = json::array();
757 json nodeTrans = json::array();
758 json nodeColors = json::array();
759
760 for (size_t i = 0; i < fNodes.size(); ++i) {
761 nodeShapeIds.push_back(fNodes[i].shapeId);
762 nodeVisibility.push_back(fNodes[i].visible);
763 for (int t = 0; t < 16; t++)
764 nodeTrans.push_back(fNodes[i].trans[t]);
765 }
766 // shape basic array
767
768 j["shapeVertices"] = shapeVertexArr;
769
770 // shape basic indices array
771 j["shapeIndices"] = shapeIndexArr;
772
773 // shape poly offset array
774 j["shapeIndicesOff"] = shapePolyOffArr;
775 j["shapeIndicesSize"] = shapePolySizeArr;
776
777 j["nodeShapeIds"] = nodeShapeIds;
778 j["nodeTrans"] = nodeTrans;
779 j["nodeVisibility"] = nodeVisibility;
780 j["fSecondarySelect"] = fAlwaysSecSelect;
781
782
783
784
785 // ship bounding box info
787 TGeoVolume *vol = top->GetVolume();
788 TGeoShape *shape = vol->GetShape();
789 shape->ComputeBBox();
790 TGeoBBox *box = dynamic_cast<TGeoBBox *>(shape);
791 if (box) {
792 const Double_t *origin = box->GetOrigin();
793
794 printf("BBox center: (%f, %f, %f)\n", origin[0], origin[1], origin[2]);
795 //printf("origin lengths: (%f, %f, %f)\n", origin[0], origin[1], origin[2]);
796
797 auto jbb = json::array();
798 jbb.push_back(origin[0] - box->GetDX());
799 jbb.push_back(origin[0] + box->GetDX());
800 jbb.push_back(origin[1] - box->GetDY());
801 jbb.push_back(origin[1] + box->GetDY());
802 jbb.push_back(origin[2] - box->GetDZ());
803 jbb.push_back(origin[2] + box->GetDZ());
804 j["bbox"] = jbb;
805 }
806 // std::cout << "Write Core json " << j.dump(1) << "\n";
807 return ret;
808}
809
811{
812 if (fGeoData) {
815 }
816}
817
818void REveGeoTopNodeViz::GetIndicesFromBrowserStack(const std::vector<int> &stack, std::set<int> &res)
819{
821 TGeoIterator it(top->GetVolume());
822 std::vector<int> nodeStack;
823 int cnt = 0;
824 TGeoNode *node;
825
826 while ((node = it.Next())) {
827 int level = it.GetLevel();
828
829 bool accept = AcceptNode(it, false);
830
831 nodeStack.resize(level);
832 if (level > 0)
833 nodeStack[level - 1] = it.GetIndex(level);
834
835 bool inside = nodeStack.size() >= stack.size() && std::equal(stack.begin(), stack.end(), nodeStack.begin());
836 if (inside) {
837 res.insert(cnt);
838 } // rnr flags
839 if (accept) cnt++;
840 } // while it
841
842 printf("GetIndicesFromBrowserStack stack size %zu res size %zu\n", stack.size(), res.size());
843}
844
846{
847 // function argument is full stack, we remove the apex path
848 size_t apexDepth = fGeoData->RefDescription().GetApexPath().size();
849 std::vector<int> stack(iStack.begin() + apexDepth, iStack.end());
850
851 // PrintStackPath(stack);
852
854 TGeoIterator it(top->GetVolume());
855 std::vector<int> nodeStack;
856 int cnt = 0;
857 TGeoNode *node;
858 while ((node = it.Next())) {
859
860 int level = it.GetLevel();
861 if (!AcceptNode(it))
862 continue;
863
864 nodeStack.resize(level);
865 if (level > 0)
866 nodeStack[level - 1] = it.GetIndex(level);
867
869 /*
870 printf("nODEcompare ======\n");
871 PrintStackPath(stack);
872 PrintStackPath(nodeStack);
873 */
874 if (nodeStack == stack) {
875 fNodes[cnt].visible = on;
876
877 break;
878 }
879 } else {
880 bool inside = nodeStack.size() >= stack.size() && std::equal(stack.begin(), stack.end(), nodeStack.begin());
881 if (inside) {
882 fNodes[cnt].visible = on;
883 }
884 } // rnr flags
885 cnt++;
886 } // while it
888}
nlohmann::json json
thread_local ElementId_t gSelId
#define R__LOG_ERROR(...)
Definition RLogger.hxx:357
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
char Text_t
General string (char)
Definition RtypesCore.h:76
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void on
Option_t Option_t TPoint TPoint const char mode
virtual Int_t WriteCoreJson(nlohmann::json &cj, Int_t rnr_offset)
Write core json.
std::unique_ptr< REveRenderData > fRenderData
! Vertex / normal / triangle index information for rendering.
ElementId_t GetElementId() const
REveRenderData * GetRenderData() const
REveGeoManagerHolder Exception-safe global variable holders.
Definition REveUtil.hxx:90
void SetChannel(unsigned connid, int chid)
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Fill core part of JSON representation.
void InitPath(const std::string &path)
void VisibilityChanged(bool on, REveGeomDescription::ERnrFlags flag, const std::vector< int > &path)
REveGeoTopNodeData(const REveGeoTopNodeData &)=delete
void ProcessSignal(const std::string &)
void GetIndicesFromBrowserStack(const std::vector< int > &stack, std::set< int > &outStack)
std::string GetHighlightTooltip(const std::set< int > &secondary_idcs) const override
REveGeoTopNodeViz(const REveGeoTopNodeViz &)=delete
void CollectNodes(TGeoVolume *volume, std::vector< BNode > &bnl, std::vector< BShape > &browsables)
void BuildRenderData() override
Write transformation Matrix to render data.
void VisibilityChanged(bool on, REveGeomDescription::ERnrFlags flag, const std::vector< int > &path)
Int_t WriteCoreJson(nlohmann::json &j, Int_t rnr_offset) override
Write core json.
void SetGeoData(REveGeoTopNodeData *d, bool rebuild=true)
bool AcceptNode(TGeoIterator &it, bool skip=true) const
void CollectShapes(TGeoNode *node, std::set< TGeoShape * > &shapes, std::vector< BShape > &browsables)
void SetFromPath(std::vector< std::string > absPath)
TGeoNode * LocateNodeWithPath(const std::vector< std::string > &path) const
void InitPath(const std::vector< std::string > &path)
bool ChangeEveVisibility(const std::vector< int > &stack, ERnrFlags rnrFlag, bool on)
std::vector< RGeomNodeVisibility > fVisibilitySelf
virtual void RefineGeoItem(ROOT::RGeoItem &item, const std::vector< int > &stack) override
Method which allows to add/modify information in RGeoItem which will be provided to client - like tit...
bool GetVisiblityForStack(const std::vector< int > &stack)
const std::vector< std::string > & GetApexPath() const
std::vector< RGeomNodeVisibility > fVisibilityRec
virtual void WebWindowCallback(unsigned connid, const std::string &kind) override
Table signal handling.
REveSelection * GetSelection() const
std::shared_ptr< ROOT::RWebWindow > GetWebWindow() const
void NewElementPicked(ElementId_t id, bool multi, bool secondary, const std::set< int > &secondary_idcs={})
Called from GUI when user picks or un-picks an element.
static void ColorFromIdx(Color_t ci, UChar_t col[4], Bool_t alpha=kTRUE)
Fill col with RGBA values corresponding to index ci.
Definition REveUtil.cxx:138
Representation of single item in the geometry browser.
Definition RGeomData.hxx:88
std::vector< RGeomNode > fDesc
! converted description, send to client
void SetVisLevel(int lvl=3)
Set maximal visible level.
void IssueSignal(const void *handler, const std::string &kind)
Issue signal, which distributed on all handlers - excluding source handler.
int GetVisLevel() const
Returns maximal visible level.
std::vector< int > MakeStackByPath(const std::vector< std::string > &path)
Produce stack based on string path Used to highlight geo volumes by browser hover event.
void Build(TGeoManager *mgr, const std::string &volname="")
Collect information about geometry hierarchy into flat list like it done in JSROOT ClonedNodes....
std::shared_ptr< RWebWindow > fWebWindow
! web window to show geometry
virtual void WebWindowCallback(unsigned connid, const std::string &arg)
Process data from client.
RGeomDescription & fDesc
! geometry description, shared with external
const_iterator begin() const
const_iterator end() const
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:36
Box class.
Definition TGeoBBox.h:18
Composite shapes are Boolean combinations of two or more shape components.
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:459
A geometry iterator.
Definition TGeoNode.h:249
const TGeoMatrix * GetCurrentMatrix() const
Returns global matrix for current node.
Int_t GetLevel() const
Definition TGeoNode.h:295
TGeoNode * GetNode(Int_t level) const
Returns current node at a given level.
Int_t GetIndex(Int_t i) const
Definition TGeoNode.h:294
TGeoNode * Next()
Returns next node.
void Skip()
Stop iterating the current branch.
The manager class for any TGeo geometry.
Definition TGeoManager.h:46
static TGeoManager * Import(const char *filename, const char *name="", Option_t *option="")
static function Import a geometry from a gdml or ROOT file
Geometrical transformation package.
Definition TGeoMatrix.h:39
A node represent a volume positioned inside another.They store links to both volumes and to the TGeoM...
Definition TGeoNode.h:39
TGeoVolume * GetVolume() const
Definition TGeoNode.h:100
Int_t GetNdaughters() const
Definition TGeoNode.h:92
Base abstract class for all shapes.
Definition TGeoShape.h:25
virtual void ComputeBBox()=0
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
Int_t GetNdaughters() const
Definition TGeoVolume.h:363
TGeoNode * GetNode(const char *name) const
get the pointer to a daughter node
TGeoShape * GetShape() const
Definition TGeoVolume.h:191
TGeoNode * FindNode(const char *name) const
search a daughter inside the list of nodes
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Stopwatch class.
Definition TStopwatch.h:28
Basic string class.
Definition TString.h:138
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
const Int_t n
Definition legend1.C:16
Namespace for ROOT features in testing.
Definition TROOT.h:100
ROOT::RLogChannel & REveLog()
Log channel for Eve diagnostics.
Definition REveTypes.cxx:51
R__EXTERN REveManager * gEve
TMarker m
Definition textangle.C:8