1*bebf13c0SMatthew G. Knepley #include <petscsys.h> 2*bebf13c0SMatthew G. Knepley #include <petscblaslapack.h> 3*bebf13c0SMatthew G. Knepley 4*bebf13c0SMatthew G. Knepley /*@C 5*bebf13c0SMatthew G. Knepley PetscLinearRegression - Gives the best least-squares linear fit to some x-y data points 6*bebf13c0SMatthew G. Knepley 7*bebf13c0SMatthew G. Knepley Input Parameters: 8*bebf13c0SMatthew G. Knepley + n - The number of points 9*bebf13c0SMatthew G. Knepley . x - The x-values 10*bebf13c0SMatthew G. Knepley - y - The y-values 11*bebf13c0SMatthew G. Knepley 12*bebf13c0SMatthew G. Knepley Output Parameters: 13*bebf13c0SMatthew G. Knepley + slope - The slope of the best-fit line 14*bebf13c0SMatthew G. Knepley - intercept - The y-intercept of the best-fit line 15*bebf13c0SMatthew G. Knepley 16*bebf13c0SMatthew G. Knepley Level: intermediate 17*bebf13c0SMatthew G. Knepley 18*bebf13c0SMatthew G. Knepley .seealso: PetscConvEstGetConvRate() 19*bebf13c0SMatthew G. Knepley @*/ 20*bebf13c0SMatthew G. Knepley PetscErrorCode PetscLinearRegression(PetscInt n, const PetscReal x[], const PetscReal y[], PetscReal *slope, PetscReal *intercept) 21*bebf13c0SMatthew G. Knepley { 22*bebf13c0SMatthew G. Knepley PetscScalar H[4]; 23*bebf13c0SMatthew G. Knepley PetscReal *X, *Y, beta[2]; 24*bebf13c0SMatthew G. Knepley PetscInt i, j, k; 25*bebf13c0SMatthew G. Knepley PetscErrorCode ierr; 26*bebf13c0SMatthew G. Knepley 27*bebf13c0SMatthew G. Knepley PetscFunctionBegin; 28*bebf13c0SMatthew G. Knepley *slope = *intercept = 0.0; 29*bebf13c0SMatthew G. Knepley ierr = PetscMalloc2(n*2, &X, n*2, &Y);CHKERRQ(ierr); 30*bebf13c0SMatthew G. Knepley for (k = 0; k < n; ++k) { 31*bebf13c0SMatthew G. Knepley /* X[n,2] = [1, x] */ 32*bebf13c0SMatthew G. Knepley X[k*2+0] = 1.0; 33*bebf13c0SMatthew G. Knepley X[k*2+1] = x[k]; 34*bebf13c0SMatthew G. Knepley } 35*bebf13c0SMatthew G. Knepley /* H = X^T X */ 36*bebf13c0SMatthew G. Knepley for (i = 0; i < 2; ++i) { 37*bebf13c0SMatthew G. Knepley for (j = 0; j < 2; ++j) { 38*bebf13c0SMatthew G. Knepley H[i*2+j] = 0.0; 39*bebf13c0SMatthew G. Knepley for (k = 0; k < n; ++k) { 40*bebf13c0SMatthew G. Knepley H[i*2+j] += X[k*2+i] * X[k*2+j]; 41*bebf13c0SMatthew G. Knepley } 42*bebf13c0SMatthew G. Knepley } 43*bebf13c0SMatthew G. Knepley } 44*bebf13c0SMatthew G. Knepley /* H = (X^T X)^{-1} */ 45*bebf13c0SMatthew G. Knepley { 46*bebf13c0SMatthew G. Knepley PetscBLASInt two = 2, ipiv[2], info; 47*bebf13c0SMatthew G. Knepley PetscScalar work[2]; 48*bebf13c0SMatthew G. Knepley 49*bebf13c0SMatthew G. Knepley ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 50*bebf13c0SMatthew G. Knepley PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&two, &two, H, &two, ipiv, &info)); 51*bebf13c0SMatthew G. Knepley PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&two, H, &two, ipiv, work, &two, &info)); 52*bebf13c0SMatthew G. Knepley ierr = PetscFPTrapPop();CHKERRQ(ierr); 53*bebf13c0SMatthew G. Knepley } 54*bebf13c0SMatthew G. Knepley /* Y = H X^T */ 55*bebf13c0SMatthew G. Knepley for (i = 0; i < 2; ++i) { 56*bebf13c0SMatthew G. Knepley for (k = 0; k < n; ++k) { 57*bebf13c0SMatthew G. Knepley Y[i*n+k] = 0.0; 58*bebf13c0SMatthew G. Knepley for (j = 0; j < 2; ++j) { 59*bebf13c0SMatthew G. Knepley Y[i*n+k] += PetscRealPart(H[i*2+j]) * X[k*2+j]; 60*bebf13c0SMatthew G. Knepley } 61*bebf13c0SMatthew G. Knepley } 62*bebf13c0SMatthew G. Knepley } 63*bebf13c0SMatthew G. Knepley /* beta = Y error = [y-intercept, slope] */ 64*bebf13c0SMatthew G. Knepley for (i = 0; i < 2; ++i) { 65*bebf13c0SMatthew G. Knepley beta[i] = 0.0; 66*bebf13c0SMatthew G. Knepley for (k = 0; k < n; ++k) { 67*bebf13c0SMatthew G. Knepley beta[i] += Y[i*n+k] * y[k]; 68*bebf13c0SMatthew G. Knepley } 69*bebf13c0SMatthew G. Knepley } 70*bebf13c0SMatthew G. Knepley ierr = PetscFree2(X, Y);CHKERRQ(ierr); 71*bebf13c0SMatthew G. Knepley *intercept = beta[0]; 72*bebf13c0SMatthew G. Knepley *slope = beta[1]; 73*bebf13c0SMatthew G. Knepley PetscFunctionReturn(0); 74*bebf13c0SMatthew G. Knepley } 75