| 1 | /* ======================================================================== *\
|
|---|
| 2 | !
|
|---|
| 3 | ! *
|
|---|
| 4 | ! * This file is part of MARS, the MAGIC 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 appear 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 12/2000 <mailto:tbretz@uni-sw.gwdg.de>
|
|---|
| 19 | ! Author(s): Harald Kornmayer 1/2001
|
|---|
| 20 | ! Author(s): Rudolf Bock 10/2001 <mailto:Rudolf.Bock@cern.ch>
|
|---|
| 21 | ! Author(s): Wolfgang Wittek 06/2002 <mailto:wittek@mppmu.mpg.de>
|
|---|
| 22 | !
|
|---|
| 23 | ! Copyright: MAGIC Software Development, 2000-2002
|
|---|
| 24 | !
|
|---|
| 25 | !
|
|---|
| 26 | \* ======================================================================== */
|
|---|
| 27 |
|
|---|
| 28 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 29 | //
|
|---|
| 30 | // MHillasSrc
|
|---|
| 31 | //
|
|---|
| 32 | // Storage Container for image parameters
|
|---|
| 33 | //
|
|---|
| 34 | // source-dependent image parameters
|
|---|
| 35 | //
|
|---|
| 36 | //
|
|---|
| 37 | // Version 1:
|
|---|
| 38 | // ----------
|
|---|
| 39 | // fAlpha angle between major axis and line source-to-center
|
|---|
| 40 | // fDist distance from source to center of ellipse
|
|---|
| 41 | //
|
|---|
| 42 | //
|
|---|
| 43 | // Version 2:
|
|---|
| 44 | // ----------
|
|---|
| 45 | // fHeadTail added
|
|---|
| 46 | //
|
|---|
| 47 | //
|
|---|
| 48 | // Version 3:
|
|---|
| 49 | // ----------
|
|---|
| 50 | // fCosDeltaAlpha cosine of angle between d and a, where
|
|---|
| 51 | // - d is the vector from the source position to the
|
|---|
| 52 | // center of the ellipse
|
|---|
| 53 | // - a is a vector along the main axis of the ellipse,
|
|---|
| 54 | // defined with positive x-component
|
|---|
| 55 | //
|
|---|
| 56 | //
|
|---|
| 57 | // Version 4:
|
|---|
| 58 | // ----------
|
|---|
| 59 | //
|
|---|
| 60 | // fHeadTail removed
|
|---|
| 61 | //
|
|---|
| 62 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 63 | #include "MHillasSrc.h"
|
|---|
| 64 |
|
|---|
| 65 | #include <fstream>
|
|---|
| 66 | #include <TArrayF.h>
|
|---|
| 67 |
|
|---|
| 68 | #include "MLog.h"
|
|---|
| 69 | #include "MLogManip.h"
|
|---|
| 70 |
|
|---|
| 71 | #include "MSrcPosCam.h"
|
|---|
| 72 |
|
|---|
| 73 | ClassImp(MHillasSrc);
|
|---|
| 74 |
|
|---|
| 75 | using namespace std;
|
|---|
| 76 |
|
|---|
| 77 | // --------------------------------------------------------------------------
|
|---|
| 78 | //
|
|---|
| 79 | // Default constructor.
|
|---|
| 80 | //
|
|---|
| 81 | MHillasSrc::MHillasSrc(const char *name, const char *title)
|
|---|
| 82 | {
|
|---|
| 83 | fName = name ? name : "MHillasSrc";
|
|---|
| 84 | fTitle = title ? title : "Parameters depending in source position";
|
|---|
| 85 | }
|
|---|
| 86 |
|
|---|
| 87 | void MHillasSrc::Reset()
|
|---|
| 88 | {
|
|---|
| 89 | fDist = -1;
|
|---|
| 90 | fAlpha = 0;
|
|---|
| 91 | fCosDeltaAlpha = 0;
|
|---|
| 92 | }
|
|---|
| 93 |
|
|---|
| 94 | // --------------------------------------------------------------------------
|
|---|
| 95 | //
|
|---|
| 96 | // Calculation of source-dependent parameters
|
|---|
| 97 | // In case you don't call Calc from within an eventloop make sure, that
|
|---|
| 98 | // you call the Reset member function before.
|
|---|
| 99 | //
|
|---|
| 100 | Bool_t MHillasSrc::Calc(const MHillas *hillas)
|
|---|
| 101 | {
|
|---|
| 102 | const Double_t mx = hillas->GetMeanX(); // [mm]
|
|---|
| 103 | const Double_t my = hillas->GetMeanY(); // [mm]
|
|---|
| 104 |
|
|---|
| 105 | const Double_t sx = mx - fSrcPos->GetX(); // [mm]
|
|---|
| 106 | const Double_t sy = my - fSrcPos->GetY(); // [mm]
|
|---|
| 107 |
|
|---|
| 108 | const Double_t sd = hillas->GetSinDelta(); // [1]
|
|---|
| 109 | const Double_t cd = hillas->GetCosDelta(); // [1]
|
|---|
| 110 |
|
|---|
| 111 | //
|
|---|
| 112 | // Distance from source position to center of ellipse.
|
|---|
| 113 | // If the distance is 0 distance, Alpha is not specified.
|
|---|
| 114 | // The calculation has failed and returnes kFALSE.
|
|---|
| 115 | //
|
|---|
| 116 | const Double_t dist = sqrt(sx*sx + sy*sy); // [mm]
|
|---|
| 117 | if (dist==0)
|
|---|
| 118 | return kFALSE;
|
|---|
| 119 |
|
|---|
| 120 | //
|
|---|
| 121 | // Calculate Alpha and Cosda = cos(d,a)
|
|---|
| 122 | // The sign of Cosda will be used for quantities containing
|
|---|
| 123 | // a head-tail information
|
|---|
| 124 | //
|
|---|
| 125 | // *OLD* const Double_t arg = (sy-tand*sx) / (dist*sqrt(tand*tand+1));
|
|---|
| 126 | // *OLD* fAlpha = asin(arg)*kRad2Deg;
|
|---|
| 127 |
|
|---|
| 128 | const Double_t arg1 = cd*sy-sd*sx; // [mm]
|
|---|
| 129 | const Double_t arg2 = cd*sx+sd*sy; // [mm]
|
|---|
| 130 |
|
|---|
| 131 | fAlpha = asin(arg1/dist)*kRad2Deg; // [deg]
|
|---|
| 132 | fCosDeltaAlpha = arg2/dist; // [1]
|
|---|
| 133 | fDist = dist; // [mm]
|
|---|
| 134 |
|
|---|
| 135 | SetReadyToSave();
|
|---|
| 136 |
|
|---|
| 137 | return kTRUE;
|
|---|
| 138 | }
|
|---|
| 139 |
|
|---|
| 140 | void MHillasSrc::Print(Option_t *) const
|
|---|
| 141 | {
|
|---|
| 142 | *fLog << all;
|
|---|
| 143 | *fLog << "Source dependant Image Parameters (" << GetName() << ")" << endl;
|
|---|
| 144 | *fLog << " - Dist [mm] = " << fDist << endl;
|
|---|
| 145 | *fLog << " - Alpha [deg] = " << fAlpha << endl;
|
|---|
| 146 | *fLog << " - CosDeltaAlpha = " << fCosDeltaAlpha << endl;
|
|---|
| 147 | }
|
|---|
| 148 | // --------------------------------------------------------------------------
|
|---|
| 149 | //
|
|---|
| 150 | // This function is ment for special usage, please never try to set
|
|---|
| 151 | // values via this function
|
|---|
| 152 | //
|
|---|
| 153 | void MHillasSrc::Set(const TArrayF &arr)
|
|---|
| 154 | {
|
|---|
| 155 | if (arr.GetSize() != 3)
|
|---|
| 156 | return;
|
|---|
| 157 |
|
|---|
| 158 | fAlpha = arr.At(0); // [deg] angle of major axis with vector to src
|
|---|
| 159 | fDist = arr.At(1); // [mm] distance between src and center of ellipse
|
|---|
| 160 | fCosDeltaAlpha = arr.At(2); // [1] cosine of angle between d and a
|
|---|
| 161 | }
|
|---|