trkfit.cc 50.7 KB
Newer Older
1 2 3 4 5 6
/*  PAON4 analysis software 
    classes and functions to read in and perform array geometry determination 
    using satellites and celestial sources tracks  
    R. Ansari, Fevrier 2019                                             */


7 8
#include <iomanip>

9 10 11
#include "pexceptions.h"
#include "trkfit.h"
#include "datacards.h"
12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28
#include "array.h"

#include "acbeam.h"
#include "gacfit.h"
#include "gcxfit.h"
#include "gcxfitbaseline.h"

#include "p4autils.h"


//------------------- Print Level for this file --------------------------
static int _prtlevel_ =0;
void TrkFit_SetPrintLevel(int lev) 
{ 
  _prtlevel_=lev; 
  return;
}
29

30 31 32 33 34 35 36 37 38 39 40 41
void TrkFit_FitLibInfo() 
{
  cout << "============================================================================"<<endl;
#ifndef TKF_AVEC_MINUIT
  cout << "============ Classe TkF_Fitter : Fitting with Sophya GeneralFit ============"<<endl;
#else 
  cout << "============= Classe TkF_Fitter : Fitting with Minuit MnMigrad ============="<<endl;
#endif
  cout << "============================================================================"<<endl;
  return;
}

42 43
//------------------- TrkInputDataSet -------------------------------------

44 45

TrkInputDataSet::TrkInputDataSet(string dcfilename, string inp_path)
46 47
  : zenang(0.) , theta_0(0.) , phi_0(0.)
{
48
  setInputBasePath(inp_path);
49 50 51 52 53 54 55 56 57
  ReadDatacardFile(dcfilename);
}


static vector<string> * dataflnm_p_ = NULL;
static vector<double> * tstart_p_ = NULL;
static vector<double> * tend_p_ = NULL;
static vector<double> * v_freqs_p_ = NULL;
static vector<string> * trkflnm_p_ = NULL;
58 59
static vector<bool> * v_noAC_p_ = NULL;
static vector<bool> * v_noCx_p_ = NULL;
60 61 62 63 64 65 66 67 68 69 70 71 72
static size_t trk_cnt = 0;

static int decode_trkcard(string const& key, string const& toks)
{
  if (key != "trk") {  // CA NE DEVRAIT PAS ARRIVER 
    cout << "decode_trkcard/ERROR  BAD key = " << key << " ( <> trk"<<endl;
    return 1;
  }
  if (! dataflnm_p_ ) { // CA NE DEVRAIT PAS ARRIVER
    cout << "decode_trkcard/ERROR  dataflnm_p_ = NULL !"<<endl;
    return 1;
  }
  char flnmdata[256], flnmtrk[256];
73
  char sflags[64];
74
  double ts,te,freq;
75
  sscanf(toks.c_str(),"%s %lg,%lg %lg %s %s",flnmdata,&ts,&te,&freq,flnmtrk,sflags);
76 77 78 79 80 81

  dataflnm_p_->push_back(flnmdata);
  tstart_p_->push_back(ts*60.);
  tend_p_->push_back(te*60.);
  v_freqs_p_->push_back(freq);
  trkflnm_p_->push_back(flnmtrk);
82 83 84 85 86 87 88 89 90 91 92
  size_t ll=strlen(sflags);
  bool noAC=false;
  bool noCx=false;
  if (ll>0) {
    for(size_t l=0; l<ll; l++)  sflags[l]=toupper(sflags[l]);
    string sflg=sflags;
    if ((sflg == "NOAC")||(sflg=="NOACCX"))  noAC=true;
    if ((sflg == "NOCX")||(sflg=="NOACCX"))  noCx=true;
  }
  v_noAC_p_->push_back(noAC);
  v_noCx_p_->push_back(noCx);
93 94 95 96 97
  trk_cnt++;
  return 0;
}


98 99 100 101 102 103
void TrkInputDataSet::setInputBasePath(string inp_path)
{
  if (inp_path.length()>0)  input_base_path=inp_path;
  return;
}

104 105 106
size_t TrkInputDataSet::ReadDatacardFile(string dcfilename)
{
  DataCards dc;
107
  string match="trk";
108 109 110 111 112 113 114 115 116 117 118 119 120
  dc.AddProcF(decode_trkcard, match);

  zenang=0.; theta_0=0.;  phi_0=0.;
  dataflnm.clear();
  tstart.clear();
  tend.clear();
  v_freqs.clear();
  trkflnm.clear();
  dataflnm_p_ = &dataflnm;
  tstart_p_ = &tstart;
  tend_p_ = &tend;
  v_freqs_p_ = &v_freqs;
  trkflnm_p_ = &trkflnm;
121 122
  v_noAC_p_=&v_noAC;
  v_noCx_p_=&v_noCx;
123
  trk_cnt = 0;
124
  // @trk visiDataTableFile tstart,tend freq TrackFileName [FLAG]
125
  //  tstart , tend in minutes freq in MHz
126 127 128
  //  optional FLAG   = NOAC  NOCX   NOACCX   
  //  NOAC : don't use for Auto-correlation fit ;  NOCX : don't use for cross-cor fits 
  //  NOACCX : don't use for Auto-correlation or cross-cor fits 
129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169
  dc.ReadFile(dcfilename);
  if (dc.HasKey("inpath"))   {   // @inpath  InputFilesDirectoryPath    
    input_base_path = dc.SParam("inpath",0,"");
  }
  if (dc.HasKey("zenang"))   {   // @zenang  Zenith Angle in degree   
    zenang = dc.DParam("zenang",0,0.);
    if (zenang<0.) {
      theta_0 = Angle(-zenang, Angle::Degree).ToRadian();  phi_0 = Angle::PioTwoCst()+Angle::OnePiCst();
    }
    else {
      theta_0 = Angle(+zenang, Angle::Degree).ToRadian();  phi_0 = Angle::PioTwoCst();
    }
  }

  dataflnm_p_ = NULL;
  tstart_p_ = NULL;
  tend_p_ = NULL;
  v_freqs_p_ = NULL;
  trkflnm_p_ = NULL;

  if (trk_cnt != trkflnm.size()) {  // ca ne devrait pas arriver
    cout << " TrkInputDataSet::ReadDatacardFile()/BUG  trk_cnt != trkflnm.size()"<<endl;
    throw PError("TrkInputDataSet::ReadDatacardFile() trk_cnt != trkflnm.size()");
  }
  trk_cnt=0;
  dcfilename_ = dcfilename;
  return trkflnm.size();
}

ostream & TrkInputDataSet::Print(ostream & os) const
{
  os << "TrkInputDataSet(dcfilename="<<dcfilename_<<")/Info:  dec-shift(zenithAngle)= "<<zenang<<" NbTrk="<<NbTrk()<<endl;
  os << "...InputBaseDirectoryPath="<<input_base_path<<endl;
  for(size_t i=0; i<NbTrk(); i++)  {
    os <<"["<<i<<"] data= "<< dataflnm[i]<<"  ts,te(min)= "<<tstart[i]/60.<<","<<tend[i]/60.<<" freg(MHz)= "<<v_freqs[i]
       <<" TrkFile="<<trkflnm[i]<<endl;
  }
  return os;
}


170 171 172 173 174 175 176 177 178 179 180 181

//------------------------ ACxDataSet -------------------------------------

AcxDataSet::AcxDataSet(TrkInputDataSet & tkds)
  : tot_npoints(0),zenang(0.),theta_0(0.),phi_0(0.)
{
  ReadData(tkds);
}

AcxDataSet::AcxDataSet(AcxDataSet const & a)
  : v_time_data(a.v_time_data), vv_data(a.vv_data), vv_err(a.vv_err), 
    v_min_data(a.v_min_data), v_max_data(a.v_max_data),
182 183
    vv_cxdata(a.vv_cxdata), vv_cxerr(a.vv_cxerr),
    v_min_cxdata(a.v_min_cxdata), v_max_cxdata(a.v_max_cxdata), 
184
    tot_npoints(a.tot_npoints), v_freqs(a.v_freqs), v_noAC(a.v_noAC), v_noCx(a.v_noCx), 
185 186
    zenang(a.zenang), theta_0(a.theta_0), phi_0(a.phi_0),
    v_acbeams(a.v_acbeams), v_cxbeams(a.v_cxbeams),
