/* ======================================================================== *\ ! ! * ! * This file is part of MARS, the MAGIC Analysis and Reconstruction ! * Software. It is distributed to you in the hope that it can be a useful ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. ! * It is distributed WITHOUT ANY WARRANTY. ! * ! * Permission to use, copy, modify and distribute this software and its ! * documentation for any purpose is hereby granted without fee, ! * provided that the above copyright notice appear in all copies and ! * that both that copyright notice and this permission notice appear ! * in supporting documentation. It is provided "as is" without express ! * or implied warranty. ! * ! ! ! Author(s): Thomas Bretz 12/2000 ! Author(s): Harald Kornmayer 1/2001 ! Author(s): Rudolf Bock 10/2001 ! Author(s): Wolfgang Wittek 06/2002 ! ! Copyright: MAGIC Software Development, 2000-2002 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MHillasSrc // // Storage Container for image parameters // // source-dependent image parameters // // Version 1: // ---------- // fAlpha angle between major axis and line source-to-center // fDist distance from source to center of ellipse // // Version 2: // ---------- // fHeadTail // // Version 3: // ---------- // fCosDeltaLenth cosine of angle between d and a, where // - d is the vector from the source position to the // center of the ellipse // - a is a vector along the main axis of the ellipse // ///////////////////////////////////////////////////////////////////////////// #include "MHillasSrc.h" #include #include #include "MLog.h" #include "MLogManip.h" #include "MSrcPosCam.h" ClassImp(MHillasSrc); // -------------------------------------------------------------------------- // // Default constructor. // MHillasSrc::MHillasSrc(const char *name, const char *title) { fName = name ? name : "MHillasSrc"; fTitle = title ? title : "Parameters depending in source position"; } void MHillasSrc::Reset() { fDist = -1; fAlpha = 0; fHeadTail = 0; fCosDeltaAlpha = 0; } // -------------------------------------------------------------------------- // // Calculation of source-dependent parameters // In case you don't call Calc from within an eventloop make sure, that // you call the Reset member function before. // Bool_t MHillasSrc::Calc(const MHillas *hillas) { fHillas = hillas; const Double_t mx = GetMeanX(); // [mm] const Double_t my = GetMeanY(); // [mm] const Double_t sx = mx - fSrcPos->GetX(); // [mm] const Double_t sy = my - fSrcPos->GetY(); // [mm] const Double_t sd = sin(GetDelta()); // [1] const Double_t cd = cos(GetDelta()); // [1] const Double_t tand = tan(GetDelta()); // [1] fHeadTail = cd*sx + sd*sy; // [mm] // // Distance from source position to center of ellipse. // If the distance is 0 distance, Alpha is not specified. // The calculation has failed and returnes kFALSE. // Double_t dist = sqrt(sx*sx + sy*sy); // [mm] if (dist==0) { //*fLog << warn << GetDescriptor() << ": Event has Dist==0... skipped." << endl; return kFALSE; } // // Calculate Alpha and Cosda = cos(d,a) // The sign of Cosda will be used for quantities containing // a head-tail information // const Double_t arg = (sy-tand*sx) / (dist*sqrt(tand*tand+1)); fAlpha = asin(arg)*kRad2Deg; // [deg] fCosDeltaAlpha = fHeadTail/dist; // [1] fDist = dist; // [mm] SetReadyToSave(); return kTRUE; } void MHillasSrc::Print(Option_t *) const { *fLog << all; *fLog << "Source dependant Image Parameters (" << GetName() << ")" << endl; *fLog << " - Dist [mm] = " << fDist << endl; *fLog << " - Alpha [deg] = " << fAlpha << endl; *fLog << " - HeadTail [mm] = " << fHeadTail << endl; *fLog << " - CosDeltaAlpha = " << fCosDeltaAlpha << endl; } // -------------------------------------------------------------------------- // // This function is ment for special usage, please never try to set // values via this function // void MHillasSrc::Set(const TArrayF &arr) { if (arr.GetSize() != 4) return; fAlpha = arr.At(0); // [deg] angle of major axis with vector to src fDist = arr.At(1); // [mm] distance between src and center of ellipse fHeadTail = arr.At(2); // [mm] fCosDeltaAlpha = arr.At(3); // [1] cosine of angle between d and a } // ----------------------------------------------------------------------- // // overloaded MParContainer to read MHillasSrc from an ascii file // /* void MHillasSrc::AsciiRead(ifstream &fin) { fin >> fAlpha; fin >> fDist; fin >> fHeadTail; } */ // ----------------------------------------------------------------------- // // overloaded MParContainer to write MHillasSrc to an ascii file /* void MHillasSrc::AsciiWrite(ofstream &fout) const { fout << fAlpha << " " << fDist; } */