Ejemplo n.º 1
0
int
tend_estimMain(int argc, char **argv, char *me, hestParm *hparm) {
  int pret;
  hestOpt *hopt = NULL;
  char *perr, *err;
  airArray *mop;

  Nrrd **nin, *nin4d, *nbmat, *nterr, *nB0, *nout;
  char *outS, *terrS, *bmatS, *eb0S;
  float soft, scale, sigma;
  int dwiax, EE, knownB0, oldstuff, estmeth, verbose, fixneg;
  unsigned int ninLen, axmap[4], wlsi, *skip, skipNum, skipIdx;
  double valueMin, thresh;

  Nrrd *ngradKVP=NULL, *nbmatKVP=NULL;
  double bKVP, bval;

  tenEstimateContext *tec;

  hestOptAdd(&hopt, "old", NULL, airTypeInt, 0, 0, &oldstuff, NULL,
             "instead of the new tenEstimateContext code, use "
             "the old tenEstimateLinear code");
  hestOptAdd(&hopt, "sigma", "sigma", airTypeFloat, 1, 1, &sigma, "nan",
             "Rician noise parameter");
  hestOptAdd(&hopt, "v", "verbose", airTypeInt, 1, 1, &verbose, "0",
             "verbosity level");
  hestOptAdd(&hopt, "est", "estimate method", airTypeEnum, 1, 1, &estmeth,
             "lls",
             "estimation method to use. \"lls\": linear-least squares",
             NULL, tenEstimate1Method);
  hestOptAdd(&hopt, "wlsi", "WLS iters", airTypeUInt, 1, 1, &wlsi, "1",
             "when using weighted-least-squares (\"-est wls\"), how "
             "many iterations to do after the initial weighted fit.");
  hestOptAdd(&hopt, "fixneg", NULL, airTypeInt, 0, 0, &fixneg, NULL,
             "after estimating the tensor, ensure that there are no negative "
             "eigenvalues by adding (to all eigenvalues) the amount by which "
             "the smallest is negative (corresponding to increasing the "
             "non-DWI image value).");
  hestOptAdd(&hopt, "ee", "filename", airTypeString, 1, 1, &terrS, "",
             "Giving a filename here allows you to save out the tensor "
             "estimation error: a value which measures how much error there "
             "is between the tensor model and the given diffusion weighted "
             "measurements for each sample.  By default, no such error "
             "calculation is saved.");
  hestOptAdd(&hopt, "eb", "filename", airTypeString, 1, 1, &eb0S, "",
             "In those cases where there is no B=0 reference image given "
             "(\"-knownB0 false\"), "
             "giving a filename here allows you to save out the B=0 image "
             "which is estimated from the data.  By default, this image value "
             "is estimated but not saved.");
  hestOptAdd(&hopt, "t", "thresh", airTypeDouble, 1, 1, &thresh, "nan",
             "value at which to threshold the mean DWI value per pixel "
             "in order to generate the \"confidence\" mask.  By default, "
             "the threshold value is calculated automatically, based on "
             "histogram analysis.");
  hestOptAdd(&hopt, "soft", "soft", airTypeFloat, 1, 1, &soft, "0",
             "how fuzzy the confidence boundary should be.  By default, "
             "confidence boundary is perfectly sharp");
  hestOptAdd(&hopt, "scale", "scale", airTypeFloat, 1, 1, &scale, "1",
             "After estimating the tensor, scale all of its elements "
             "(but not the confidence value) by this amount.  Can help with "
             "downstream numerical precision if values are very large "
             "or small.");
  hestOptAdd(&hopt, "mv", "min val", airTypeDouble, 1, 1, &valueMin, "1.0",
             "minimum plausible value (especially important for linear "
             "least squares estimation)");
  hestOptAdd(&hopt, "B", "B-list", airTypeString, 1, 1, &bmatS, NULL,
             "6-by-N list of B-matrices characterizing "
             "the diffusion weighting for each "
             "image.  \"tend bmat\" is one source for such a matrix; see "
             "its usage info for specifics on how the coefficients of "
             "the B-matrix are ordered. "
             "An unadorned plain text file is a great way to "
             "specify the B-matrix.\n  **OR**\n "
             "Can say just \"-B kvp\" to try to learn B matrices from "
             "key/value pair information in input images.");
  hestOptAdd(&hopt, "b", "b", airTypeDouble, 1, 1, &bval, "nan",
             "\"b\" diffusion-weighting factor (units of sec/mm^2)");
  hestOptAdd(&hopt, "knownB0", "bool", airTypeBool, 1, 1, &knownB0, NULL,
             "Determines of the B=0 non-diffusion-weighted reference image "
             "is known, or if it has to be estimated along with the tensor "
             "elements.\n "
             "\b\bo if \"true\": in the given list of diffusion gradients or "
             "B-matrices, there are one or more with zero norm, which are "
             "simply averaged to find the B=0 reference image value\n "
             "\b\bo if \"false\": there may or may not be diffusion-weighted "
             "images among the input; the B=0 image value is going to be "
             "estimated along with the diffusion model");
  hestOptAdd(&hopt, "i", "dwi0 dwi1", airTypeOther, 1, -1, &nin, "-",
             "all the diffusion-weighted images (DWIs), as seperate 3D nrrds, "
             "**OR**: One 4D nrrd of all DWIs stacked along axis 0",
             &ninLen, NULL, nrrdHestNrrd);
  hestOptAdd(&hopt, "o", "nout", airTypeString, 1, 1, &outS, "-",
             "output tensor volume");

  mop = airMopNew();
  airMopAdd(mop, hopt, (airMopper)hestOptFree, airMopAlways);
  USAGE(_tend_estimInfoL);
  JUSTPARSE();
  airMopAdd(mop, hopt, (airMopper)hestParseFree, airMopAlways);

  nout = nrrdNew();
  airMopAdd(mop, nout, (airMopper)nrrdNuke, airMopAlways);
  nbmat = nrrdNew();
  airMopAdd(mop, nbmat, (airMopper)nrrdNuke, airMopAlways);

  /* figure out B-matrix */
  if (strcmp("kvp", airToLower(bmatS))) {
    /* its NOT coming from key/value pairs */
    if (!AIR_EXISTS(bval)) {
      fprintf(stderr, "%s: need to specify scalar b-value\n", me);
      airMopError(mop); return 1;
    }
    if (nrrdLoad(nbmat, bmatS, NULL)) {
      airMopAdd(mop, err=biffGetDone(NRRD), airFree, airMopAlways);
      fprintf(stderr, "%s: trouble loading B-matrix:\n%s\n", me, err);
      airMopError(mop); return 1;
    }
    nin4d = nin[0];
    skip = NULL;
    skipNum = 0;
  } else {
    /* it IS coming from key/value pairs */
    if (1 != ninLen) {
      fprintf(stderr, "%s: require a single 4-D DWI volume for "
              "key/value pair based calculation of B-matrix\n", me);
      airMopError(mop); return 1;
    }
    if (oldstuff) {
      if (knownB0) {
        fprintf(stderr, "%s: sorry, key/value-based DWI info not compatible "
                "with older implementation of knownB0\n", me);
        airMopError(mop); return 1;
      }
    }
    if (tenDWMRIKeyValueParse(&ngradKVP, &nbmatKVP, &bKVP,
                              &skip, &skipNum, nin[0])) {
      airMopAdd(mop, err=biffGetDone(TEN), airFree, airMopAlways);
      fprintf(stderr, "%s: trouble parsing DWI info:\n%s\n", me, err);
      airMopError(mop); return 1;
    }
    if (AIR_EXISTS(bval)) {
      fprintf(stderr, "%s: WARNING: key/value pair derived b-value %g "
              "over-riding %g from command-line", me, bKVP, bval);
    }
    bval = bKVP;
    if (ngradKVP) {
      airMopAdd(mop, ngradKVP, (airMopper)nrrdNuke, airMopAlways);
      if (tenBMatrixCalc(nbmat, ngradKVP)) {
        airMopAdd(mop, err=biffGetDone(TEN), airFree, airMopAlways);
        fprintf(stderr, "%s: trouble finding B-matrix:\n%s\n", me, err);
        airMopError(mop); return 1;
      }
    } else {
      airMopAdd(mop, nbmatKVP, (airMopper)nrrdNuke, airMopAlways);
      if (nrrdConvert(nbmat, nbmatKVP, nrrdTypeDouble)) {
        airMopAdd(mop, err=biffGetDone(NRRD), airFree, airMopAlways);
        fprintf(stderr, "%s: trouble converting B-matrix:\n%s\n", me, err);
        airMopError(mop); return 1;
      }
    }
    /* this will work because of the impositions of tenDWMRIKeyValueParse */
    dwiax = ((nrrdKindList == nin[0]->axis[0].kind ||
              nrrdKindVector == nin[0]->axis[0].kind)
             ? 0
             : ((nrrdKindList == nin[0]->axis[1].kind ||
                 nrrdKindVector == nin[0]->axis[1].kind)
                ? 1
                : ((nrrdKindList == nin[0]->axis[2].kind ||
                    nrrdKindVector == nin[0]->axis[2].kind)
                   ? 2
                   : 3)));
    if (0 == dwiax) {
      nin4d = nin[0];
    } else {
      axmap[0] = dwiax;
      axmap[1] = 1 > dwiax ? 1 : 0;
      axmap[2] = 2 > dwiax ? 2 : 1;
      axmap[3] = 3 > dwiax ? 3 : 2;
      nin4d = nrrdNew();
      airMopAdd(mop, nin4d, (airMopper)nrrdNuke, airMopAlways);
      if (nrrdAxesPermute(nin4d, nin[0], axmap)) {
        airMopAdd(mop, err=biffGetDone(NRRD), airFree, airMopAlways);
        fprintf(stderr, "%s: trouble creating DWI volume:\n%s\n", me, err);
        airMopError(mop); return 1;
      }
    }
  }

  nterr = NULL;
  nB0 = NULL;
  if (!oldstuff) {
    if (1 != ninLen) {
      fprintf(stderr, "%s: sorry, currently need single 4D volume "
              "for new implementation\n", me);
      airMopError(mop); return 1;
    }
    if (!AIR_EXISTS(thresh)) {
      if (tend_estimThresholdFind(&thresh, nbmat, nin4d)) {
        airMopAdd(mop, err=biffGetDone(TEN), airFree, airMopAlways);
        fprintf(stderr, "%s: trouble finding threshold:\n%s\n", me, err);
        airMopError(mop); return 1;
      }
      /* HACK to lower threshold a titch */
      thresh *= 0.93;
      fprintf(stderr, "%s: using mean DWI threshold %g\n", me, thresh);
    }
    tec = tenEstimateContextNew();
    tec->progress = AIR_TRUE;
    airMopAdd(mop, tec, (airMopper)tenEstimateContextNix, airMopAlways);
    EE = 0;
    if (!EE) tenEstimateVerboseSet(tec, verbose);
    if (!EE) tenEstimateNegEvalShiftSet(tec, fixneg);
    if (!EE) EE |= tenEstimateMethodSet(tec, estmeth);
    if (!EE) EE |= tenEstimateBMatricesSet(tec, nbmat, bval, !knownB0);
    if (!EE) EE |= tenEstimateValueMinSet(tec, valueMin);
    for (skipIdx=0; skipIdx<skipNum; skipIdx++) {
      /* fprintf(stderr, "%s: skipping %u\n", me, skip[skipIdx]); */
      if (!EE) EE |= tenEstimateSkipSet(tec, skip[skipIdx], AIR_TRUE);
    }
    switch(estmeth) {
    case tenEstimate1MethodLLS:
      if (airStrlen(terrS)) {
        tec->recordErrorLogDwi = AIR_TRUE;
        /* tec->recordErrorDwi = AIR_TRUE; */
      }
      break;
    case tenEstimate1MethodNLS:
      if (airStrlen(terrS)) {
        tec->recordErrorDwi = AIR_TRUE;
      }
      break;
    case tenEstimate1MethodWLS:
      if (!EE) tec->WLSIterNum = wlsi;
      if (airStrlen(terrS)) {
        tec->recordErrorDwi = AIR_TRUE;
      }
      break;
    case tenEstimate1MethodMLE:
      if (!(AIR_EXISTS(sigma) && sigma > 0.0)) {
        fprintf(stderr, "%s: can't do %s w/out sigma > 0 (not %g)\n",
                me, airEnumStr(tenEstimate1Method, tenEstimate1MethodMLE),
                sigma);
        airMopError(mop); return 1;
      }
      if (!EE) EE |= tenEstimateSigmaSet(tec, sigma);
      if (airStrlen(terrS)) {
        tec->recordLikelihoodDwi = AIR_TRUE;
      }
      break;
    }
    if (!EE) EE |= tenEstimateThresholdSet(tec, thresh, soft);
    if (!EE) EE |= tenEstimateUpdate(tec);
    if (EE) {
      airMopAdd(mop, err=biffGetDone(TEN), airFree, airMopAlways);
      fprintf(stderr, "%s: trouble setting up estimation:\n%s\n", me, err);
      airMopError(mop); return 1;
    }
    if (tenEstimate1TensorVolume4D(tec, nout, &nB0,
                                   airStrlen(terrS) 
                                   ? &nterr 
                                   : NULL, 
                                   nin4d, nrrdTypeFloat)) {
      airMopAdd(mop, err=biffGetDone(TEN), airFree, airMopAlways);
      fprintf(stderr, "%s: trouble doing estimation:\n%s\n", me, err);
      airMopError(mop); return 1;
    }
    if (airStrlen(terrS)) {
      airMopAdd(mop, nterr, (airMopper)nrrdNuke, airMopAlways);
    }
  } else {
    EE = 0;
    if (1 == ninLen) {
      EE = tenEstimateLinear4D(nout, airStrlen(terrS) ? &nterr : NULL, &nB0,
                               nin4d, nbmat, knownB0, thresh, soft, bval);
    } else {
      EE = tenEstimateLinear3D(nout, airStrlen(terrS) ? &nterr : NULL, &nB0,
                               (const Nrrd**)nin, ninLen, nbmat,
                               knownB0, thresh, soft, bval);
    }
    if (EE) {
      airMopAdd(mop, err=biffGetDone(TEN), airFree, airMopAlways);
      fprintf(stderr, "%s: trouble making tensor volume:\n%s\n", me, err);
      airMopError(mop); return 1;
    }
  }
  if (nterr) {
    /* it was allocated by tenEstimate*, we have to clean it up */
    airMopAdd(mop, nterr, (airMopper)nrrdNuke, airMopAlways);
  }
  if (nB0) {
    /* it was allocated by tenEstimate*, we have to clean it up */
    airMopAdd(mop, nB0, (airMopper)nrrdNuke, airMopAlways);
  }
  if (1 != scale) {
    if (tenSizeScale(nout, nout, scale)) {
      airMopAdd(mop, err=biffGetDone(TEN), airFree, airMopAlways);
      fprintf(stderr, "%s: trouble doing scaling:\n%s\n", me, err);
      airMopError(mop); return 1;
    }
  }
  if (nterr) {
    if (nrrdSave(terrS, nterr, NULL)) {
      airMopAdd(mop, err=biffGetDone(NRRD), airFree, airMopAlways);
      fprintf(stderr, "%s: trouble writing error image:\n%s\n", me, err);
      airMopError(mop); return 1;
    }
  }
  if (!knownB0 && airStrlen(eb0S)) {
    if (nrrdSave(eb0S, nB0, NULL)) {
      airMopAdd(mop, err=biffGetDone(NRRD), airFree, airMopAlways);
      fprintf(stderr, "%s: trouble writing estimated B=0 image:\n%s\n",
              me, err);
      airMopError(mop); return 1;
    }
  }

  if (nrrdSave(outS, nout, NULL)) {
    airMopAdd(mop, err=biffGetDone(NRRD), airFree, airMopAlways);
    fprintf(stderr, "%s: trouble writing:\n%s\n", me, err);
    airMopError(mop); return 1;
  }

  airMopOkay(mop);
  return 0;
}
Ejemplo n.º 2
0
void
_tenDwiGageAnswer(gageContext *ctx, gagePerVolume *pvl) {
  char me[]="_tenDwiGageAnswer";
  unsigned int dwiIdx;
  tenDwiGageKindData *kindData;
  tenDwiGagePvlData *pvlData;
  double *dwiAll, dwiMean=0, tentmp[7];

  kindData = AIR_CAST(tenDwiGageKindData *, pvl->kind->data);
  pvlData = AIR_CAST(tenDwiGagePvlData *, pvl->data);

  dwiAll = pvl->directAnswer[tenDwiGageAll];
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageAll)) {
    /* done if doV */
    if (ctx->verbose) {
      for (dwiIdx=0; dwiIdx<pvl->kind->valLen; dwiIdx++) {
        fprintf(stderr, "%s(%d+%g,%d+%g,%d+%g): dwi[%u] = %g\n", me,
                ctx->point.xi, ctx->point.xf,
                ctx->point.yi, ctx->point.yf,
                ctx->point.zi, ctx->point.zf,
                dwiIdx, dwiAll[dwiIdx]);
      }
      fprintf(stderr, "%s: type(ngrad) = %d = %s\n", me,
              kindData->ngrad->type,
              airEnumStr(nrrdType, kindData->ngrad->type));
    }
  }

  /*
    if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageB0)) {
    if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageJustDWI)) {
    done if doV
    }
  */
  /* HEY this isn't valid for multiple b-values */
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageADC)) {
    double logdwi, logb0;
    logb0 = log(AIR_MAX(kindData->valueMin,
                        pvl->directAnswer[tenDwiGageB0][0]));
    for (dwiIdx=1; dwiIdx<pvl->kind->valLen; dwiIdx++) {
      logdwi = log(AIR_MAX(kindData->valueMin,
                           pvl->directAnswer[tenDwiGageJustDWI][dwiIdx-1]));
      pvl->directAnswer[tenDwiGageADC][dwiIdx-1]
        = (logb0 - logdwi)/kindData->bval;
    }
  }
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageMeanDWIValue)) {
    dwiMean = 0;
    for (dwiIdx=1; dwiIdx<pvl->kind->valLen; dwiIdx++) {
      dwiMean += dwiAll[dwiIdx];
    }
    dwiMean /= pvl->kind->valLen;
    pvl->directAnswer[tenDwiGageMeanDWIValue][0] = dwiMean;
  }

  /* note: the gage interface to tenEstimate functionality 
     allows you exactly one kind of tensor estimation (per kind),
     so the function call to do the estimation is actually
     repeated over and over again; the copy into the answer
     buffer is what changes... */
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensorLLS)) {
    tenEstimate1TensorSingle_d(pvlData->tec1, tentmp, dwiAll);
    TEN_T_COPY(pvl->directAnswer[tenDwiGageTensorLLS], tentmp);
  }
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensorLLSError)) {
    pvl->directAnswer[tenDwiGageTensorLLSError][0] = pvlData->tec1->errorDwi;
  }  
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensorLLSErrorLog)) {
    pvl->directAnswer[tenDwiGageTensorLLSErrorLog][0] 
      = pvlData->tec1->errorLogDwi;
  }  
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensorWLS)) {
    tenEstimate1TensorSingle_d(pvlData->tec1, tentmp, dwiAll);
    TEN_T_COPY(pvl->directAnswer[tenDwiGageTensorWLS], tentmp);
  }
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensorNLS)) {
    tenEstimate1TensorSingle_d(pvlData->tec1, tentmp, dwiAll);
    TEN_T_COPY(pvl->directAnswer[tenDwiGageTensorNLS], tentmp);
  }
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensorMLE)) {
    tenEstimate1TensorSingle_d(pvlData->tec1, tentmp, dwiAll);
    TEN_T_COPY(pvl->directAnswer[tenDwiGageTensorMLE], tentmp);
  }
  /* HEY: have to implement all the different kinds of errors */

  /* BEGIN sneakiness ........ */
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensor)) {
    gageItemEntry *item;
    item = pvl->kind->table + tenDwiGageTensor;
    TEN_T_COPY(pvl->directAnswer[tenDwiGageTensor],
               pvl->directAnswer[item->prereq[0]]);
  }
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensorError)) {
    gageItemEntry *item;
    item = pvl->kind->table + tenDwiGageTensorError;
    pvl->directAnswer[tenDwiGageTensorError][0]
      = pvl->directAnswer[item->prereq[0]][0];
  }
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensorErrorLog)) {
    gageItemEntry *item;
    item = pvl->kind->table + tenDwiGageTensorErrorLog;
    pvl->directAnswer[tenDwiGageTensorErrorLog][0]
      = pvl->directAnswer[item->prereq[0]][0];
  }
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensorLikelihood)) {
    gageItemEntry *item;
    item = pvl->kind->table + tenDwiGageTensorLikelihood;
    pvl->directAnswer[tenDwiGageTensorLikelihood][0]
      = pvl->directAnswer[item->prereq[0]][0];
  }
  /* END sneakiness ........ */

  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageFA)) {
    pvl->directAnswer[tenDwiGageFA][0]
      = pvl->directAnswer[tenDwiGageTensor][0]
      * tenAnisoTen_d(pvl->directAnswer[tenDwiGageTensor],
                      tenAniso_FA);
  }

  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGageTensorAllDWIError)) {
    const double *grads;
    int gradcount;
    double *ten, d;
    int i;
    
    /* HEY: should switch to tenEstimate-based DWI simulation */
    ten = pvl->directAnswer[tenDwiGageTensor];
    gradcount = pvl->kind->valLen -1; /* Dont count b0 */
    grads = ((const double*) kindData->ngrad->data) +3; /* Ignore b0 grad */
    for( i=0; i < gradcount; i++ ) {
      d = dwiAll[0]*exp(- pvlData->tec1->bValue 
                        * TEN_T3V_CONTR(ten, grads + 3*i));
      pvl->directAnswer[tenDwiGageTensorAllDWIError][i] = dwiAll[i+1] - d;
    }
  }

  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGage2TensorQSeg)) {
    const double *grads;
    int gradcount;
    double *twoten;
    unsigned int valIdx, E;
    
    twoten = pvl->directAnswer[tenDwiGage2TensorQSeg];
    
    gradcount = pvl->kind->valLen -1; /* Dont count b0 */
    grads = ((const double*) kindData->ngrad->data) +3; /* Ignore b0 grad */
    if (dwiAll[0] != 0) { /*  S0 = 0 */
      _tenQball(pvlData->tec2->bValue, gradcount, dwiAll, grads,
                pvlData->qvals);
      _tenQvals2points(gradcount, pvlData->qvals, grads, pvlData->qpoints);
      _tenSegsamp2(gradcount, pvlData->qvals, grads, pvlData->qpoints,
                   pvlData->wght + 1, pvlData->dists );
    } else {
      /* stupid; should really return right here since data is garbage */
      for (valIdx=1; valIdx < AIR_CAST(unsigned int, gradcount+1); valIdx++) {
        pvlData->wght[valIdx] = valIdx % 2;
      }
    }
    
    E = 0;
    for (valIdx=1; valIdx<pvl->kind->valLen; valIdx++) {
      if (!E) E |= tenEstimateSkipSet(pvlData->tec2, valIdx,
                                      pvlData->wght[valIdx]);
    }
    if (!E) E |= tenEstimateUpdate(pvlData->tec2);
    if (!E) E |= tenEstimate1TensorSingle_d(pvlData->tec2,
                                            twoten + 0, dwiAll);
    for (valIdx=1; valIdx<pvl->kind->valLen; valIdx++) {
      if (!E) E |= tenEstimateSkipSet(pvlData->tec2, valIdx,
                                      1 - pvlData->wght[valIdx]);
    }
    if (!E) E |= tenEstimateUpdate(pvlData->tec2);
    if (!E) E |= tenEstimate1TensorSingle_d(pvlData->tec2,
                                            twoten + 7, dwiAll);
    if (E) {
      fprintf(stderr, "!%s: (trouble) %s\n", me, biffGetDone(TEN));
    }
    
    /* hack: confidence for two-tensor fit */
    twoten[0] = (twoten[0] + twoten[7])/2;
    twoten[7] = 0.5; /* fraction that is the first tensor (initial value) */
    /* twoten[1 .. 6] = first tensor */
    /* twoten[8 .. 13] = second tensor */
    
    /* Compute fraction between tensors if not garbage in this voxel */
    if (twoten[0] > 0.5) {
      double exp0,exp1,d,e=0,g=0, a=0,b=0;
      int i;
      
      for( i=0; i < gradcount; i++ ) {
        exp0 = exp(-pvlData->tec2->bValue * TEN_T3V_CONTR(twoten + 0,
                                                          grads + 3*i));
        exp1 = exp(-pvlData->tec2->bValue * TEN_T3V_CONTR(twoten + 7,
                                                          grads + 3*i));
        
        d = dwiAll[i+1] / dwiAll[0];
        e = exp0 - exp1;
        g = d - exp1;
        
        a += .5*e*e;
        b += e*g;
      }
      
      twoten[7] = AIR_CLAMP(0, 0.5*(b/a), 1);
    }
  }
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGage2TensorQSegError)) {
    const double *grads;
    int gradcount;
    double *twoten, d;
    int i;
    
    /* HEY: should switch to tenEstimate-based DWI simulation */
    if (dwiAll[0] != 0) { /* S0 = 0 */
      twoten = pvl->directAnswer[tenDwiGage2TensorQSeg];
      gradcount = pvl->kind->valLen -1; /* Dont count b0 */
      grads = ((const double*) kindData->ngrad->data) +3; /* Ignore b0 grad */
      
      pvl->directAnswer[tenDwiGage2TensorQSegError][0] = 0;
      for( i=0; i < gradcount; i++ ) {
        d = twoten[7]*exp(-pvlData->tec2->bValue * TEN_T3V_CONTR(twoten + 0,
                                                                 grads + 3*i));
        d += (1 - twoten[7])*exp(-pvlData->tec2->bValue 
                                 *TEN_T3V_CONTR(twoten + 7, grads + 3*i));
        d = dwiAll[i+1]/dwiAll[0] - d;
        pvl->directAnswer[tenDwiGage2TensorQSegError][0] += d*d;
      }
      pvl->directAnswer[tenDwiGage2TensorQSegError][0] = 
        sqrt( pvl->directAnswer[tenDwiGage2TensorQSegError][0] );
    } else {
      /* HEY: COMPLETELY WRONG!! An error is not defined! */
      pvl->directAnswer[tenDwiGage2TensorQSegError][0] = 0;
    }
    /* printf("%f\n",pvl->directAnswer[tenDwiGage2TensorQSegError][0]); */
  }
  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGage2TensorQSegAndError)) {
    double *twoten, *err, *twotenerr;
    
    twoten = pvl->directAnswer[tenDwiGage2TensorQSeg];
    err = pvl->directAnswer[tenDwiGage2TensorQSegError];
    twotenerr = pvl->directAnswer[tenDwiGage2TensorQSegAndError];
    TEN_T_COPY(twotenerr + 0, twoten + 0);
    TEN_T_COPY(twotenerr + 7, twoten + 7);
    twotenerr[14] = err[0];
  }

  if (GAGE_QUERY_ITEM_TEST(pvl->query, tenDwiGage2TensorPeled)) {
#if TEEM_LEVMAR
#define PARAMS 4
    double *twoTen, Cp /* , residual, AICSingFit, AICTwoFit */;
    /* Vars for the NLLS */
    double guess[PARAMS], loBnd[PARAMS], upBnd[PARAMS],
      opts[LM_OPTS_SZ], *grad, *egrad, tenA[7], tenB[7],
      matA[9], matB[9], matTmp[9], rott[9];
    unsigned int gi;
    int lmret;
    
    /* Pointer to the location where the two tensor will be written */
    twoTen = pvl->directAnswer[tenDwiGage2TensorPeled];
    /* Estimate the DWI error, error is given as standard deviation */
    pvlData->tec2->recordErrorDwi = AIR_FALSE;
    /* Estimate the single tensor */
    tenEstimate1TensorSingle_d(pvlData->tec2, pvlData->ten1, dwiAll);
    /* Get the eigenValues and eigen vectors for this tensor */
    tenEigensolve_d(pvlData->ten1Eval, pvlData->ten1Evec, pvlData->ten1);
    /* Get westins Cp */
    Cp = tenAnisoEval_d(pvlData->ten1Eval, tenAniso_Cp1);
    
    /* Calculate the residual, need the variance to sqr it */
    /* residual = pvlData->tec2->errorDwi*pvlData->tec2->errorDwi; */
    /* Calculate the AIC for single tensor fit */
    /* AICSingFit = _tenComputeAIC(residual, pvlData->tec2->dwiNum, 6); */

    /* the CP-based test is gone; caller's responsibility */
      
    /* rotate DW gradients by inverse of eigenvector column matrix
       and place into pvlData->nten1EigenGrads (which has been
       allocated by _tenDwiGagePvlDataNew()) */
    grad = AIR_CAST(double *, kindData->ngrad->data);
    egrad = AIR_CAST(double *, pvlData->nten1EigenGrads->data);
    for (gi=0; gi<kindData->ngrad->axis[1].size; gi++) {
      /* yes, this is also transforming some zero-length (B0) gradients;
         that's harmless */
      ELL_3MV_MUL(egrad, pvlData->ten1Evec, grad);
      grad += 3;
      egrad += 3;
    }
    
    /* Lower and upper bounds for the NLLS routine */
    loBnd[0] = 0.0;
    loBnd[1] = 0.0;       
    loBnd[2] = -AIR_PI/2;
    loBnd[3] = -AIR_PI/2;
    upBnd[0] = pvlData->ten1Eval[0]*5;
    upBnd[1] = 1.0;
    upBnd[2] = AIR_PI/2;
    upBnd[3] = AIR_PI/2;
    /* Starting point for the NLLS */
    guess[0] = pvlData->ten1Eval[0];
    guess[1] = 0.5;

    guess[2] = AIR_PI/4;
    guess[3] = -AIR_PI/4;
    /*
    guess[2] = AIR_AFFINE(0, airDrandMT_r(pvlData->randState), 1,
                          AIR_PI/6, AIR_PI/3);
    guess[3] = AIR_AFFINE(0, airDrandMT_r(pvlData->randState), 1,
                          -AIR_PI/6, -AIR_PI/3);
    */
    /* Fill in the constraints for the LM optimization, 
       the threshold of error difference */
    opts[0] = pvlData->levmarTau;
    opts[1] = pvlData->levmarEps1;
    opts[2] = pvlData->levmarEps2;
    opts[3] = pvlData->levmarEps3;
    /* Very imp to set this opt, note that only forward
       differences are used to approx Jacobian */
    opts[4] = pvlData->levmarDelta;
    
    /* run NLLS, results are stored back into guess[] */
    pvlData->levmarUseFastExp = AIR_FALSE;
    lmret = dlevmar_bc_dif(_tenLevmarPeledCB, guess, pvlData->tec2->dwi,
                           PARAMS, pvlData->tec2->dwiNum, loBnd, upBnd,
                           pvlData->levmarMaxIter, opts,
                           pvlData->levmarInfo,
                           NULL, NULL, pvlData);
    if (-1 == lmret) {
      ctx->errNum = 1;
      sprintf(ctx->errStr, "%s: dlevmar_bc_dif() failed!", me);
    } else {
      /* Get the AIC for the two tensor fit, use the levmarinfo
         to get the residual */
      /*
        residual = pvlData->levmarInfo[1]/pvlData->tec2->dwiNum;
        AICTwoFit = _tenComputeAIC(residual, pvlData->tec2->dwiNum, 12);
      */
      /* Form the tensors using the estimated pp, returned in guess */
      _tenPeledRotate2D(tenA, guess[0], pvlData->ten1Eval[2], guess[2]);
      _tenPeledRotate2D(tenB, guess[0], pvlData->ten1Eval[2], guess[3]);
      TEN_T2M(matA, tenA);
      TEN_T2M(matB, tenB);
      
      ELL_3M_TRANSPOSE(rott, pvlData->ten1Evec);
      ELL_3M_MUL(matTmp, matA, pvlData->ten1Evec);
      ELL_3M_MUL(matA, rott, matTmp);
      ELL_3M_MUL(matTmp, matB, pvlData->ten1Evec);
      ELL_3M_MUL(matB, rott, matTmp);
      
      /* Copy two two tensors */
      /* guess[1] is population fraction of first tensor */
      if (guess[1] > 0.5) {
        twoTen[7] = guess[1];
        TEN_M2T(twoTen + 0, matA);
        TEN_M2T(twoTen + 7, matB);
      } else {
        twoTen[7] = 1 - guess[1];
        TEN_M2T(twoTen + 0, matB);
        TEN_M2T(twoTen + 7, matA);
      }
      twoTen[0] = 1;
    }
#undef PARAMS
#else
    double *twoTen;
    twoTen = pvl->directAnswer[tenDwiGage2TensorPeled];
    TEN_T_SET(twoTen + 0, AIR_NAN, AIR_NAN, AIR_NAN, AIR_NAN,
              AIR_NAN, AIR_NAN, AIR_NAN);
    TEN_T_SET(twoTen + 7, AIR_NAN, AIR_NAN, AIR_NAN, AIR_NAN,
              AIR_NAN, AIR_NAN, AIR_NAN);
    fprintf(stderr, "%s: sorry, not compiled with TEEM_LEVMAR\n", me);
#endif
  }