187 188 189 190 191 192 193
    v_RcFit_ac(a.v_RcFit_ac), v_xi2red_ac(a.v_xi2red_ac),
    v_Ddish(a.v_Ddish), v_thetaant(a.v_thetaant), v_phiant(a.v_phiant),
    v_err_Ddish(a.v_err_Ddish), v_err_thetaant(a.v_err_thetaant), v_err_phiant(a.v_err_phiant),
    v_RcFit_cx(a.v_RcFit_cx), v_xi2red_cx(a.v_xi2red_cx),
    v_phase(a.v_phase), v_phi_0(a.v_phi_0), v_a_phi(a.v_a_phi),
    v_err_phi_0(a.v_err_phi_0), v_err_a_phi(a.v_err_a_phi),
    v_Acx(a.v_Acx), v_Bcx(a.v_Bcx)
194 195 196 197 198 199 200
{
}

AcxDataSet & AcxDataSet::operator = (AcxDataSet const & a)
{
  v_time_data=a.v_time_data; vv_data=a.vv_data; vv_err=a.vv_err; 
  v_min_data=a.v_min_data;   v_max_data=a.v_max_data;
201 202
  vv_cxdata=a.vv_cxdata;   vv_cxerr=a.vv_cxerr;
  v_min_cxdata=a.v_min_cxdata;  v_max_cxdata=a.v_max_cxdata; 
203
  tot_npoints=a.tot_npoints; v_freqs=a.v_freqs;   v_noAC=a.v_noAC;  v_noCx=a.v_noCx;
204 205
  zenang=a.zenang;  theta_0=a.theta_0;  phi_0=a.phi_0;
  v_acbeams=a.v_acbeams;  v_cxbeams=a.v_cxbeams;
206 207
  v_RcFit_ac=a.v_RcFit_ac;  v_xi2red_ac=a.v_xi2red_ac;
  v_Ddish=a.v_Ddish;  v_thetaant=a.v_thetaant;  v_phiant=a.v_phiant;
208
  v_err_Ddish=a.v_err_Ddish;  v_err_thetaant=a.v_err_thetaant;  v_err_phiant=a.v_err_phiant;
209 210 211 212
  v_RcFit_cx=a.v_RcFit_cx; v_xi2red_cx=a.v_xi2red_cx;
  v_phase=a.v_phase; v_phi_0=a.v_phi_0;  v_a_phi=a.v_a_phi;
  v_err_phi_0=a.v_err_phi_0; v_err_a_phi=a.v_err_a_phi; 
  v_Acx=a.v_Acx;  v_Bcx=a.v_Bcx;
213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231
  return (*this);
}

size_t AcxDataSet::ReadData(TrkInputDataSet & tkds)    
{
  cout << "---- AcxDataSet::AcxDataSet() reading 4 PAON4 auto-correlation & 6 Cross-cor signals/DataTables for"
       <<tkds.NbTrk()<<" tracks ..."<<endl;

  if (tkds.NbTrk() != v_time_data.size()) {
    v_time_data.resize(tkds.NbTrk());
    vv_data.resize(tkds.NbTrk());
    vv_err.resize(tkds.NbTrk());
    v_min_data.resize(tkds.NbTrk());
    v_max_data.resize(tkds.NbTrk());
    vv_cxdata.resize(tkds.NbTrk());
    vv_cxerr.resize(tkds.NbTrk());
    v_min_cxdata.resize(tkds.NbTrk());
    v_max_cxdata.resize(tkds.NbTrk());    
  }
232
  v_freqs=tkds.v_freqs;  v_noAC=tkds.v_noAC;  v_noCx=tkds.v_noCx;
233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310
  zenang=tkds.zenang;   theta_0=tkds.theta_0;    phi_0=tkds.phi_0;
  size_t NB_ANTENNES=getNbAutoCor();   // nombre d'antennes 
  size_t NB_CXCORS=getNbCrossCor();
  tot_npoints = 0;   // total number of points for fit 
  const char * acname[4]={"V11","V22","V33","V44"};
  const char * cxname[6]={"V12","V13","V14","V23","V24","V34"};
  
  for(size_t j=0; j<tkds.dataflnm.size(); j++) {
    string flnm = tkds.input_base_path+tkds.dataflnm[j]+".ppf";
    cout << "1."<<j+1<<" Extracting data from data file DataTable: " << flnm<<endl
	 << " ... For time interval (Trk"<<j+1<<") "<<tkds.tstart[j]<<" < t < "<<tkds.tend[j]<<endl;
    DataTable dt_data;
    PInPersist pin(flnm);
    pin >> dt_data;
    dt_data.SetShowMinMaxFlag(true);
    size_t ktime = dt_data.IndexNom("timesec");
    vector<double> vtm;
    dt_data.GetColumn(ktime, vtm);
    vector< vector<double> > v_vac(NB_ANTENNES);
    for(size_t ii=0; ii<NB_ANTENNES; ii++) {   // 4 auto-correlations
      size_t kac = dt_data.IndexNom(acname[ii]);
      dt_data.GetColumn(kac, v_vac[ii]);
      vector<double> vtmp, vetmp;
      vv_data[j].push_back(vtmp);
      vv_err[j].push_back(vetmp);
      v_min_data[j].push_back(9.e19);
      v_max_data[j].push_back(-9.e19);
    }
    vector< vector <complex<double> > > v_vcx(NB_CXCORS);
    for(size_t ii=0; ii<NB_CXCORS; ii++) {   // 6 cross-correlations
      size_t kac = dt_data.IndexNom(cxname[ii]);
      dt_data.GetColumn(kac, v_vcx[ii]);
      vector< complex<double> > vtmp;
      vector<double> vetmp;
      vv_cxdata[j].push_back(vtmp);
      vv_cxerr[j].push_back(vetmp);
      v_min_cxdata[j].push_back(9.e19);
      v_max_cxdata[j].push_back(-9.e19);
    }
    
    vector< vector<double> > & v_data = vv_data[j];
    vector< vector<double> > & v_err = vv_err[j];
    vector< vector< complex<double> > > & v_cxdata = vv_cxdata[j];
    vector< vector<double> > & v_cxerr = vv_cxerr[j];

    for(size_t k=0; k<vtm.size(); k++) {
      if ((vtm[k]<tkds.tstart[j])||(vtm[k]>tkds.tend[j]))  continue;
      v_time_data[j].push_back(vtm[k]);
      for(size_t ii=0; ii<NB_ANTENNES; ii++) {
	vector<double> & vac = v_vac[ii];
	v_data[ii].push_back(vac[k]);
	v_err[ii].push_back(0.1*sqrt(fabs(vac[k])));   // calcul d'erreur, a affiner 
	if (vac[k]<v_min_data[j][ii])  v_min_data[j][ii]=vac[k];
	if (vac[k]>v_max_data[j][ii])  v_max_data[j][ii]=vac[k];
      }
      for(size_t ii=0; ii<NB_CXCORS; ii++) {   // 6 cross-correlations
	vector< complex<double> > & vcx = v_vcx[ii];
	v_cxdata[ii].push_back(vcx[k]);
	double acx=std::abs(vcx[k]);
	v_cxerr[ii].push_back(0.1*sqrt(acx));
	if (acx<v_min_cxdata[j][ii])  v_min_cxdata[j][ii]=acx;
	if (acx>v_max_cxdata[j][ii])  v_max_cxdata[j][ii]=acx;
      }
    }
    
    tot_npoints += v_time_data[j].size();   // total number of points for fit 
    cout << " ... Done for " << j+1 << " data size="<<v_time_data[j].size()<<endl;
    cout << "  Data-AutoCor Min,Max[A1...A4]="; 
    for(size_t ii=0; ii<NB_ANTENNES; ii++)
      cout<<setw(10)<<v_min_data[j][ii]<<","<<setw(10)<<v_max_data[j][ii]<<" ; ";   cout << endl;
    cout << "  Data-CxCorr (abs) Min,Max[Cx1...Cx6]="; 
    for(size_t ii=0; ii<NB_ANTENNES; ii++)
      cout<<setw(10)<<v_min_cxdata[j][ii]<<","<<setw(10)<<v_max_cxdata[j][ii]<<" ; ";   cout << endl;

  }
  return tot_npoints;
}

