Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TGLScenePad.cxx
Go to the documentation of this file.
1 // @(#)root/gl:$Id$
2 // Author: Matevz Tadel, Jun 2007
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2004, 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 "TGLScenePad.h"
13 
14 #include "TGLViewer.h"
15 #include "TGLLogicalShape.h"
16 #include "TGLPhysicalShape.h"
17 #include "TGLObject.h"
18 #include "TGLStopwatch.h"
19 #include "TBuffer3D.h"
20 #include "TBuffer3DTypes.h"
21 #include "TPolyMarker3D.h"
22 #include "TColor.h"
23 #include "TROOT.h"
24 #include "TH3.h"
25 
26 #include "TGLFaceSet.h"
27 #include "TGLPolyLine.h"
28 #include "TGLPolyMarker.h"
29 #include "TGLCylinder.h"
30 #include "TGLSphere.h"
31 
32 #include "TVirtualPad.h"
33 #include "TAtt3D.h"
34 #include "TClass.h"
35 #include "TList.h"
36 #include "TMath.h"
37 
38 #include "TGLPlot3D.h"
39 
40 
41 /** \class TGLScenePad
42 \ingroup opengl
43 Implements VirtualViewer3D interface and fills the base-class
44 visualization structures from pad contents.
45 */
46 
47 ClassImp(TGLScenePad);
48 
49 ////////////////////////////////////////////////////////////////////////////////
50 
51 TGLScenePad::TGLScenePad(TVirtualPad* pad) :
52  TVirtualViewer3D(),
53  TGLScene(),
54 
55  fPad (pad),
56  fInternalPIDs (kFALSE),
57  fNextInternalPID (1), // 0 reserved
58  fLastPID (0), // 0 reserved
59  fAcceptedPhysicals (0),
60  fComposite (0),
61  fCSLevel (0),
62  fSmartRefresh (kFALSE)
63 {
64  // Constructor.
65 }
66 
67 
68 /******************************************************************************/
69 // Histo import and Sub-pad traversal
70 /******************************************************************************/
71 
72 ////////////////////////////////////////////////////////////////////////////////
73 /// Scale and rotate a histo object to mimic placement in canvas.
74 
75 void TGLScenePad::AddHistoPhysical(TGLLogicalShape* log, const Float_t *histoColor)
76 {
77  Double_t how = ((Double_t) gPad->GetWh()) / gPad->GetWw();
78 
79  Double_t lw = gPad->GetAbsWNDC();
80  Double_t lh = gPad->GetAbsHNDC() * how;
81  Double_t lm = TMath::Min(lw, lh);
82 
83  const TGLBoundingBox& bb = log->BoundingBox();
84 
85  // Timur always packs histos in a square: let's just take x-diff.
86  Double_t size = TMath::Sqrt(3) * (bb.XMax() - bb.XMin());
87  Double_t scale = lm / size;
88  TGLVector3 scaleVec(scale, scale, scale);
89 
90  Double_t tx = gPad->GetAbsXlowNDC() + lw;
91  Double_t ty = gPad->GetAbsYlowNDC() * how + lh;
92  TGLVector3 transVec(0, ty, tx); // For viewer convention (starts looking along -x).
93 
94  // XXXX plots no longer centered at 0. Or they never were?
95  // Impossible to translate and scale them as they should be, it
96  // seems. This requires further investigation, eventually.
97  //
98  // bb.Dump();
99  // printf("lm=%f, size=%f, scale=%f, tx=%f, ty=%f\n",
100  // lm, size, scale, tx, ty);
101  //
102  // TGLVector3 c(bb.Center().Arr());
103  // c.Negate();
104  // c.Dump();
105  // mat.Translate(c);
106 
107  TGLMatrix mat;
108  mat.Scale(scaleVec);
109  mat.Translate(transVec);
110  mat.RotateLF(3, 2, TMath::PiOver2());
111  mat.RotateLF(1, 3, TMath::DegToRad()*gPad->GetTheta());
112  mat.RotateLF(1, 2, TMath::DegToRad()*(gPad->GetPhi() - 90));
113  Float_t rgba[4] = {1.f, 1.f, 1.f, 1.f};
114  if (histoColor) {
115  rgba[0] = histoColor[0];
116  rgba[1] = histoColor[1];
117  rgba[2] = histoColor[2];
118  rgba[3] = histoColor[3];
119  }
120  TGLPhysicalShape* phys = new TGLPhysicalShape(fNextInternalPID++, *log, mat, false, rgba);
121  AdoptPhysical(*phys);
122 
123  // Part of XXXX above.
124  // phys->BoundingBox().Dump();
125 }
126 
127 namespace {
128 
129 ////////////////////////////////////////////////////////////////////////////////
130 ///TTree::Draw can create polymarker + empty TH3 (to draw as a frame around marker).
131 ///Unfortunately, this is not good for GL - this will be two unrelated
132 ///objects in two unrelated coordinate systems.
133 ///So, this function checks list contents, and if it founds empty TH3 and polymarker,
134 ///the must be combined as one object.
135 ///Later we'll reconsider the design.
136 
137 Bool_t HasPolymarkerAndFrame(const TList *lst)
138 {
139  Bool_t gotEmptyTH3 = kFALSE;
140  Bool_t gotMarker = kFALSE;
141 
142  TObjOptLink *lnk = lst ? (TObjOptLink*)lst->FirstLink() : 0;
143  for (; lnk; lnk = (TObjOptLink*)lnk->Next()) {
144  const TObject *obj = lnk->GetObject();
145  if (const TH3 *th3 = dynamic_cast<const TH3*>(obj)) {
146  if(!th3->GetEntries())
147  gotEmptyTH3 = kTRUE;
148  } else if (dynamic_cast<const TPolyMarker3D *>(obj))
149  gotMarker = kTRUE;
150  }
151 
152  return gotMarker && gotEmptyTH3;
153 }
154 
155 }
156 
157 
158 
159 ////////////////////////////////////////////////////////////////////////////////
160 /// Iterate over pad-primitives and import them.
161 
162 void TGLScenePad::SubPadPaint(TVirtualPad* pad)
163 {
164  TVirtualPad *padsav = gPad;
165  TVirtualViewer3D *vv3dsav = pad->GetViewer3D();
166  gPad = pad;
167  pad->SetViewer3D(this);
168 
169  TList *prims = pad->GetListOfPrimitives();
170 
171  if (HasPolymarkerAndFrame(prims)) {
172  ComposePolymarker(prims);
173  } else {
174  TObjOptLink *lnk = (prims) ? (TObjOptLink*)prims->FirstLink() : 0;
175  for (; lnk; lnk = (TObjOptLink*)lnk->Next())
176  ObjectPaint(lnk->GetObject(), lnk->GetOption());
177  }
178 
179  pad->SetViewer3D(vv3dsav);
180  gPad = padsav;
181 }
182 
183 
184 ////////////////////////////////////////////////////////////////////////////////
185 /// Override of virtual TVirtualViewer3D::ObjectPaint().
186 /// Special handling of 2D/3D histograms to activate Timur's
187 /// histo-painters.
188 
189 void TGLScenePad::ObjectPaint(TObject* obj, Option_t* opt)
190 {
191  TGLPlot3D* log = TGLPlot3D::CreatePlot(obj, opt, gPad);
192  if (log)
193  {
194  AdoptLogical(*log);
195  AddHistoPhysical(log);
196  }
197  else if (obj->InheritsFrom(TAtt3D::Class()))
198  {
199  // Handle 3D primitives here.
200  obj->Paint(opt);
201  }
202  else if (obj->InheritsFrom(TVirtualPad::Class()))
203  {
204  SubPadPaint(dynamic_cast<TVirtualPad*>(obj));
205  }
206  else
207  {
208  // Handle 2D primitives here.
209  obj->Paint(opt);
210  }
211 }
212 
213 ////////////////////////////////////////////////////////////////////////////////
214 /// Entry point for requesting update of scene's contents from
215 /// gl-viewer.
216 
217 void TGLScenePad::PadPaintFromViewer(TGLViewer* viewer)
218 {
219  Bool_t sr = fSmartRefresh;
220  fSmartRefresh = viewer->GetSmartRefresh();
221 
222  PadPaint(fPad);
223 
224  fSmartRefresh = sr;
225 }
226 
227 ////////////////////////////////////////////////////////////////////////////////
228 /// Entry point for updating scene contents via VirtualViewer3D
229 /// interface.
230 /// For now this is handled by TGLViewer as it remains
231 /// the 'Viewer3D' of given pad.
232 
233 void TGLScenePad::PadPaint(TVirtualPad* pad)
234 {
235  if (pad != fPad)
236  {
237  Error("TGLScenePad::PadPaint", "Mismatch between pad argument and data-member!");
238  return;
239  }
240 
241  BeginScene();
242  SubPadPaint(fPad);
243  EndScene();
244 }
245 
246 
247 //==============================================================================
248 // VV3D
249 //==============================================================================
250 
251 ////////////////////////////////////////////////////////////////////////////////
252 /// Start building of the scene.
253 /// Old contents is dropped, unless smart-refresh is in active. Then
254 /// the object supporting it are kept in a cache and possibly reused.
255 ///
256 /// TVirtualViewer3D interface overload - see base/src/TVirtualViewer3D.cxx
257 /// for description of viewer architecture.
258 
259 void TGLScenePad::BeginScene()
260 {
261  if (gDebug>2) {
262  Info("TGLScenePad::BeginScene", "entering.");
263  }
264 
265  if ( ! BeginUpdate()) {
266  Error("TGLScenePad::BeginScene", "could not take scene lock.");
267  return;
268  }
269 
270  UInt_t destroyedLogicals = 0;
271  UInt_t destroyedPhysicals = 0;
272 
273  TGLStopwatch stopwatch;
274  if (gDebug > 2) {
275  stopwatch.Start();
276  }
277 
278  // Rebuilds can potentially invalidate all logical and
279  // physical shapes.
280  // Physicals must be removed first.
281  destroyedPhysicals = DestroyPhysicals();
282  if (fSmartRefresh) {
283  destroyedLogicals = BeginSmartRefresh();
284  } else {
285  destroyedLogicals = DestroyLogicals();
286  }
287 
288  // Potentially using external physical IDs
289  fInternalPIDs = kFALSE;
290 
291  // Reset internal physical ID counter
292  fNextInternalPID = 1;
293  fLastPID = 0;
294 
295  // Reset tracing info
296  fAcceptedPhysicals = 0;
297 
298  if (gDebug > 2) {
299  Info("TGLScenePad::BeginScene", "destroyed %d physicals %d logicals in %f msec",
300  destroyedPhysicals, destroyedLogicals, stopwatch.End());
301  DumpMapSizes();
302  }
303 }
304 
305 ////////////////////////////////////////////////////////////////////////////////
306 /// End building of the scene.
307 /// TVirtualViewer3D interface overload - see base/src/TVirtualViewer3D.cxx
308 /// for description of viewer architecture
309 
310 void TGLScenePad::EndScene()
311 {
312  if (fSmartRefresh) {
313  EndSmartRefresh();
314  }
315 
316  EndUpdate();
317 
318  if (gDebug > 2) {
319  Info("TGLScenePad::EndScene", "Accepted %d physicals", fAcceptedPhysicals);
320  DumpMapSizes();
321  }
322 }
323 
324 ////////////////////////////////////////////////////////////////////////////////
325 /// Add an object to the viewer, using internal physical IDs
326 /// TVirtualViewer3D interface overload - see base/src/TVirtualViewer3D.cxx
327 /// for description of viewer architecture
328 
329 Int_t TGLScenePad::AddObject(const TBuffer3D& buffer, Bool_t* addChildren)
330 {
331  // If this is called we are generating internal physical IDs
332  fInternalPIDs = kTRUE;
333  Int_t sections = AddObject(fNextInternalPID, buffer, addChildren);
334  return sections;
335 }
336 
337 ////////////////////////////////////////////////////////////////////////////////
338 /// Add an object to the scene, using an external physical ID
339 /// TVirtualViewer3D interface overload - see base/src/TVirtualViewer3D.cxx
340 /// for description of viewer architecture
341 
342 Int_t TGLScenePad::AddObject(UInt_t physicalID, const TBuffer3D& buffer, Bool_t* addChildren)
343 {
344  // TODO: Break this up and make easier to understand. This is
345  // pretty convoluted due to the large number of cases it has to
346  // deal with:
347  // i) existing physical and/or logical;
348  // ii) external provider may or may not supply bounding box;
349  // iii) local/global reference frame;
350  // iv) deferred filling of some sections of the buffer;
351  // v) internal or external physical IDs;
352  // vi) composite components as special case.
353  //
354  // The buffer filling means the function is re-entrant which adds
355  // to complication.
356 
357  if (physicalID == 0) {
358  Error("TGLScenePad::AddObject", "0 physical ID reserved");
359  return TBuffer3D::kNone;
360  }
361 
362  // Internal and external physical IDs cannot be mixed in a scene build
363  if (fInternalPIDs && physicalID != fNextInternalPID) {
364  Error("TGLScenePad::AddObject", "invalid next physical ID - mix of internal + external IDs?");
365  return TBuffer3D::kNone;
366  }
367 
368  // We always take all children ... interest is viewer dependent.
369  if (addChildren)
370  *addChildren = kTRUE;
371 
372  // Scene should be modify locked
373  if (CurrentLock() != kModifyLock) {
374  Error("TGLScenePad::AddObject", "expected scene to be modify-locked.");
375  return TBuffer3D::kNone;
376  }
377 
378  // Note that 'object' here is really a physical/logical pair described
379  // in buffer + physical ID.
380 
381  // If adding component to a current partial composite do this now
382  if (fComposite) {
383  RootCsg::TBaseMesh *newMesh = RootCsg::ConvertToMesh(buffer);
384  // Solaris CC can't create stl pair with enumerate type
385  fCSTokens.push_back(std::make_pair(static_cast<UInt_t>(TBuffer3D::kCSNoOp), newMesh));
386  return TBuffer3D::kNone;
387  }
388 
389  // TODO: Could be a data member - save possible double lookup?
390  TGLPhysicalShape *physical = FindPhysical(physicalID);
391  TGLLogicalShape *logical = 0;
392 
393  // If we have a valid (non-zero) ID, see if the logical is already cached.
394  // If it is not, try to create a direct renderer object.
395  if (buffer.fID)
396  {
397  logical = FindLogical(buffer.fID);
398  if (!logical)
399  logical = AttemptDirectRenderer(buffer.fID);
400  }
401 
402  // First attempt to add this physical.
403  if (physicalID != fLastPID)
404  {
405  // Existing physical.
406  // MT comment: I don't think this should ever happen.
407  if (physical)
408  {
409  // If we have physical we should have logical cached, too.
410  if (!logical) {
411  Error("TGLScenePad::AddObject", "cached physical with no associated cached logical");
412  }
413 
414  // Since we already have logical no need for further checks.
415  // Done ... prepare for next object.
416  if (fInternalPIDs)
417  ++fNextInternalPID;
418 
419  return TBuffer3D::kNone;
420  }
421 
422  // Need any extra sections in buffer?
423  Bool_t includeRaw = (logical == 0);
424  Int_t extraSections = ValidateObjectBuffer(buffer, includeRaw);
425  if (extraSections != TBuffer3D::kNone)
426  return extraSections;
427 
428  fLastPID = physicalID;
429  }
430 
431  if (fLastPID != physicalID) {
432  Error("TGLScenePad::AddObject", "internal physical ID tracking error?");
433  }
434 
435  // Being here means we need to add a physical, maybe logical as well.
436  if (physical) {
437  Error("TGLScenePad::AddObject", "expecting to require physical");
438  return TBuffer3D::kNone;
439  }
440 
441  // Create logical if required.
442  if (!logical)
443  {
444  logical = CreateNewLogical(buffer);
445  if (!logical) {
446  Error("TGLScenePad::AddObject", "failed to create logical");
447  return TBuffer3D::kNone;
448  }
449  // Add logical to scene
450  AdoptLogical(*logical);
451  }
452 
453  // Create the physical, bind it to the logical and add it to the scene.
454  physical = CreateNewPhysical(physicalID, buffer, *logical);
455 
456  if (physical)
457  {
458  AdoptPhysical(*physical);
459  buffer.fPhysicalID = physicalID;
460  ++fAcceptedPhysicals;
461  if (gDebug>3 && fAcceptedPhysicals%1000 == 0) {
462  Info("TGLScenePad::AddObject", "added %d physicals", fAcceptedPhysicals);
463  }
464  }
465  else
466  {
467  Error("TGLScenePad::AddObject", "failed to create physical");
468  }
469 
470  // Done ... prepare for next object.
471  if (fInternalPIDs)
472  fNextInternalPID++;
473 
474  return TBuffer3D::kNone;
475 }
476 
477 ////////////////////////////////////////////////////////////////////////////////
478 /// Open new composite container.
479 /// TVirtualViewer3D interface overload - see base/src/TVirtualViewer3D.cxx
480 /// for description of viewer architecture.
481 
482 Bool_t TGLScenePad::OpenComposite(const TBuffer3D& buffer, Bool_t* addChildren)
483 {
484  if (fComposite) {
485  Error("TGLScenePad::OpenComposite", "composite already open");
486  return kFALSE;
487  }
488  UInt_t extraSections = AddObject(buffer, addChildren);
489  if (extraSections != TBuffer3D::kNone) {
490  Error("TGLScenePad::OpenComposite", "expected top level composite to not require extra buffer sections");
491  }
492 
493  // If composite was created it is of interest - we want the rest of the
494  // child components
495  if (fComposite) {
496  return kTRUE;
497  } else {
498  return kFALSE;
499  }
500 }
501 
502 ////////////////////////////////////////////////////////////////////////////////
503 /// Close composite container
504 /// TVirtualViewer3D interface overload - see base/src/TVirtualViewer3D.cxx
505 /// for description of viewer architecture
506 
507 void TGLScenePad::CloseComposite()
508 {
509  // If we have a partially complete composite build it now
510  if (fComposite) {
511  // TODO: Why is this member and here - only used in BuildComposite()
512  fCSLevel = 0;
513 
514  RootCsg::TBaseMesh *resultMesh = BuildComposite();
515  fComposite->SetFromMesh(resultMesh);
516  delete resultMesh;
517  for (UInt_t i = 0; i < fCSTokens.size(); ++i) delete fCSTokens[i].second;
518  fCSTokens.clear();
519  fComposite = 0;
520  }
521 }
522 
523 ////////////////////////////////////////////////////////////////////////////////
524 /// Add composite operation used to combine objects added via AddObject
525 /// TVirtualViewer3D interface overload - see base/src/TVirtualViewer3D.cxx
526 /// for description of viewer architecture
527 
528 void TGLScenePad::AddCompositeOp(UInt_t operation)
529 {
530  fCSTokens.push_back(std::make_pair(operation, (RootCsg::TBaseMesh *)0));
531 }
532 
533 
534 // Protected methods
535 
536 ////////////////////////////////////////////////////////////////////////////////
537 /// Validate if the passed 'buffer' contains all sections we require to add object.
538 /// Returns Int_t combination of TBuffer::ESection flags still required - or
539 /// TBuffer3D::kNone if buffer is valid.
540 /// If 'includeRaw' is kTRUE check for kRaw/kRawSizes - skip otherwise.
541 /// See base/src/TVirtualViewer3D.cxx for description of viewer architecture
542 
543 Int_t TGLScenePad::ValidateObjectBuffer(const TBuffer3D& buffer, Bool_t includeRaw) const
544 {
545  // kCore: Should always be filled
546  if (!buffer.SectionsValid(TBuffer3D::kCore)) {
547  Error("TGLScenePad::ValidateObjectBuffer", "kCore section of buffer should be filled always");
548  return TBuffer3D::kNone;
549  }
550 
551  // Need to check raw (kRaw/kRawSizes)?
552  if (!includeRaw) {
553  return TBuffer3D::kNone;
554  }
555 
556  // kRawSizes / kRaw: These are on demand based on shape type
557  Bool_t needRaw = kFALSE;
558 
559  // We need raw tesselation in these cases:
560  //
561  // 1. Shape type is NOT kSphere / kTube / kTubeSeg / kCutTube / kComposite
562  if (buffer.Type() != TBuffer3DTypes::kSphere &&
563  buffer.Type() != TBuffer3DTypes::kTube &&
564  buffer.Type() != TBuffer3DTypes::kTubeSeg &&
565  buffer.Type() != TBuffer3DTypes::kCutTube &&
566  buffer.Type() != TBuffer3DTypes::kComposite)
567  {
568  needRaw = kTRUE;
569  }
570  // 2. Sphere type is kSPHE, but the sphere is hollow and/or cut - we
571  // do not support native drawing of these currently
572  else if (buffer.Type() == TBuffer3DTypes::kSphere)
573  {
574  const TBuffer3DSphere * sphereBuffer = dynamic_cast<const TBuffer3DSphere *>(&buffer);
575  if (sphereBuffer) {
576  if (!sphereBuffer->IsSolidUncut()) {
577  needRaw = kTRUE;
578  }
579  } else {
580  Error("TGLScenePad::ValidateObjectBuffer", "failed to cast buffer of type 'kSphere' to TBuffer3DSphere");
581  return TBuffer3D::kNone;
582  }
583  }
584  // 3. kBoundingBox is not filled - we generate a bounding box from
585  else if (!buffer.SectionsValid(TBuffer3D::kBoundingBox))
586  {
587  needRaw = kTRUE;
588  }
589  // 4. kShapeSpecific is not filled - except in case of top level composite
590  else if (!buffer.SectionsValid(TBuffer3D::kShapeSpecific) &&
591  buffer.Type() != TBuffer3DTypes::kComposite)
592  {
593  needRaw = kTRUE;
594  }
595  // 5. We are a component (not the top level) of a composite shape
596  else if (fComposite)
597  {
598  needRaw = kTRUE;
599  }
600 
601  if (needRaw && !buffer.SectionsValid(TBuffer3D::kRawSizes|TBuffer3D::kRaw)) {
602  return TBuffer3D::kRawSizes|TBuffer3D::kRaw;
603  } else {
604  return TBuffer3D::kNone;
605  }
606 }
607 
608 ////////////////////////////////////////////////////////////////////////////////
609 /// Create and return a new TGLLogicalShape from the supplied buffer
610 
611 TGLLogicalShape* TGLScenePad::CreateNewLogical(const TBuffer3D& buffer) const
612 {
613  TGLLogicalShape * newLogical = 0;
614 
615  if (buffer.fColor == 1) // black -> light-brown; std behaviour for geom
616  const_cast<TBuffer3D&>(buffer).fColor = 42;
617 
618  switch (buffer.Type())
619  {
620  case TBuffer3DTypes::kLine:
621  newLogical = new TGLPolyLine(buffer);
622  break;
623  case TBuffer3DTypes::kMarker:
624  newLogical = new TGLPolyMarker(buffer);
625  break;
626  case TBuffer3DTypes::kSphere:
627  {
628  const TBuffer3DSphere * sphereBuffer = dynamic_cast<const TBuffer3DSphere *>(&buffer);
629  if (sphereBuffer)
630  {
631  // We can only draw solid uncut spheres natively at present.
632  // If somebody already passed the raw buffer, they probably want us to use it.
633  if (sphereBuffer->IsSolidUncut() && !buffer.SectionsValid(TBuffer3D::kRawSizes|TBuffer3D::kRaw))
634  {
635  newLogical = new TGLSphere(*sphereBuffer);
636  } else {
637  newLogical = new TGLFaceSet(buffer);
638  }
639  } else {
640  Error("TGLScenePad::CreateNewLogical", "failed to cast buffer of type 'kSphere' to TBuffer3DSphere");
641  }
642  break;
643  }
644  case TBuffer3DTypes::kTube:
645  case TBuffer3DTypes::kTubeSeg:
646  case TBuffer3DTypes::kCutTube:
647  {
648  const TBuffer3DTube * tubeBuffer = dynamic_cast<const TBuffer3DTube *>(&buffer);
649  if (tubeBuffer)
650  {
651  // If somebody already passed the raw buffer, they probably want us to use it.
652  if (!buffer.SectionsValid(TBuffer3D::kRawSizes|TBuffer3D::kRaw)) {
653  newLogical = new TGLCylinder(*tubeBuffer);
654  } else {
655  newLogical = new TGLFaceSet(buffer);
656  }
657  } else {
658  Error("TGLScenePad::CreateNewLogical", "failed to cast buffer of type 'kTube/kTubeSeg/kCutTube' to TBuffer3DTube");
659  }
660  break;
661  }
662  case TBuffer3DTypes::kComposite:
663  {
664  // Create empty faceset and record partial complete composite object
665  // Will be populated with mesh in CloseComposite()
666  if (fComposite)
667  {
668  Error("TGLScenePad::CreateNewLogical", "composite already open");
669  }
670  fComposite = new TGLFaceSet(buffer);
671  newLogical = fComposite;
672  break;
673  }
674  default:
675  newLogical = new TGLFaceSet(buffer);
676  break;
677  }
678 
679  return newLogical;
680 }
681 
682 ////////////////////////////////////////////////////////////////////////////////
683 /// Create and return a new TGLPhysicalShape with id 'ID', using
684 /// 'buffer' placement information (translation etc), and bound to
685 /// suppled 'logical'
686 
687 TGLPhysicalShape*
688 TGLScenePad::CreateNewPhysical(UInt_t ID, const TBuffer3D& buffer,
689  const TGLLogicalShape& logical) const
690 {
691  // Extract indexed color from buffer
692  // TODO: Still required? Better use proper color triplet in buffer?
693  Int_t colorIndex = buffer.fColor;
694  if (colorIndex < 0) colorIndex = 42;
695  Float_t rgba[4];
696  TGLScene::RGBAFromColorIdx(rgba, colorIndex, buffer.fTransparency);
697  return new TGLPhysicalShape(ID, logical, buffer.fLocalMaster,
698  buffer.fReflection, rgba);
699 }
700 
701 
702 ////////////////////////////////////////////////////////////////////////////////
703 
704 void TGLScenePad::ComposePolymarker(const TList *lst)
705 {
706  TPolyMarker3D *pm = 0;
707  TH3 *th3 = 0;
708  TObjOptLink *lnk = (TObjOptLink*)lst->FirstLink();
709  for (; lnk; lnk = (TObjOptLink*)lnk->Next()) {
710  TObject *obj = lnk->GetObject();
711  if (TPolyMarker3D *dPm = dynamic_cast<TPolyMarker3D*>(obj)) {
712  if(!pm)
713  pm = dPm;
714  } else if (TH3 *dTH3 = dynamic_cast<TH3*>(obj)) {
715  if(!th3 && !dTH3->GetEntries())
716  th3 = dTH3;
717  } else
718  ObjectPaint(obj, lnk->GetOption());
719 
720  if (pm && th3) {
721  //Create a new TH3 plot, containing polymarker.
722  TGLPlot3D* log = TGLPlot3D::CreatePlot(th3, pm);
723  AdoptLogical(*log);
724  //Try to extract polymarker's color and
725  //create a physical shape with correct color.
726  const Color_t cInd = pm->GetMarkerColor();
727  if (TColor *c = gROOT->GetColor(cInd)) {
728  Float_t rgba[4] = {0.f, 0.f, 0.f, 1.};
729  c->GetRGB(rgba[0], rgba[1], rgba[2]);
730  AddHistoPhysical(log, rgba);
731  } else
732  AddHistoPhysical(log);
733 
734  //Composition was added into gl-viewer.
735  pm = 0;
736  th3 = 0;
737  }
738  }
739 }
740 
741 ////////////////////////////////////////////////////////////////////////////////
742 /// Build and return composite shape mesh
743 
744 RootCsg::TBaseMesh* TGLScenePad::BuildComposite()
745 {
746  const CSPart_t &currToken = fCSTokens[fCSLevel];
747  UInt_t opCode = currToken.first;
748 
749  if (opCode != TBuffer3D::kCSNoOp) {
750  ++fCSLevel;
751  RootCsg::TBaseMesh *left = BuildComposite();
752  RootCsg::TBaseMesh *right = BuildComposite();
753  //RootCsg::TBaseMesh *result = 0;
754  switch (opCode) {
755  case TBuffer3D::kCSUnion:
756  return RootCsg::BuildUnion(left, right);
757  case TBuffer3D::kCSIntersection:
758  return RootCsg::BuildIntersection(left, right);
759  case TBuffer3D::kCSDifference:
760  return RootCsg::BuildDifference(left, right);
761  default:
762  Error("BuildComposite", "Wrong operation code %d\n", opCode);
763  return 0;
764  }
765  } else return fCSTokens[fCSLevel++].second;
766 }
767 
768 ////////////////////////////////////////////////////////////////////////////////
769 /// Try to construct an appropriate logical-shape sub-class based
770 /// on id'class, following convention that SomeClassGL is a suitable
771 /// renderer for class SomeClass.
772 
773 TGLLogicalShape* TGLScenePad::AttemptDirectRenderer(TObject* id)
774 {
775  TClass* cls = TGLObject::GetGLRenderer(id->IsA());
776  if (cls == 0)
777  return 0;
778 
779  TGLObject* rnr = reinterpret_cast<TGLObject*>(cls->New());
780  if (rnr) {
781  Bool_t status;
782  try
783  {
784  status = rnr->SetModel(id);
785  }
786  catch (std::exception&)
787  {
788  status = kFALSE;
789  }
790  if (!status)
791  {
792  Warning("TGLScenePad::AttemptDirectRenderer", "failed initializing direct rendering.");
793  delete rnr;
794  return 0;
795  }
796  rnr->SetBBox();
797  AdoptLogical(*rnr);
798  }
799  return rnr;
800 }