| 1 | /* ======================================================================== *\
|
|---|
| 2 | !
|
|---|
| 3 | ! *
|
|---|
| 4 | ! * This file is part of CheObs, the Modular Analysis and Reconstruction
|
|---|
| 5 | ! * Software. It is distributed to you in the hope that it can be a useful
|
|---|
| 6 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
|
|---|
| 7 | ! * It is distributed WITHOUT ANY WARRANTY.
|
|---|
| 8 | ! *
|
|---|
| 9 | ! * Permission to use, copy, modify and distribute this software and its
|
|---|
| 10 | ! * documentation for any purpose is hereby granted without fee,
|
|---|
| 11 | ! * provided that the above copyright notice appears in all copies and
|
|---|
| 12 | ! * that both that copyright notice and this permission notice appear
|
|---|
| 13 | ! * in supporting documentation. It is provided "as is" without express
|
|---|
| 14 | ! * or implied warranty.
|
|---|
| 15 | ! *
|
|---|
| 16 | !
|
|---|
| 17 | !
|
|---|
| 18 | ! Author(s): Thomas Bretz, 1/2009 <mailto:thomas.bretz@epfl.ch>
|
|---|
| 19 | !
|
|---|
| 20 | ! Copyright: CheObs Software Development, 2000-2010
|
|---|
| 21 | !
|
|---|
| 22 | !
|
|---|
| 23 | \* ======================================================================== */
|
|---|
| 24 |
|
|---|
| 25 | //////////////////////////////////////////////////////////////////////////////
|
|---|
| 26 | //
|
|---|
| 27 | // MSimRays
|
|---|
| 28 | //
|
|---|
| 29 | // Task to produce rays from a light source at either infinity or a given
|
|---|
| 30 | // height from a given local sky position.
|
|---|
| 31 | //
|
|---|
| 32 | // The sky position is defined by an MPointingPos object in the parameter
|
|---|
| 33 | // list (if none exists, the source is at the reflector axis). Its
|
|---|
| 34 | // default name is "MPointingPos".
|
|---|
| 35 | //
|
|---|
| 36 | // The height of the light/point source is set by SetHeight in units of km.
|
|---|
| 37 | // A value <= 0 means infinity.
|
|---|
| 38 | //
|
|---|
| 39 | // The number of rays produced per event is defined by SetNumPhotons(n).
|
|---|
| 40 | // The default is 1000.
|
|---|
| 41 | //
|
|---|
| 42 | //////////////////////////////////////////////////////////////////////////////
|
|---|
| 43 | #include "MSimRays.h"
|
|---|
| 44 |
|
|---|
| 45 | #include <TMath.h> // root >=5.20
|
|---|
| 46 | #include <TRandom.h>
|
|---|
| 47 | #include <TRotation.h>
|
|---|
| 48 |
|
|---|
| 49 | #include "MLog.h"
|
|---|
| 50 | #include "MLogManip.h"
|
|---|
| 51 |
|
|---|
| 52 | #include "MParList.h"
|
|---|
| 53 |
|
|---|
| 54 | #include "MSpline3.h"
|
|---|
| 55 | #include "MQuaternion.h"
|
|---|
| 56 |
|
|---|
| 57 | #include "MPhotonEvent.h"
|
|---|
| 58 | #include "MPhotonData.h"
|
|---|
| 59 |
|
|---|
| 60 | #include "MOptics.h"
|
|---|
| 61 | #include "MPointingPos.h"
|
|---|
| 62 |
|
|---|
| 63 | ClassImp(MSimRays);
|
|---|
| 64 |
|
|---|
| 65 | using namespace std;
|
|---|
| 66 |
|
|---|
| 67 | // --------------------------------------------------------------------------
|
|---|
| 68 | //
|
|---|
| 69 | // Default Constructor.
|
|---|
| 70 | //
|
|---|
| 71 | MSimRays::MSimRays(const char* name, const char *title)
|
|---|
| 72 | : fEvt(0), fReflector(0), fPointPos(0), fSource(0),
|
|---|
| 73 | fNumPhotons(1000), fHeight(-1), fWavelengthMin(-1), fWavelengthMax(-1),
|
|---|
| 74 | fRandomDist(0), fNameReflector("MReflector"), fNamePointPos("MPointingPos"),
|
|---|
| 75 | fNameSource("Source")
|
|---|
| 76 | {
|
|---|
| 77 | fName = name ? name : "MSimRays";
|
|---|
| 78 | fTitle = title ? title : "Task to calculate reflection os a mirror";
|
|---|
| 79 | }
|
|---|
| 80 |
|
|---|
| 81 | MSimRays::~MSimRays()
|
|---|
| 82 | {
|
|---|
| 83 | delete fRandomDist;
|
|---|
| 84 | }
|
|---|
| 85 |
|
|---|
| 86 | // --------------------------------------------------------------------------
|
|---|
| 87 | //
|
|---|
| 88 | // Search for the necessary parameter containers.
|
|---|
| 89 | //
|
|---|
| 90 | Int_t MSimRays::PreProcess(MParList *pList)
|
|---|
| 91 | {
|
|---|
| 92 | fEvt = (MPhotonEvent*)pList->FindCreateObj("MPhotonEvent");
|
|---|
| 93 | if (!fEvt)
|
|---|
| 94 | return kFALSE;
|
|---|
| 95 |
|
|---|
| 96 | if (!pList->FindCreateObj("MCorsikaEvtHeader"))
|
|---|
| 97 | return kFALSE;
|
|---|
| 98 |
|
|---|
| 99 | fReflector = (MOptics*)pList->FindObject(fNameReflector, "MOptics");
|
|---|
| 100 | if (!fReflector)
|
|---|
| 101 | {
|
|---|
| 102 | *fLog << inf << fNameReflector << " [MOptics] not found..." << endl;
|
|---|
| 103 | return kFALSE;
|
|---|
| 104 | }
|
|---|
| 105 |
|
|---|
| 106 | fSource = (MPointingPos*)pList->FindObject(fNameSource, "MPointingPos");
|
|---|
| 107 | if (!fSource)
|
|---|
| 108 | {
|
|---|
| 109 | // *fLog << inf << fNameSource << " [MPointingPos] not found..." << endl;
|
|---|
| 110 | // return kFALSE;
|
|---|
| 111 | }
|
|---|
| 112 |
|
|---|
| 113 | fPointPos = (MPointingPos*)pList->FindObject(fNamePointPos, "MPointingPos");
|
|---|
| 114 | if (!fPointPos)
|
|---|
| 115 | {
|
|---|
| 116 | *fLog << inf << fNamePointPos << " [MPointingPos] not found..." << endl;
|
|---|
| 117 | return kFALSE;
|
|---|
| 118 | }
|
|---|
| 119 |
|
|---|
| 120 | return kTRUE;
|
|---|
| 121 | }
|
|---|
| 122 |
|
|---|
| 123 | // --------------------------------------------------------------------------
|
|---|
| 124 | //
|
|---|
| 125 | // Converts the photons into the telscope coordinate frame using the
|
|---|
| 126 | // pointing position from MPointingPos.
|
|---|
| 127 | //
|
|---|
| 128 | // Reflects all photons on all mirrors and stores the final photons on
|
|---|
| 129 | // the focal plane. Also intermediate photons are stored for debugging.
|
|---|
| 130 | //
|
|---|
| 131 | Int_t MSimRays::Process()
|
|---|
| 132 | {
|
|---|
| 133 | // Get arrays from event container
|
|---|
| 134 | fEvt->Resize(fNumPhotons);
|
|---|
| 135 |
|
|---|
| 136 | TClonesArray &arr = fEvt->GetArray();
|
|---|
| 137 |
|
|---|
| 138 | const Int_t num = arr.GetEntriesFast();
|
|---|
| 139 |
|
|---|
| 140 | const Double_t maxr = fReflector->GetMaxR();
|
|---|
| 141 |
|
|---|
| 142 | const Double_t deltazd = fSource ? fSource->GetZdRad() : 0;
|
|---|
| 143 | const Double_t deltaaz = fSource ? fSource->GetAzRad() : 0;
|
|---|
| 144 |
|
|---|
| 145 | const Double_t zd = fPointPos->GetZdRad() + deltazd;
|
|---|
| 146 | const Double_t az = fPointPos->GetAzRad() + deltaaz;
|
|---|
| 147 |
|
|---|
| 148 |
|
|---|
| 149 | // cm -> m
|
|---|
| 150 | // s -> ns
|
|---|
| 151 | // length -> time
|
|---|
| 152 | const Double_t conv = 1./(TMath::C()*100/1e9);
|
|---|
| 153 |
|
|---|
| 154 | // Local sky coordinates (direction of telescope axis)
|
|---|
| 155 | //const Double_t zd = fPointing->GetZdRad(); // x==north
|
|---|
| 156 | //const Double_t az = fPointing->GetAzRad();
|
|---|
| 157 |
|
|---|
| 158 | // Height of point source [cm] (0 means infinity)
|
|---|
| 159 | const Double_t h = fHeight * 100000;
|
|---|
| 160 |
|
|---|
| 161 | // Rotation matrix to derotate sky
|
|---|
| 162 | // For the new coordinate system see the Wiki
|
|---|
| 163 | TRotation rot; // The signs are positive because we align the incident point on ground to the telescope axis
|
|---|
| 164 | rot.RotateX( zd); // Rotate point on ground to align it with the telescope axis
|
|---|
| 165 | rot.RotateZ(-az); // tilt the point from ground to make it parallel to the mirror plane
|
|---|
| 166 |
|
|---|
| 167 | Int_t idx = 0;
|
|---|
| 168 | while (idx<num)
|
|---|
| 169 | {
|
|---|
| 170 | MPhotonData &dat = *static_cast<MPhotonData*>(arr.UncheckedAt(idx));
|
|---|
| 171 |
|
|---|
| 172 | Double_t x, y;
|
|---|
| 173 | if (fHeight<0)
|
|---|
| 174 | {
|
|---|
| 175 | // Parallel light
|
|---|
| 176 | // --------------
|
|---|
| 177 | const Double_t r = gRandom->Uniform();
|
|---|
| 178 | gRandom->Circle(x, y, maxr*TMath::Sqrt(r));
|
|---|
| 179 | }
|
|---|
| 180 | else
|
|---|
| 181 | {
|
|---|
| 182 | // Point source
|
|---|
| 183 | // ------------
|
|---|
| 184 | // Adapted from: http://mathworld.wolfram.com/SpherePointPicking.html
|
|---|
| 185 | // Note that theta and phi is exchanged!
|
|---|
| 186 |
|
|---|
| 187 | // The maximum zenith angle is theta=atan(maxr/h)
|
|---|
| 188 | // cos(theta) = cos(atan(maxr/h)) = 1/sqrt(1+maxr^2/h^2)
|
|---|
| 189 | const double min_cost = 1./TMath::Sqrt(1.+maxr*maxr/h/h);
|
|---|
| 190 | const double cos_theta = gRandom->Uniform(min_cost, 1);
|
|---|
| 191 |
|
|---|
| 192 | gRandom->Circle(x, y, h*TMath::Sqrt(1./cos_theta/cos_theta - 1));
|
|---|
| 193 |
|
|---|
| 194 | // const double cos_theta = gRandom->Uniform(ct, 1);
|
|---|
| 195 | // const double sin_theta = TMath::Sqrt(1.-cos_theta*cos_theta);
|
|---|
| 196 | // gRandom->Circle(x, y, h*sin_theta/cos_theta);
|
|---|
| 197 |
|
|---|
| 198 | // Homogeneous on a sphere
|
|---|
| 199 | // const double phi = TMath::TwoPi() * gRandom->Uniform();
|
|---|
| 200 | // x = sin_theta * cos(phi);
|
|---|
| 201 | // y = sin_theta * sin(phi);
|
|---|
| 202 | // z = cos_theta;
|
|---|
| 203 |
|
|---|
| 204 | // Project the photons to a plane at z=1
|
|---|
| 205 | // x /= cos_theta;
|
|---|
| 206 | // y /= cos_theta;
|
|---|
| 207 | // z /= cos_theta; // z = 1
|
|---|
| 208 |
|
|---|
| 209 | // The radius of the sphere is h
|
|---|
| 210 | // x *= h;
|
|---|
| 211 | // y *= h;
|
|---|
| 212 | // z *= h; // z = h
|
|---|
| 213 | }
|
|---|
| 214 |
|
|---|
| 215 | // The is the incident direction of the photon
|
|---|
| 216 | // h==0 means infinitiy
|
|---|
| 217 | const TVector3 u = fHeight>0 ? TVector3(x, y, -h).Unit() : TVector3(0, 0, -1);
|
|---|
| 218 |
|
|---|
| 219 | // w is pointing away from the direction the photon comes from
|
|---|
| 220 | // CORSIKA-orig: x(north), y(west), z(up), t(time)
|
|---|
| 221 | // NOW: x(east), y(north), z(up), t(time)
|
|---|
| 222 | MQuaternion p(TVector3(x, y, 0), fHeight>0 ? TMath::Sqrt(x*x + y*y + h*h): 0);
|
|---|
| 223 | MQuaternion w(u, conv);
|
|---|
| 224 |
|
|---|
| 225 | // Rotate the coordinates into the reflector's coordinate system.
|
|---|
| 226 | // It is assumed that the z-plane is parallel to the focal plane.
|
|---|
| 227 | // (The reflector coordinate system is defined by the telescope orientation)
|
|---|
| 228 | p *= rot;
|
|---|
| 229 | w *= rot;
|
|---|
| 230 |
|
|---|
| 231 | // Now propagate the photon to the z-plane in the new coordinate system
|
|---|
| 232 | p.PropagateZ0(w);
|
|---|
| 233 |
|
|---|
| 234 | // Shift the coordinate system to the telescope. Corsika's
|
|---|
| 235 | // coordinate system is always w.r.t. to the particle axis
|
|---|
| 236 | //p += impact;
|
|---|
| 237 |
|
|---|
| 238 | // Store new position and direction in the reflector's coordinate frame
|
|---|
| 239 | dat.SetPosition(p);
|
|---|
| 240 | dat.SetDirection(w);
|
|---|
| 241 |
|
|---|
| 242 | if (fRandomDist)
|
|---|
| 243 | {
|
|---|
| 244 | dat.SetWavelength(fRandomDist->FindX(gRandom->Uniform()));
|
|---|
| 245 | }
|
|---|
| 246 | else
|
|---|
| 247 | {
|
|---|
| 248 | if (fWavelengthMin>0 && fWavelengthMax>0)
|
|---|
| 249 | dat.SimWavelength(fWavelengthMin, fWavelengthMax);
|
|---|
| 250 | }
|
|---|
| 251 |
|
|---|
| 252 | idx++;
|
|---|
| 253 | }
|
|---|
| 254 |
|
|---|
| 255 | // Doesn't seem to be too time consuming. But we could also sort later!
|
|---|
| 256 | // (after cones, inside the camera)
|
|---|
| 257 | // fEvt->Sort(kTRUE);
|
|---|
| 258 |
|
|---|
| 259 | return kTRUE;
|
|---|
| 260 | }
|
|---|
| 261 |
|
|---|
| 262 | void MSimRays::SetDistribution(const MSpline3 &s)
|
|---|
| 263 | {
|
|---|
| 264 | delete fRandomDist;
|
|---|
| 265 | fRandomDist = new MSpline3(s.GetIntegralSpline());
|
|---|
| 266 | }
|
|---|
| 267 |
|
|---|
| 268 | bool MSimRays::ReadDistribution(const char *filename, const char *fmt)
|
|---|
| 269 | {
|
|---|
| 270 | TGraph g(filename, fmt);
|
|---|
| 271 | if (g.GetN()<2)
|
|---|
| 272 | return false;
|
|---|
| 273 |
|
|---|
| 274 | SetDistribution(MSpline3(g));
|
|---|
| 275 | return true;
|
|---|
| 276 | }
|
|---|
| 277 |
|
|---|
| 278 | // --------------------------------------------------------------------------
|
|---|
| 279 | //
|
|---|
| 280 | // Height: -1
|
|---|
| 281 | // NumPhotons: 1000
|
|---|
| 282 | //
|
|---|
| 283 | Int_t MSimRays::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
|
|---|
| 284 | {
|
|---|
| 285 | Bool_t rc = kFALSE;
|
|---|
| 286 | if (IsEnvDefined(env, prefix, "Height", print))
|
|---|
| 287 | {
|
|---|
| 288 | rc = kTRUE;
|
|---|
| 289 | fHeight = GetEnvValue(env, prefix, "Height", fHeight);
|
|---|
| 290 | }
|
|---|
| 291 | if (IsEnvDefined(env, prefix, "NumPhotons", print))
|
|---|
| 292 | {
|
|---|
| 293 | rc = kTRUE;
|
|---|
| 294 | fNumPhotons = GetEnvValue(env, prefix, "NumPhotons", (Int_t)fNumPhotons);
|
|---|
| 295 | }
|
|---|
| 296 |
|
|---|
| 297 | return rc;
|
|---|
| 298 | }
|
|---|