xref: /petsc/src/sys/utils/mathfit.c (revision 792fecdfe9134cce4d631112660ddd34f063bc17)
15f80ce2aSJacob Faibussowitsch #include <petsc/private/petscimpl.h>
2bebf13c0SMatthew G. Knepley #include <petscblaslapack.h>
3bebf13c0SMatthew G. Knepley 
4bebf13c0SMatthew G. Knepley /*@C
5bebf13c0SMatthew G. Knepley   PetscLinearRegression - Gives the best least-squares linear fit to some x-y data points
6bebf13c0SMatthew G. Knepley 
7bebf13c0SMatthew G. Knepley   Input Parameters:
8bebf13c0SMatthew G. Knepley + n - The number of points
9bebf13c0SMatthew G. Knepley . x - The x-values
10bebf13c0SMatthew G. Knepley - y - The y-values
11bebf13c0SMatthew G. Knepley 
12bebf13c0SMatthew G. Knepley   Output Parameters:
13bebf13c0SMatthew G. Knepley + slope     - The slope of the best-fit line
14bebf13c0SMatthew G. Knepley - intercept - The y-intercept of the best-fit line
15bebf13c0SMatthew G. Knepley 
16bebf13c0SMatthew G. Knepley   Level: intermediate
17bebf13c0SMatthew G. Knepley 
18db781477SPatrick Sanan .seealso: `PetscConvEstGetConvRate()`
19bebf13c0SMatthew G. Knepley @*/
20bebf13c0SMatthew G. Knepley PetscErrorCode PetscLinearRegression(PetscInt n, const PetscReal x[], const PetscReal y[], PetscReal *slope, PetscReal *intercept)
21bebf13c0SMatthew G. Knepley {
22bebf13c0SMatthew G. Knepley   PetscScalar  H[4];
23bebf13c0SMatthew G. Knepley   PetscReal   *X, *Y, beta[2];
24bebf13c0SMatthew G. Knepley 
25bebf13c0SMatthew G. Knepley   PetscFunctionBegin;
265f80ce2aSJacob Faibussowitsch   if (n) {
275f80ce2aSJacob Faibussowitsch     PetscValidRealPointer(x,2);
285f80ce2aSJacob Faibussowitsch     PetscValidRealPointer(y,3);
295f80ce2aSJacob Faibussowitsch   }
305f80ce2aSJacob Faibussowitsch   PetscValidRealPointer(slope,4);
315f80ce2aSJacob Faibussowitsch   PetscValidRealPointer(intercept,5);
329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n*2, &X, n*2, &Y));
335f80ce2aSJacob Faibussowitsch   for (PetscInt k = 0; k < n; ++k) {
34bebf13c0SMatthew G. Knepley     /* X[n,2] = [1, x] */
35bebf13c0SMatthew G. Knepley     X[k*2+0] = 1.0;
36bebf13c0SMatthew G. Knepley     X[k*2+1] = x[k];
37bebf13c0SMatthew G. Knepley   }
38bebf13c0SMatthew G. Knepley   /* H = X^T X */
395f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < 2; ++i) {
405f80ce2aSJacob Faibussowitsch     for (PetscInt j = 0; j < 2; ++j) {
41bebf13c0SMatthew G. Knepley       H[i*2+j] = 0.0;
425f80ce2aSJacob Faibussowitsch       for (PetscInt k = 0; k < n; ++k) H[i*2+j] += X[k*2+i] * X[k*2+j];
43bebf13c0SMatthew G. Knepley     }
44bebf13c0SMatthew G. Knepley   }
45bebf13c0SMatthew G. Knepley   /* H = (X^T X)^{-1} */
46bebf13c0SMatthew G. Knepley   {
47bebf13c0SMatthew G. Knepley     PetscBLASInt two = 2, ipiv[2], info;
48bebf13c0SMatthew G. Knepley     PetscScalar  work[2];
49bebf13c0SMatthew G. Knepley 
509566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
51*792fecdfSBarry Smith     PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&two, &two, H, &two, ipiv, &info));
52*792fecdfSBarry Smith     PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&two, H, &two, ipiv, work, &two, &info));
539566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPop());
54bebf13c0SMatthew G. Knepley   }
55bebf13c0SMatthew G. Knepley     /* Y = H X^T */
565f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < 2; ++i) {
575f80ce2aSJacob Faibussowitsch     for (PetscInt k = 0; k < n; ++k) {
58bebf13c0SMatthew G. Knepley       Y[i*n+k] = 0.0;
595f80ce2aSJacob Faibussowitsch       for (PetscInt j = 0; j < 2; ++j) Y[i*n+k] += PetscRealPart(H[i*2+j]) * X[k*2+j];
60bebf13c0SMatthew G. Knepley     }
61bebf13c0SMatthew G. Knepley   }
62bebf13c0SMatthew G. Knepley   /* beta = Y error = [y-intercept, slope] */
635f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < 2; ++i) {
64bebf13c0SMatthew G. Knepley     beta[i] = 0.0;
655f80ce2aSJacob Faibussowitsch     for (PetscInt k = 0; k < n; ++k) beta[i] += Y[i*n+k] * y[k];
66bebf13c0SMatthew G. Knepley   }
679566063dSJacob Faibussowitsch   PetscCall(PetscFree2(X, Y));
68bebf13c0SMatthew G. Knepley   *intercept = beta[0];
69bebf13c0SMatthew G. Knepley   *slope     = beta[1];
70bebf13c0SMatthew G. Knepley   PetscFunctionReturn(0);
71bebf13c0SMatthew G. Knepley }
72