Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RGeomData.cxx
Go to the documentation of this file.
1// @(#)root/eve7:$Id$
2// Author: Sergey Linev, 14.12.2018
3
4/*************************************************************************
5 * Copyright (C) 1995-2019, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include <ROOT/RGeomData.hxx>
13
16#include <ROOT/RLogger.hxx>
17#include "CsgOps.h"
18
19#include "TMath.h"
20#include "TColor.h"
21#include "TROOT.h"
22#include "TGeoNode.h"
23#include "TGeoVolume.h"
24#include "TGeoBBox.h"
25#include "TGeoManager.h"
26#include "TGeoMatrix.h"
27#include "TGeoMedium.h"
28#include "TGeoMaterial.h"
29#include "TGeoBoolNode.h"
30#include "TGeoCompositeShape.h"
31#include "TBuffer3D.h"
32#include "TBufferJSON.h"
33
34#include <algorithm>
35
36namespace ROOT {
37namespace Experimental {
38
39
41{
42 static RLogChannel sLog("ROOT.Geom");
43 return sLog;
44}
45
46
47/** Iterator of hierarchical geometry structures */
48
50
52 int fParentId{-1};
53 unsigned fChild{0};
54 int fNodeId{0};
55
56 std::vector<int> fStackParents;
57 std::vector<int> fStackChilds;
58
59public:
60
62
63 const std::string &GetName() const { return fDesc.fDesc[fNodeId].name; }
64
65 bool IsValid() const { return fNodeId >= 0; }
66
67 int GetNodeId() const { return fNodeId; }
68
69 bool HasChilds() const { return (fNodeId < 0) ? true : fDesc.fDesc[fNodeId].chlds.size() > 0; }
70
71 int NumChilds() const { return (fNodeId < 0) ? 1 : fDesc.fDesc[fNodeId].chlds.size(); }
72
73 bool Enter()
74 {
75 if (fNodeId < 0) {
76 Reset();
77 fNodeId = 0;
78 return true;
79 }
80
81 if (fNodeId >= (int) fDesc.fDesc.size())
82 return false;
83
84 auto &node = fDesc.fDesc[fNodeId];
85 if (node.chlds.size() == 0) return false;
86 fStackParents.emplace_back(fParentId);
87 fStackChilds.emplace_back(fChild);
89 fChild = 0;
90 fNodeId = node.chlds[fChild];
91 return true;
92 }
93
94 bool Leave()
95 {
96 if (fStackParents.size() == 0) {
97 fNodeId = -1;
98 return false;
99 }
100 fParentId = fStackParents.back();
101 fChild = fStackChilds.back();
102
103 fStackParents.pop_back();
104 fStackChilds.pop_back();
105
106 if (fParentId < 0) {
107 fNodeId = 0;
108 } else {
110 }
111 return true;
112 }
113
114 bool Next()
115 {
116 // does not have parents
117 if ((fNodeId <= 0) || (fParentId < 0)) {
118 Reset();
119 return false;
120 }
121
122 auto &prnt = fDesc.fDesc[fParentId];
123 if (++fChild >= prnt.chlds.size()) {
124 fNodeId = -1; // not valid node, only Leave can be called
125 return false;
126 }
127
128 fNodeId = prnt.chlds[fChild];
129 return true;
130 }
131
132 bool Reset()
133 {
134 fParentId = -1;
135 fNodeId = -1;
136 fChild = 0;
137 fStackParents.clear();
138 fStackChilds.clear();
139
140 return true;
141 }
142
143 bool NextNode()
144 {
145 if (Enter()) return true;
146
147 if (Next()) return true;
148
149 while (Leave()) {
150 if (Next()) return true;
151 }
152
153 return false;
154 }
155
156 /** Navigate to specified path - path specified as string and should start with "/" */
157 bool Navigate(const std::string &path)
158 {
159 size_t pos = path.find("/");
160 if (pos != 0) return false;
161
162 Reset(); // set to the top of element
163
164 while (++pos < path.length()) {
165 auto last = pos;
166
167 pos = path.find("/", last);
168
169 if (pos == std::string::npos) pos = path.length();
170
171 std::string folder = path.substr(last, pos-last);
172
173 if (!Enter()) return false;
174
175 bool find = false;
176
177 do {
178 find = (folder.compare(GetName()) == 0);
179 } while (!find && Next());
180
181 if (!find) return false;
182 }
183
184 return true;
185 }
186
187 /** Navigate to specified path */
188 bool Navigate(const std::vector<std::string> &path)
189 {
190 Reset(); // set to the top of element
191
192 for (auto &folder : path) {
193
194 if (!Enter()) return false;
195
196 bool find = false;
197
198 do {
199 find = (folder.compare(GetName()) == 0);
200 } while (!find && Next());
201
202 if (!find) return false;
203 }
204
205 return true;
206 }
207
208
209 /// Returns array of ids to currently selected node
210 std::vector<int> CurrentIds() const
211 {
212 std::vector<int> res;
213 if (IsValid()) {
214 for (unsigned n=1;n<fStackParents.size();++n)
215 res.emplace_back(fStackParents[n]);
216 if (fParentId >= 0) res.emplace_back(fParentId);
217 res.emplace_back(fNodeId);
218 }
219 return res;
220 }
221
222};
223} // namespace Experimental
224} // namespace ROOT
225
226
227using namespace ROOT::Experimental;
228
229using namespace std::string_literals;
230
231/////////////////////////////////////////////////////////////////////
232/// Pack matrix into vector, which can be send to client
233/// Following sizes can be used for vector:
234/// 0 - Identity matrix
235/// 3 - Translation
236/// 4 - Scale (last element always 1)
237/// 9 - Rotation
238/// 16 - Full size
239
240void RGeomDescription::PackMatrix(std::vector<float> &vect, TGeoMatrix *matr)
241{
242 vect.clear();
243
244 if (!matr || matr->IsIdentity())
245 return;
246
247 auto trans = matr->GetTranslation();
248 auto scale = matr->GetScale();
249 auto rotate = matr->GetRotationMatrix();
250
251 bool is_translate = matr->IsA() == TGeoTranslation::Class(),
252 is_scale = matr->IsA() == TGeoScale::Class(),
253 is_rotate = matr->IsA() == TGeoRotation::Class();
254
255 if (!is_translate && !is_scale && !is_rotate) {
256 // check if trivial matrix
257
258 auto test = [](double val, double chk) { return (val==chk) || (TMath::Abs(val-chk) < 1e-20); };
259
260 bool no_scale = test(scale[0],1) && test(scale[1],1) && test(scale[2],1);
261 bool no_trans = test(trans[0],0) && test(trans[1],0) && test(trans[2],0);
262 bool no_rotate = test(rotate[0],1) && test(rotate[1],0) && test(rotate[2],0) &&
263 test(rotate[3],0) && test(rotate[4],1) && test(rotate[5],0) &&
264 test(rotate[6],0) && test(rotate[7],0) && test(rotate[8],1);
265
266 if (no_scale && no_trans && no_rotate)
267 return;
268
269 if (no_scale && no_trans && !no_rotate) {
270 is_rotate = true;
271 } else if (no_scale && !no_trans && no_rotate) {
272 is_translate = true;
273 } else if (!no_scale && no_trans && no_rotate) {
274 is_scale = true;
275 }
276 }
277
278 if (is_translate) {
279 vect.resize(3);
280 vect[0] = trans[0];
281 vect[1] = trans[1];
282 vect[2] = trans[2];
283 return;
284 }
285
286 if (is_scale) {
287 vect.resize(4);
288 vect[0] = scale[0];
289 vect[1] = scale[1];
290 vect[2] = scale[2];
291 vect[3] = 1;
292 return;
293 }
294
295 if (is_rotate) {
296 vect.resize(9);
297 for (int n=0;n<9;++n)
298 vect[n] = rotate[n];
299 return;
300 }
301
302 vect.resize(16);
303 vect[0] = rotate[0]; vect[4] = rotate[1]; vect[8] = rotate[2]; vect[12] = trans[0];
304 vect[1] = rotate[3]; vect[5] = rotate[4]; vect[9] = rotate[5]; vect[13] = trans[1];
305 vect[2] = rotate[6]; vect[6] = rotate[7]; vect[10] = rotate[8]; vect[14] = trans[2];
306 vect[3] = 0; vect[7] = 0; vect[11] = 0; vect[15] = 1;
307}
308
309/////////////////////////////////////////////////////////////////////
310/// Collect information about geometry hierarchy into flat list
311/// like it done in JSROOT ClonedNodes.createClones
312
313void RGeomDescription::Build(TGeoManager *mgr, const std::string &volname)
314{
316 if (!mgr) return;
317
318 auto topnode = mgr->GetTopNode();
319
320 if (!volname.empty()) {
321 auto vol = mgr->GetVolume(volname.c_str());
322 if (vol) {
323 TGeoNode *node;
324 TGeoIterator next(mgr->GetTopVolume());
325 while ((node = next()) != nullptr) {
326 if (node->GetVolume() == vol) break;
327 }
328 if (node) topnode = node;
329 }
330 }
331
332 // by top node visibility always enabled and harm logic
333 // later visibility can be controlled by other means
334 // mgr->GetTopNode()->GetVolume()->SetVisibility(kFALSE);
335
336 int maxnodes = mgr->GetMaxVisNodes();
337
339 SetVisLevel(mgr->GetVisLevel());
340 SetMaxVisNodes(maxnodes);
341 SetMaxVisFaces( (maxnodes > 5000 ? 5000 : (maxnodes < 1000 ? 1000 : maxnodes)) * 100);
342
343 BuildDescription(topnode, topnode->GetVolume());
344}
345
346/////////////////////////////////////////////////////////////////////
347/// Collect information about geometry from single volume
348/// like it done in JSROOT ClonedNodes.createClones
349
351{
353
354 if (!vol) return;
355
356 fDrawVolume = vol;
357
359}
360
361/////////////////////////////////////////////////////////////////////
362/// Clear geometry description
363
365{
366 fDesc.clear();
367 fNodes.clear();
368 fSortMap.clear();
370 fDrawIdCut = 0;
371 fDrawVolume = nullptr;
372}
373
374/////////////////////////////////////////////////////////////////////
375/// Build geometry description
376
378{
379 // vector to remember numbers
380 std::vector<int> numbers;
381 int offset = 1000000000;
382
383 // try to build flat list of all nodes
384 TGeoNode *snode = topnode;
385 TGeoIterator iter(topvolume);
386 do {
387 if (!snode) {
388 numbers.emplace_back(offset);
389 fNodes.emplace_back(nullptr);
390 } else if (snode->GetNumber() >= offset) {
391 // artificial offset already applied, used as identifier
392 iter.Skip(); // no need to look inside
393 } else {
394 numbers.emplace_back(snode->GetNumber());
395 snode->SetNumber(offset + fNodes.size()); // use id with shift 1e9
396 fNodes.emplace_back(snode);
397 }
398 } while ((snode = iter()) != nullptr);
399
400 fDesc.reserve(fNodes.size());
401 numbers.reserve(fNodes.size());
402 fSortMap.reserve(fNodes.size());
403
404 // array for sorting
405 std::vector<RGeomNode *> sortarr;
406 sortarr.reserve(fNodes.size());
407
408 // create vector of desc and childs
409 int cnt = 0;
410 for (auto node : fNodes) {
411
412 fDesc.emplace_back(node ? node->GetNumber() - offset : 0);
413 TGeoVolume *vol = node ? node->GetVolume() : topvolume;
414
415 auto &desc = fDesc[cnt++];
416
417 sortarr.emplace_back(&desc);
418
419 desc.name = node ? node->GetName() : vol->GetName();
420
421 auto shape = dynamic_cast<TGeoBBox *>(vol->GetShape());
422 if (shape) {
423 desc.vol = shape->GetDX() * shape->GetDY() * shape->GetDZ();
424 desc.nfaces = 12; // TODO: get better value for each shape - excluding composite
425 }
426
427 CopyMaterialProperties(vol, desc);
428
429 auto chlds = node ? node->GetNodes() : vol->GetNodes();
430
431 PackMatrix(desc.matr, node ? node->GetMatrix() : nullptr);
432
433 if (chlds)
434 for (int n = 0; n <= chlds->GetLast(); ++n) {
435 auto chld = dynamic_cast<TGeoNode *> (chlds->At(n));
436 desc.chlds.emplace_back(chld->GetNumber() - offset);
437 }
438 }
439
440 // recover numbers
441 cnt = 0;
442 for (auto node : fNodes) {
443 auto number = numbers[cnt++];
444 if (node) node->SetNumber(number);
445 }
446
447 // sort in volume descent order
448 std::sort(sortarr.begin(), sortarr.end(), [](RGeomNode *a, RGeomNode * b) { return a->vol > b->vol; });
449
450 cnt = 0;
451 for (auto &elem: sortarr) {
452 fSortMap.emplace_back(elem->id);
453 elem->sortid = cnt++; // keep place in sorted array to correctly apply cut
454 }
455
456 MarkVisible(); // set visibility flags
457
459}
460
461/////////////////////////////////////////////////////////////////////
462/// Get volume for specified nodeid
463/// If specific volume was configured, it will be returned for nodeid==0
464
466{
467 auto node = fNodes[nodeid];
468 if (node) return node->GetVolume();
469 return nodeid == 0 ? fDrawVolume : nullptr;
470}
471
472/////////////////////////////////////////////////////////////////////
473/// Set visibility flag for each nodes
474
476{
477 int res = 0;
478 for (int nodeid = 0; nodeid < (int) fNodes.size(); nodeid++) {
479
480 auto node = fNodes[nodeid];
481 auto vol = GetVolume(nodeid);
482 auto &desc = fDesc[nodeid];
483 desc.vis = 0;
484 desc.nochlds = false;
485
486 if (on_screen) {
487 if (!node || node->IsOnScreen())
488 desc.vis = 99;
489 } else {
490 if (vol->IsVisible() && !vol->TestAttBit(TGeoAtt::kVisNone))
491 desc.vis = 99;
492
493 if (node && !node->IsVisDaughters())
494 desc.nochlds = true;
495
496 if ((desc.vis > 0) && (desc.chlds.size() > 0) && !desc.nochlds)
497 desc.vis = 1;
498 }
499
500 if (desc.IsVisible() && desc.CanDisplay()) res++;
501 }
502
503 return res;
504}
505
506/////////////////////////////////////////////////////////////////////
507/// Count total number of visible childs under each node
508
510{
511 for (auto &node : fDesc)
512 node.idshift = -1;
513
514 using ScanFunc_t = std::function<int(RGeomNode &)>;
515
516 ScanFunc_t scan_func = [&, this](RGeomNode &node) {
517 if (node.idshift < 0) {
518 node.idshift = 0;
519 for(auto id : node.chlds)
520 node.idshift += scan_func(fDesc[id]);
521 }
522
523 return node.idshift + 1;
524 };
525
526 if (fDesc.size() > 0)
527 scan_func(fDesc[0]);
528}
529
530/////////////////////////////////////////////////////////////////////
531/// Iterate over all nodes and call function for visible
532
533int RGeomDescription::ScanNodes(bool only_visible, int maxlvl, RGeomScanFunc_t func)
534{
535 if (fDesc.empty()) return 0;
536
537 std::vector<int> stack;
538 stack.reserve(25); // reserve enough space for most use-cases
539 int counter = 0;
540
541 using ScanFunc_t = std::function<int(int, int)>;
542
543 ScanFunc_t scan_func = [&, this](int nodeid, int lvl) {
544 auto &desc = fDesc[nodeid];
545 int res = 0;
546
547 if (desc.nochlds && (lvl > 0)) lvl = 0;
548
549 // same logic as in JSROOT ClonedNodes.scanVisible
550 bool is_visible = (lvl >= 0) && (desc.vis > lvl) && desc.CanDisplay();
551
552 if (is_visible || !only_visible)
553 if (func(desc, stack, is_visible, counter))
554 res++;
555
556 counter++; // count sequence id of current position in scan, will be used later for merging drawing lists
557
558 if ((desc.chlds.size() > 0) && ((lvl > 0) || !only_visible)) {
559 auto pos = stack.size();
560 stack.emplace_back(0);
561 for (unsigned k = 0; k < desc.chlds.size(); ++k) {
562 stack[pos] = k; // stack provides index in list of childs
563 res += scan_func(desc.chlds[k], lvl - 1);
564 }
565 stack.pop_back();
566 } else {
567 counter += desc.idshift;
568 }
569
570 return res;
571 };
572
573 if (!maxlvl && (GetVisLevel() > 0)) maxlvl = GetVisLevel();
574 if (!maxlvl) maxlvl = 4;
575 if (maxlvl > 97) maxlvl = 97; // check while vis property of node is 99 normally
576
577 return scan_func(0, maxlvl);
578}
579
580/////////////////////////////////////////////////////////////////////
581/// Collect nodes which are used in visibles
582
584{
585 // TODO: for now reset all flags, later can be kept longer
586 for (auto &node : fDesc)
587 node.useflag = false;
588
589 drawing.cfg = &fCfg;
590
591 drawing.numnodes = fDesc.size();
592
593 for (auto &item : drawing.visibles) {
594 int nodeid = 0;
595 for (auto &chindx : item.stack) {
596 auto &node = fDesc[nodeid];
597 if (!node.useflag) {
598 node.useflag = true;
599 drawing.nodes.emplace_back(&node);
600 }
601 if (chindx >= (int)node.chlds.size())
602 break;
603 nodeid = node.chlds[chindx];
604 }
605
606 auto &node = fDesc[nodeid];
607 if (!node.useflag) {
608 node.useflag = true;
609 drawing.nodes.emplace_back(&node);
610 }
611 }
612
613 // printf("SELECT NODES %d\n", (int) drawing.nodes.size());
614}
615
616/////////////////////////////////////////////////////////////////////
617/// Find description object for requested shape
618/// If not exists - will be created
619
620std::string RGeomDescription::ProcessBrowserRequest(const std::string &msg)
621{
622 std::string res;
623
624 auto request = TBufferJSON::FromJSON<RBrowserRequest>(msg);
625
626 if (msg.empty()) {
627 request = std::make_unique<RBrowserRequest>();
628 request->first = 0;
629 request->number = 100;
630 }
631
632 if (!request)
633 return res;
634
635 if (request->path.empty() && (request->first == 0) && (GetNumNodes() < (IsPreferredOffline() ? 1000000 : 1000))) {
636
637 std::vector<RGeomNodeBase *> vect(fDesc.size(), nullptr);
638
639 int cnt = 0;
640 for (auto &item : fDesc)
641 vect[cnt++]= &item;
642
643 res = "DESCR:"s + TBufferJSON::ToJSON(&vect,GetJsonComp()).Data();
644 } else {
645 std::vector<Browsable::RItem> temp_nodes;
646 bool toplevel = request->path.empty();
647
648 // create temporary object for the short time
649 RBrowserReply reply;
650 reply.path = request->path;
651 reply.first = request->first;
652
653 RGeomBrowserIter iter(*this);
654 if (iter.Navigate(request->path)) {
655
656 reply.nchilds = iter.NumChilds();
657 // scan childs of selected nodes
658 if (iter.Enter()) {
659
660 while ((request->first > 0) && iter.Next()) {
661 request->first--;
662 }
663
664 while (iter.IsValid() && (request->number > 0)) {
665 temp_nodes.emplace_back(iter.GetName(), iter.NumChilds());
666 if (toplevel) temp_nodes.back().SetExpanded(true);
667 request->number--;
668 if (!iter.Next()) break;
669 }
670 }
671 }
672
673 for (auto &n : temp_nodes)
674 reply.nodes.emplace_back(&n);
675
676 res = "BREPL:"s + TBufferJSON::ToJSON(&reply, GetJsonComp()).Data();
677 }
678
679 return res;
680}
681
682/////////////////////////////////////////////////////////////////////
683/// Find description object for requested shape
684/// If not exists - will be created
685
687{
688 for (auto &descr : fShapes)
689 if (descr.fShape == shape)
690 return descr;
691
692 fShapes.emplace_back(shape);
693 auto &elem = fShapes.back();
694 elem.id = fShapes.size() - 1;
695 return elem;
696}
697
698
699////////////////////////////////////////////////////////////////////////
700/// Function produces mesh for provided shape, applying matrix to the result
701
702std::unique_ptr<RootCsg::TBaseMesh> MakeGeoMesh(TGeoMatrix *matr, TGeoShape *shape)
703{
704 TGeoCompositeShape *comp = dynamic_cast<TGeoCompositeShape *> (shape);
705
706 std::unique_ptr<RootCsg::TBaseMesh> res;
707
708 if (!comp) {
709 std::unique_ptr<TBuffer3D> b3d(shape->MakeBuffer3D());
710
711 if (matr) {
712 Double_t *v = b3d->fPnts;
713 Double_t buf[3];
714 for (UInt_t i = 0; i < b3d->NbPnts(); ++i) {
715 buf[0] = v[i*3];
716 buf[1] = v[i*3+1];
717 buf[2] = v[i*3+2];
718 matr->LocalToMaster(buf, &v[i*3]);
719 }
720 }
721
722 res.reset(RootCsg::ConvertToMesh(*b3d.get()));
723 } else {
724 auto node = comp->GetBoolNode();
725
726 TGeoHMatrix mleft, mright;
727 if (matr) { mleft = *matr; mright = *matr; }
728
729 mleft.Multiply(node->GetLeftMatrix());
730 auto left = MakeGeoMesh(&mleft, node->GetLeftShape());
731
732 mright.Multiply(node->GetRightMatrix());
733 auto right = MakeGeoMesh(&mright, node->GetRightShape());
734
735 if (node->IsA() == TGeoUnion::Class()) res.reset(RootCsg::BuildUnion(left.get(), right.get()));
736 if (node->IsA() == TGeoIntersection::Class()) res.reset(RootCsg::BuildIntersection(left.get(), right.get()));
737 if (node->IsA() == TGeoSubtraction::Class()) res.reset(RootCsg::BuildDifference(left.get(), right.get()));
738 }
739
740 return res;
741}
742
743
744/////////////////////////////////////////////////////////////////////
745/// Find description object and create render information
746
748{
749 auto &elem = FindShapeDescr(shape);
750
751 if (elem.nfaces == 0) {
752
753 int boundary = 3; //
754 if (shape->IsComposite()) {
755 // composite is most complex for client, therefore by default build on server
756 boundary = 1;
757 } else if (!shape->IsCylType()) {
758 // simple box geometry is compact and can be delivered as raw
759 boundary = 2;
760 }
761
762 if (IsBuildShapes() < boundary) {
763 elem.nfaces = 1;
764 elem.fShapeInfo.shape = shape;
765 } else {
766
767 int old_nsegm = -1;
768 if (fCfg.nsegm > 0 && gGeoManager) {
769 old_nsegm = gGeoManager->GetNsegments();
771 }
772
773 auto mesh = MakeGeoMesh(nullptr, shape);
774
775 if (old_nsegm > 0 && gGeoManager)
776 gGeoManager->SetNsegments(old_nsegm);
777
778
779 Int_t num_vertices = mesh->NumberOfVertices(), num_polynoms = 0;
780
781 for (unsigned polyIndex = 0; polyIndex < mesh->NumberOfPolys(); ++polyIndex) {
782
783 auto size_of_polygon = mesh->SizeOfPoly(polyIndex);
784
785 if (size_of_polygon >= 3)
786 num_polynoms += (size_of_polygon - 2);
787 }
788
789 Int_t index_buffer_size = num_polynoms * 3, // triangle indexes
790 vertex_buffer_size = num_vertices * 3; // X,Y,Z array
791
792 elem.nfaces = num_polynoms;
793
794 std::vector<float> vertices(vertex_buffer_size);
795
796 for (Int_t i = 0; i < num_vertices; ++i) {
797 auto v = mesh->GetVertex(i);
798 vertices[i*3] = v[0];
799 vertices[i*3+1] = v[1];
800 vertices[i*3+2] = v[2];
801 }
802
803 elem.fRawInfo.raw.resize(vertices.size() * sizeof(float));
804
805 memcpy(reinterpret_cast<char *>(elem.fRawInfo.raw.data()), vertices.data(), vertices.size() * sizeof(float));
806
807 auto &indexes = elem.fRawInfo.idx;
808
809 indexes.resize(index_buffer_size);
810 int pos = 0;
811
812 for (unsigned polyIndex = 0; polyIndex < mesh->NumberOfPolys(); ++polyIndex) {
813 auto size_of_polygon = mesh->SizeOfPoly(polyIndex);
814
815 // add first triangle
816 if (size_of_polygon >= 3)
817 for (int i = 0; i < 3; ++i)
818 indexes[pos++] = mesh->GetVertexIndex(polyIndex, i);
819
820 // add following triangles
821 if (size_of_polygon > 3)
822 for (unsigned vertex = 3; vertex < size_of_polygon; vertex++) {
823 indexes[pos++] = mesh->GetVertexIndex(polyIndex, 0);
824 indexes[pos++] = mesh->GetVertexIndex(polyIndex, vertex-1);
825 indexes[pos++] = mesh->GetVertexIndex(polyIndex, vertex);
826 }
827 }
828
829 }
830 }
831
832 return elem;
833}
834
835/////////////////////////////////////////////////////////////////////
836/// Copy material properties
837
839{
840 if (!volume) return;
841
842 TColor *col{nullptr};
843
844 if ((volume->GetFillColor() > 1) && (volume->GetLineColor() == 1))
845 col = gROOT->GetColor(volume->GetFillColor());
846 else if (volume->GetLineColor() >= 0)
847 col = gROOT->GetColor(volume->GetLineColor());
848
849 if (volume->GetMedium() && (volume->GetMedium() != TGeoVolume::DummyMedium()) && volume->GetMedium()->GetMaterial()) {
850 auto material = volume->GetMedium()->GetMaterial();
851
852 auto fillstyle = material->GetFillStyle();
853 if ((fillstyle>=3000) && (fillstyle<=3100)) node.opacity = (3100 - fillstyle) / 100.;
854 if (!col) col = gROOT->GetColor(material->GetFillColor());
855 }
856
857 if (col) {
858 node.color = std::to_string((int)(col->GetRed()*255)) + "," +
859 std::to_string((int)(col->GetGreen()*255)) + "," +
860 std::to_string((int)(col->GetBlue()*255));
861 if (node.opacity == 1.)
862 node.opacity = col->GetAlpha();
863 } else {
864 node.color.clear();
865 }
866}
867
868/////////////////////////////////////////////////////////////////////
869/// Reset shape info, which used to pack binary data
870
872{
873 for (auto &s: fShapes)
874 s.reset();
875}
876
877/////////////////////////////////////////////////////////////////////
878/// Produce JSON string which can be directly used with `build`
879/// function from JSROOT to create three.js model of configured geometry
880///
881/// Collect all information required to draw geometry on the client
882/// This includes list of each visible nodes, meshes and matrixes
883///
884/// Example of usage:
885///
886/// void geom() {
887/// auto f = TFile::Open("file_name.root");
888/// auto vol = f->Get<TGeoVolume>("object_name");
889/// ROOT::Experimental::RGeomDescription desc;
890/// desc.Build(vol);
891/// std::ofstream fout("geom.json");
892/// fout << desc.ProduceJson();
893/// }
894///
895/// In JSROOT one loads data from JSON file and call `build` function to
896/// produce three.js model
897
899{
900 std::vector<int> viscnt(fDesc.size(), 0);
901
902 int level = GetVisLevel();
903
904 // first count how many times each individual node appears
905 int numnodes = ScanNodes(true, level, [&viscnt](RGeomNode &node, std::vector<int> &, bool, int) {
906 viscnt[node.id]++;
907 return true;
908 });
909
910 if (GetMaxVisNodes() > 0) {
911 while ((numnodes > GetMaxVisNodes()) && (level > 1)) {
912 level--;
913 viscnt.assign(viscnt.size(), 0);
914 numnodes = ScanNodes(true, level, [&viscnt](RGeomNode &node, std::vector<int> &, bool, int) {
915 viscnt[node.id]++;
916 return true;
917 });
918 }
919 }
920
921 fActualLevel = level;
922 fDrawIdCut = 0;
923
924 int totalnumfaces = 0, totalnumnodes = 0;
925
926 //for (auto &node : fDesc)
927 // node.SetDisplayed(false);
928
929 // build all shapes in volume decreasing order
930 for (auto &sid: fSortMap) {
931 fDrawIdCut++; //
932 auto &desc = fDesc[sid];
933
934 if ((viscnt[sid] <= 0) || (desc.vol <= 0)) continue;
935
936 auto shape = GetVolume(sid)->GetShape();
937 if (!shape) continue;
938
939 // now we need to create TEveGeoPolyShape, which can provide all rendering data
940 auto &shape_descr = MakeShapeDescr(shape);
941
942 // should not happen, but just in case
943 if (shape_descr.nfaces <= 0) {
944 R__LOG_ERROR(RGeomLog()) << "No faces for the shape " << shape->GetName() << " class " << shape->ClassName();
945 continue;
946 }
947
948 // check how many faces are created
949 totalnumfaces += shape_descr.nfaces * viscnt[sid];
950 if ((GetMaxVisFaces() > 0) && (totalnumfaces > GetMaxVisFaces())) break;
951
952 // also avoid too many nodes
953 totalnumnodes += viscnt[sid];
954 if ((GetMaxVisNodes() > 0) && (totalnumnodes > GetMaxVisNodes())) break;
955
956 // desc.SetDisplayed(true);
957 }
958
959 // finally we should create data for streaming to the client
960 // it includes list of visible nodes and rawdata
961
962 RGeomDrawing drawing;
964 bool has_shape = false;
965
966 ScanNodes(true, level, [&, this](RGeomNode &node, std::vector<int> &stack, bool, int seqid) {
967 if (node.sortid < fDrawIdCut) {
968 drawing.visibles.emplace_back(node.id, seqid, stack);
969
970 auto &item = drawing.visibles.back();
971 item.color = node.color;
972 item.opacity = node.opacity;
973
974 auto volume = GetVolume(node.id);
975
976 auto &sd = MakeShapeDescr(volume->GetShape());
977
978 item.ri = sd.rndr_info();
979 if (sd.has_shape()) has_shape = true;
980 }
981 return true;
982 });
983
984 CollectNodes(drawing);
985
986 return MakeDrawingJson(drawing, has_shape);
987}
988
989/////////////////////////////////////////////////////////////////////
990/// Collect all information required to draw geometry on the client
991/// This includes list of each visible nodes, meshes and matrixes
992
994{
995 fDrawJson = "GDRAW:"s + ProduceJson();
996}
997
998/////////////////////////////////////////////////////////////////////
999/// Clear raw data. Will be rebuild when next connection will be established
1000
1002{
1003 fDrawJson.clear();
1004}
1005
1006/////////////////////////////////////////////////////////////////////
1007/// return true when node used in main geometry drawing and does not have childs
1008/// for such nodes one could provide optimize toggling of visibility flags
1009
1011{
1012 if ((nodeid < 0) || (nodeid >= (int)fDesc.size()))
1013 return false;
1014
1015 auto &desc = fDesc[nodeid];
1016
1017 return (desc.sortid < fDrawIdCut) && desc.IsVisible() && desc.CanDisplay() && (desc.chlds.size()==0);
1018}
1019
1020
1021/////////////////////////////////////////////////////////////////////
1022/// Search visible nodes for provided name
1023/// If number of found elements less than 100, create description and shapes for them
1024/// Returns number of match elements
1025
1026int RGeomDescription::SearchVisibles(const std::string &find, std::string &hjson, std::string &json)
1027{
1028 hjson.clear();
1029 json.clear();
1030
1031 if (find.empty()) {
1032 hjson = "FOUND:RESET";
1033 return 0;
1034 }
1035
1036 std::vector<int> nodescnt(fDesc.size(), 0), viscnt(fDesc.size(), 0);
1037
1038 int nmatches = 0;
1039
1040 auto match_func = [&find](RGeomNode &node) {
1041 return (node.vol > 0) && (node.name.compare(0, find.length(), find) == 0);
1042 };
1043
1044 // first count how many times each individual node appears
1045 ScanNodes(false, 0, [&nodescnt,&viscnt,&match_func,&nmatches](RGeomNode &node, std::vector<int> &, bool is_vis, int) {
1046
1047 if (match_func(node)) {
1048 nmatches++;
1049 nodescnt[node.id]++;
1050 if (is_vis) viscnt[node.id]++;
1051 };
1052 return true;
1053 });
1054
1055 // do not send too much data, limit could be made configurable later
1056 if (nmatches==0) {
1057 hjson = "FOUND:NO";
1058 return nmatches;
1059 }
1060
1061 if ((GetMaxVisNodes() > 0) && (nmatches > 10 * GetMaxVisNodes())) {
1062 hjson = "FOUND:Too many " + std::to_string(nmatches);
1063 return nmatches;
1064 }
1065
1066 // now build all necessary shapes and check number of faces - not too many
1067
1068 int totalnumfaces = 0, totalnumnodes = 0, scnt = 0;
1069 bool send_rawdata = true;
1070
1071 // build all shapes in volume decreasing order
1072 for (auto &sid: fSortMap) {
1073 if (scnt++ < fDrawIdCut) continue; // no need to send most significant shapes
1074
1075 if (viscnt[sid] == 0) continue; // this node is not used at all
1076
1077 auto &desc = fDesc[sid];
1078 if ((viscnt[sid] <= 0) && (desc.vol <= 0)) continue;
1079
1080 auto shape = GetVolume(sid)->GetShape();
1081 if (!shape) continue;
1082
1083 // create shape raw data
1084 auto &shape_descr = MakeShapeDescr(shape);
1085
1086 // should not happen, but just in case
1087 if (shape_descr.nfaces <= 0) {
1088 R__LOG_ERROR(RGeomLog()) << "No faces for the shape " << shape->GetName() << " class " << shape->ClassName();
1089 continue;
1090 }
1091
1092 // check how many faces are created
1093 totalnumfaces += shape_descr.nfaces * viscnt[sid];
1094 if ((GetMaxVisFaces() > 0) && (totalnumfaces > GetMaxVisFaces())) { send_rawdata = false; break; }
1095
1096 // also avoid too many nodes
1097 totalnumnodes += viscnt[sid];
1098 if ((GetMaxVisNodes() > 0) && (totalnumnodes > GetMaxVisNodes())) { send_rawdata = false; break; }
1099 }
1100
1101 // only for debug purposes - remove later
1102 // send_rawdata = false;
1103
1104 // finally we should create data for streaming to the client
1105 // it includes list of visible nodes and rawdata (if there is enough space)
1106
1107
1108 std::vector<RGeomNodeBase> found_desc; ///<! hierarchy of nodes, used for search
1109 std::vector<int> found_map(fDesc.size(), -1); ///<! mapping between nodeid - > foundid
1110
1111 // these are only selected nodes to produce hierarchy
1112
1113 found_desc.emplace_back(0);
1114 found_desc[0].vis = fDesc[0].vis;
1115 found_desc[0].name = fDesc[0].name;
1116 found_desc[0].color = fDesc[0].color;
1117 found_map[0] = 0;
1118
1120
1121 RGeomDrawing drawing;
1122 bool has_shape = true;
1123
1124 ScanNodes(false, 0, [&, this](RGeomNode &node, std::vector<int> &stack, bool is_vis, int seqid) {
1125 // select only nodes which should match
1126 if (!match_func(node))
1127 return true;
1128
1129 // add entries into hierarchy of found elements
1130 int prntid = 0;
1131 for (auto &s : stack) {
1132 int chldid = fDesc[prntid].chlds[s];
1133 if (found_map[chldid] <= 0) {
1134 int newid = found_desc.size();
1135 found_desc.emplace_back(newid); // potentially original id can be used here
1136 found_map[chldid] = newid; // re-map into reduced hierarchy
1137
1138 found_desc.back().vis = fDesc[chldid].vis;
1139 found_desc.back().name = fDesc[chldid].name;
1140 found_desc.back().color = fDesc[chldid].color;
1141 }
1142
1143 auto pid = found_map[prntid];
1144 auto cid = found_map[chldid];
1145
1146 // now add entry into childs lists
1147 auto &pchlds = found_desc[pid].chlds;
1148 if (std::find(pchlds.begin(), pchlds.end(), cid) == pchlds.end())
1149 pchlds.emplace_back(cid);
1150
1151 prntid = chldid;
1152 }
1153
1154 // no need to add visibles
1155 if (!is_vis) return true;
1156
1157 drawing.visibles.emplace_back(node.id, seqid, stack);
1158
1159 // no need to transfer shape if it provided with main drawing list
1160 // also no binary will be transported when too many matches are there
1161 if (!send_rawdata || (node.sortid < fDrawIdCut)) {
1162 // do not include render data
1163 return true;
1164 }
1165
1166 auto &item = drawing.visibles.back();
1167 auto volume = GetVolume(node.id);
1168
1169 item.color = node.color;
1170 item.opacity = node.opacity;
1171
1172 auto &sd = MakeShapeDescr(volume->GetShape());
1173
1174 item.ri = sd.rndr_info();
1175 if (sd.has_shape()) has_shape = true;
1176 return true;
1177 });
1178
1179 hjson = "FESCR:"s + TBufferJSON::ToJSON(&found_desc, GetJsonComp()).Data();
1180
1181 CollectNodes(drawing);
1182
1183 json = "FDRAW:"s + MakeDrawingJson(drawing, has_shape);
1184
1185 return nmatches;
1186}
1187
1188/////////////////////////////////////////////////////////////////////////////////
1189/// Returns nodeid for given stack array, returns -1 in case of failure
1190
1191int RGeomDescription::FindNodeId(const std::vector<int> &stack)
1192{
1193 int nodeid = 0;
1194
1195 for (auto &chindx: stack) {
1196 auto &node = fDesc[nodeid];
1197 if (chindx >= (int) node.chlds.size()) return -1;
1198 nodeid = node.chlds[chindx];
1199 }
1200
1201 return nodeid;
1202}
1203
1204/////////////////////////////////////////////////////////////////////////////////
1205/// Creates stack for given array of ids, first element always should be 0
1206
1207std::vector<int> RGeomDescription::MakeStackByIds(const std::vector<int> &ids)
1208{
1209 std::vector<int> stack;
1210
1211 if (ids[0] != 0) {
1212 printf("Wrong first id\n");
1213 return stack;
1214 }
1215
1216 int nodeid = 0;
1217
1218 for (unsigned k = 1; k < ids.size(); ++k) {
1219
1220 int prntid = nodeid;
1221 nodeid = ids[k];
1222
1223 if (nodeid >= (int) fDesc.size()) {
1224 printf("Wrong node id %d\n", nodeid);
1225 stack.clear();
1226 return stack;
1227 }
1228 auto &chlds = fDesc[prntid].chlds;
1229 auto pos = std::find(chlds.begin(), chlds.end(), nodeid);
1230 if (pos == chlds.end()) {
1231 printf("Wrong id %d not a child of %d - fail to find stack num %d\n", nodeid, prntid, (int) chlds.size());
1232 stack.clear();
1233 return stack;
1234 }
1235
1236 stack.emplace_back(std::distance(chlds.begin(), pos));
1237 }
1238
1239 return stack;
1240}
1241
1242/////////////////////////////////////////////////////////////////////////////////
1243/// Produce stack based on string path
1244/// Used to highlight geo volumes by browser hover event
1245
1246std::vector<int> RGeomDescription::MakeStackByPath(const std::vector<std::string> &path)
1247{
1248 std::vector<int> res;
1249
1250 RGeomBrowserIter iter(*this);
1251
1252 if (iter.Navigate(path))
1253 res = MakeStackByIds(iter.CurrentIds());
1254
1255 return res;
1256}
1257
1258/////////////////////////////////////////////////////////////////////////////////
1259/// Produce list of node ids for given stack
1260/// If found nodes preselected - use their ids
1261
1262std::vector<int> RGeomDescription::MakeIdsByStack(const std::vector<int> &stack)
1263{
1264 std::vector<int> ids;
1265
1266 ids.emplace_back(0);
1267 int nodeid = 0;
1268 bool failure = false;
1269
1270 for (auto s : stack) {
1271 auto &chlds = fDesc[nodeid].chlds;
1272 if (s >= (int) chlds.size()) { failure = true; break; }
1273
1274 ids.emplace_back(chlds[s]);
1275
1276 nodeid = chlds[s];
1277 }
1278
1279 if (failure) {
1280 printf("Fail to convert stack into list of nodes\n");
1281 ids.clear();
1282 }
1283
1284 return ids;
1285}
1286
1287/////////////////////////////////////////////////////////////////////////////////
1288/// Returns path string for provided stack
1289
1290std::vector<std::string> RGeomDescription::MakePathByStack(const std::vector<int> &stack)
1291{
1292 std::vector<std::string> path;
1293
1294 auto ids = MakeIdsByStack(stack);
1295 for (auto &id : ids)
1296 path.emplace_back(fDesc[id].name);
1297
1298 return path;
1299}
1300
1301
1302/////////////////////////////////////////////////////////////////////////////////
1303/// Return string with only part of nodes description which were modified
1304/// Checks also volume
1305
1307{
1308 std::vector<RGeomNodeBase *> nodes;
1309 auto vol = GetVolume(nodeid);
1310
1311 // we take not only single node, but all there same volume is referenced
1312 // nodes.push_back(&fDesc[nodeid]);
1313
1314 int id = 0;
1315 for (auto &desc : fDesc)
1316 if (GetVolume(id++) == vol)
1317 nodes.emplace_back(&desc);
1318
1319 return "MODIF:"s + TBufferJSON::ToJSON(&nodes, GetJsonComp()).Data();
1320}
1321
1322
1323/////////////////////////////////////////////////////////////////////////////////
1324/// Produce shape rendering data for given stack
1325/// All nodes, which are referencing same shape will be transferred
1326/// Returns true if new render information provided
1327
1328bool RGeomDescription::ProduceDrawingFor(int nodeid, std::string &json, bool check_volume)
1329{
1330 // only this shape is interesting
1331
1332 TGeoVolume *vol = (nodeid < 0) ? nullptr : GetVolume(nodeid);
1333
1334 if (!vol || !vol->GetShape()) {
1335 json.append("NO");
1336 return false;
1337 }
1338
1339 RGeomDrawing drawing;
1340
1341 ScanNodes(true, 0, [&, this](RGeomNode &node, std::vector<int> &stack, bool, int seq_id) {
1342 // select only nodes which reference same shape
1343
1344 if (check_volume) {
1345 if (GetVolume(node.id) != vol) return true;
1346 } else {
1347 if (node.id != nodeid) return true;
1348 }
1349
1350 drawing.visibles.emplace_back(node.id, seq_id, stack);
1351
1352 auto &item = drawing.visibles.back();
1353
1354 item.color = node.color;
1355 item.opacity = node.opacity;
1356 return true;
1357 });
1358
1359 // no any visible nodes were done
1360 if (drawing.visibles.size()==0) {
1361 json.append("NO");
1362 return false;
1363 }
1364
1366
1367 bool has_shape = false, has_raw = false;
1368
1369 auto &sd = MakeShapeDescr(vol->GetShape());
1370
1371 // assign shape data
1372 for (auto &item : drawing.visibles) {
1373 item.ri = sd.rndr_info();
1374 if (sd.has_shape()) has_shape = true;
1375 if (sd.has_raw()) has_raw = true;
1376 }
1377
1378 CollectNodes(drawing);
1379
1380 json.append(MakeDrawingJson(drawing, has_shape));
1381
1382 return has_raw || has_shape;
1383}
1384
1385/////////////////////////////////////////////////////////////////////////////////
1386/// Produce JSON for the drawing
1387/// If TGeoShape appears in the drawing, one has to keep typeinfo
1388/// But in this case one can exclude several classes which are not interesting,
1389/// but appears very often
1390
1391std::string RGeomDescription::MakeDrawingJson(RGeomDrawing &drawing, bool has_shapes)
1392{
1393 int comp = GetJsonComp();
1394
1395 if (!has_shapes || (comp < TBufferJSON::kSkipTypeInfo))
1396 return TBufferJSON::ToJSON(&drawing, comp).Data();
1397
1398 comp = comp % TBufferJSON::kSkipTypeInfo; // no typeinfo skipping
1399
1401 json.SetCompact(comp);
1402 json.SetSkipClassInfo(TClass::GetClass<RGeomDrawing>());
1403 json.SetSkipClassInfo(TClass::GetClass<RGeomNode>());
1404 json.SetSkipClassInfo(TClass::GetClass<RGeomVisible>());
1405 json.SetSkipClassInfo(TClass::GetClass<RGeomShapeRenderInfo>());
1406 json.SetSkipClassInfo(TClass::GetClass<RGeomRawRenderInfo>());
1407
1408 return json.StoreObject(&drawing, TClass::GetClass<RGeomDrawing>()).Data();
1409}
1410
1411/////////////////////////////////////////////////////////////////////////////////
1412/// Change visibility for specified element
1413/// Returns true if changes was performed
1414
1415bool RGeomDescription::ChangeNodeVisibility(int nodeid, bool selected)
1416{
1417 auto &dnode = fDesc[nodeid];
1418
1419 auto vol = GetVolume(nodeid);
1420
1421 // nothing changed
1422 if (vol->IsVisible() == selected)
1423 return false;
1424
1425 dnode.vis = selected ? 99 : 0;
1426 vol->SetVisibility(selected);
1427 if (dnode.chlds.size() > 0) {
1428 if (selected) dnode.vis = 1; // visibility disabled when any child
1429 vol->SetVisDaughters(selected);
1430 }
1431
1432 int id = 0;
1433 for (auto &desc: fDesc)
1434 if (GetVolume(id++) == vol)
1435 desc.vis = dnode.vis;
1436
1437 ClearDrawData(); // after change raw data is no longer valid
1438
1439 return true;
1440}
1441
1442/////////////////////////////////////////////////////////////////////////////////
1443/// Change visibility for specified element
1444/// Returns true if changes was performed
1445
1446std::unique_ptr<RGeomNodeInfo> RGeomDescription::MakeNodeInfo(const std::vector<std::string> &path)
1447{
1448 std::unique_ptr<RGeomNodeInfo> res;
1449
1450 RGeomBrowserIter iter(*this);
1451
1452 if (iter.Navigate(path)) {
1453
1454 auto node = fNodes[iter.GetNodeId()];
1455
1456 auto &desc = fDesc[iter.GetNodeId()];
1457
1458 res = std::make_unique<RGeomNodeInfo>();
1459
1460 res->path = path;
1461 res->node_name = node ? node->GetName() : "node_name";
1462 res->node_type = node ? node->ClassName() : "no class";
1463
1464 auto vol = GetVolume(iter.GetNodeId());
1465
1466 TGeoShape *shape = vol ? vol->GetShape() : nullptr;
1467
1468 if (shape) {
1469 res->shape_name = shape->GetName();
1470 res->shape_type = shape->ClassName();
1471 }
1472
1473 if (shape && desc.CanDisplay()) {
1474
1475 auto &shape_descr = MakeShapeDescr(shape);
1476
1477 res->ri = shape_descr.rndr_info(); // temporary pointer, can be used preserved for short time
1478 }
1479 }
1480
1481 return res;
1482}
1483
1484/////////////////////////////////////////////////////////////////////////////////
1485/// Change configuration by client
1486/// Returns true if any parameter was really changed
1487
1489{
1490 auto cfg = TBufferJSON::FromJSON<RGeomConfig>(json);
1491 if (!cfg) return false;
1492
1493 auto json1 = TBufferJSON::ToJSON(cfg.get());
1494 auto json2 = TBufferJSON::ToJSON(&fCfg);
1495
1496 if (json1 == json2)
1497 return false;
1498
1499 fCfg = *cfg; // use assign
1500
1501 ClearDrawData();
1502
1503 return true;
1504}
1505
1506
1507
1508
1509
nlohmann::json json
std::unique_ptr< RootCsg::TBaseMesh > MakeGeoMesh(TGeoMatrix *matr, TGeoShape *shape)
Function produces mesh for provided shape, applying matrix to the result.
#define R__LOG_ERROR(...)
Definition RLogger.hxx:362
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
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 Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
R__EXTERN TGeoManager * gGeoManager
#define gROOT
Definition TROOT.h:405
Reply on browser request.
std::vector< const Browsable::RItem * > nodes
list of pointers, no ownership!
std::vector< std::string > path
reply path
int nchilds
total number of childs in the node
int first
first node in returned list
Iterator of hierarchical geometry structures.
Definition RGeomData.cxx:49
RGeomBrowserIter(RGeomDescription &desc)
Definition RGeomData.cxx:61
bool Navigate(const std::string &path)
Navigate to specified path - path specified as string and should start with "/".
bool Navigate(const std::vector< std::string > &path)
Navigate to specified path
const std::string & GetName() const
Definition RGeomData.cxx:63
std::vector< int > CurrentIds() const
Returns array of ids to currently selected node.
int nsegm
number of segments for cylindrical shapes
std::vector< int > fSortMap
! nodes in order large -> smaller volume
std::string ProcessBrowserRequest(const std::string &req="")
Find description object for requested shape If not exists - will be created.
void PackMatrix(std::vector< float > &arr, TGeoMatrix *matr)
Pack matrix into vector, which can be send to client Following sizes can be used for vector: 0 - Iden...
bool ProduceDrawingFor(int nodeid, std::string &json, bool check_volume=false)
Produce shape rendering data for given stack All nodes, which are referencing same shape will be tran...
std::vector< int > MakeIdsByStack(const std::vector< int > &stack)
Produce list of node ids for given stack If found nodes preselected - use their ids.
int fActualLevel
! level can be reduced when selecting nodes
RGeomConfig fCfg
! configuration parameter editable from GUI
int MarkVisible(bool on_screen=false)
Set visibility flag for each nodes.
int GetJsonComp() const
Returns JSON compression level for data transfer.
bool IsPrincipalEndNode(int nodeid)
return true when node used in main geometry drawing and does not have childs for such nodes one could...
int GetMaxVisFaces() const
Returns maximal visible number of faces, ignored when non-positive.
TGeoVolume * fDrawVolume
! select volume independent from TGeoMaanger
int GetNumNodes() const
Number of unique nodes in the geometry.
void SetMaxVisNodes(int cnt)
Set maximal number of nodes which should be selected for drawing.
void SetVisLevel(int lvl=3)
Set maximal visible level.
int GetMaxVisNodes() const
Returns maximal visible number of nodes, ignored when non-positive.
void ClearDescription()
Clear geometry description.
std::vector< int > MakeStackByIds(const std::vector< int > &ids)
Creates stack for given array of ids, first element always should be 0.
int GetVisLevel() const
Returns maximal visible level.
std::vector< RGeomNode > fDesc
converted description, send to client
void SetMaxVisFaces(int cnt)
Set maximal number of faces which should be selected for drawing.
std::vector< ShapeDescr > fShapes
! shapes with created descriptions
std::string MakeDrawingJson(RGeomDrawing &drawing, bool has_shapes=false)
Produce JSON for the drawing If TGeoShape appears in the drawing, one has to keep typeinfo But in thi...
TGeoVolume * GetVolume(int nodeid)
Get volume for specified nodeid If specific volume was configured, it will be returned for nodeid==0.
std::unique_ptr< RGeomNodeInfo > MakeNodeInfo(const std::vector< std::string > &path)
Change visibility for specified element Returns true if changes was performed.
void ProduceDrawData()
Collect all information required to draw geometry on the client This includes list of each visible no...
void SetNSegments(int n=0)
Set number of segments for cylindrical shapes, if 0 - default value will be used.
void BuildDescription(TGeoNode *topnode, TGeoVolume *topvolume)
Build geometry description.
int SearchVisibles(const std::string &find, std::string &hjson, std::string &json)
Search visible nodes for provided name If number of found elements less than 100, create description ...
void ClearDrawData()
Clear raw data. Will be rebuild when next connection will be established.
std::vector< std::string > MakePathByStack(const std::vector< int > &stack)
Returns path string for provided stack.
bool ChangeConfiguration(const std::string &json)
Change configuration by client Returns true if any parameter was really changed.
void CopyMaterialProperties(TGeoVolume *vol, RGeomNode &node)
Copy material properties.
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.
int ScanNodes(bool only_visible, int maxlvl, RGeomScanFunc_t func)
Iterate over all nodes and call function for visible.
bool IsPreferredOffline() const
Is offline operations preferred.
std::vector< TGeoNode * > fNodes
! flat list of all nodes
ShapeDescr & MakeShapeDescr(TGeoShape *shape)
Find description object and create render information.
int IsBuildShapes() const
Returns true if binary 3D model build already by C++ server (default)
ShapeDescr & FindShapeDescr(TGeoShape *shape)
Find description object for requested shape If not exists - will be created.
std::string fDrawJson
! JSON with main nodes drawn by client
std::string ProduceJson()
Produce JSON string which can be directly used with build function from JSROOT to create three....
int FindNodeId(const std::vector< int > &stack)
Returns nodeid for given stack array, returns -1 in case of failure.
bool ChangeNodeVisibility(int nodeid, bool selected)
Change visibility for specified element Returns true if changes was performed.
void CollectNodes(RGeomDrawing &drawing)
Collect nodes which are used in visibles.
void ProduceIdShifts()
Count total number of visible childs under each node.
std::string ProduceModifyReply(int nodeid)
Return string with only part of nodes description which were modified Checks also volume.
void ResetRndrInfos()
Reset shape info, which used to pack binary data.
int fDrawIdCut
! sortid used for selection of most-significant nodes
void Build(TGeoManager *mgr, const std::string &volname="")
Collect information about geometry hierarchy into flat list like it done in JSROOT ClonedNodes....
Object with full description for drawing geometry It includes list of visible items and list of nodes...
RGeomConfig * cfg
current configurations
std::vector< RGeomVisible > visibles
all visible items
std::vector< RGeomNode * > nodes
all used nodes to display visible items and not known for client
int numnodes
total number of nodes in description
int sortid
! place in sorted array, to check cuts, or id of original node when used search structures
Definition RGeomData.hxx:45
int id
node id, index in array
Definition RGeomData.hxx:38
std::string color
rgb code without rgb() prefix
Definition RGeomData.hxx:44
Full node description including matrices and other attributes.
Definition RGeomData.hxx:54
float opacity
! opacity of the color
Definition RGeomData.hxx:61
A log configuration for a channel, e.g.
Definition RLogger.hxx:101
virtual Color_t GetFillColor() const
Return the fill area color.
Definition TAttFill.h:30
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition TAttFill.h:31
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:33
Class for serializing object to and from JavaScript Object Notation (JSON) format.
Definition TBufferJSON.h:30
static TString ToJSON(const T *obj, Int_t compact=0, const char *member_name=nullptr)
Definition TBufferJSON.h:75
@ kSkipTypeInfo
do not store typenames in JSON
Definition TBufferJSON.h:48
void SetCompact(int level)
Set level of space/newline/array compression Lower digit of compact parameter define formatting rules...
The color creation and management class.
Definition TColor.h:19
static Int_t GetColor(const char *hexcolor)
Static method returning color number for color specified by hex color string of form: "#rrggbb",...
Definition TColor.cxx:1823
@ kVisNone
Definition TGeoAtt.h:26
Box class.
Definition TGeoBBox.h:18
virtual Double_t GetDX() const
Definition TGeoBBox.h:74
Composite shapes are Boolean combinations of two or more shape components.
TGeoBoolNode * GetBoolNode() const
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:421
void Multiply(const TGeoMatrix *right)
multiply to the right with an other transformation if right is identity matrix, just return
static TClass * Class()
A geometry iterator.
Definition TGeoNode.h:245
void Skip()
Stop iterating the current branch.
The manager class for any TGeo geometry.
Definition TGeoManager.h:45
TGeoVolume * GetVolume(const char *name) const
Search for a named volume. All trailing blanks stripped.
Int_t GetMaxVisNodes() const
TGeoNode * GetTopNode() const
void SetNsegments(Int_t nseg)
Set number of segments for approximating circles in drawing.
Int_t GetVisLevel() const
Returns current depth to which geometry is drawn.
Int_t GetNsegments() const
Get number of segments approximating circles.
TGeoVolume * GetTopVolume() const
Geometrical transformation package.
Definition TGeoMatrix.h:41
virtual const Double_t * GetTranslation() const =0
virtual TClass * IsA() const
Definition TGeoMatrix.h:109
virtual void LocalToMaster(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
virtual const Double_t * GetScale() const =0
Bool_t IsIdentity() const
Definition TGeoMatrix.h:66
virtual const Double_t * GetRotationMatrix() const =0
TGeoMaterial * GetMaterial() const
Definition TGeoMedium.h:52
A node represent a volume positioned inside another.They store links to both volumes and to the TGeoM...
Definition TGeoNode.h:41
TGeoVolume * GetVolume() const
Definition TGeoNode.h:97
Int_t GetNumber() const
Definition TGeoNode.h:95
void SetNumber(Int_t number)
Definition TGeoNode.h:116
static TClass * Class()
static TClass * Class()
Base abstract class for all shapes.
Definition TGeoShape.h:26
virtual Bool_t IsComposite() const
Definition TGeoShape.h:130
virtual Bool_t IsCylType() const =0
virtual const char * GetName() const
Get the shape name.
virtual TBuffer3D * MakeBuffer3D() const
Definition TGeoShape.h:143
static TClass * Class()
static TClass * Class()
static TClass * Class()
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:174
TObjArray * GetNodes()
Definition TGeoVolume.h:168
static TGeoMedium * DummyMedium()
TGeoShape * GetShape() const
Definition TGeoVolume.h:189
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
const char * Data() const
Definition TString.h:380
const Int_t n
Definition legend1.C:16
std::function< bool(RGeomNode &, std::vector< int > &, bool, int)> RGeomScanFunc_t
RLogChannel & RGeomLog()
Log channel for Eve diagnostics.
Definition RGeomData.cxx:40
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
Definition test.py:1