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

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