Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
REveGeomData.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/REveGeomData.hxx>
13 
14 #include <ROOT/RBrowserItem.hxx>
16 #include <ROOT/REveUtil.hxx>
17 #include <ROOT/RLogger.hxx>
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 "TGeoCompositeShape.h"
30 #include "TObjArray.h"
31 #include "TBuffer3D.h"
32 #include "TBufferJSON.h"
33 
34 #include <algorithm>
35 
36 using namespace std::string_literals;
37 
38 
39 /** Base class for iterating of hierarchical structure */
40 
41 namespace ROOT {
42 namespace Experimental {
43 
44 
45 class RGeomBrowserIter {
46 
47  REveGeomDescription &fDesc;
48  int fParentId{-1};
49  unsigned fChild{0};
50  int fNodeId{0};
51 
52  std::vector<int> fStackParents;
53  std::vector<int> fStackChilds;
54 
55 public:
56 
57  RGeomBrowserIter(REveGeomDescription &desc) : fDesc(desc) {}
58 
59  const std::string &GetName() const { return fDesc.fDesc[fNodeId].name; }
60 
61  bool IsValid() const { return fNodeId >= 0; }
62 
63  int GetNodeId() const { return fNodeId; }
64 
65  bool HasChilds() const { return (fNodeId < 0) ? true : fDesc.fDesc[fNodeId].chlds.size() > 0; }
66 
67  int NumChilds() const { return (fNodeId < 0) ? 1 : fDesc.fDesc[fNodeId].chlds.size(); }
68 
69  bool Enter()
70  {
71  if (fNodeId < 0) {
72  Reset();
73  fNodeId = 0;
74  return true;
75  }
76 
77  auto &node = fDesc.fDesc[fNodeId];
78  if (node.chlds.size() == 0) return false;
79  fStackParents.emplace_back(fParentId);
80  fStackChilds.emplace_back(fChild);
81  fParentId = fNodeId;
82  fChild = 0;
83  fNodeId = node.chlds[fChild];
84  return true;
85  }
86 
87  bool Leave()
88  {
89  if (fStackParents.size() == 0) {
90  fNodeId = -1;
91  return false;
92  }
93  fParentId = fStackParents.back();
94  fChild = fStackChilds.back();
95 
96  fStackParents.pop_back();
97  fStackChilds.pop_back();
98 
99  if (fParentId < 0) {
100  fNodeId = 0;
101  } else {
102  fNodeId = fDesc.fDesc[fParentId].chlds[fChild];
103  }
104  return true;
105  }
106 
107  bool Next()
108  {
109  // does not have parents
110  if ((fNodeId <= 0) || (fParentId < 0)) {
111  Reset();
112  return false;
113  }
114 
115  auto &prnt = fDesc.fDesc[fParentId];
116  if (++fChild >= prnt.chlds.size()) {
117  fNodeId = -1; // not valid node, only Leave can be called
118  return false;
119  }
120 
121  fNodeId = prnt.chlds[fChild];
122  return true;
123  }
124 
125  bool Reset()
126  {
127  fParentId = -1;
128  fNodeId = -1;
129  fChild = 0;
130  fStackParents.clear();
131  fStackChilds.clear();
132 
133  return true;
134  }
135 
136  bool NextNode()
137  {
138  if (Enter()) return true;
139 
140  if (Next()) return true;
141 
142  while (Leave()) {
143  if (Next()) return true;
144  }
145 
146  return false;
147  }
148 
149  /** Navigate to specified path. For now path should start from '/' */
150 
151  bool Navigate(const std::string &path)
152  {
153  size_t pos = path.find("/");
154  if (pos != 0) return false;
155 
156  Reset(); // set to the top of element
157 
158  while (++pos < path.length()) {
159  auto last = pos;
160 
161  pos = path.find("/", last);
162 
163  if (pos == std::string::npos) pos = path.length();
164 
165  std::string folder = path.substr(last, pos-last);
166 
167  if (!Enter()) return false;
168 
169  bool find = false;
170 
171  do {
172  find = (folder.compare(GetName()) == 0);
173  } while (!find && Next());
174 
175  if (!find) return false;
176  }
177 
178  return true;
179  }
180 
181  /// Returns array of ids to currently selected node
182  std::vector<int> CurrentIds() const
183  {
184  std::vector<int> res;
185  if (IsValid()) {
186  for (unsigned n=1;n<fStackParents.size();++n)
187  res.emplace_back(fStackParents[n]);
188  if (fParentId >= 0) res.emplace_back(fParentId);
189  res.emplace_back(fNodeId);
190  }
191  return res;
192  }
193 
194 };
195 } // namespace Experimental
196 } // namespace ROOT
197 
198 /////////////////////////////////////////////////////////////////////
199 /// Pack matrix into vector, which can be send to client
200 /// Following sizes can be used for vector:
201 /// 0 - Identity matrix
202 /// 3 - Translation
203 /// 4 - Scale (last element always 1)
204 /// 9 - Rotation
205 /// 16 - Full size
206 
207 void ROOT::Experimental::REveGeomDescription::PackMatrix(std::vector<float> &vect, TGeoMatrix *matr)
208 {
209  vect.clear();
210 
211  if (!matr || matr->IsIdentity()) {
212  return;
213  }
214 
215  auto trans = matr->GetTranslation();
216  auto scale = matr->GetScale();
217  auto rotate = matr->GetRotationMatrix();
218 
219  bool is_translate = matr->IsA() == TGeoTranslation::Class(),
220  is_scale = matr->IsA() == TGeoScale::Class(),
221  is_rotate = matr->IsA() == TGeoRotation::Class();
222 
223  if (!is_translate && !is_scale && !is_rotate) {
224  // check if trivial matrix
225 
226  auto test = [](double val, double chk) { return (val==chk) || (TMath::Abs(val-chk) < 1e-20); };
227 
228  bool no_scale = test(scale[0],1) && test(scale[1],1) && test(scale[2],1);
229  bool no_trans = test(trans[0],0) && test(trans[1],0) && test(trans[2],0);
230  bool no_rotate = test(rotate[0],1) && test(rotate[1],0) && test(rotate[2],0) &&
231  test(rotate[3],0) && test(rotate[4],1) && test(rotate[5],0) &&
232  test(rotate[6],0) && test(rotate[7],0) && test(rotate[8],1);
233 
234  if (no_scale && no_trans && no_rotate)
235  return;
236 
237  if (no_scale && no_trans && !no_rotate) {
238  is_rotate = true;
239  } else if (no_scale && !no_trans && no_rotate) {
240  is_translate = true;
241  } else if (!no_scale && no_trans && no_rotate) {
242  is_scale = true;
243  }
244  }
245 
246  if (is_translate) {
247  vect.resize(3);
248  vect[0] = trans[0];
249  vect[1] = trans[1];
250  vect[2] = trans[2];
251  return;
252  }
253 
254  if (is_scale) {
255  vect.resize(4);
256  vect[0] = scale[0];
257  vect[1] = scale[1];
258  vect[2] = scale[2];
259  vect[3] = 1;
260  return;
261  }
262 
263  if (is_rotate) {
264  vect.resize(9);
265  for (int n=0;n<9;++n)
266  vect[n] = rotate[n];
267  return;
268  }
269 
270  vect.resize(16);
271  vect[0] = rotate[0]; vect[4] = rotate[1]; vect[8] = rotate[2]; vect[12] = trans[0];
272  vect[1] = rotate[3]; vect[5] = rotate[4]; vect[9] = rotate[5]; vect[13] = trans[1];
273  vect[2] = rotate[6]; vect[6] = rotate[7]; vect[10] = rotate[8]; vect[14] = trans[2];
274  vect[3] = 0; vect[7] = 0; vect[11] = 0; vect[15] = 1;
275 }
276 
277 /////////////////////////////////////////////////////////////////////
278 /// Collect information about geometry hierarchy into flat list
279 /// like it done JSROOT.GEO.ClonedNodes.prototype.CreateClones
280 
281 void ROOT::Experimental::REveGeomDescription::Build(TGeoManager *mgr, const std::string &volname)
282 {
283  fDesc.clear();
284  fNodes.clear();
285  fSortMap.clear();
286  ClearDrawData();
287  fDrawIdCut = 0;
288 
289  if (!mgr) return;
290 
291  auto topnode = mgr->GetTopNode();
292  if (!volname.empty()) {
293  auto vol = mgr->GetVolume(volname.c_str());
294  if (vol) {
295  TGeoNode *node;
296  TGeoIterator next(mgr->GetTopVolume());
297  while ((node=next())) {
298  if (node->GetVolume() == vol) break;
299  }
300  if (node) { topnode = node; printf("Find node with volume\n"); }
301  }
302  }
303 
304  // by top node visibility always enabled and harm logic
305  // later visibility can be controlled by other means
306  // mgr->GetTopNode()->GetVolume()->SetVisibility(kFALSE);
307 
308  int maxnodes = mgr->GetMaxVisNodes();
309 
310  SetNSegments(mgr->GetNsegments());
311  SetVisLevel(mgr->GetVisLevel());
312  SetMaxVisNodes(maxnodes);
313  SetMaxVisFaces( (maxnodes > 5000 ? 5000 : (maxnodes < 1000 ? 1000 : maxnodes)) * 100);
314 
315  // vector to remember numbers
316  std::vector<int> numbers;
317  int offset = 1000000000;
318 
319  // try to build flat list of all nodes
320  TGeoNode *snode = topnode;
321  TGeoIterator iter(topnode->GetVolume());
322  do {
323  // artificial offset, used as identifier
324  if (snode->GetNumber() >= offset) {
325  iter.Skip(); // no need to look inside
326  } else {
327  numbers.emplace_back(snode->GetNumber());
328  snode->SetNumber(offset + fNodes.size()); // use id with shift 1e9
329  fNodes.emplace_back(snode);
330  }
331  } while ((snode = iter()) != nullptr);
332 
333  fDesc.reserve(fNodes.size());
334  numbers.reserve(fNodes.size());
335  fSortMap.reserve(fNodes.size());
336 
337  // array for sorting
338  std::vector<REveGeomNode *> sortarr;
339  sortarr.reserve(fNodes.size());
340 
341  // create vector of desc and childs
342  int cnt = 0;
343  for (auto &node: fNodes) {
344 
345  fDesc.emplace_back(node->GetNumber() - offset);
346  auto &desc = fDesc[cnt++];
347 
348  sortarr.emplace_back(&desc);
349 
350  desc.name = node->GetName();
351 
352  auto shape = dynamic_cast<TGeoBBox *>(node->GetVolume()->GetShape());
353  if (shape) {
354  desc.vol = shape->GetDX()*shape->GetDY()*shape->GetDZ();
355  desc.nfaces = 12; // TODO: get better value for each shape - excluding composite
356  }
357 
358  CopyMaterialProperties(node->GetVolume(), desc);
359 
360  auto chlds = node->GetNodes();
361 
362  PackMatrix(desc.matr, node->GetMatrix());
363 
364  if (chlds)
365  for (int n = 0; n <= chlds->GetLast(); ++n) {
366  auto chld = dynamic_cast<TGeoNode *> (chlds->At(n));
367  desc.chlds.emplace_back(chld->GetNumber()-offset);
368  }
369  }
370 
371  // recover numbers
372  cnt = 0;
373  for (auto &node: fNodes)
374  node->SetNumber(numbers[cnt++]);
375 
376  // sort in volume descent order
377  std::sort(sortarr.begin(), sortarr.end(), [](REveGeomNode *a, REveGeomNode * b) { return a->vol > b->vol; });
378 
379  cnt = 0;
380  for (auto &elem: sortarr) {
381  fSortMap.emplace_back(elem->id);
382  elem->sortid = cnt++; // keep place in sorted array to correctly apply cut
383  }
384 
385  MarkVisible(); // set visibility flags
386 
387  ProduceIdShifts();
388 }
389 
390 /////////////////////////////////////////////////////////////////////
391 /// Set visibility flag for each nodes
392 
393 int ROOT::Experimental::REveGeomDescription::MarkVisible(bool on_screen)
394 {
395  int res = 0, cnt = 0;
396  for (auto &node: fNodes) {
397  auto &desc = fDesc[cnt++];
398  desc.vis = 0;
399  desc.nochlds = false;
400 
401  if (on_screen) {
402  if (node->IsOnScreen())
403  desc.vis = 99;
404  } else {
405  auto vol = node->GetVolume();
406 
407  if (vol->IsVisible() && !vol->TestAttBit(TGeoAtt::kVisNone))
408  desc.vis = 99;
409 
410  if (!node->IsVisDaughters())
411  desc.nochlds = true;
412 
413  if ((desc.vis > 0) && (desc.chlds.size() > 0) && !desc.nochlds)
414  desc.vis = 1;
415  }
416 
417  if (desc.IsVisible() && desc.CanDisplay()) res++;
418  }
419 
420  return res;
421 }
422 
423 /////////////////////////////////////////////////////////////////////
424 /// Count total number of visible childs under each node
425 
426 void ROOT::Experimental::REveGeomDescription::ProduceIdShifts()
427 {
428  for (auto &node : fDesc)
429  node.idshift = -1;
430 
431  using ScanFunc_t = std::function<int(REveGeomNode &)>;
432 
433  ScanFunc_t scan_func = [&, this](REveGeomNode &node) {
434  if (node.idshift < 0) {
435  node.idshift = 0;
436  for(auto id : node.chlds)
437  node.idshift += scan_func(fDesc[id]);
438  }
439 
440  return node.idshift + 1;
441  };
442 
443  if (fDesc.size() > 0)
444  scan_func(fDesc[0]);
445 }
446 
447 /////////////////////////////////////////////////////////////////////
448 /// Iterate over all nodes and call function for visible
449 
450 int ROOT::Experimental::REveGeomDescription::ScanNodes(bool only_visible, int maxlvl, REveGeomScanFunc_t func)
451 {
452  std::vector<int> stack;
453  stack.reserve(25); // reserve enough space for most use-cases
454  int counter{0};
455 
456  using ScanFunc_t = std::function<int(int, int)>;
457 
458  ScanFunc_t scan_func = [&, this](int nodeid, int lvl) {
459  auto &desc = fDesc[nodeid];
460  int res = 0;
461 
462  if (desc.nochlds && (lvl > 0)) lvl = 0;
463 
464  // same logic as in JSROOT.GEO.ClonedNodes.prototype.ScanVisible
465  bool is_visible = (lvl >= 0) && (desc.vis > lvl) && desc.CanDisplay();
466 
467  if (is_visible || !only_visible)
468  if (func(desc, stack, is_visible, counter))
469  res++;
470 
471  counter++; // count sequence id of current position in scan, will be used later for merging drawing lists
472 
473  if ((desc.chlds.size() > 0) && ((lvl > 0) || !only_visible)) {
474  auto pos = stack.size();
475  stack.emplace_back(0);
476  for (unsigned k = 0; k < desc.chlds.size(); ++k) {
477  stack[pos] = k; // stack provides index in list of chdils
478  res += scan_func(desc.chlds[k], lvl - 1);
479  }
480  stack.pop_back();
481  } else {
482  counter += desc.idshift;
483  }
484 
485  return res;
486  };
487 
488  if (!maxlvl && (GetVisLevel() > 0)) maxlvl = GetVisLevel();
489  if (!maxlvl) maxlvl = 4;
490  if (maxlvl > 97) maxlvl = 97; // check while vis property of node is 99 normally
491 
492  return scan_func(0, maxlvl);
493 }
494 
495 /////////////////////////////////////////////////////////////////////
496 /// Collect nodes which are used in visibles
497 
498 void ROOT::Experimental::REveGeomDescription::CollectNodes(REveGeomDrawing &drawing)
499 {
500  // TODO: for now reset all flags, later can be kept longer
501  for (auto &node : fDesc)
502  node.useflag = false;
503 
504  drawing.cfg = &fCfg;
505 
506  drawing.numnodes = fDesc.size();
507 
508  for (auto &item : drawing.visibles) {
509  int nodeid{0};
510  for (auto &chindx : item.stack) {
511  auto &node = fDesc[nodeid];
512  if (!node.useflag) {
513  node.useflag = true;
514  drawing.nodes.emplace_back(&node);
515  }
516  if (chindx >= (int)node.chlds.size())
517  break;
518  nodeid = node.chlds[chindx];
519  }
520 
521  auto &node = fDesc[nodeid];
522  if (!node.useflag) {
523  node.useflag = true;
524  drawing.nodes.emplace_back(&node);
525  }
526  }
527 
528  printf("SELECT NODES %d\n", (int) drawing.nodes.size());
529 }
530 
531 /////////////////////////////////////////////////////////////////////
532 /// Find description object for requested shape
533 /// If not exists - will be created
534 
535 std::string ROOT::Experimental::REveGeomDescription::ProcessBrowserRequest(const std::string &msg)
536 {
537  std::string res;
538 
539  auto request = TBufferJSON::FromJSON<RBrowserRequest>(msg);
540 
541  if (msg.empty()) {
542  request = std::make_unique<RBrowserRequest>();
543  request->path = "/";
544  request->first = 0;
545  request->number = 100;
546  }
547 
548  if (!request)
549  return res;
550 
551  if ((request->path.compare("/") == 0) && (request->first == 0) && (GetNumNodes() < (IsPreferredOffline() ? 1000000 : 1000))) {
552 
553  std::vector<REveGeomNodeBase *> vect(fDesc.size(), nullptr);
554 
555  int cnt = 0;
556  for (auto &item : fDesc)
557  vect[cnt++]= &item;
558 
559  res = "DESCR:"s + TBufferJSON::ToJSON(&vect,GetJsonComp()).Data();
560 
561  // example how iterator can be used
562  RGeomBrowserIter iter(*this);
563  int nelements = 0;
564  while (iter.NextNode())
565  nelements++;
566  printf("Total number of valid nodes %d\n", nelements);
567 
568  } else {
569  std::vector<RBrowserItem> temp_nodes;
570  bool toplevel = (request->path.compare("/") == 0);
571 
572  // create temporary object for the short time
573  RBrowserReply reply;
574  reply.path = request->path;
575  reply.first = request->first;
576 
577  RGeomBrowserIter iter(*this);
578  if (iter.Navigate(request->path)) {
579 
580  reply.nchilds = iter.NumChilds();
581  // scan childs of selected nodes
582  if (iter.Enter()) {
583 
584  while ((request->first > 0) && iter.Next()) {
585  request->first--;
586  }
587 
588  while (iter.IsValid() && (request->number > 0)) {
589  temp_nodes.emplace_back(iter.GetName(), iter.NumChilds());
590  if (toplevel) temp_nodes.back().SetExpanded(true);
591  request->number--;
592  if (!iter.Next()) break;
593  }
594  }
595  }
596 
597  for (auto &n : temp_nodes)
598  reply.nodes.emplace_back(&n);
599 
600  res = "BREPL:"s + TBufferJSON::ToJSON(&reply, GetJsonComp()).Data();
601  }
602 
603  return res;
604 }
605 
606 /////////////////////////////////////////////////////////////////////
607 /// Find description object for requested shape
608 /// If not exists - will be created
609 
610 ROOT::Experimental::REveGeomDescription::ShapeDescr &ROOT::Experimental::REveGeomDescription::FindShapeDescr(TGeoShape *shape)
611 {
612  for (auto &descr : fShapes)
613  if (descr.fShape == shape)
614  return descr;
615 
616  fShapes.emplace_back(shape);
617  auto &elem = fShapes.back();
618  elem.id = fShapes.size() - 1;
619  return elem;
620 }
621 
622 /////////////////////////////////////////////////////////////////////
623 /// Find description object and create render information
624 
625 ROOT::Experimental::REveGeomDescription::ShapeDescr &
626 ROOT::Experimental::REveGeomDescription::MakeShapeDescr(TGeoShape *shape)
627 {
628  auto &elem = FindShapeDescr(shape);
629 
630  if (elem.nfaces == 0) {
631 
632  TGeoCompositeShape *comp = nullptr;
633 
634  int boundary = 3; //
635  if (shape->IsComposite()) {
636  comp = dynamic_cast<TGeoCompositeShape *>(shape);
637  // composite is most complex for client, therefore by default build on server
638  boundary = 1;
639  } else if (!shape->IsCylType()) {
640  // simple box geometry is compact and can be delivered as raw
641  boundary = 2;
642  }
643 
644  if (IsBuildShapes() < boundary) {
645  elem.nfaces = 1;
646  elem.fShapeInfo.shape = shape;
647  } else {
648 
649  auto poly = std::make_unique<REveGeoPolyShape>();
650 
651  if (comp) {
652  poly->BuildFromComposite(comp, GetNSegments());
653  } else {
654  poly->BuildFromShape(shape, GetNSegments());
655  }
656 
657  REveRenderData rd;
658 
659  poly->FillRenderData(rd);
660 
661  elem.nfaces = poly->GetNumFaces();
662 
663  elem.fRawInfo.raw.resize(rd.GetBinarySize());
664  rd.Write( reinterpret_cast<char *>(elem.fRawInfo.raw.data()), elem.fRawInfo.raw.size() );
665  elem.fRawInfo.sz[0] = rd.SizeV();
666  elem.fRawInfo.sz[1] = rd.SizeN();
667  elem.fRawInfo.sz[2] = rd.SizeI();
668 
669  }
670  }
671 
672  return elem;
673 }
674 
675 /////////////////////////////////////////////////////////////////////
676 /// Copy material properties
677 
678 void ROOT::Experimental::REveGeomDescription::CopyMaterialProperties(TGeoVolume *volume, REveGeomNode &node)
679 {
680  if (!volume) return;
681 
682  TColor *col{nullptr};
683 
684  if ((volume->GetFillColor() > 1) && (volume->GetLineColor() == 1))
685  col = gROOT->GetColor(volume->GetFillColor());
686  else if (volume->GetLineColor() >= 0)
687  col = gROOT->GetColor(volume->GetLineColor());
688 
689  if (volume->GetMedium() && (volume->GetMedium() != TGeoVolume::DummyMedium()) && volume->GetMedium()->GetMaterial()) {
690  auto material = volume->GetMedium()->GetMaterial();
691 
692  auto fillstyle = material->GetFillStyle();
693  if ((fillstyle>=3000) && (fillstyle<=3100)) node.opacity = (3100 - fillstyle) / 100.;
694  if (!col) col = gROOT->GetColor(material->GetFillColor());
695  }
696 
697  if (col) {
698  node.color = std::to_string((int)(col->GetRed()*255)) + "," +
699  std::to_string((int)(col->GetGreen()*255)) + "," +
700  std::to_string((int)(col->GetBlue()*255));
701  if (node.opacity == 1.)
702  node.opacity = col->GetAlpha();
703  } else {
704  node.color.clear();
705  }
706 }
707 
708 /////////////////////////////////////////////////////////////////////
709 /// Reset shape info, which used to pack binary data
710 
711 void ROOT::Experimental::REveGeomDescription::ResetRndrInfos()
712 {
713  for (auto &s: fShapes)
714  s.reset();
715 }
716 
717 /////////////////////////////////////////////////////////////////////
718 /// Collect all information required to draw geometry on the client
719 /// This includes list of each visible nodes, meshes and matrixes
720 
721 bool ROOT::Experimental::REveGeomDescription::CollectVisibles()
722 {
723  std::vector<int> viscnt(fDesc.size(), 0);
724 
725  int level = GetVisLevel();
726 
727  // first count how many times each individual node appears
728  int numnodes = ScanNodes(true, level, [&viscnt](REveGeomNode &node, std::vector<int> &, bool, int) {
729  viscnt[node.id]++;
730  return true;
731  });
732 
733  if (GetMaxVisNodes() > 0) {
734  while ((numnodes > GetMaxVisNodes()) && (level > 1)) {
735  level--;
736  viscnt.assign(viscnt.size(), 0);
737  numnodes = ScanNodes(true, level, [&viscnt](REveGeomNode &node, std::vector<int> &, bool, int) {
738  viscnt[node.id]++;
739  return true;
740  });
741  }
742  }
743 
744  fActualLevel = level;
745  fDrawIdCut = 0;
746 
747  int totalnumfaces{0}, totalnumnodes{0};
748 
749  //for (auto &node : fDesc)
750  // node.SetDisplayed(false);
751 
752  // build all shapes in volume decreasing order
753  for (auto &sid: fSortMap) {
754  fDrawIdCut++; //
755  auto &desc = fDesc[sid];
756 
757  if ((viscnt[sid] <= 0) || (desc.vol <= 0)) continue;
758 
759  auto shape = fNodes[sid]->GetVolume()->GetShape();
760  if (!shape) continue;
761 
762  // now we need to create TEveGeoPolyShape, which can provide all rendering data
763  auto &shape_descr = MakeShapeDescr(shape);
764 
765  // should not happen, but just in case
766  if (shape_descr.nfaces <= 0) {
767  R__ERROR_HERE("webeve") << "No faces for the shape " << shape->GetName() << " class " << shape->ClassName();
768  continue;
769  }
770 
771  // check how many faces are created
772  totalnumfaces += shape_descr.nfaces * viscnt[sid];
773  if ((GetMaxVisFaces() > 0) && (totalnumfaces > GetMaxVisFaces())) break;
774 
775  // also avoid too many nodes
776  totalnumnodes += viscnt[sid];
777  if ((GetMaxVisNodes() > 0) && (totalnumnodes > GetMaxVisNodes())) break;
778 
779  // desc.SetDisplayed(true);
780  }
781 
782  // finally we should create data for streaming to the client
783  // it includes list of visible nodes and rawdata
784 
785  REveGeomDrawing drawing;
786  ResetRndrInfos();
787  bool has_shape = false;
788 
789  ScanNodes(true, level, [&, this](REveGeomNode &node, std::vector<int> &stack, bool, int seqid) {
790  if (node.sortid < fDrawIdCut) {
791  drawing.visibles.emplace_back(node.id, seqid, stack);
792 
793  auto &item = drawing.visibles.back();
794  item.color = node.color;
795  item.opacity = node.opacity;
796 
797  auto volume = fNodes[node.id]->GetVolume();
798 
799  auto &sd = MakeShapeDescr(volume->GetShape());
800 
801  item.ri = sd.rndr_info();
802  if (sd.has_shape()) has_shape = true;
803  }
804  return true;
805  });
806 
807  CollectNodes(drawing);
808 
809  fDrawJson = "GDRAW:"s + MakeDrawingJson(drawing, has_shape);
810 
811  return true;
812 }
813 
814 /////////////////////////////////////////////////////////////////////
815 /// Clear raw data. Will be rebuild when next connection will be established
816 
817 void ROOT::Experimental::REveGeomDescription::ClearDrawData()
818 {
819  fDrawJson.clear();
820 }
821 
822 /////////////////////////////////////////////////////////////////////
823 /// return true when node used in main geometry drawing and does not have childs
824 /// for such nodes one could provide optimize toggling of visibility flags
825 
826 bool ROOT::Experimental::REveGeomDescription::IsPrincipalEndNode(int nodeid)
827 {
828  if ((nodeid < 0) || (nodeid >= (int)fDesc.size()))
829  return false;
830 
831  auto &desc = fDesc[nodeid];
832 
833  return (desc.sortid < fDrawIdCut) && desc.IsVisible() && desc.CanDisplay() && (desc.chlds.size()==0);
834 }
835 
836 
837 /////////////////////////////////////////////////////////////////////
838 /// Search visible nodes for provided name
839 /// If number of found elements less than 100, create description and shapes for them
840 /// Returns number of match elements
841 
842 int ROOT::Experimental::REveGeomDescription::SearchVisibles(const std::string &find, std::string &hjson, std::string &json)
843 {
844  hjson.clear();
845  json.clear();
846 
847  if (find.empty()) {
848  hjson = "FOUND:RESET";
849  return 0;
850  }
851 
852  std::vector<int> nodescnt(fDesc.size(), 0), viscnt(fDesc.size(), 0);
853 
854  int nmatches{0};
855 
856  auto match_func = [&find](REveGeomNode &node) {
857  return (node.vol > 0) && (node.name.compare(0, find.length(), find) == 0);
858  };
859 
860  // first count how many times each individual node appears
861  ScanNodes(false, 0, [&nodescnt,&viscnt,&match_func,&nmatches](REveGeomNode &node, std::vector<int> &, bool is_vis, int) {
862 
863  if (match_func(node)) {
864  nmatches++;
865  nodescnt[node.id]++;
866  if (is_vis) viscnt[node.id]++;
867  };
868  return true;
869  });
870 
871  // do not send too much data, limit could be made configurable later
872  if (nmatches==0) {
873  hjson = "FOUND:NO";
874  return nmatches;
875  }
876 
877  if (nmatches > 10 * GetMaxVisNodes()) {
878  hjson = "FOUND:Too many " + std::to_string(nmatches);
879  return nmatches;
880  }
881 
882  // now build all necessary shapes and check number of faces - not too many
883 
884  int totalnumfaces{0}, totalnumnodes{0}, scnt{0};
885  bool send_rawdata{true};
886 
887  // build all shapes in volume decreasing order
888  for (auto &sid: fSortMap) {
889  if (scnt++ < fDrawIdCut) continue; // no need to send most significant shapes
890 
891  if (viscnt[sid] == 0) continue; // this node is not used at all
892 
893  auto &desc = fDesc[sid];
894  if ((viscnt[sid] <= 0) && (desc.vol <= 0)) continue;
895 
896  auto shape = fNodes[sid]->GetVolume()->GetShape();
897  if (!shape) continue;
898 
899  // create shape raw data
900  auto &shape_descr = MakeShapeDescr(shape);
901 
902  // should not happen, but just in case
903  if (shape_descr.nfaces <= 0) {
904  R__ERROR_HERE("webeve") << "No faces for the shape " << shape->GetName() << " class " << shape->ClassName();
905  continue;
906  }
907 
908  // check how many faces are created
909  totalnumfaces += shape_descr.nfaces * viscnt[sid];
910  if ((GetMaxVisFaces() > 0) && (totalnumfaces > GetMaxVisFaces())) { send_rawdata = false; break; }
911 
912  // also avoid too many nodes
913  totalnumnodes += viscnt[sid];
914  if ((GetMaxVisNodes() > 0) && (totalnumnodes > GetMaxVisNodes())) { send_rawdata = false; break; }
915  }
916 
917  // only for debug purposes - remove later
918  // send_rawdata = false;
919 
920  // finally we should create data for streaming to the client
921  // it includes list of visible nodes and rawdata (if there is enough space)
922 
923 
924  std::vector<REveGeomNodeBase> found_desc; ///<! hierarchy of nodes, used for search
925  std::vector<int> found_map(fDesc.size(), -1); ///<! mapping between nodeid - > foundid
926 
927  // these are only selected nodes to produce hierarchy
928 
929  found_desc.emplace_back(0);
930  found_desc[0].vis = fDesc[0].vis;
931  found_desc[0].name = fDesc[0].name;
932  found_desc[0].color = fDesc[0].color;
933  found_map[0] = 0;
934 
935  ResetRndrInfos();
936 
937  REveGeomDrawing drawing;
938  bool has_shape = true;
939 
940  ScanNodes(false, 0, [&, this](REveGeomNode &node, std::vector<int> &stack, bool is_vis, int seqid) {
941  // select only nodes which should match
942  if (!match_func(node))
943  return true;
944 
945  // add entries into hierarchy of found elements
946  int prntid = 0;
947  for (auto &s : stack) {
948  int chldid = fDesc[prntid].chlds[s];
949  if (found_map[chldid] <= 0) {
950  int newid = found_desc.size();
951  found_desc.emplace_back(newid); // potentially original id can be used here
952  found_map[chldid] = newid; // re-map into reduced hierarchy
953 
954  found_desc.back().vis = fDesc[chldid].vis;
955  found_desc.back().name = fDesc[chldid].name;
956  found_desc.back().color = fDesc[chldid].color;
957  }
958 
959  auto pid = found_map[prntid];
960  auto cid = found_map[chldid];
961 
962  // now add entry into childs lists
963  auto &pchlds = found_desc[pid].chlds;
964  if (std::find(pchlds.begin(), pchlds.end(), cid) == pchlds.end())
965  pchlds.emplace_back(cid);
966 
967  prntid = chldid;
968  }
969 
970  // no need to add visibles
971  if (!is_vis) return true;
972 
973  drawing.visibles.emplace_back(node.id, seqid, stack);
974 
975  // no need to transfer shape if it provided with main drawing list
976  // also no binary will be transported when too many matches are there
977  if (!send_rawdata || (node.sortid < fDrawIdCut)) {
978  // do not include render data
979  return true;
980  }
981 
982  auto &item = drawing.visibles.back();
983  auto volume = fNodes[node.id]->GetVolume();
984 
985  item.color = node.color;
986  item.opacity = node.opacity;
987 
988  auto &sd = MakeShapeDescr(volume->GetShape());
989 
990  item.ri = sd.rndr_info();
991  if (sd.has_shape()) has_shape = true;
992  return true;
993  });
994 
995  hjson = "FESCR:"s + TBufferJSON::ToJSON(&found_desc, GetJsonComp()).Data();
996 
997  CollectNodes(drawing);
998 
999  json = "FDRAW:"s + MakeDrawingJson(drawing, has_shape);
1000 
1001  return nmatches;
1002 }
1003 
1004 /////////////////////////////////////////////////////////////////////////////////
1005 /// Returns nodeid for given stack array, returns -1 in case of failure
1006 
1007 int ROOT::Experimental::REveGeomDescription::FindNodeId(const std::vector<int> &stack)
1008 {
1009  int nodeid{0};
1010 
1011  for (auto &chindx: stack) {
1012  auto &node = fDesc[nodeid];
1013  if (chindx >= (int) node.chlds.size()) return -1;
1014  nodeid = node.chlds[chindx];
1015  }
1016 
1017  return nodeid;
1018 }
1019 
1020 /////////////////////////////////////////////////////////////////////////////////
1021 /// Creates stack for given array of ids, first element always should be 0
1022 
1023 std::vector<int> ROOT::Experimental::REveGeomDescription::MakeStackByIds(const std::vector<int> &ids)
1024 {
1025  std::vector<int> stack;
1026 
1027  if (ids[0] != 0) {
1028  printf("Wrong first id\n");
1029  return stack;
1030  }
1031 
1032  int nodeid = 0;
1033 
1034  for (unsigned k = 1; k < ids.size(); ++k) {
1035 
1036  int prntid = nodeid;
1037  nodeid = ids[k];
1038 
1039  if (nodeid >= (int) fDesc.size()) {
1040  printf("Wrong node id %d\n", nodeid);
1041  stack.clear();
1042  return stack;
1043  }
1044  auto &chlds = fDesc[prntid].chlds;
1045  auto pos = std::find(chlds.begin(), chlds.end(), nodeid);
1046  if (pos == chlds.end()) {
1047  printf("Wrong id %d not a child of %d - fail to find stack num %d\n", nodeid, prntid, (int) chlds.size());
1048  stack.clear();
1049  return stack;
1050  }
1051 
1052  stack.emplace_back(std::distance(chlds.begin(), pos));
1053  }
1054 
1055  return stack;
1056 }
1057 
1058 /////////////////////////////////////////////////////////////////////////////////
1059 /// Produce stack based on string path
1060 /// Used to highlight geo volumes by browser hover event
1061 
1062 std::vector<int> ROOT::Experimental::REveGeomDescription::MakeStackByPath(const std::string &path)
1063 {
1064  std::vector<int> res;
1065 
1066  RGeomBrowserIter iter(*this);
1067 
1068  if (iter.Navigate(path)) {
1069 // auto ids = iter.CurrentIds();
1070 // printf("path %s ", path.c_str());
1071 // for (auto &id: ids)
1072 // printf("%d ", id);
1073 // printf("\n");
1074  res = MakeStackByIds(iter.CurrentIds());
1075  }
1076 
1077  return res;
1078 }
1079 
1080 /////////////////////////////////////////////////////////////////////////////////
1081 /// Produce list of node ids for given stack
1082 /// If found nodes preselected - use their ids
1083 
1084 std::vector<int> ROOT::Experimental::REveGeomDescription::MakeIdsByStack(const std::vector<int> &stack)
1085 {
1086  std::vector<int> ids;
1087 
1088  ids.emplace_back(0);
1089  int nodeid = 0;
1090  bool failure = false;
1091 
1092  for (auto s : stack) {
1093  auto &chlds = fDesc[nodeid].chlds;
1094  if (s >= (int) chlds.size()) { failure = true; break; }
1095 
1096  ids.emplace_back(chlds[s]);
1097 
1098  nodeid = chlds[s];
1099  }
1100 
1101  if (failure) {
1102  printf("Fail to convert stack into list of nodes\n");
1103  ids.clear();
1104  }
1105 
1106  return ids;
1107 }
1108 
1109 /////////////////////////////////////////////////////////////////////////////////
1110 /// Returns path string for provided stack
1111 
1112 std::string ROOT::Experimental::REveGeomDescription::MakePathByStack(const std::vector<int> &stack)
1113 {
1114  std::string path;
1115 
1116  auto ids = MakeIdsByStack(stack);
1117  if (ids.size() > 0) {
1118  path = "/";
1119  for (auto &id : ids) {
1120  path.append(fDesc[id].name);
1121  path.append("/");
1122  }
1123  }
1124 
1125  return path;
1126 }
1127 
1128 
1129 /////////////////////////////////////////////////////////////////////////////////
1130 /// Return string with only part of nodes description which were modified
1131 /// Checks also volume
1132 
1133 std::string ROOT::Experimental::REveGeomDescription::ProduceModifyReply(int nodeid)
1134 {
1135  std::vector<REveGeomNodeBase *> nodes;
1136  auto vol = fNodes[nodeid]->GetVolume();
1137 
1138  // we take not only single node, but all there same volume is referenced
1139  // nodes.push_back(&fDesc[nodeid]);
1140 
1141  int id{0};
1142  for (auto &desc : fDesc)
1143  if (fNodes[id++]->GetVolume() == vol)
1144  nodes.emplace_back(&desc);
1145 
1146  return "MODIF:"s + TBufferJSON::ToJSON(&nodes, GetJsonComp()).Data();
1147 }
1148 
1149 
1150 /////////////////////////////////////////////////////////////////////////////////
1151 /// Produce shape rendering data for given stack
1152 /// All nodes, which are referencing same shape will be transferred
1153 /// Returns true if new render information provided
1154 
1155 bool ROOT::Experimental::REveGeomDescription::ProduceDrawingFor(int nodeid, std::string &json, bool check_volume)
1156 {
1157  // only this shape is interesting
1158 
1159  TGeoVolume *vol = (nodeid < 0) ? nullptr : fNodes[nodeid]->GetVolume();
1160 
1161  if (!vol || !vol->GetShape()) {
1162  json.append("NO");
1163  return false;
1164  }
1165 
1166  REveGeomDrawing drawing;
1167 
1168  ScanNodes(true, 0, [&, this](REveGeomNode &node, std::vector<int> &stack, bool, int seq_id) {
1169  // select only nodes which reference same shape
1170 
1171  if (check_volume) {
1172  if (fNodes[node.id]->GetVolume() != vol) return true;
1173  } else {
1174  if (node.id != nodeid) return true;
1175  }
1176 
1177  drawing.visibles.emplace_back(node.id, seq_id, stack);
1178 
1179  auto &item = drawing.visibles.back();
1180 
1181  item.color = node.color;
1182  item.opacity = node.opacity;
1183  return true;
1184  });
1185 
1186  // no any visible nodes were done
1187  if (drawing.visibles.size()==0) {
1188  json.append("NO");
1189  return false;
1190  }
1191 
1192  ResetRndrInfos();
1193 
1194  bool has_shape = false, has_raw = false;
1195 
1196  auto &sd = MakeShapeDescr(vol->GetShape());
1197 
1198  // assign shape data
1199  for (auto &item : drawing.visibles) {
1200  item.ri = sd.rndr_info();
1201  if (sd.has_shape()) has_shape = true;
1202  if (sd.has_raw()) has_raw = true;
1203  }
1204 
1205  CollectNodes(drawing);
1206 
1207  json.append(MakeDrawingJson(drawing, has_shape));
1208 
1209  return has_raw || has_shape;
1210 }
1211 
1212 /////////////////////////////////////////////////////////////////////////////////
1213 /// Produce JSON for the drawing
1214 /// If TGeoShape appears in the drawing, one has to keep typeinfo
1215 /// But in this case one can exclude several classes which are not interesting,
1216 /// but appears very often
1217 
1218 std::string ROOT::Experimental::REveGeomDescription::MakeDrawingJson(REveGeomDrawing &drawing, bool has_shapes)
1219 {
1220  int comp = GetJsonComp();
1221 
1222  if (!has_shapes || (comp < TBufferJSON::kSkipTypeInfo))
1223  return TBufferJSON::ToJSON(&drawing, comp).Data();
1224 
1225  comp = comp % TBufferJSON::kSkipTypeInfo; // no typeingo skipping
1226 
1227  TBufferJSON json;
1228  json.SetCompact(comp);
1229  json.SetSkipClassInfo(TClass::GetClass<REveGeomDrawing>());
1230  json.SetSkipClassInfo(TClass::GetClass<REveGeomNode>());
1231  json.SetSkipClassInfo(TClass::GetClass<REveGeomVisible>());
1232  json.SetSkipClassInfo(TClass::GetClass<RGeomShapeRenderInfo>());
1233  json.SetSkipClassInfo(TClass::GetClass<RGeomRawRenderInfo>());
1234 
1235  return json.StoreObject(&drawing, TClass::GetClass<REveGeomDrawing>()).Data();
1236 }
1237 
1238 /////////////////////////////////////////////////////////////////////////////////
1239 /// Change visibility for specified element
1240 /// Returns true if changes was performed
1241 
1242 bool ROOT::Experimental::REveGeomDescription::ChangeNodeVisibility(int nodeid, bool selected)
1243 {
1244  auto &dnode = fDesc[nodeid];
1245 
1246  auto vol = fNodes[nodeid]->GetVolume();
1247 
1248  // nothing changed
1249  if (vol->IsVisible() == selected)
1250  return false;
1251 
1252  dnode.vis = selected ? 99 : 0;
1253  vol->SetVisibility(selected);
1254  if (dnode.chlds.size() > 0) {
1255  if (selected) dnode.vis = 1; // visibility disabled when any child
1256  vol->SetVisDaughters(selected);
1257  }
1258 
1259  int id{0};
1260  for (auto &desc: fDesc)
1261  if (fNodes[id++]->GetVolume() == vol)
1262  desc.vis = dnode.vis;
1263 
1264  ClearDrawData(); // after change raw data is no longer valid
1265 
1266  return true;
1267 }
1268 
1269 /////////////////////////////////////////////////////////////////////////////////
1270 /// Change visibility for specified element
1271 /// Returns true if changes was performed
1272 
1273 std::unique_ptr<ROOT::Experimental::REveGeomNodeInfo> ROOT::Experimental::REveGeomDescription::MakeNodeInfo(const std::string &path)
1274 {
1275  std::unique_ptr<REveGeomNodeInfo> res;
1276 
1277  RGeomBrowserIter iter(*this);
1278 
1279  if (iter.Navigate(path)) {
1280 
1281  auto node = fNodes[iter.GetNodeId()];
1282 
1283  auto &desc = fDesc[iter.GetNodeId()];
1284 
1285  res = std::make_unique<REveGeomNodeInfo>();
1286 
1287  res->fullpath = path;
1288  res->node_name = node->GetName();
1289  res->node_type = node->ClassName();
1290 
1291  TGeoShape *shape = node->GetVolume() ? node->GetVolume()->GetShape() : nullptr;
1292 
1293  if (shape) {
1294  res->shape_name = shape->GetName();
1295  res->shape_type = shape->ClassName();
1296  }
1297 
1298  if (shape && desc.CanDisplay()) {
1299 
1300  auto &shape_descr = MakeShapeDescr(shape);
1301 
1302  res->ri = shape_descr.rndr_info(); // temporary pointer, can be used preserved for short time
1303  }
1304  }
1305 
1306  return res;
1307 }
1308 
1309 /////////////////////////////////////////////////////////////////////////////////
1310 /// Change configuration by client
1311 /// Returns true if any parameter was really changed
1312 
1313 bool ROOT::Experimental::REveGeomDescription::ChangeConfiguration(const std::string &json)
1314 {
1315  auto cfg = TBufferJSON::FromJSON<REveGeomConfig>(json);
1316  if (!cfg) return false;
1317 
1318  auto json1 = TBufferJSON::ToJSON(cfg.get());
1319  auto json2 = TBufferJSON::ToJSON(&fCfg);
1320 
1321  if (json1 == json2)
1322  return false;
1323 
1324  fCfg = *cfg; // use assign
1325 
1326  ClearDrawData();
1327 
1328  return true;
1329 }
1330 
1331