Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
REveElement.cxx
Go to the documentation of this file.
1 // @(#)root/eve7:$Id$
2 // Authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007, 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/REveElement.hxx>
13 #include <ROOT/REveUtil.hxx>
14 #include <ROOT/REveScene.hxx>
15 #include <ROOT/REveCompound.hxx>
16 #include <ROOT/REveTrans.hxx>
17 #include <ROOT/REveManager.hxx>
18 #include <ROOT/REveSelection.hxx>
21 #include <ROOT/REveRenderData.hxx>
22 
23 #include "TGeoMatrix.h"
24 
25 #include "TClass.h"
26 #include "TPRegexp.h"
27 #include "TROOT.h"
28 #include "TColor.h"
29 
30 #include "json.hpp"
31 #include <cassert>
32 
33 
34 #include <algorithm>
35 
36 using namespace ROOT::Experimental;
37 namespace REX = ROOT::Experimental;
38 
39 /** \class REveElement
40 \ingroup REve
41 Base class for REveUtil visualization elements, providing hierarchy
42 management, rendering control and list-tree item management.
43 
44 Class of acceptable children can be limited by setting the
45 fChildClass member.
46 */
47 
48 ////////////////////////////////////////////////////////////////////////////////
49 /// Default constructor.
50 
51 REveElement::REveElement(const std::string& name, const std::string& title) :
52  fName (name),
53  fTitle (title)
54 {
55 }
56 
57 ////////////////////////////////////////////////////////////////////////////////
58 /// Copy constructor. Does shallow copy.
59 /// For deep-cloning and children-cloning, see:
60 /// ~~~ {.cpp}
61 /// REveElement* CloneElementRecurse(Int_t level)
62 /// void CloneChildrenRecurse(REveElement* dest, Int_t level)
63 /// ~~~
64 /// 'void* UserData' is NOT copied.
65 /// If the element is projectable, its projections are NOT copied.
66 ///
67 /// Not implemented for most sub-classes, let us know.
68 /// Note that sub-classes of REveProjected are NOT and will NOT be copyable.
69 
70 REveElement::REveElement(const REveElement& e) :
71  fName (e.fName),
72  fTitle (e.fTitle),
73  fChildClass (e.fChildClass),
74  fVizTag (e.fVizTag),
75  fDestroyOnZeroRefCnt (e.fDestroyOnZeroRefCnt),
76  fRnrSelf (e.fRnrSelf),
77  fRnrChildren (e.fRnrChildren),
78  fCanEditMainColor (e.fCanEditMainColor),
79  fCanEditMainTransparency(e.fCanEditMainTransparency),
80  fCanEditMainTrans (e.fCanEditMainTrans),
81  fMainTransparency (e.fMainTransparency),
82  fPickable (e.fPickable),
83  fCSCBits (e.fCSCBits)
84 {
85  SetVizModel(e.fVizModel);
86  // FIXME: from Sergey: one have to use other way to referencing main color
87  if (e.fMainColorPtr)
88  fMainColorPtr = (Color_t*)((const char*) this + ((const char*) e.fMainColorPtr - (const char*) &e));
89  if (e.fMainTrans)
90  fMainTrans = std::make_unique<REveTrans>(*e.fMainTrans.get());
91 }
92 
93 ////////////////////////////////////////////////////////////////////////////////
94 /// Destructor. Do not call this method directly, either call Destroy() or
95 /// Annihilate(). See also DestroyElements() and AnnihilateElements() if you
96 /// need to delete all children of an element.
97 
98 REveElement::~REveElement()
99 {
100  if (fDestructing != kAnnihilate)
101  {
102  fDestructing = kStandard;
103  RemoveElementsInternal();
104 
105  if (fMother) {
106  fMother->RemoveElementLocal(this);
107  fMother->fChildren.remove(this);
108  }
109 
110  if (fScene) {
111  fScene->SceneElementRemoved( fElementId);
112  }
113 
114  for (auto &au : fAunts)
115  {
116  au->RemoveNieceInternal(this);
117  }
118  }
119 }
120 
121 ElementId_t REveElement::get_mother_id() const
122 {
123  return fMother ? fMother->GetElementId() : 0;
124 }
125 
126 ElementId_t REveElement::get_scene_id() const
127 {
128  return fScene ? fScene->GetElementId() : 0;
129 }
130 
131 void REveElement::assign_element_id_recurisvely()
132 {
133  assert(fElementId == 0);
134 
135  REX::gEve->AssignElementId(this);
136  for (auto &c : fChildren)
137  c->assign_element_id_recurisvely();
138 }
139 
140 void REveElement::assign_scene_recursively(REveScene* s)
141 {
142  assert(fScene == nullptr);
143 
144  fScene = s;
145 
146  // XXX MT -- Why do we have fDestructing here? Can this really happen?
147  // If yes, shouldn't we block it in AddElement() already?
148  if (fDestructing == kNone && fScene && fScene->IsAcceptingChanges())
149  {
150  StampElementAdded();
151  }
152  for (auto &c : fChildren)
153  c->assign_scene_recursively(s);
154 }
155 
156 ////////////////////////////////////////////////////////////////////////////////
157 /// Called before the element is deleted, thus offering the last chance
158 /// to detach from acquired resources and from the framework itself.
159 /// Here the request is just passed to REveManager.
160 /// If you override it, make sure to call base-class version.
161 
162 void REveElement::PreDeleteElement()
163 {
164  if (fElementId != 0)
165  {
166  REX::gEve->PreDeleteElement(this);
167  }
168 }
169 
170 ////////////////////////////////////////////////////////////////////////////////
171 /// Clone the element via copy constructor.
172 /// Should be implemented for all classes that require cloning support.
173 
174 REveElement* REveElement::CloneElement() const
175 {
176  return new REveElement(*this);
177 }
178 
179 ////////////////////////////////////////////////////////////////////////////////
180 /// Clone elements and recurse 'level' deep over children.
181 /// - If level == 0, only the element itself is cloned (default).
182 /// - If level == -1, all the hierarchy is cloned.
183 
184 REveElement* REveElement::CloneElementRecurse(Int_t level) const
185 {
186  REveElement* el = CloneElement();
187  if (level--)
188  {
189  CloneChildrenRecurse(el, level);
190  }
191  return el;
192 }
193 
194 ////////////////////////////////////////////////////////////////////////////////
195 /// Clone children and attach them to the dest element.
196 /// If level == 0, only the direct descendants are cloned (default).
197 /// If level == -1, all the hierarchy is cloned.
198 
199 void REveElement::CloneChildrenRecurse(REveElement* dest, Int_t level) const
200 {
201  for (auto &c: fChildren)
202  dest->AddElement(c->CloneElementRecurse(level));
203 }
204 
205 
206 //==============================================================================
207 
208 ////////////////////////////////////////////////////////////////////////////////
209 /// Set name of an element.
210 
211 void REveElement::SetName(const std::string& name)
212 {
213  fName = name;
214  NameTitleChanged();
215 }
216 
217 ////////////////////////////////////////////////////////////////////////////////
218 /// Set title of an element.
219 
220 void REveElement::SetTitle(const std::string& title)
221 {
222  fTitle = title;
223  NameTitleChanged();
224 }
225 
226 ////////////////////////////////////////////////////////////////////////////////
227 /// Set name and title of an element.
228 /// Here we attempt to cast the assigned object into TNamed and call
229 /// SetNameTitle() there.
230 /// If you override this call NameTitleChanged() from there.
231 
232 void REveElement::SetNameTitle(const std::string& name, const std::string& title)
233 {
234  fName = name;
235  fTitle = title;
236  NameTitleChanged();
237 }
238 
239 ////////////////////////////////////////////////////////////////////////////////
240 /// Virtual function called when a name or title of the element has
241 /// been changed.
242 /// If you override this, call also the version of your direct base-class.
243 
244 void REveElement::NameTitleChanged()
245 {
246  // Should send out some message. Need a new stamp type?
247 }
248 
249 ////////////////////////////////////////////////////////////////////////////////
250 /// Set visualization-parameter model element.
251 /// Calling of this function from outside of EVE should in principle
252 /// be avoided as it can lead to dis-synchronization of viz-tag and
253 /// viz-model.
254 
255 void REveElement::SetVizModel(REveElement* model)
256 {
257  fVizModel = model;
258 }
259 
260 ////////////////////////////////////////////////////////////////////////////////
261 /// Find model element in VizDB that corresponds to previously
262 /// assigned fVizTag and set fVizModel accordingly.
263 /// If the tag is not found in VizDB, the old model-element is kept
264 /// and false is returned.
265 
266 Bool_t REveElement::SetVizModelByTag()
267 {
268  REveElement* model = REX::gEve->FindVizDBEntry(fVizTag);
269  if (model)
270  {
271  SetVizModel(model);
272  return kTRUE;
273  }
274  else
275  {
276  return kFALSE;
277  }
278 }
279 
280 ////////////////////////////////////////////////////////////////////////////////
281 /// Set the VizTag, find model-element from the VizDB and copy
282 /// visualization-parameters from it. If the model is not found and
283 /// fallback_tag is non-null, search for it is attempted as well.
284 /// For example: ApplyVizTag("TPC Clusters", "Clusters");
285 ///
286 /// If the model-element can not be found a warning is printed and
287 /// false is returned.
288 
289 Bool_t REveElement::ApplyVizTag(const TString& tag, const TString& fallback_tag)
290 {
291  REveElement* model;
292 
293  if ((model = REX::gEve->FindVizDBEntry(tag)) != nullptr)
294  {
295  SetVizTag(tag);
296  }
297  else if ( ! fallback_tag.IsNull() && (model = REX::gEve->FindVizDBEntry(fallback_tag)) != nullptr)
298  {
299  SetVizTag(fallback_tag);
300  }
301 
302  if (model)
303  {
304  SetVizModel(model);
305  CopyVizParamsFromDB();
306  return true;
307  }
308  Warning("REveElement::ApplyVizTag", "entry for tag '%s' not found in VizDB.", tag.Data());
309  return false;
310 }
311 
312 ////////////////////////////////////////////////////////////////////////////////
313 /// Propagate visualization parameters to dependent elements.
314 ///
315 /// MainColor is propagated independently in SetMainColor().
316 /// In this case, as fMainColor is a pointer to Color_t, it should
317 /// be set in TProperClass::CopyVizParams().
318 ///
319 /// Render state is not propagated. Maybe it should be, at least optionally.
320 
321 void REveElement::PropagateVizParamsToProjecteds()
322 {
323  REveProjectable* pable = dynamic_cast<REveProjectable*>(this);
324  if (pable && pable->HasProjecteds())
325  {
326  pable->PropagateVizParams();
327  }
328 }
329 
330 ////////////////////////////////////////////////////////////////////////////////
331 /// Propagate visualization parameters from element el (defaulting
332 /// to this) to all children.
333 ///
334 /// The primary use of this is for model-elements from
335 /// visualization-parameter database.
336 
337 void REveElement::PropagateVizParamsToChildren(REveElement* el)
338 {
339  if (!el) el = this;
340 
341  for (auto &c : fChildren)
342  {
343  c->CopyVizParams(el);
344  }
345 }
346 
347 ////////////////////////////////////////////////////////////////////////////////
348 /// Copy visualization parameters from element el.
349 /// This method needs to be overriden by any class that introduces
350 /// new parameters.
351 
352 void REveElement::CopyVizParams(const REveElement* el)
353 {
354  fCanEditMainColor = el->fCanEditMainColor;
355  fCanEditMainTransparency = el->fCanEditMainTransparency;
356  fMainTransparency = el->fMainTransparency;
357  if (fMainColorPtr == & fDefaultColor) fDefaultColor = el->GetMainColor();
358 
359  AddStamp(kCBColorSelection | kCBObjProps);
360 }
361 
362 ////////////////////////////////////////////////////////////////////////////////
363 /// Copy visualization parameters from the model-element fVizModel.
364 /// A warning is printed if the model-element fVizModel is not set.
365 
366 void REveElement::CopyVizParamsFromDB()
367 {
368  if (fVizModel)
369  {
370  CopyVizParams(fVizModel);
371  }
372  else
373  {
374  Warning("REveElement::CopyVizParamsFromDB", "VizModel has not been set.");
375  }
376 }
377 
378 ////////////////////////////////////////////////////////////////////////////////
379 /// Save visualization parameters for this element with given tag.
380 ///
381 /// This function creates the instantiation code, calls virtual
382 /// WriteVizParams() and, at the end, writes out the code for
383 /// registration of the model into the VizDB.
384 
385 void REveElement::SaveVizParams(std::ostream& out, const TString& tag, const TString& var)
386 {
387  static const REveException eh("REveElement::GetObject ");
388 
389  TString t = " ";
390  TString cls(IsA()->GetName());
391 
392  out << "\n";
393 
394  TString intro = " TAG='" + tag + "', CLASS='" + cls + "'";
395  out << " //" << intro << "\n";
396  out << " //" << TString('-', intro.Length()) << "\n";
397  out << t << cls << "* " << var <<" = new " << cls << ";\n";
398 
399  WriteVizParams(out, var);
400 
401  out << t << "REX::gEve->InsertVizDBEntry(\"" << tag << "\", "<< var <<");\n";
402 }
403 
404 ////////////////////////////////////////////////////////////////////////////////
405 /// Write-out visual parameters for this object.
406 /// This is a virtual function and all sub-classes are required to
407 /// first call the base-element version.
408 /// The name of the element pointer is 'x%03d', due to cint limitations.
409 /// Three spaces should be used for indentation, same as in
410 /// SavePrimitive() methods.
411 
412 void REveElement::WriteVizParams(std::ostream& out, const TString& var)
413 {
414  TString t = " " + var + "->";
415 
416  out << t << "SetElementName(\"" << fName << "\");\n";
417  out << t << "SetElementTitle(\"" << fTitle << "\");\n";
418  out << t << "SetEditMainColor(" << fCanEditMainColor << ");\n";
419  out << t << "SetEditMainTransparency(" << fCanEditMainTransparency << ");\n";
420  out << t << "SetMainTransparency(" << fMainTransparency << ");\n";
421 }
422 
423 ////////////////////////////////////////////////////////////////////////////////
424 /// Set visual parameters for this object for given tag.
425 
426 void REveElement::VizDB_Apply(const std::string& tag)
427 {
428  if (ApplyVizTag(tag))
429  {
430  PropagateVizParamsToProjecteds();
431  REX::gEve->Redraw3D();
432  }
433 }
434 
435 ////////////////////////////////////////////////////////////////////////////////
436 /// Reset visual parameters for this object from VizDB.
437 /// The model object must be already set.
438 
439 void REveElement::VizDB_Reapply()
440 {
441  if (fVizModel)
442  {
443  CopyVizParamsFromDB();
444  PropagateVizParamsToProjecteds();
445  REX::gEve->Redraw3D();
446  }
447 }
448 
449 ////////////////////////////////////////////////////////////////////////////////
450 /// Copy visual parameters from this element to viz-db model.
451 /// If update is set, all clients of the model will be updated to
452 /// the new value.
453 /// A warning is printed if the model-element fVizModel is not set.
454 
455 void REveElement::VizDB_UpdateModel(Bool_t update)
456 {
457  if (fVizModel)
458  {
459  fVizModel->CopyVizParams(this);
460  if (update)
461  {
462  // XXX Back references from vizdb templates have been removed in Eve7.
463  // XXX We could traverse all scenes and elementes and reset those that
464  // XXX have a matching fVizModel. Or something.
465  Error("VizDB_UpdateModel", "update from vizdb -> elements not implemented.");
466  // fVizModel->PropagateVizParamsToElements(fVizModel);
467  // REX::gEve->Redraw3D();
468  }
469  }
470  else
471  {
472  Warning("VizDB_UpdateModel", "VizModel has not been set.");
473  }
474 }
475 
476 ////////////////////////////////////////////////////////////////////////////////
477 /// Create a replica of element and insert it into VizDB with given tag.
478 /// If replace is true an existing element with the same tag will be replaced.
479 /// If update is true, existing client of tag will be updated.
480 
481 void REveElement::VizDB_Insert(const std::string& tag, Bool_t replace, Bool_t update)
482 {
483  static const REveException eh("REveElement::GetObject ");
484 
485  TClass* cls = IsA();
486  REveElement* el = reinterpret_cast<REveElement*>(cls->New());
487  if (!el) {
488  Error("VizDB_Insert", "Creation of replica failed.");
489  return;
490  }
491  el->CopyVizParams(this);
492  Bool_t succ = REX::gEve->InsertVizDBEntry(tag, el, replace, update);
493  if (succ && update)
494  REX::gEve->Redraw3D();
495 }
496 
497 ////////////////////////////////////////////////////////////////////////////////
498 /// Add el into the list aunts.
499 ///
500 /// Adding aunt is subordinate to adding a niece.
501 /// This is an internal function.
502 
503 void REveElement::AddAunt(REveAunt *au)
504 {
505  assert(au != nullptr);
506 
507  fAunts.emplace_back(au);
508 }
509 
510 ////////////////////////////////////////////////////////////////////////////////
511 /// Remove el from the list of aunts.
512 /// Removing aunt is subordinate to removing a niece.
513 /// This is an internal function.
514 
515 void REveElement::RemoveAunt(REveAunt *au)
516 {
517  assert(au != nullptr);
518 
519  fAunts.remove(au);
520 }
521 
522 /******************************************************************************/
523 
524 ////////////////////////////////////////////////////////////////////////////////
525 /// Check external references to this and eventually auto-destruct
526 /// the render-element.
527 
528 void REveElement::CheckReferenceCount(const std::string& from)
529 {
530  if (fDestructing != kNone)
531  return;
532 
533  if (fMother == nullptr && fDestroyOnZeroRefCnt && fDenyDestroy <= 0)
534  {
535  if (gDebug > 0)
536  Info("REveElement::CheckReferenceCount", "(called from %s) auto-destructing '%s' on zero reference count.",
537  from.c_str(), GetCName());
538 
539  PreDeleteElement();
540  delete this;
541  }
542 }
543 
544 ////////////////////////////////////////////////////////////////////////////////
545 /// Return class for this element
546 
547 TClass *REveElement::IsA() const
548 {
549  return TClass::GetClass(typeid(*this), kTRUE, kTRUE);
550 }
551 
552 ////////////////////////////////////////////////////////////////////////////////
553 /// Export render-element to CINT with variable name var_name.
554 
555 void REveElement::ExportToCINT(const char *var_name)
556 {
557  const char* cname = IsA()->GetName();
558  gROOT->ProcessLine(TString::Format("%s* %s = (%s*)0x%lx;", cname, var_name, cname, (ULong_t)this));
559 }
560 
561 ////////////////////////////////////////////////////////////////////////////////
562 /// Set render state of this element, i.e. if it will be published
563 /// on next scene update pass.
564 /// Returns true if the state has changed.
565 
566 Bool_t REveElement::SetRnrSelf(Bool_t rnr)
567 {
568  if (SingleRnrState())
569  {
570  return SetRnrState(rnr);
571  }
572 
573  if (rnr != fRnrSelf)
574  {
575  fRnrSelf = rnr;
576  StampVisibility();
577  PropagateRnrStateToProjecteds();
578 
579  return kTRUE;
580  }
581  return kFALSE;
582 }
583 
584 ////////////////////////////////////////////////////////////////////////////////
585 /// Set render state of this element's children, i.e. if they will
586 /// be published on next scene update pass.
587 /// Returns true if the state has changed.
588 
589 Bool_t REveElement::SetRnrChildren(Bool_t rnr)
590 {
591  if (SingleRnrState())
592  {
593  return SetRnrState(rnr);
594  }
595 
596  if (rnr != fRnrChildren)
597  {
598  fRnrChildren = rnr;
599  StampVisibility();
600  PropagateRnrStateToProjecteds();
601  return kTRUE;
602  }
603  return kFALSE;
604 }
605 
606 ////////////////////////////////////////////////////////////////////////////////
607 /// Set state for rendering of this element and its children.
608 /// Returns true if the state has changed.
609 
610 Bool_t REveElement::SetRnrSelfChildren(Bool_t rnr_self, Bool_t rnr_children)
611 {
612  if (SingleRnrState())
613  {
614  return SetRnrState(rnr_self);
615  }
616 
617  if (fRnrSelf != rnr_self || fRnrChildren != rnr_children)
618  {
619  fRnrSelf = rnr_self;
620  fRnrChildren = rnr_children;
621  StampVisibility();
622  PropagateRnrStateToProjecteds();
623  return kTRUE;
624  }
625  return kFALSE;
626 }
627 
628 ////////////////////////////////////////////////////////////////////////////////
629 /// Set render state of this element and of its children to the same
630 /// value.
631 /// Returns true if the state has changed.
632 
633 Bool_t REveElement::SetRnrState(Bool_t rnr)
634 {
635  if (fRnrSelf != rnr || fRnrChildren != rnr)
636  {
637  fRnrSelf = fRnrChildren = rnr;
638  StampVisibility();
639  PropagateRnrStateToProjecteds();
640  return kTRUE;
641  }
642  return kFALSE;
643 }
644 
645 ////////////////////////////////////////////////////////////////////////////////
646 /// Propagate render state to the projected replicas of this element.
647 /// Maybe this should be optional on REX::gEve/element level.
648 
649 void REveElement::PropagateRnrStateToProjecteds()
650 {
651  REveProjectable *pable = dynamic_cast<REveProjectable*>(this);
652  if (pable && pable->HasProjecteds())
653  {
654  pable->PropagateRenderState(fRnrSelf, fRnrChildren);
655  }
656 }
657 
658 ////////////////////////////////////////////////////////////////////////////////
659 /// Set up element to use built-in main color and set flags allowing editing
660 /// of main color and transparency.
661 
662 void REveElement::SetupDefaultColorAndTransparency(Color_t col, Bool_t can_edit_color, Bool_t can_edit_transparency)
663 {
664  fMainColorPtr = & fDefaultColor;
665  fDefaultColor = col;
666  fCanEditMainColor = can_edit_color;
667  fCanEditMainTransparency = can_edit_transparency;
668 }
669 
670 ////////////////////////////////////////////////////////////////////////////////
671 /// Set main color of the element.
672 
673 void REveElement::SetMainColor(Color_t color)
674 {
675  Color_t old_color = GetMainColor();
676 
677  if (fMainColorPtr)
678  {
679  *fMainColorPtr = color;
680  StampColorSelection();
681  }
682 
683  PropagateMainColorToProjecteds(color, old_color);
684 }
685 
686 ////////////////////////////////////////////////////////////////////////////////
687 /// Convert pixel to Color_t and call SetMainColor().
688 
689 void REveElement::SetMainColorPixel(Pixel_t pixel)
690 {
691  SetMainColor(TColor::GetColor(pixel));
692 }
693 
694 ////////////////////////////////////////////////////////////////////////////////
695 /// Convert RGB values to Color_t and call SetMainColor.
696 
697 void REveElement::SetMainColorRGB(UChar_t r, UChar_t g, UChar_t b)
698 {
699  SetMainColor(TColor::GetColor(r, g, b));
700 }
701 
702 ////////////////////////////////////////////////////////////////////////////////
703 /// Convert RGB values to Color_t and call SetMainColor.
704 
705 void REveElement::SetMainColorRGB(Float_t r, Float_t g, Float_t b)
706 {
707  SetMainColor(TColor::GetColor(r, g, b));
708 }
709 
710 ////////////////////////////////////////////////////////////////////////////////
711 /// Propagate color to projected elements.
712 
713 void REveElement::PropagateMainColorToProjecteds(Color_t color, Color_t old_color)
714 {
715  REveProjectable* pable = dynamic_cast<REveProjectable*>(this);
716  if (pable && pable->HasProjecteds())
717  {
718  pable->PropagateMainColor(color, old_color);
719  }
720 }
721 
722 ////////////////////////////////////////////////////////////////////////////////
723 /// Set main-transparency.
724 /// Transparency is clamped to [0, 100].
725 
726 void REveElement::SetMainTransparency(Char_t t)
727 {
728  Char_t old_t = GetMainTransparency();
729 
730  if (t > 100) t = 100;
731  fMainTransparency = t;
732  StampColorSelection();
733 
734  PropagateMainTransparencyToProjecteds(t, old_t);
735 }
736 
737 ////////////////////////////////////////////////////////////////////////////////
738 /// Set main-transparency via float alpha variable.
739 /// Value of alpha is clamped t0 [0, 1].
740 
741 void REveElement::SetMainAlpha(Float_t alpha)
742 {
743  if (alpha < 0) alpha = 0;
744  if (alpha > 1) alpha = 1;
745  SetMainTransparency((Char_t) (100.0f*(1.0f - alpha)));
746 }
747 
748 ////////////////////////////////////////////////////////////////////////////////
749 /// Propagate transparency to projected elements.
750 
751 void REveElement::PropagateMainTransparencyToProjecteds(Char_t t, Char_t old_t)
752 {
753  REveProjectable* pable = dynamic_cast<REveProjectable*>(this);
754  if (pable && pable->HasProjecteds())
755  {
756  pable->PropagateMainTransparency(t, old_t);
757  }
758 }
759 
760 ////////////////////////////////////////////////////////////////////////////////
761 /// Return pointer to main transformation. If 'create' flag is set (default)
762 /// it is created if not yet existing.
763 
764 REveTrans *REveElement::PtrMainTrans(Bool_t create)
765 {
766  if (!fMainTrans && create)
767  InitMainTrans();
768 
769  return fMainTrans.get();
770 }
771 
772 ////////////////////////////////////////////////////////////////////////////////
773 /// Return reference to main transformation. It is created if not yet
774 /// existing.
775 
776 REveTrans &REveElement::RefMainTrans()
777 {
778  if (!fMainTrans)
779  InitMainTrans();
780 
781  return *fMainTrans.get();
782 }
783 
784 ////////////////////////////////////////////////////////////////////////////////
785 /// Initialize the main transformation to identity matrix.
786 /// If can_edit is true (default), the user will be able to edit the
787 /// transformation parameters via GUI.
788 
789 void REveElement::InitMainTrans(Bool_t can_edit)
790 {
791  if (fMainTrans)
792  fMainTrans->UnitTrans();
793  else
794  fMainTrans = std::make_unique<REveTrans>();
795  fCanEditMainTrans = can_edit;
796 }
797 
798 ////////////////////////////////////////////////////////////////////////////////
799 /// Destroy the main transformation matrix, it will always be taken
800 /// as identity. Editing of transformation parameters is disabled.
801 
802 void REveElement::DestroyMainTrans()
803 {
804  fMainTrans.reset(nullptr);
805  fCanEditMainTrans = kFALSE;
806 }
807 
808 ////////////////////////////////////////////////////////////////////////////////
809 /// Set transformation matrix from column-major array.
810 
811 void REveElement::SetTransMatrix(Double_t* carr)
812 {
813  RefMainTrans().SetFrom(carr);
814 }
815 
816 ////////////////////////////////////////////////////////////////////////////////
817 /// Set transformation matrix from TGeo's matrix.
818 
819 void REveElement::SetTransMatrix(const TGeoMatrix& mat)
820 {
821  RefMainTrans().SetFrom(mat);
822 }
823 
824 ////////////////////////////////////////////////////////////////////////////////
825 /// Check if el can be added to this element.
826 ///
827 /// Here we make sure the new child is not equal to this and, if fChildClass
828 /// is set, that it is inherited from it.
829 
830 Bool_t REveElement::AcceptElement(REveElement* el)
831 {
832  if (el == this)
833  return kFALSE;
834  if (fChildClass && ! el->IsA()->InheritsFrom(fChildClass))
835  return kFALSE;
836  return kTRUE;
837 }
838 
839 ////////////////////////////////////////////////////////////////////////////////
840 /// Add el to the list of children.
841 
842 void REveElement::AddElement(REveElement *el)
843 {
844  static const REveException eh("REveElement::AddElement ");
845 
846  if (!el) throw eh + "called with nullptr argument.";
847  if ( ! AcceptElement(el)) throw eh + Form("parent '%s' rejects '%s'.", GetCName(), el->GetCName());
848  if (el->fElementId) throw eh + "element already has an id.";
849  // if (el->fScene) throw eh + "element already has a Scene.";
850  if (el->fMother) throw eh + "element already has a Mother.";
851 
852  // XXX Implement reparent --> MoveElement() ????
853  // Actually, better to do new = old.Clone(), RemoveElement(old), AddElement(new);
854  // Or do magick with Inc/DecDenyDestroy().
855  // PITA with existing children !!!! Need to re-scene them.
856 
857  if (fElementId) el->assign_element_id_recurisvely();
858  if (fScene && ! el->fScene) el->assign_scene_recursively(fScene);
859 
860  el->fMother = this;
861 
862  fChildren.emplace_back(el);
863 
864  // XXXX This should be element added. Also, should be different for
865  // "full (re)construction". Scenes should manage that and have
866  // state like: none - constructing - clearing - nominal - updating.
867  // I recon this means an element should have a ptr to its scene.
868  //
869  // ElementChanged();
870 }
871 
872 ////////////////////////////////////////////////////////////////////////////////
873 /// Remove el from the list of children.
874 
875 void REveElement::RemoveElement(REveElement* el)
876 {
877  static const REveException eh("REveElement::RemoveElement ");
878 
879  if (!el) throw eh + "called with nullptr argument.";
880  if (el->fMother != this) throw eh + "this element is not mother of el.";
881 
882  RemoveElementLocal(el);
883 
884  el->fScene->SceneElementRemoved(fElementId);
885  el->fMother = nullptr;
886  el->fScene = nullptr;
887 
888  el->CheckReferenceCount();
889 
890  fChildren.remove(el);
891 
892  // XXXX This should be ElementRemoved(). Also, think about recursion, deletion etc.
893  // Also, this seems to be done above, in the call to fScene.
894  //
895  // ElementChanged();
896 }
897 
898 ////////////////////////////////////////////////////////////////////////////////
899 /// Perform additional local removal of el.
900 /// Called from RemoveElement() which does whole untangling.
901 /// Put into special function as framework-related handling of
902 /// element removal should really be common to all classes and
903 /// clearing of local structures happens in between removal
904 /// of list-tree-items and final removal.
905 /// If you override this, you should also override
906 /// RemoveElementsLocal().
907 
908 void REveElement::RemoveElementLocal(REveElement* /*el*/)
909 {
910 }
911 
912 ////////////////////////////////////////////////////////////////////////////////
913 /// Remove all elements. This assumes removing of all elements can
914 /// be done more efficiently then looping over them and removing one
915 /// by one. This protected function performs the removal on the
916 /// level of REveElement.
917 
918 void REveElement::RemoveElementsInternal()
919 {
920  RemoveElementsLocal();
921 
922  for (auto &c : fChildren)
923  {
924  c->fScene->SceneElementRemoved(c->fElementId);
925  c->fMother = nullptr;
926  c->fScene = nullptr;
927 
928  c->CheckReferenceCount();
929  }
930 
931  fChildren.clear();
932 }
933 
934 ////////////////////////////////////////////////////////////////////////////////
935 /// Remove all elements. This assumes removing of all elements can
936 /// be done more efficiently then looping over them and removing
937 /// them one by one.
938 
939 void REveElement::RemoveElements()
940 {
941  if (HasChildren())
942  {
943  RemoveElementsInternal();
944  // ElementChanged();
945  }
946 }
947 
948 ////////////////////////////////////////////////////////////////////////////////
949 /// Perform additional local removal of all elements.
950 /// See comment to RemoveElementLocal(REveElement*).
951 
952 void REveElement::RemoveElementsLocal()
953 {
954 }
955 
956 ////////////////////////////////////////////////////////////////////////////////
957 /// If this is a projectable, loop over all projected replicas and
958 /// add the projected image of child 'el' there. This is supposed to
959 /// be called after you add a child to a projectable after it has
960 /// already been projected.
961 /// You might also want to call RecheckImpliedSelections() on this
962 /// element or 'el'.
963 ///
964 /// If 'same_depth' flag is true, the same depth as for parent object
965 /// is used in every projection. Otherwise current depth of each
966 /// relevant projection-manager is used.
967 
968 void REveElement::ProjectChild(REveElement* el, Bool_t same_depth)
969 {
970  REveProjectable* pable = dynamic_cast<REveProjectable*>(this);
971  if (pable && HasChild(el))
972  {
973  for (auto &pp: pable->RefProjecteds())
974  {
975  auto pmgr = pp->GetManager();
976  Float_t cd = pmgr->GetCurrentDepth();
977  if (same_depth) pmgr->SetCurrentDepth(pp->GetDepth());
978 
979  pmgr->SubImportElements(el, pp->GetProjectedAsElement());
980 
981  if (same_depth) pmgr->SetCurrentDepth(cd);
982  }
983  }
984 }
985 
986 ////////////////////////////////////////////////////////////////////////////////
987 /// If this is a projectable, loop over all projected replicas and
988 /// add the projected image of all children there. This is supposed
989 /// to be called after you destroy all children and then add new
990 /// ones after this element has already been projected.
991 /// You might also want to call RecheckImpliedSelections() on this
992 /// element.
993 ///
994 /// If 'same_depth' flag is true, the same depth as for the
995 /// projected element is used in every projection. Otherwise current
996 /// depth of each relevant projection-manager is used.
997 
998 void REveElement::ProjectAllChildren(Bool_t same_depth)
999 {
1000  REveProjectable* pable = dynamic_cast<REveProjectable*>(this);
1001  if (pable)
1002  {
1003  for (auto &pp: pable->RefProjecteds())
1004  {
1005  REveProjectionManager *pmgr = pp->GetManager();
1006  Float_t cd = pmgr->GetCurrentDepth();
1007  if (same_depth) pmgr->SetCurrentDepth(pp->GetDepth());
1008 
1009  pmgr->SubImportChildren(this, pp->GetProjectedAsElement());
1010 
1011  if (same_depth) pmgr->SetCurrentDepth(cd);
1012  }
1013  }
1014 }
1015 
1016 ////////////////////////////////////////////////////////////////////////////////
1017 /// Check if element el is a child of this element.
1018 
1019 Bool_t REveElement::HasChild(REveElement* el)
1020 {
1021  return (std::find(fChildren.begin(), fChildren.end(), el) != fChildren.end());
1022 }
1023 
1024 ////////////////////////////////////////////////////////////////////////////////
1025 /// Find the first child with given name. If cls is specified (non
1026 /// 0), it is also checked.
1027 ///
1028 /// Returns nullptr if not found.
1029 
1030 REveElement *REveElement::FindChild(const TString& name, const TClass* cls)
1031 {
1032  for (auto &c: fChildren)
1033  {
1034  if (name.CompareTo(c->GetCName()) == 0)
1035  {
1036  if (!cls || c->IsA()->InheritsFrom(cls))
1037  return c;
1038  }
1039  }
1040  return nullptr;
1041 }
1042 
1043 ////////////////////////////////////////////////////////////////////////////////
1044 /// Find the first child whose name matches regexp. If cls is
1045 /// specified (non 0), it is also checked.
1046 ///
1047 /// Returns nullptr if not found.
1048 
1049 REveElement* REveElement::FindChild(TPRegexp &regexp, const TClass *cls)
1050 {
1051  for (auto &c: fChildren)
1052  {
1053  if (regexp.MatchB(c->GetName()))
1054  {
1055  if (!cls || c->IsA()->InheritsFrom(cls))
1056  return c;
1057  }
1058  }
1059  return nullptr;
1060 }
1061 
1062 ////////////////////////////////////////////////////////////////////////////////
1063 /// Find all children with given name and append them to matches
1064 /// list. If class is specified (non 0), it is also checked.
1065 ///
1066 /// Returns number of elements added to the list.
1067 
1068 Int_t REveElement::FindChildren(List_t& matches,
1069  const TString& name, const TClass* cls)
1070 {
1071  Int_t count = 0;
1072  for (auto &c: fChildren)
1073  {
1074  if (name.CompareTo(c->GetCName()) == 0)
1075  {
1076  if (!cls || c->IsA()->InheritsFrom(cls))
1077  {
1078  matches.push_back(c);
1079  ++count;
1080  }
1081  }
1082  }
1083  return count;
1084 }
1085 
1086 ////////////////////////////////////////////////////////////////////////////////
1087 /// Find all children whose name matches regexp and append them to
1088 /// matches list.
1089 ///
1090 /// Returns number of elements added to the list.
1091 
1092 Int_t REveElement::FindChildren(List_t &matches,
1093  TPRegexp &regexp, const TClass *cls)
1094 {
1095  Int_t count = 0;
1096  for (auto &c : fChildren)
1097  {
1098  if (regexp.MatchB(c->GetCName()))
1099  {
1100  if (!cls || c->IsA()->InheritsFrom(cls))
1101  {
1102  matches.push_back(c);
1103  ++count;
1104  }
1105  }
1106  }
1107  return count;
1108 }
1109 
1110 ////////////////////////////////////////////////////////////////////////////////
1111 /// Returns the first child element or 0 if the list is empty.
1112 
1113 REveElement* REveElement::FirstChild() const
1114 {
1115  return HasChildren() ? fChildren.front() : nullptr;
1116 }
1117 
1118 ////////////////////////////////////////////////////////////////////////////////
1119 /// Returns the last child element or 0 if the list is empty.
1120 
1121 REveElement* REveElement::LastChild () const
1122 {
1123  return HasChildren() ? fChildren.back() : nullptr;
1124 }
1125 
1126 ////////////////////////////////////////////////////////////////////////////////
1127 /// Enable rendering of children and their list contents.
1128 /// Arguments control how to set self/child rendering.
1129 
1130 void REveElement::EnableListElements(Bool_t rnr_self, Bool_t rnr_children)
1131 {
1132  for (auto &c: fChildren)
1133  c->SetRnrSelfChildren(rnr_self, rnr_children);
1134 }
1135 
1136 ////////////////////////////////////////////////////////////////////////////////
1137 /// Disable rendering of children and their list contents.
1138 /// Arguments control how to set self/child rendering.
1139 ///
1140 /// Same as above function, but default arguments are different. This
1141 /// is convenient for calls via context menu.
1142 
1143 void REveElement::DisableListElements(Bool_t rnr_self, Bool_t rnr_children)
1144 {
1145  for (auto &c: fChildren)
1146  c->SetRnrSelfChildren(rnr_self, rnr_children);
1147 }
1148 
1149 ////////////////////////////////////////////////////////////////////////////////
1150 /// Protected member function called from REveElement::Annihilate().
1151 
1152 void REveElement::AnnihilateRecursively()
1153 {
1154  // projected were already destroyed in REveElement::Anihilate(), now only clear its list
1155  REveProjectable* pable = dynamic_cast<REveProjectable*>(this);
1156  if (pable && pable->HasProjecteds())
1157  pable->ClearProjectedList();
1158 
1159  // same as REveElement::RemoveElementsInternal(), except parents are ignored
1160  RemoveElementsLocal();
1161  for (auto &c : fChildren)
1162  c->AnnihilateRecursively();
1163 
1164  fChildren.clear();
1165 
1166  fDestructing = kAnnihilate;
1167  PreDeleteElement();
1168 
1169  delete this;
1170 }
1171 
1172 ////////////////////////////////////////////////////////////////////////////////
1173 /// Optimized destruction without check of reference-count.
1174 /// Parents are not notified about child destruction.
1175 /// The method should only be used when an element does not have
1176 /// more than one parent -- otherwise an exception is thrown.
1177 
1178 void REveElement::Annihilate()
1179 {
1180  static const REveException eh("REveElement::Annihilate ");
1181 
1182  fDestructing = kAnnihilate;
1183 
1184  // recursive annihilation of projecteds
1185  REveProjectable* pable = dynamic_cast<REveProjectable*>(this);
1186  if (pable && pable->HasProjecteds())
1187  {
1188  pable->AnnihilateProjecteds();
1189  }
1190 
1191  // detach from the parent
1192  if (fMother)
1193  {
1194  fMother->RemoveElement(this);
1195  }
1196 
1197  // XXXX wont the above already start off the destruction cascade ?????
1198 
1199  AnnihilateRecursively();
1200 
1201  // XXXX ????? Annihilate flag ???? Is it different than regular remove ????
1202  // REX::gEve->Redraw3D();
1203 }
1204 
1205 ////////////////////////////////////////////////////////////////////////////////
1206 /// Annihilate elements.
1207 
1208 void REveElement::AnnihilateElements()
1209 {
1210  while (!fChildren.empty())
1211  {
1212  auto c = fChildren.front();
1213  c->Annihilate();
1214  }
1215 }
1216 
1217 ////////////////////////////////////////////////////////////////////////////////
1218 /// Destroy this element. Throws an exception if deny-destroy is in force.
1219 /// This method should be called instead of a destructor.
1220 /// Note that an exception will be thrown if the element has been
1221 /// protected against destruction with IncDenyDestroy().
1222 
1223 void REveElement::Destroy()
1224 {
1225  static const REveException eh("REveElement::Destroy ");
1226 
1227  if (fDenyDestroy > 0)
1228  throw eh + TString::Format("element '%s' (%s*) %p is protected against destruction.",
1229  GetCName(), IsA()->GetName(), this);
1230 
1231  PreDeleteElement();
1232  delete this;
1233  REX::gEve->Redraw3D();
1234 }
1235 
1236 ////////////////////////////////////////////////////////////////////////////////
1237 /// Destroy this element. Prints a warning if deny-destroy is in force.
1238 
1239 void REveElement::DestroyOrWarn()
1240 {
1241  try
1242  {
1243  Destroy();
1244  }
1245  catch (REveException &exc)
1246  {
1247  ::Warning("REveElement::DestroyOrWarn", "Error while destroy element %p : %s", this, exc.what());
1248  }
1249 }
1250 
1251 ////////////////////////////////////////////////////////////////////////////////
1252 /// Destroy all children of this element.
1253 
1254 void REveElement::DestroyElements()
1255 {
1256  while (HasChildren())
1257  {
1258  auto c = fChildren.front();
1259  if (c->fDenyDestroy <= 0)
1260  {
1261  try {
1262  c->Destroy();
1263  }
1264  catch (REveException &exc) {
1265  ::Warning("REveElement::DestroyElements", "element destruction failed: '%s'.", exc.what());
1266  RemoveElement(c);
1267  }
1268  }
1269  else
1270  {
1271  if (gDebug > 0)
1272  ::Info("REveElement::DestroyElements", "element '%s' is protected against destruction, removing locally.", c->GetCName());
1273  RemoveElement(c);
1274  }
1275  }
1276 
1277  REX::gEve->Redraw3D();
1278 }
1279 
1280 ////////////////////////////////////////////////////////////////////////////////
1281 /// Returns state of flag determining if the element will be
1282 /// destroyed when reference count reaches zero.
1283 /// This is true by default.
1284 
1285 Bool_t REveElement::GetDestroyOnZeroRefCnt() const
1286 {
1287  return fDestroyOnZeroRefCnt;
1288 }
1289 
1290 ////////////////////////////////////////////////////////////////////////////////
1291 /// Sets the state of flag determining if the element will be
1292 /// destroyed when reference count reaches zero.
1293 /// This is true by default.
1294 
1295 void REveElement::SetDestroyOnZeroRefCnt(Bool_t d)
1296 {
1297  fDestroyOnZeroRefCnt = d;
1298 }
1299 
1300 ////////////////////////////////////////////////////////////////////////////////
1301 /// Returns the number of times deny-destroy has been requested on
1302 /// the element.
1303 
1304 Int_t REveElement::GetDenyDestroy() const
1305 {
1306  return fDenyDestroy;
1307 }
1308 
1309 ////////////////////////////////////////////////////////////////////////////////
1310 /// Increases the deny-destroy count of the element.
1311 /// Call this if you store an external pointer to the element.
1312 
1313 void REveElement::IncDenyDestroy()
1314 {
1315  ++fDenyDestroy;
1316 }
1317 
1318 ////////////////////////////////////////////////////////////////////////////////
1319 /// Decreases the deny-destroy count of the element.
1320 /// Call this after releasing an external pointer to the element.
1321 
1322 void REveElement::DecDenyDestroy()
1323 {
1324  if (--fDenyDestroy <= 0)
1325  CheckReferenceCount("REveElement::DecDenyDestroy ");
1326 }
1327 
1328 ////////////////////////////////////////////////////////////////////////////////
1329 /// Set pickable state on the element and all its children.
1330 
1331 void REveElement::SetPickableRecursively(Bool_t p)
1332 {
1333  fPickable = p;
1334  for (auto &c: fChildren)
1335  c->SetPickableRecursively(p);
1336 }
1337 
1338 ////////////////////////////////////////////////////////////////////////////////
1339 /// Returns the master element - that is:
1340 /// - master of projectable, if this is a projected;
1341 /// - master of compound, if fCompound is set;
1342 /// - master of mother, if kSCBTakeMotherAsMaster bit is set;
1343 /// If non of the above is true, *this* is returned.
1344 
1345 REveElement *REveElement::GetSelectionMaster()
1346 {
1347  if (fSelectionMaster) return fSelectionMaster;
1348 
1349  REveProjected* proj = dynamic_cast<REveProjected*>(this);
1350  if (proj)
1351  {
1352  return dynamic_cast<REveElement*>(proj->GetProjectable())->GetSelectionMaster();
1353  }
1354  if (fCompound)
1355  {
1356  return fCompound->GetSelectionMaster();
1357  }
1358  if (TestCSCBits(kCSCBTakeMotherAsMaster) && fMother)
1359  {
1360  return fMother->GetSelectionMaster();
1361  }
1362  return this;
1363 }
1364 
1365 ////////////////////////////////////////////////////////////////////////////////
1366 /// Populate set impSelSet with derived / dependant elements.
1367 ///
1368 /// If this is a REveProjectable, the projected replicas are added
1369 /// to the set. Thus it does not have to be reimplemented for each
1370 /// sub-class of REveProjected.
1371 ///
1372 /// Note that this also takes care of projections of REveCompound
1373 /// class, which is also a projectable.
1374 
1375 void REveElement::FillImpliedSelectedSet(Set_t &impSelSet)
1376 {
1377  REveProjectable* p = dynamic_cast<REveProjectable*>(this);
1378  if (p)
1379  p->AddProjectedsToSet(impSelSet);
1380 }
1381 
1382 ////////////////////////////////////////////////////////////////////////////////
1383 /// Call this if it is possible that implied-selection or highlight
1384 /// has changed for this element or for implied-selection this
1385 /// element is member of and you want to maintain consistent
1386 /// selection state.
1387 /// This can happen if you add elements into compounds in response
1388 /// to user-interaction.
1389 
1390 void REveElement::RecheckImpliedSelections()
1391 {
1392  // XXXX MT 2019-01 --- RecheckImpliedSelections
1393  //
1394  // With removal of selection state from this class there might be some
1395  // corner cases requiring checking of implied-selected state in
1396  // selection/highlight objects.
1397  //
1398  // This could be done as part of begin / end changes on the EveManager level.
1399  //
1400  // See also those functions in TEveSelection.
1401 
1402  // if (fSelected || fImpliedSelected)
1403  // REX::gEve->GetSelection()->RecheckImpliedSetForElement(this);
1404 
1405  // if (fHighlighted || fImpliedHighlighted)
1406  // REX::gEve->GetHighlight()->RecheckImpliedSetForElement(this);
1407 }
1408 
1409 ////////////////////////////////////////////////////////////////////////////////
1410 /// Add (bitwise or) given stamps to fChangeBits.
1411 /// Register this element to REX::gEve as stamped.
1412 /// This method is virtual so that sub-classes can add additional
1413 /// actions. The base-class method should still be called (or replicated).
1414 
1415 void REveElement::AddStamp(UChar_t bits)
1416 {
1417  if (fDestructing == kNone && fScene && fScene->IsAcceptingChanges())
1418  {
1419  if (gDebug > 0)
1420  ::Info(Form("%s::AddStamp", GetCName()), "%d + (%d) -> %d", fChangeBits, bits, fChangeBits | bits);
1421 
1422  if (fChangeBits == 0)
1423  {
1424  fScene->SceneElementChanged(this);
1425  }
1426 
1427  fChangeBits |= bits;
1428  }
1429 }
1430 
1431 ////////////////////////////////////////////////////////////////////////////////
1432 /// Write transformation Matrix to render data
1433 ////////////////////////////////////////////////////////////////////////////////
1434 
1435 void REveElement::BuildRenderData()
1436 {
1437  if (fMainTrans.get())
1438  {
1439  fRenderData->SetMatrix(fMainTrans->Array());
1440  }
1441 }
1442 
1443 ////////////////////////////////////////////////////////////////////////////////
1444 /// Convert Bool_t to string - kTRUE or kFALSE.
1445 /// Needed in WriteVizParams().
1446 
1447 const std::string& REveElement::ToString(Bool_t b)
1448 {
1449  static const std::string true_str ("kTRUE");
1450  static const std::string false_str("kFALSE");
1451 
1452  return b ? true_str : false_str;
1453 }
1454 
1455 ////////////////////////////////////////////////////////////////////////////////
1456 /// Write core json. If rnr_offset is negative, render data shall not be
1457 /// written.
1458 /// Returns number of bytes written into binary render data.
1459 
1460 Int_t REveElement::WriteCoreJson(nlohmann::json &j, Int_t rnr_offset)
1461 {
1462  j["_typename"] = IsA()->GetName();
1463  j["fName"] = fName;
1464  j["fTitle"] = fTitle;
1465  j["fElementId"] = GetElementId();
1466  j["fMotherId"] = get_mother_id();
1467  j["fSceneId"] = get_scene_id();
1468  j["fMasterId"] = GetSelectionMaster()->GetElementId();
1469 
1470  j["fRnrSelf"] = GetRnrSelf();
1471  j["fRnrChildren"] = GetRnrChildren();
1472 
1473  j["fMainColor"] = GetMainColor();
1474  j["fMainTransparency"] = GetMainTransparency();
1475  j["fPickable"] = fPickable;
1476 
1477  Int_t ret = 0;
1478 
1479  if (rnr_offset >= 0) {
1480  BuildRenderData();
1481 
1482  if (fRenderData) {
1483  nlohmann::json rd = {};
1484 
1485  rd["rnr_offset"] = rnr_offset;
1486  rd["rnr_func"] = fRenderData->GetRnrFunc();
1487  rd["vert_size"] = fRenderData->SizeV();
1488  rd["norm_size"] = fRenderData->SizeN();
1489  rd["index_size"] = fRenderData->SizeI();
1490  rd["trans_size"] = fRenderData->SizeT();
1491 
1492  j["render_data"] = rd;
1493 
1494  ret = fRenderData->GetBinarySize();
1495  }
1496  }
1497 
1498  return ret;
1499 }