testRSQLManager.cc
Go to the documentation of this file.
1 // Test giving out infos on Stations
2 
3 #include <iostream>
4 #include <sstream>
5 #include <string>
6 #include <vector>
7 #include <det/Detector.h>
8 #include <rdet/RDetector.h>
9 #include <rdet/Station.h>
10 #include <rdet/Channel.h>
11 #include <sdet/SDetector.h>
12 #include <sdet/Station.h>
13 #include <fwk/CentralConfig.h>
14 #include <utl/TimeStamp.h>
15 #include <utl/UTCDateTime.h>
16 #include <cppunit/extensions/HelperMacros.h>
17 #include <tst/Verify.h>
18 #include <utl/Point.h>
19 #include <utl/UTMPoint.h>
20 #include <utl/ReferenceEllipsoid.h>
21 #include <fwk/CoordinateSystemRegistry.h>
22 #include <fwk/LocalCoordinateSystem.h>
23 #include <TROOT.h>
24 #include <TH2F.h>
25 #include <TFile.h>
26 #include <TCanvas.h>
27 #include <TPaveText.h>
28 #include <TAttText.h>
29 
30 using namespace det;
31 using namespace std;
32 using namespace utl;
33 using namespace tst;
34 using namespace rdet;
35 
36 
37 #define ASSERT_CLOSE(x, y, eps) CPPUNIT_ASSERT(Verify<CloseTo>(x, y, eps))
38 #define ASSERT_EQUAL(x, y) CPPUNIT_ASSERT(Verify<Equal>(x, y))
39 
40 
41 class RSQLManagerTest : public CppUnit::TestFixture {
42 
43  CPPUNIT_TEST_SUITE(RSQLManagerTest);
44  CPPUNIT_TEST(testStationListManager);
45  CPPUNIT_TEST(testBeaconManager);
46  CPPUNIT_TEST_SUITE_END();
47 
48 public:
49  void setUp()
50  {
51  ErrorLogger::GetInstance().SetVerbosity(Verbosity::eVerbose);
52 
53  // CentralConfig must be created first.
54  //
55  /*CentralConfig* fCentralConfig =*/
56  fwk::CentralConfig::GetInstance (BOOTSTRAPFILE);
57 
58  // Create a detector
59  //
60  /*Detector& theDetector =*/
61  Detector::GetInstance();
62  }
63 
65  {
66  Detector& detector = Detector::GetInstance();
67  TimeStamp detTime = UTCDateTime(2010, 10, 1, 16, 45, 12).GetTimeStamp();
68  detector.Update(detTime);
69 
70  // try to get the name of the RDetector.Station id=1
71  const RDetector& RDetector = detector.GetRDetector();
72  const Station& RStation = RDetector.GetStation(101);
73  cout << "Name of rdet::Station with Id 101: " << RStation.GetName() << endl;
74 
75  utl::Point pos1 = RDetector.GetStation(101).GetPosition();
76  utl::Point pos2 = RDetector.GetStation(110).GetPosition();
77  utl::Point pos3 = RDetector.GetStation(112).GetPosition();
78 
79  const UTMPoint utmpos1(pos1, ReferenceEllipsoid::GetWGS84());
80  const UTMPoint utmpos2(pos2, ReferenceEllipsoid::GetWGS84());
81  const UTMPoint utmpos3(pos3, ReferenceEllipsoid::GetWGS84());
82 
83  ASSERT_CLOSE(utmpos1.GetEasting(), 450885.827, 1e-8);
84  ASSERT_CLOSE(utmpos1.GetNorthing(), 6114721.732, 1e-8);
85  ASSERT_CLOSE(utmpos1.GetHeight(), 1566.993, 1e-6);
86 
87  ASSERT_CLOSE(utmpos2.GetEasting(), 450755.970, 1e-8);
88  ASSERT_CLOSE(utmpos2.GetNorthing(), 6114509.317, 1e-8);
89  ASSERT_CLOSE(utmpos2.GetHeight(), 1567.973, 1e-6);
90 
91  ASSERT_CLOSE(utmpos3.GetEasting(), 451251.865, 1e-8);
92  ASSERT_CLOSE(utmpos3.GetNorthing(), 6114503.188, 1e-8);
93  ASSERT_CLOSE(utmpos3.GetHeight(), 1562.512, 1e-6);
94 
95  ASSERT_CLOSE(RDetector.GetStation(101).GetChannel(1).GetSamplingFrequency(), 0.2, 1e-4);
96  ASSERT_CLOSE(RDetector.GetStation(110).GetChannel(1).GetSamplingFrequency(), 0.2, 1e-4);
97 
98  cout << "update detector" << endl;
99  detTime = UTCDateTime(2011, 5, 1, 16, 45, 12).GetTimeStamp();
100  detector.Update(detTime);
101  pos1 = RDetector.GetStation(101).GetPosition();
102 
103  cout << "update detector" << endl;
104  detTime = UTCDateTime(2012, 5, 1, 16, 45, 12).GetTimeStamp();
105  detector.Update(detTime);
106 
107  pos1 = RDetector.GetStation(101).GetPosition();
108  pos2 = RDetector.GetStation(110).GetPosition();
109  pos3 = RDetector.GetStation(112).GetPosition();
110 
111  const UTMPoint utmpos4(pos1, ReferenceEllipsoid::GetWGS84());
112  const UTMPoint utmpos5(pos2, ReferenceEllipsoid::GetWGS84());
113  const UTMPoint utmpos6(pos3, ReferenceEllipsoid::GetWGS84());
114 
115  ASSERT_CLOSE(utmpos4.GetEasting(), 450885.827, 1e-8);
116  ASSERT_CLOSE(utmpos4.GetNorthing(), 6114721.732, 1e-8);
117  ASSERT_CLOSE(utmpos4.GetHeight(), 1566.993, 1e-6);
118 
119  ASSERT_CLOSE(utmpos5.GetEasting(), 450755.970, 1e-8);
120  ASSERT_CLOSE(utmpos5.GetNorthing(), 6114509.317, 1e-8);
121  ASSERT_CLOSE(utmpos5.GetHeight(), 1567.973, 1e-6);
122 
123  ASSERT_CLOSE(utmpos6.GetEasting(), 451251.865, 1e-8);
124  ASSERT_CLOSE(utmpos6.GetNorthing(), 6114503.188, 1e-8);
125  ASSERT_CLOSE(utmpos6.GetHeight(), 1562.512, 1e-6);
126 
127  ASSERT_CLOSE(RDetector.GetStation(101).GetChannel(1).GetSamplingFrequency(), 0.18, 1e-4);
128  ASSERT_CLOSE(RDetector.GetStation(110).GetChannel(1).GetSamplingFrequency(), 0.18, 1e-4);
129  }
130 
132  {
133  Detector& detector = Detector::GetInstance();
134  const TimeStamp detTime = UTCDateTime(2013, 11, 1, 16, 45, 12).GetTimeStamp();
135  detector.Update(detTime);
136  const RDetector& rDetector = detector.GetRDetector();
137  const std::vector<double>& beaconFrequencies = rDetector.GetBeaconFrequencies();
138  std::cout << "RDetector::GetBeaconFrequencies() returns frequencies (in GHz):" << std::endl;
139  for (const auto& f : beaconFrequencies)
140  std::cout << f << std::endl;
141  // test Rdetector::GetBeaconReferencePhase
142  const double beaconRefPhase = rDetector.GetBeaconReferencePhase(101, 0.061523);
143  ASSERT_CLOSE(beaconRefPhase, 0.261232, 1e-5);
144  }
145 
146 };
147 
148 
void Update(const utl::TimeStamp &time, const bool invData=true, const bool invComp=true, const bool forceRadio=false)
Update detector: deletes currently constructed stations and sets new time.
Definition: Detector.cc:179
double GetBeaconReferencePhase(const int stationId, const double beaconFreq) const
Get beacon reference phases for one station.
Definition: RDetector.cc:197
Point object.
Definition: Point.h:32
Detector description interface for Station-related data.
Class to hold and convert a point in geodetic coordinates.
Definition: UTMPoint.h:40
Traditional name.
Definition: Verbosity.h:17
double GetNorthing() const
Get the northing.
Definition: UTMPoint.h:206
CPPUNIT_TEST_SUITE_REGISTRATION(testAiresShowerFile)
const std::vector< double > & GetBeaconFrequencies() const
Get vector of Beacon Frequencies for given detector-time from Manager.
Definition: RDetector.cc:183
A TimeStamp holds GPS second and nanosecond for some event.
Definition: TimeStamp.h:110
Detector description interface for RDetector-related data.
Definition: RDetector.h:46
double GetHeight() const
Get the height.
Definition: UTMPoint.h:212
Top of the hierarchy of the detector description interface.
Definition: Detector.h:81
double GetEasting() const
Get the easting.
Definition: UTMPoint.h:209
const Channel & GetChannel(const int id) const
Get specified Channel by id.
#define ASSERT_CLOSE(x, y, eps)
static CentralConfig * GetInstance()
Use this the first time you get an instance of central configuration.
const rdet::RDetector & GetRDetector() const
Definition: Detector.cc:143
void testStationListManager()
double GetSamplingFrequency() const
Get sampling Frequency of ADC (unit?)
utl::Point GetPosition() const
Tank position in Site Cartesian Coordinates.
const Station & GetStation(const int stationId) const
Get station by Station Id.
Definition: RDetector.cc:141
TimeStamp GetTimeStamp() const
Definition: UTCDateTime.cc:115
std::string GetName() const
Station name.

, generated on Tue Sep 26 2023.