311 312 313 314 315 316 317 318 319 320 321 322
ostream & AcxDataSet::PrintACFitSummary(ostream & os)
{
  const char* acnames[4]={"AC-1","AC-2","AC-3","AC-4"}; 
  os << "--------- Fitted Parameters and errors from AutoCorrelations (D-dish, Theta,Phi Antennes) "<<endl;
  for(size_t i=0; i<getNbAutoCor(); i++) {
    double thetaant=v_thetaant[i];  double err_thetaant=v_err_thetaant[i];
    double phiant=v_phiant[i];  double err_phiant=v_err_phiant[i];
    double elevdeg=90.-Angle(thetaant).ToDegree();
    double err_elevdeg=Angle(err_thetaant).ToDegree();
    double azimdeg=90.-Angle(phiant).ToDegree();
    if (azimdeg<0.)  azimdeg += 360.;
    double err_azimdeg=Angle(err_phiant).ToDegree();
323
    os<<acnames[i]<<" D(m)= "<<setw(8)<<v_Ddish[i]<<" +/- "<<setw(8)<<v_err_Ddish[i];
324 325
    os<<" Elev(deg)= "<<setw(8)<<elevdeg<<" +/- "<<setw(8)<<err_elevdeg;
    os<<" Azim(deg)= "<<setw(8)<<azimdeg<<" +/- "<<setw(8)<<err_azimdeg;
326
    os<<"  RcFit="<<setw(6)<<v_RcFit_ac[i]<<" Xi2Red="<<setw(8)<<v_xi2red_ac[i]<<endl;
327 328 329 330
  }
  return os;
}

331 332 333
ostream & AcxDataSet::PrintCxPhaseFitSummary(ostream & os)
{
  const char* cxnames[6]={"Cx-1x2","Cx-1x3","Cx-1x4","Cx-2x3","Cx-2x4","Cx-3x4"};
334
  os << "--------- Cx-Fitted phases @1300 MHz ";
335 336 337 338 339 340 341 342
  for(size_t i=0; i<getNbCrossCor(); i++) os<<setw(8)<<Angle(v_phase[i]).ToDegree()<<" ; ";   cout<<endl;
  double dphi23=v_phase[1]-v_phase[0];   if (dphi23<0.) dphi23+=(2.*M_PI);
  double dphi24=v_phase[2]-v_phase[0];   if (dphi24<0.) dphi24+=(2.*M_PI);
  double dphi34=v_phase[2]-v_phase[1];   if (dphi34<0.) dphi34+=(2.*M_PI);
  os << "--- Compatibility of fitted phases (@1300 MHz) among the six Baselines "<<endl;
  os<<" Cx-2x3: Phi3-Phi2= "<<setw(6)<<Angle(dphi23).ToDegree()<<" EqualTo? Phi23= "<<setw(6)<<Angle(v_phase[3]).ToDegree()<<endl;
  os<<" Cx-2x4: Phi4-Phi2= "<<setw(6)<<Angle(dphi24).ToDegree()<<" EqualTo? Phi24= "<<setw(6)<<Angle(v_phase[4]).ToDegree()<<endl;
  os<<" Cx-3x4: Phi4-Phi3= "<<setw(6)<<Angle(dphi34).ToDegree()<<" EqualTo? Phi34= "<<setw(6)<<Angle(v_phase[5]).ToDegree()<<endl;
343
  os << "--------- Cx-Fitted Phase(freq) parameters and errors  Phi(freq)=phi0+a_phi*(freq-1250.)/250. "<<endl;
344 345
  for(size_t i=0; i<getNbCrossCor(); i++) {
    os<<cxnames[i]<<" phi0= "<<setw(8)<<Angle(v_phi_0[i]).ToDegree()<<" +/- "<<setw(8)<<Angle(v_err_phi_0[i]).ToDegree()
346
      <<" a_phi= "<<setw(8)<<Angle(v_a_phi[i]).ToDegree()<<" +/- "<<setw(8)<<Angle(v_err_a_phi[i]).ToDegree()
347
      <<"  RcFit="<<setw(6)<<v_RcFit_cx[i]<<" Xi2Red="<<setw(8)<<v_xi2red_cx[i]<<endl;
348 349 350 351
  }
  return os;
}

352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371

//------------------------ TrackSet -------------------------------------
TrackSet::TrackSet(TrackSet const & a)
  : v_time_sat(a.v_time_sat), v_sat_elev(a.v_sat_elev), v_sat_azim(a.v_sat_azim),
    v_interp_elev(a.v_interp_elev), v_interp_sazim(a.v_interp_sazim)							   
{
}

TrackSet & TrackSet::operator = (TrackSet const & a)
{
  v_time_sat=a.v_time_sat;  v_sat_elev=a.v_sat_elev;  v_sat_azim=a.v_sat_azim;
  v_interp_elev=a.v_interp_elev;  v_interp_sazim=a.v_interp_sazim; 
  return *this;
}

TrackSet::TrackSet(TrkInputDataSet & tkds)
{
  ReadData(tkds);
}

372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426
size_t TrackSet::ReadTrackFile(string flnm, vector<double> & tims, vector<double> & elevs, vector<double> & azims, SLinInterp1D & li_elev, SLinInterp1D & li_sazim)
{
  cout <<"TrackSet::ReadTrackFile() Extracting data from source/satellite track DataTables: Filename= " << flnm << endl;
  DataTable dt_sat;
  PInPersist pin(flnm);
  pin >> dt_sat;
  dt_sat.SetShowMinMaxFlag(true);
  size_t ktime = dt_sat.IndexNom("timesec");
  dt_sat.GetColumn(ktime, tims);
  size_t kelev = dt_sat.IndexNom("elevation");
  dt_sat.GetColumn(kelev, elevs);
  size_t kazim = dt_sat.IndexNom("azimuth");
  dt_sat.GetColumn(kazim, azims);
  li_elev.DefinePoints(tims, elevs);
  double last_azim=azims[0];
  //    vector<double> cazim(v_sat_azim[j].size());
  // azimuth values, shifted possibly +360 +720 deg ... to avoid jumping from 360 deg to 0 deg  
  vector<double> shifted_azim(azims.size());   
  double azim_offset=0.;
  double min_azim_offset=0.;
  bool fgneg_azim_offset=false;
  for(size_t k=0; k<azims.size(); k++)  {
    double azim=azims[k];
    if ((k>0)&&(azim<last_azim)) {
      if ((last_azim>300.)&&(azim<60.))  {
	azim_offset += 360.;
	if (_prtlevel_>0) 
	  cout << "TrackSet::ReadTrackFile()/Info-Warning: 360 to 0 deg. Jump k="<<k<<" last_azim="<<last_azim<<" azimuth= "<<azim<<" Offset->"<<azim_offset<<endl;
      }
    }
    else if ((k>0)&&(azim>last_azim)) {
      if ((last_azim<60)&&(azim>300.))  {
	azim_offset -= 360.;
	if (_prtlevel_>0) 
	  cout << "TrackSet::ReadTrackFile()/Info-Warning: 0 to 360 deg. Jump: k="<<k<<" last_azim="<<last_azim<<" azimuth= "<<azim<<" Offset->"<<azim_offset<<endl;
      }
    }
    if (azim_offset<min_azim_offset)  min_azim_offset=azim_offset;
    last_azim = azim;
    shifted_azim[k]=azim+azim_offset;
    /*
      double phisrcdeg=90.-v_sat_azim[j][k];
      if (phisrcdeg<0.)  phisrcdeg+=360.;
      double phisrc=Angle(phisrcdeg,Angle::Degree).ToRadian();
      cazim[k]=cos(phisrc);
    */
  }
  if (min_azim_offset < -300.) {
    cout << "TrackSet::ReadTrackFile()/Info-Warning: - correcting for negative azim_offset -> Adding " << -min_azim_offset <<" deg."<<endl;
    for(size_t k=0; k<shifted_azim.size(); k++)   shifted_azim[k] -= min_azim_offset;
  }
  li_sazim.DefinePoints(tims, shifted_azim);
  return tims.size();
}

