Logo ROOT   6.16/01
Reference Guide
TFoamCell.cxx
Go to the documentation of this file.
1// @(#)root/foam:$Id$
2// Author: S. Jadach <mailto:Stanislaw.jadach@ifj.edu.pl>, P.Sawicki <mailto:Pawel.Sawicki@ifj.edu.pl>
3
4//_________________________________________________________________________________
5//
6// Class TFoamCell used in TFoam
7// ==============================
8// Objects of this class are hyper-rectangular cells organized in the binary tree.
9// Special algorithm for encoding relative positioning of the cells
10// allow to save total memory allocation needed for the system of cells.
11//
12//_________________________________________________________________________________
13
14#include "Riostream.h"
15#include "TFoamCell.h"
16#include "TFoamVect.h"
17
18
20
21////////////////////////////////////////////////////////////////////////////////
22/// Default constructor for streamer
23
25{
26 fParent = 0;
27 fDaught0 = 0;
28 fDaught1 = 0;
29}
30
31////////////////////////////////////////////////////////////////////////////////
32/// User constructor allocating single empty Cell
33
35{
36 if ( kDim >0) {
37 //---------=========----------
38 fDim = kDim;
39 fSerial = 0;
40 fStatus = 1;
41 fParent = 0;
42 fDaught0 = 0;
43 fDaught1 = 0;
44 fXdiv = 0.0;
45 fBest = 0;
46 fVolume = 0.0;
47 fIntegral = 0.0;
48 fDrive = 0.0;
49 fPrimary = 0.0;
50 } else
51 Error("TFoamCell","Dimension has to be >0 \n ");
52}
53
54////////////////////////////////////////////////////////////////////////////////
55/// Copy constructor (not tested!)
56
58{
59 Error("TFoamCell", "+++++ NEVER USE Copy constructor for TFoamCell \n");
60 fStatus = From.fStatus;
61 fParent = From.fParent;
62 fDaught0 = From.fDaught0;
63 fDaught1 = From.fDaught1;
64 fXdiv = From.fXdiv;
65 fBest = From.fBest;
66 fVolume = From.fVolume;
67 fIntegral = From.fIntegral;
68 fDrive = From.fDrive;
69 fPrimary = From.fPrimary;
70}
71
72////////////////////////////////////////////////////////////////////////////////
73/// Destructor
74
76{
77}
78
79////////////////////////////////////////////////////////////////////////////////
80/// Substitution operator = (never used)
81
83{
84 Info("TFoamCell", "operator=\n ");
85 if (&From == this) return *this;
86 fStatus = From.fStatus;
87 fParent = From.fParent;
88 fDaught0 = From.fDaught0;
89 fDaught1 = From.fDaught1;
90 fXdiv = From.fXdiv;
91 fBest = From.fBest;
92 fVolume = From.fVolume;
93 fIntegral = From.fIntegral;
94 fDrive = From.fDrive;
95 fPrimary = From.fPrimary;
96 return *this;
97}
98
99
100////////////////////////////////////////////////////////////////////////////////
101/// Fills in certain data into newly allocated cell
102
103void TFoamCell::Fill(Int_t Status, TFoamCell *Parent, TFoamCell *Daugh1, TFoamCell *Daugh2)
104{
105 fStatus = Status;
106 fParent = Parent;
107 fDaught0 = Daugh1;
108 fDaught1 = Daugh2;
109}
110
111////////////////////////////////////////////////////////////////////////////////
112// GETTERS/SETTERS
113////////////////////////////////////////////////////////////////////////////////
114
115////////////////////////////////////////////////////////////////////////////////
116/// Provides size and position of the cell
117/// These parameter are calculated by analyzing information in all parents
118/// cells up to the root cell. It takes time but saves memory.
119
120void TFoamCell::GetHcub( TFoamVect &cellPosi, TFoamVect &cellSize) const
121{
122 if(fDim<1) return;
123 const TFoamCell *pCell,*dCell;
124 cellPosi = 0.0; cellSize=1.0; // load all components
125 dCell = this;
126 while(dCell != 0) {
127 pCell = dCell->GetPare();
128 if( pCell== 0) break;
129 Int_t kDiv = pCell->fBest;
130 Double_t xDivi = pCell->fXdiv;
131 if(dCell == pCell->GetDau0() ) {
132 cellSize[kDiv] *=xDivi;
133 cellPosi[kDiv] *=xDivi;
134 } else if( dCell == pCell->GetDau1() ) {
135 cellSize[kDiv] *=(1.0-xDivi);
136 cellPosi[kDiv] =cellPosi[kDiv]*(1.0-xDivi)+xDivi;
137 } else {
138 Error("GetHcub ","Something wrong with linked tree \n");
139 }
140 dCell=pCell;
141 }//while
142}//GetHcub
143
144////////////////////////////////////////////////////////////////////////////////
145/// Provides size of the cell
146/// Size parameters are calculated by analyzing information in all parents
147/// cells up to the root cell. It takes time but saves memory.
148
149void TFoamCell::GetHSize( TFoamVect &cellSize) const
150{
151 if(fDim<1) return;
152 const TFoamCell *pCell,*dCell;
153 cellSize=1.0; // load all components
154 dCell = this;
155 while(dCell != 0) {
156 pCell = dCell->GetPare();
157 if( pCell== 0) break;
158 Int_t kDiv = pCell->fBest;
159 Double_t xDivi = pCell->fXdiv;
160 if(dCell == pCell->GetDau0() ) {
161 cellSize[kDiv]=cellSize[kDiv]*xDivi;
162 } else if(dCell == pCell->GetDau1() ) {
163 cellSize[kDiv]=cellSize[kDiv]*(1.0-xDivi);
164 } else {
165 Error("GetHSize ","Something wrong with linked tree \n");
166 }
167 dCell=pCell;
168 }//while
169}//GetHSize
170
171////////////////////////////////////////////////////////////////////////////////
172/// Calculates volume of the cell using size params which are calculated
173
175{
176 Int_t k;
177 Double_t volu=1.0;
178 if(fDim>0) { // h-cubical subspace
179 TFoamVect cellSize(fDim);
180 GetHSize(cellSize);
181 for(k=0; k<fDim; k++) volu *= cellSize[k];
182 }
183 fVolume =volu;
184}
185
186////////////////////////////////////////////////////////////////////////////////
187/// Printout of the cell geometry parameters for the debug purpose
188
189void TFoamCell::Print(Option_t *option) const
190{
191 if(!option) Error("Print", "No option set\n");
192
193 std::cout << " Status= "<< fStatus <<",";
194 std::cout << " Volume= "<< fVolume <<",";
195 std::cout << " TrueInteg= " << fIntegral <<",";
196 std::cout << " DriveInteg= "<< fDrive <<",";
197 std::cout << " PrimInteg= " << fPrimary <<",";
198 std::cout<< std::endl;
199 std::cout << " Xdiv= "<<fXdiv<<",";
200 std::cout << " Best= "<<fBest<<",";
201 std::cout << " Parent= {"<< (GetPare() ? GetPare()->GetSerial() : -1) <<"} "; // extra DEBUG
202 std::cout << " Daught0= {"<< (GetDau0() ? GetDau0()->GetSerial() : -1 )<<"} "; // extra DEBUG
203 std::cout << " Daught1= {"<< (GetDau1() ? GetDau1()->GetSerial() : -1 )<<"} "; // extra DEBUG
204 std::cout<< std::endl;
205 //
206 //
207 if(fDim>0 ) {
208 TFoamVect cellPosi(fDim); TFoamVect cellSize(fDim);
209 GetHcub(cellPosi,cellSize);
210 std::cout <<" Posi= "; cellPosi.Print("1"); std::cout<<","<< std::endl;
211 std::cout <<" Size= "; cellSize.Print("1"); std::cout<<","<< std::endl;
212 }
213}
214///////////////////////////////////////////////////////////////////
215// End of class TFoamCell //
216///////////////////////////////////////////////////////////////////
int Int_t
Definition: RtypesCore.h:41
double Double_t
Definition: RtypesCore.h:55
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:363
const Int_t kDiv
Definition: TPCON.h:30
Double_t fXdiv
Definition: TFoamCell.h:34
virtual ~TFoamCell()
Destructor.
Definition: TFoamCell.cxx:75
TFoamCell * GetDau1() const
Definition: TFoamCell.h:72
void CalcVolume()
Calculates volume of the cell using size params which are calculated.
Definition: TFoamCell.cxx:174
void Print(Option_t *option) const
Printout of the cell geometry parameters for the debug purpose.
Definition: TFoamCell.cxx:189
TRef fDaught1
Definition: TFoamCell.h:31
Int_t GetSerial() const
Definition: TFoamCell.h:76
Double_t fIntegral
Definition: TFoamCell.h:38
Int_t fStatus
Definition: TFoamCell.h:28
TFoamCell * GetDau0() const
Definition: TFoamCell.h:71
Double_t fVolume
Definition: TFoamCell.h:37
TRef fParent
Definition: TFoamCell.h:29
TRef fDaught0
Definition: TFoamCell.h:30
void GetHcub(TFoamVect &, TFoamVect &) const
Provides size and position of the cell These parameter are calculated by analyzing information in all...
Definition: TFoamCell.cxx:120
TFoamCell * GetPare() const
Definition: TFoamCell.h:70
Double_t fDrive
Definition: TFoamCell.h:39
Short_t fDim
Definition: TFoamCell.h:23
Double_t fPrimary
Definition: TFoamCell.h:40
TFoamCell()
Default constructor for streamer.
Definition: TFoamCell.cxx:24
TFoamCell & operator=(const TFoamCell &)
Substitution operator = (never used)
Definition: TFoamCell.cxx:82
void Fill(Int_t, TFoamCell *, TFoamCell *, TFoamCell *)
Fills in certain data into newly allocated cell.
Definition: TFoamCell.cxx:103
Int_t fBest
Definition: TFoamCell.h:35
Int_t fSerial
Definition: TFoamCell.h:27
void GetHSize(TFoamVect &) const
Provides size of the cell Size parameters are calculated by analyzing information in all parents cell...
Definition: TFoamCell.cxx:149
void Print(Option_t *option) const
Printout of all vector components on "std::cout".
Definition: TFoamVect.cxx:204
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:854