313 {
314 double *x = NULL, *y = NULL, *sy = NULL, *sx = NULL, *diff = NULL, xOffset,
315 xScaleFactor;
316 double *xOrig = NULL, *yOrig = NULL, *sxOrig, *syOrig, *sy0 = NULL;
317 long terms, normTerm, ySigmasValid;
318 int64_t i, j, points, pointsOrig;
319 long symmetry, chebyshev, autoOffset, copyParameters = 0;
320 double sigmas;
321 long sigmasMode, sparseInterval;
322 unsigned long flags;
323 double *coef, *coefSigma;
324 double chi, xLow, xHigh, rmsResidual;
325 char *xName, *yName, *xSigmaName, *ySigmaName;
326 char *input, *output, **coefUnits;
328 long isFit, iArg, modifySigmas;
329 long generateSigmas, verbose, ignoreSigmas;
330
331 long invalid = 0;
332 int32_t *order;
333 SCANNED_ARG *s_arg;
334 double xMin, xMax, revpowThreshold, revpowCompleteThres, goodEnoughChi;
335 long rangeFitOnly = 0;
336 double rms_average(double *d_x, int64_t d_n);
337 char *fitLabelFormat = "%g";
338 static char rpnSeqBuffer[SDDS_MAXLINE];
339 unsigned long pipeFlags, reviseOrders, majorOrderFlag;
341 short columnMajorOrder = -1;
342
343 sxOrig = syOrig = NULL;
344 rmsResidual = 0;
345
347 argc =
scanargs(&s_arg, argc, argv);
348 if (argc < 2 || argc > (3 + N_OPTIONS)) {
349 fprintf(stderr, "usage: %s%s%s\n", USAGE, additional_help1,
350 additional_help2);
351 exit(EXIT_FAILURE);
352 }
353
354 input = output = NULL;
355 xName = yName = xSigmaName = ySigmaName = NULL;
356 modifySigmas = reviseOrders = chebyshev = 0;
357 order = NULL;
358 symmetry = NO_SYMMETRY;
359 xMin = xMax = 0;
360 autoOffset = 0;
361 generateSigmas = 0;
362 sigmasMode = -1;
363 sigmas = 1;
364 sparseInterval = 1;
365 terms = 2;
366 verbose = ignoreSigmas = 0;
367 normTerm = -1;
368 xOffset = 0;
369 xScaleFactor = 1;
370 coefUnits = NULL;
373 pipeFlags = 0;
374 evalParameters.file = evalParameters.valuesFile = evalParameters.valuesColumn = NULL;
375 evalParameters.initialized = evalParameters.inputInitialized = 0;
376
377 for (iArg = 1; iArg < argc; iArg++) {
378 if (s_arg[iArg].arg_type == OPTION) {
379 switch (
match_string(s_arg[iArg].list[0], option, N_OPTIONS, 0)) {
380 case CLO_MAJOR_ORDER:
381 majorOrderFlag = 0;
382 s_arg[iArg].n_items--;
383 if (s_arg[iArg].n_items > 0 &&
385 &s_arg[iArg].n_items, 0, "row", -1, NULL, 0,
386 SDDS_ROW_MAJOR_ORDER, "column", -1, NULL, 0,
387 SDDS_COLUMN_MAJOR_ORDER, NULL)))
388 SDDS_Bomb(
"invalid -majorOrder syntax/values");
389 if (majorOrderFlag & SDDS_COLUMN_MAJOR_ORDER)
390 columnMajorOrder = 1;
391 else if (majorOrderFlag & SDDS_ROW_MAJOR_ORDER)
392 columnMajorOrder = 0;
393 break;
394 case CLO_MODIFYSIGMAS:
395 modifySigmas = 1;
396 break;
397 case CLO_AUTOOFFSET:
398 autoOffset = 1;
399 break;
400 case CLO_ORDERS:
401 if (s_arg[iArg].n_items < 2)
403 order =
tmalloc(
sizeof(*order) * (terms = s_arg[iArg].n_items - 1));
404 for (i = 1; i < s_arg[iArg].n_items; i++) {
405 if (sscanf(s_arg[iArg].list[i], "%" SCNd32, order + i - 1) != 1)
406 SDDS_Bomb(
"unable to scan order from -orders list");
407 }
408 break;
409 case CLO_RANGE:
410 rangeFitOnly = 0;
411 if ((s_arg[iArg].n_items != 3 && s_arg[iArg].n_items != 4) ||
412 1 != sscanf(s_arg[iArg].list[1], "%lf", &xMin) ||
413 1 != sscanf(s_arg[iArg].list[2], "%lf", &xMax) || xMin >= xMax)
415 if (s_arg[iArg].n_items == 4) {
416 if (strncmp(
str_tolower(s_arg[iArg].list[3]),
"fitonly",
417 strlen(s_arg[iArg].list[3])) == 0) {
418 rangeFitOnly = 1;
419 } else
421 }
422 break;
423 case CLO_GENERATESIGMAS:
424 generateSigmas = FLGS_GENERATESIGMAS;
425 if (s_arg[iArg].n_items > 1) {
426 if (s_arg[iArg].n_items != 2)
427 SDDS_Bomb(
"incorrect -generateSigmas syntax");
428 if (strncmp(s_arg[iArg].list[1], "keepsmallest",
429 strlen(s_arg[iArg].list[1])) == 0)
430 generateSigmas |= FLGS_KEEPSMALLEST;
431 if (strncmp(s_arg[iArg].list[1], "keeplargest",
432 strlen(s_arg[iArg].list[1])) == 0)
433 generateSigmas |= FLGS_KEEPLARGEST;
434 if ((generateSigmas & FLGS_KEEPSMALLEST) &&
435 (generateSigmas & FLGS_KEEPLARGEST))
436 SDDS_Bomb(
"ambiguous -generateSigmas syntax");
437 }
438 break;
439 case CLO_TERMS:
440 if (s_arg[iArg].n_items != 2 ||
441 sscanf(s_arg[iArg].list[1], "%ld", &terms) != 1)
443 break;
444 case CLO_XOFFSET:
445 if (s_arg[iArg].n_items != 2 ||
446 sscanf(s_arg[iArg].list[1], "%lf", &xOffset) != 1)
448 break;
449 case CLO_SYMMETRY:
450 if (s_arg[iArg].n_items == 2) {
451 if ((symmetry =
match_string(s_arg[iArg].list[1], symmetry_options,
452 N_SYMMETRY_OPTIONS, 0)) < 0)
453 SDDS_Bomb(
"unknown option used with -symmetry");
454 } else
456 break;
457 case CLO_SIGMAS:
458 if (s_arg[iArg].n_items != 3)
460 if (sscanf(s_arg[iArg].list[1], "%lf", &sigmas) != 1)
461 SDDS_Bomb(
"couldn't scan value for -sigmas");
462 if ((sigmasMode =
match_string(s_arg[iArg].list[2], sigmas_options,
463 N_SIGMAS_OPTIONS, 0)) < 0)
465 break;
466 case CLO_SPARSE:
467 if (s_arg[iArg].n_items != 2)
469 if (sscanf(s_arg[iArg].list[1], "%ld", &sparseInterval) != 1)
470 SDDS_Bomb(
"couldn't scan value for -sparse");
471 if (sparseInterval < 1)
473 break;
474 case CLO_VERBOSE:
475 verbose = 1;
476 break;
477 case CLO_NORMALIZE:
478 normTerm = 0;
479 if (s_arg[iArg].n_items > 2 ||
480 (s_arg[iArg].n_items == 2 &&
481 sscanf(s_arg[iArg].list[1], "%ld", &normTerm) != 1) ||
482 normTerm < 0)
484 break;
485 case CLO_REVISEORDERS:
486 revpowThreshold = 0.1;
487 revpowCompleteThres = 10;
488 goodEnoughChi = 1;
489 s_arg[iArg].n_items -= 1;
491 &s_arg[iArg].n_items, 0,
493 "complete",
SDDS_DOUBLE, &revpowCompleteThres, 1, REVPOW_COMPLETE,
495 "verbose", -1, NULL, 1, REVPOW_VERBOSE, NULL) ||
496 revpowThreshold < 0 || revpowCompleteThres < 0 || goodEnoughChi < 0)
497 SDDS_Bomb(
"invalid -reviseOrders syntax");
498 reviseOrders |= REVPOW_ACTIVE;
499 break;
500 case CLO_CHEBYSHEV:
501 if (s_arg[iArg].n_items > 2 ||
502 (s_arg[iArg].n_items == 2 &&
503 strncmp(s_arg[iArg].list[1], "convert",
504 strlen(s_arg[iArg].list[1])) != 0))
506 chebyshev = s_arg[iArg].n_items;
509 break;
510 case CLO_XFACTOR:
511 if (s_arg[iArg].n_items != 2 ||
512 sscanf(s_arg[iArg].list[1], "%lf", &xScaleFactor) != 1 ||
513 xScaleFactor == 0)
515 break;
516 case CLO_COLUMNS:
517 if (s_arg[iArg].n_items < 3 || s_arg[iArg].n_items > 5)
519 xName = s_arg[iArg].list[1];
520 yName = s_arg[iArg].list[2];
521 s_arg[iArg].n_items -= 3;
522 if (!
scanItemList(&flags, s_arg[iArg].list + 3, &s_arg[iArg].n_items, 0,
523 "xsigma",
SDDS_STRING, &xSigmaName, 1, 0,
"ysigma",
526 break;
527 case CLO_FITLABELFORMAT:
528 if (s_arg[iArg].n_items != 2)
529 SDDS_Bomb(
"invalid -fitLabelFormat syntax");
530 fitLabelFormat = s_arg[iArg].list[1];
531 break;
532 case CLO_PIPE:
534 &pipeFlags))
536 break;
537 case CLO_EVALUATE:
538 if (s_arg[iArg].n_items < 2)
540 evalParameters.file = s_arg[iArg].list[1];
541 s_arg[iArg].n_items -= 2;
542 s_arg[iArg].list += 2;
543 if (!
scanItemList(&evalParameters.flags, s_arg[iArg].list,
544 &s_arg[iArg].n_items, 0,
545 "begin",
SDDS_DOUBLE, &evalParameters.begin, 1, EVAL_BEGIN_GIVEN,
546 "end",
SDDS_DOUBLE, &evalParameters.end, 1, EVAL_END_GIVEN,
547 "number",
SDDS_LONG64, &evalParameters.number, 1, EVAL_NUMBER_GIVEN,
548 "valuesfile",
SDDS_STRING, &evalParameters.valuesFile, 1, EVAL_VALUESFILE_GIVEN,
549 "valuescolumn",
SDDS_STRING, &evalParameters.valuesColumn, 1, EVAL_VALUESCOLUMN_GIVEN,
550 "reusepage", 0, NULL, 0, EVAL_REUSE_PAGE_GIVEN,
551 NULL))
553 if (evalParameters.flags & EVAL_VALUESFILE_GIVEN || evalParameters.flags & EVAL_VALUESCOLUMN_GIVEN) {
554 if (evalParameters.flags & (EVAL_BEGIN_GIVEN | EVAL_END_GIVEN | EVAL_NUMBER_GIVEN))
555 SDDS_Bomb(
"invalid -evaluate syntax: given begin/end/number or valuesFile/valuesColumn, not a mixture.");
556 if (!(evalParameters.flags & EVAL_VALUESFILE_GIVEN && evalParameters.flags & EVAL_VALUESCOLUMN_GIVEN))
557 SDDS_Bomb(
"invalid -evaluate syntax: give both valuesFile and valuesColumn, not just one");
558 }
559 evalParameters.initialized = 0;
560 break;
561 case CLO_COPY_PARAMETERS:
562 copyParameters = 1;
563 break;
564 default:
565 bomb(
"unknown switch", USAGE);
566 break;
567 }
568 } else {
569 if (input == NULL)
570 input = s_arg[iArg].list[0];
571 else if (output == NULL)
572 output = s_arg[iArg].list[0];
573 else
575 }
576 }
577
579
580 if (symmetry && order)
581 SDDS_Bomb(
"can't specify both -symmetry and -orders");
582 if (chebyshev && order)
583 SDDS_Bomb(
"can't specify both -chebyshev and -orders");
584 if (chebyshev && symmetry)
585 SDDS_Bomb(
"can't specify both -chebyshev and -symmetry");
586 if (!xName || !yName)
587 SDDS_Bomb(
"you must specify a column name for x and y");
588
589 if (modifySigmas && !xSigmaName)
590 SDDS_Bomb(
"you must specify x sigmas with -modifySigmas");
591 if (generateSigmas) {
592 if (modifySigmas)
593 SDDS_Bomb(
"you can't specify both -generateSigmas and -modifySigmas");
594 }
595 if (ySigmaName) {
596 if (sigmasMode != -1)
597 SDDS_Bomb(
"you can't specify both -sigmas and a y sigma name");
598 }
599 ySigmasValid = 0;
600 if (sigmasMode != -1 || generateSigmas || ySigmaName || modifySigmas)
601 ySigmasValid = 1;
602
603 if (normTerm >= 0 && normTerm >= terms)
604 SDDS_Bomb(
"can't normalize to that term--not that many terms");
605 if (reviseOrders && !(sigmasMode != -1 || generateSigmas || ySigmaName))
606 SDDS_Bomb(
"can't use -reviseOrders unless a y sigma or -generateSigmas is given");
607
608 if (symmetry == EVEN_SYMMETRY) {
609 order =
tmalloc(
sizeof(*order) * terms);
610 for (i = 0; i < terms; i++)
611 order[i] = 2 * i;
612 } else if (symmetry == ODD_SYMMETRY) {
613 order =
tmalloc(
sizeof(*order) * terms);
614 for (i = 0; i < terms; i++)
615 order[i] = 2 * i + 1;
616 } else if (!order) {
617 order =
tmalloc(
sizeof(*order) * terms);
618 for (i = 0; i < terms; i++)
619 order[i] = i;
620 }
621 coef =
tmalloc(
sizeof(*coef) * terms);
622 coefSigma =
tmalloc(
sizeof(*coefSigma) * terms);
623 iTerm =
tmalloc(
sizeof(*iTerm) * terms);
624 iTermSig =
tmalloc(
sizeof(*iTermSig) * terms);
625
628 checkInputFile(&SDDSin, xName, yName, xSigmaName, ySigmaName);
629 coefUnits = initializeOutputFile(&SDDSout, output, &SDDSin, input, xName,
630 yName, xSigmaName, ySigmaName, ySigmasValid,
631 order, terms, chebyshev, copyParameters);
632 if (columnMajorOrder != -1)
633 SDDSout.layout.data_mode.column_major = columnMajorOrder;
634 else
635 SDDSout.layout.data_mode.column_major =
636 SDDSin.layout.data_mode.column_major;
638
639 invalid = 0;
641 pointsOrig = 0;
642 invalid = 1;
643 isFit = 0;
644 } else {
646 fprintf(stderr, "error: unable to read column %s\n", xName);
648 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
649 }
651 fprintf(stderr, "error: unable to read column %s\n", yName);
653 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
654 }
655 sx = NULL;
657 fprintf(stderr, "error: unable to read column %s\n", xSigmaName);
659 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
660 }
661 sy0 = NULL;
663 fprintf(stderr, "error: unable to read column %s\n", ySigmaName);
665 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
666 }
667 if (!sy0)
668 sy0 =
tmalloc(
sizeof(*sy0) * points);
669
670 if (xMin != xMax || sparseInterval != 1) {
671 xOrig =
tmalloc(
sizeof(*xOrig) * points);
672 yOrig =
tmalloc(
sizeof(*yOrig) * points);
673 if (sx)
674 sxOrig =
tmalloc(
sizeof(*sxOrig) * points);
675 if (ySigmasValid)
676 syOrig =
tmalloc(
sizeof(*syOrig) * points);
677 pointsOrig = points;
678 for (i = j = 0; i < points; i++) {
679 xOrig[i] = x[i];
680 yOrig[i] = y[i];
681 if (ySigmasValid)
682 syOrig[i] = sy0[i];
683 if (sx)
684 sxOrig[i] = sx[i];
685 }
686 if (xMin != xMax) {
687 for (i = j = 0; i < points; i++) {
688 if (xOrig[i] <= xMax && xOrig[i] >= xMin) {
689 x[j] = xOrig[i];
690 y[j] = yOrig[i];
691 if (ySigmasValid)
692 sy0[j] = syOrig[i];
693 if (sx)
694 sx[j] = sxOrig[i];
695 j++;
696 }
697 }
698 points = j;
699 }
700 if (sparseInterval != 1) {
701 for (i = j = 0; i < points; i++) {
702 if (i % sparseInterval == 0) {
703 x[j] = x[i];
704 y[j] = y[i];
705 if (ySigmasValid)
706 sy0[j] = sy0[i];
707 if (sx)
708 sx[j] = sx[i];
709 j++;
710 }
711 }
712 points = j;
713 }
714 } else {
715 xOrig = x;
716 yOrig = y;
717 sxOrig = sx;
718 syOrig = sy0;
719 pointsOrig = points;
720 }
721
723
724 if (sigmasMode == ABSOLUTE_SIGMAS) {
725 for (i = 0; i < points; i++)
726 sy0[i] = sigmas;
727 if (sy0 != syOrig)
728 for (i = 0; i < pointsOrig; i++)
729 syOrig[i] = sigmas;
730 } else if (sigmasMode == FRACTIONAL_SIGMAS) {
731 for (i = 0; i < points; i++)
732 sy0[i] = sigmas * fabs(y[i]);
733 if (sy0 != syOrig)
734 for (i = 0; i < pointsOrig; i++)
735 syOrig[i] = fabs(yOrig[i]) * sigmas;
736 }
737
738 if (!ySigmasValid || generateSigmas)
739 for (i = 0; i < points; i++)
740 sy0[i] = 1;
741 else
742 for (i = 0; i < points; i++)
743 if (sy0[i] == 0)
744 SDDS_Bomb(
"y sigma = 0 for one or more points.");
745
746 diff =
tmalloc(
sizeof(*x) * points);
747 sy =
tmalloc(
sizeof(*sy) * points);
748 for (i = 0; i < points; i++)
749 sy[i] = sy0[i];
750
752 xOffset = 0;
753
756 if (chebyshev) {
757 if (xOffset) {
758
759 xScaleFactor = MAX(fabs(xHigh - xOffset), fabs(xLow - xOffset));
760 } else {
761
762 xOffset = (xHigh + xLow) / 2;
763 xScaleFactor = (xHigh - xLow) / 2;
764 }
767 }
768
769 if (generateSigmas || modifySigmas) {
770
771 isFit =
lsfg(x, y, sy, points, terms, order, coef, coefSigma, &chi,
772 diff, basis_fn);
773 if (!isFit)
775 if (verbose) {
776 fputs("initial_fit:", stdout);
777 print_coefs(stdout, xOffset, xScaleFactor, chebyshev, coef, NULL,
778 order, terms, chi, normTerm, "");
779 fprintf(stdout, "unweighted rms deviation from fit: %21.15e\n",
780 rms_average(diff, points));
781 }
782 if (modifySigmas) {
783 if (!ySigmasValid) {
784 for (i = 0; i < points; i++)
785 sy[i] =
786 fabs(
eval_sum(basis_dfn, coef, order, terms, x[i]) * sx[i]);
787 } else
788 for (i = 0; i < points; i++) {
789 sy[i] = sqrt(
790 sqr(sy0[i]) +
791 sqr(
eval_sum(basis_dfn, coef, order, terms, x[i]) * sx[i]));
792 }
793 }
794 if (generateSigmas) {
795 double sigma;
796 for (i = sigma = 0; i < points; i++) {
797 sigma += sqr(diff[i]);
798 }
799 sigma = sqrt(sigma / (points - terms));
800 for (i = 0; i < points; i++) {
801 if (generateSigmas & FLGS_KEEPSMALLEST) {
802 if (sigma < sy[i])
803 sy[i] = sigma;
804 } else if (generateSigmas & FLGS_KEEPLARGEST) {
805 if (sigma > sy[i])
806 sy[i] = sigma;
807 } else {
808 sy[i] = sigma;
809 }
810 }
811 for (i = 0; i < pointsOrig; i++) {
812 if (generateSigmas & FLGS_KEEPSMALLEST) {
813 if (sigma < sy0[i])
814 sy0[i] = sigma;
815 } else if (generateSigmas & FLGS_KEEPLARGEST) {
816 if (sigma > sy0[i])
817 sy0[i] = sigma;
818 } else {
819 sy0[i] = sigma;
820 }
821 }
822 }
823 }
824
825 if (reviseOrders & REVPOW_ACTIVE) {
826 terms = reviseFitOrders(
827 x, y, sy, points, terms, order, coef, coefSigma, diff, basis_fn,
828 reviseOrders, xOffset, xScaleFactor, normTerm, ySigmasValid,
829 chebyshev, revpowThreshold, revpowCompleteThres, goodEnoughChi);
830 reviseOrders = 0;
831 }
832
833 isFit =
lsfg(x, y, sy, points, terms, order, coef, coefSigma, &chi, diff,
834 basis_fn);
835 if (isFit) {
836 rmsResidual = rms_average(diff, points);
837 if (verbose) {
838 print_coefs(stdout, xOffset, xScaleFactor, chebyshev, coef,
839 (ySigmasValid ? coefSigma : NULL), order, terms, chi,
840 normTerm, "");
841 fprintf(stdout, "unweighted rms deviation from fit: %21.15e\n",
842 rmsResidual);
843 }
844 } else if (verbose)
845 fprintf(stdout, "fit failed.\n");
846
847 if (evalParameters.file)
848 makeEvaluationTable(&evalParameters, x, points, coef, order, terms,
849 &SDDSin, xName, yName);
850 }
851
852 if (!
SDDS_StartPage(&SDDSout, rangeFitOnly ? pointsOrig : points))
854 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
855 rpnSeqBuffer[0] = 0;
856 if (!invalid) {
857 setCoefficientData(&SDDSout, coef, (ySigmasValid ? coefSigma : NULL),
858 coefUnits, order, terms, chebyshev, fitLabelFormat,
859 rpnSeqBuffer);
860 if (rangeFitOnly) {
861 double *residual;
862 compareOriginalToFit(xOrig, yOrig, &residual, pointsOrig, &rmsResidual,
863 coef, order, terms);
864
866 pointsOrig, ix) ||
868 pointsOrig, iy) ||
870 pointsOrig, iResidual))
872 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
873 for (i = 0; i < pointsOrig; i++)
874 residual[i] = yOrig[i] - residual[i];
875
877 pointsOrig, iFit))
879 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
880
881 if (ixSigma != -1 &&
883 pointsOrig, ixSigma))
885 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
886 if (ySigmasValid && iySigma != -1 &&
888 pointsOrig, iySigma))
890 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
891 free(residual);
892 } else {
893 for (i = 0; i < points; i++)
894 diff[i] =
895 -diff[i];
897 ix) ||
899 iy) ||
901 points, iResidual))
903 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
904 for (i = 0; i < points; i++)
905 diff[i] =
906 y[i] - diff[i];
908 points, iFit))
910 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
911
912 if (ixSigma != -1 &&
914 ixSigma))
916 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
917 if (ySigmasValid && iySigma != -1 &&
919 iySigma))
921 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
922 }
923 }
924
927 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
929 &SDDSout, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iRpnSequence,
930 invalid ? "" : rpnSeqBuffer, iRmsResidual,
931 invalid ? NaN : rmsResidual, iChiSq, invalid ? NaN : chi, iTerms,
932 terms, iSigLevel,
934 invalid ? NaN : xOffset, iFactor, invalid ? NaN : xScaleFactor,
935 iFitIsValid, isFit ? 'y' : 'n', -1) ||
938 SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
939 if (!invalid) {
940 free(diff);
941 free(sy);
942 if (xOrig != x)
943 free(xOrig);
944 if (yOrig != y)
945 free(yOrig);
946 if (syOrig != sy0)
947 free(syOrig);
948 if (sxOrig != sx)
949 free(sxOrig);
950 free(x);
951 free(y);
952 free(sx);
953 free(sy0);
954 }
955 }
958 exit(EXIT_FAILURE);
959 }
960 if (evalParameters.initialized && !
SDDS_Terminate(&(evalParameters.dataset)))
962
963 free(coef);
964 free(coefSigma);
965 if (coefUnits)
966 free(coefUnits);
967 if (order)
968 free(order);
969
970 return (EXIT_SUCCESS);
971}
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.