427 428 429 430 431 432 433 434 435 436 437 438 439 440
size_t TrackSet::ReadData(TrkInputDataSet & tkds)
{
  cout << "---- TrackSet::ReadData() ; reading source (satellites, ..) for "
       <<tkds.NbTrk()<<" tracks ..."<<endl;
  if (tkds.NbTrk() != v_time_sat.size()) {
    v_time_sat.resize(tkds.NbTrk());
    v_sat_elev.resize(tkds.NbTrk());
    v_sat_azim.resize(tkds.NbTrk());
    v_interp_elev.resize(tkds.NbTrk());
    v_interp_sazim.resize(tkds.NbTrk());
  }

  for(size_t j=0; j<tkds.NbTrk(); j++) {
    string flnm = tkds.input_base_path+tkds.trkflnm[j]+".ppf";
441 442 443 444 445
    size_t npts=ReadTrackFile(flnm, v_time_sat[j], v_sat_elev[j], v_sat_azim[j], v_interp_elev[j], v_interp_sazim[j]);
    cout<<"["<<j+1<<"]  DONE timevec.size()="<<npts<<"  SLinInterp1D for elevation / azimuth created ..."<<endl;
    if (_prtlevel_>0) {
      cout << v_interp_elev[j];
      cout << v_interp_sazim[j];
446
    }
447 448 449 450 451
  }
  return 0;
}


452
//------------------------ ACxSetFitter -------------------------------------
453
ACxSetFitter::ACxSetFitter(AcxDataSet & data, TrackSet & tks)
454
  : fggaussbeam_(true), D_dish(5.), acxd_(data), tks_(tks), fit_ac_done(false), fit_cx_done(false), 
455 456 457 458 459
    v_RcFit_ac(tks.getNbAutoCor()), v_xi2red_ac(tks.getNbAutoCor()),
    v_Ddish(tks.getNbAutoCor()), v_thetaant(tks.getNbAutoCor()), 
    v_phiant(tks.getNbAutoCor()), v_A(tks.getNbAutoCor()), v_B(tks.getNbAutoCor()), 
    v_err_Ddish(tks.getNbAutoCor()), v_err_thetaant(tks.getNbAutoCor()), 
    v_err_phiant(tks.getNbAutoCor()), v_err_A(tks.getNbAutoCor()), v_err_B(tks.getNbAutoCor()), 
460
    v_acbeams(tks.getNbAutoCor()),
461
    v_RcFit_cx(tks.getNbCrossCor()), v_xi2red_cx(tks.getNbCrossCor()),
462 463 464 465
    v_phase(tks.getNbCrossCor()), v_phi_0(tks.getNbCrossCor()), v_a_phi(tks.getNbCrossCor()), 
    v_Acx(tks.getNbCrossCor()), v_Bcx(tks.getNbCrossCor()), 
    v_err_phi_0(tks.getNbCrossCor()), v_err_a_phi(tks.getNbCrossCor()), 
    v_err_Acx(tks.getNbCrossCor()), v_err_Bcx(tks.getNbCrossCor()),
466
    v_cxbeams(tks.getNbCrossCor())
467
{
468 469 470 471
  if (data.NbTrk() != tks.NbTrk())
    throw ParmError("ACxSetFitter(data, tks) NOT same number of tracks NbTrk() in data and tks");
  if (data.NbTrk() < 1)
    throw ParmError("ACxSetFitter(data, tks) 0 tracks in data data.NbTrk()<1 ");
472 473 474 475 476 477 478 479
}

