56void print_coefs(FILE *fprec,
double x_offset,
double x_scale,
long chebyshev,
double *coef,
double *s_coef,
57 int32_t *order,
long n_terms,
double chi,
long norm_term,
char *prepend);
58char **makeCoefficientUnits(
SDDS_DATASET *SDDSout,
char *xName,
char *yName, int32_t *order,
long terms);
59long setCoefficientData(
SDDS_DATASET *SDDSout,
double *coef,
double *coefSigma,
char **coefUnits,
long *order,
long terms);
61 SDDS_DATASET *SDDSin,
char *input,
char *xName,
char **yNames,
char *xSigmaName,
62 char **ySigmaNames,
long sigmasValid, int32_t *order,
long terms,
long isChebyshev,
63 long numCols,
long copyParameters);
64void checkInputFile(
SDDS_DATASET *SDDSin,
char *xName,
char **yNames,
char *xSigmaName,
char **ySigmaNames,
long numYNames);
65long coefficient_index(int32_t *order,
long terms,
long order_of_interest);
66void makeFitLabel(
char *buffer,
long bufsize,
char *fitLabelFormat,
double *coef, int32_t *order,
long terms,
long colIndex);
68char **ResolveColumnNames(
SDDS_DATASET *SDDSin,
char **wildcardList,
long length, int32_t *numYNames);
69char **GenerateYSigmaNames(
char *controlString,
char **yNames,
long numYNames);
70void RemoveElementFromStringArray(
char **array,
long index,
long length);
71char **RemoveNonNumericColumnsFromNameArray(
SDDS_DATASET *SDDSin,
char **columns, int32_t *numColumns);
72void compareOriginalToFit(
double *x,
double *y,
double **residual, int64_t points,
double *rmsResidual,
double *coef, int32_t *order,
long terms);
76double tcheby(
double x,
long n);
77double dtcheby(
double x,
long n);
78double ipower(
double x,
long n);
79double dipower(
double x,
long n);
102 CLO_SIGMAINDEPENDENT,
110char *option[N_OPTIONS] = {
138 "sddsmpfit [options] [<inputfile>] [<outputfile>]\n\n"
140 " -pipe=[input][,output] Use standard input/output.\n"
141 " -independent=<xName> Specify the independent data column.\n"
142 " -dependent=<yname1-wildcard>[,<yname2-wildcard>,...] \n"
143 " Specify dependent data columns with wildcards.\n"
144 " -sigmaIndependent=<xSigma> Specify independent sigma column.\n"
145 " -sigmaDependent=<ySigmaFormatString> Specify dependent sigma columns using a format string.\n"
146 " -terms=<number> Number of terms in the fit.\n"
147 " -symmetry={none|odd|even} Symmetry of the fit about x_offset.\n"
148 " -orders=<number>[,<number>,...] Specify the orders to use in the fit.\n"
149 " -reviseOrders[=threshold=<value>][,verbose] Revise orders to eliminate poorly-determined coefficients.\n"
150 " -chebyshev[=convert] Use Chebyshev T polynomials. 'convert' to output ordinary polynomials.\n"
151 " -xOffset=<value> Set the x offset for the fit.\n"
152 " -xFactor=<value> Set the scaling factor for x values.\n"
153 " -sigmas=<value>,{absolute|fractional} Set sigma values as absolute or fractional.\n"
154 " -minimumSigma=<value> Set the minimum sigma value.\n"
155 " -modifySigmas Modify y sigmas using x sigmas and initial fit.\n"
156 " -generateSigmas={keepLargest|keepSmallest} Generate y sigmas from RMS deviation.\n"
157 " -sparse=<interval> Sample data at the specified interval.\n"
158 " -range=<lower>,<upper>[,fitOnly] Set the range for fitting and evaluation. 'fitOnly' compares the fit to original data range.\n"
159 " -normalize[=<termNumber>] Normalize the fit to make a specified term unity.\n"
160 " -verbose Enable verbose output.\n"
161 " -evaluate=<filename>[,begin=<value>][,end=<value>][,number=<integer>] \n"
162 " Evaluate the fit over a range of points.\n"
163 " -fitLabelFormat=<sprintf-string> Set the format string for fit labels.\n"
164 " -infoFile=<filename> Specify an information file for fit coefficients and statistics.\n"
165 " -copyParameters Copy parameters from input to output.\n\n"
166 "Program by Michael Borland, revised by Brad Dolin.\n"
167 "Compiled on " __DATE__
" " __TIME__
", SVN revision: " SVN_VERSION
"\n";
169static char *additional_help =
"\n\
170sddsmpfit does fits to the form y = SUM(i){ A[i] *P(x-x_offset, i)}, where P(x,i) is the ith basis\n\
171function evaluated at x. sddsmpfit returns the A[i] and estimates of the errors in the values.\n\
172By default P(x,i) = x^i. One can also select Chebyshev T polynomials as the basis functions.\n\n\
173-independent specify name of independent data column to use.\n\
174-dependent specify names of dependent data columns to use, using wildcards,\n\
175 separated by commas.\n\
176-sigmaIndependent specify name of independent sigma values to use\n\
177-sigmaDependent specify names of dependent sigma values to use by specifying a printf-style control\n\
178 string to generate the names from the independent variable names (e.g., %sSigma)\n\
179-terms number of terms desired in fit.\n\
180-symmetry symmetry of desired fit about x_offset.\n\
181-orders orders (P[i]) to use in fitting.\n\
182-reviseOrders the orders used in the fit are modified from the specified ones\n\
183 in order eliminate poorly-determined coefficients, based on fitting\n\
184 of the first data page.\n";
185static char *additional_help2 =
"-chebyshev use Chebyshev T polynomials (xOffset is set automatically).\n\
186 Giving the `convert' option causes the fit to be written out in\n\
187 terms of ordinary polynomials.\n\
188-xOffset desired value of x to fit about.\n\
189-xFactor desired factor to multiply x values by before fitting.\n\
190-sigmas specify absolute or fractional sigma for all points.\n\
191-minimumSigma specify minimum sigma value. If the value is less than this\n\
192 it is replaced by this value.\n\
193-modifySigmas modify the y sigmas using the x sigmas and an initial fit.\n\
194-generateSigmas generate y sigmas from the rms deviation from an initial fit.\n\
195 optionally keep the sigmas from the data if larger/smaller than rms\n\
197-sparse specify integer interval at which to sample data.\n\
198-range specify range of independent variable over which to perform fit and evaluation.\n\
199 If 'fitOnly' is given, then fit is compared to data over the original range.\n\
200-normalize normalize so that specified term is unity.\n\
201-evaluate specify evaluation of fit over a selected range of\n\
202 equispaced points.\n\
203-infoFile specify name of optional information file containing coefficients and fit statistics.\n\
204-copyParameters specify that parameters from input should be copied to output.\n\
205-verbose generates extra output that may be useful.\n\n";
208#define EVEN_SYMMETRY 1
209#define ODD_SYMMETRY 2
210#define N_SYMMETRY_OPTIONS 3
211char *symmetry_options[N_SYMMETRY_OPTIONS] = {
"none",
"even",
"odd"};
213#define ABSOLUTE_SIGMAS 0
214#define FRACTIONAL_SIGMAS 1
215#define N_SIGMAS_OPTIONS 2
216char *sigmas_options[N_SIGMAS_OPTIONS] = {
"absolute",
"fractional"};
218#define FLGS_GENERATESIGMAS 1
219#define FLGS_KEEPLARGEST 2
220#define FLGS_KEEPSMALLEST 4
222#define REVPOW_ACTIVE 0x0001
223#define REVPOW_VERBOSE 0x0002
225static long *iIntercept = NULL, *iInterceptO = NULL, *iInterceptSigma = NULL, *iInterceptSigmaO = NULL;
226static long *iSlope = NULL, *iSlopeO = NULL, *iSlopeSigma = NULL, *iSlopeSigmaO = NULL;
227static long *iCurvature = NULL, *iCurvatureO = NULL, *iCurvatureSigma = NULL, *iCurvatureSigmaO = NULL;
228static long iOffset = -1, iOffsetO = -1, iFactor = -1, iFactorO = -1;
229static long *iChiSq = NULL, *iChiSqO = NULL, *iRmsResidual = NULL, *iRmsResidualO = NULL, *iSigLevel = NULL, *iSigLevelO = NULL;
230static long *iFitIsValid = NULL, *iFitIsValidO = NULL, *iFitLabel = NULL, *iFitLabelO = NULL, iTerms = -1, iTermsO = -1;
232static long ix = -1, ixSigma = -1;
233static long *iy = NULL, *iySigma = NULL;
234static long *iFit = NULL, *iResidual = NULL;
236static long iOrder = -1, *iCoefficient = NULL, *iCoefficientSigma = NULL, *iCoefficientUnits = NULL;
238static char *xSymbol, **ySymbols;
240#define EVAL_BEGIN_GIVEN 0x0001U
241#define EVAL_END_GIVEN 0x0002U
242#define EVAL_NUMBER_GIVEN 0x0004U
244#define MAX_Y_SIGMA_NAME_SIZE 1024
256void makeEvaluationTable(
EVAL_PARAMETERS *evalParameters,
double *x, int64_t points,
double *coef, int32_t *order,
257 long terms,
char *xName,
char **yName,
long yNames,
long iYName);
259static double (*basis_fn)(
double xa,
long ordera);
260static double (*basis_dfn)(
double xa,
long ordera);
262int main(
int argc,
char **argv) {
263 double **y = NULL, **sy = NULL, **diff = NULL;
264 double *x = NULL, *sx = NULL;
265 double xOffset, xScaleFactor;
266 double *xOrig = NULL, **yOrig = NULL, *sxOrig = NULL, **syOrig = NULL, **sy0 = NULL;
267 long terms, normTerm, ip, ySigmasValid;
268 int64_t i, j, points, pointsOrig;
269 long symmetry, chebyshev, termsGiven;
270 double sigmas, minimumSigma;
271 long sigmasMode, sparseInterval;
272 double **coef = NULL, **coefSigma = NULL;
273 double *chi = NULL, xLow, xHigh, *rmsResidual = NULL;
274 char *xName = NULL, *yName = NULL, **yNames = NULL, *xSigmaName = NULL;
275 char **ySigmaNames = NULL, *ySigmaControlString = NULL;
276 char *input = NULL, *output = NULL;
278 long *isFit = NULL, iArg, modifySigmas, termIndex;
279 long generateSigmas, verbose, ignoreSigmas;
280 long outputInitialized, copyParameters = 0;
281 int32_t *order = NULL;
283 double xMin, xMax, revpowThreshold;
284 double rms_average(
double *d_x, int64_t d_n);
285 char *infoFile = NULL;
286 char *fitLabelFormat =
"%g";
287 static char fitLabelBuffer[SDDS_MAXLINE];
288 unsigned long pipeFlags, reviseOrders;
290 long rangeFitOnly = 0;
293 long cloDependentIndex = -1, numDependentItems;
297 argc =
scanargs(&s_arg, argc, argv);
298 if (argc < 2 || argc > (3 + N_OPTIONS)) {
299 fprintf(stderr,
"usage: %s\n", USAGE);
300 fprintf(stderr,
"%s%s", additional_help, additional_help2);
304 input = output = NULL;
305 xName = yName = xSigmaName = ySigmaControlString = NULL;
306 yNames = ySigmaNames = NULL;
307 numDependentItems = 0;
308 modifySigmas = reviseOrders = chebyshev = 0;
310 symmetry = NO_SYMMETRY;
318 verbose = ignoreSigmas = 0;
325 evalParameters.file = NULL;
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_MODIFYSIGMAS:
337 SDDS_Bomb(
"give -order or -terms, not both");
338 if (s_arg[iArg].n_items < 2)
340 order =
tmalloc(
sizeof(*order) * (terms = s_arg[iArg].n_items - 1));
341 for (i = 1; i < s_arg[iArg].n_items; i++) {
342 if (sscanf(s_arg[iArg].list[i],
"%" SCNd32, order + i - 1) != 1)
343 SDDS_Bomb(
"unable to scan order from -orders list");
348 if ((s_arg[iArg].n_items != 3 && s_arg[iArg].n_items != 4) || 1 != sscanf(s_arg[iArg].list[1],
"%lf", &xMin) || 1 != sscanf(s_arg[iArg].list[2],
"%lf", &xMax) || xMin >= xMax)
350 if (s_arg[iArg].n_items == 4) {
351 if (strncmp(
str_tolower(s_arg[iArg].list[3]),
"fitonly", strlen(s_arg[iArg].list[3])) == 0) {
357 case CLO_GENERATESIGMAS:
358 generateSigmas = FLGS_GENERATESIGMAS;
359 if (s_arg[iArg].n_items > 1) {
360 if (s_arg[iArg].n_items != 2)
361 SDDS_Bomb(
"incorrect -generateSigmas synax");
362 if (strncmp(s_arg[iArg].list[1],
"keepsmallest", strlen(s_arg[iArg].list[1])) == 0)
363 generateSigmas |= FLGS_KEEPSMALLEST;
364 if (strncmp(s_arg[iArg].list[1],
"keeplargest", strlen(s_arg[iArg].list[1])) == 0)
365 generateSigmas |= FLGS_KEEPLARGEST;
366 if ((generateSigmas & FLGS_KEEPSMALLEST) && (generateSigmas & FLGS_KEEPLARGEST))
367 SDDS_Bomb(
"ambiguous -generateSigmas synax");
372 SDDS_Bomb(
"give -order or -terms, not both");
373 if (s_arg[iArg].n_items != 2 || sscanf(s_arg[iArg].list[1],
"%ld", &terms) != 1)
378 if (s_arg[iArg].n_items != 2 || sscanf(s_arg[iArg].list[1],
"%lf", &xOffset) != 1)
382 if (s_arg[iArg].n_items == 2) {
383 if ((symmetry =
match_string(s_arg[iArg].list[1], symmetry_options, N_SYMMETRY_OPTIONS, 0)) < 0)
384 SDDS_Bomb(
"unknown option used with -symmetry");
389 if (s_arg[iArg].n_items != 3)
391 if (sscanf(s_arg[iArg].list[1],
"%lf", &sigmas) != 1)
392 SDDS_Bomb(
"couldn't scan value for -sigmas");
393 if ((sigmasMode =
match_string(s_arg[iArg].list[2], sigmas_options, N_SIGMAS_OPTIONS, 0)) < 0)
397 if (s_arg[iArg].n_items != 2)
398 SDDS_Bomb(
"incorrect -minimumSigma syntax");
399 if (sscanf(s_arg[iArg].list[1],
"%lf", &minimumSigma) != 1)
400 SDDS_Bomb(
"couldn't scan value for -minimumSigma");
403 if (s_arg[iArg].n_items != 2)
405 if (sscanf(s_arg[iArg].list[1],
"%ld", &sparseInterval) != 1)
406 SDDS_Bomb(
"couldn't scan value for -sparse");
407 if (sparseInterval < 1)
415 if (s_arg[iArg].n_items > 2 ||
416 (s_arg[iArg].n_items == 2 && sscanf(s_arg[iArg].list[1],
"%ld", &normTerm) != 1) ||
420 case CLO_REVISEORDERS:
421 revpowThreshold = 0.1;
422 s_arg[iArg].n_items -= 1;
423 if (!
scanItemList(&reviseOrders, s_arg[iArg].list + 1, &s_arg[iArg].n_items, 0,
425 "verbose", -1, NULL, 1, REVPOW_VERBOSE, NULL))
426 SDDS_Bomb(
"invalid -reviseOrders syntax");
427 reviseOrders |= REVPOW_ACTIVE;
428 revpowThreshold = fabs(revpowThreshold);
431 if (s_arg[iArg].n_items > 2 ||
432 (s_arg[iArg].n_items == 2 && strncmp(s_arg[iArg].list[1],
"convert", strlen(s_arg[iArg].list[1])) != 0))
434 chebyshev = s_arg[iArg].n_items;
439 if (s_arg[iArg].n_items != 2 ||
440 sscanf(s_arg[iArg].list[1],
"%lf", &xScaleFactor) != 1 || xScaleFactor == 0)
443 case CLO_INDEPENDENT:
444 if (s_arg[iArg].n_items != 2)
445 SDDS_Bomb(
"invalid -independent syntax");
446 xName = s_arg[iArg].list[1];
449 numDependentItems = s_arg[iArg].n_items - 1;
450 cloDependentIndex = iArg;
451 if (numDependentItems < 1)
454 case CLO_SIGMAINDEPENDENT:
455 if (s_arg[iArg].n_items != 2)
456 SDDS_Bomb(
"invalid -sigmaIndependent syntax");
457 xSigmaName = s_arg[iArg].list[1];
459 case CLO_SIGMADEPENDENT:
460 if (s_arg[iArg].n_items != 2)
461 SDDS_Bomb(
"invalid -sigmaDependent syntax");
462 ySigmaControlString = s_arg[iArg].list[1];
464 case CLO_FITLABELFORMAT:
465 if (s_arg[iArg].n_items != 2)
466 SDDS_Bomb(
"invalid -fitLabelFormat syntax");
467 fitLabelFormat = s_arg[iArg].list[1];
470 if (!
processPipeOption(s_arg[iArg].list + 1, s_arg[iArg].n_items - 1, &pipeFlags))
474 if (s_arg[iArg].n_items != 2)
476 infoFile = s_arg[iArg].list[1];
479 if (s_arg[iArg].n_items < 2)
481 evalParameters.file = s_arg[iArg].list[1];
482 s_arg[iArg].n_items -= 2;
483 s_arg[iArg].list += 2;
484 if (!
scanItemList(&evalParameters.flags, s_arg[iArg].list, &s_arg[iArg].n_items, 0,
485 "begin",
SDDS_DOUBLE, &evalParameters.begin, 1, EVAL_BEGIN_GIVEN,
486 "end",
SDDS_DOUBLE, &evalParameters.end, 1, EVAL_END_GIVEN,
487 "number",
SDDS_LONG64, &evalParameters.number, 1, EVAL_NUMBER_GIVEN, NULL))
490 case CLO_COPYPARAMETERS:
494 bomb(
"unknown switch", USAGE);
499 input = s_arg[iArg].list[0];
500 else if (output == NULL)
501 output = s_arg[iArg].list[0];
509 if (symmetry && order)
510 SDDS_Bomb(
"can't specify both -symmetry and -orders");
511 if (!xName || !numDependentItems)
512 SDDS_Bomb(
"you must specify a column name for x and y");
513 if (modifySigmas && !xSigmaName)
514 SDDS_Bomb(
"you must specify x sigmas with -modifySigmas");
515 if (generateSigmas) {
517 SDDS_Bomb(
"you can't specify both -generateSigmas and -modifySigmas");
519 if (ySigmaControlString) {
520 if (sigmasMode != -1)
521 SDDS_Bomb(
"you can't specify both -sigmas and a y sigma name");
524 if (sigmasMode != -1 || generateSigmas || ySigmaControlString || modifySigmas)
527 if (normTerm >= 0 && normTerm >= terms)
528 SDDS_Bomb(
"can't normalize to that term--not that many terms");
529 if (reviseOrders && !(sigmasMode != -1 || generateSigmas || ySigmaNames))
530 SDDS_Bomb(
"can't use -reviseOrders unless a y sigma or -generateSigmas is given");
532 if (symmetry == EVEN_SYMMETRY) {
533 order =
tmalloc(
sizeof(*order) * terms);
534 for (i = 0; i < terms; i++)
536 }
else if (symmetry == ODD_SYMMETRY) {
537 order =
tmalloc(
sizeof(*order) * terms);
538 for (i = 0; i < terms; i++)
539 order[i] = 2 * i + 1;
541 order =
tmalloc(
sizeof(*order) * terms);
542 for (i = 0; i < terms; i++)
548 outputInitialized = 0;
549 yNames = ResolveColumnNames(&SDDSin, s_arg[cloDependentIndex].list + 1, numDependentItems, &numYNames);
550 if (ySigmaControlString != NULL)
551 ySigmaNames = GenerateYSigmaNames(ySigmaControlString, yNames, numYNames);
553 checkInputFile(&SDDSin, xName, yNames, xSigmaName, ySigmaNames, numYNames);
554 sy0 =
tmalloc(
sizeof(
double *) * numYNames);
555 y =
tmalloc(
sizeof(
double *) * numYNames);
556 sy =
tmalloc(
sizeof(
double *) * numYNames);
557 isFit =
tmalloc(
sizeof(
long) * numYNames);
558 chi =
tmalloc(
sizeof(
double) * numYNames);
559 coef =
tmalloc(
sizeof(
double *) * numYNames);
560 coefSigma =
tmalloc(
sizeof(
double *) * numYNames);
561 for (colIndex = 0; colIndex < numYNames; colIndex++) {
562 coef[colIndex] =
tmalloc(
sizeof(
double) * terms);
563 coefSigma[colIndex] =
tmalloc(
sizeof(
double) * terms);
565 iCoefficient =
tmalloc(
sizeof(
long) * numYNames);
566 iCoefficientSigma =
tmalloc(
sizeof(
long) * numYNames);
567 iCoefficientUnits =
tmalloc(
sizeof(
long) * numYNames);
575 fprintf(stderr,
"error: unable to read column %s\n", xName);
578 for (i = 0; i < numYNames; i++) {
580 fprintf(stderr,
"error: unable to read column %s\n", yNames[i]);
586 fprintf(stderr,
"error: unable to read column %s\n", xSigmaName);
589 for (colIndex = 0; colIndex < numYNames; colIndex++)
590 sy0[colIndex] =
tmalloc(
sizeof(
double) * points);
592 for (i = 0; i < numYNames; i++) {
594 fprintf(stderr,
"error: unable to read column %s\n", ySigmaNames[i]);
600 if (minimumSigma > 0) {
602 for (i = 0; i < numYNames; i++) {
603 for (j = 0; j < points; j++)
604 if (sy0[i][j] < minimumSigma)
605 sy0[i][j] = minimumSigma;
609 if (xMin != xMax || sparseInterval != 1) {
610 xOrig =
tmalloc(
sizeof(*xOrig) * points);
611 yOrig =
tmalloc(
sizeof(*yOrig) * numYNames);
612 for (colIndex = 0; colIndex < numYNames; colIndex++)
613 yOrig[colIndex] =
tmalloc(
sizeof(
double) * points);
615 sxOrig =
tmalloc(
sizeof(*sxOrig) * points);
617 syOrig =
tmalloc(
sizeof(*syOrig) * numYNames);
618 for (colIndex = 0; colIndex < numYNames; colIndex++)
619 syOrig[colIndex] =
tmalloc(
sizeof(
double) * points);
622 for (i = j = 0; i < points; i++) {
626 for (colIndex = 0; colIndex < numYNames; colIndex++) {
627 yOrig[colIndex][i] = y[colIndex][i];
629 syOrig[colIndex][i] = sy0[colIndex][i];
633 for (i = j = 0; i < points; i++) {
634 if (xOrig[i] <= xMax && xOrig[i] >= xMin) {
636 for (colIndex = 0; colIndex < numYNames; colIndex++) {
637 y[colIndex][j] = yOrig[colIndex][i];
639 sy0[colIndex][j] = syOrig[colIndex][i];
648 if (sparseInterval != 1) {
649 for (i = j = 0; i < points; i++) {
650 if (i % sparseInterval == 0) {
652 for (colIndex = 0; colIndex < numYNames; colIndex++) {
653 y[colIndex][j] = y[colIndex][i];
655 sy0[colIndex][j] = sy0[colIndex][i];
674 if (sigmasMode == ABSOLUTE_SIGMAS) {
675 for (colIndex = 0; colIndex < numYNames; colIndex++) {
676 for (i = 0; i < points; i++)
677 sy0[colIndex][i] = sigmas;
678 if (sy0[colIndex] != syOrig[colIndex])
679 for (i = 0; i < pointsOrig; i++)
680 syOrig[colIndex][i] = sigmas;
682 }
else if (sigmasMode == FRACTIONAL_SIGMAS) {
683 for (colIndex = 0; colIndex < numYNames; colIndex++) {
684 for (i = 0; i < points; i++)
685 sy0[colIndex][i] = sigmas * fabs(y[colIndex][i]);
686 if (sy0[colIndex] != syOrig[colIndex])
687 for (i = 0; i < pointsOrig; i++)
688 syOrig[colIndex][i] = fabs(yOrig[colIndex][i]) * sigmas;
692 for (i = 0; i < numYNames; i++) {
693 if (minimumSigma > 0) {
695 for (j = 0; j < points; j++)
696 if (sy0[i][j] < minimumSigma)
697 sy0[i][j] = minimumSigma;
701 if (!ySigmasValid || generateSigmas)
702 for (colIndex = 0; colIndex < numYNames; colIndex++) {
703 for (i = 0; i < points; i++)
704 sy0[colIndex][i] = 1;
707 for (i = 0; i < points; i++)
708 for (colIndex = 0; colIndex < numYNames; colIndex++) {
709 if (sy0[colIndex][i] == 0)
710 SDDS_Bomb(
"y sigma = 0 for one or more points.");
713 diff =
tmalloc(
sizeof(*diff) * numYNames);
714 sy =
tmalloc(
sizeof(*sy) * numYNames);
715 for (colIndex = 0; colIndex < numYNames; colIndex++) {
716 diff[colIndex] =
tmalloc(
sizeof(
double) * points);
717 sy[colIndex] =
tmalloc(
sizeof(
double) * points);
720 for (i = 0; i < points; i++) {
721 for (colIndex = 0; colIndex < numYNames; colIndex++)
722 sy[colIndex][i] = sy0[colIndex][i];
728 xOffset = (xHigh + xLow) / 2;
733 if (generateSigmas || modifySigmas) {
735 for (colIndex = 0; colIndex < numYNames; colIndex++) {
736 isFit[colIndex] =
lsfg(x, y[colIndex], sy[colIndex], points, terms, order, coef[colIndex], coefSigma[colIndex], &chi[colIndex], diff[colIndex], basis_fn);
737 if (!isFit[colIndex]) {
738 fprintf(stderr,
"Column %s: ", yNames[colIndex]);
742 fprintf(stderr,
"Column %s: ", yNames[colIndex]);
743 fputs(
"initial_fit:", stderr);
744 print_coefs(stderr, xOffset, xScaleFactor, chebyshev, coef[colIndex], NULL, order, terms, chi[colIndex], normTerm,
"");
745 fprintf(stderr,
"unweighted rms deviation from fit: %21.15le\n", rms_average(diff[colIndex], points));
749 for (i = 0; i < points; i++)
750 sy[colIndex][i] = fabs(
eval_sum(basis_dfn, coef[colIndex], order, terms, x[i]) * sx[i]);
752 for (i = 0; i < points; i++) {
753 sy[colIndex][i] = sqrt(sqr(sy0[colIndex][i]) + sqr(
eval_sum(basis_dfn, coef[colIndex], order, terms, x[i]) * sx[i]));
756 if (generateSigmas) {
758 for (i = sigma = 0; i < points; i++) {
759 sigma += sqr(diff[colIndex][i]);
761 sigma = sqrt(sigma / (points - terms));
762 for (i = 0; i < points; i++) {
763 if (generateSigmas & FLGS_KEEPSMALLEST) {
764 if (sigma < sy[colIndex][i])
765 sy[colIndex][i] = sigma;
766 }
else if (generateSigmas & FLGS_KEEPLARGEST) {
767 if (sigma > sy[colIndex][i])
768 sy[colIndex][i] = sigma;
770 sy[colIndex][i] = sigma;
773 for (i = 0; i < pointsOrig; i++) {
774 if (generateSigmas & FLGS_KEEPSMALLEST) {
775 if (sigma < sy0[colIndex][i])
776 sy0[colIndex][i] = sigma;
777 }
else if (generateSigmas & FLGS_KEEPLARGEST) {
778 if (sigma > sy0[colIndex][i])
779 sy0[colIndex][i] = sigma;
781 sy0[colIndex][i] = sigma;
788 if (reviseOrders & REVPOW_ACTIVE) {
790 long bestTerms, newBest;
794 bestOrder =
tmalloc(
sizeof(*bestOrder) * bestTerms);
795 for (ip = 0; ip < terms; ip++)
796 bestOrder[ip] = order[ip];
798 for (colIndex = 0; colIndex < numYNames; colIndex++) {
799 isFit[colIndex] =
lsfg(x, y[colIndex], sy[colIndex], points, bestTerms, bestOrder, coef[colIndex], coefSigma[colIndex], &bestChi, diff[colIndex], basis_fn);
800 if (!isFit[colIndex]) {
801 fprintf(stderr,
"Column %s: ", yNames[colIndex]);
803 if (reviseOrders & REVPOW_VERBOSE) {
804 fprintf(stderr,
"Column %s: ", yNames[colIndex]);
805 fputs(
"fit to revise orders:", stderr);
806 print_coefs(stderr, xOffset, xScaleFactor, chebyshev, coef[colIndex], (ySigmasValid ? coefSigma[colIndex] : NULL), bestOrder, bestTerms, bestChi, normTerm,
"");
807 fprintf(stderr,
"unweighted rms deviation from fit: %21.15le\n", rms_average(diff[colIndex], points));
813 terms = bestTerms - 1;
814 for (ip = bestTerms - 1; ip >= 0; ip--) {
815 for (i = j = 0; i < bestTerms; i++)
817 order[j++] = bestOrder[i];
818 isFit[colIndex] =
lsfg(x, y[colIndex], sy[colIndex], points, terms, order, coef[colIndex], coefSigma[colIndex], &chi[colIndex], diff[colIndex], basis_fn);
819 if (!isFit[colIndex]) {
820 fprintf(stderr,
"Column %s: ", yNames[colIndex]);
823 if (reviseOrders & REVPOW_VERBOSE) {
824 fprintf(stderr,
"Column %s: ", yNames[colIndex]);
825 fputs(
"new trial fit:", stderr);
826 print_coefs(stderr, xOffset, xScaleFactor, chebyshev, coef[colIndex], (ySigmasValid ? coefSigma[colIndex] : NULL), order, terms, chi[colIndex], normTerm,
"");
827 fprintf(stderr,
"unweighted rms deviation from fit: %21.15le\n", rms_average(diff[colIndex], points));
829 if (chi[colIndex] - bestChi < revpowThreshold) {
830 bestChi = chi[colIndex];
833 for (i = 0; i < terms; i++)
834 bestOrder[i] = order[i];
835 if (reviseOrders & REVPOW_VERBOSE) {
836 fputs(
"new best fit:", stderr);
837 print_coefs(stderr, xOffset, xScaleFactor, chebyshev, coef[colIndex], (ySigmasValid ? coefSigma[colIndex] : NULL), bestOrder, bestTerms, bestChi, normTerm,
"");
838 fprintf(stderr,
"unweighted rms deviation from fit: %21.15le\n", rms_average(diff[colIndex], points));
847 for (ip = 0; ip < terms; ip++)
848 order[ip] = bestOrder[ip];
854 if (!outputInitialized) {
855 initializeOutputFile(&SDDSout, &SDDSoutInfo, output, infoFile, &SDDSin, input, xName, yNames, xSigmaName, ySigmaNames, ySigmasValid, order, terms, chebyshev, numYNames, copyParameters);
857 outputInitialized = 1;
859 if (evalParameters.file)
860 setupEvaluationFile(&evalParameters, xName, yNames, numYNames, &SDDSin);
862 rmsResidual =
tmalloc(
sizeof(
double) * numYNames);
863 for (colIndex = 0; colIndex < numYNames; colIndex++) {
864 isFit[colIndex] =
lsfg(x, y[colIndex], sy[colIndex], points, terms, order, coef[colIndex], coefSigma[colIndex], &chi[colIndex], diff[colIndex], basis_fn);
865 if (isFit[colIndex]) {
866 rmsResidual[colIndex] = rms_average(diff[colIndex], points);
868 fprintf(stderr,
"Column: %s\n", yNames[colIndex]);
869 print_coefs(stderr, xOffset, xScaleFactor, chebyshev, coef[colIndex], (ySigmasValid ? coefSigma[colIndex] : NULL), order, terms, chi[colIndex], normTerm,
"");
870 fprintf(stderr,
"unweighted rms deviation from fit: %21.15le\n", rmsResidual[colIndex]);
873 fprintf(stderr,
"fit failed for %s.\n", yNames[colIndex]);
875 if (evalParameters.file)
876 makeEvaluationTable(&evalParameters, x, points, coef[colIndex], order, terms, xName, yNames, numYNames, colIndex);
879 if (outputInitialized) {
880 if (!
SDDS_StartPage(&SDDSout, rangeFitOnly ? pointsOrig : points) ||
883 if (copyParameters) {
889 if (!
SDDS_SetColumnFromDoubles(&SDDSout, SDDS_SET_BY_INDEX, rangeFitOnly ? xOrig : x, rangeFitOnly ? pointsOrig : points, ix) ||
892 for (colIndex = 0; colIndex < numYNames; colIndex++) {
894 double *residual, rmsResidual0;
895 compareOriginalToFit(xOrig, yOrig[colIndex], &residual, pointsOrig, &rmsResidual0, coef[colIndex], order, terms);
899 for (i = 0; i < pointsOrig; i++)
900 residual[i] = yOrig[colIndex][i] - residual[i];
905 for (i = 0; i < points; i++)
906 diff[colIndex][i] = -diff[colIndex][i];
910 for (i = 0; i < points; i++)
911 diff[colIndex][i] = y[colIndex][i] - diff[colIndex][i];
917 !
SDDS_SetColumnFromDoubles(&SDDSout, SDDS_SET_BY_INDEX, rangeFitOnly ? sxOrig : sx, rangeFitOnly ? pointsOrig : points, ixSigma))
919 for (colIndex = 0; colIndex < numYNames; colIndex++) {
920 if (ySigmasValid && iySigma[colIndex] != -1 &&
921 !
SDDS_SetColumnFromDoubles(&SDDSout, SDDS_SET_BY_INDEX, rangeFitOnly ? syOrig[colIndex] : sy[colIndex], rangeFitOnly ? pointsOrig : points, iySigma[colIndex]))
925 termIndex = coefficient_index(order, terms, 0);
926 if (iIntercept[colIndex] != -1 &&
927 !
SDDS_SetParameters(&SDDSoutInfo, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iIntercept[colIndex], coef[colIndex][termIndex], -1))
929 if (iInterceptSigma[colIndex] != -1 &&
930 !
SDDS_SetParameters(&SDDSoutInfo, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iInterceptSigma[colIndex], coefSigma[colIndex][termIndex], -1))
933 termIndex = coefficient_index(order, terms, 1);
934 if (iSlope[colIndex] != -1 &&
935 !
SDDS_SetParameters(&SDDSoutInfo, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iSlope[colIndex], coef[colIndex][termIndex], -1))
937 if (iSlopeSigma[colIndex] != -1 &&
938 !
SDDS_SetParameters(&SDDSoutInfo, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iSlopeSigma[colIndex], coefSigma[colIndex][termIndex], -1))
941 termIndex = coefficient_index(order, terms, 2);
942 if (iCurvature[colIndex] != -1 &&
943 !
SDDS_SetParameters(&SDDSoutInfo, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iCurvature[colIndex], coef[colIndex][termIndex], -1))
945 if (iCurvatureSigma[colIndex] != -1 &&
946 !
SDDS_SetParameters(&SDDSoutInfo, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iCurvatureSigma[colIndex], coefSigma[colIndex][termIndex], -1))
948 if (iFitLabel[colIndex] != -1) {
949 makeFitLabel(fitLabelBuffer, SDDS_MAXLINE, fitLabelFormat, coef[colIndex], order, terms, colIndex);
950 if (!
SDDS_SetParameters(&SDDSoutInfo, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iFitLabel[colIndex], fitLabelBuffer, -1))
958 iRmsResidual[colIndex], rmsResidual[colIndex], iChiSq[colIndex], chi[colIndex],
959 iTerms, terms, iSigLevel[colIndex],
ChiSqrSigLevel(chi[colIndex], points - terms),
960 iOffset, xOffset, iFactor, xScaleFactor, iFitIsValid[colIndex], isFit[colIndex] ?
'y' :
'n', -1))
964 termIndex = coefficient_index(order, terms, 0);
965 if (iInterceptO[colIndex] != -1 &&
966 !
SDDS_SetParameters(&SDDSout, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iInterceptO[colIndex], coef[colIndex][termIndex], -1))
968 if (iInterceptSigmaO[colIndex] != -1 &&
969 !
SDDS_SetParameters(&SDDSout, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iInterceptSigmaO[colIndex], coefSigma[colIndex][termIndex], -1))
972 termIndex = coefficient_index(order, terms, 1);
973 if (iSlopeO[colIndex] != -1 &&
974 !
SDDS_SetParameters(&SDDSout, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iSlopeO[colIndex], coef[colIndex][termIndex], -1))
976 if (iSlopeSigmaO[colIndex] != -1 &&
977 !
SDDS_SetParameters(&SDDSout, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iSlopeSigmaO[colIndex], coefSigma[colIndex][termIndex], -1))
980 termIndex = coefficient_index(order, terms, 2);
981 if (iCurvatureO[colIndex] != -1 &&
982 !
SDDS_SetParameters(&SDDSout, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iCurvatureO[colIndex], coef[colIndex][termIndex], -1))
984 if (iCurvatureSigmaO[colIndex] != -1 &&
985 !
SDDS_SetParameters(&SDDSout, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iCurvatureSigmaO[colIndex], coefSigma[colIndex][termIndex], -1))
987 if (iFitLabelO[colIndex] != -1) {
988 makeFitLabel(fitLabelBuffer, SDDS_MAXLINE, fitLabelFormat, coef[colIndex], order, terms, colIndex);
989 if (!
SDDS_SetParameters(&SDDSout, SDDS_SET_BY_INDEX | SDDS_PASS_BY_VALUE, iFitLabelO[colIndex], fitLabelBuffer, -1))
993 iRmsResidualO[colIndex], rmsResidual[colIndex], iChiSqO[colIndex], chi[colIndex],
994 iTermsO, terms, iSigLevelO[colIndex],
ChiSqrSigLevel(chi[colIndex], points - terms),
995 iOffsetO, xOffset, iFactorO, xScaleFactor, iFitIsValidO[colIndex], isFit[colIndex] ?
'y' :
'n', -1))
1007 for (colIndex = 0; colIndex < numYNames; colIndex++) {
1008 free(diff[colIndex]);
1010 if (yOrig[colIndex] != y[colIndex])
1011 free(yOrig[colIndex]);
1012 if (syOrig && sy0 && syOrig[colIndex] != sy0[colIndex])
1013 free(syOrig[colIndex]);
1015 if (sy0 && sy0[colIndex])
1016 free(sy0[colIndex]);
1019 return (EXIT_SUCCESS);
1022void print_coefs(FILE * fpo,
double xOffset,
double xScaleFactor,
long chebyshev,
double *coef,
double *coefSigma,
1023 int32_t *order,
long terms,
double chi,
long normTerm,
char *prepend) {
1027 fprintf(fpo,
"%s%ld-term Chebyshev T polynomial least-squares fit about x=%21.15le, scaled by %21.15le:\n", prepend, terms, xOffset, xScaleFactor);
1029 fprintf(fpo,
"%s%ld-term polynomial least-squares fit about x=%21.15le:\n", prepend, terms, xOffset);
1030 if (normTerm >= 0 && terms > normTerm) {
1031 if (coef[normTerm] != 0)
1032 fprintf(fpo,
"%s coefficients are normalized with factor %21.15le to make a[%ld]==1\n", prepend, coef[normTerm], (order ? order[normTerm] : normTerm));
1034 fprintf(fpo,
"%s can't normalize coefficients as requested: a[%ld]==0\n", prepend, (order ? order[normTerm] : normTerm));
1040 for (i = 0; i < terms; i++) {
1041 fprintf(fpo,
"%sa[%ld] = %21.15le ", prepend, (order ? order[i] : i), (normTerm < 0 ? coef[i] : coef[i] / coef[normTerm]));
1043 fprintf(fpo,
"+/- %21.15le\n", (normTerm < 0 ? coefSigma[i] : coefSigma[i] / fabs(coef[normTerm])));
1048 fprintf(fpo,
"%sreduced chi-squared = %21.15le\n", prepend, chi);
1051void RemoveElementFromStringArray(
char **array,
long index,
long length) {
1054 for (lh = index; lh < length - 1; lh++)
1055 array[lh] = array[lh + 1];
1058char **RemoveNonNumericColumnsFromNameArray(
SDDS_DATASET * SDDSin,
char **columns, int32_t *numColumns) {
1059 long i, numNumericColumns = *numColumns;
1061 for (i = 0; i < *numColumns; i++) {
1063 printf(
"Removing %s because not a numeric type.\n", columns[i]);
1064 RemoveElementFromStringArray(columns, i, *numColumns);
1065 numNumericColumns--;
1069 *numColumns = numNumericColumns;
1073char **ResolveColumnNames(
SDDS_DATASET * SDDSin,
char **wildcardList,
long length, int32_t *numYNames) {
1079 for (i = 0; i < length; i++) {
1084 bomb(
"Error matching columns in ResolveColumnNames: No matches.", NULL);
1086 result = RemoveNonNumericColumnsFromNameArray(SDDSin, result, numYNames);
1090char **GenerateYSigmaNames(
char *controlString,
char **yNames,
long numYNames) {
1092 char **result, sigmaName[MAX_Y_SIGMA_NAME_SIZE];
1094 result =
tmalloc(
sizeof(
char *) * numYNames);
1095 for (i = 0; i < numYNames; i++) {
1096 sprintf(sigmaName, controlString, yNames[i]);
1097 nameLength = strlen(sigmaName);
1098 result[i] =
tmalloc(
sizeof(
char) * (nameLength + 1));
1099 strcpy(result[i], sigmaName);
1104void makeFitLabel(
char *buffer,
long bufsize,
char *fitLabelFormat,
double *coef, int32_t *order,
long terms,
long colIndex) {
1106 static char buffer1[SDDS_MAXLINE], buffer2[SDDS_MAXLINE];
1108 sprintf(buffer,
"%s = ", ySymbols[colIndex]);
1109 for (i = 0; i < terms; i++) {
1111 sprintf(buffer1, fitLabelFormat, coef[i]);
1114 strcpy(buffer1,
" +");
1115 sprintf(buffer1 + 2, fitLabelFormat, coef[i]);
1117 sprintf(buffer1, fitLabelFormat, coef[i]);
1118 strcat(buffer1,
"*");
1119 strcat(buffer1, xSymbol);
1121 sprintf(buffer2,
"$a%" PRId32
"$n", order[i]);
1122 strcat(buffer1, buffer2);
1125 if ((
long)(strlen(buffer) + strlen(buffer1)) > (
long)(0.95 * bufsize)) {
1126 fprintf(stderr,
"buffer overflow making fit label!\n");
1129 strcat(buffer, buffer1);
1133double rms_average(
double *x, int64_t n) {
1137 for (i = sum2 = 0; i < n; i++)
1140 return (sqrt(sum2 / n));
1143long coefficient_index(int32_t * order,
long terms,
long order_of_interest) {
1145 for (i = 0; i < terms; i++)
1146 if (order[i] == order_of_interest)
1151void checkInputFile(
SDDS_DATASET * SDDSin,
char *xName,
char **yNames,
char *xSigmaName,
char **ySigmaNames,
long numYNames) {
1156 SDDS_Bomb(
"x column doesn't exist or is nonnumeric");
1162 if (xSigmaName && !(ptr =
SDDS_FindColumn(SDDSin, FIND_NUMERIC_TYPE, xSigmaName, NULL)))
1163 SDDS_Bomb(
"x sigma column doesn't exist or is nonnumeric");
1168 for (i = 0; i < numYNames; i++) {
1170 if (!(ptr =
SDDS_FindColumn(SDDSin, FIND_NUMERIC_TYPE, ySigmaNames[i], NULL)))
1171 SDDS_Bomb(
"y sigma column doesn't exist or is nonnumeric");
1179 SDDS_DATASET *SDDSin,
char *input,
char *xName,
char **yNames,
char *xSigmaName,
1180 char **ySigmaNames,
long sigmasValid, int32_t *order,
long terms,
long isChebyshev,
1181 long numCols,
long copyParameters) {
1182 char buffer[SDDS_MAXLINE], buffer1[SDDS_MAXLINE], buffer2[SDDS_MAXLINE], buffer3[SDDS_MAXLINE];
1183 char *xUnits, *yUnits, ***coefUnits;
1188 coefUnits =
tmalloc(
sizeof(
char **) * numCols);
1189 ySymbols =
tmalloc(
sizeof(
char *) * numCols);
1190 iIntercept =
tmalloc(
sizeof(
long) * numCols);
1191 iInterceptO =
tmalloc(
sizeof(
long) * numCols);
1192 iInterceptSigma =
tmalloc(
sizeof(
long) * numCols);
1193 iInterceptSigmaO =
tmalloc(
sizeof(
long) * numCols);
1194 iSlope =
tmalloc(
sizeof(
long) * numCols);
1195 iSlopeO =
tmalloc(
sizeof(
long) * numCols);
1196 iSlopeSigma =
tmalloc(
sizeof(
long) * numCols);
1197 iSlopeSigmaO =
tmalloc(
sizeof(
long) * numCols);
1198 iCurvature =
tmalloc(
sizeof(
long) * numCols);
1199 iCurvatureO =
tmalloc(
sizeof(
long) * numCols);
1200 iCurvatureSigma =
tmalloc(
sizeof(
long) * numCols);
1201 iCurvatureSigmaO =
tmalloc(
sizeof(
long) * numCols);
1202 iChiSq =
tmalloc(
sizeof(
long) * numCols);
1203 iChiSqO =
tmalloc(
sizeof(
long) * numCols);
1204 iRmsResidual =
tmalloc(
sizeof(
long) * numCols);
1205 iRmsResidualO =
tmalloc(
sizeof(
long) * numCols);
1206 iSigLevel =
tmalloc(
sizeof(
long) * numCols);
1207 iSigLevelO =
tmalloc(
sizeof(
long) * numCols);
1208 iFitIsValid =
tmalloc(
sizeof(
long) * numCols);
1209 iFitIsValidO =
tmalloc(
sizeof(
long) * numCols);
1210 iFitLabel =
tmalloc(
sizeof(
long) * numCols);
1211 iFitLabelO =
tmalloc(
sizeof(
long) * numCols);
1212 iy =
tmalloc(
sizeof(
long) * numCols);
1213 iySigma =
tmalloc(
sizeof(
long) * numCols);
1214 iFit =
tmalloc(
sizeof(
long) * numCols);
1215 iResidual =
tmalloc(
sizeof(
long) * numCols);
1217 for (colIndex = 0; colIndex < numCols; colIndex++) {
1218 ySymbols[colIndex] = NULL;
1219 coefUnits[colIndex] =
tmalloc(
sizeof(
char *) * terms);
1220 iInterceptSigma[colIndex] = -1;
1221 iInterceptSigmaO[colIndex] = -1;
1222 iIntercept[colIndex] = -1;
1223 iInterceptO[colIndex] = -1;
1224 iInterceptSigma[colIndex] = -1;
1225 iInterceptSigmaO[colIndex] = -1;
1226 iSlope[colIndex] = -1;
1227 iSlopeO[colIndex] = -1;
1228 iSlopeSigma[colIndex] = -1;
1229 iSlopeSigmaO[colIndex] = -1;
1230 iCurvature[colIndex] = -1;
1231 iCurvatureO[colIndex] = -1;
1232 iCurvatureSigma[colIndex] = -1;
1233 iCurvatureSigmaO[colIndex] = -1;
1234 iChiSq[colIndex] = -1;
1235 iChiSqO[colIndex] = -1;
1236 iRmsResidual[colIndex] = -1;
1237 iRmsResidualO[colIndex] = -1;
1238 iSigLevel[colIndex] = -1;
1239 iSigLevelO[colIndex] = -1;
1240 iFitIsValid[colIndex] = -1;
1241 iFitIsValidO[colIndex] = -1;
1242 iFitLabel[colIndex] = -1;
1243 iFitLabelO[colIndex] = -1;
1245 iySigma[colIndex] = -1;
1246 iFit[colIndex] = -1;
1247 iResidual[colIndex] = -1;
1250 if (!
SDDS_InitializeOutput(SDDSout, SDDS_BINARY, 0, NULL,
"sddsmpfit output: fitted data", output) ||
1254 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1256 for (colIndex = 0; colIndex < numCols; colIndex++) {
1260 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1264 for (colIndex = 0; colIndex < numCols; colIndex++)
1266 ySymbols[colIndex] = yNames[colIndex];
1268 for (colIndex = 0; colIndex < numCols; colIndex++) {
1276 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1278 for (colIndex = 0; colIndex < numCols; colIndex++) {
1279 sprintf(buffer,
"%sFit", yNames[colIndex]);
1280 sprintf(buffer1,
"Fit[%s]", ySymbols[colIndex]);
1283 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1285 SDDS_Bomb(
"unable to get index of just-defined fit output column");
1287 sprintf(buffer,
"%sResidual", yNames[colIndex]);
1288 sprintf(buffer1,
"Residual[%s]", ySymbols[colIndex]);
1291 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1293 SDDS_Bomb(
"unable to get index of just-defined residual output column");
1295 if (sigmasValid && !ySigmaNames) {
1296 sprintf(buffer,
"%sSigma", yNames[colIndex]);
1298 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1301 sprintf(buffer1,
"Sigma[%s]", ySymbols[colIndex]);
1303 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1307 if (!(coefUnits[colIndex] = makeCoefficientUnits(SDDSout, xName, yNames[colIndex], order, terms)))
1308 SDDS_Bomb(
"unable to make coefficient units");
1311 if (outputInfo && !
SDDS_InitializeOutput(SDDSoutInfo, SDDS_BINARY, 0, NULL,
"sddsmpfit output: fit information", outputInfo))
1312 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1316 SDDS_DefineParameter(SDDSoutInfo,
"Basis", NULL, NULL,
"Function basis for fit", NULL,
SDDS_STRING, isChebyshev ?
"Chebyshev T polynomials" :
"ordinary polynomials") < 0)
1317 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1320 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1323 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1324 sprintf(buffer,
"%sOffset", xName);
1325 sprintf(buffer1,
"Offset of %s for fit", xName);
1327 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1328 sprintf(buffer,
"%sScale", xName);
1329 sprintf(buffer1,
"Scale factor of %s for fit", xName);
1331 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1333 for (colIndex = 0; colIndex < numCols; colIndex++) {
1335 sprintf(buffer1,
"%sCoefficient", yNames[colIndex]);
1336 sprintf(buffer2,
"%sCoefficientSigma", yNames[colIndex]);
1337 sprintf(buffer3,
"%sCoefficientUnits", yNames[colIndex]);
1340 (sigmasValid &&
SDDS_DefineColumn(SDDSoutInfo, buffer2,
"$gs$r$ba$n",
"[CoefficientUnits]",
"sigma of coefficient of term in fit", NULL,
SDDS_DOUBLE, 0) < 0))
1341 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1348 sprintf(buffer1,
"%sReducedChiSquared", yNames[colIndex]);
1349 sprintf(buffer2,
"%sRmsResidual", yNames[colIndex]);
1350 sprintf(buffer3,
"%sSignificanceLevel", yNames[colIndex]);
1353 "Reduced chi-squared of fit",
1356 (iRmsResidual[colIndex] =
1358 (iSigLevel[colIndex] =
SDDS_DefineParameter(SDDSoutInfo, buffer3, NULL, NULL,
"Probability that data is from fit function", NULL,
SDDS_DOUBLE, NULL)) < 0)
1359 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1363 sprintf(buffer,
"%sFitIsValid", yNames[colIndex]);
1365 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1369 sprintf(buffer,
"%sSddsmpfitlabel", yNames[colIndex]);
1371 if ((i = coefficient_index(order, terms, 0)) >= 0) {
1372 sprintf(buffer,
"%sIntercept", yNames[colIndex]);
1374 sprintf(buffer,
"%sInterceptSigma", yNames[colIndex]);
1376 iInterceptSigma[colIndex] =
SDDS_DefineParameter(SDDSoutInfo, buffer, buffer, coefUnits[colIndex][i],
"Sigma of intercept of fit", NULL,
SDDS_DOUBLE, NULL);
1378 sprintf(buffer,
"%sSlope", yNames[colIndex]);
1379 if ((i = coefficient_index(order, terms, 1)) >= 0) {
1382 sprintf(buffer,
"%sSlopeSigma", yNames[colIndex]);
1386 if ((i = coefficient_index(order, terms, 2)) >= 0) {
1387 sprintf(buffer,
"%sCurvature", yNames[colIndex]);
1390 sprintf(buffer,
"%sCurvatureSigma", yNames[colIndex]);
1391 iCurvatureSigma[colIndex] =
SDDS_DefineParameter(SDDSoutInfo, buffer, buffer, coefUnits[colIndex][i],
"Sigma of curvature of fit", NULL,
SDDS_DOUBLE, NULL);
1395 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1400 if (
SDDS_DefineParameter(SDDSout,
"Basis", NULL, NULL,
"Function basis for fit", NULL,
SDDS_STRING, isChebyshev ?
"Chebyshev T polynomials" :
"ordinary polynomials") < 0)
1401 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1404 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1407 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1408 sprintf(buffer,
"%sOffset", xName);
1409 sprintf(buffer1,
"Offset of %s for fit", xName);
1411 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1412 sprintf(buffer,
"%sScale", xName);
1413 sprintf(buffer1,
"Scale factor of %s for fit", xName);
1415 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1417 for (colIndex = 0; colIndex < numCols; colIndex++) {
1419 sprintf(buffer1,
"%sReducedChiSquared", yNames[colIndex]);
1420 sprintf(buffer2,
"%sRmsResidual", yNames[colIndex]);
1421 sprintf(buffer3,
"%sSignificanceLevel", yNames[colIndex]);
1424 "Reduced chi-squared of fit",
1427 (iRmsResidualO[colIndex] =
1429 (iSigLevelO[colIndex] =
SDDS_DefineParameter(SDDSout, buffer3, NULL, NULL,
"Probability that data is from fit function", NULL,
SDDS_DOUBLE, NULL)) < 0)
1430 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1434 sprintf(buffer,
"%sFitIsValid", yNames[colIndex]);
1436 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1439 sprintf(buffer,
"%sSddsmpfitlabel", yNames[colIndex]);
1441 if ((i = coefficient_index(order, terms, 0)) >= 0) {
1442 sprintf(buffer,
"%sIntercept", yNames[colIndex]);
1444 sprintf(buffer,
"%sInterceptSigma", yNames[colIndex]);
1446 iInterceptSigmaO[colIndex] =
SDDS_DefineParameter(SDDSout, buffer, buffer, coefUnits[colIndex][i],
"Sigma of intercept of fit", NULL,
SDDS_DOUBLE, NULL);
1448 sprintf(buffer,
"%sSlope", yNames[colIndex]);
1449 if ((i = coefficient_index(order, terms, 1)) >= 0) {
1452 sprintf(buffer,
"%sSlopeSigma", yNames[colIndex]);
1456 if ((i = coefficient_index(order, terms, 2)) >= 0) {
1457 sprintf(buffer,
"%sCurvature", yNames[colIndex]);
1460 sprintf(buffer,
"%sCurvatureSigma", yNames[colIndex]);
1461 iCurvatureSigmaO[colIndex] =
SDDS_DefineParameter(SDDSout, buffer, buffer, coefUnits[colIndex][i],
"Sigma of curvature of fit", NULL,
SDDS_DOUBLE, NULL);
1465 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1469 if (copyParameters) {
1471 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1473 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1477 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1482char **makeCoefficientUnits(
SDDS_DATASET * SDDSout,
char *xName,
char *yName, int32_t *order,
long terms) {
1483 char *xUnits, *yUnits, buffer[SDDS_MAXLINE];
1484 char **coefUnits = NULL;
1489 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1491 coefUnits =
tmalloc(
sizeof(*coefUnits) * terms);
1495 for (i = 0; i < terms; i++)
1496 coefUnits[i] = yUnits;
1500 for (i = 0; i < terms; i++) {
1501 if (order[i] == 0) {
1502 if (strcmp(yUnits,
"1") != 0)
1506 }
else if (strcmp(xUnits, yUnits) == 0) {
1508 sprintf(buffer,
"1/%s$a%" PRId32
"$n", xUnits, order[i] - 1);
1514 sprintf(buffer,
"%s/%s$a%" PRId32
"$n", yUnits, xUnits, order[i]);
1516 sprintf(buffer,
"%s/%s", yUnits, xUnits);
1527 SDDSout = &evalParameters->dataset;
1528 if (!
SDDS_InitializeOutput(SDDSout, SDDS_BINARY, 0, NULL,
"sddsmpfit output: evaluation of fits", evalParameters->file) ||
1530 SDDS_Bomb(
"Problem setting up evaluation file");
1531 for (i = 0; i < yNames; i++)
1533 SDDS_Bomb(
"Problem setting up evaluation file");
1535 SDDS_Bomb(
"Problem setting up evaluation file");
1538void makeEvaluationTable(
EVAL_PARAMETERS * evalParameters,
double *x, int64_t points,
double *coef, int32_t *order,
long terms,
1539 char *xName,
char **yName,
long yNames,
long iYName) {
1540 static double *xEval = NULL, *yEval = NULL;
1541 static int64_t maxEvals = 0;
1545 if (!(evalParameters->flags & EVAL_BEGIN_GIVEN) || !(evalParameters->flags & EVAL_END_GIVEN)) {
1548 if (!(evalParameters->flags & EVAL_BEGIN_GIVEN))
1549 evalParameters->begin = min;
1550 if (!(evalParameters->flags & EVAL_END_GIVEN))
1551 evalParameters->end = max;
1553 if (!(evalParameters->flags & EVAL_NUMBER_GIVEN))
1554 evalParameters->number = points;
1555 if (evalParameters->number > 1)
1556 delta = (evalParameters->end - evalParameters->begin) / (evalParameters->number - 1);
1560 if (!xEval || maxEvals < evalParameters->number) {
1561 if (!(xEval = (
double *)
SDDS_Realloc(xEval,
sizeof(*xEval) * evalParameters->number)) ||
1562 !(yEval = (
double *)
SDDS_Realloc(yEval,
sizeof(*yEval) * evalParameters->number)))
1564 maxEvals = evalParameters->number;
1567 for (i = 0; i < evalParameters->number; i++) {
1568 xEval[i] = evalParameters->begin + i * delta;
1569 yEval[i] =
eval_sum(basis_fn, coef, order, terms, xEval[i]);
1573 !
SDDS_StartPage(&evalParameters->dataset, evalParameters->number)) ||
1576 (iYName == yNames - 1 && !
SDDS_WritePage(&evalParameters->dataset)))
1577 SDDS_PrintErrors(stderr, SDDS_VERBOSE_PrintErrors | SDDS_EXIT_PrintErrors);
1580void compareOriginalToFit(
double *x,
double *y,
double **residual, int64_t points,
double *rmsResidual,
double *coef, int32_t *order,
long terms) {
1582 double residualSum2, fit;
1584 *residual =
tmalloc(
sizeof(**residual) * points);
1587 for (i = 0; i < points; i++) {
1588 fit =
eval_sum(basis_fn, coef, order, terms, x[i]);
1589 (*residual)[i] = y[i] - fit;
1590 residualSum2 += sqr((*residual)[i]);
1592 *rmsResidual = sqrt(residualSum2 / points);
SDDS (Self Describing Data Set) Data Types Definitions and Function Prototypes.
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_SetColumnFromLongs(SDDS_DATASET *SDDS_dataset, int32_t mode, int32_t *data, int64_t rows,...)
Sets the values for a single data column using long integer numbers.
int32_t SDDS_ChangeColumnInformation(SDDS_DATASET *SDDS_dataset, char *field_name, void *memory, int32_t mode,...)
Modifies a specific field in a column definition within the SDDS dataset.
int32_t SDDS_GetColumnInformation(SDDS_DATASET *SDDS_dataset, char *field_name, void *memory, int32_t mode,...)
Retrieves information about a specified column in the SDDS dataset.
int32_t SDDS_InitializeOutput(SDDS_DATASET *SDDS_dataset, int32_t data_mode, int32_t lines_per_row, const char *description, const char *contents, const char *filename)
Initializes the SDDS output dataset.
int32_t SDDS_WritePage(SDDS_DATASET *SDDS_dataset)
Writes the current data table to the output file.
int32_t SDDS_DefineColumn(SDDS_DATASET *SDDS_dataset, const char *name, const char *symbol, const char *units, const char *description, const char *format_string, int32_t type, int32_t field_length)
Defines a data column within the SDDS dataset.
int32_t SDDS_WriteLayout(SDDS_DATASET *SDDS_dataset)
Writes the SDDS layout header to the output file.
int32_t SDDS_DefineParameter(SDDS_DATASET *SDDS_dataset, const char *name, const char *symbol, const char *units, const char *description, const char *format_string, int32_t type, char *fixed_value)
Defines a data parameter with a fixed string value.
int32_t SDDS_TransferColumnDefinition(SDDS_DATASET *target, SDDS_DATASET *source, char *name, char *newName)
Transfers a column definition from a source dataset to a target dataset.
int32_t SDDS_TransferAllParameterDefinitions(SDDS_DATASET *SDDS_target, SDDS_DATASET *SDDS_source, uint32_t mode)
Transfers all parameter definitions from a source dataset to a target dataset.
int32_t SDDS_GetColumnIndex(SDDS_DATASET *SDDS_dataset, char *name)
Retrieves the index of a named column in the SDDS dataset.
int32_t SDDS_CheckColumn(SDDS_DATASET *SDDS_dataset, char *name, char *units, int32_t type, FILE *fp_message)
Checks if a column exists in the SDDS dataset with the specified name, units, and type.
char ** SDDS_GetColumnNames(SDDS_DATASET *SDDS_dataset, int32_t *number)
Retrieves the names of all columns in the SDDS dataset.
void SDDS_PrintErrors(FILE *fp, int32_t mode)
Prints recorded error messages to a specified file stream.
void SDDS_RegisterProgramName(const char *name)
Registers the executable program name for use in error messages.
int32_t SDDS_NumberOfErrors()
Retrieves the number of errors recorded by SDDS library routines.
int32_t SDDS_StringIsBlank(char *s)
Checks if a string is blank (contains only whitespace characters).
char * SDDS_FindColumn(SDDS_DATASET *SDDS_dataset, int32_t mode,...)
Finds the first column in the SDDS dataset that matches the specified criteria.
void SDDS_Bomb(char *message)
Terminates the program after printing an error message and recorded errors.
int32_t SDDS_CopyString(char **target, const char *source)
Copies a source string to a target string with memory allocation.
void * SDDS_Realloc(void *old_ptr, size_t new_size)
Reallocates memory to a new size.
#define SDDS_STRING
Identifier for the string data type.
#define SDDS_LONG
Identifier for the signed 32-bit integer data type.
#define SDDS_CHARACTER
Identifier for the character data type.
#define SDDS_ANY_NUMERIC_TYPE
Special identifier used by SDDS_Check*() routines to accept any numeric type.
#define SDDS_DOUBLE
Identifier for the double data type.
#define SDDS_LONG64
Identifier for the signed 64-bit integer data type.
void * tmalloc(uint64_t size_of_block)
Allocates a memory block of the specified size with zero initialization.
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.
double eval_sum(double(*fn)(double x, long ord), double *coef, int32_t *order, long n_coefs, double x0)
Evaluate a sum of basis functions.
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.
double dtcheby(double x, long n)
Evaluate the derivative of the Chebyshev polynomial T_n(x).
void set_argument_offset(double offset)
Set the offset applied to the input argument of basis functions.
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.