Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
schroedinger_hydrogen.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_graphics
3/// \notebook
4/// \preview Plot the Amplitude of a Hydrogen Atom.
5///
6/// Visualize the Amplitude of a Hydrogen Atom in the n = 2, l = 0, m = 0 state.
7/// Demonstrates how TH2F can be used in Quantum Mechanics.
8///
9/// The formula for Hydrogen in this energy state is \f$ \psi_{200} = \frac{1}{4\sqrt{2\pi}a_0
10/// ^{\frac{3}{2}}}(2-\frac{\sqrt{x^2+y^2}}{a_0})e^{-\frac{\sqrt{x^2+y^2}}{2a_0}} \f$
11///
12/// \macro_image
13/// \macro_code
14///
15/// \author Advait Dhingra
16
17#include <cmath>
18
19double WaveFunction(double x, double y)
20{
21 double r = sqrt(x * x + y * y);
22
23 double w = (1 / pow((4 * sqrt(2 * TMath::Pi()) * 1), 1.5)) *
24 (2 - (r / 1) * pow(TMath::E(), (-1 * r) / 2)); // Wavefunction formula for psi 2,0,0
25
26 return w * w; // Amplitude
27}
28
29void schroedinger_hydrogen()
30{
31 TH2F *h2D = new TH2F("Hydrogen Atom",
32 "Hydrogen in n = 2, l = 0, m = 0 state; Position in x direction; Position in y direction", 200,
33 -10, 10, 200, -10, 10);
34
35 for (float i = -10; i < 10; i += 0.01) {
36 for (float j = -10; j < 10; j += 0.01) {
37 h2D->Fill(i, j, WaveFunction(i, j));
38 }
39 }
40
41 gStyle->SetPalette(kCividis);
42 gStyle->SetOptStat(0);
43
44 TCanvas *c1 = new TCanvas("c1", "Schroedinger's Hydrogen Atom", 750, 1500);
45 c1->Divide(1, 2);
46
47 auto c1_1 = c1->cd(1);
48 c1_1->SetRightMargin(0.14);
49 h2D->GetXaxis()->SetLabelSize(0.03);
50 h2D->GetYaxis()->SetLabelSize(0.03);
51 h2D->GetZaxis()->SetLabelSize(0.03);
52 h2D->SetContour(50);
53 h2D->Draw("colz");
54
55 TLatex *l = new TLatex(
56 -10, -12.43,
57 "The Electron is more likely to be found in the yellow areas and less likely to be found in the blue areas.");
58 l->SetTextFont(42);
59 l->SetTextSize(0.02);
60 l->Draw();
61
62 auto c1_2 = c1->cd(2);
63 c1_2->SetTheta(42.);
64
65 TH2D *h2Dc = (TH2D *)h2D->Clone();
66 h2Dc->SetTitle("3D view of probability amplitude;;");
67 h2Dc->Draw("surf2");
68}
ROOT::R::TRInterface & r
Definition Object.C:4
@ kCividis
Definition TColor.h:139
externTStyle * gStyle
Definition TStyle.h:442
virtual void SetLabelSize(Float_t size=0.04)
Set size of axis labels.
Definition TAttAxis.cxx:184
The Canvas class.
Definition TCanvas.h:23
TAxis * GetZaxis()
Definition TH1.h:573
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6836
TAxis * GetXaxis()
Definition TH1.h:571
TAxis * GetYaxis()
Definition TH1.h:572
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3097
virtual void SetContour(Int_t nlevels, const Double_t *levels=nullptr)
Set the number and values of contour levels.
Definition TH1.cxx:8620
TObject * Clone(const char *newname="") const override
Make a complete copy of the underlying object.
Definition TH1.cxx:2786
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:400
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:345
Int_t Fill(Double_t) override
Invalid Fill method.
Definition TH2.cxx:363
RVec< PromoteTypes< T0, T1 > > pow(const T0 &x, const RVec< T1 > &v)
Definition RVec.hxx:1842
return c1
Definition legend1.C:41
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
constexpr Double_t E()
Base of natural log: .
Definition TMath.h:96
constexpr Double_t Pi()
Definition TMath.h:40
TLine l
Definition textangle.C:4