int ACxSetFitter::doACfit(string outfilename)
{
  cout << "======================================================================================"<<endl;
  cout << "---- ACxSetFitter::doACfit() ; Performing antenna pointing fit ..."<<endl;
  ofstream ofr(outfilename.c_str());
  ofr << "#### Pointing/dish diameter fit on autocorrelation (ACxSetFitter::doACfit() "<<endl
480
      << "## NumAntenna RcFit Xi2red  Deff err_Deff  Elevation err_Elev  Azimuth err_Ezim  A0 err_A0 B0 err_B0 A1 err_A1 B1 err_B1 ..."<<endl;
481 482 483 484 485 486 487 488
  size_t NB_ANTENNES = acxd_.getNbAutoCor();
  size_t NTRK = acxd_.NbTrk();

  for(size_t ii=0; ii<NB_ANTENNES; ii++)  { 
    v_A[ii].resize(NTRK);     v_B[ii].resize(NTRK); 
    v_err_A[ii].resize(NTRK);     v_err_B[ii].resize(NTRK); 
  }
  int tot_npoints_fit = 0;
489 490 491 492
  for(size_t j=0; j<NTRK; j++) {
    if (acxd_.v_noAC[j])  continue;
    tot_npoints_fit += acxd_.v_time_data[j].size();
  }
493 494 495 496
  for(size_t ii=0; ii<NB_ANTENNES; ii++) {
    cout << "-------- doACfit() 1."<<ii+1<<" Creating General Fit for AutoCor Antenna= " << ii+1 << endl;
    GeneralFitData gdata(1, tot_npoints_fit);
    for(size_t j=0; j<NTRK; j++) {
497
      if (acxd_.v_noAC[j])  continue;
498 499 500 501 502 503
      vector< vector<double> > & v_data = acxd_.vv_data[j];
      vector< vector<double> > & v_err = acxd_.vv_err[j];
      for(size_t k=0; k<acxd_.v_time_data[j].size(); k++) {
	gdata.AddData1(acxd_.v_time_data[j][k],v_data[ii][k],v_err[ii][k]); // Fill x, y and error on y     
      }
    }
504 505 506 507
    TkF_ACXi2 gxi2( acxd_.v_time_data, acxd_.vv_data, acxd_.vv_err, acxd_.v_freqs, acxd_.v_noAC, 
		    tks_.v_interp_elev, tks_.v_interp_sazim, ii, fggaussbeam_);  // MyACGenXi2
    //    GeneralFit mFit(&gxi2);
    TkF_Fitter mFit(gxi2);
508
    mFit.SetData(&gdata);        // connect data to the fitter , here the data is unused - gxi2 includes its data 
509
    mFit.SetMaxStep(5000);
510
    // SetParam(int n,double value, double step,double min=1., double max=-1.);
511
    mFit.SetParam(0,"D_dish",D_dish,0.1,D_dish*0.7,D_dish*1.4);
512 513 514 515 516 517 518 519 520 521 522 523 524
    // mFit.SetFix(0, D_dish);
    
    double thetaAntenne=0., phiAntenne=0.;
    if (fabs(acxd_.zenang)>1.e-6) {
      if (acxd_.zenang<0)  {
	thetaAntenne=Angle(-acxd_.zenang,Angle::Degree).ToRadian();
	phiAntenne=Angle(270.,Angle::Degree).ToRadian();
      }
      else {
	thetaAntenne=Angle(acxd_.zenang,Angle::Degree).ToRadian();
	phiAntenne=Angle(90.,Angle::Degree).ToRadian();
      }
    }
525
    mFit.SetParam(1,"ThetaAntenne",thetaAntenne,M_PI/1440,0.,M_PI/4.); // thetaAntenne+M_PI/30.); // 
526 527 528
    mFit.SetParam(2,"PhiAntenne",phiAntenne,M_PI/180.,0.,2.*M_PI);
    // mFit.SetFix(1, thetaAntenne);
    // mFit.SetFix(2, phiAntenne);
529

530
    //DEL    size_t jj=0;
531
    for(size_t j=0; j<NTRK; j++) {
532 533
      double A = acxd_.v_max_data[j][ii];
      double B = acxd_.v_min_data[j][ii];
534 535 536 537 538
      A -= B;   
      if (A<1.e-9)  { 
	cout << " doACfit()/Warning NumAnt/ii="<<ii<<" NumTrk/j="<<j<<" Negative A , A="<<A<<" B="<<B<<" A->"<<0.1*B<<endl;
	A=0.1*B;
      }
539 540
      v_A[ii][j]=A;   v_err_A[ii][j]=0.;
      v_B[ii][j]=B;   v_err_B[ii][j]=0.;
541
      //DEL if (acxd_.v_noAC[j])  continue;
542 543
      char pname[32];
      sprintf(pname,"A%d",(int)(j+1));
544
      mFit.SetParam(2*j+3,pname,A,A/10.,A/20,A*5);
545
      sprintf(pname,"B%d",(int)(j+1));
546 547 548 549 550 551 552
      mFit.SetParam(2*j+4,pname,B,B/10.,B/20,B*5);
      // mFit.SetFix(2*jj+4, B);
      if (acxd_.v_noAC[j]) {
	mFit.SetFix(2*j+3, A);
	mFit.SetFix(2*j+4, B);
      }
      //DEL      jj++;
553 554 555
    }
    //DBG mFit.PrintFit();
    //    cout << "do_p4_trkfit 2."<<ii+1<<" Performing the fit for AutoCor Antenna= " << ii+1 << endl;
556
    int rcfit = mFit.doFit();
557
    double xi2red=mFit.GetChi2Red();
558
    if (_prtlevel_>1) mFit.PrintFit();
559
    v_RcFit_ac[ii]=rcfit;  v_xi2red_ac[ii]=xi2red;
560 561 562 563 564 565 566 567 568
    if(rcfit>0) { 
      cout<< "------- Fit result for Antenna No="<<ii+1<<" Reduce_Chisquare = " << mFit.GetChi2Red()
	  << " nstep="<<mFit.GetNStep() << " rc="<<rcfit<<endl;
    }
    else {
      cout << "---Fit failed for "<<ii+1<<"--- Fit_Error, rc = " << rcfit << "  nstep="<<mFit.GetNStep()<<endl;
      ofr <<setw(4)<<ii+1<<" ERROR FIT RC="<<rcfit<<"  nstep="<<mFit.GetNStep()<<endl;
      if (_prtlevel_>0) mFit.PrintFitErr(rcfit);
    } 
569 570 571 572 573 574 575 576 577

    ofr <<setw(4)<<ii+1<<" "<<setw(8)<<mFit.GetChi2Red()<<" "; 
    double Dfit=mFit.GetParm(0);   double err_Dfit=mFit.GetParmErr(0);
    cout <<setw(16)<<"DishDiameter= "<<setw(10)<<Dfit<<" +/- "<<setw(10)<<err_Dfit<<" m."<<endl;
    ofr <<setw(5)<<rcfit<<" "<<setw(8)<<Dfit<<" "<<setw(8)<<err_Dfit<<"  "; 
    v_Ddish[ii]=Dfit;
    v_err_Ddish[ii]=err_Dfit;
    double thetaant=mFit.GetParm(1);   double err_thetaant=mFit.GetParmErr(1);
    v_thetaant[ii]=thetaant;
578
    v_err_thetaant[ii]=err_thetaant;
579 580 581 582 583 584 585 586 587 588 589
    double elevdeg=90.-Angle(thetaant).ToDegree();
    double err_elevdeg=Angle(err_thetaant).ToDegree();
    cout <<setw(16)<<"ThetaAntenne= "<<setw(12)<<Angle(thetaant).ToDegree()<< " +/- "
	 <<setw(12)<<Angle(err_thetaant).ToDegree()<<" (elevation="
	 <<setw(8)<<elevdeg<<" +/- "<<setw(8)<<err_elevdeg<<") deg."<<endl;
    ofr <<setw(8)<<elevdeg<<" "<<setw(8)<<err_elevdeg<<"  "; 
    double phiant=mFit.GetParm(2);   double err_phiant=mFit.GetParmErr(2);
    double azimdeg=90.-Angle(phiant).ToDegree();
    if (azimdeg<0.)  azimdeg += 360.;
    double err_azimdeg=Angle(err_phiant).ToDegree();
    v_phiant[ii]=phiant;
590
    v_err_phiant[ii]=err_phiant;
591 592 593
    cout <<setw(16)<<"PhiAntenne= "<<setw(12)<<Angle(phiant).ToDegree()<< " +/- "
	 <<setw(12)<<Angle(err_phiant).ToDegree()<<" (azimuth  ="
	 <<setw(8)<<azimdeg<<" +/- "<<setw(8)<<err_azimdeg<<" ) deg."<<endl;
594
    ofr <<setw(8)<<azimdeg<<" "<<setw(8)<<err_azimdeg<<"  ";
595
    //DEL    jj=0;
596
    for(size_t j=0; j<NTRK; j++) {
597 598 599
      double A=1.,B=0.,err_A=0.,err_B=0.;
      A=mFit.GetParm(3+2*j);    err_A=mFit.GetParmErr(3+2*j);
      B=mFit.GetParm(4+2*j);    err_B=mFit.GetParmErr(4+2*j);
600
      v_A[ii][j]=A;  v_err_A[ii][j]=err_A;  v_B[ii][j]=B;  v_err_B[ii][j]=err_B;
601
      cout << "  Trk/Sat["<<j<<"] -> A= "<<A<<" +/- "<<err_A<<"  B= "<<B<<" +/- "<<err_B<<(acxd_.v_noAC[j]?" FIXED":"")<<endl;
602
      if (acxd_.v_noAC[j])
603
	ofr <<setw(8)<<A<<" "<<setw(8)<<" NOFIT "<<" "<<setw(8)<<B<<" "<<setw(8)<<" FIXED "<<" ";
604 605
      else 
	ofr <<setw(8)<<A<<" "<<setw(8)<<err_A<<" "<<setw(8)<<B<<" "<<setw(8)<<err_B<<" ";
606 607 608 609 610 611 612 613 614 615 616
    }
    ofr << endl;
    double clight = PhysQty::c().SIValue();
    double lambda = clight/(acxd_.v_freqs[0]*1.e6);
    ACBeam acb1(Dfit, thetaant, phiant, lambda);
    acb1.setGaussianLobe(fggaussbeam_);
    ACBeam acb2(Dfit, thetaant, phiant, lambda);
    acb2.setGaussianLobe(fggaussbeam_);
    Vector3d baseline0(0.,0.,0.);
    v_acbeams[ii]=CxBeam(acb1, acb2, baseline0);
    
617 618
  }
  
619 620
  fit_ac_done=true;
  acxd_.v_acbeams=v_acbeams;
621 622 623 624
  acxd_.v_RcFit_ac=v_RcFit_ac;
  acxd_.v_xi2red_ac=v_xi2red_ac;
  acxd_.v_Ddish=v_Ddish;
  acxd_.v_thetaant=v_thetaant;
625
  acxd_.v_phiant=v_phiant;
626 627 628
  acxd_.v_err_Ddish=v_err_Ddish;
  acxd_.v_err_thetaant=v_err_thetaant;
  acxd_.v_err_phiant=v_err_phiant;
629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644
  return 0;
}

int ACxSetFitter::saveExpectedAC(string outcheckfilename)
{
  if (outcheckfilename.length()<1)  return 1;
  cout << "-----ACxSetFitter::saveExpectedAC() : computing expected signal for fitted params , will be saved to file "
       <<outcheckfilename<<endl;
  POutPersist pos(outcheckfilename);
  size_t NB_ANTENNES = acxd_.getNbAutoCor();
  size_t NTRK = acxd_.NbTrk();

  for(size_t ii=0; ii<NB_ANTENNES; ii++)     {
    if (_prtlevel_>1) 
      cout << "... Computing DataSignal & Expected Signal for fitted params and dish "<<ii+1<<endl;
    
645
    MyACSignal macs(acxd_.v_time_data, acxd_.vv_data, acxd_.vv_err, acxd_.v_freqs, acxd_.v_noAC, 
646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661
		    tks_.v_interp_elev, tks_.v_interp_sazim, ii, fggaussbeam_);
      
    double Ddishfit=v_Ddish[ii];
    double thetafit=v_thetaant[ii];
    double phifit=v_phiant[ii];
    
    char oname[32];
    for(size_t j=0; j<NTRK; j++)  {
      double A = v_A[ii][j];
      double B = v_B[ii][j];
      Vector signal = macs.getDataSignal(j);
      sprintf(oname,"ac_%d_%d",(int)ii+1,(int)j+1);
      pos << PPFNameTag(oname)<<signal;
      Vector expsignal = macs.getExpectedSignal(j, Ddishfit, thetafit, phifit, A, B);
      sprintf(oname,"simac_%d_%d",(int)ii+1,(int)j+1);      
      pos << PPFNameTag(oname)<<expsignal;
662 663 664 665 666
      if (ii==0)  {
	Vector tmvec = macs.getTimeVec(j);
	sprintf(oname,"tim_%d",(int)j+1);
	pos << PPFNameTag(oname)<<tmvec;
      }
667 668 669 670 671 672
    }
  } 
  return 0;
}


