xref: /petsc/include/petscregressor.h (revision 10632adc56aeb24e1cd352dce800f73ee42de22b)
134b254c5SRichard Tran Mills #pragma once
234b254c5SRichard Tran Mills 
334b254c5SRichard Tran Mills #include <petsctao.h>
434b254c5SRichard Tran Mills 
534b254c5SRichard Tran Mills /* MANSEC = ML */
634b254c5SRichard Tran Mills /* SUBMANSEC = PetscRegressor */
734b254c5SRichard Tran Mills 
834b254c5SRichard Tran Mills /*S
934b254c5SRichard Tran Mills    PetscRegressor - Abstract PETSc object that manages regression and classification problems
1034b254c5SRichard Tran Mills 
1134b254c5SRichard Tran Mills    Level: beginner
1234b254c5SRichard Tran Mills 
13789736e1SBarry Smith    Notes:
14789736e1SBarry Smith    For linear problems `PetscRegressor` supports ordinary least squares, lasso, and ridge regression using the `PetscRegressorType` of `PETSCREGRESSORLINEAR`
15789736e1SBarry Smith    and `PetscRegressorLinearType` of `REGRESSOR_LINEAR_OLS`, `REGRESSOR_LINEAR_LASSO`, and `REGRESSOR_LINEAR_RIDGE`.
16789736e1SBarry Smith 
1734b254c5SRichard Tran Mills    We have slightly abused the term "regressor" in the naming of this component of PETSc.
1834b254c5SRichard Tran Mills    Statisticians would say that we are doing "regression", and a "regressor", in this context, strictly means an
1934b254c5SRichard Tran Mills    independent (or "predictor") variable in the regression analysis. However, "regressor" has taken on an informal
2034b254c5SRichard Tran Mills    meaning in the machine-learning community of something along the lines of "algorithm or implementation used to fit
2134b254c5SRichard Tran Mills    a regression model". Examples are `MLPRegressor` (multi-layer perceptron regressor) or `RandomForestRegressor`
2234b254c5SRichard Tran Mills    from the scikit-learn toolkit (which is itself not consistent about the use of the term "regressor", since it has a
2334b254c5SRichard Tran Mills    `LinearRegression` component instead of a `LinearRegressor` component).
2434b254c5SRichard Tran Mills 
25789736e1SBarry Smith .seealso: `PetscRegressorCreate()`, `PetscRegressorLinearType`, `PetscRegressorSetType()`, `PetscRegressorType`, `PetscRegressorDestroy()`,
26789736e1SBarry Smith           `PETSCREGRESSORLINEAR`, `PetscRegressorLinearType`, `REGRESSOR_LINEAR_OLS`, `REGRESSOR_LINEAR_LASSO`, `REGRESSOR_LINEAR_RIDGE`.
2734b254c5SRichard Tran Mills S*/
2834b254c5SRichard Tran Mills typedef struct _p_PetscRegressor *PetscRegressor;
2934b254c5SRichard Tran Mills 
3034b254c5SRichard Tran Mills /*J
3134b254c5SRichard Tran Mills   PetscRegressorType - String with the name of a PETSc regression method.
3234b254c5SRichard Tran Mills 
3334b254c5SRichard Tran Mills   Level: beginner
3434b254c5SRichard Tran Mills 
35789736e1SBarry Smith .seealso: [](ch_regressor), `PetscRegressorSetType()`, `PetscRegressor`, `PetscRegressorRegister()`, `PetscRegressorCreate()`, `PetscRegressorSetFromOptions()`,
36789736e1SBarry Smith           `PETSCREGRESSORLINEAR`
3734b254c5SRichard Tran Mills J*/
3834b254c5SRichard Tran Mills typedef const char *PetscRegressorType;
3934b254c5SRichard Tran Mills #define PETSCREGRESSORLINEAR "linear"
4034b254c5SRichard Tran Mills 
4134b254c5SRichard Tran Mills /*E
4234b254c5SRichard Tran Mills   PetscRegressorLinearType - Type of linear regression
4334b254c5SRichard Tran Mills 
4434b254c5SRichard Tran Mills   Values:
45789736e1SBarry Smith +  `REGRESSOR_LINEAR_OLS`    - ordinary least squares
4634b254c5SRichard Tran Mills .  `REGRESSOR_LINEAR_LASSO`  - lasso
4734b254c5SRichard Tran Mills -  `REGRESSOR_LINEAR_RIDGE`  - ridge
4834b254c5SRichard Tran Mills 
4934b254c5SRichard Tran Mills   Level: advanced
5034b254c5SRichard Tran Mills 
51*0664cd31SRichard Tran Mills   Note:
52*0664cd31SRichard Tran Mills   One can perform binary classification using the ridge regressor type by converting labels into the
53*0664cd31SRichard Tran Mills   values -1 and +1, corresponding to the two classes, and then performing a ridge regression.
54*0664cd31SRichard Tran Mills   Observations with a negative prediction value are then placed in the -1 class, while those with positive values
55*0664cd31SRichard Tran Mills   are placed in the +1 class.
56*0664cd31SRichard Tran Mills   This is the approach used in the RidgeClassifer implementation provided by the scikit-learn library.
57*0664cd31SRichard Tran Mills 
5834b254c5SRichard Tran Mills .seealso: `PetscRegressor`, `PETSCREGRESSORLINEAR`
5934b254c5SRichard Tran Mills E*/
6034b254c5SRichard Tran Mills 
6134b254c5SRichard Tran Mills typedef enum {
6234b254c5SRichard Tran Mills   REGRESSOR_LINEAR_OLS,
6334b254c5SRichard Tran Mills   REGRESSOR_LINEAR_LASSO,
6434b254c5SRichard Tran Mills   REGRESSOR_LINEAR_RIDGE
6534b254c5SRichard Tran Mills } PetscRegressorLinearType;
6634b254c5SRichard Tran Mills PETSC_EXTERN const char *const PetscRegressorLinearTypes[];
6734b254c5SRichard Tran Mills 
6834b254c5SRichard Tran Mills PETSC_EXTERN PetscFunctionList PetscRegressorList;
6934b254c5SRichard Tran Mills PETSC_EXTERN PetscClassId      PETSCREGRESSOR_CLASSID;
7034b254c5SRichard Tran Mills 
7134b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorInitializePackage(void);
7234b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorFinalizePackage(void);
7334b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorRegister(const char[], PetscErrorCode (*)(PetscRegressor));
7434b254c5SRichard Tran Mills 
7534b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorCreate(MPI_Comm, PetscRegressor *);
7634b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorReset(PetscRegressor);
7734b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorDestroy(PetscRegressor *);
7834b254c5SRichard Tran Mills 
7934b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorSetOptionsPrefix(PetscRegressor, const char[]);
8034b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorAppendOptionsPrefix(PetscRegressor, const char[]);
8134b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorGetOptionsPrefix(PetscRegressor, const char *[]);
8234b254c5SRichard Tran Mills 
8334b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorSetType(PetscRegressor, PetscRegressorType);
8434b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorGetType(PetscRegressor, PetscRegressorType *);
8534b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorSetRegularizerWeight(PetscRegressor, PetscReal);
8634b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorSetUp(PetscRegressor);
8734b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorSetFromOptions(PetscRegressor);
8834b254c5SRichard Tran Mills 
8934b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorView(PetscRegressor, PetscViewer);
9034b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorViewFromOptions(PetscRegressor, PetscObject, const char[]);
9134b254c5SRichard Tran Mills 
9234b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorFit(PetscRegressor, Mat, Vec);
9334b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorPredict(PetscRegressor, Mat, Vec);
9434b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorGetTao(PetscRegressor, Tao *);
9534b254c5SRichard Tran Mills 
9634b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearSetFitIntercept(PetscRegressor, PetscBool);
9734b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearSetUseKSP(PetscRegressor, PetscBool);
9834b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearGetKSP(PetscRegressor, KSP *);
9934b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearGetCoefficients(PetscRegressor, Vec *);
10034b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearGetIntercept(PetscRegressor, PetscScalar *);
10134b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearSetType(PetscRegressor, PetscRegressorLinearType);
10234b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearGetType(PetscRegressor, PetscRegressorLinearType *);
103