///
/// Return a metric in <i>naturalized</i> coordinates.
/// Frequency coordinates of spindown order \f$s\f$ are scaled by
/// \f[ \frac{2\pi}{(s+1)!} \left(\frac{\overline{\Delta T}}{2}\right)^{s+1} \f]
/// where \f$\overline{\Delta T}\equiv\sum_{k}^{N} \Delta T_k\f$ is the average segment-length
/// over all \f$N\f$ segments.
///
/// Sky coordinates are scaled by Holgers' units, see Eq.(44) in PRD82,042002(2010),
/// without the equatorial rotation in alpha:
/// \f[ \frac{2\pi f R_{E}  \cos(\delta_{D})}{c} \f]
/// where \f$f\f$ is the frequency and
/// \f$R_{E}\f$ the Earth radius, and \f$\delta_{D}\f$ is the detectors latitude.
///
/// If \c p_transform is non-NULL, return the naturalization transform in \c *p_transform.
/// If \c p_gpr_ij is non-NULL, apply the transform to the metric \c *p_gpr_ij.
///
/// \note \c *p_gpr_ij and \c *p_transform will be allocated if \c NULL. \c *p_gpr_ij and \c g_ij
/// may point to the same matrix.
///
int
XLALNaturalizeMetric(
  gsl_matrix **p_gpr_ij,			///< [in,out,optional] Pointer to transformed matrix \f$g^{\prime}_{ij}\f$
  gsl_matrix **p_transform,			///< [in,out,optional] Pointer to naturalization transform
  const gsl_matrix *g_ij,			///< [in] Matrix to transform \f$g_{ij}\f$
  const DopplerMetricParams *metricParams	///< [in] Input parameters used to calculate naturalization transform
  )
{

  // Check input
  XLAL_CHECK( g_ij != NULL, XLAL_EFAULT );
  XLAL_CHECK( g_ij->size1 == g_ij->size2, XLAL_ESIZE );
  if ( p_transform != NULL ) {
    if ( *p_transform != NULL ) {
      XLAL_CHECK( (*p_transform)->size1 == g_ij->size1, XLAL_ESIZE );
      XLAL_CHECK( (*p_transform)->size2 == g_ij->size2, XLAL_ESIZE );
    } else {
      GAMAT( *p_transform, g_ij->size1, g_ij->size2 );
    }
  }
  XLAL_CHECK( metricParams != NULL, XLAL_EFAULT );

  // Compute average segment length
  XLAL_CHECK ( XLALSegListIsInitialized ( &(metricParams->segmentList) ), XLAL_EINVAL, "Passed un-initialzied segment list 'metricParams->segmentList'\n");
  UINT4 Nseg = metricParams->segmentList.length;
  REAL8 sumTseg = 0;
  for ( UINT4 k=0; k < Nseg; k ++ )
    {
      LIGOTimeGPS *startTime_k = &(metricParams->segmentList.segs[k].start);
      LIGOTimeGPS *endTime_k   = &(metricParams->segmentList.segs[k].end);
      sumTseg += XLALGPSDiff( endTime_k, startTime_k );
    }
  REAL8 avgTseg = sumTseg / Nseg;

  // Compute naturalization transform
  gsl_matrix *GAMAT( transform, g_ij->size1, g_ij->size2 );
  gsl_matrix_set_zero( transform );
  for (size_t i = 0; i < g_ij->size1; ++i) {
    const DopplerCoordinateID coordID = metricParams->coordSys.coordIDs[i];
    const double Freq = metricParams->signalParams.Doppler.fkdot[0];
    const double T = avgTseg;
    double scale = 0;
    switch (coordID) {
    case DOPPLERCOORD_NONE:
      scale = 1;
      break;

    case DOPPLERCOORD_FREQ:
    case DOPPLERCOORD_GC_NU0:
      scale = LAL_TWOPI * LAL_FACT_INV[1] * (0.5 * T);
      break;

    case DOPPLERCOORD_F1DOT:
    case DOPPLERCOORD_GC_NU1:
      scale = LAL_TWOPI * LAL_FACT_INV[2] * POW2(0.5 * T);
      break;

    case DOPPLERCOORD_F2DOT:
    case DOPPLERCOORD_GC_NU2:
      scale = LAL_TWOPI * LAL_FACT_INV[3] * POW3(0.5 * T);
      break;

    case DOPPLERCOORD_F3DOT:
    case DOPPLERCOORD_GC_NU3:
      scale = LAL_TWOPI * LAL_FACT_INV[4] * POW4(0.5 * T);
      break;

    case DOPPLERCOORD_ALPHA:
    case DOPPLERCOORD_DELTA:
    case DOPPLERCOORD_N2X_EQU:
    case DOPPLERCOORD_N2Y_EQU:
    case DOPPLERCOORD_N2X_ECL:
    case DOPPLERCOORD_N2Y_ECL:
    case DOPPLERCOORD_N3X_EQU:
    case DOPPLERCOORD_N3Y_EQU:
    case DOPPLERCOORD_N3Z_EQU:
    case DOPPLERCOORD_N3X_ECL:
    case DOPPLERCOORD_N3Y_ECL:
    case DOPPLERCOORD_N3Z_ECL:
    case DOPPLERCOORD_N3SX_EQU:
    case DOPPLERCOORD_N3SY_EQU:
    case DOPPLERCOORD_N3OX_ECL:
    case DOPPLERCOORD_N3OY_ECL:
      {
        const REAL8 rEarth_c = (LAL_REARTH_SI / LAL_C_SI);
        REAL8 cosdD = cos ( metricParams->multiIFO.sites[0].frDetector.vertexLatitudeRadians );
        scale = LAL_TWOPI * Freq * rEarth_c * cosdD;
      }
      break;

    default:
      scale = 1;
      break;
    } // switch(coordID)
    gsl_matrix_set( transform, i, i, 1.0 / scale );
  }

  // Apply transform
  if ( p_gpr_ij != NULL ) {
    XLAL_CHECK( XLALTransformMetric( p_gpr_ij, transform, g_ij ) == XLAL_SUCCESS, XLAL_EFUNC );
  }

  // Return transform
  if ( p_transform != NULL ) {
    gsl_matrix_memcpy( *p_transform, transform );
  }

  // Cleanup
  GFMAT( transform );

  return XLAL_SUCCESS;

} // XLALNaturalizeMetric()
Exemple #2
0
int
XLALOutputDopplerMetric ( FILE *fp, const DopplerPhaseMetric *Pmetric, const DopplerFstatMetric *Fmetric, const ResultHistory_t *history )
{
  UINT4 i;
  REAL8 A, B, C, D;

  // ----- input sanity checks
  XLAL_CHECK ( fp != NULL, XLAL_EFAULT );
  XLAL_CHECK ( Pmetric != NULL || Fmetric != NULL, XLAL_EFAULT );
  const DopplerMetricParams *meta = (Pmetric != NULL) ? &(Pmetric->meta) : &(Fmetric->meta);
  XLAL_CHECK ( XLALSegListIsInitialized ( &(meta->segmentList) ), XLAL_EINVAL, "Got un-initialized segment list in 'metric->meta.segmentList'\n" );
  UINT4 Nseg = meta->segmentList.length;
  XLAL_CHECK ( Nseg >= 1, XLAL_EDOM, "Got invalid zero-length segment list 'metric->meta.segmentList'\n" );

  /* useful shortcuts */
  const PulsarDopplerParams *doppler = &(meta->signalParams.Doppler);
  const PulsarAmplitudeParams *Amp = &(meta->signalParams.Amp);

  /* output history info */
  if ( history )
    {
      if ( history->app_name ) fprintf (fp, "%%%% app_name: %s\n", history->app_name );
      if ( history->cmdline) fprintf (fp, "%%%% commandline: %s\n", history->cmdline );
      if ( history->VCSInfoString ) fprintf (fp, "%%%% Code Version: %s\n", history->VCSInfoString );
    }

  fprintf ( fp, "DopplerCoordinates = { " );
  for ( i=0; i < meta->coordSys.dim; i ++ )
    {
      if ( i > 0 ) fprintf ( fp, ", " );
      fprintf ( fp, "\"%s\"", XLALDopplerCoordinateName(meta->coordSys.coordIDs[i]));
    }
  fprintf ( fp, "};\n");

  { /* output projection info */
    const char *pname;
    if ( meta->projectCoord < 0 )
      pname = "None";
    else
      pname = XLALDopplerCoordinateName ( meta->coordSys.coordIDs[meta->projectCoord] );

    fprintf ( fp, "%%%% Projection onto subspace orthogonal to coordinate: '%s'\n", pname);
  }

  fprintf ( fp, "%%%% DetectorMotionType = '%s'\n", XLALDetectorMotionName(meta->detMotionType) );
  fprintf ( fp, "h0 = %g;\ncosi = %g;\npsi = %g;\nphi0 = %g;\n", Amp->h0, Amp->cosi, Amp->psi, Amp->phi0 );
  fprintf ( fp, "%%%% DopplerPoint = {\n");
  fprintf ( fp, "refTime = %.1f;\n", XLALGPSGetREAL8 ( &doppler->refTime ) );
  fprintf ( fp, "Alpha   = %f;\nDelta = %f;\n", doppler->Alpha, doppler->Delta );
  fprintf ( fp, "fkdot   = [%f, %g, %g, %g ];\n", doppler->fkdot[0], doppler->fkdot[1], doppler->fkdot[2], doppler->fkdot[3] );
  if ( doppler->asini > 0 )
    {
      fprintf ( fp, "%%%% 	   orbit = { \n");
      fprintf ( fp, "%%%% 		tp = {%d, %d}\n", doppler->tp.gpsSeconds, doppler->tp.gpsNanoSeconds );
      fprintf ( fp, "%%%% 		argp  = %g\n", doppler->argp );
      fprintf ( fp, "%%%% 		asini = %g\n", doppler->asini );
      fprintf ( fp, "%%%% 		ecc = %g\n", doppler->ecc );
      fprintf ( fp, "%%%% 		period = %g\n", doppler->period );
      fprintf ( fp, "%%%% 	   }\n");
    } /* if doppler->orbit */
  fprintf ( fp, "%%%% }\n");

  LIGOTimeGPS *tStart = &(meta->segmentList.segs[0].start);
  LIGOTimeGPS *tEnd   = &(meta->segmentList.segs[Nseg-1].end);
  REAL8 Tspan = XLALGPSDiff ( tEnd, tStart );
  fprintf ( fp, "startTime = %.1f;\n", XLALGPSGetREAL8 ( tStart ) );
  fprintf ( fp, "Tspan     = %.1f;\n", Tspan );
  fprintf ( fp, "Nseg      = %d;\n", Nseg );
  fprintf ( fp, "detectors = {");
  for ( i=0; i < meta->multiIFO.length; i ++ )
    {
      if ( i > 0 ) fprintf ( fp, ", ");
      fprintf ( fp, "\"%s\"", meta->multiIFO.sites[i].frDetector.name );
    }
  fprintf ( fp, "};\n");
  fprintf ( fp, "detectorWeights = [");
  for ( i=0; i < meta->multiNoiseFloor.length; i ++ )
    {
      if ( i > 0 ) fprintf ( fp, ", ");
      fprintf ( fp, "%f", meta->multiNoiseFloor.sqrtSn[i] );
    }
  fprintf ( fp, "];\n");

  /* ----- output phase metric ---------- */
  if ( Pmetric != NULL )
    {
      fprintf ( fp, "\ng_ij = \\\n" ); XLALfprintfGSLmatrix ( fp, METRIC_FORMAT,  Pmetric->g_ij );
      fprintf ( fp, "maxrelerr_gPh = %.2e;\n", Pmetric->maxrelerr );

      gsl_matrix *gN_ij = NULL;
      if ( XLALNaturalizeMetric ( &gN_ij, NULL, Pmetric->g_ij, meta ) != XLAL_SUCCESS ) {
        XLALPrintError ("%s: something failed Naturalizing phase metric g_ij!\n", __func__ );
        XLAL_ERROR ( XLAL_EFUNC );
      }
      fprintf ( fp, "\ngN_ij = \\\n" ); XLALfprintfGSLmatrix ( fp, METRIC_FORMAT,  gN_ij );
      gsl_matrix_free ( gN_ij );

      gsl_matrix *gDN_ij = NULL;
      if ( XLALDiagNormalizeMetric ( &gDN_ij, NULL, Pmetric->g_ij ) != XLAL_SUCCESS ) {
        XLALPrintError ("%s: something failed NormDiagonalizing phase metric g_ij!\n", __func__ );
        XLAL_ERROR ( XLAL_EFUNC );
      }
      fprintf ( fp, "\ngDN_ij = \\\n" ); XLALfprintfGSLmatrix ( fp, METRIC_FORMAT,  gDN_ij );
      gsl_matrix_free ( gDN_ij );
    }

  /* ----- output F-metric (and related matrices ---------- */
  if ( Fmetric != NULL )
    {
      fprintf ( fp, "\ngF_ij = \\\n" );   XLALfprintfGSLmatrix ( fp, METRIC_FORMAT,  Fmetric->gF_ij );
      fprintf ( fp, "\ngFav_ij = \\\n" ); XLALfprintfGSLmatrix ( fp, METRIC_FORMAT,  Fmetric->gFav_ij );
      fprintf ( fp, "\nm1_ij = \\\n" );   XLALfprintfGSLmatrix ( fp, METRIC_FORMAT,  Fmetric->m1_ij );
      fprintf ( fp, "\nm2_ij = \\\n" );   XLALfprintfGSLmatrix ( fp, METRIC_FORMAT,  Fmetric->m2_ij );
      fprintf ( fp, "\nm3_ij = \\\n" );   XLALfprintfGSLmatrix ( fp, METRIC_FORMAT,  Fmetric->m3_ij );
      fprintf ( fp, "maxrelerr_gF = %.2e;\n", Fmetric->maxrelerr );
    }

  /*  ----- output Fisher matrix ---------- */
  if ( Fmetric != NULL && Fmetric->Fisher_ab != NULL )
    {
      A = gsl_matrix_get ( Fmetric->Fisher_ab, 0, 0 );
      B = gsl_matrix_get ( Fmetric->Fisher_ab, 1, 1 );
      C = gsl_matrix_get ( Fmetric->Fisher_ab, 0, 1 );

      D = A * B - C * C;

      fprintf ( fp, "\nA = %.16g;\nB = %.16g;\nC = %.16g;\nD = %.16g;\n", A, B, C, D );
      fprintf ( fp, "\nrho2 = %.16g;\n", Fmetric->rho2 );

      fprintf (fp, "\nFisher_ab = \\\n" ); XLALfprintfGSLmatrix ( fp, METRIC_FORMAT,  Fmetric->Fisher_ab );
    }

  // ---------- output segment list at the end, as this can potentially become quite long and distracting
  char *seglist_octave;
  XLAL_CHECK ( (seglist_octave = XLALSegList2String ( &(meta->segmentList) )) != NULL, XLAL_EFUNC, "XLALSegList2String() with xlalErrno = %d\n", xlalErrno );
  fprintf ( fp, "\n\nsegmentList = %s;\n", seglist_octave );
  XLALFree ( seglist_octave );


  return XLAL_SUCCESS;

} /* XLALOutputDopplerMetric() */