source: trunk/MagicSoft/Mars/mtemp/MTelAxisFromStars.cc@ 4796

Last change on this file since 4796 was 4796, checked in by wittek, 20 years ago
*** empty log message ***
File size: 38.9 KB
Line 
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! Author(s): Wolfgang Wittek 07/2004 <mailto:wittek@mppmu.mpg.de>
18!
19! Copyright: MAGIC Software Development, 2000-2004
20!
21!
22\* ======================================================================== */
23
24/////////////////////////////////////////////////////////////////////////////
25//
26// MTelAxisFromStars
27//
28// This task
29// - determines the transformation from expected positions of stars
30// in the camera to measured positions of these stars in the camera
31// - applies this transformation to expected positions of other objects
32// to obtain the estimated positions of these objects in the camera
33// - puts the estimated positions into the relevant containers
34//
35// Input Containers :
36// MStarCam[MStarCam], MStarCamSource[MStarCam]
37//
38// Output Containers :
39// MSkyCamTrans, MSrcPosCam
40//
41/////////////////////////////////////////////////////////////////////////////
42#include <TMath.h>
43#include <TList.h>
44#include <TSystem.h>
45
46#include <fstream>
47
48#include "MTelAxisFromStars.h"
49
50#include "MParList.h"
51
52
53#include "MLog.h"
54#include "MLogManip.h"
55
56#include "MReportDrive.h"
57#include "MPointingPos.h"
58#include "MSrcPosCam.h"
59
60#include "MStarCam.h"
61#include "MStarPos.h"
62#include "MSkyCamTrans.h"
63#include "MStarCamTrans.h"
64
65ClassImp(MTelAxisFromStars);
66
67using namespace std;
68
69// --------------------------------------------------------------------------
70//
71// Constructor
72//
73MTelAxisFromStars::MTelAxisFromStars(const char *name, const char *title)
74{
75 fName = name ? name : "MTelAxisFromStars";
76 fTitle = title ? title : "Calculate source position from star positions";
77
78 // if scale factor fLambda should NOT be fixed set fFixdScaleFactor to
79 // -1.0; otherwise set it to the value requested
80 fFixedScaleFactor = 1.0;
81
82 // if rotation angle fAlfa should NOT be fixed set fFixdRotationAngle to
83 // -1.0; otherwise set it to the requested value
84 fFixedRotationAngle = 0.0;
85
86 // default type of input is : the result of the Gauss fit
87 // type 0 : result from the weighted average
88 // type 1 : result from the Gauss fit
89 fInputType = 1;
90
91 // default value of fAberr
92 // the value 1.07 is valid if the expected position (with aberration)
93 // in the camera is calculated as the average of the positions obtained
94 // from the reflections at the individual mirrors
95 fAberr = 1.07;
96
97}
98
99// --------------------------------------------------------------------------
100//
101// Destructor
102//
103MTelAxisFromStars::~MTelAxisFromStars()
104{
105 delete fStarCamTrans;
106}
107
108// --------------------------------------------------------------------------
109//
110// Set links to containers
111//
112
113Int_t MTelAxisFromStars::PreProcess(MParList *pList)
114{
115 fDrive = (MReportDrive*)pList->FindObject(AddSerialNumber("MReportDrive"));
116 if (!fDrive)
117 {
118 *fLog << err << AddSerialNumber("MReportDrive")
119 << " not found... aborting." << endl;
120 return kFALSE;
121 }
122
123
124 fStarCam = (MStarCam*)pList->FindObject("MStarCam", "MStarCam");
125 if (!fStarCam)
126 {
127 *fLog << err << "MStarCam not found... aborting." << endl;
128 return kFALSE;
129 }
130
131
132 fSourceCam = (MStarCam*)pList->FindObject("MSourceCam", "MStarCam");
133 if (!fSourceCam)
134 {
135 *fLog << warn << "MSourceCam[MStarCam] not found... continue " << endl;
136 }
137
138
139 //-----------------------------------------------------------------
140 fSrcPos = (MSrcPosCam*)pList->FindCreateObj(AddSerialNumber("MSrcPosCam"), "MSrcPosCam");
141 if (!fSrcPos)
142 return kFALSE;
143
144 fPntPos = (MSrcPosCam*)pList->FindCreateObj(AddSerialNumber("MSrcPosCam"), "MPntPosCam");
145 if (!fPntPos)
146 return kFALSE;
147
148 fPointPos = (MPointingPos*)pList->FindCreateObj(AddSerialNumber("MPointingPos"), "MPointingPos");
149 if (!fPointPos)
150 return kFALSE;
151
152 fSourcePos = (MPointingPos*)pList->FindCreateObj(AddSerialNumber("MPointingPos"), "MSourcePos");
153 if (!fSourcePos)
154 return kFALSE;
155
156 fSkyCamTrans = (MSkyCamTrans*)pList->FindCreateObj(AddSerialNumber("MSkyCamTrans"));
157 if (!fSkyCamTrans)
158 return kFALSE;
159
160
161 //-----------------------------------------------------------------
162 // book an MStarCamTrans object
163 // this is needed when calling one of the member functions of MStarCamTrans
164
165 MGeomCam *geom = (MGeomCam*)pList->FindObject(AddSerialNumber("MGeomCam"));
166 if (!geom)
167 {
168 *fLog << err << AddSerialNumber("MGeomCam")
169 << " not found... aborting." << endl;
170 return kFALSE;
171 }
172
173 MObservatory *obs = (MObservatory*)pList->FindObject(AddSerialNumber("MObservatory"));
174 if (!obs)
175 {
176 *fLog << err << AddSerialNumber("MObservatory")
177 << " not found... aborting." << endl;
178 return kFALSE;
179 }
180
181 //-----------------------------------------------------------------
182 fStarCamTrans = new MStarCamTrans(*geom, *obs);
183
184 *fLog << all << "MTelAxisFromStars::Preprocess; the optical aberration factor is set equal to : "
185 << fAberr ;
186
187 *fLog << all << "MTelAxisFromStars::Preprocess; input type is set equal to : "
188 << fInputType ;
189 if (fInputType == 0)
190 *fLog << " (calculated star positions)" << endl;
191 else
192 *fLog << " (fitted star positions)" << endl;
193
194 *fLog << all << "MTelAxisFromStars::Preprocess; scale factor will be fixed at : "
195 << fFixedScaleFactor << endl;
196
197 *fLog << all << "MTelAxisFromStars::Preprocess; rotation angle will be fixed at : "
198 << fFixedRotationAngle << endl;
199
200 return kTRUE;
201}
202
203// --------------------------------------------------------------------------
204//
205// Set optical aberration factor
206//
207// fAberr is the ratio between
208// the distance from the camera center with optical aberration and
209// the distance from the camera center with an ideal imaging
210//
211// fAberr = r_real/r_ideal
212//
213void MTelAxisFromStars::SetOpticalAberr(Double_t aberr)
214{
215 fAberr = aberr;
216}
217
218// --------------------------------------------------------------------------
219//
220// Set the type of the input
221//
222// type = 0 calculated star positions (by averaging)
223// type = 1 fitted star positions (by Gauss fit)
224//
225void MTelAxisFromStars::SetInputType(Int_t type)
226{
227 fInputType = type;
228}
229
230// --------------------------------------------------------------------------
231//
232// Fix the scale factor fLambda
233//
234//
235void MTelAxisFromStars::FixScaleFactorAt(Double_t lambda)
236{
237 fFixedScaleFactor = lambda;
238}
239
240
241// --------------------------------------------------------------------------
242//
243// Fix rotation angle fAlfa
244//
245//
246void MTelAxisFromStars::FixRotationAngleAt(Double_t alfa)
247{
248 fFixedRotationAngle = alfa; // [degrees]
249}
250
251
252// --------------------------------------------------------------------------
253//
254// Process
255//
256// call FindSkyCamTrans to find the Sky-Camera transformation
257// call TransSkyCam to transform some sky directions
258// into the camera system
259// put the estimated source position into MSrcPosCam
260//
261
262Int_t MTelAxisFromStars::Process()
263{
264 //Int_t run = fRun->GetRunNumber();
265 //*fLog << "MTelAxisFromStars::Process; run = " << run << endl;
266
267 //--------------------------------------
268 // Define the input for FindSkyCamTrans
269 //
270
271 // get the expected (axy[0], axy[1]) and the measured positions
272 // (bxy[0], bxy[1]) of stars in the camera from MStarCam
273 Int_t fNumStars = fStarCam->GetNumStars();
274
275 if (fNumStars <= 0)
276 return kTRUE;
277
278 TArrayD axy[2];
279 axy[0].Set(fNumStars);
280 axy[1].Set(fNumStars);
281
282 TArrayD bxy[2];
283 bxy[0].Set(fNumStars);
284 bxy[1].Set(fNumStars);
285
286 // error matrix of bxy
287 TArrayD exy[2][2];
288 exy[0][0].Set(fNumStars);
289 exy[0][1].Set(fNumStars);
290 exy[1][0].Set(fNumStars);
291 exy[1][1].Set(fNumStars);
292
293 // transformation parameters
294 Double_t fLambda;
295 Double_t fAlfa;
296 Double_t fA[2][2];
297 Double_t fD[2];
298 Double_t fErrD[2][2];
299 Int_t fNumIter;
300 Int_t fNdof;
301 Double_t fChi2;
302 Double_t fChi2Prob;
303
304 MStarPos *star = 0;
305 TIter next(fStarCam->GetList());
306 Int_t ix = 0;
307
308 // loop over all stars
309 while ( (star = (MStarPos*)next()) )
310 {
311 axy[0][ix] = star->GetXExp();
312 axy[1][ix] = star->GetYExp();
313
314 if (fInputType == 0)
315 {
316 // values from averaging
317 bxy[0][ix] = star->GetMeanXCalc();
318 bxy[1][ix] = star->GetMeanYCalc();
319
320 // this is the error matrix for (MeanXCalc, MeanYCalc);
321 // this is the error matrix which should be used
322 exy[0][0][ix] = star->GetXXErrCalc();
323 exy[0][1][ix] = star->GetXYErrCalc();
324 exy[1][0][ix] = star->GetXYErrCalc();
325 exy[1][1][ix] = star->GetYYErrCalc();
326
327 //exy[0][0][ix] = star->GetSigmaXCalc()*star->GetSigmaXCalc();
328 //exy[0][1][ix] = 0.0;
329 //exy[1][0][ix] = 0.0;
330 //exy[1][1][ix] = star->GetSigmaYCalc()*star->GetSigmaYCalc();
331 }
332
333 else if (fInputType == 1)
334 {
335 // values from Gauss fit
336 bxy[0][ix] = star->GetMeanXFit();
337 bxy[1][ix] = star->GetMeanYFit();
338
339 // this is the error matrix for (MeanXFit, MeanYFit);
340 // this is the error matrix which should be used
341 exy[0][0][ix] = star->GetXXErrFit();
342 exy[0][1][ix] = star->GetXYErrFit();
343 exy[1][0][ix] = star->GetXYErrFit();
344 exy[1][1][ix] = star->GetYYErrFit();
345
346 // this is the error matrix constructed from SigmaXFit and SigmaYFit;
347 // it is used because the errors above are too small, at present
348 //exy[0][0][ix] = star->GetSigmaXFit() * star->GetSigmaXFit();
349 //exy[0][1][ix] = star->GetCorrXYFit() *
350 // star->GetSigmaXFit() * star->GetSigmaYFit();
351 //exy[1][0][ix] = exy[0][1][ix];
352 //exy[1][1][ix] = star->GetSigmaYFit() * star->GetSigmaYFit();
353 }
354
355 else
356 {
357 *fLog << err << "MTelAxisFromStars::Process; type of input is not defined"
358 << endl;
359 return kFALSE;
360 }
361
362 // don't include stars with undefined error
363 Double_t deter = exy[0][0][ix]*exy[1][1][ix]
364 - exy[0][1][ix]*exy[1][0][ix];
365
366 //*fLog << "ix ,deter, xx, xy, yy = " << ix << ": "
367 // << deter << ", " << exy[0][0][ix] << ", "
368 // << exy[0][1][ix] << ", " << exy[1][1][ix] << endl;
369 if (deter <= 0.0)
370 continue;
371
372 //*fLog << "MTelAxisFromStars : " << endl;
373 //*fLog << " ix, XExp, YExp, XFit, YFit, SigmaX2, SigmaXY, SigmaY2 = "
374 // << ix << " : "
375 // << axy[0][ix] << ", " << axy[1][ix] << ", "
376 // << bxy[0][ix] << ", " << bxy[1][ix] << ", "
377 // << exy[0][0][ix] << ", " << exy[0][1][ix] << ", "
378 // << exy[1][1][ix] << endl;
379
380 ix++;
381 }
382
383 //--------------------------------------
384 // Find the transformation from expected positions (axy[1], axy[2])
385 // to measured positions (bxy[1], bxy[2]) in the camera
386
387 Int_t fNStars = ix;
388
389 if (ix < fNumStars)
390 {
391 // reset the sizes of the arrays
392 Int_t fNStars = ix;
393 axy[0].Set(fNStars);
394 axy[1].Set(fNStars);
395
396 bxy[0].Set(fNStars);
397 bxy[1].Set(fNStars);
398
399 exy[0][0].Set(fNStars);
400 exy[0][1].Set(fNStars);
401 exy[1][0].Set(fNStars);
402 exy[1][1].Set(fNStars);
403 }
404
405 Bool_t fitOK;
406 if (fNStars < 1)
407 {
408 //*fLog << "MTelAxisFromStars::Process; no star for MTelAxisFromStars"
409 // << endl;
410 fitOK = kFALSE;
411 }
412 else
413 {
414 fitOK = FindSkyCamTrans(axy, bxy, exy,
415 fFixedRotationAngle, fFixedScaleFactor, fLambda,
416 fAlfa , fA, fD, fErrD,
417 fNumIter, fNdof, fChi2, fChi2Prob);
418 }
419
420 if (!fitOK && fNStars >= 1)
421 {
422 *fLog << err
423 << "MTelAxisFromStars::Process; Fit to find transformation from star to camera system failed"
424 << endl;
425
426 //if (fNStars > 0)
427 //{
428 // *fLog << err
429 // << " fNumIter, fNdof, fChi2, fChi2Prob = " << fNumIter
430 // << ", " << fNdof << ", " << fChi2 << ", " << fChi2Prob << endl;
431 //}
432
433 return kTRUE;
434 }
435
436
437 //--------------------------------------
438 // Put the transformation parameters into the MSkyCamTrans container
439
440 fSkyCamTrans->SetParameters(fLambda, fAlfa, fA, fD, fErrD,
441 fNumStars, fNumIter, fNdof, fChi2, fChi2Prob);
442 fSkyCamTrans->SetReadyToSave();
443
444 //--------------------------------------
445 // Put the camera position (X, Y)
446 // obtained by transforming the camera center (0, 0)
447 // into MPntPosCam[MSrcPosCam]
448
449 fPntPos->SetXY(fD[0], fD[1]);
450 fPntPos->SetReadyToSave();
451
452
453 //--------------------------------------
454 // Put the corrected pointing position into MPointingPos
455 //
456 SetPointingPosition(fStarCamTrans, fDrive, fSkyCamTrans, fPointPos);
457
458
459 //--------------------------------------
460 // Put the estimated position of the source into SrcPosCam
461 //
462 // get the source direction from MReportDrive
463 // Note : this has to be changed for the wobble mode, where the source
464 // isn't in the center of the camera
465 Double_t decsource = fDrive->GetDec();
466 Double_t rasource = fDrive->GetRa();
467 //
468 Double_t radrive = fDrive->GetRa();
469 Double_t hdrive = fDrive->GetHa();
470 Double_t hsource = hdrive + radrive - rasource;
471 fSourcePos->SetSkyPosition(rasource, decsource, hsource);
472
473 SetSourcePosition(fStarCamTrans, fPointPos, fSourcePos, fSrcPos);
474
475 //*fLog << "after calling SetSourcePosition : , X, Y ,fD = "
476 // << fSrcPos->GetX() << ", " << fSrcPos->GetY() << ", "
477 // << fD[0] << ", " << fD[1] << endl;
478
479 //--------------------------------------
480 // Apply the transformation to some expected positions (asxy[1], asxy[2])
481 // to obtain estimated positions (bsxy[1], bsxy[2]) in the camera
482 // and their error matrices esxy[2][2]
483
484 // get the expected positions (asxy[1], asxy[2]) from another MStarCam
485 // container (with the name "MSourceCam")
486 Int_t fNumStarsSource = 0;
487
488 if (fSourceCam)
489 fNumStarsSource = fSourceCam->GetNumStars();
490
491 //*fLog << "MTelAxisFromStars::Process; fNumStarsSource = "
492 // << fNumStarsSource << endl;
493
494 if (fNumStarsSource > 0)
495 {
496 TArrayD asxy[2];
497 asxy[0].Set(fNumStarsSource);
498 asxy[1].Set(fNumStarsSource);
499
500 TArrayD bsxy[2];
501 bsxy[0].Set(fNumStarsSource);
502 bsxy[1].Set(fNumStarsSource);
503
504 TArrayD esxy[2][2];
505 esxy[0][0].Set(fNumStarsSource);
506 esxy[0][1].Set(fNumStarsSource);
507 esxy[1][0].Set(fNumStarsSource);
508 esxy[1][1].Set(fNumStarsSource);
509
510 MStarPos *starSource = 0;
511 TIter nextSource(fSourceCam->GetList());
512 ix = 0;
513 while ( (starSource = (MStarPos*)nextSource()) )
514 {
515 asxy[0][ix] = starSource->GetXExp();
516 asxy[1][ix] = starSource->GetYExp();
517
518 ix++;
519 }
520
521 TransSkyCam(fLambda, fA, fD, fErrD, asxy, bsxy, esxy);
522
523 // put the estimated positions into the MStarCam container
524 // with name "MSourceCam"
525 TIter setnextSource(fSourceCam->GetList());
526 ix = 0;
527 while ( (starSource = (MStarPos*)setnextSource()) )
528 {
529 Double_t corr = esxy[0][1][ix] /
530 TMath::Sqrt( esxy[0][0][ix] * esxy[1][1][ix] );
531 starSource->SetFitValues(100.0, 100.0, bsxy[0][ix], bsxy[1][ix],
532 TMath::Sqrt(esxy[0][0][ix]), TMath::Sqrt(esxy[1][1][ix]), corr,
533 esxy[0][0][ix], esxy[0][1][ix], esxy[1][1][ix],
534 fChi2, fNdof);
535
536 ix++;
537 }
538
539 }
540
541 //--------------------------------------
542
543 return kTRUE;
544}
545
546
547
548//---------------------------------------------------------------------------
549//
550// FindSkyCamTrans
551//
552// This routine determines the transformation
553//
554// ( cos(alfa) -sin(alfa) )
555// b = lambda * A * a + d A = ( )
556// ^ ^ ^ ( sin(alfa) cos(alfa) )
557// | | |
558// scale rotation shift
559// factor matrix
560//
561// from sky coordinates 'a' (projected onto the camera) to camera
562// coordinates 'b', using the positions of known stars in the camera.
563// The latter positions may have been determined by analysing the
564// DC currents in the different pixels.
565//
566// Input : a[2] x and y coordinates of stars projected onto the camera;
567// they were obtained from (RA, dec) of the stars and
568// (ThetaTel, PhiTel) and the time of observation;
569// these are the 'expected positions' of stars in the camera
570// b[2] 'measured positions' of these stars in the camera;
571// they may have been obtained from the DC currents
572// e[2][2] error matrix of b[2]
573// fixedrotationangle value [in degrees] at which rotation angle
574// alfa should be fixed; -1 means don't fix
575// fixedscalefactor value at which scale factor lambda
576// should be fixed; -1 means don't fix
577//
578// Output : lambda, alfadeg, A[2][2], d[2] fit results;
579// parameters describing the transformation
580// from 'expected positions' to the 'measured
581// positions' in the camera
582// errd[2][2] error matrix of d[2]
583// fNumIter number of iterations
584// fNdoF number of degrees of freedom
585// fChi2 chi-square value
586// fChi2Prob chi-square probability
587//
588// The units are assumed to be
589// [degrees] for alfadeg
590// [mm] for a, b, d
591// [1] for lambda
592
593Bool_t MTelAxisFromStars::FindSkyCamTrans(
594 TArrayD a[2], TArrayD b[2], TArrayD e[2][2],
595 Double_t &fixedrotationang, Double_t &fixedscalefac, Double_t &lambda,
596 Double_t &alfadeg, Double_t A[2][2], Double_t d[2], Double_t errd[2][2],
597 Int_t &fNumIter, Int_t &fNdof, Double_t &fChi2, Double_t &fChi2Prob)
598{
599 Int_t fNumStars = a[0].GetSize();
600
601 //*fLog << "MTelAxisFromStars::FindSkyCamTrans; expected and measured positions :"
602 // << endl;
603 for (Int_t ix=0; ix<fNumStars; ix++)
604 {
605 //*fLog << " ix, a[0], a[1], b[0], b[1], errxx, errxy, erryy = "
606 // << ix << " : "
607 // << a[0][ix] << ", " << a[1][ix] << "; "
608 // << b[0][ix] << ", " << b[1][ix] << "; "
609 // << e[0][0][ix] << ", " << e[0][1][ix] << ", "
610 // << e[1][1][ix] << endl;
611 }
612
613
614 //-------------------------------------------
615 // fix some parameters if the number of degrees of freedom is too low
616 // (<= 0.0)
617
618 Double_t fixedscalefactor = fixedscalefac;
619 Double_t fixedrotationangle = fixedrotationang;
620
621 // calculate number of degrees of freedom
622 fNdof = 2 * fNumStars - 4;
623 if (fixedscalefactor != -1.0)
624 fNdof += 1;
625 if (fixedrotationangle != -1.0)
626 fNdof += 1;
627
628 // if there is only 1 star fix both rotation angle and scale factor
629 if (fNumStars == 1)
630 {
631 if (fixedscalefactor == -1.0)
632 {
633 fixedscalefactor = 1.0;
634 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; scale factor is fixed at "
635 << fixedscalefactor << endl;
636 }
637 if (fixedrotationangle == -1.0)
638 {
639 fixedrotationangle = 0.0;
640 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; rotation angle is fixed at "
641 << fixedrotationangle << endl;
642 }
643 }
644 // otherwise fix only 1 parameter if possible
645 else if (fNdof < 0)
646 {
647 if (fNdof < -2)
648 {
649 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; number of degrees of freedom is too low : "
650 << fNdof << "; fNumStars = " << fNumStars << endl;
651 return kFALSE;
652 }
653 else if (fNdof == -2)
654 {
655 if (fixedscalefactor == -1.0 && fixedrotationangle == -1.0)
656 {
657 fixedscalefactor = 1.0;
658 fixedrotationangle = 0.0;
659 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; scale factor and rotation angle are fixed at "
660 << fixedscalefactor << " and " << fixedrotationangle
661 << " respectively" << endl;
662 }
663 else
664 {
665 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; number of degrees of freedom is too low : "
666 << fNdof << "; fNumStars = " << fNumStars << endl;
667 return kFALSE;
668 }
669 }
670 else if (fNdof == -1)
671 {
672 if (fixedrotationangle == -1.0)
673 {
674 fixedrotationangle = 0.0;
675 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; rotation angle is fixed at "
676 << fixedrotationangle << endl;
677 }
678 else if (fixedscalefactor == -1.0)
679 {
680 fixedscalefactor = 1.0;
681 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; scale factor is fixed at "
682 << fixedscalefactor << endl;
683 }
684 else
685 {
686 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; number of degrees of freedom is too low : "
687 << fNdof << "; fNumStars = " << fNumStars<< endl;
688 return kFALSE;
689 }
690 }
691 }
692
693 // recalculate number of degrees of freedom
694 fNdof = 2 * fNumStars - 4;
695 if (fixedscalefactor != -1.0)
696 fNdof += 1;
697 if (fixedrotationangle != -1.0)
698 fNdof += 1;
699
700 if (fNdof < 0)
701 return kFALSE;
702 //-------------------------------------------
703
704
705 // get first approximation of scaling factor
706 if (fixedscalefactor != -1.0)
707 lambda = fixedscalefactor;
708 else
709 lambda = 1.0;
710
711 Double_t lambdaold = lambda;
712 Double_t dlambda = 0.0;
713
714 // get first approximation of rotation angle
715 Double_t alfa = 0.0;
716 if (fixedrotationangle != -1.0)
717 alfa = fixedrotationangle / TMath::RadToDeg();
718
719 Double_t alfaold = alfa;
720 // maximum allowed change of alfa in 1 iteration step (5 degrees)
721 Double_t dalfamax = 5.0 / TMath::RadToDeg();
722 Double_t dalfa = 0.0;
723
724 Double_t cosal = TMath::Cos(alfa);
725 Double_t sinal = TMath::Sin(alfa);
726
727 A[0][0] = cosal;
728 A[0][1] = -sinal;
729 A[1][0] = sinal;
730 A[1][1] = cosal;
731
732
733 Double_t absdold2 = 10000.0;
734 Double_t fChangeofd2 = 10000.0;
735
736
737 TArrayD Aa[2];
738 Aa[0].Set(fNumStars);
739 Aa[1].Set(fNumStars);
740
741
742 Double_t sumEbminlamAa[2];
743
744 TArrayD Ebminlambracd[2];
745 Ebminlambracd[0].Set(fNumStars);
746 Ebminlambracd[1].Set(fNumStars);
747
748 TArrayD EAa[2];
749 EAa[0].Set(fNumStars);
750 EAa[1].Set(fNumStars);
751
752 // invert the error matrices
753 TArrayD c[2][2];
754 c[0][0].Set(fNumStars);
755 c[0][1].Set(fNumStars);
756 c[1][0].Set(fNumStars);
757 c[1][1].Set(fNumStars);
758
759 for (Int_t ix=0; ix<fNumStars; ix++)
760 {
761 Double_t XX = e[0][0][ix];
762 Double_t XY = e[0][1][ix];
763 Double_t YY = e[1][1][ix];
764
765 // get inverse of error matrix
766 Double_t determ = XX*YY - XY*XY;
767 c[0][0][ix] = YY / determ;
768 c[0][1][ix] = -XY / determ;
769 c[1][0][ix] = -XY / determ;
770 c[1][1][ix] = XX / determ;
771 }
772
773
774
775 // calculate sum of inverted error matrices
776 Double_t determsumc;
777 Double_t sumc[2][2];
778 sumc[0][0] = 0.0;
779 sumc[0][1] = 0.0;
780 sumc[1][0] = 0.0;
781 sumc[1][1] = 0.0;
782
783 for (Int_t ix=0; ix<fNumStars; ix++)
784 {
785 sumc[0][0] += c[0][0][ix];
786 sumc[0][1] += c[0][1][ix];
787 sumc[1][0] += c[1][0][ix];
788 sumc[1][1] += c[1][1][ix];
789 }
790 determsumc = sumc[0][0]*sumc[1][1] - sumc[0][1]*sumc[1][0];
791
792 // calculate inverse of sum of inverted error matrices
793 Double_t sumcinv[2][2];
794 sumcinv[0][0] = sumc[1][1] / determsumc;
795 sumcinv[0][1] = -sumc[0][1] / determsumc;
796 sumcinv[1][0] = -sumc[1][0] / determsumc;
797 sumcinv[1][1] = sumc[0][0] / determsumc;
798
799 //*fLog << "sumcinv = " << sumcinv[0][0] << ", " << sumcinv[0][1]
800 // << ", " << sumcinv[1][1] << endl;
801
802
803 // minimize chi2 by iteration ***** start **********************
804
805 // stop iteration when change in |d|*|d| is less than 'told2'
806 // and change in alfa is less than 'toldalfa'
807 // and change in lambda is less than 'toldlambda'
808 // or chi2 is less than 'tolchi2'
809 Double_t told2 = 0.3*0.3; // [mm*mm]; 1/100 of an inner pixel diameter
810 Double_t toldalfa = 0.01 / TMath::RadToDeg(); // 0.01 degrees
811 Double_t toldlambda = 0.00006; // uncertainty of 1 mm of distance
812 // between camera and reflector
813 Double_t tolchi2 = 1.e-5;
814
815 Int_t fNumIterMax = 100;
816 fNumIter = 0;
817
818 for (Int_t i=0; i<fNumIterMax; i++)
819 {
820 fNumIter++;
821
822 // get next approximation of d ------------------
823 for (Int_t ix=0; ix<fNumStars; ix++)
824 {
825 Aa[0][ix] = A[0][0] * a[0][ix] + A[0][1]*a[1][ix];
826 Aa[1][ix] = A[1][0] * a[0][ix] + A[1][1]*a[1][ix];
827
828 //*fLog << "ix, Aa = " << ix << " : " << Aa[0][ix] << ", "
829 // << Aa[1][ix] << endl;
830 }
831
832 sumEbminlamAa[0] = 0.0;
833 sumEbminlamAa[1] = 0.0;
834
835 for (Int_t ix=0; ix<fNumStars; ix++)
836 {
837 sumEbminlamAa[0] += c[0][0][ix] * (b[0][ix] - lambda*Aa[0][ix])
838 + c[0][1][ix] * (b[1][ix] - lambda*Aa[1][ix]);
839
840 sumEbminlamAa[1] += c[1][0][ix] * (b[0][ix] - lambda*Aa[0][ix])
841 + c[1][1][ix] * (b[1][ix] - lambda*Aa[1][ix]);
842 }
843
844 //*fLog << "sumEbminlamAa = " << sumEbminlamAa[0] << ", "
845 // << sumEbminlamAa[1] << endl;
846
847 d[0] = sumcinv[0][0] * sumEbminlamAa[0]
848 + sumcinv[0][1] * sumEbminlamAa[1] ;
849
850 d[1] = sumcinv[1][0] * sumEbminlamAa[0]
851 + sumcinv[1][1] * sumEbminlamAa[1] ;
852
853 Double_t absdnew2 = d[0]*d[0] + d[1]*d[1];
854 fChangeofd2 = absdnew2 - absdold2;
855
856 //*fLog << "fNumIter : " << fNumIter
857 // << "; alfa, lambda, d[0], d[1], absdold2, absdnew2 = " << endl;
858 //*fLog << alfa << ", " << lambda << ", " << d[0] << ", " << d[1]
859 // << ", " << absdold2 << ", " << absdnew2 << endl;
860
861
862 if ( fabs(fChangeofd2) < told2 && fabs(dalfa) < toldalfa &&
863 fabs(dlambda) < toldlambda )
864 {
865 //*fLog << "Iteration stopped because of small changes : fChangeofd2, dalfa, dlambda = "
866 // << fChangeofd2 << ", " << dalfa << ", " << dlambda << endl;
867 break;
868 }
869 absdold2 = absdnew2;
870
871 // get next approximation of matrix A ----------------
872 if (fFixedRotationAngle == -1.0)
873 {
874 for (Int_t ix=0; ix<fNumStars; ix++)
875 {
876 Ebminlambracd[0][ix] =
877 c[0][0][ix] * ( b[0][ix] - lambda*Aa[0][ix] - d[0] )
878 + c[0][1][ix] * ( b[1][ix] - lambda*Aa[1][ix] - d[1] );
879
880 Ebminlambracd[1][ix] =
881 c[1][0][ix] * ( b[0][ix] - lambda*Aa[0][ix] - d[0] )
882 + c[1][1][ix] * ( b[1][ix] - lambda*Aa[1][ix] - d[1] );
883
884 //*fLog << "ix, Ebminlambracd = " << ix << " : "
885 // << Ebminlambracd[0][ix] << ", "
886 // << Ebminlambracd[1][ix] << endl;
887 }
888
889 // stop iteration if fChi2 is small enough
890 fChi2 = 0.0;
891 for (Int_t ix=0; ix<fNumStars; ix++)
892 {
893 fChi2 += (b[0][ix]-lambda*Aa[0][ix]-d[0] ) * Ebminlambracd[0][ix]
894 + (b[1][ix]-lambda*Aa[1][ix]-d[1] ) * Ebminlambracd[1][ix];
895 }
896 if ( fChi2 < tolchi2 )
897 {
898 //*fLog << "iteration stopped because of small fChi2 : "
899 // << fChi2 << endl;
900 break;
901 }
902
903
904 Double_t dchi2dA[2][2];
905 dchi2dA[0][0] = 0.0;
906 dchi2dA[0][1] = 0.0;
907 dchi2dA[1][0] = 0.0;
908 dchi2dA[1][1] = 0.0;
909
910 for (Int_t ix=0; ix<fNumStars; ix++)
911 {
912 dchi2dA[0][0] += Ebminlambracd[0][ix] * a[0][ix];
913 dchi2dA[0][1] += Ebminlambracd[0][ix] * a[1][ix];
914 dchi2dA[1][0] += Ebminlambracd[1][ix] * a[0][ix];
915 dchi2dA[1][1] += Ebminlambracd[1][ix] * a[1][ix];
916 }
917
918 //*fLog << "dchi2dA = " << dchi2dA[0][0] << ", " << dchi2dA[0][1]
919 // << ", " << dchi2dA[1][0] << ", " << dchi2dA[1][1] << endl;
920
921 // ********* 1st derivative (d chi2) / (d alfa) ************
922 Double_t dchi2dalfa = -2.0*lambda *
923 ( - sinal*(dchi2dA[0][0]+dchi2dA[1][1])
924 + cosal*(dchi2dA[1][0]-dchi2dA[0][1]) );
925
926
927 //Double_t dalfa1st = - fChi2 / dchi2dalfa;
928
929 //*fLog << "fChi2, dchi2dalfa = " << fChi2 << ", "
930 // << dchi2dalfa << endl;
931 //*fLog << "proposed change of alfa using 1st derivative = "
932 // << dalfa1st << endl;
933
934 // ********* 2nd derivative (d2 chi2) / (d alfa2) ******
935 Double_t term1 = 0.0;
936 Double_t term2 = 0.0;
937 Double_t term3 = 0.0;
938 Double_t term4 = 0.0;
939
940 for (Int_t ix=0; ix<fNumStars; ix++)
941 {
942 term1 += a[0][ix]*c[0][0][ix]*a[0][ix] + a[1][ix]*c[1][0][ix]*a[0][ix]
943 + a[0][ix]*c[0][1][ix]*a[1][ix] + a[1][ix]*c[1][1][ix]*a[1][ix];
944
945 term2 += a[0][ix]*c[1][0][ix]*a[0][ix] - a[1][ix]*c[0][0][ix]*a[0][ix]
946 + a[0][ix]*c[1][1][ix]*a[1][ix] - a[1][ix]*c[0][1][ix]*a[1][ix];
947
948 term3 = a[0][ix]*c[0][0][ix]*a[1][ix] + a[1][ix]*c[1][0][ix]*a[1][ix]
949 - a[0][ix]*c[0][1][ix]*a[0][ix] - a[1][ix]*c[1][1][ix]*a[0][ix];
950
951 term4 += a[0][ix]*c[1][0][ix]*a[1][ix] - a[1][ix]*c[0][0][ix]*a[1][ix]
952 - a[0][ix]*c[1][1][ix]*a[0][ix] + a[1][ix]*c[0][1][ix]*a[0][ix];
953 }
954
955 Double_t d2chi2dalfa2 =
956 - 2.0*lambda * ( - cosal*(dchi2dA[0][0]+dchi2dA[1][1])
957 - sinal*(dchi2dA[1][0]-dchi2dA[0][1]) )
958 + 2.0*lambda*lambda * ( sinal*sinal * term1 - sinal*cosal * term2
959 + sinal*cosal * term3 - cosal*cosal * term4);
960
961 // Gauss-Newton step
962 Double_t dalfa2nd = - dchi2dalfa / d2chi2dalfa2;
963
964 //*fLog << "proposed change of alfa using 2st derivative = "
965 // << dalfa2nd << endl;
966
967 //dalfa = dalfa1st;
968 dalfa = dalfa2nd;
969
970 // ******************************************
971
972
973 // restrict change of alfa
974 if ( fabs(dalfa) > dalfamax )
975 {
976 dalfa = TMath::Sign( dalfamax, dalfa );
977 }
978 alfa = alfaold + dalfa;
979
980 if ( alfa < -5.0/TMath::RadToDeg() )
981 alfa = -5.0/TMath::RadToDeg();
982 else if ( alfa > 5.0/TMath::RadToDeg() )
983 alfa = 5.0/TMath::RadToDeg();
984
985 dalfa = alfa - alfaold;
986
987 alfaold = alfa;
988
989 sinal = TMath::Sin(alfa);
990 cosal = TMath::Cos(alfa);
991
992 A[0][0] = cosal;
993 A[0][1] = -sinal;
994 A[1][0] = sinal;
995 A[1][1] = cosal;
996
997 //*fLog << "alfa-alfaold = " << dalfa << endl;
998 //*fLog << "new alfa = " << alfa << endl;
999 }
1000
1001
1002 // get next approximation of lambda ----------------
1003 if (fFixedScaleFactor == -1.0)
1004 {
1005 for (Int_t ix=0; ix<fNumStars; ix++)
1006 {
1007 Aa[0][ix] = A[0][0]*a[0][ix] + A[0][1]*a[1][ix];
1008 Aa[1][ix] = A[1][0]*a[0][ix] + A[1][1]*a[1][ix];
1009
1010 EAa[0][ix] =
1011 c[0][0][ix] * Aa[0][ix] + c[0][1][ix] * Aa[1][ix];
1012 EAa[1][ix] =
1013 c[1][0][ix] * Aa[0][ix] + c[1][1][ix] * Aa[1][ix];
1014
1015 //*fLog << "ix, Aa = " << ix << " : " << Aa[0][ix] << ", "
1016 // << Aa[1][ix] << endl;
1017
1018 //*fLog << "ix, EAa = " << ix << " : " << EAa[0][ix] << ", "
1019 // << EAa[1][ix] << endl;
1020 }
1021
1022 Double_t num = 0.0;
1023 Double_t denom = 0.0;
1024 for (Int_t ix=0; ix<fNumStars; ix++)
1025 {
1026 num += (b[0][ix]-d[0]) * EAa[0][ix]
1027 + (b[1][ix]-d[1]) * EAa[1][ix];
1028
1029 denom += Aa[0][ix] * EAa[0][ix]
1030 + Aa[1][ix] * EAa[1][ix];
1031
1032 //*fLog << "ix : b-d = " << ix << " : " << b[0][ix]-d[0]
1033 // << ", " << b[1][ix]-d[1] << endl;
1034
1035 //*fLog << "ix : Aa = " << ix << " : " << Aa[0][ix]
1036 // << ", " << Aa[1][ix] << endl;
1037 }
1038
1039 lambda = num / denom;
1040
1041 if ( lambda < 0.9 )
1042 lambda = 0.9;
1043 else if ( lambda > 1.1 )
1044 lambda = 1.1;
1045
1046 dlambda = lambda - lambdaold;
1047 lambdaold = lambda;
1048
1049 //*fLog << "num, denom, lambda, dlambda = " << num
1050 // << ", " << denom << ", " << lambda << ", "
1051 // << dlambda << endl;
1052 }
1053
1054 }
1055 //------- end of iteration *****************************************
1056
1057 alfadeg = alfa * TMath::RadToDeg();
1058
1059 // calculate error matrix of d[2]
1060 errd[0][0] = sumcinv[0][0];
1061 errd[0][1] = sumcinv[0][1];
1062 errd[1][0] = sumcinv[1][0];
1063 errd[1][1] = sumcinv[1][1];
1064
1065 // evaluate quality of fit
1066
1067 // calculate chi2
1068 for (Int_t ix=0; ix<fNumStars; ix++)
1069 {
1070 Ebminlambracd[0][ix] =
1071 c[0][0][ix] * ( b[0][ix] - lambda*Aa[0][ix] - d[0] )
1072 + c[0][1][ix] * ( b[1][ix] - lambda*Aa[1][ix] - d[1] );
1073
1074 Ebminlambracd[1][ix] =
1075 c[1][0][ix] * (b[0][ix] - lambda*Aa[0][ix] - d[0] )
1076 + c[1][1][ix] * (b[1][ix] - lambda*Aa[1][ix] - d[1] );
1077 }
1078
1079 fChi2 = 0.0;
1080 for (Int_t ix=0; ix<fNumStars; ix++)
1081 {
1082 fChi2 += (b[0][ix] - lambda*Aa[0][ix] - d[0] ) * Ebminlambracd[0][ix]
1083 + (b[1][ix] - lambda*Aa[1][ix] - d[1] ) * Ebminlambracd[1][ix];
1084 }
1085
1086 fChi2Prob = TMath::Prob(fChi2, fNdof);
1087
1088 //*fLog << "MTelAxisFromStars::FindSkyCamTrans :" << endl;
1089 //*fLog << " fNumStars, fChi2, fNdof, fChi2Prob, fNumIter, fChangeofd2, dalfa, dlambda = "
1090 // << fNumStars << ", " << fChi2 << ", " << fNdof << ", "
1091 // << fChi2Prob << ", "
1092 // << fNumIter << ", " << fChangeofd2 << ", " << dalfa << ", "
1093 // << dlambda << endl;
1094 //*fLog << " lambda, alfadeg, d[0], d[1] = " << lambda << ", "
1095 // << alfadeg << ", " << d[0] << ", " << d[1] << endl;
1096
1097 return kTRUE;
1098}
1099
1100// --------------------------------------------------------------------------
1101//
1102// Apply transformation (lambda, A, d)
1103// to the expected positions (a[1], a[2])
1104// to obtain the estimated positions (b[1], b[2])
1105//
1106// e[2][2] is the error matrix of b[2]
1107
1108void MTelAxisFromStars::TransSkyCam(
1109 Double_t &lambda, Double_t A[2][2], Double_t d[2], Double_t errd[2][2],
1110 TArrayD a[2], TArrayD b[2], TArrayD e[2][2])
1111{
1112 Int_t numpos = a[0].GetSize();
1113 if (numpos <= 0)
1114 return;
1115
1116 //*fLog << "MTelAxisFromStars::TransSkyCam; expected and estimated positions :"
1117 // << endl;
1118
1119 for (Int_t ix=0; ix<numpos; ix++)
1120 {
1121 //*fLog << "MTelAxisFromStars; ix = " << ix << endl;
1122
1123 b[0][ix] = lambda * (A[0][0]*a[0][ix] + A[0][1]*a[1][ix]) + d[0];
1124 b[1][ix] = lambda * (A[1][0]*a[0][ix] + A[1][1]*a[1][ix]) + d[1];
1125
1126 e[0][0][ix] = errd[0][0];
1127 e[0][1][ix] = errd[0][1];
1128 e[1][0][ix] = errd[1][0];
1129 e[1][1][ix] = errd[1][1];
1130
1131 //*fLog << " ix, a[0], a[1], b[0], b[1], errxx, errxy, erryy = "
1132 // << ix << " : "
1133 // << a[0][ix] << ", " << a[1][ix] << "; "
1134 // << b[0][ix] << ", " << b[1][ix] << "; "
1135 // << e[0][0][ix] << ", " << e[0][1][ix] << ", "
1136 // << e[1][1][ix] << endl;
1137 }
1138}
1139
1140// --------------------------------------------------------------------------
1141//
1142// SetPointingPosition
1143//
1144// put the corrected pointing direction into MPointingPos[MPointingPos];
1145// this direction corresponds to the position (0,0) in the camera
1146//
1147
1148Bool_t MTelAxisFromStars::SetPointingPosition(MStarCamTrans *fstarcamtrans,
1149 MReportDrive *fdrive, MSkyCamTrans *ftrans, MPointingPos *fpointpos)
1150{
1151 Double_t decdrive = fdrive->GetDec();
1152 Double_t hdrive = fdrive->GetHa();
1153 Double_t radrive = fdrive->GetRa();
1154
1155 // this is the estimated position (with optical aberration) in the camera
1156 // corresponding to the direction in MReportDrive
1157 Double_t Xpoint = (ftrans->GetShiftD())[0];
1158 Double_t Ypoint = (ftrans->GetShiftD())[1];
1159
1160 // get the sky direction corresponding to the position (0,0) in the camera
1161 Double_t decpoint = 0.0;
1162 Double_t hpoint = 0.0;
1163 fstarcamtrans->CelCamToCel0(decdrive, hdrive,
1164 Xpoint/fAberr, Ypoint/fAberr, decpoint, hpoint);
1165 Double_t rapoint = radrive - hpoint + hdrive;
1166 fpointpos->SetSkyPosition(rapoint, decpoint, hpoint);
1167
1168 // get the local direction corresponding to the position (0,0) in the camera
1169 Double_t thetadrive = fdrive->GetNominalZd();
1170 Double_t phidrive = fdrive->GetNominalAz();
1171 Double_t thetapoint = 0.0;
1172 Double_t phipoint = 0.0;
1173 fstarcamtrans->LocCamToLoc0(thetadrive, phidrive,
1174 Xpoint/fAberr, Ypoint/fAberr, thetapoint, phipoint);
1175 fpointpos->SetLocalPosition(thetapoint, phipoint);
1176 fpointpos->SetReadyToSave();
1177
1178 //*fLog << "SetPointingPosition : decdrive, hdrive, radrive Xpoint, Ypoint = "
1179 // << decdrive << ", " << hdrive << ", " << radrive << ", "
1180 // << Xpoint << ", " << Ypoint << endl;
1181
1182 //*fLog << "SetPointingPosition : thetadrive, phidrive, thetapoint, phipoint = "
1183 // << thetadrive << ", " << phidrive << ", " << thetapoint << ", "
1184 // << phipoint << endl;
1185
1186 return kTRUE;
1187}
1188
1189// --------------------------------------------------------------------------
1190//
1191// SetSourcePosition
1192//
1193// put the estimated position of the source in the camera into
1194// MSrcPosCam[MSrcPosCam]
1195//
1196// and the estimated local direction of the source into
1197// MSourcePos[MPointingPos]
1198//
1199
1200Bool_t MTelAxisFromStars::SetSourcePosition(MStarCamTrans *fstarcamtrans,
1201 MPointingPos *fpointpos, MPointingPos *fsourcepos, MSrcPosCam *fsrcpos)
1202{
1203 // get the corrected pointing direction
1204 // corresponding to the position (0,0) in the camera
1205 Double_t decpoint = fpointpos->GetDec();
1206 Double_t hpoint = fpointpos->GetHa();
1207
1208 // get the sky direction of the source
1209 Double_t decsource = fsourcepos->GetDec();
1210 Double_t hsource = fsourcepos->GetHa();
1211
1212 // get the estimated position (Xsource, Ysource) of the source in the camera;
1213 // this is a position for an ideal imaging, without optical aberration
1214 Double_t Xsource = 0.0;
1215 Double_t Ysource = 0.0;
1216 fstarcamtrans->Cel0CelToCam(decpoint, hpoint,
1217 decsource, hsource, Xsource, Ysource);
1218 fsrcpos->SetXY(Xsource*fAberr, Ysource*fAberr);
1219 fsrcpos->SetReadyToSave();
1220
1221 // get the estimated local direction of the source
1222 Double_t thetapoint = fpointpos->GetZd();
1223 Double_t phipoint = fpointpos->GetAz();
1224 Double_t thetasource = 0.0;
1225 Double_t phisource = 0.0;
1226 fstarcamtrans->Loc0CamToLoc(thetapoint, phipoint,
1227 Xsource, Ysource, thetasource, phisource);
1228 fsourcepos->SetLocalPosition(thetasource, phisource);
1229 fsourcepos->SetReadyToSave();
1230
1231 //*fLog << "SetSourcePosition : decpoint, hpoint, decsource, hsource, Xsource, Ysource = "
1232 // << decpoint << ", " << hpoint << ", " << decsource << ", "
1233 // << hsource << ", " << Xsource << ", " << Ysource << endl;
1234 //*fLog << "SetSourcePosition : thetapoint, phipoint, thetasource, phisource = "
1235 // << thetapoint << ", " << phipoint << ", " << thetasource << ", "
1236 // << phisource << endl;
1237
1238 return kTRUE;
1239}
1240
1241// --------------------------------------------------------------------------
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
Note: See TracBrowser for help on using the repository browser.