xref: /petsc/include/petscregressor.h (revision 34b254c57d2aa195261fbc0db2d1455fb6d091da)
1*34b254c5SRichard Tran Mills #pragma once
2*34b254c5SRichard Tran Mills 
3*34b254c5SRichard Tran Mills #include <petsctao.h>
4*34b254c5SRichard Tran Mills 
5*34b254c5SRichard Tran Mills /* MANSEC = ML */
6*34b254c5SRichard Tran Mills /* SUBMANSEC = PetscRegressor */
7*34b254c5SRichard Tran Mills 
8*34b254c5SRichard Tran Mills /*S
9*34b254c5SRichard Tran Mills    PetscRegressor - Abstract PETSc object that manages regression and classification problems
10*34b254c5SRichard Tran Mills 
11*34b254c5SRichard Tran Mills    Level: beginner
12*34b254c5SRichard Tran Mills 
13*34b254c5SRichard Tran Mills    Note:
14*34b254c5SRichard Tran Mills    We have slightly abused the term "regressor" in the naming of this component of PETSc.
15*34b254c5SRichard Tran Mills    Statisticians would say that we are doing "regression", and a "regressor", in this context, strictly means an
16*34b254c5SRichard Tran Mills    independent (or "predictor") variable in the regression analysis. However, "regressor" has taken on an informal
17*34b254c5SRichard Tran Mills    meaning in the machine-learning community of something along the lines of "algorithm or implementation used to fit
18*34b254c5SRichard Tran Mills    a regression model". Examples are `MLPRegressor` (multi-layer perceptron regressor) or `RandomForestRegressor`
19*34b254c5SRichard Tran Mills    from the scikit-learn toolkit (which is itself not consistent about the use of the term "regressor", since it has a
20*34b254c5SRichard Tran Mills    `LinearRegression` component instead of a `LinearRegressor` component).
21*34b254c5SRichard Tran Mills 
22*34b254c5SRichard Tran Mills .seealso: `PetscRegressorCreate()`, `PetscRegressorSetType()`, `PetscRegressorType`, `PetscRegressorDestroy()`
23*34b254c5SRichard Tran Mills S*/
24*34b254c5SRichard Tran Mills typedef struct _p_PetscRegressor *PetscRegressor;
25*34b254c5SRichard Tran Mills 
26*34b254c5SRichard Tran Mills /*J
27*34b254c5SRichard Tran Mills   PetscRegressorType - String with the name of a PETSc regression method.
28*34b254c5SRichard Tran Mills 
29*34b254c5SRichard Tran Mills   Level: beginner
30*34b254c5SRichard Tran Mills 
31*34b254c5SRichard Tran Mills .seealso: [](ch_regressor), `PetscRegressorSetType()`, `PetscRegressor`, `PetscRegressorRegister()`, `PetscRegressorCreate()`, `PetscRegressorSetFromOptions()`
32*34b254c5SRichard Tran Mills J*/
33*34b254c5SRichard Tran Mills typedef const char *PetscRegressorType;
34*34b254c5SRichard Tran Mills #define PETSCREGRESSORLINEAR "linear"
35*34b254c5SRichard Tran Mills 
36*34b254c5SRichard Tran Mills /*E
37*34b254c5SRichard Tran Mills   PetscRegressorLinearType - Type of linear regression
38*34b254c5SRichard Tran Mills 
39*34b254c5SRichard Tran Mills   Values:
40*34b254c5SRichard Tran Mills +  `REGRESSOR_LINEAR_OLS`    - ordinary least square
41*34b254c5SRichard Tran Mills .  `REGRESSOR_LINEAR_LASSO`  - lasso
42*34b254c5SRichard Tran Mills -  `REGRESSOR_LINEAR_RIDGE`  - ridge
43*34b254c5SRichard Tran Mills 
44*34b254c5SRichard Tran Mills   Level: advanced
45*34b254c5SRichard Tran Mills 
46*34b254c5SRichard Tran Mills .seealso: `PetscRegressor`, `PETSCREGRESSORLINEAR`
47*34b254c5SRichard Tran Mills E*/
48*34b254c5SRichard Tran Mills 
49*34b254c5SRichard Tran Mills typedef enum {
50*34b254c5SRichard Tran Mills   REGRESSOR_LINEAR_OLS,
51*34b254c5SRichard Tran Mills   REGRESSOR_LINEAR_LASSO,
52*34b254c5SRichard Tran Mills   REGRESSOR_LINEAR_RIDGE
53*34b254c5SRichard Tran Mills } PetscRegressorLinearType;
54*34b254c5SRichard Tran Mills PETSC_EXTERN const char *const PetscRegressorLinearTypes[];
55*34b254c5SRichard Tran Mills 
56*34b254c5SRichard Tran Mills PETSC_EXTERN PetscFunctionList PetscRegressorList;
57*34b254c5SRichard Tran Mills PETSC_EXTERN PetscClassId      PETSCREGRESSOR_CLASSID;
58*34b254c5SRichard Tran Mills 
59*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorInitializePackage(void);
60*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorFinalizePackage(void);
61*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorRegister(const char[], PetscErrorCode (*)(PetscRegressor));
62*34b254c5SRichard Tran Mills 
63*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorCreate(MPI_Comm, PetscRegressor *);
64*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorReset(PetscRegressor);
65*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorDestroy(PetscRegressor *);
66*34b254c5SRichard Tran Mills 
67*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorSetOptionsPrefix(PetscRegressor, const char[]);
68*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorAppendOptionsPrefix(PetscRegressor, const char[]);
69*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorGetOptionsPrefix(PetscRegressor, const char *[]);
70*34b254c5SRichard Tran Mills 
71*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorSetType(PetscRegressor, PetscRegressorType);
72*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorGetType(PetscRegressor, PetscRegressorType *);
73*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorSetRegularizerWeight(PetscRegressor, PetscReal);
74*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorSetUp(PetscRegressor);
75*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorSetFromOptions(PetscRegressor);
76*34b254c5SRichard Tran Mills 
77*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorView(PetscRegressor, PetscViewer);
78*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorViewFromOptions(PetscRegressor, PetscObject, const char[]);
79*34b254c5SRichard Tran Mills 
80*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorFit(PetscRegressor, Mat, Vec);
81*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorPredict(PetscRegressor, Mat, Vec);
82*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorGetTao(PetscRegressor, Tao *);
83*34b254c5SRichard Tran Mills 
84*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearSetFitIntercept(PetscRegressor, PetscBool);
85*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearSetUseKSP(PetscRegressor, PetscBool);
86*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearGetKSP(PetscRegressor, KSP *);
87*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearGetCoefficients(PetscRegressor, Vec *);
88*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearGetIntercept(PetscRegressor, PetscScalar *);
89*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearSetType(PetscRegressor, PetscRegressorLinearType);
90*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearGetType(PetscRegressor, PetscRegressorLinearType *);
91