673
int ACxSetFitter::doCxfit(string outfilenamecx, bool useAac, bool fg_B0, bool fgphi0only)
674 675 676 677 678
{
  size_t NB_ANTENNES=acxd_.getNbAutoCor();   // nombre d'antennes 
  size_t NB_CXCORS=acxd_.getNbCrossCor();
  size_t NTRK = acxd_.NbTrk();

679

680
  //------ Valeurs de phases et pente de Phi(freq) = Phi0 + aPhi * (freq-1250.)/250.   (en degres) 
681 682 683
  double phi0deg_I[6] = {250. , 110, 60., 220., 170., 310.};
  double aphideg_I[6] = {-71. , 383., 449., 454., 520., 66.};

684
  cout << "======================================================================================"<<endl;
685
  cout << "---------- ACxSetFitter::doCxfit() ; Performing cross-cor phase fit for NTrk="<<NTRK<<endl;
686 687 688
  if (useAac) cout << " ... Using Amplitude from auto-correlations fit for initial fit parameter value..."<<endl; 
  ofstream ofr(outfilenamecx.c_str());
  ofr << "#### cross-cor phase fit (ACxSetFitter::doCxfit() ) "<<endl
689
      << "## NumCxCor RcFit Xi2red Phi0 err_Phi0 a_Phi err_a_Phi (deg) A0 err_A0 B0 errB0 A1 err_A1  ..."<<endl;
690 691
  int tot_npoints_fit = 0;
  for(size_t j=0; j<NTRK; j++) tot_npoints_fit += 2*(acxd_.v_time_data[j].size());
692 693
  cout << " Total number of data points for fit="<< tot_npoints_fit<<endl;

694 695 696 697
  size_t Anum1[6]={0,0,0,1,1,2};
  size_t Anum2[6]={1,2,3,2,3,3};
  for(size_t ii=0; ii<NB_CXCORS; ii++) {
    v_Acx[ii].resize(NTRK);   
698 699 700
    v_Bcx[ii].resize(NTRK);  
    v_err_Acx[ii].resize(NTRK);   
    v_err_Bcx[ii].resize(NTRK);  
701 702
    for(size_t j=0; j<NTRK; j++) {
      v_Acx[ii][j]=1.;   v_Bcx[ii][j]=complex<double>(0.,0.);
703
      v_err_Acx[ii][j]=1.;   v_err_Bcx[ii][j]=complex<double>(0.,0.);
704
    }
705 706 707 708 709 710 711 712 713 714 715
    Vector3d baseline=P4Coords::getBaseline(Anum1[ii]+1,Anum2[ii]+1);
    cout << "--------- 1."<<ii+1<<" doCxfit() Doing fit for CrossCor= " << ii << " FxF= " 
	 << Anum1[ii]+1<<"x"<<Anum2[ii]+1<<" Baseline="<<baseline<<endl;
    GeneralFitData gdata(1, tot_npoints_fit);
    for(size_t j=0; j<NTRK; j++) {
      vector< vector< complex<double> > > & v_cxdata = acxd_.vv_cxdata[j];
      vector< vector<double> > & v_cxerr = acxd_.vv_cxerr[j];
      for(size_t k=0; k<acxd_.v_time_data[j].size(); k++) {
	gdata.AddData1(acxd_.v_time_data[j][k],v_cxdata[ii][k].real(),v_cxerr[ii][k]); // Fill x, y and error on y
	gdata.AddData1(acxd_.v_time_data[j][k],v_cxdata[ii][k].imag(),v_cxerr[ii][k]); // Fill x, y and error on y     
      }
716
    }
717 718 719 720 721 722 723
    double clight = PhysQty::c().SIValue();
    double lambda = clight/(acxd_.v_freqs[0]*1.e6);
    ACBeam acb1(v_Ddish[Anum1[ii]], v_thetaant[Anum1[ii]], v_phiant[Anum1[ii]], lambda);
    acb1.setGaussianLobe(fggaussbeam_);
    ACBeam acb2(v_Ddish[Anum2[ii]], v_thetaant[Anum2[ii]], v_phiant[Anum2[ii]], lambda);
    acb2.setGaussianLobe(fggaussbeam_);
    CxBeam cxbeam(acb1, acb2, baseline);
724
    v_cxbeams[ii]=cxbeam;
725

726
    TkF_CxXi2 gxi2( acxd_.v_time_data, acxd_.vv_cxdata, acxd_.vv_cxerr, acxd_.v_freqs, acxd_.v_noCx, 
727 728 729
		    tks_.v_interp_elev, tks_.v_interp_sazim, cxbeam, ii);  // MyCxGenXi2
    //    GeneralFit mFit(&gxi2);
    TkF_Fitter mFit(gxi2);
730
    mFit.SetData(&gdata);        // connect data to the fitter , here the data is unused - gxi2 includes its data 
731
    mFit.SetMaxStep(3000);
732
    // SetParam(int n,double value, double step,double min=1., double max=-1.);
733 734
    mFit.SetParam(0,"Phi_0",Angle(phi0deg_I[ii],Angle::Degree).ToRadian(),M_PI/180.,-0.5*M_PI,3*M_PI);
    mFit.SetParam(1,"a_phi",Angle(aphideg_I[ii],Angle::Degree).ToRadian(),0.05,-15.,15.);
735 736 737 738 739
    if (fgphi0only) {
      cout << " ACxSetFitter::doCxfit() Fitting Phi0 Only (frequency independent phase)"<<endl;
      mFit.SetFix(1,0.);
    }
    else cout << " ACxSetFitter::doCxfit() Fitting  Phase(freq) = Phi0 + a_Phi * (freq-1250.)/250. "<<endl;
740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756

    char oname[32];
    vector<double> v_amp(NTRK);
    for(size_t j=0; j<NTRK; j++) {
      double A=1.; // v_max_cxdata[j][ii]; 
      TVector< complex<double> >  signal = gxi2.getDataSignal(j);
      Vector asig = SOPHYA::abs(signal);
      double mins, maxs;
      asig.MinMax(mins, maxs);
      TVector< complex<double> >  expsignal = gxi2.getExpectedSignal(j, 0., A);
      Vector aexpsig = SOPHYA::abs(expsignal);
      double mine, maxe;
      aexpsig.MinMax(mine, maxe);
      A=maxs/maxe;
      v_amp[j]=A; 
    }

757
    double fparm[100];  fparm[0]=0.;
758 759
    fparm[0]=Angle(phi0deg_I[ii], Angle::Degree).ToRadian();  // Angle(phi0deg_I[ii],Angle::Degree).ToRadian()
    if (fgphi0only)  fparm[1]=0.;
760
    else  fparm[1]=Angle(aphideg_I[ii], Angle::Degree).ToRadian();;
761
    
762 763 764 765
    double bestxi2 = 9.e19;
    double bestphase=0.;
    int bestnpts,npts;
    int bestafact;
766
    double afact[12]={0.15,0.3,0.5,0.75,1.0,1.25,1.5,1.75,2.0,2.4,2.8,3.2};
767
    bool fg_ph_I=false;   // if true , phase value from phi0deg_I aphideg_I 
768 769 770
    for(int ia=0; ia<12; ia++) {
      for(size_t j=0; j<NTRK; j++) {
	double Aac=sqrt(v_A[Anum1[ii]][j] * v_A[Anum2[ii]][j]);
771 772
	//DBG	cout << " *DBG* j="<<j<<" ia="<<ia<<" vA="<<v_A[Anum1[ii]][j]<<" x "<<v_A[Anum2[ii]][j]
	//     <<"  -> "<<Aac<<endl;
773
	fparm[2+3*j]=(useAac?Aac:v_amp[j]);
774
	fparm[2+3*j]*=afact[ia];   fparm[3+3*j]=fparm[4+3*j]=0.;
775
      }
776 777 778 779
      for(int jp=-1; jp<180; jp++) {
	double ph = jp*2.;
	if (jp==-1) fparm[0]=Angle(phi0deg_I[ii], Angle::Degree).ToRadian(); 
	else  fparm[0]=Angle(ph, Angle::Degree).ToRadian();
780
	double xi2 = gxi2.getXi2(fparm, npts);
781
	//DBG	cout << " *DBG* ia="<<ia<<" afact="<<afact[ia]<<" ph="<<ph<<" xi2="<<xi2<<endl;
782
	if (xi2 < bestxi2) {
783
	  if (jp==-1) fg_ph_I=true;  else fg_ph_I=false; 
784 785 786 787
	  bestxi2 = xi2; bestphase=fparm[0]; bestnpts=npts;  bestafact=afact[ia];
	}
      }
    }
788
    mFit.SetParam(0,"Phi_0",bestphase,M_PI/720.,-0.5*M_PI,2.5*M_PI);
789 790 791
    cout << "2."<<ii+1<<" Scan param bestxi2_red="<<bestxi2/(double)(tot_npoints_fit-(2+NTRK))
	 <<" bestphase="<<Angle(bestphase).ToDegree()<<(fg_ph_I?" (Phase from phi0deg_I)":" ")
	 <<" bestnpts="<<bestnpts<<" bestafact="<<bestafact<< " A= ";  
792
    v_phi_0[ii]=bestphase;
793 794
    for(size_t j=0; j<NTRK; j++)  {
      cout << v_amp[j] << " , ";  
795
      double Aac=sqrt(v_A[Anum1[ii]][j] * v_A[Anum2[ii]][j]);
796 797 798 799 800 801 802 803 804
      v_Acx[ii][j]=(useAac?Aac:v_amp[j]);
    }
    cout << endl;
    for(size_t j=0; j<NTRK; j++) {
      char pname[32];
      sprintf(pname,"A%d",(int)(j+1));
      double Aac=sqrt(v_A[Anum1[ii]][j] * v_A[Anum2[ii]][j]);
      double A=(useAac?Aac:v_amp[j]);
      //DBG      cout << "*DBG* j="<<j<<" Aac= "<<Aac<<" v_amp="<<v_amp[j]<<"  A= "<<A<<"  A1="<<v_A[Anum1[ii]][j]<<" A2="<<v_A[Anum2[ii]][j]<<endl;
805
      mFit.SetParam(2+3*j,pname,A,A/10.,A/4,A*4);
806 807
      sprintf(pname,"Bre%d",(int)(j+1));
      mFit.SetParam(3+3*j,pname,0.,A/25.,-A/5,A/5.);
808 809
      sprintf(pname,"Bim%d",(int)(j+1));
      mFit.SetParam(4+3*j,pname,0.,A/25.,-A/5,A/5.);
810 811
      if (acxd_.v_noCx[j]) {
	mFit.SetFix(2+3*j,A);
812 813 814
	mFit.SetFix(3+3*j,0.);
	mFit.SetFix(4+3*j,0.);
      }
815 816 817 818 819 820
      else {
	if (fg_B0) {
	  mFit.SetFix(3+3*j,0.);
	  mFit.SetFix(4+3*j,0.);
	}
      }
821 822 823 824
    }
    //DBG mFit.PrintFit();
    if (_prtlevel_>1)    
      cout << " 3."<<ii+1<<" Performing the fit for CrossCor " << ii << " FxF= " << Anum1[ii]+1<<"x"<<Anum2[ii]+1<<endl;
825
    int rcfit = mFit.doFit();
826
    v_RcFit_cx[ii]=rcfit;   v_xi2red_cx[ii]=mFit.GetChi2Red();
827 828 829
    if (_prtlevel_>1) mFit.PrintFit();
    if(rcfit>0) { 
      //      cout<< "-------------------------- Result for Cross No " << ii << endl; 
830
      cout<< "------ Fit result for Cross No "<<ii+1<<" Reduce_Chisquare = " << mFit.GetChi2Red()
831 832 833 834 835 836 837
	  << " nstep="<<mFit.GetNStep() << " rc="<<rcfit<<endl;
    }
    else {
      cout << "---Fit failed for "<<ii<<" Fit_Error, rc = " << rcfit << "  nstep="<<mFit.GetNStep()<<endl;
      ofr <<setw(4)<<ii+1<<" ERROR FIT RC="<<rcfit<<"  nstep="<<mFit.GetNStep()<<endl;
      if (_prtlevel_>0) mFit.PrintFitErr(rcfit);
    }
838 839

    ofr <<setw(4)<<ii+1<<" "<<setw(5)<<rcfit<<setw(8)<<mFit.GetChi2Red()<<" "; 
840 841 842 843 844 845 846 847 848 849 850 851
    double phi0=mFit.GetParm(0);   double err_phi0=mFit.GetParmErr(0);
    double aphi=mFit.GetParm(1);   double err_aphi=mFit.GetParmErr(1);
    // on calcule la phase ajustee pour la frequence de reference 1300 MHz 
    double phase=gxi2.getPhase4Freq(phi0,aphi,1300.);
    while (phase<0.) phase += 2.*M_PI;
    while (phase>2.*M_PI) phase -= 2.*M_PI;
    cout <<"Phase(@1300MHz)= "<<setw(10)<<Angle(phase).ToDegree()<<"  phi_0= "<<setw(10)
	 <<Angle(phi0).ToDegree()<<" +/- "<<setw(10)<<Angle(err_phi0).ToDegree()<<" deg."
	 <<" a_phi= "<<setw(8)<<Angle(aphi).ToDegree()<<" +/- "<<setw(10)
	 <<Angle(err_aphi).ToDegree()<<" deg/250 MHz"<<endl;
    ofr <<setw(8)<<Angle(phi0).ToDegree()<<" "<<setw(8)<<Angle(err_phi0).ToDegree()<<"  "
	<<setw(8)<<Angle(aphi).ToDegree()<<" "<<setw(8)<<Angle(err_aphi).ToDegree()<<"  ";
852
    v_phase[ii]=phase;
853 854 855 856
    v_phi_0[ii]=phi0;
    v_err_phi_0[ii]=err_phi0;
    v_a_phi[ii]=aphi;
    v_err_a_phi[ii]=err_aphi;
857 858 859
    for(size_t j=0; j<NTRK; j++) {
      double Aac=sqrt(v_A[Anum1[ii]][j] * v_A[Anum2[ii]][j]);
      double Ai=(useAac?Aac:v_amp[j]);
860
      double A=mFit.GetParm(2+3*j);   double err_A=mFit.GetParmErr(2+3*j);
861 862 863 864
      complex<double> B(mFit.GetParm(3+3*j), mFit.GetParm(4+3*j));
      complex<double> err_B(mFit.GetParmErr(3+3*j), mFit.GetParmErr(4+3*j));
      cout << "  Trk["<<j<<"]  A= "<<A<<" +/- "<<err_A<<"  (A/Ai="<<A/Ai<<")"<<
	" B= "<<B<<" +/- "<<err_B<<(acxd_.v_noCx[j]?" NoFIT":" ")<<endl; 
865
      v_Acx[ii][j]=A;  
866
      v_Bcx[ii][j]=B;
867
      v_err_Acx[ii][j]=err_A; 
868
      ofr <<setw(8)<<A<<" "<<setw(8)<<err_A<<" "<<setw(14)<<B<<" "<<setw(12)<<err_B<<" ";
869
    }
870
    ofr << endl; 
871
  }
872 873 874 875 876 877 878 879 880 881 882
  cout << " --- Fitted phases: ";
  for(size_t i=0; i<NB_CXCORS; i++) cout<<setw(6)<<Angle(v_phase[i]).ToDegree()<<" ; ";   cout<<endl;
  double dphi23=v_phase[1]-v_phase[0];   if (dphi23<0.) dphi23+=(2.*M_PI);
  double dphi24=v_phase[2]-v_phase[0];   if (dphi24<0.) dphi24+=(2.*M_PI);
  double dphi34=v_phase[2]-v_phase[1];   if (dphi34<0.) dphi34+=(2.*M_PI);
  cout<<" Cx-2x3: "<<setw(6)<<Angle(dphi23).ToDegree()<<" ==? "<<setw(6)<<Angle(v_phase[3]).ToDegree()<<endl;
  cout<<" Cx-2x4: "<<setw(6)<<Angle(dphi24).ToDegree()<<" ==? "<<setw(6)<<Angle(v_phase[4]).ToDegree()<<endl;
  cout<<" Cx-3x4: "<<setw(6)<<Angle(dphi34).ToDegree()<<" ==? "<<setw(6)<<Angle(v_phase[5]).ToDegree()<<endl;
  ofr<<"# Cx-2x3: "<<setw(6)<<Angle(dphi23).ToDegree()<<" ==? "<<setw(6)<<Angle(v_phase[3]).ToDegree()<<endl;
  ofr<<"# Cx-2x4: "<<setw(6)<<Angle(dphi24).ToDegree()<<" ==? "<<setw(6)<<Angle(v_phase[4]).ToDegree()<<endl;
  ofr<<"# Cx-3x4: "<<setw(6)<<Angle(dphi34).ToDegree()<<" ==? "<<setw(6)<<Angle(v_phase[5]).ToDegree()<<endl;
883 884
  fit_cx_done=true;
  acxd_.v_cxbeams=v_cxbeams;
885 886
  acxd_.v_RcFit_cx=v_RcFit_cx;
  acxd_.v_xi2red_cx=v_xi2red_cx;
887
  acxd_.v_phase=v_phase;
888 889
  acxd_.v_phi_0=v_phi_0;
  acxd_.v_a_phi=v_a_phi;
890 891
  acxd_.v_err_phi_0=v_err_phi_0;
  acxd_.v_err_a_phi=v_err_a_phi;
892 893
  acxd_.v_Acx=v_Acx;
  acxd_.v_Bcx=v_Bcx;
894 895
  return 0;
} 
896

