266 {
267 double *x = NULL, *y = NULL, *sy = NULL, *sx = NULL, *diff = NULL, xOffset,
268 xScaleFactor;
269 double *xOrig = NULL, *yOrig = NULL, *sxOrig, *syOrig, *sy0 = NULL;
270 long terms, normTerm, ySigmasValid;
271 int64_t i, j, points, pointsOrig;
272 long symmetry, chebyshev, autoOffset, copyParameters = 0;
273 double sigmas;
274 long sigmasMode, sparseInterval;
275 unsigned long flags;
276 double *coef, *coefSigma;
277 double chi, xLow, xHigh, rmsResidual;
278 char *xName, *yName, *xSigmaName, *ySigmaName;
279 char *input, *output, **coefUnits;
281 long isFit, iArg, modifySigmas;
282 long generateSigmas, verbose, ignoreSigmas;
283 long npages = 0, invalid = 0;
284 int32_t *order;
285 SCANNED_ARG *s_arg;
286 double xMin, xMax, revpowThreshold, revpowCompleteThres, goodEnoughChi;
287 long rangeFitOnly = 0;
288 double rms_average(double *d_x, int64_t d_n);
289 char *fitLabelFormat = "%g";
290 static char rpnSeqBuffer[SDDS_MAXLINE];
291 unsigned long pipeFlags, reviseOrders, majorOrderFlag;
293 short columnMajorOrder = -1;
294
295 sxOrig = syOrig = NULL;
296 rmsResidual = 0;
297
299 argc =
scanargs(&s_arg, argc, argv);
300 if (argc < 2 || argc > (3 + N_OPTIONS)) {
301 fprintf(stderr, "usage: %s%s%s\n", USAGE, additional_help1,
302 additional_help2);
303 exit(EXIT_FAILURE);
304 }
305
306 input = output = NULL;
307 xName = yName = xSigmaName = ySigmaName = NULL;
308 modifySigmas = reviseOrders = chebyshev = 0;
309 order = NULL;
310 symmetry = NO_SYMMETRY;
311 xMin = xMax = 0;
312 autoOffset = 0;
313 generateSigmas = 0;
314 sigmasMode = -1;
315 sigmas = 1;
316 sparseInterval = 1;
317 terms = 2;
318 verbose = ignoreSigmas = 0;
319 normTerm = -1;
320 xOffset = 0;
321 xScaleFactor = 1;
322 coefUnits = NULL;
325 pipeFlags = 0;
326 evalParameters.file = evalParameters.valuesFile = evalParameters.valuesColumn = NULL;
327 evalParameters.initialized = evalParameters.inputInitialized = 0;
328
329 for (iArg = 1; iArg < argc; iArg++) {
330 if (s_arg[iArg].arg_type == OPTION) {
331 switch (
match_string(s_arg[iArg].list[0], option, N_OPTIONS, 0)) {
332 case CLO_MAJOR_ORDER:
333 majorOrderFlag = 0;
334 s_arg[iArg].n_items--;
335 if (s_arg[iArg].n_items > 0 &&
337 &s_arg[iArg].n_items, 0, "row", -1, NULL, 0,
338 SDDS_ROW_MAJOR_ORDER, "column", -1, NULL, 0,
339 SDDS_COLUMN_MAJOR_ORDER, NULL)))
340 SDDS_Bomb(
"invalid -majorOrder syntax/values");
341 if (majorOrderFlag & SDDS_COLUMN_MAJOR_ORDER)
342 columnMajorOrder = 1;
343 else if (majorOrderFlag & SDDS_ROW_MAJOR_ORDER)
344 columnMajorOrder = 0;
345 break;
346 case CLO_MODIFYSIGMAS:
347 modifySigmas = 1;
348 break;
349 case CLO_AUTOOFFSET:
350 autoOffset = 1;
351 break;
352 case CLO_ORDERS:
353 if (s_arg[iArg].n_items < 2)
355 order =
tmalloc(
sizeof(*order) * (terms = s_arg[iArg].n_items - 1));
356 for (i = 1; i < s_arg[iArg].n_items; i++) {
357 if (sscanf(s_arg[iArg].list[i], "%" SCNd32, order + i - 1) != 1)
358 SDDS_Bomb(
"unable to scan order from -orders list");
359 }
360 break;
361 case CLO_RANGE:
362 rangeFitOnly = 0;
363 if ((s_arg[iArg].n_items != 3 && s_arg[iArg].n_items != 4) ||
364 1 != sscanf(s_arg[iArg].list[1], "%lf", &xMin) ||
365 1 != sscanf(s_arg[iArg].list[2], "%lf", &xMax) || xMin >= xMax)
367 if (s_arg[iArg].n_items == 4) {
368 if (strncmp(
str_tolower(s_arg[iArg].list[3]),
"fitonly",
369 strlen(s_arg[iArg].list[3])) == 0) {
370 rangeFitOnly = 1;
371 } else
373 }
374 break;
375 case CLO_GENERATESIGMAS:
376 generateSigmas = FLGS_GENERATESIGMAS;
377 if (s_arg[iArg].n_items > 1) {
378 if (s_arg[iArg].n_items != 2)
379 SDDS_Bomb(
"incorrect -generateSigmas syntax");
380 if (strncmp(s_arg[iArg].list[1], "keepsmallest",
381 strlen(s_arg[iArg].list[1])) == 0)
382 generateSigmas |= FLGS_KEEPSMALLEST;
383 if (strncmp(s_arg[iArg].list[1], "keeplargest",
384 strlen(s_arg[iArg].list[1])) == 0)
385 generateSigmas |= FLGS_KEEPLARGEST;
386 if ((generateSigmas & FLGS_KEEPSMALLEST) &&
387 (generateSigmas & FLGS_KEEPLARGEST))
388 SDDS_Bomb(
"ambiguous -generateSigmas syntax");
389 }
390 break;
391 case CLO_TERMS:
392 if (s_arg[iArg].n_items != 2 ||
393 sscanf(s_arg[iArg].list[1], "%ld", &terms) != 1)
395 break;
396 case CLO_XOFFSET:
397 if (s_arg[iArg].n_items != 2 ||
398 sscanf(s_arg[iArg].list[1], "%lf", &xOffset) != 1)
400 break;
401 case CLO_SYMMETRY:
402 if (s_arg[iArg].n_items == 2) {
403 if ((symmetry =
match_string(s_arg[iArg].list[1], symmetry_options,
404 N_SYMMETRY_OPTIONS, 0)) < 0)
405 SDDS_Bomb(
"unknown option used with -symmetry");
406 } else
408 break;
409 case CLO_SIGMAS:
410 if (s_arg[iArg].n_items != 3)
412 if (sscanf(s_arg[iArg].list[1], "%lf", &sigmas) != 1)
413 SDDS_Bomb(
"couldn't scan value for -sigmas");
414 if ((sigmasMode =
match_string(s_arg[iArg].list[2], sigmas_options,
415 N_SIGMAS_OPTIONS, 0)) < 0)
417 break;
418 case CLO_SPARSE:
419 if (s_arg[iArg].n_items != 2)
421 if (sscanf(s_arg[iArg].list[1], "%ld", &sparseInterval) != 1)
422 SDDS_Bomb(
"couldn't scan value for -sparse");
423 if (sparseInterval < 1)
425 break;
426 case CLO_VERBOSE:
427 verbose = 1;
428 break;
429 case CLO_NORMALIZE:
430 normTerm = 0;
431 if (s_arg[iArg].n_items > 2 ||
432 (s_arg[iArg].n_items == 2 &&
433 sscanf(s_arg[iArg].list[1], "%ld", &normTerm) != 1) ||
434 normTerm < 0)
436 break;
437 case CLO_REVISEORDERS:
438 revpowThreshold = 0.1;
439 revpowCompleteThres = 10;
440 goodEnoughChi = 1;
441 s_arg[iArg].n_items -= 1;
443 &s_arg[iArg].n_items, 0,
445 "complete",
SDDS_DOUBLE, &revpowCompleteThres, 1, REVPOW_COMPLETE,
447 "verbose", -1, NULL, 1, REVPOW_VERBOSE, NULL) ||
448 revpowThreshold < 0 || revpowCompleteThres < 0 || goodEnoughChi < 0)
449 SDDS_Bomb(
"invalid -reviseOrders syntax");
450 reviseOrders |= REVPOW_ACTIVE;
451 break;
452 case CLO_CHEBYSHEV:
453 if (s_arg[iArg].n_items > 2 ||
454 (s_arg[iArg].n_items == 2 &&
455 strncmp(s_arg[iArg].list[1], "convert",
456 strlen(s_arg[iArg].list[1])) != 0))
458 chebyshev = s_arg[iArg].n_items;
461 break;
462 case CLO_XFACTOR:
463 if (s_arg[iArg].n_items != 2 ||
464 sscanf(s_arg[iArg].list[1], "%lf", &xScaleFactor) != 1 ||
465 xScaleFactor == 0)
467 break;
468 case CLO_COLUMNS:
469 if (s_arg[iArg].n_items < 3 || s_arg[iArg].n_items > 5)
471 xName = s_arg[iArg].list[1];
472 yName = s_arg[iArg].list[2];
473 s_arg[iArg].n_items -= 3;
474 if (!
scanItemList(&flags, s_arg[iArg].list + 3, &s_arg[iArg].n_items, 0,
475 "xsigma",
SDDS_STRING, &xSigmaName, 1, 0,
"ysigma",
478 break;
479 case CLO_FITLABELFORMAT:
480 if (s_arg[iArg].n_items != 2)
481 SDDS_Bomb(
"invalid -fitLabelFormat syntax");
482 fitLabelFormat = s_arg[iArg].list[1];
483 break;
484 case CLO_PIPE:
486 &pipeFlags))
488 break;
489 case CLO_EVALUATE:
490 if (s_arg[iArg].n_items < 2)
492 evalParameters.file = s_arg[iArg].list[1];
493 s_arg[iArg].n_items -= 2;
494 s_arg[iArg].list += 2;
495 if (!
scanItemList(&evalParameters.flags, s_arg[iArg].list,
496 &s_arg[iArg].n_items, 0,
497 "begin",
SDDS_DOUBLE, &evalParameters.begin, 1, EVAL_BEGIN_GIVEN,
498 "end",
SDDS_DOUBLE, &evalParameters.end, 1, EVAL_END_GIVEN,
499 "number",
SDDS_LONG64, &evalParameters.number, 1, EVAL_NUMBER_GIVEN,
500 "valuesfile",
SDDS_STRING, &evalParameters.valuesFile, 1, EVAL_VALUESFILE_GIVEN,
501 "valuescolumn",
SDDS_STRING, &evalParameters.valuesColumn, 1, EVAL_VALUESCOLUMN_GIVEN,
502 "reusepage", 0, NULL, 0, EVAL_REUSE_PAGE_GIVEN,
503 NULL))
505 if (evalParameters.flags & EVAL_VALUESFILE_GIVEN || evalParameters.flags & EVAL_VALUESCOLUMN_GIVEN) {
506 if (evalParameters.flags & (EVAL_BEGIN_GIVEN | EVAL_END_GIVEN | EVAL_NUMBER_GIVEN))
507 SDDS_Bomb(
"invalid -evaluate syntax: given begin/end/number or valuesFile/valuesColumn, not a mixture.");
508 if (!(evalParameters.flags & EVAL_VALUESFILE_GIVEN && evalParameters.flags & EVAL_VALUESCOLUMN_GIVEN))
509 SDDS_Bomb(
"invalid -evaluate syntax: give both valuesFile and valuesColumn, not just one");
510 }
511 evalParameters.initialized = 0;
512 break;
513 case CLO_COPY_PARAMETERS:
514 copyParameters = 1;
515 break;
516 default:
517 bomb(
"unknown switch", USAGE);
518 break;
519 }
520 } else {
521 if (input == NULL)
522 input = s_arg[iArg].list[0];
523 else if (output == NULL)
524 output = s_arg[iArg].list[0];
525 else
527 }
528 }
529
531
532 if (symmetry && order)
533 SDDS_Bomb(
"can't specify both -symmetry and -orders");
534 if (chebyshev && order)
535 SDDS_Bomb(
"can't specify both -chebyshev and -orders");
536 if (chebyshev && symmetry)
537 SDDS_Bomb(
"can't specify both -chebyshev and -symmetry");
538 if (!xName || !yName)
539 SDDS_Bomb(
"you must specify a column name for x and y");
540
541 if (modifySigmas && !xSigmaName)
542 SDDS_Bomb(
"you must specify x sigmas with -modifySigmas");
543 if (generateSigmas) {
544 if (modifySigmas)
545 SDDS_Bomb(
"you can't specify both -generateSigmas and -modifySigmas");
546 }
547 if (ySigmaName) {
548 if (sigmasMode != -1)
549 SDDS_Bomb(
"you can't specify both -sigmas and a y sigma name");
550 }
551 ySigmasValid = 0;
552 if (sigmasMode != -1 || generateSigmas || ySigmaName || modifySigmas)
553 ySigmasValid = 1;
554
555 if (normTerm >= 0 && normTerm >= terms)
556 SDDS_Bomb(
"can't normalize to that term--not that many terms");
557 if (reviseOrders && !(sigmasMode != -1 || generateSigmas || ySigmaName))
559 "can't use -reviseOrders unless a y sigma or -generateSigmas is given");
560
561 if (symmetry == EVEN_SYMMETRY) {
562 order =
tmalloc(
sizeof(*order) * terms);
563 for (i = 0; i < terms; i++)
564 order[i] = 2 * i;
565 } else if (symmetry == ODD_SYMMETRY) {
566 order =
tmalloc(
sizeof(*order) * terms);
567 for (i = 0; i < terms; i++)
568 order[i] = 2 * i + 1;
569 } else if (!order) {
570 order =
tmalloc(
sizeof(*order) * terms);
571 for (i = 0; i < terms; i++)
572 order[i] = i;
573 }
574 coef =
tmalloc(
sizeof(*coef) * terms);
575 coefSigma =
tmalloc(
sizeof(*coefSigma) * terms);
576 iTerm =
tmalloc(
sizeof(*iTerm) * terms);
577 iTermSig =
tmalloc(
sizeof(*iTermSig) * terms);
578
581 checkInputFile(&SDDSin, xName, yName, xSigmaName, ySigmaName);
582 coefUnits = initializeOutputFile(&SDDSout, output, &SDDSin, input, xName,
583 yName, xSigmaName, ySigmaName, ySigmasValid,
584 order, terms, chebyshev, copyParameters);
585 if (columnMajorOrder != -1)
586 SDDSout.layout.data_mode.column_major = columnMajorOrder;
587 else
588 SDDSout.layout.data_mode.column_major =
589 SDDSin.layout.data_mode.column_major;
591 npages++;
592 invalid = 0;
594 pointsOrig = 0;
595 invalid = 1;
596 isFit = 0;
597 } else {
599 fprintf(stderr, "error: unable to read column %s\n", xName);
601 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
602 }
604 fprintf(stderr, "error: unable to read column %s\n", yName);
606 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
607 }
608 sx = NULL;
610 fprintf(stderr, "error: unable to read column %s\n", xSigmaName);
612 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
613 }
614 sy0 = NULL;
616 fprintf(stderr, "error: unable to read column %s\n", ySigmaName);
618 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
619 }
620 if (!sy0)
621 sy0 =
tmalloc(
sizeof(*sy0) * points);
622
623 if (xMin != xMax || sparseInterval != 1) {
624 xOrig =
tmalloc(
sizeof(*xOrig) * points);
625 yOrig =
tmalloc(
sizeof(*yOrig) * points);
626 if (sx)
627 sxOrig =
tmalloc(
sizeof(*sxOrig) * points);
628 if (ySigmasValid)
629 syOrig =
tmalloc(
sizeof(*syOrig) * points);
630 pointsOrig = points;
631 for (i = j = 0; i < points; i++) {
632 xOrig[i] = x[i];
633 yOrig[i] = y[i];
634 if (ySigmasValid)
635 syOrig[i] = sy0[i];
636 if (sx)
637 sxOrig[i] = sx[i];
638 }
639 if (xMin != xMax) {
640 for (i = j = 0; i < points; i++) {
641 if (xOrig[i] <= xMax && xOrig[i] >= xMin) {
642 x[j] = xOrig[i];
643 y[j] = yOrig[i];
644 if (ySigmasValid)
645 sy0[j] = syOrig[i];
646 if (sx)
647 sx[j] = sxOrig[i];
648 j++;
649 }
650 }
651 points = j;
652 }
653 if (sparseInterval != 1) {
654 for (i = j = 0; i < points; i++) {
655 if (i % sparseInterval == 0) {
656 x[j] = x[i];
657 y[j] = y[i];
658 if (ySigmasValid)
659 sy0[j] = sy0[i];
660 if (sx)
661 sx[j] = sx[i];
662 j++;
663 }
664 }
665 points = j;
666 }
667 } else {
668 xOrig = x;
669 yOrig = y;
670 sxOrig = sx;
671 syOrig = sy0;
672 pointsOrig = points;
673 }
674
676
677 if (sigmasMode == ABSOLUTE_SIGMAS) {
678 for (i = 0; i < points; i++)
679 sy0[i] = sigmas;
680 if (sy0 != syOrig)
681 for (i = 0; i < pointsOrig; i++)
682 syOrig[i] = sigmas;
683 } else if (sigmasMode == FRACTIONAL_SIGMAS) {
684 for (i = 0; i < points; i++)
685 sy0[i] = sigmas * fabs(y[i]);
686 if (sy0 != syOrig)
687 for (i = 0; i < pointsOrig; i++)
688 syOrig[i] = fabs(yOrig[i]) * sigmas;
689 }
690
691 if (!ySigmasValid || generateSigmas)
692 for (i = 0; i < points; i++)
693 sy0[i] = 1;
694 else
695 for (i = 0; i < points; i++)
696 if (sy0[i] == 0)
697 SDDS_Bomb(
"y sigma = 0 for one or more points.");
698
699 diff =
tmalloc(
sizeof(*x) * points);
700 sy =
tmalloc(
sizeof(*sy) * points);
701 for (i = 0; i < points; i++)
702 sy[i] = sy0[i];
703
705 xOffset = 0;
706
709 if (chebyshev) {
710 if (xOffset) {
711
712 xScaleFactor = MAX(fabs(xHigh - xOffset), fabs(xLow - xOffset));
713 } else {
714
715 xOffset = (xHigh + xLow) / 2;
716 xScaleFactor = (xHigh - xLow) / 2;
717 }
720 }
721
722 if (generateSigmas || modifySigmas) {
723
724 isFit =
lsfg(x, y, sy, points, terms, order, coef, coefSigma, &chi,
725 diff, basis_fn);
726 if (!isFit)
728 if (verbose) {
729 fputs("initial_fit:", stdout);
730 print_coefs(stdout, xOffset, xScaleFactor, chebyshev, coef, NULL,
731 order, terms, chi, normTerm, "");
732 fprintf(stdout, "unweighted rms deviation from fit: %21.15e\n",
733 rms_average(diff, points));
734 }
735 if (modifySigmas) {
736 if (!ySigmasValid) {
737 for (i = 0; i < points; i++)
738 sy[i] =
739 fabs(
eval_sum(basis_dfn, coef, order, terms, x[i]) * sx[i]);
740 } else
741 for (i = 0; i < points; i++) {
742 sy[i] = sqrt(
743 sqr(sy0[i]) +
744 sqr(
eval_sum(basis_dfn, coef, order, terms, x[i]) * sx[i]));
745 }
746 }
747 if (generateSigmas) {
748 double sigma;
749 for (i = sigma = 0; i < points; i++) {
750 sigma += sqr(diff[i]);
751 }
752 sigma = sqrt(sigma / (points - terms));
753 for (i = 0; i < points; i++) {
754 if (generateSigmas & FLGS_KEEPSMALLEST) {
755 if (sigma < sy[i])
756 sy[i] = sigma;
757 } else if (generateSigmas & FLGS_KEEPLARGEST) {
758 if (sigma > sy[i])
759 sy[i] = sigma;
760 } else {
761 sy[i] = sigma;
762 }
763 }
764 for (i = 0; i < pointsOrig; i++) {
765 if (generateSigmas & FLGS_KEEPSMALLEST) {
766 if (sigma < sy0[i])
767 sy0[i] = sigma;
768 } else if (generateSigmas & FLGS_KEEPLARGEST) {
769 if (sigma > sy0[i])
770 sy0[i] = sigma;
771 } else {
772 sy0[i] = sigma;
773 }
774 }
775 }
776 }
777
778 if (reviseOrders & REVPOW_ACTIVE) {
779 terms = reviseFitOrders(
780 x, y, sy, points, terms, order, coef, coefSigma, diff, basis_fn,
781 reviseOrders, xOffset, xScaleFactor, normTerm, ySigmasValid,
782 chebyshev, revpowThreshold, revpowCompleteThres, goodEnoughChi);
783 reviseOrders = 0;
784 }
785
786 isFit =
lsfg(x, y, sy, points, terms, order, coef, coefSigma, &chi, diff,
787 basis_fn);
788 if (isFit) {
789 rmsResidual = rms_average(diff, points);
790 if (verbose) {
791 print_coefs(stdout, xOffset, xScaleFactor, chebyshev, coef,
792 (ySigmasValid ? coefSigma : NULL), order, terms, chi,
793 normTerm, "");
794 fprintf(stdout, "unweighted rms deviation from fit: %21.15e\n",
795 rmsResidual);
796 }
797 } else if (verbose)
798 fprintf(stdout, "fit failed.\n");
799
800 if (evalParameters.file)
801 makeEvaluationTable(&evalParameters, x, points, coef, order, terms,
802 &SDDSin, xName, yName);
803 }
804
805 if (!
SDDS_StartPage(&SDDSout, rangeFitOnly ? pointsOrig : points))
807 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
808 rpnSeqBuffer[0] = 0;
809 if (!invalid) {
810 setCoefficientData(&SDDSout, coef, (ySigmasValid ? coefSigma : NULL),
811 coefUnits, order, terms, chebyshev, fitLabelFormat,
812 rpnSeqBuffer);
813 if (rangeFitOnly) {
814 double *residual;
815 compareOriginalToFit(xOrig, yOrig, &residual, pointsOrig, &rmsResidual,
816 coef, order, terms);
817
819 pointsOrig, ix) ||
821 pointsOrig, iy) ||
823 pointsOrig, iResidual))
825 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
826 for (i = 0; i < pointsOrig; i++)
827 residual[i] = yOrig[i] - residual[i];
828
830 pointsOrig, iFit))
832 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
833
834 if (ixSigma != -1 &&
836 pointsOrig, ixSigma))
838 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
839 if (ySigmasValid && iySigma != -1 &&
841 pointsOrig, iySigma))
843 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
844 free(residual);
845 } else {
846 for (i = 0; i < points; i++)
847 diff[i] =
848 -diff[i];
850 ix) ||
852 iy) ||
854 points, iResidual))
856 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
857 for (i = 0; i < points; i++)
858 diff[i] =
859 y[i] - diff[i];
861 points, iFit))
863 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
864
865 if (ixSigma != -1 &&
867 ixSigma))
869 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
870 if (ySigmasValid && iySigma != -1 &&
872 iySigma))
874 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
875 }
876 }
877
880 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
882 &SDDSout, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iRpnSequence,
883 invalid ? "" : rpnSeqBuffer, iRmsResidual,
884 invalid ? NaN : rmsResidual, iChiSq, invalid ? NaN : chi, iTerms,
885 terms, iSigLevel,
887 invalid ? NaN : xOffset, iFactor, invalid ? NaN : xScaleFactor,
888 iFitIsValid, isFit ? 'y' : 'n', -1) ||
891 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
892 if (!invalid) {
893 free(diff);
894 free(sy);
895 if (xOrig != x)
896 free(xOrig);
897 if (yOrig != y)
898 free(yOrig);
899 if (syOrig != sy0)
900 free(syOrig);
901 if (sxOrig != sx)
902 free(sxOrig);
903 free(x);
904 free(y);
905 free(sx);
906 free(sy0);
907 }
908 }
911 exit(EXIT_FAILURE);
912 }
913 if (evalParameters.initialized && !
SDDS_Terminate(&(evalParameters.dataset)))
915
916 free(coef);
917 free(coefSigma);
918 if (coefUnits)
919 free(coefUnits);
920 if (order)
921 free(order);
922
923 return (EXIT_SUCCESS);
924}
int32_t SDDS_CopyParameters(SDDS_DATASET *SDDS_target, SDDS_DATASET *SDDS_source)
int32_t SDDS_StartPage(SDDS_DATASET *SDDS_dataset, int64_t expected_n_rows)
int32_t SDDS_SetParameters(SDDS_DATASET *SDDS_dataset, int32_t mode,...)
int32_t SDDS_SetColumnFromDoubles(SDDS_DATASET *SDDS_dataset, int32_t mode, double *data, int64_t rows,...)
Sets the values for a single data column using double-precision floating-point numbers.
int32_t SDDS_WritePage(SDDS_DATASET *SDDS_dataset)
Writes the current data table to the output file.
void SDDS_RegisterProgramName(const char *name)
Registers the executable program name for use in error messages.
#define SDDS_LONG64
Identifier for the signed 64-bit integer data type.
void bomb(char *error, char *usage)
Reports error messages to the terminal and aborts the program.
int find_min_max(double *min, double *max, double *list, int64_t n)
Finds the minimum and maximum values in a list of doubles.
long lsfg(double *xd, double *yd, double *sy, long n_pts, long n_terms, int32_t *order, double *coef, double *s_coef, double *chi, double *diff, double(*fn)(double x, long ord))
Computes generalized least squares fits using a function passed by the caller.
long match_string(char *string, char **option, long n_options, long mode)
Matches a given string against an array of option strings based on specified modes.
int scanargs(SCANNED_ARG **scanned, int argc, char **argv)
long processPipeOption(char **item, long items, unsigned long *flags)
void processFilenames(char *programName, char **input, char **output, unsigned long pipeFlags, long noWarnings, long *tmpOutputUsed)
long scanItemList(unsigned long *flags, char **item, long *items, unsigned long mode,...)
Scans a list of items and assigns values based on provided keywords and types.
void set_argument_scale(double scale)
Set the scale factor applied to the input argument of basis functions.
void set_argument_offset(double offset)
Set the offset applied to the input argument of basis functions.
double dtcheby(double x, long n)
Evaluate the derivative of the Chebyshev polynomial T_n(x).
double ipower(double x, long n)
Evaluate a power function x^n.
double dipower(double x, long n)
Evaluate the derivative of x^n.
double tcheby(double x, long n)
Evaluate the Chebyshev polynomial of the first kind T_n(x).
double ChiSqrSigLevel(double ChiSquared0, long nu)
Computes the probability that a chi-squared variable exceeds a given value.
char * str_tolower(char *s)
Convert a string to lower case.