1 /* Example inspired by the toy example in https://www.r-bloggers.com/2020/06/understanding-lasso-and-ridge-regression-2/ 2 * blog post by Dr. Atakan Ekiz. 3 * Here we wish to predict the number of shark attacks (that is, this number is our response variable), 4 * using the following predictor variables: 5 * - percentage of swimmers who watched the movie Jaws 6 * - the number of swimmers in the water 7 * - the average temperature of the day 8 * - the price of your favorite tech stock of the day (totally uncorrelated variable) */ 9 10 static char help[] = "Tests basic creation and destruction of PetscRegressor objects.\n\n"; 11 12 #include <petscregressor.h> 13 14 int main(int argc, char **args) 15 { 16 PetscRegressor regressor; 17 PetscMPIInt rank; 18 Mat X; 19 Vec y, y_predicted, coefficients; 20 PetscScalar intercept; 21 22 PetscScalar y_array[20] = {98, 53, 39, 127, 73, 42, 71, 61, 83, 74, 85, 82, 62, 60, 43, 69, 67, 69, 85, 3}; // Number of shark attacks 23 24 PetscScalar X_array[80] = {37.92934, 513, 92.89899, 137.2139, // % watched Jaws, #swimmers, temperature, stock price 25 52.77429, 451, 87.86271, 145.7987, // 26 60.84441, 456, 88.28927, 149.7299, // 27 26.54302, 546, 89.43875, 147.1180, // 28 54.29125, 431, 88.01132, 124.3068, // 29 55.06056, 355, 88.06297, 114.1730, // 30 44.25260, 557, 87.78536, 112.5773, // 31 44.53368, 398, 87.49603, 125.1628, // 32 44.35548, 498, 88.95234, 124.8483, // 33 41.09962, 406, 89.00630, 115.9223, // 34 45.22807, 610, 86.38794, 148.1111, // 35 40.01614, 452, 88.83585, 131.7050, // 36 42.23746, 429, 87.78222, 106.3717, // 37 50.64459, 450, 87.97008, 121.1523, // 38 59.59494, 337, 89.67538, 145.7158, // 39 48.89715, 383, 91.12611, 123.3896, // 40 44.88990, 282, 93.29563, 145.4085, // 41 40.88805, 366, 88.45329, 129.8872, // 42 41.62828, 471, 93.21182, 131.5871, // 43 74.15835, 453, 87.68438, 143.4579}; 44 45 PetscInt rows_ix[20] = {0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19}; 46 PetscInt cols_ix[4] = {0, 1, 2, 3}; 47 48 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 49 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 50 51 PetscCall(VecCreate(PETSC_COMM_WORLD, &y)); 52 PetscCall(VecSetSizes(y, PETSC_DECIDE, 20)); 53 PetscCall(VecSetFromOptions(y)); 54 PetscCall(VecDuplicate(y, &y_predicted)); 55 PetscCall(MatCreate(PETSC_COMM_WORLD, &X)); 56 PetscCall(MatSetSizes(X, PETSC_DECIDE, PETSC_DECIDE, 20, 4)); 57 PetscCall(MatSetFromOptions(X)); 58 PetscCall(MatSetUp(X)); 59 60 if (!rank) { 61 PetscCall(VecSetValues(y, 20, rows_ix, y_array, INSERT_VALUES)); 62 PetscCall(MatSetValues(X, 20, rows_ix, 4, cols_ix, X_array, ADD_VALUES)); 63 } 64 PetscCall(VecAssemblyBegin(y)); 65 PetscCall(VecAssemblyEnd(y)); 66 PetscCall(MatAssemblyBegin(X, MAT_FINAL_ASSEMBLY)); 67 PetscCall(MatAssemblyEnd(X, MAT_FINAL_ASSEMBLY)); 68 69 PetscCall(PetscRegressorCreate(PETSC_COMM_WORLD, ®ressor)); 70 PetscCall(PetscRegressorSetType(regressor, PETSCREGRESSORLINEAR)); 71 PetscRegressorSetFromOptions(regressor); 72 PetscCall(PetscRegressorFit(regressor, X, y)); 73 PetscCall(PetscRegressorPredict(regressor, X, y_predicted)); 74 PetscCall(PetscRegressorLinearGetIntercept(regressor, &intercept)); 75 PetscCall(PetscRegressorLinearGetCoefficients(regressor, &coefficients)); 76 77 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Intercept is %lf\n", intercept)); 78 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Coefficients are\n")); 79 PetscCall(VecView(coefficients, PETSC_VIEWER_STDOUT_WORLD)); 80 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Predicted values are\n")); 81 PetscCall(VecView(y_predicted, PETSC_VIEWER_STDOUT_WORLD)); 82 83 PetscCall(PetscRegressorDestroy(®ressor)); 84 85 PetscCall(PetscFinalize()); 86 return 0; 87 } 88