897

898 899 900 901 902 903 904 905 906 907 908
int ACxSetFitter::saveExpectedCx(string outcheckfilename)
{
  cout << "ACxSetFitter::saveExpectedCx() saving expected cross-cor (and visi-data) to file "<<outcheckfilename<<endl;
  POutPersist pox(outcheckfilename);
  size_t NB_CXCORS=acxd_.getNbCrossCor();
  size_t NTRK = acxd_.NbTrk();

  char oname[32];

  for(size_t ii=0; ii<NB_CXCORS; ii++) {
    CxBeam cxbeam=v_cxbeams[ii];
909
    MyCxSignal cxsig( acxd_.v_time_data, acxd_.vv_cxdata, acxd_.vv_cxerr, acxd_.v_freqs, acxd_.v_noCx, 
910 911 912 913 914
		      tks_.v_interp_elev, tks_.v_interp_sazim, cxbeam, ii);
    for(size_t j=0; j<NTRK; j++) {
      TVector< complex<double> >  signal = cxsig.getDataSignal(j);
      sprintf(oname,"cx_%d_%d",(int)ii+1,(int)j+1);
      pox << PPFNameTag(oname)<<signal;
915 916 917
      //DBG      cout << " *DBG* getPhase4Freq() phi0="<<acxd_.v_phi_0[ii]<<" a_phi="<<acxd_.v_a_phi[ii]<<" freq="<<acxd_.v_freqs[j]<<endl;
      double phase=cxsig.getPhase4Freq(acxd_.v_phi_0[ii],acxd_.v_a_phi[ii],acxd_.v_freqs[j]);
      TVector< complex<double> >  expsignal = cxsig.getExpectedSignal(j, phase, v_Acx[ii][j]);
918 919 920 921 922 923 924 925 926 927 928
      sprintf(oname,"simcx_%d_%d",(int)ii+1,(int)j+1);
      pox << PPFNameTag(oname)<<expsignal;
      if (ii==0)  {
	Vector tmvec = cxsig.getTimeVec(j);
	sprintf(oname,"tim_%d",(int)j+1);
	pox << PPFNameTag(oname)<<tmvec;
      }
    }
  }
  return 0;
}  
929 930 931

