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

Last change on this file since 5200 was 4798, checked in by wittek, 20 years ago
*** empty log message ***
File size: 39.4 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 if (fInputType == 1)
532 {
533 starSource->SetFitValues(100.0, 100.0, bsxy[0][ix], bsxy[1][ix],
534 TMath::Sqrt(esxy[0][0][ix]), TMath::Sqrt(esxy[1][1][ix]), corr,
535 esxy[0][0][ix], esxy[0][1][ix], esxy[1][1][ix],
536 fChi2, fNdof);
537 }
538 else
539 {
540 starSource->SetCalcValues(100.0, 100.0, bsxy[0][ix], bsxy[1][ix],
541 TMath::Sqrt(esxy[0][0][ix]), TMath::Sqrt(esxy[1][1][ix]), corr,
542 esxy[0][0][ix], esxy[0][1][ix], esxy[1][1][ix]);
543 }
544
545 ix++;
546 }
547
548 }
549
550 //--------------------------------------
551
552 return kTRUE;
553}
554
555
556
557//---------------------------------------------------------------------------
558//
559// FindSkyCamTrans
560//
561// This routine determines the transformation
562//
563// ( cos(alfa) -sin(alfa) )
564// b = lambda * A * a + d A = ( )
565// ^ ^ ^ ( sin(alfa) cos(alfa) )
566// | | |
567// scale rotation shift
568// factor matrix
569//
570// from sky coordinates 'a' (projected onto the camera) to camera
571// coordinates 'b', using the positions of known stars in the camera.
572// The latter positions may have been determined by analysing the
573// DC currents in the different pixels.
574//
575// Input : a[2] x and y coordinates of stars projected onto the camera;
576// they were obtained from (RA, dec) of the stars and
577// (ThetaTel, PhiTel) and the time of observation;
578// these are the 'expected positions' of stars in the camera
579// b[2] 'measured positions' of these stars in the camera;
580// they may have been obtained from the DC currents
581// e[2][2] error matrix of b[2]
582// fixedrotationangle value [in degrees] at which rotation angle
583// alfa should be fixed; -1 means don't fix
584// fixedscalefactor value at which scale factor lambda
585// should be fixed; -1 means don't fix
586//
587// Output : lambda, alfadeg, A[2][2], d[2] fit results;
588// parameters describing the transformation
589// from 'expected positions' to the 'measured
590// positions' in the camera
591// errd[2][2] error matrix of d[2]
592// fNumIter number of iterations
593// fNdoF number of degrees of freedom
594// fChi2 chi-square value
595// fChi2Prob chi-square probability
596//
597// The units are assumed to be
598// [degrees] for alfadeg
599// [mm] for a, b, d
600// [1] for lambda
601
602Bool_t MTelAxisFromStars::FindSkyCamTrans(
603 TArrayD a[2], TArrayD b[2], TArrayD e[2][2],
604 Double_t &fixedrotationang, Double_t &fixedscalefac, Double_t &lambda,
605 Double_t &alfadeg, Double_t A[2][2], Double_t d[2], Double_t errd[2][2],
606 Int_t &fNumIter, Int_t &fNdof, Double_t &fChi2, Double_t &fChi2Prob)
607{
608 Int_t fNumStars = a[0].GetSize();
609
610 //*fLog << "MTelAxisFromStars::FindSkyCamTrans; expected and measured positions :"
611 // << endl;
612 for (Int_t ix=0; ix<fNumStars; ix++)
613 {
614 //*fLog << " ix, a[0], a[1], b[0], b[1], errxx, errxy, erryy = "
615 // << ix << " : "
616 // << a[0][ix] << ", " << a[1][ix] << "; "
617 // << b[0][ix] << ", " << b[1][ix] << "; "
618 // << e[0][0][ix] << ", " << e[0][1][ix] << ", "
619 // << e[1][1][ix] << endl;
620 }
621
622
623 //-------------------------------------------
624 // fix some parameters if the number of degrees of freedom is too low
625 // (<= 0.0)
626
627 Double_t fixedscalefactor = fixedscalefac;
628 Double_t fixedrotationangle = fixedrotationang;
629
630 // calculate number of degrees of freedom
631 fNdof = 2 * fNumStars - 4;
632 if (fixedscalefactor != -1.0)
633 fNdof += 1;
634 if (fixedrotationangle != -1.0)
635 fNdof += 1;
636
637 // if there is only 1 star fix both rotation angle and scale factor
638 if (fNumStars == 1)
639 {
640 if (fixedscalefactor == -1.0)
641 {
642 fixedscalefactor = 1.0;
643 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; scale factor is fixed at "
644 << fixedscalefactor << endl;
645 }
646 if (fixedrotationangle == -1.0)
647 {
648 fixedrotationangle = 0.0;
649 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; rotation angle is fixed at "
650 << fixedrotationangle << endl;
651 }
652 }
653 // otherwise fix only 1 parameter if possible
654 else if (fNdof < 0)
655 {
656 if (fNdof < -2)
657 {
658 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; number of degrees of freedom is too low : "
659 << fNdof << "; fNumStars = " << fNumStars << endl;
660 return kFALSE;
661 }
662 else if (fNdof == -2)
663 {
664 if (fixedscalefactor == -1.0 && fixedrotationangle == -1.0)
665 {
666 fixedscalefactor = 1.0;
667 fixedrotationangle = 0.0;
668 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; scale factor and rotation angle are fixed at "
669 << fixedscalefactor << " and " << fixedrotationangle
670 << " respectively" << endl;
671 }
672 else
673 {
674 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; number of degrees of freedom is too low : "
675 << fNdof << "; fNumStars = " << fNumStars << endl;
676 return kFALSE;
677 }
678 }
679 else if (fNdof == -1)
680 {
681 if (fixedrotationangle == -1.0)
682 {
683 fixedrotationangle = 0.0;
684 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; rotation angle is fixed at "
685 << fixedrotationangle << endl;
686 }
687 else if (fixedscalefactor == -1.0)
688 {
689 fixedscalefactor = 1.0;
690 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; scale factor is fixed at "
691 << fixedscalefactor << endl;
692 }
693 else
694 {
695 *fLog << warn << "MTelAxisFromStars::FindSkyCamTrans; number of degrees of freedom is too low : "
696 << fNdof << "; fNumStars = " << fNumStars<< endl;
697 return kFALSE;
698 }
699 }
700 }
701
702 // recalculate number of degrees of freedom
703 fNdof = 2 * fNumStars - 4;
704 if (fixedscalefactor != -1.0)
705 fNdof += 1;
706 if (fixedrotationangle != -1.0)
707 fNdof += 1;
708
709 if (fNdof < 0)
710 return kFALSE;
711 //-------------------------------------------
712
713
714 // get first approximation of scaling factor
715 if (fixedscalefactor != -1.0)
716 lambda = fixedscalefactor;
717 else
718 lambda = 1.0;
719
720 Double_t lambdaold = lambda;
721 Double_t dlambda = 0.0;
722
723 // get first approximation of rotation angle
724 Double_t alfa = 0.0;
725 if (fixedrotationangle != -1.0)
726 alfa = fixedrotationangle / TMath::RadToDeg();
727
728 Double_t alfaold = alfa;
729 // maximum allowed change of alfa in 1 iteration step (5 degrees)
730 Double_t dalfamax = 5.0 / TMath::RadToDeg();
731 Double_t dalfa = 0.0;
732
733 Double_t cosal = TMath::Cos(alfa);
734 Double_t sinal = TMath::Sin(alfa);
735
736 A[0][0] = cosal;
737 A[0][1] = -sinal;
738 A[1][0] = sinal;
739 A[1][1] = cosal;
740
741
742 Double_t absdold2 = 10000.0;
743 Double_t fChangeofd2 = 10000.0;
744
745
746 TArrayD Aa[2];
747 Aa[0].Set(fNumStars);
748 Aa[1].Set(fNumStars);
749
750
751 Double_t sumEbminlamAa[2];
752
753 TArrayD Ebminlambracd[2];
754 Ebminlambracd[0].Set(fNumStars);
755 Ebminlambracd[1].Set(fNumStars);
756
757 TArrayD EAa[2];
758 EAa[0].Set(fNumStars);
759 EAa[1].Set(fNumStars);
760
761 // invert the error matrices
762 TArrayD c[2][2];
763 c[0][0].Set(fNumStars);
764 c[0][1].Set(fNumStars);
765 c[1][0].Set(fNumStars);
766 c[1][1].Set(fNumStars);
767
768 for (Int_t ix=0; ix<fNumStars; ix++)
769 {
770 Double_t XX = e[0][0][ix];
771 Double_t XY = e[0][1][ix];
772 Double_t YY = e[1][1][ix];
773
774 // get inverse of error matrix
775 Double_t determ = XX*YY - XY*XY;
776 c[0][0][ix] = YY / determ;
777 c[0][1][ix] = -XY / determ;
778 c[1][0][ix] = -XY / determ;
779 c[1][1][ix] = XX / determ;
780 }
781
782
783
784 // calculate sum of inverted error matrices
785 Double_t determsumc;
786 Double_t sumc[2][2];
787 sumc[0][0] = 0.0;
788 sumc[0][1] = 0.0;
789 sumc[1][0] = 0.0;
790 sumc[1][1] = 0.0;
791
792 for (Int_t ix=0; ix<fNumStars; ix++)
793 {
794 sumc[0][0] += c[0][0][ix];
795 sumc[0][1] += c[0][1][ix];
796 sumc[1][0] += c[1][0][ix];
797 sumc[1][1] += c[1][1][ix];
798 }
799 determsumc = sumc[0][0]*sumc[1][1] - sumc[0][1]*sumc[1][0];
800
801 // calculate inverse of sum of inverted error matrices
802 Double_t sumcinv[2][2];
803 sumcinv[0][0] = sumc[1][1] / determsumc;
804 sumcinv[0][1] = -sumc[0][1] / determsumc;
805 sumcinv[1][0] = -sumc[1][0] / determsumc;
806 sumcinv[1][1] = sumc[0][0] / determsumc;
807
808 //*fLog << "sumcinv = " << sumcinv[0][0] << ", " << sumcinv[0][1]
809 // << ", " << sumcinv[1][1] << endl;
810
811
812 // minimize chi2 by iteration ***** start **********************
813
814 // stop iteration when change in |d|*|d| is less than 'told2'
815 // and change in alfa is less than 'toldalfa'
816 // and change in lambda is less than 'toldlambda'
817 // or chi2 is less than 'tolchi2'
818 Double_t told2 = 0.3*0.3; // [mm*mm]; 1/100 of an inner pixel diameter
819 Double_t toldalfa = 0.01 / TMath::RadToDeg(); // 0.01 degrees
820 Double_t toldlambda = 0.00006; // uncertainty of 1 mm of distance
821 // between camera and reflector
822 Double_t tolchi2 = 1.e-5;
823
824 Int_t fNumIterMax = 100;
825 fNumIter = 0;
826
827 for (Int_t i=0; i<fNumIterMax; i++)
828 {
829 fNumIter++;
830
831 // get next approximation of d ------------------
832 for (Int_t ix=0; ix<fNumStars; ix++)
833 {
834 Aa[0][ix] = A[0][0] * a[0][ix] + A[0][1]*a[1][ix];
835 Aa[1][ix] = A[1][0] * a[0][ix] + A[1][1]*a[1][ix];
836
837 //*fLog << "ix, Aa = " << ix << " : " << Aa[0][ix] << ", "
838 // << Aa[1][ix] << endl;
839 }
840
841 sumEbminlamAa[0] = 0.0;
842 sumEbminlamAa[1] = 0.0;
843
844 for (Int_t ix=0; ix<fNumStars; ix++)
845 {
846 sumEbminlamAa[0] += c[0][0][ix] * (b[0][ix] - lambda*Aa[0][ix])
847 + c[0][1][ix] * (b[1][ix] - lambda*Aa[1][ix]);
848
849 sumEbminlamAa[1] += c[1][0][ix] * (b[0][ix] - lambda*Aa[0][ix])
850 + c[1][1][ix] * (b[1][ix] - lambda*Aa[1][ix]);
851 }
852
853 //*fLog << "sumEbminlamAa = " << sumEbminlamAa[0] << ", "
854 // << sumEbminlamAa[1] << endl;
855
856 d[0] = sumcinv[0][0] * sumEbminlamAa[0]
857 + sumcinv[0][1] * sumEbminlamAa[1] ;
858
859 d[1] = sumcinv[1][0] * sumEbminlamAa[0]
860 + sumcinv[1][1] * sumEbminlamAa[1] ;
861
862 Double_t absdnew2 = d[0]*d[0] + d[1]*d[1];
863 fChangeofd2 = absdnew2 - absdold2;
864
865 //*fLog << "fNumIter : " << fNumIter
866 // << "; alfa, lambda, d[0], d[1], absdold2, absdnew2 = " << endl;
867 //*fLog << alfa << ", " << lambda << ", " << d[0] << ", " << d[1]
868 // << ", " << absdold2 << ", " << absdnew2 << endl;
869
870
871 if ( fabs(fChangeofd2) < told2 && fabs(dalfa) < toldalfa &&
872 fabs(dlambda) < toldlambda )
873 {
874 //*fLog << "Iteration stopped because of small changes : fChangeofd2, dalfa, dlambda = "
875 // << fChangeofd2 << ", " << dalfa << ", " << dlambda << endl;
876 break;
877 }
878 absdold2 = absdnew2;
879
880 // get next approximation of matrix A ----------------
881 if (fFixedRotationAngle == -1.0)
882 {
883 for (Int_t ix=0; ix<fNumStars; ix++)
884 {
885 Ebminlambracd[0][ix] =
886 c[0][0][ix] * ( b[0][ix] - lambda*Aa[0][ix] - d[0] )
887 + c[0][1][ix] * ( b[1][ix] - lambda*Aa[1][ix] - d[1] );
888
889 Ebminlambracd[1][ix] =
890 c[1][0][ix] * ( b[0][ix] - lambda*Aa[0][ix] - d[0] )
891 + c[1][1][ix] * ( b[1][ix] - lambda*Aa[1][ix] - d[1] );
892
893 //*fLog << "ix, Ebminlambracd = " << ix << " : "
894 // << Ebminlambracd[0][ix] << ", "
895 // << Ebminlambracd[1][ix] << endl;
896 }
897
898 // stop iteration if fChi2 is small enough
899 fChi2 = 0.0;
900 for (Int_t ix=0; ix<fNumStars; ix++)
901 {
902 fChi2 += (b[0][ix]-lambda*Aa[0][ix]-d[0] ) * Ebminlambracd[0][ix]
903 + (b[1][ix]-lambda*Aa[1][ix]-d[1] ) * Ebminlambracd[1][ix];
904 }
905 if ( fChi2 < tolchi2 )
906 {
907 //*fLog << "iteration stopped because of small fChi2 : "
908 // << fChi2 << endl;
909 break;
910 }
911
912
913 Double_t dchi2dA[2][2];
914 dchi2dA[0][0] = 0.0;
915 dchi2dA[0][1] = 0.0;
916 dchi2dA[1][0] = 0.0;
917 dchi2dA[1][1] = 0.0;
918
919 for (Int_t ix=0; ix<fNumStars; ix++)
920 {
921 dchi2dA[0][0] += Ebminlambracd[0][ix] * a[0][ix];
922 dchi2dA[0][1] += Ebminlambracd[0][ix] * a[1][ix];
923 dchi2dA[1][0] += Ebminlambracd[1][ix] * a[0][ix];
924 dchi2dA[1][1] += Ebminlambracd[1][ix] * a[1][ix];
925 }
926
927 //*fLog << "dchi2dA = " << dchi2dA[0][0] << ", " << dchi2dA[0][1]
928 // << ", " << dchi2dA[1][0] << ", " << dchi2dA[1][1] << endl;
929
930 // ********* 1st derivative (d chi2) / (d alfa) ************
931 Double_t dchi2dalfa = -2.0*lambda *
932 ( - sinal*(dchi2dA[0][0]+dchi2dA[1][1])
933 + cosal*(dchi2dA[1][0]-dchi2dA[0][1]) );
934
935
936 //Double_t dalfa1st = - fChi2 / dchi2dalfa;
937
938 //*fLog << "fChi2, dchi2dalfa = " << fChi2 << ", "
939 // << dchi2dalfa << endl;
940 //*fLog << "proposed change of alfa using 1st derivative = "
941 // << dalfa1st << endl;
942
943 // ********* 2nd derivative (d2 chi2) / (d alfa2) ******
944 Double_t term1 = 0.0;
945 Double_t term2 = 0.0;
946 Double_t term3 = 0.0;
947 Double_t term4 = 0.0;
948
949 for (Int_t ix=0; ix<fNumStars; ix++)
950 {
951 term1 += a[0][ix]*c[0][0][ix]*a[0][ix] + a[1][ix]*c[1][0][ix]*a[0][ix]
952 + a[0][ix]*c[0][1][ix]*a[1][ix] + a[1][ix]*c[1][1][ix]*a[1][ix];
953
954 term2 += a[0][ix]*c[1][0][ix]*a[0][ix] - a[1][ix]*c[0][0][ix]*a[0][ix]
955 + a[0][ix]*c[1][1][ix]*a[1][ix] - a[1][ix]*c[0][1][ix]*a[1][ix];
956
957 term3 = a[0][ix]*c[0][0][ix]*a[1][ix] + a[1][ix]*c[1][0][ix]*a[1][ix]
958 - a[0][ix]*c[0][1][ix]*a[0][ix] - a[1][ix]*c[1][1][ix]*a[0][ix];
959
960 term4 += a[0][ix]*c[1][0][ix]*a[1][ix] - a[1][ix]*c[0][0][ix]*a[1][ix]
961 - a[0][ix]*c[1][1][ix]*a[0][ix] + a[1][ix]*c[0][1][ix]*a[0][ix];
962 }
963
964 Double_t d2chi2dalfa2 =
965 - 2.0*lambda * ( - cosal*(dchi2dA[0][0]+dchi2dA[1][1])
966 - sinal*(dchi2dA[1][0]-dchi2dA[0][1]) )
967 + 2.0*lambda*lambda * ( sinal*sinal * term1 - sinal*cosal * term2
968 + sinal*cosal * term3 - cosal*cosal * term4);
969
970 // Gauss-Newton step
971 Double_t dalfa2nd = - dchi2dalfa / d2chi2dalfa2;
972
973 //*fLog << "proposed change of alfa using 2st derivative = "
974 // << dalfa2nd << endl;
975
976 //dalfa = dalfa1st;
977 dalfa = dalfa2nd;
978
979 // ******************************************
980
981
982 // restrict change of alfa
983 if ( fabs(dalfa) > dalfamax )
984 {
985 dalfa = TMath::Sign( dalfamax, dalfa );
986 }
987 alfa = alfaold + dalfa;
988
989 if ( alfa < -5.0/TMath::RadToDeg() )
990 alfa = -5.0/TMath::RadToDeg();
991 else if ( alfa > 5.0/TMath::RadToDeg() )
992 alfa = 5.0/TMath::RadToDeg();
993
994 dalfa = alfa - alfaold;
995
996 alfaold = alfa;
997
998 sinal = TMath::Sin(alfa);
999 cosal = TMath::Cos(alfa);
1000
1001 A[0][0] = cosal;
1002 A[0][1] = -sinal;
1003 A[1][0] = sinal;
1004 A[1][1] = cosal;
1005
1006 //*fLog << "alfa-alfaold = " << dalfa << endl;
1007 //*fLog << "new alfa = " << alfa << endl;
1008 }
1009
1010
1011 // get next approximation of lambda ----------------
1012 if (fFixedScaleFactor == -1.0)
1013 {
1014 for (Int_t ix=0; ix<fNumStars; ix++)
1015 {
1016 Aa[0][ix] = A[0][0]*a[0][ix] + A[0][1]*a[1][ix];
1017 Aa[1][ix] = A[1][0]*a[0][ix] + A[1][1]*a[1][ix];
1018
1019 EAa[0][ix] =
1020 c[0][0][ix] * Aa[0][ix] + c[0][1][ix] * Aa[1][ix];
1021 EAa[1][ix] =
1022 c[1][0][ix] * Aa[0][ix] + c[1][1][ix] * Aa[1][ix];
1023
1024 //*fLog << "ix, Aa = " << ix << " : " << Aa[0][ix] << ", "
1025 // << Aa[1][ix] << endl;
1026
1027 //*fLog << "ix, EAa = " << ix << " : " << EAa[0][ix] << ", "
1028 // << EAa[1][ix] << endl;
1029 }
1030
1031 Double_t num = 0.0;
1032 Double_t denom = 0.0;
1033 for (Int_t ix=0; ix<fNumStars; ix++)
1034 {
1035 num += (b[0][ix]-d[0]) * EAa[0][ix]
1036 + (b[1][ix]-d[1]) * EAa[1][ix];
1037
1038 denom += Aa[0][ix] * EAa[0][ix]
1039 + Aa[1][ix] * EAa[1][ix];
1040
1041 //*fLog << "ix : b-d = " << ix << " : " << b[0][ix]-d[0]
1042 // << ", " << b[1][ix]-d[1] << endl;
1043
1044 //*fLog << "ix : Aa = " << ix << " : " << Aa[0][ix]
1045 // << ", " << Aa[1][ix] << endl;
1046 }
1047
1048 lambda = num / denom;
1049
1050 if ( lambda < 0.9 )
1051 lambda = 0.9;
1052 else if ( lambda > 1.1 )
1053 lambda = 1.1;
1054
1055 dlambda = lambda - lambdaold;
1056 lambdaold = lambda;
1057
1058 //*fLog << "num, denom, lambda, dlambda = " << num
1059 // << ", " << denom << ", " << lambda << ", "
1060 // << dlambda << endl;
1061 }
1062
1063 }
1064 //------- end of iteration *****************************************
1065
1066 alfadeg = alfa * TMath::RadToDeg();
1067
1068 // calculate error matrix of d[2]
1069 errd[0][0] = sumcinv[0][0];
1070 errd[0][1] = sumcinv[0][1];
1071 errd[1][0] = sumcinv[1][0];
1072 errd[1][1] = sumcinv[1][1];
1073
1074 // evaluate quality of fit
1075
1076 // calculate chi2
1077 for (Int_t ix=0; ix<fNumStars; ix++)
1078 {
1079 Ebminlambracd[0][ix] =
1080 c[0][0][ix] * ( b[0][ix] - lambda*Aa[0][ix] - d[0] )
1081 + c[0][1][ix] * ( b[1][ix] - lambda*Aa[1][ix] - d[1] );
1082
1083 Ebminlambracd[1][ix] =
1084 c[1][0][ix] * (b[0][ix] - lambda*Aa[0][ix] - d[0] )
1085 + c[1][1][ix] * (b[1][ix] - lambda*Aa[1][ix] - d[1] );
1086 }
1087
1088 fChi2 = 0.0;
1089 for (Int_t ix=0; ix<fNumStars; ix++)
1090 {
1091 fChi2 += (b[0][ix] - lambda*Aa[0][ix] - d[0] ) * Ebminlambracd[0][ix]
1092 + (b[1][ix] - lambda*Aa[1][ix] - d[1] ) * Ebminlambracd[1][ix];
1093 }
1094
1095 fChi2Prob = TMath::Prob(fChi2, fNdof);
1096
1097 //*fLog << "MTelAxisFromStars::FindSkyCamTrans :" << endl;
1098 //*fLog << " fNumStars, fChi2, fNdof, fChi2Prob, fNumIter, fChangeofd2, dalfa, dlambda = "
1099 // << fNumStars << ", " << fChi2 << ", " << fNdof << ", "
1100 // << fChi2Prob << ", "
1101 // << fNumIter << ", " << fChangeofd2 << ", " << dalfa << ", "
1102 // << dlambda << endl;
1103 //*fLog << " lambda, alfadeg, d[0], d[1] = " << lambda << ", "
1104 // << alfadeg << ", " << d[0] << ", " << d[1] << endl;
1105
1106 return kTRUE;
1107}
1108
1109// --------------------------------------------------------------------------
1110//
1111// Apply transformation (lambda, A, d)
1112// to the expected positions (a[1], a[2])
1113// to obtain the estimated positions (b[1], b[2])
1114//
1115// e[2][2] is the error matrix of b[2]
1116
1117void MTelAxisFromStars::TransSkyCam(
1118 Double_t &lambda, Double_t A[2][2], Double_t d[2], Double_t errd[2][2],
1119 TArrayD a[2], TArrayD b[2], TArrayD e[2][2])
1120{
1121 Int_t numpos = a[0].GetSize();
1122 if (numpos <= 0)
1123 return;
1124
1125 //*fLog << "MTelAxisFromStars::TransSkyCam; lambda, A, d = "
1126 // << lambda << "; " << A[0][0] << ", " << A[0][1] << ", "
1127 // << A[1][1] << "; " << d[0] << ", " << d[1] << endl;
1128
1129 //*fLog << "MTelAxisFromStars::TransSkyCam; expected and estimated positions :"
1130 // << endl;
1131
1132 for (Int_t ix=0; ix<numpos; ix++)
1133 {
1134 // *fLog << "MTelAxisFromStars; ix = " << ix << endl;
1135
1136 b[0][ix] = lambda * (A[0][0]*a[0][ix] + A[0][1]*a[1][ix]) + d[0];
1137 b[1][ix] = lambda * (A[1][0]*a[0][ix] + A[1][1]*a[1][ix]) + d[1];
1138
1139 e[0][0][ix] = errd[0][0];
1140 e[0][1][ix] = errd[0][1];
1141 e[1][0][ix] = errd[1][0];
1142 e[1][1][ix] = errd[1][1];
1143
1144 // *fLog << " ix, a[0], a[1], b[0], b[1], errxx, errxy, erryy = "
1145 // << ix << " : "
1146 // << a[0][ix] << ", " << a[1][ix] << "; "
1147 // << b[0][ix] << ", " << b[1][ix] << "; "
1148 // << e[0][0][ix] << ", " << e[0][1][ix] << ", "
1149 // << e[1][1][ix] << endl;
1150 }
1151}
1152
1153// --------------------------------------------------------------------------
1154//
1155// SetPointingPosition
1156//
1157// put the corrected pointing direction into MPointingPos[MPointingPos];
1158// this direction corresponds to the position (0,0) in the camera
1159//
1160
1161Bool_t MTelAxisFromStars::SetPointingPosition(MStarCamTrans *fstarcamtrans,
1162 MReportDrive *fdrive, MSkyCamTrans *ftrans, MPointingPos *fpointpos)
1163{
1164 Double_t decdrive = fdrive->GetDec();
1165 Double_t hdrive = fdrive->GetHa();
1166 Double_t radrive = fdrive->GetRa();
1167
1168 // this is the estimated position (with optical aberration) in the camera
1169 // corresponding to the direction in MReportDrive
1170 Double_t Xpoint = (ftrans->GetShiftD())[0];
1171 Double_t Ypoint = (ftrans->GetShiftD())[1];
1172
1173 // get the sky direction corresponding to the position (0,0) in the camera
1174 Double_t decpoint = 0.0;
1175 Double_t hpoint = 0.0;
1176 fstarcamtrans->CelCamToCel0(decdrive, hdrive,
1177 Xpoint/fAberr, Ypoint/fAberr, decpoint, hpoint);
1178 Double_t rapoint = radrive - hpoint + hdrive;
1179 fpointpos->SetSkyPosition(rapoint, decpoint, hpoint);
1180
1181 // get the local direction corresponding to the position (0,0) in the camera
1182 Double_t thetadrive = fdrive->GetNominalZd();
1183 Double_t phidrive = fdrive->GetNominalAz();
1184 Double_t thetapoint = 0.0;
1185 Double_t phipoint = 0.0;
1186 fstarcamtrans->LocCamToLoc0(thetadrive, phidrive,
1187 Xpoint/fAberr, Ypoint/fAberr, thetapoint, phipoint);
1188 fpointpos->SetLocalPosition(thetapoint, phipoint);
1189 fpointpos->SetReadyToSave();
1190
1191 //*fLog << "SetPointingPosition : decdrive, hdrive, radrive Xpoint, Ypoint = "
1192 // << decdrive << ", " << hdrive << ", " << radrive << ", "
1193 // << Xpoint << ", " << Ypoint << endl;
1194
1195 //*fLog << "SetPointingPosition : thetadrive, phidrive, thetapoint, phipoint = "
1196 // << thetadrive << ", " << phidrive << ", " << thetapoint << ", "
1197 // << phipoint << endl;
1198
1199 return kTRUE;
1200}
1201
1202// --------------------------------------------------------------------------
1203//
1204// SetSourcePosition
1205//
1206// put the estimated position of the source in the camera into
1207// MSrcPosCam[MSrcPosCam]
1208//
1209// and the estimated local direction of the source into
1210// MSourcePos[MPointingPos]
1211//
1212
1213Bool_t MTelAxisFromStars::SetSourcePosition(MStarCamTrans *fstarcamtrans,
1214 MPointingPos *fpointpos, MPointingPos *fsourcepos, MSrcPosCam *fsrcpos)
1215{
1216 // get the corrected pointing direction
1217 // corresponding to the position (0,0) in the camera
1218 Double_t decpoint = fpointpos->GetDec();
1219 Double_t hpoint = fpointpos->GetHa();
1220
1221 // get the sky direction of the source
1222 Double_t decsource = fsourcepos->GetDec();
1223 Double_t hsource = fsourcepos->GetHa();
1224
1225 // get the estimated position (Xsource, Ysource) of the source in the camera;
1226 // this is a position for an ideal imaging, without optical aberration
1227 Double_t Xsource = 0.0;
1228 Double_t Ysource = 0.0;
1229 fstarcamtrans->Cel0CelToCam(decpoint, hpoint,
1230 decsource, hsource, Xsource, Ysource);
1231 fsrcpos->SetXY(Xsource*fAberr, Ysource*fAberr);
1232 fsrcpos->SetReadyToSave();
1233
1234 // get the estimated local direction of the source
1235 Double_t thetapoint = fpointpos->GetZd();
1236 Double_t phipoint = fpointpos->GetAz();
1237 Double_t thetasource = 0.0;
1238 Double_t phisource = 0.0;
1239 fstarcamtrans->Loc0CamToLoc(thetapoint, phipoint,
1240 Xsource, Ysource, thetasource, phisource);
1241 fsourcepos->SetLocalPosition(thetasource, phisource);
1242 fsourcepos->SetReadyToSave();
1243
1244 //*fLog << "SetSourcePosition : decpoint, hpoint, decsource, hsource, Xsource, Ysource = "
1245 // << decpoint << ", " << hpoint << ", " << decsource << ", "
1246 // << hsource << ", " << Xsource << ", " << Ysource << endl;
1247 //*fLog << "SetSourcePosition : thetapoint, phipoint, thetasource, phisource = "
1248 // << thetapoint << ", " << phipoint << ", " << thetasource << ", "
1249 // << phisource << endl;
1250
1251 return kTRUE;
1252}
1253
1254// --------------------------------------------------------------------------
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
Note: See TracBrowser for help on using the repository browser.