Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
glvox2.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_gl
3/// This macro demonstrates how to use "glcol" option for TH3
4/// and how to create user defined TRANSFER FUNCTION:
5/// transfer function maps bin value to voxel's opacity.
6/// codomain is [0, 1] (1. - non-transparent, 0.5 is semitransparent, etc.)
7/// To pass transparency function into painting algorithm, you have to:
8/// 1. Create TF1 object (with symbolic expression like "0.5 * (sin(x) + 1)":
9/// ~~~{.cpp}
10/// ...
11/// TF1 * tf = new TF1("TransferFunction", "0.5 * (sin(x) + 1)", -10., 10.);
12/// ...
13/// ~~~
14/// IMPORTANT, the name of TF1 object MUST be "TransferFunction".
15/// 2. Add this function into a hist's list of functions:
16/// ~~~{.cpp}
17/// ...
18/// TList * lof = hist->GetListOfFunctions();
19/// if (lof) lof->Add(tf);
20/// ...
21/// ~~~
22/// It's also possible to use your own function and pass it into TF1, please read
23/// TF1 documentation to learn how.
24///
25/// This macro is to be compiled: TF1 is extremely slow with interpreted function
26/// as an argument.
27///
28/// \macro_image(nobatch)
29/// \macro_code
30///
31/// \author Timur Pocheptsov
32
33#include "TStyle.h"
34#include "TList.h"
35#include "TH3.h"
36#include "TF1.h"
37
38namespace {
39
40Double_t my_transfer_function(const Double_t *x, const Double_t * /*param*/)
41{
42 // Bin values in our example range from -2 to 1.
43 // Let's make values from -2. to -1.5 more transparent:
44 if (*x < -1.5)
45 return 0.008;
46
47 if (*x < -0.5)
48 return 0.015;
49
50 if (*x < 0.)
51 return 0.02;
52
53 if (*x < 0.5)
54 return 0.03;
55
56 if (*x < 0.8)
57 return 0.04;
58
59 return 0.05;
60}
61
62} // namespace
63
64void glvox2()
65{
66 // Create and fill TH3.
67 const UInt_t nX = 30;
68 const Double_t xMin = -1., xMax = 1., xStep = (xMax - xMin) / (nX - 1);
69
70 const UInt_t nY = 30;
71 const Double_t yMin = -1., yMax = 1., yStep = (yMax - yMin) / (nY - 1);
72
73 const UInt_t nZ = 30;
74 const Double_t zMin = -1., zMax = 1., zStep = (zMax - zMin) / (nZ - 1);
75
76 TH3F *hist = new TH3F("glvoxel", "glvoxel", nX, -1., 1., nY, -1., 1., nZ, -1., 1.);
77
78 // Fill the histogram to create a "sphere".
79 for (UInt_t i = 0; i < nZ; ++i) {
80 const Double_t z = zMin + i * zStep;
81
82 for (UInt_t j = 0; j < nY; ++j) {
83 const Double_t y = yMin + j * yStep;
84
85 for (UInt_t k = 0; k < nX; ++k) {
86 const Double_t x = xMin + k * xStep;
87
88 const Double_t val = 1. - (x * x + y * y + z * z);
89 hist->SetBinContent(k + 1, j + 1, i + 1, val);
90 }
91 }
92 }
93
94 // Now, specify the transfer function.
95 TList *lf = hist->GetListOfFunctions();
96 if (lf) {
97 TF1 *tf = new TF1("TransferFunction", my_transfer_function, -2., 1.);
98 lf->Add(tf);
99 }
100
102
103 hist->Draw("glcol");
104}
unsigned int UInt_t
Definition RtypesCore.h:46
double Double_t
Definition RtypesCore.h:59
R__EXTERN TStyle * gStyle
Definition TStyle.h:436
1-Dim function class
Definition TF1.h:233
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3068
TList * GetListOfFunctions() const
Definition TH1.h:245
3-D histogram with a float per channel (see TH1 documentation)
Definition TH3.h:318
void SetBinContent(Int_t bin, Double_t content) override
Set bin content.
Definition TH3.cxx:3465
A doubly linked list.
Definition TList.h:38
void Add(TObject *obj) override
Definition TList.h:81
void SetCanvasPreferGL(Bool_t prefer=kTRUE)
Definition TStyle.h:341
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17