Logo ROOT   6.12/07
Reference Guide
TVolumeView.cxx
Go to the documentation of this file.
1 // @(#)root/table:$Id$
2 // Author: Valery Fine(fine@bnl.gov) 25/12/98
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, 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 
13 #include "Riostream.h"
14 #include <assert.h>
15 #include <stdlib.h>
16 
17 #include "TCanvas.h"
18 #include "TPad.h"
19 #include "TCernLib.h"
20 #include "TBrowser.h"
21 #include "TVolumeView.h"
22 #include "TVolumeViewIter.h"
23 #include "TVolumePosition.h"
24 #include "TROOT.h"
25 #include "TView.h"
26 #include "TTablePadView3D.h"
27 #include "TGeometry.h"
28 #include "TVirtualPad.h"
29 #include "TObjArray.h"
30 #include "TVirtualViewer3D.h"
31 #include "TBuffer3D.h"
32 
33 //////////////////////////////////////////////////////////////////////////
34 // //
35 // TVolumeView //
36 // //
37 // TVolumeView class is a special kind of TDataSet with one extra //
38 // pointer to wrap any TObject onto TDataSet object //
39 // //
40 // BE CAREFUL !!! //
41 // One has to use it carefully no control over that extra object //
42 // is performed. This means: the object m_Obj data-member points to can//
43 // be destroyed with no this kbject notifying. //
44 // There is no tool /protection to check whether m_Obj is till alive. //
45 // It is one's code responsilitiy //
46 // //
47 //////////////////////////////////////////////////////////////////////////
48 
50 
51 ////////////////////////////////////////////////////////////////////////////////
52 
54  : TObjectSet(viewNode->GetName(),(TObject *)nodePosition),fListOfShapes(0)
55  // ,fListOfAttributes(0)
56 {
57  //
58  // This ctor creates a TVolumeView structure from the "marked" nodes
59  // of the "viewNode" input structure
60  // It re-calculates all positions according of the new topology
61  // All new TVolume became UNMARKED though
62  //
63  if (!gGeometry) new TGeometry;
64  if (viewNode) {
65  SetTitle(viewNode->GetTitle());
66  EDataSetPass mode = kContinue;
67  TVolumeViewIter next(viewNode,0);
68  TVolumeView *nextView = 0;
69  while ( (nextView = (TVolumeView *)next(mode)) ){
70  mode = kContinue;
71  if (nextView->IsMarked()) {
72  TVolumePosition *position =next[0];
73  if (!position->GetNode()) {
74  Error("TVolumeView ctor","%s %s ",GetName(),nextView->GetName());
75  }
76  Add(new TVolumeView(nextView,position));
77  mode = kPrune;
78  }
79  }
80  }
81 }
82 
83 ////////////////////////////////////////////////////////////////////////////////
84 
86  : TObjectSet(viewNode->GetName(),(TObject *)0),fListOfShapes(0)
87  // ,fListOfAttributes(0)
88 {
89  //
90  // This ctor creates a TVolumeView structure containing:
91  //
92  // - viewNode on the top
93  // - skip ALL node from the original viewNode untill topNode found
94  // - include all "marked" node below "topNode" if any
95  // topNode is always included
96  //
97  // It re-calculates all positions according of the new topology
98  //
99  if (!gGeometry) new TGeometry;
100  if (viewNode && topNode) {
101  SetTitle(viewNode->GetTitle());
102  // define the depth of the "top" Node
103  EDataSetPass mode = kContinue;
104  TVolumeViewIter next(viewNode,0);
105  TVolumeView *nextView = 0;
106  while ( (nextView = (TVolumeView *)next(mode)) ){
107  mode = kContinue;
108  // Skip till "top Node" found
109  if (topNode != nextView) continue;
110  TVolumePosition *position = next[0];
111  if (!position->GetNode()) {
112  Error("TVolumeView ctor","%s %s ",GetName(),nextView->GetName());
113  }
114  Add(new TVolumeView(nextView,position));
115  break;
116  }
117  }
118 }
119 
120 ////////////////////////////////////////////////////////////////////////////////
121 
122 TVolumeView::TVolumeView(TVolumeView *viewNode,const Char_t *nodeName1,const Char_t *nodeName2)
123  : TObjectSet(viewNode->GetName(),(TObject *)0),fListOfShapes(0)
124  // ,fListOfAttributes(0)
125 {
126  //
127  // This ctor creates a TVolumeView structure containing:
128  //
129  // - viewNode on the top
130  // - skip ALL node from the original viewNode untill topNodeName found
131  // - include all "marked" node below "topNodename" if any
132  // topNodeName is always included
133  //
134  // It re-calculates all positions according of the new topology
135  //
136  const Char_t *foundName[2] = {nodeName1, nodeName2};
137  Bool_t found = kFALSE;
138  if (!gGeometry) new TGeometry;
139  if (viewNode && nodeName1 && nodeName1[0]) {
140  SetTitle(viewNode->GetTitle());
141  // define the depth of the "top" Node
142  EDataSetPass mode = kContinue;
143  TVolumeViewIter next(viewNode,0);
144  TVolumeView *nextView = 0;
145  while ( (nextView = (TVolumeView *)next(mode)) ){
146  mode = kContinue;
147  // Skip till "top Node" found
148  Int_t i = 0;
149  found = kFALSE;
150  for (i=0;i<2;i++) {
151  if (foundName[i]) {
152  if (strcmp(nextView->GetName(),foundName[i])) continue;
153  foundName[i] = 0;
154  found = kTRUE;
155  break;
156  }
157  }
158  if (!found) continue;
159  TVolumePosition *position = next[0];
160  if (!position->GetNode()) {
161  Error("TVolumeView ctor","%s %s ",GetName(),nextView->GetName());
162  }
163  Add(new TVolumeView(nextView,position));
164  mode = kPrune;
165  }
166  }
167 }
168 
169 ////////////////////////////////////////////////////////////////////////////////
170 
171 TVolumeView::TVolumeView(TVolumeView *viewNode,const TVolumeView *node1,const TVolumeView *node2)
172  : TObjectSet(viewNode->GetName(),(TObject *)0),fListOfShapes(0)
173  // ,fListOfAttributes(0)
174 {
175  //
176  // This ctor creates a TVolumeView structure containing:
177  //
178  // - viewNode on the top
179  // - skip ALL node from the original viewNode untill topNodeName found
180  // - include all "marked" node below "topNodename" if any
181  // topNodeName is always included
182  //
183  // It re-calculates all positions according of the new topology
184  //
185  const TVolumeView *foundView[2] = {node1, node2};
186  const Int_t nViews = sizeof(foundView)/sizeof(const TVolumeView *);
187  Bool_t found = kFALSE;
188  if (!gGeometry) new TGeometry;
189  if (viewNode) {
190  SetTitle(viewNode->GetTitle());
191  // define the depth of the "top" Node
192  EDataSetPass mode = kContinue;
193  TVolumeViewIter next(viewNode,0);
194  TVolumeView *nextView = 0;
195  while ( (nextView = (TVolumeView *)next(mode)) ){
196  mode = kContinue;
197  // Skip till "top Node" found
198  Int_t i = 0;
199  found = kFALSE;
200  for (i=0;i<nViews;i++) {
201  if (foundView[i]) {
202  if (nextView != foundView[i]) continue;
203  foundView[i] = 0;
204  found = kTRUE;
205  break;
206  }
207  }
208  if (!found) continue;
209  TVolumePosition *position = next[0];
210  if (!position->GetNode()) {
211  Error("TVolumeView ctor","%s %s ",GetName(),nextView->GetName());
212  }
213  Add(new TVolumeView(nextView,position));
214  mode = kPrune;
215  }
216  }
217 }
218 
219 ////////////////////////////////////////////////////////////////////////////////
220 ///
221 /// Creates TVolumeView (view) with a topology similar with TVolume *pattern
222 ///
223 /// Parameters:
224 /// -----------
225 /// pattern - the pattern dataset
226 /// iopt = kStruct - clone only my structural links
227 /// kAll - clone all links
228 /// kRefs - clone only refs
229 /// kMarked - clone marked (not implemented yet) only
230 ///
231 /// All new-created sets become the structural ones anyway.
232 ///
233 /// std::cout << "ctor for " << GetName() << " - " << GetTitle() << std::endl;
234 
236  const TVolumePosition *nodePosition,EDataSetPass iopt, TVolumeView *rootVolume)
237  : TObjectSet(pattern.GetName(),(TObject *)nodePosition),fListOfShapes(0)
238 {
239  if (!gGeometry) new TGeometry;
240  if (!nodePosition) {
241  // Create the trivial position if any
242  nodePosition = new TVolumePosition(&pattern);
243  SetObject((TObject*)nodePosition);
244  }
245  if (!rootVolume) {
246  rootVolume = this;
247  nodePosition = 0;
248  }
249  SetTitle(pattern.GetTitle());
250  if ( pattern.IsMarked() ) Mark();
251  TVolumePosition *position = 0;
252  const TList *list = pattern.GetListOfPositions();
253  if (!list || maxDepLevel == 1 || maxDepLevel < 0) return;
254 
255  TIter nextPosition(list);
256  Bool_t optSel = (iopt == kStruct);
257 // Bool_t optAll = (iopt == kAll);
258  Bool_t optMarked = (iopt == kMarked);
259 
260  const TRotMatrix *thisMatrix = 0;
261  Double_t thisTranslation[3] = {0,0,0};
262  if (nodePosition ) {
263  thisMatrix = nodePosition->GetMatrix();
264  for (int i =0; i< 3; i++) thisTranslation[i]= nodePosition->GetX(i);
265  }
266  while ( (position = (TVolumePosition *)nextPosition()) ) {
267  // define the the related TVolume
268  TVolume *node = position->GetNode();
269  Double_t *positionMatrix = ((TRotMatrix *)position->GetMatrix())->GetMatrix();
270  if (node) {
271  UInt_t positionId = position->GetId();
272  Double_t newTranslation[3] = {position->GetX(),position->GetY(),position->GetZ()};
273  Double_t newMatrix[9];
274  TRotMatrix currentMatrix;
275 
276  if (nodePosition) {
277  if (positionMatrix) {
278  TGeometry::UpdateTempMatrix(thisTranslation,thisMatrix?((TRotMatrix *)thisMatrix)->GetMatrix():0
279  ,position->GetX(),position->GetY(),position->GetZ(),positionMatrix
280  ,newTranslation,newMatrix);
281  currentMatrix.SetMatrix(newMatrix);
282  } else {
283  TCL::vadd(thisTranslation, newTranslation,newTranslation,3);
284  currentMatrix.SetMatrix(((TRotMatrix *)thisMatrix)->GetMatrix());
285  }
286  } else {
287  if (positionMatrix)
288  currentMatrix.SetMatrix(positionMatrix);
289  else {
290  TCL::ucopy(thisTranslation,newTranslation,3);
291  currentMatrix.SetMatrix(TVolume::GetIdentity()->GetMatrix());
292  }
293  }
294  TVolumePosition nextPos(node,newTranslation[0],newTranslation[1],
295  newTranslation[2], &currentMatrix);
296  nextPos.SetId(positionId);
297  if (optMarked && !node->IsMarked()) {
298  TVolumeView fakeView(*node,maxDepLevel,&nextPos,iopt,rootVolume);
299  fakeView.DoOwner(kFALSE);
300  continue;
301  }
302 
303  if (optSel) {
304  TDataSet *parent = node->GetParent();
305  if ( parent && (parent != (TDataSet *)&pattern) ) continue;
306  }
307  TRotMatrix *newRotation = new TRotMatrix();
308  newRotation->SetMatrix(currentMatrix.GetMatrix());
309  TVolumePosition *nP = new TVolumePosition(node,newTranslation[0],newTranslation[1],
310  newTranslation[2], newRotation);
311  nP->SetId(positionId);
312  rootVolume->Add(new TVolumeView(*node,maxDepLevel?maxDepLevel-1:0,nP,iopt));
313  } else
314  Error("TVolumeView ctor","Position with NO node attached has been supplied");
315 
316  }
317 }
318 ////////////////////////////////////////////////////////////////////////////////
319 ///to be documented
320 
322  TObjectSet(viewNode.GetName(),(TObject *)viewNode.GetPosition())
323  ,TAtt3D()
324  ,fListOfShapes(viewNode.GetListOfShapes())
325 {
326  if (viewNode.IsOwner()) {
327  viewNode.DoOwner(kFALSE); DoOwner();
328  }
329 }
330 
331 ////////////////////////////////////////////////////////////////////////////////
332 
333 TVolumeView::TVolumeView(Double_t *translate, Double_t *rotate, UInt_t positionId, TVolume *topNode,
334  const Char_t *thisNodePath, const Char_t *matrixName, Int_t matrixType)
335  // : fListOfAttributes(0)
336 {
337  // Special ctor to back TVolumeView::SavePrimitive() method
338  if (!gGeometry) new TGeometry;
339  fListOfShapes = 0;
340  TVolume *thisNode = 0;
341  Double_t thisX = translate[0];
342  Double_t thisY = translate[1];
343  Double_t thisZ = translate[2];
344 
345  // Find TVolume by path;
346  if (topNode) {
347  thisNode = (TVolume *)topNode->Find(thisNodePath);
348  if (!thisNode->InheritsFrom(TVolume::Class())) {
349  Error("TVolumeView","wrong node <%s> on path: \"%s\"",thisNode->GetName(),thisNodePath);
350  thisNode = 0;
351  }
352  }
353 
354  TRotMatrix *thisRotMatrix = 0;
355  if (matrixName && strlen(matrixName)) thisRotMatrix = gGeometry->GetRotMatrix(matrixName);
356  TVolumePosition *thisPosition = 0;
357  if (thisRotMatrix)
358  thisPosition = new TVolumePosition(thisNode,thisX, thisY, thisZ, matrixName);
359  else if (matrixType==2)
360  thisPosition = new TVolumePosition(thisNode,thisX, thisY, thisZ);
361  else if (rotate) {
362  const Char_t *title = "rotation";
363  thisRotMatrix = new TRotMatrix((char *)matrixName,(char *)title,rotate);
364  thisPosition = new TVolumePosition(thisNode,thisX, thisY, thisZ, thisRotMatrix);
365  } else
366  Error("TVolumeView"," No rotation matrix is defined");
367  if (thisPosition) thisPosition->SetId(positionId);
368  SetObject(thisPosition);
369  if (thisNode) {
370  SetName(thisNode->GetName());
371  SetTitle(thisNode->GetTitle());
372  }
373 }
374 
375 ////////////////////////////////////////////////////////////////////////////////
376 ///to be documented
377 
379  : TObjectSet(thisNode?thisNode->GetName():"",(TObject *)nodePosition),fListOfShapes(0)
380 {
381  if (!gGeometry) new TGeometry;
383  if (thisNode)
384  SetTitle(thisNode->GetTitle());
385 }
386 
387 ////////////////////////////////////////////////////////////////////////////////
388 /// default dtor (empty for this class)
389 
391 {
392 }
393 
394 ////////////////////////////////////////////////////////////////////////////////
395 /// Add the TVolume in the Tnode data-structure refered
396 /// by this TVolumeView object
397 /// Return TVolume * the input TVolume * was attached to
398 
400 {
401  TVolume *closedNode = 0;
402  TVolumePosition *pos ;
403  if ( node && (pos = GetPosition() ) && (closedNode = pos->GetNode()) )
404  closedNode->Add(node);
405  return closedNode;
406 }
407 
408 ////////////////////////////////////////////////////////////////////////////////
409 ///to be documented
410 
411 void TVolumeView::Add(TShape *shape, Bool_t IsMaster)
412 {
413  if (!shape) return;
414  if (!fListOfShapes) fListOfShapes = new TList;
415  if (IsMaster)
416  fListOfShapes->AddFirst(shape);
417  else
418  fListOfShapes->Add(shape);
419 }
420 
421 ////////////////////////////////////////////////////////////////////////////////
422 ///to be documented
423 
425 {
427 // TVolumePosition *pos = GetPosition();
428 // if (pos) pos->Browse(b);
429 // b->Add(pos);
430 }
431 
432 ////////////////////////////////////////////////////////////////////////////////
433 /// Compute distance from point px,py to a TVolumeView.
434 ///
435 /// Compute the closest distance of approach from point px,py to the position of
436 /// this node.
437 /// The distance is computed in pixels units.
438 ///
439 /// It is restricted by 2 levels of TVolumes.
440 
442 {
443  const Int_t big = 9999;
444  const Int_t inaxis = 7;
445  const Int_t maxdist = 5;
446 
447  Int_t dist = big;
448 
449  Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
450  Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
451  Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
452  Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
453 
454 //*-*- return if point is not in the user area
455  if (px < puxmin - inaxis) return big;
456  if (py > puymin + inaxis) return big;
457  if (px > puxmax + inaxis) return big;
458  if (py < puymax - inaxis) return big;
459 
460  TView *view =gPad->GetView();
461  if (!view) return big;
462 
463  TVolumePosition *position = GetPosition();
464  TVolume *thisNode = 0;
465  TShape *thisShape = 0;
466  if (position) {
467  thisNode = position->GetNode();
468  position->UpdatePosition();
469  if (thisNode) {
470  thisShape = thisNode->GetShape();
471  if (!(thisNode->GetVisibility() & TVolume::kThisUnvisible) &&
472  thisShape && thisShape->GetVisibility()) {
473  dist = thisShape->DistancetoPrimitive(px,py);
474  if (dist < maxdist) {
475  gPad->SetSelected(this);
476  return 0;
477  }
478  }
479  }
480  }
481 
482 // if ( TestBit(kSonsInvisible) ) return dist;
483 
484 //*-*- Loop on all sons
485  TSeqCollection *fNodes = GetCollection();
486  Int_t nsons = fNodes?fNodes->GetSize():0;
487  Int_t dnode = dist;
488  if (nsons) {
489  gGeometry->PushLevel();
490  TVolume *node;
491  TIter next(fNodes);
492  while ((node = (TVolume *)next())) {
493  dnode = node->DistancetoPrimitive(px,py);
494  if (dnode <= 0) break;
495  if (dnode < dist) dist = dnode;
496  if (gGeometry->GeomLevel() > 2) break;
497  }
498  gGeometry->PopLevel();
499  }
500 
501  if (gGeometry->GeomLevel()==0 && dnode > maxdist) {
502  gPad->SetSelected(view);
503  return 0;
504  } else
505  return dnode;
506 }
507 
508 ////////////////////////////////////////////////////////////////////////////////
509 /// Draw Referenced node with current parameters.
510 
512 {
513  TString opt = option;
514  opt.ToLower();
515 //*-*- Clear pad if option "same" not given
516  if (!gPad) {
517  gROOT->MakeDefCanvas();
518  }
519  if (!opt.Contains("same")) gPad->Clear();
520 
521 //*-*- Draw Referenced node
524 
525  // Check geometry level
526 
527  Int_t iopt = atoi(option);
528  TDataSet *parent = 0;
529  char buffer[12];
530  if (iopt < 0) {
531  snprintf(buffer,12,"%d",-iopt);
532  option = buffer;
533  // select parent to draw
534  parent = this;
535  do parent = parent->GetParent();
536  while (parent && ++iopt);
537  }
538  if (parent) parent->AppendPad(option);
539  else AppendPad(option);
540 
541 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,03,05)
542  // the new (4.03/05) way to active 3D viewer
543  // Create a 3-D view
544  TView *view = gPad->GetView();
545  if (!view) {
546  view = TView::CreateView(1,0,0);
547  // Set the view to perform a first autorange (frame) draw.
548  // TViewer3DPad will revert view to normal painting after this
549  view->SetAutoRange(kTRUE);
550  }
551 
552  // Create a 3D viewer to draw us
553  gPad->GetViewer3D(option);
554 #else
555  Paint(option);
556 #endif
557 }
558 
559 ////////////////////////////////////////////////////////////////////////////////
560 ///to be documented
561 
563 {
564  TVolumePosition *pos = GetPosition();
565  if (pos)
566  return pos->GetNode();
567  return 0;
568 }
569 
570 ////////////////////////////////////////////////////////////////////////////////
571 ///
572 /// Calculate the position of the vertrex of the outlined cube in repect
573 /// of the given TVolumeView object
574 ///
575 
576 Int_t TVolumeView::GetGlobalRange(const TVolumeView *rootNode,Float_t *globalMin,Float_t *globalMax)
577 {
578  if (rootNode) {
579  SetTitle(rootNode->GetTitle());
580  EDataSetPass mode = kContinue;
581  TVolumeViewIter next((TVolumeView *)rootNode,0);
582  TVolumeView *nextView = 0;
583  // Find itself.
584  while ( (nextView = (TVolumeView *)next(mode)) && nextView != this ){}
585  if (nextView == this) {
586  TVolumePosition *position = next[0];
587  if (!position->GetNode()) {
588  Error("TVolumeView ctor","%s %s ",GetName(),nextView->GetName());
589  }
590  // Calculate the range of the outlined cube verteces.
591  GetLocalRange(globalMin,globalMax);
592  Float_t offSet[3] = {Float_t(position->GetX()),Float_t(position->GetY()),Float_t(position->GetZ())};
593  for (Int_t i=0;i<3;i++) {
594  globalMin[i] += offSet[i];
595  globalMax[i] += offSet[i];
596  }
597  }
598  return next.GetDepth();
599  }
600  else return -1;
601 }
602 
603 ////////////////////////////////////////////////////////////////////////////////
604 /// GetRange
605 ///
606 /// Calculates the size of 3 box the node occupies.
607 /// Return:
608 /// two floating point arrays with the bound of box
609 /// surroundind all shapes of this TModeView
610 ///
611 
613 {
614  TVirtualPad *savePad = gPad;
615  // Create a dummy TPad;
616  TCanvas dummyPad("--Dumm--","dum",1,1);
617  // Assing 3D TView
618  TView *view = TView::CreateView(1,0,0);
619 
622  view->SetAutoRange(kTRUE);
623  Paint("range");
624  view->GetRange(&min[0],&max[0]);
625  delete view;
626  // restore "current pad"
627  if (savePad) savePad->cd();
628 }
629 
630 ////////////////////////////////////////////////////////////////////////////////
631 ///to be documented
632 
634 {
635  if (!gPad) return 0;
636  static char info[512];
637  Double_t x[3] = {0,0,0.5};
638  ((TPad *)gPad)->AbsPixeltoXY(px,py,x[0],x[1]);
639  TView *view =gPad->GetView();
640  if (view) {
641  Double_t min[3], max[3];
642  view->GetRange(min,max);
643  for (int i =0; i<3;i++) min[i] = (max[i]+min[i])/2;
644  view->WCtoNDC(min,max);
645  min[0] = x[0]; min[1] = x[1];
646  min[2] = max[2];
647  view->NDCtoWC(min, x);
648  }
649  TShape *shape = GetShape();
650  if (shape)
651  snprintf(info,512,"%6.2f/%6.2f/%6.2f: %s/%s, shape=%s/%s",x[0],x[1],x[2],GetName(),GetTitle(),shape->GetName(),shape->ClassName());
652  else
653  snprintf(info,512,"%6.2f/%6.2f/%6.2f: %s/%s",x[0],x[1],x[2],GetName(),GetTitle());
654  return info;
655 }
656 
657 ////////////////////////////////////////////////////////////////////////////////
658 ///to be documented
659 
660 TVolumePosition *TVolumeView::Local2Master(const Char_t *localName, const Char_t *masterName)
661 {
662  TVolumeView *masterNode = this;
663  TVolumePosition *position = 0;
664  if (masterName && masterName[0]) masterNode = (TVolumeView *)Find(masterName);
665  if (masterNode) {
666  TVolumeViewIter transform(masterNode,0);
667  if (transform(localName)) position = transform[0];
668  }
669  return position;
670 }
671 
672 ////////////////////////////////////////////////////////////////////////////////
673 ///to be documented
674 
676 {
677  TVolumePosition *position = 0;
678  if (!masterNode) masterNode = this;
679  if (masterNode && localNode) {
680  TVolumeViewIter transform((TVolumeView *)masterNode,0);
681  TVolumeView *nextNode = 0;
682  while ((nextNode = (TVolumeView *)transform()) && nextNode != localNode) { }
683  if (nextNode) position = transform[0];
684  }
685  return position;
686 }
687 
688 ////////////////////////////////////////////////////////////////////////////////
689 ///
690 /// calculate transformation master = (M-local->master )*local + (T-local->master )
691 /// where
692 /// M-local->master - rotation matrix 3 x 3 from the master node to the local node
693 /// T-local->master - trasport vector 3 from the master node to the local node
694 ///
695 /// returns a "master" pointer if transformation has been found
696 /// otherwise 0;
697 ///
698 
700  const Char_t *localName, const Char_t *masterName, Int_t nVector)
701 {
702  Float_t *trans = 0;
703  TVolumePosition *position = 0;
704  TVolumeView *masterNode = this;
705  if (masterName && masterName[0]) masterNode = (TVolumeView *)Find(masterName);
706  if (masterNode) {
707  TVolumeViewIter transform(masterNode,0);
708  if (transform(localName) && (position = (TVolumePosition *) transform.GetPosition()) )
709  trans = position->Local2Master(local,master,nVector);
710  }
711  return trans;
712 }
713 
714 ////////////////////////////////////////////////////////////////////////////////
715 ///
716 /// calculate transformation master = (M-local->master )*local + (T-local->master )
717 /// where
718 /// M-local->master - rotation matrix 3 x 3 from the master node to the local node
719 /// T-local->master - trasport vector 3 from the master node to the local node
720 ///
721 /// returns a "master" pointer if transformation has been found
722 /// otherwise 0;
723 ///
724 
726  const TVolumeView *localNode,
727  const TVolumeView *masterNode, Int_t nVector)
728 {
729  Float_t *trans = 0;
730  TVolumePosition *position = 0;
731  if (!masterNode) masterNode = this;
732  if (masterNode && localNode) {
733  TVolumeViewIter transform((TVolumeView *)masterNode,0);
734  TVolumeView *nextNode = 0;
735  while ((nextNode = (TVolumeView *)transform()) && nextNode != localNode) { }
736  if (nextNode && (position = (TVolumePosition *) transform.GetPosition()) )
737  trans = position->Local2Master(local,master,nVector);
738  }
739  return trans;
740 }
741 
742 ////////////////////////////////////////////////////////////////////////////////
743 /// Paint Referenced node with current parameters.
744 ///
745 /// - vis = 1 (default) shape is drawn
746 /// - vis = 0 shape is not drawn but its sons may be not drawn
747 /// - vis = -1 shape is not drawn. Its sons are not drawn
748 /// - vis = -2 shape is drawn. Its sons are not drawn
749 ///
750 /// It draw the TVolumeView layers from the iFirst one (form the zero) till
751 /// iLast one reached.
752 ///
753 /// restrict the levels for "range" option
754 
756 {
757  Int_t level = gGeometry->GeomLevel();
758  if (!option) return;
759  if (option[0]=='r' && level > 3 ) return;
760 
761  Int_t iFirst = atoi(option);
762  Int_t iLast = 0;
763  const char *delim = strpbrk( option,":-,");
764  if (delim) iLast = atoi(delim+1);
765  if (iLast < iFirst) {
766  iLast = iFirst-1;
767  iFirst = 0;
768  }
769 
770  if ( (0 < iLast) && (iLast < level) ) return;
771 
772  TTablePadView3D *view3D = (TTablePadView3D*)gPad->GetView3D();
773 
774  TVolume *thisNode = 0;
775  TVolumePosition *position = GetPosition();
776 
777  // UpdatePosition does change the current matrix and it MUST be called FIRST !!!
778  if (position) {
779  thisNode = position->GetNode();
780  position->UpdatePosition(option);
781  }
782 
783  // if (option[0] !='r' ) printf(" Level %d first = %d iLast %d \n",level, iFirst, iLast);
784  if (level >= iFirst) {
785  PaintShape(option);
786  if (thisNode) thisNode->PaintShape(option);
787  }
788  // if ( thisNode->TestBit(kSonsInvisible) ) return;
789 
790 //*-*- Paint all sons
791  TSeqCollection *nodes = GetCollection();
792  Int_t nsons = nodes?nodes->GetSize():0;
793 
794  if(!nsons) return;
795 
796  gGeometry->PushLevel();
797  TVolumeView *node;
798  TIter next(nodes);
799  while ((node = (TVolumeView *)next())) {
800  if (view3D) view3D->PushMatrix();
801 
802  node->Paint(option);
803 
804  if (view3D) view3D->PopMatrix();
805  }
806  gGeometry->PopLevel();
807 }
808 
809 ////////////////////////////////////////////////////////////////////////////////
810 /// Paint shape of the node
811 /// To be called from the TObject::Paint method only
812 
814 {
815  Bool_t rangeView = option && option[0]=='r';
816 
817  TIter nextShape(fListOfShapes);
818  TShape *shape = 0;
819  while( (shape = (TShape *)nextShape()) ) {
820  if (!shape->GetVisibility()) continue;
821  if (!rangeView) {
822  TTablePadView3D *view3D = (TTablePadView3D*)gPad->GetView3D();
823  if (view3D)
824  view3D->SetLineAttr(shape->GetLineColor(),shape->GetLineWidth(),option);
825  }
826 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,03,05)
827  // It MUST be the TShape::PAint method:
828  Bool_t viewerWantsSons = kTRUE;
829  TVirtualViewer3D * viewer3D = gPad->GetViewer3D();
830  if (viewer3D) {
831  // We only provide master frame positions in these shapes
832  // so don't ask viewer preference
833 
834  // Ask all shapes for kCore/kBoundingBox/kShapeSpecific
835  // Not all will support the last two - which is fine
836  const TBuffer3D & buffer =
838 
839  // TShape sets buffer id based on TNode * gNode
840  // As we not using TNode we need to override this
841  const_cast<TBuffer3D &>(buffer).fID = this;
842 
843  Int_t reqSections = viewer3D->AddObject(buffer, &viewerWantsSons);
844  if (reqSections != TBuffer3D::kNone) {
845  shape->GetBuffer3D(reqSections);
846  viewer3D->AddObject(buffer);
847  }
848  }
849 #else
850  shape->Paint(option);
851 #endif
852  }
853 }
854 
855 ////////////////////////////////////////////////////////////////////////////////
856 /// return the full path of this data set
857 
859 {
860  TString str;
861  TVolumeView *parent = (TVolumeView *)GetParent();
862  if (parent) {
863  str = parent->PathP();
864  str += "/";
865  }
866  str += GetName();
868  if (p) {
869  char buffer[10];
870  snprintf(buffer,10,";%d",p->GetId());
871  str += buffer;
872  }
873  return str;
874 }
875 
876 ////////////////////////////////////////////////////////////////////////////////
877 ///to be documented
878 
879 void TVolumeView::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
880 {
881  const Char_t *sceleton[] = {
882  "TVolumeView *CreateNodeView(TVolume *topNode) {"
883  ," TString thisNodePath = "
884  ," UInt_t thisPositionId = "
885  ," Double_t thisTranslate[3] = "
886  ," "
887  ," TString matrixName = "
888  ," Int_t matrixType = "
889  ," Double_t thisMatrix[] = { "
890  ," "
891  ," "
892  ," };"
893  ," return = new TVolumeView(thisTranslate, thisMatrix, thisPositionId, topNode,"
894  ," thisNodePath.Data(),matrixName.Data(), matrixType);"
895  ,"}"
896  };
897 //------------------- end of sceleton ---------------------
898  Int_t sceletonSize = sizeof(sceleton)/sizeof(const Char_t*);
899  TVolumePosition *thisPosition = GetPosition();
900  TVolume *thisFullNode = GetNode();
901  TString thisNodePath = thisFullNode ? thisFullNode->Path() : TString("");
902  // Define position
903  UInt_t thisPositionId = thisPosition ? thisPosition->GetId():0;
904  Double_t thisX = thisPosition ? thisPosition->GetX():0;
905  Double_t thisY = thisPosition ? thisPosition->GetY():0;
906  Double_t thisZ = thisPosition ? thisPosition->GetZ():0;
907 
908  const TRotMatrix *matrix = thisPosition ? thisPosition->GetMatrix():0;
909  Int_t matrixType = 2;
910  TString matrixName = " ";
911  Double_t thisMatrix[] = { 0,0,0, 0,0,0, 0,0,0 };
912  if (matrix) {
913  matrixName = matrix->GetName();
914  memcpy(thisMatrix,((TRotMatrix *)matrix)->GetMatrix(),9*sizeof(Double_t));
915  matrixType = matrix->GetType();
916  }
917  Int_t im = 0;
918  for (Int_t lineNumber =0; lineNumber < sceletonSize; lineNumber++) {
919  out << sceleton[lineNumber]; // std::cout << lineNumber << ". " << sceleton[lineNumber];
920  switch (lineNumber) {
921  case 1: out << "\"" << thisNodePath.Data() << "\";" ; // std::cout << "\"" << thisNodePath.Data() << "\";" ;
922  break;
923  case 2: out << thisPositionId << ";" ; // std::cout << "\"" << thisNodePath.Data() << "\";" ;
924  break;
925  case 3: out << "{" << thisX << ", " << thisY << ", "<< thisZ << "};"; // std::cout << thisX << ";" ;
926  break;
927  case 5: out << "\"" << matrixName << "\";" ; // std::cout << "\"" << matrixName << "\";" ;
928  break;
929  case 6: out << matrixType << ";" ; // std::cout << matrixType << ";" ;
930  break;
931  case 7: out << thisMatrix[im++] << ", "; out << thisMatrix[im++] << ", "; out << thisMatrix[im++] << ", ";
932  break;
933  case 8: out << thisMatrix[im++] << ", "; out << thisMatrix[im++] << ", "; out << thisMatrix[im++] << ", ";
934  break;
935  case 9: out << thisMatrix[im++] << ", "; out << thisMatrix[im++] << ", "; out << thisMatrix[im++];
936  break;
937  default:
938  break;
939  };
940 // std::cout << " " << std::endl;
941  out << " " << std::endl;
942  }
943 }
944 
945 ////////////////////////////////////////////////////////////////////////////////
946 ///to be documented
947 
949 {
950  TVolume *thisNode = GetNode();
951  if (thisNode) thisNode->SetLineAttributes();
952 }
953 
954 ////////////////////////////////////////////////////////////////////////////////
955 ///to be documented
956 
958 {
959  TVolume *node = GetNode();
960  if (node) node->SetVisibility(TVolume::ENodeSEEN(vis));
961 }
962 
963 ////////////////////////////////////////////////////////////////////////////////
964 /// Return total size of this 3-D Node with its attributes.
965 
967 {
968  if (GetListOfShapes()) {
969  TIter nextShape(GetListOfShapes());
970  TShape *shape = 0;
971  while( (shape = (TShape *)nextShape()) ) {
972  if (shape->GetVisibility()) shape->Sizeof3D();
973  }
974  }
975 
976  TVolume *thisNode = GetNode();
977  if (thisNode && !(thisNode->GetVisibility()&TVolume::kThisUnvisible) ) {
978  TIter nextShape(thisNode->GetListOfShapes());
979  TShape *shape = 0;
980  while( (shape = (TShape *)nextShape()) ) {
981  if (shape->GetVisibility()) shape->Sizeof3D();
982  }
983  }
984 
985 // if ( TestBit(kSonsInvisible) ) return;
986 
987  TVolumeView *node;
988  TDataSetIter next((TVolumeView *)this);
989  while ((node = (TVolumeView *)next())) node->Sizeof3D();
990 }
virtual Double_t GetY() const
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual UInt_t GetId() const
std::string GetName(const std::string &scope_name)
Definition: Cppyy.cxx:145
virtual Double_t * GetMatrix()
Definition: TRotMatrix.h:54
virtual Int_t GetDepth() const
Definition: TDataSetIter.h:70
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
virtual ENodeSEEN GetVisibility() const
Definition: TVolume.h:82
virtual void SetLineAttributes()
Invoke the DialogCanvas Line attributes.
Definition: TAttLine.cxx:280
TRotMatrix * GetRotMatrix(const char *name) const
Return pointer to RotMatrix with name.
Definition: TGeometry.cxx:356
virtual void Sizeof3D() const
Return total size of this 3-D Node with its attributes.
virtual void WCtoNDC(const Float_t *pw, Float_t *pn)=0
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Computes distance from point (px,py) to the object.
Definition: TObject.cxx:186
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
ENodeSEEN
Definition: TVolume.h:38
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition: TNamed.cxx:140
virtual void SetAutoRange(Bool_t autorange=kTRUE)=0
virtual Bool_t DoOwner(Bool_t done=kTRUE)
Set / Reset the ownerships and returns the previous status of the ownerships.
Definition: TObjectSet.cxx:85
virtual void AddFirst(TObject *obj)
Add object at the beginning of the list.
Definition: TList.cxx:97
virtual TVolume * GetNode() const
to be documented
See TView3D.
Definition: TView.h:25
Use this attribute class when an object should have 3D capabilities.
Definition: TAtt3D.h:19
TList * fListOfShapes
Definition: TVolumeView.h:28
#define gROOT
Definition: TROOT.h:402
Basic string class.
Definition: TString.h:125
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1099
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
virtual Double_t GetZ() const
TList * GetListOfShapes() const
Definition: TVolume.h:80
virtual void SetMatrix(const Double_t *matrix)
copy predefined 3x3 matrix into TRotMatrix object
Definition: TRotMatrix.cxx:216
virtual void NDCtoWC(const Float_t *pn, Float_t *pw)=0
virtual TVolumePosition * GetPosition() const
Definition: TVolumeView.h:52
virtual TVirtualPad * cd(Int_t subpadnumber=0)=0
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:35
virtual void PaintShape(Option_t *option="")
Paint shape of the volume To be called from the TObject::Paint method only.
Definition: TVolume.cxx:635
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections) const
Stub to avoid forcing implementation at this stage.
Definition: TShape.cxx:253
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:105
const TRotMatrix * GetMatrix() const
#define SafeDelete(p)
Definition: RConfig.h:509
Sequenceable collection abstract base class.
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:128
Double_t x[n]
Definition: legend1.C:17
virtual void Add(TDataSet *dataset)
Definition: TVolume.h:97
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a TVolumeView.
void Class()
Definition: Class.C:29
Int_t GeomLevel() const
Definition: TGeometry.h:74
virtual void SetVisibility(ENodeSEEN vis=TVolume::kBothVisible)
Set visibility for this volume and its sons.
Definition: TVolume.cxx:747
virtual TShape * GetShape() const
Definition: TVolumeView.h:83
virtual TList * GetListOfPositions()
Definition: TVolume.h:83
virtual Int_t GetType() const
Definition: TRotMatrix.h:56
virtual TVolumePosition * Local2Master(const TVolumeView *localNode, const TVolumeView *masterNode=0)
to be documented
virtual Int_t PushLevel()
Definition: TGeometry.h:99
EDataSetPass
Definition: TDataSet.h:40
Abstract 3D shapes viewer.
virtual void Paint(Option_t *option="")
Paint Referenced node with current parameters.
static TRotMatrix * GetIdentity()
Return a pointer the "identity" matrix.
Definition: TVolume.cxx:501
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition: TVirtualPad.h:49
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
to be documented
Definition: TVolume.cxx:339
A doubly linked list.
Definition: TList.h:44
TGeometry description.
Definition: TGeometry.h:39
virtual void GetLocalRange(Float_t *min, Float_t *max)
GetRange.
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
virtual Bool_t IsMarked() const
Definition: TVolume.h:98
virtual Double_t GetX(Int_t indx=0) const
virtual void SetObject(TObject *obj)
The depricated method (left here for the sake of the backward compatibility)
Definition: TObjectSet.h:59
This is the base class for all geometry shapes.
Definition: TShape.h:35
Manages a detector rotation matrix.
Definition: TRotMatrix.h:28
virtual TVolume * AddNode(TVolume *node)
Add the TVolume in the Tnode data-structure refered by this TVolumeView object Return TVolume * the i...
virtual Double_t * Local2Master(const Double_t *local, Double_t *master, Int_t nPoints=1) const
Convert one point from local system to master reference system.
virtual char * GetObjectInfo(Int_t px, Int_t py) const
to be documented
virtual void Add(TDataSet *dataset)
Definition: TVolumeView.h:79
static int * ucopy(const int *a, int *b, int n)
Definition: TCernLib.h:325
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
virtual void SetLineAttributes()
to be documented
virtual void UpdatePosition(Option_t *option="")
to be documented
virtual Int_t AddObject(const TBuffer3D &buffer, Bool_t *addChildren=0)=0
unsigned int UInt_t
Definition: RtypesCore.h:42
The most important graphics class in the ROOT system.
Definition: TPad.h:29
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
virtual TDataSet * Find(const char *path) const
Full description see: TDataSetIter::Find.
Definition: TDataSet.cxx:362
Generic 3D primitive description class.
Definition: TBuffer3D.h:17
virtual TString Path() const
return the full path of this data set
Definition: TDataSet.cxx:626
const Bool_t kFALSE
Definition: RtypesCore.h:88
virtual void SetId(UInt_t id)
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
virtual TList * GetListOfShapes() const
Definition: TVolumeView.h:82
virtual void Paint(Option_t *option="")
This method is used only when a shape is painted outside a TNode.
Definition: TShape.cxx:143
The Canvas class.
Definition: TCanvas.h:31
virtual void Browse(TBrowser *b)
to be documented
virtual void PaintShape(Option_t *option)
Paint shape of the node To be called from the TObject::Paint method only.
#define ClassImp(name)
Definition: Rtypes.h:359
double Double_t
Definition: RtypesCore.h:55
virtual void Browse(TBrowser *b)
Browse this dataset (called by TBrowser).
Definition: TObjectSet.cxx:65
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:570
virtual void SetLineAttr(Color_t color, Int_t width, Option_t *opt="")
static float * vadd(const float *b, const float *c, float *a, int n)
Definition: TCernLib.h:380
virtual void SetGeomLevel(Int_t level=0)
Definition: TGeometry.h:104
static TView * CreateView(Int_t system=1, const Double_t *rmin=0, const Double_t *rmax=0)
Create a concrete default 3-d view via the plug-in manager.
Definition: TView.cxx:39
virtual TVolume * GetNode() const
Mother of all ROOT objects.
Definition: TObject.h:37
virtual ~TVolumeView()
default dtor (empty for this class)
virtual void UpdateTempMatrix(Double_t x=0, Double_t y=0, Double_t z=0, TRotMatrix *matrix=0)
Update temp matrix.
Definition: TGeometry.cxx:661
char Char_t
Definition: RtypesCore.h:29
virtual const TVolumePosition * GetPosition(Int_t level=0) const
to be documented
virtual void Draw(Option_t *depth="3")
Draw Referenced node with current parameters.
virtual Bool_t IsMarked() const
Definition: TVolumeView.h:81
TShape * GetShape() const
Definition: TVolume.h:79
virtual void Add(TObject *obj)
Definition: TList.h:87
virtual void GetRange(Float_t *min, Float_t *max)=0
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
#define snprintf
Definition: civetweb.c:822
virtual Bool_t IsOwner() const
Definition: TObjectSet.h:57
#define gPad
Definition: TVirtualPad.h:285
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
to be documented
virtual Int_t GetGlobalRange(const TVolumeView *rootNode, Float_t *min, Float_t *max)
Calculate the position of the vertrex of the outlined cube in repect of the given TVolumeView object...
virtual void Sizeof3D() const
Set total size of this 3D object (used by X3D interface).
Definition: TAtt3D.cxx:27
void Mark()
Definition: TDataSet.h:156
virtual void PopMatrix()
virtual Int_t GetSize() const
Definition: TCollection.h:180
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
R__EXTERN TGeometry * gGeometry
Definition: TGeometry.h:158
const Bool_t kTRUE
Definition: RtypesCore.h:87
virtual Int_t PopLevel()
Definition: TGeometry.h:100
Int_t GetVisibility() const
Definition: TShape.h:58
virtual TSeqCollection * GetCollection() const
Definition: TDataSet.h:105
virtual void PushMatrix()
virtual void SetVisibility(Int_t vis=1)
to be documented
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual TDataSet * GetParent() const
Definition: TDataSet.h:111
const char * Data() const
Definition: TString.h:345
virtual TString PathP() const
return the full path of this data set