//------------------------ CxBaselineFitter -------------------------------------
CxBaselineFitter::CxBaselineFitter(vector<AcxDataSet> & v_data, vector<TrackSet> & v_tks)
932 933
  : v_acxd(v_data), v_trks(v_tks), tot_ntrks(0), fit_done(false), simplex_done(false), 
    xi2red(-9.e9), bestfitparam(NULL), err_bestfitparam(NULL)
934 935 936 937 938
{
  if (v_acxd.size() != v_trks.size())
    throw ParmError("CxBaselineFitter::CxBaselineFitter(v_data, v_tks) NOT same size v_data,v_tks ");
  if (v_acxd.size() < 1)
    throw ParmError("CxBaselineFitter::CxBaselineFitter(v_data, v_tks) v_data.size()<1 ");
939
  
940 941 942 943 944
  tot_ntrks=0;
  for(size_t i=0; i<v_acxd.size(); i++) tot_ntrks+=v_acxd[i].NbTrk();
  if (tot_ntrks<1)
    throw ParmError("CxBaselineFitter::CxBaselineFitter(v_data, v_tks) 0 tracks ! tot_ntrks<1 ");

945
  size_t nparam = 5*(v_acxd[0].getNbAutoCor()-1);  // 5 param / antenne , phi0, aphi, dX,dY,dZ
946 947
  bestfitparam = new double[nparam];
  err_bestfitparam = new double[nparam];
948 949

  initFitParams();
950 951
}

952 953 954 955 956
CxBaselineFitter::~CxBaselineFitter()
{
  if (bestfitparam) delete[] bestfitparam;
  if (err_bestfitparam) delete[] err_bestfitparam;
}
957

958 959
void CxBaselineFitter::initFitParams()
{
960 961
  // Positions en z des 3 antennes 2,3,4 / 1 , en metres 
  double z_antennes[3]={-0.10,+0.05,-0.02};
962
  //DBG  cout << " *DBG* CxBaselineFitter::initFitParams() v_acxd[0].v_phase.size()="<<v_acxd[0].v_phase.size()<<endl;
963 964 965 966
  size_t NB_ANTENNES=v_acxd[0].getNbAutoCor();   // nombre d'antennes 
  size_t NB_CXCORS=v_acxd[0].getNbCrossCor();
  if (NB_ANTENNES != 4)
    throw PError("CxBaselineFitter::initFitParams() NB_ANTENNES != 4  Current version works only for 4 antenna");
967 968 969 970
  v_phi_0.resize(v_acxd[0].getNbAutoCor()-1);
  v_err_phi_0.resize(NB_ANTENNES-1);
  v_a_phi.resize(v_acxd[0].getNbAutoCor()-1);
  v_err_a_phi.resize(NB_ANTENNES-1);
971 972 973
  v_baselineshits.resize(NB_ANTENNES-1);
  v_err_baselineshits.resize(NB_ANTENNES-1);
  for(size_t i=0; i<(NB_ANTENNES-1); i++) {
974 975
    v_phi_0[i]=v_acxd[0].v_phi_0[i];   v_err_phi_0[i]=0.;
    v_a_phi[i]=v_acxd[0].v_a_phi[i];   v_err_a_phi[i]=0.;
976 977
    v_baselineshits[i]=Vector3d(0.,0.,0.);
    v_err_baselineshits[i]=Vector3d(0.,0.,0.);
978 979 980 981
    bestfitparam[2*i]=v_phi_0[i];
    err_bestfitparam[2*i]=0.;
    bestfitparam[2*i+1]=v_a_phi[i];
    err_bestfitparam[2*i]=0.;
982 983 984 985
    for(size_t j=0; j<3; j++) {
      bestfitparam[3*(i+1)+j]=err_bestfitparam[3*(i+1)+j]=0.;
    }
  }
986 987
  //DBG  cout << " *DBG* DONE **** CxBaselineFitter::initFitParams()"<<endl;

988 989
}

990
int CxBaselineFitter::dofit(string outfilename, int fgfixbaseline, bool fgphi0only)
991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006
{
  size_t NB_ANTENNES=v_acxd[0].getNbAutoCor();   // nombre d'antennes 
  size_t NB_CXCORS=v_acxd[0].getNbCrossCor();
  cout << "======================================================================================"<<endl;
  cout << "------- CxBaselineFitter::dofit()  Performing baseline/phase fit on the 6 cross-cors "<<" TotNbTracks="<<tot_ntrks<<endl;
  
  ofstream ofr(outfilename.c_str());
  ofr << "####  Fitted phases and baseline-shifts (CxBaselineFitter::dofit() ) "<<endl
      << "## NumAntenna  Phase BaselineShiftX  BaselineShiftY BaselineShiftZ  (Phase in degree, BaselineShift in meter) "<<endl;

  int tot_npoints_fit = 0;
  for(size_t i=0; i<v_acxd.size(); i++)
    for(size_t j=0; j<v_acxd[i].NbTrk(); j++)
      tot_npoints_fit += 2*(v_acxd[i].v_time_data[j].size())*NB_CXCORS;
  cout << " Total number of data points for fit="<< tot_npoints_fit<<endl;
  GeneralFitData gdata(1, tot_npoints_fit);
1007
  int npoints2=