1*34b254c5SRichard Tran Mills #include <../src/ml/regressor/impls/linear/linearimpl.h> /*I "petscregressor.h" I*/ 2*34b254c5SRichard Tran Mills #include <../src/tao/leastsquares/impls/brgn/brgn.h> /*I "petsctao.h" I*/ 3*34b254c5SRichard Tran Mills 4*34b254c5SRichard Tran Mills const char *const PetscRegressorLinearTypes[] = {"ols", "lasso", "ridge", "RegressorLinearType", "REGRESSOR_LINEAR_", NULL}; 5*34b254c5SRichard Tran Mills 6*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorLinearSetFitIntercept_Linear(PetscRegressor regressor, PetscBool flg) 7*34b254c5SRichard Tran Mills { 8*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 9*34b254c5SRichard Tran Mills 10*34b254c5SRichard Tran Mills PetscFunctionBegin; 11*34b254c5SRichard Tran Mills linear->fit_intercept = flg; 12*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 13*34b254c5SRichard Tran Mills } 14*34b254c5SRichard Tran Mills 15*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorLinearSetType_Linear(PetscRegressor regressor, PetscRegressorLinearType type) 16*34b254c5SRichard Tran Mills { 17*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 18*34b254c5SRichard Tran Mills 19*34b254c5SRichard Tran Mills PetscFunctionBegin; 20*34b254c5SRichard Tran Mills linear->type = type; 21*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 22*34b254c5SRichard Tran Mills } 23*34b254c5SRichard Tran Mills 24*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorLinearGetType_Linear(PetscRegressor regressor, PetscRegressorLinearType *type) 25*34b254c5SRichard Tran Mills { 26*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 27*34b254c5SRichard Tran Mills 28*34b254c5SRichard Tran Mills PetscFunctionBegin; 29*34b254c5SRichard Tran Mills *type = linear->type; 30*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 31*34b254c5SRichard Tran Mills } 32*34b254c5SRichard Tran Mills 33*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorLinearGetIntercept_Linear(PetscRegressor regressor, PetscScalar *intercept) 34*34b254c5SRichard Tran Mills { 35*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 36*34b254c5SRichard Tran Mills 37*34b254c5SRichard Tran Mills PetscFunctionBegin; 38*34b254c5SRichard Tran Mills *intercept = linear->intercept; 39*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 40*34b254c5SRichard Tran Mills } 41*34b254c5SRichard Tran Mills 42*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorLinearGetCoefficients_Linear(PetscRegressor regressor, Vec *coefficients) 43*34b254c5SRichard Tran Mills { 44*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 45*34b254c5SRichard Tran Mills 46*34b254c5SRichard Tran Mills PetscFunctionBegin; 47*34b254c5SRichard Tran Mills *coefficients = linear->coefficients; 48*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 49*34b254c5SRichard Tran Mills } 50*34b254c5SRichard Tran Mills 51*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorLinearGetKSP_Linear(PetscRegressor regressor, KSP *ksp) 52*34b254c5SRichard Tran Mills { 53*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 54*34b254c5SRichard Tran Mills 55*34b254c5SRichard Tran Mills PetscFunctionBegin; 56*34b254c5SRichard Tran Mills if (!linear->ksp) { 57*34b254c5SRichard Tran Mills PetscCall(KSPCreate(PetscObjectComm((PetscObject)regressor), &linear->ksp)); 58*34b254c5SRichard Tran Mills PetscCall(PetscObjectIncrementTabLevel((PetscObject)linear->ksp, (PetscObject)regressor, 1)); 59*34b254c5SRichard Tran Mills PetscCall(PetscObjectSetOptions((PetscObject)linear->ksp, ((PetscObject)regressor)->options)); 60*34b254c5SRichard Tran Mills } 61*34b254c5SRichard Tran Mills *ksp = linear->ksp; 62*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 63*34b254c5SRichard Tran Mills } 64*34b254c5SRichard Tran Mills 65*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorLinearSetUseKSP_Linear(PetscRegressor regressor, PetscBool flg) 66*34b254c5SRichard Tran Mills { 67*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 68*34b254c5SRichard Tran Mills 69*34b254c5SRichard Tran Mills PetscFunctionBegin; 70*34b254c5SRichard Tran Mills linear->use_ksp = flg; 71*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 72*34b254c5SRichard Tran Mills } 73*34b254c5SRichard Tran Mills 74*34b254c5SRichard Tran Mills static PetscErrorCode EvaluateResidual(Tao tao, Vec x, Vec f, void *ptr) 75*34b254c5SRichard Tran Mills { 76*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)ptr; 77*34b254c5SRichard Tran Mills 78*34b254c5SRichard Tran Mills PetscFunctionBegin; 79*34b254c5SRichard Tran Mills /* Evaluate f = A * x - b */ 80*34b254c5SRichard Tran Mills PetscCall(MatMult(linear->X, x, f)); 81*34b254c5SRichard Tran Mills PetscCall(VecAXPY(f, -1.0, linear->rhs)); 82*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 83*34b254c5SRichard Tran Mills } 84*34b254c5SRichard Tran Mills 85*34b254c5SRichard Tran Mills static PetscErrorCode EvaluateJacobian(Tao tao, Vec x, Mat J, Mat Jpre, void *ptr) 86*34b254c5SRichard Tran Mills { 87*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)ptr; 88*34b254c5SRichard Tran Mills 89*34b254c5SRichard Tran Mills PetscFunctionBegin; 90*34b254c5SRichard Tran Mills J = linear->X; 91*34b254c5SRichard Tran Mills Jpre = linear->X; 92*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 93*34b254c5SRichard Tran Mills } 94*34b254c5SRichard Tran Mills 95*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorSetUp_Linear(PetscRegressor regressor) 96*34b254c5SRichard Tran Mills { 97*34b254c5SRichard Tran Mills PetscInt M, N; 98*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 99*34b254c5SRichard Tran Mills KSP ksp; 100*34b254c5SRichard Tran Mills Tao tao; 101*34b254c5SRichard Tran Mills 102*34b254c5SRichard Tran Mills PetscFunctionBegin; 103*34b254c5SRichard Tran Mills PetscCall(MatGetSize(regressor->training, &M, &N)); 104*34b254c5SRichard Tran Mills 105*34b254c5SRichard Tran Mills if (linear->fit_intercept) { 106*34b254c5SRichard Tran Mills /* If we are fitting the intercept, we need to make A a composite matrix using MATCENTERING to preserve sparsity. 107*34b254c5SRichard Tran Mills * Though there might be some cases we don't want to do this for, depending on what kind of matrix is passed in. (Probably bad idea for dense?) 108*34b254c5SRichard Tran Mills * We will also need to ensure that the right-hand side passed to the KSP is also mean-centered, since we 109*34b254c5SRichard Tran Mills * intend to compute the intercept separately from regression coefficients (that is, we will not be adding a 110*34b254c5SRichard Tran Mills * column of all 1s to our design matrix). */ 111*34b254c5SRichard Tran Mills PetscCall(MatCreateCentering(PetscObjectComm((PetscObject)regressor), PETSC_DECIDE, M, &linear->C)); 112*34b254c5SRichard Tran Mills PetscCall(MatCreate(PetscObjectComm((PetscObject)regressor), &linear->X)); 113*34b254c5SRichard Tran Mills PetscCall(MatSetSizes(linear->X, PETSC_DECIDE, PETSC_DECIDE, M, N)); 114*34b254c5SRichard Tran Mills PetscCall(MatSetType(linear->X, MATCOMPOSITE)); 115*34b254c5SRichard Tran Mills PetscCall(MatCompositeSetType(linear->X, MAT_COMPOSITE_MULTIPLICATIVE)); 116*34b254c5SRichard Tran Mills PetscCall(MatCompositeAddMat(linear->X, regressor->training)); 117*34b254c5SRichard Tran Mills PetscCall(MatCompositeAddMat(linear->X, linear->C)); 118*34b254c5SRichard Tran Mills PetscCall(VecDuplicate(regressor->target, &linear->rhs)); 119*34b254c5SRichard Tran Mills PetscCall(MatMult(linear->C, regressor->target, linear->rhs)); 120*34b254c5SRichard Tran Mills } else { 121*34b254c5SRichard Tran Mills // When not fitting intercept, we assume that the input data are already centered. 122*34b254c5SRichard Tran Mills linear->X = regressor->training; 123*34b254c5SRichard Tran Mills linear->rhs = regressor->target; 124*34b254c5SRichard Tran Mills 125*34b254c5SRichard Tran Mills PetscCall(PetscObjectReference((PetscObject)linear->X)); 126*34b254c5SRichard Tran Mills PetscCall(PetscObjectReference((PetscObject)linear->rhs)); 127*34b254c5SRichard Tran Mills } 128*34b254c5SRichard Tran Mills 129*34b254c5SRichard Tran Mills if (linear->coefficients) PetscCall(VecDestroy(&linear->coefficients)); 130*34b254c5SRichard Tran Mills 131*34b254c5SRichard Tran Mills if (linear->use_ksp) { 132*34b254c5SRichard Tran Mills PetscCheck(linear->type == REGRESSOR_LINEAR_OLS, PetscObjectComm((PetscObject)regressor), PETSC_ERR_ARG_WRONGSTATE, "KSP can be used to fit a linear regressor only when its type is OLS"); 133*34b254c5SRichard Tran Mills 134*34b254c5SRichard Tran Mills if (!linear->ksp) PetscCall(PetscRegressorLinearGetKSP(regressor, &linear->ksp)); 135*34b254c5SRichard Tran Mills ksp = linear->ksp; 136*34b254c5SRichard Tran Mills 137*34b254c5SRichard Tran Mills PetscCall(MatCreateVecs(linear->X, &linear->coefficients, NULL)); 138*34b254c5SRichard Tran Mills /* Set up the KSP to solve the least squares problem (without solving for intercept, as this is done separately) using KSPLSQR. */ 139*34b254c5SRichard Tran Mills PetscCall(MatCreateNormal(linear->X, &linear->XtX)); 140*34b254c5SRichard Tran Mills PetscCall(KSPSetType(ksp, KSPLSQR)); 141*34b254c5SRichard Tran Mills PetscCall(KSPSetOperators(ksp, linear->X, linear->XtX)); 142*34b254c5SRichard Tran Mills PetscCall(KSPSetOptionsPrefix(ksp, ((PetscObject)regressor)->prefix)); 143*34b254c5SRichard Tran Mills PetscCall(KSPAppendOptionsPrefix(ksp, "regressor_linear_")); 144*34b254c5SRichard Tran Mills PetscCall(KSPSetFromOptions(ksp)); 145*34b254c5SRichard Tran Mills } else { 146*34b254c5SRichard Tran Mills /* Note: Currently implementation creates TAO inside of implementations. 147*34b254c5SRichard Tran Mills * Thus, all the prefix jobs are done inside implementations, not in interface */ 148*34b254c5SRichard Tran Mills const char *prefix; 149*34b254c5SRichard Tran Mills 150*34b254c5SRichard Tran Mills if (!regressor->tao) PetscCall(PetscRegressorGetTao(regressor, &tao)); 151*34b254c5SRichard Tran Mills 152*34b254c5SRichard Tran Mills PetscCall(MatCreateVecs(linear->X, &linear->coefficients, &linear->residual)); 153*34b254c5SRichard Tran Mills /* Set up the TAO object to solve the (regularized) least squares problem (without solving for intercept, which is done separately) using TAOBRGN. */ 154*34b254c5SRichard Tran Mills PetscCall(TaoSetType(tao, TAOBRGN)); 155*34b254c5SRichard Tran Mills PetscCall(TaoSetSolution(tao, linear->coefficients)); 156*34b254c5SRichard Tran Mills PetscCall(TaoSetResidualRoutine(tao, linear->residual, EvaluateResidual, linear)); 157*34b254c5SRichard Tran Mills PetscCall(TaoSetJacobianResidualRoutine(tao, linear->X, linear->X, EvaluateJacobian, linear)); 158*34b254c5SRichard Tran Mills if (!linear->use_ksp) PetscCall(TaoBRGNSetRegularizerWeight(tao, regressor->regularizer_weight)); 159*34b254c5SRichard Tran Mills // Set the regularization type and weight for the BRGN as linear->type dictates: 160*34b254c5SRichard Tran Mills // TODO BRGN needs to be BRGNSetRegularizationType 161*34b254c5SRichard Tran Mills // PetscOptionsSetValue no longer works due to functioning prefix system 162*34b254c5SRichard Tran Mills PetscCall(PetscRegressorGetOptionsPrefix(regressor, &prefix)); 163*34b254c5SRichard Tran Mills PetscCall(TaoSetOptionsPrefix(regressor->tao, prefix)); 164*34b254c5SRichard Tran Mills PetscCall(TaoAppendOptionsPrefix(tao, "regressor_linear_")); 165*34b254c5SRichard Tran Mills { 166*34b254c5SRichard Tran Mills TAO_BRGN *gn = (TAO_BRGN *)regressor->tao->data; 167*34b254c5SRichard Tran Mills 168*34b254c5SRichard Tran Mills switch (linear->type) { 169*34b254c5SRichard Tran Mills case REGRESSOR_LINEAR_OLS: 170*34b254c5SRichard Tran Mills regressor->regularizer_weight = 0.0; // OLS, by definition, uses a regularizer weight of 0 171*34b254c5SRichard Tran Mills break; 172*34b254c5SRichard Tran Mills case REGRESSOR_LINEAR_LASSO: 173*34b254c5SRichard Tran Mills gn->reg_type = BRGN_REGULARIZATION_L1DICT; 174*34b254c5SRichard Tran Mills break; 175*34b254c5SRichard Tran Mills case REGRESSOR_LINEAR_RIDGE: 176*34b254c5SRichard Tran Mills gn->reg_type = BRGN_REGULARIZATION_L2PURE; 177*34b254c5SRichard Tran Mills break; 178*34b254c5SRichard Tran Mills default: 179*34b254c5SRichard Tran Mills break; 180*34b254c5SRichard Tran Mills } 181*34b254c5SRichard Tran Mills } 182*34b254c5SRichard Tran Mills PetscCall(TaoSetFromOptions(tao)); 183*34b254c5SRichard Tran Mills } 184*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 185*34b254c5SRichard Tran Mills } 186*34b254c5SRichard Tran Mills 187*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorReset_Linear(PetscRegressor regressor) 188*34b254c5SRichard Tran Mills { 189*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 190*34b254c5SRichard Tran Mills 191*34b254c5SRichard Tran Mills PetscFunctionBegin; 192*34b254c5SRichard Tran Mills /* Destroy the PETSc objects associated with the linear regressor implementation. */ 193*34b254c5SRichard Tran Mills linear->ksp_its = 0; 194*34b254c5SRichard Tran Mills linear->ksp_tot_its = 0; 195*34b254c5SRichard Tran Mills 196*34b254c5SRichard Tran Mills PetscCall(MatDestroy(&linear->X)); 197*34b254c5SRichard Tran Mills PetscCall(MatDestroy(&linear->XtX)); 198*34b254c5SRichard Tran Mills PetscCall(MatDestroy(&linear->C)); 199*34b254c5SRichard Tran Mills PetscCall(KSPDestroy(&linear->ksp)); 200*34b254c5SRichard Tran Mills PetscCall(VecDestroy(&linear->coefficients)); 201*34b254c5SRichard Tran Mills PetscCall(VecDestroy(&linear->rhs)); 202*34b254c5SRichard Tran Mills PetscCall(VecDestroy(&linear->residual)); 203*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 204*34b254c5SRichard Tran Mills } 205*34b254c5SRichard Tran Mills 206*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorDestroy_Linear(PetscRegressor regressor) 207*34b254c5SRichard Tran Mills { 208*34b254c5SRichard Tran Mills PetscFunctionBegin; 209*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearSetFitIntercept_C", NULL)); 210*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearSetUseKSP_C", NULL)); 211*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearGetKSP_C", NULL)); 212*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearGetCoefficients_C", NULL)); 213*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearGetIntercept_C", NULL)); 214*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearSetType_C", NULL)); 215*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearGetType_C", NULL)); 216*34b254c5SRichard Tran Mills PetscCall(PetscRegressorReset_Linear(regressor)); 217*34b254c5SRichard Tran Mills PetscCall(PetscFree(regressor->data)); 218*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 219*34b254c5SRichard Tran Mills } 220*34b254c5SRichard Tran Mills 221*34b254c5SRichard Tran Mills /*@ 222*34b254c5SRichard Tran Mills PetscRegressorLinearSetFitIntercept - Set a flag to indicate that the intercept (also known as the "bias" or "offset") should 223*34b254c5SRichard Tran Mills be calculated; data are assumed to be mean-centered if false. 224*34b254c5SRichard Tran Mills 225*34b254c5SRichard Tran Mills Logically Collective 226*34b254c5SRichard Tran Mills 227*34b254c5SRichard Tran Mills Input Parameters: 228*34b254c5SRichard Tran Mills + regressor - the `PetscRegressor` context 229*34b254c5SRichard Tran Mills - flg - `PETSC_TRUE` to calculate the intercept, `PETSC_FALSE` to assume mean-centered data (default is `PETSC_TRUE`) 230*34b254c5SRichard Tran Mills 231*34b254c5SRichard Tran Mills Level: intermediate 232*34b254c5SRichard Tran Mills 233*34b254c5SRichard Tran Mills Options Database Key: 234*34b254c5SRichard Tran Mills . regressor_linear_fit_intercept <true,false> - fit the intercept 235*34b254c5SRichard Tran Mills 236*34b254c5SRichard Tran Mills Note: 237*34b254c5SRichard Tran Mills If the user indicates that the intercept should not be calculated, the intercept will be set to zero. 238*34b254c5SRichard Tran Mills 239*34b254c5SRichard Tran Mills .seealso: `PetscRegressor`, `PetscRegressorFit()` 240*34b254c5SRichard Tran Mills @*/ 241*34b254c5SRichard Tran Mills PetscErrorCode PetscRegressorLinearSetFitIntercept(PetscRegressor regressor, PetscBool flg) 242*34b254c5SRichard Tran Mills { 243*34b254c5SRichard Tran Mills PetscFunctionBegin; 244*34b254c5SRichard Tran Mills /* TODO: Add companion PetscRegressorLinearGetFitIntercept(), and put it in the .seealso: */ 245*34b254c5SRichard Tran Mills PetscValidHeaderSpecific(regressor, PETSCREGRESSOR_CLASSID, 1); 246*34b254c5SRichard Tran Mills PetscValidLogicalCollectiveBool(regressor, flg, 2); 247*34b254c5SRichard Tran Mills PetscTryMethod(regressor, "PetscRegressorLinearSetFitIntercept_C", (PetscRegressor, PetscBool), (regressor, flg)); 248*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 249*34b254c5SRichard Tran Mills } 250*34b254c5SRichard Tran Mills 251*34b254c5SRichard Tran Mills /*@ 252*34b254c5SRichard Tran Mills PetscRegressorLinearSetUseKSP - Set a flag to indicate that a `KSP` object, instead of a `Tao` one, should be used 253*34b254c5SRichard Tran Mills to fit the regressor 254*34b254c5SRichard Tran Mills 255*34b254c5SRichard Tran Mills Logically Collective 256*34b254c5SRichard Tran Mills 257*34b254c5SRichard Tran Mills Input Parameters: 258*34b254c5SRichard Tran Mills + regressor - the `PetscRegressor` context 259*34b254c5SRichard Tran Mills - flg - `PETSC_TRUE` to use a `KSP`, `PETSC_FALSE` to use a `Tao` object (default is false) 260*34b254c5SRichard Tran Mills 261*34b254c5SRichard Tran Mills Options Database Key: 262*34b254c5SRichard Tran Mills . regressor_linear_use_ksp <true,false> - use `KSP` 263*34b254c5SRichard Tran Mills 264*34b254c5SRichard Tran Mills Level: intermediate 265*34b254c5SRichard Tran Mills 266*34b254c5SRichard Tran Mills .seealso: `PetscRegressor`, `PetscRegressorLinearGetKSP()` 267*34b254c5SRichard Tran Mills @*/ 268*34b254c5SRichard Tran Mills PetscErrorCode PetscRegressorLinearSetUseKSP(PetscRegressor regressor, PetscBool flg) 269*34b254c5SRichard Tran Mills { 270*34b254c5SRichard Tran Mills PetscFunctionBegin; 271*34b254c5SRichard Tran Mills /* TODO: Add companion PetscRegressorLinearGetUseKSP(), and put it in the .seealso: */ 272*34b254c5SRichard Tran Mills PetscValidHeaderSpecific(regressor, PETSCREGRESSOR_CLASSID, 1); 273*34b254c5SRichard Tran Mills PetscValidLogicalCollectiveBool(regressor, flg, 2); 274*34b254c5SRichard Tran Mills PetscTryMethod(regressor, "PetscRegressorLinearSetUseKSP_C", (PetscRegressor, PetscBool), (regressor, flg)); 275*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 276*34b254c5SRichard Tran Mills } 277*34b254c5SRichard Tran Mills 278*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorSetFromOptions_Linear(PetscRegressor regressor, PetscOptionItems PetscOptionsObject) 279*34b254c5SRichard Tran Mills { 280*34b254c5SRichard Tran Mills PetscBool set, flg = PETSC_FALSE; 281*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 282*34b254c5SRichard Tran Mills 283*34b254c5SRichard Tran Mills PetscFunctionBegin; 284*34b254c5SRichard Tran Mills PetscOptionsHeadBegin(PetscOptionsObject, "PetscRegressor options for linear regressors"); 285*34b254c5SRichard Tran Mills PetscCall(PetscOptionsBool("-regressor_linear_fit_intercept", "Calculate intercept for linear model", "PetscRegressorLinearSetFitIntercept", flg, &flg, &set)); 286*34b254c5SRichard Tran Mills if (set) PetscCall(PetscRegressorLinearSetFitIntercept(regressor, flg)); 287*34b254c5SRichard Tran Mills PetscCall(PetscOptionsBool("-regressor_linear_use_ksp", "Use KSP instead of TAO for linear model fitting problem", "PetscRegressorLinearSetFitIntercept", flg, &flg, &set)); 288*34b254c5SRichard Tran Mills if (set) PetscCall(PetscRegressorLinearSetUseKSP(regressor, flg)); 289*34b254c5SRichard Tran Mills PetscCall(PetscOptionsEnum("-regressor_linear_type", "Linear regression method", "PetscRegressorLinearTypes", PetscRegressorLinearTypes, (PetscEnum)linear->type, (PetscEnum *)&linear->type, NULL)); 290*34b254c5SRichard Tran Mills PetscOptionsHeadEnd(); 291*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 292*34b254c5SRichard Tran Mills } 293*34b254c5SRichard Tran Mills 294*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorView_Linear(PetscRegressor regressor, PetscViewer viewer) 295*34b254c5SRichard Tran Mills { 296*34b254c5SRichard Tran Mills PetscBool isascii; 297*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 298*34b254c5SRichard Tran Mills 299*34b254c5SRichard Tran Mills PetscFunctionBegin; 300*34b254c5SRichard Tran Mills PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 301*34b254c5SRichard Tran Mills if (isascii) { 302*34b254c5SRichard Tran Mills PetscCall(PetscViewerASCIIPushTab(viewer)); 303*34b254c5SRichard Tran Mills PetscCall(PetscViewerASCIIPrintf(viewer, "PetscRegressor Linear Type: %s\n", PetscRegressorLinearTypes[linear->type])); 304*34b254c5SRichard Tran Mills if (linear->ksp) { 305*34b254c5SRichard Tran Mills PetscCall(KSPView(linear->ksp, viewer)); 306*34b254c5SRichard Tran Mills PetscCall(PetscViewerASCIIPrintf(viewer, "total KSP iterations: %" PetscInt_FMT "\n", linear->ksp_tot_its)); 307*34b254c5SRichard Tran Mills } 308*34b254c5SRichard Tran Mills if (linear->fit_intercept) PetscCall(PetscViewerASCIIPrintf(viewer, "Intercept=%g\n", (double)linear->intercept)); 309*34b254c5SRichard Tran Mills PetscCall(PetscViewerASCIIPopTab(viewer)); 310*34b254c5SRichard Tran Mills } 311*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 312*34b254c5SRichard Tran Mills } 313*34b254c5SRichard Tran Mills 314*34b254c5SRichard Tran Mills /*@ 315*34b254c5SRichard Tran Mills PetscRegressorLinearGetKSP - Returns the `KSP` context for a `PETSCREGRESSORLINEAR` object. 316*34b254c5SRichard Tran Mills 317*34b254c5SRichard Tran Mills Not Collective, but if the `PetscRegressor` is parallel, then the `KSP` object is parallel 318*34b254c5SRichard Tran Mills 319*34b254c5SRichard Tran Mills Input Parameter: 320*34b254c5SRichard Tran Mills . regressor - the `PetscRegressor` context 321*34b254c5SRichard Tran Mills 322*34b254c5SRichard Tran Mills Output Parameter: 323*34b254c5SRichard Tran Mills . ksp - the `KSP` context 324*34b254c5SRichard Tran Mills 325*34b254c5SRichard Tran Mills Level: beginner 326*34b254c5SRichard Tran Mills 327*34b254c5SRichard Tran Mills Note: 328*34b254c5SRichard Tran Mills This routine will always return a `KSP`, but, depending on the type of the linear regressor and the options that are set, the regressor may actually use a `Tao` object instead of this `KSP`. 329*34b254c5SRichard Tran Mills 330*34b254c5SRichard Tran Mills .seealso: `PetscRegressorGetTao()` 331*34b254c5SRichard Tran Mills @*/ 332*34b254c5SRichard Tran Mills PetscErrorCode PetscRegressorLinearGetKSP(PetscRegressor regressor, KSP *ksp) 333*34b254c5SRichard Tran Mills { 334*34b254c5SRichard Tran Mills PetscFunctionBegin; 335*34b254c5SRichard Tran Mills PetscValidHeaderSpecific(regressor, PETSCREGRESSOR_CLASSID, 1); 336*34b254c5SRichard Tran Mills PetscAssertPointer(ksp, 2); 337*34b254c5SRichard Tran Mills PetscUseMethod(regressor, "PetscRegressorLinearGetKSP_C", (PetscRegressor, KSP *), (regressor, ksp)); 338*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 339*34b254c5SRichard Tran Mills } 340*34b254c5SRichard Tran Mills 341*34b254c5SRichard Tran Mills /*@ 342*34b254c5SRichard Tran Mills PetscRegressorLinearGetCoefficients - Get a vector of the fitted coefficients from a linear regression model 343*34b254c5SRichard Tran Mills 344*34b254c5SRichard Tran Mills Not Collective but the vector is parallel 345*34b254c5SRichard Tran Mills 346*34b254c5SRichard Tran Mills Input Parameter: 347*34b254c5SRichard Tran Mills . regressor - the `PetscRegressor` context 348*34b254c5SRichard Tran Mills 349*34b254c5SRichard Tran Mills Output Parameter: 350*34b254c5SRichard Tran Mills . coefficients - the vector of the coefficients 351*34b254c5SRichard Tran Mills 352*34b254c5SRichard Tran Mills Level: beginner 353*34b254c5SRichard Tran Mills 354*34b254c5SRichard Tran Mills .seealso: `PetscRegressor`, `PetscRegressorLinearGetIntercept()`, `PETSCREGRESSORLINEAR`, `Vec` 355*34b254c5SRichard Tran Mills @*/ 356*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearGetCoefficients(PetscRegressor regressor, Vec *coefficients) 357*34b254c5SRichard Tran Mills { 358*34b254c5SRichard Tran Mills PetscFunctionBegin; 359*34b254c5SRichard Tran Mills PetscValidHeaderSpecific(regressor, PETSCREGRESSOR_CLASSID, 1); 360*34b254c5SRichard Tran Mills PetscAssertPointer(coefficients, 2); 361*34b254c5SRichard Tran Mills PetscUseMethod(regressor, "PetscRegressorLinearGetCoefficients_C", (PetscRegressor, Vec *), (regressor, coefficients)); 362*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 363*34b254c5SRichard Tran Mills } 364*34b254c5SRichard Tran Mills 365*34b254c5SRichard Tran Mills /*@ 366*34b254c5SRichard Tran Mills PetscRegressorLinearGetIntercept - Get the intercept from a linear regression model 367*34b254c5SRichard Tran Mills 368*34b254c5SRichard Tran Mills Not Collective 369*34b254c5SRichard Tran Mills 370*34b254c5SRichard Tran Mills Input Parameter: 371*34b254c5SRichard Tran Mills . regressor - the `PetscRegressor` context 372*34b254c5SRichard Tran Mills 373*34b254c5SRichard Tran Mills Output Parameter: 374*34b254c5SRichard Tran Mills . intercept - the intercept 375*34b254c5SRichard Tran Mills 376*34b254c5SRichard Tran Mills Level: beginner 377*34b254c5SRichard Tran Mills 378*34b254c5SRichard Tran Mills .seealso: `PetscRegressor`, `PetscRegressorLinearSetFitIntercept()`, `PetscRegressorLinearGetCoefficients()`, `PETSCREGRESSORLINEAR` 379*34b254c5SRichard Tran Mills @*/ 380*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorLinearGetIntercept(PetscRegressor regressor, PetscScalar *intercept) 381*34b254c5SRichard Tran Mills { 382*34b254c5SRichard Tran Mills PetscFunctionBegin; 383*34b254c5SRichard Tran Mills PetscValidHeaderSpecific(regressor, PETSCREGRESSOR_CLASSID, 1); 384*34b254c5SRichard Tran Mills PetscAssertPointer(intercept, 2); 385*34b254c5SRichard Tran Mills PetscUseMethod(regressor, "PetscRegressorLinearGetIntercept_C", (PetscRegressor, PetscScalar *), (regressor, intercept)); 386*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 387*34b254c5SRichard Tran Mills } 388*34b254c5SRichard Tran Mills 389*34b254c5SRichard Tran Mills /*@C 390*34b254c5SRichard Tran Mills PetscRegressorLinearSetType - Sets the type of linear regression to be performed 391*34b254c5SRichard Tran Mills 392*34b254c5SRichard Tran Mills Logically Collective 393*34b254c5SRichard Tran Mills 394*34b254c5SRichard Tran Mills Input Parameters: 395*34b254c5SRichard Tran Mills + regressor - the `PetscRegressor` context (should be of type `PETSCREGRESSORLINEAR`) 396*34b254c5SRichard Tran Mills - type - a known linear regression method 397*34b254c5SRichard Tran Mills 398*34b254c5SRichard Tran Mills Options Database Key: 399*34b254c5SRichard Tran Mills . -regressor_linear_type - Sets the linear regression method; use -help for a list of available methods 400*34b254c5SRichard Tran Mills (for instance "-regressor_linear_type ols" or "-regressor_linear_type lasso") 401*34b254c5SRichard Tran Mills 402*34b254c5SRichard Tran Mills Level: intermediate 403*34b254c5SRichard Tran Mills 404*34b254c5SRichard Tran Mills .seealso: `PetscRegressorLinearGetType()`, `PetscRegressorLinearType`, `PetscRegressorSetType()` 405*34b254c5SRichard Tran Mills @*/ 406*34b254c5SRichard Tran Mills PetscErrorCode PetscRegressorLinearSetType(PetscRegressor regressor, PetscRegressorLinearType type) 407*34b254c5SRichard Tran Mills { 408*34b254c5SRichard Tran Mills PetscFunctionBegin; 409*34b254c5SRichard Tran Mills PetscValidHeaderSpecific(regressor, PETSCREGRESSOR_CLASSID, 1); 410*34b254c5SRichard Tran Mills PetscValidLogicalCollectiveEnum(regressor, type, 2); 411*34b254c5SRichard Tran Mills PetscTryMethod(regressor, "PetscRegressorLinearSetType_C", (PetscRegressor, PetscRegressorLinearType), (regressor, type)); 412*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 413*34b254c5SRichard Tran Mills } 414*34b254c5SRichard Tran Mills 415*34b254c5SRichard Tran Mills /*@ 416*34b254c5SRichard Tran Mills PetscRegressorLinearGetType - Return the type for the `PETSCREGRESSORLINEAR` solver 417*34b254c5SRichard Tran Mills 418*34b254c5SRichard Tran Mills Input Parameter: 419*34b254c5SRichard Tran Mills . regressor - the `PetscRegressor` solver context 420*34b254c5SRichard Tran Mills 421*34b254c5SRichard Tran Mills Output Parameter: 422*34b254c5SRichard Tran Mills . type - `PETSCREGRESSORLINEAR` type 423*34b254c5SRichard Tran Mills 424*34b254c5SRichard Tran Mills Level: advanced 425*34b254c5SRichard Tran Mills 426*34b254c5SRichard Tran Mills .seealso: `PetscRegressor`, `PETSCREGRESSORLINEAR`, `PetscRegressorLinearSetType()`, `PetscRegressorLinearType` 427*34b254c5SRichard Tran Mills @*/ 428*34b254c5SRichard Tran Mills PetscErrorCode PetscRegressorLinearGetType(PetscRegressor regressor, PetscRegressorLinearType *type) 429*34b254c5SRichard Tran Mills { 430*34b254c5SRichard Tran Mills PetscFunctionBegin; 431*34b254c5SRichard Tran Mills PetscValidHeaderSpecific(regressor, PETSCREGRESSOR_CLASSID, 1); 432*34b254c5SRichard Tran Mills PetscAssertPointer(type, 2); 433*34b254c5SRichard Tran Mills PetscUseMethod(regressor, "PetscRegressorLinearGetType_C", (PetscRegressor, PetscRegressorLinearType *), (regressor, type)); 434*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 435*34b254c5SRichard Tran Mills } 436*34b254c5SRichard Tran Mills 437*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorFit_Linear(PetscRegressor regressor) 438*34b254c5SRichard Tran Mills { 439*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 440*34b254c5SRichard Tran Mills KSP ksp; 441*34b254c5SRichard Tran Mills PetscScalar target_mean, *column_means_global, *column_means_local, column_means_dot_coefficients; 442*34b254c5SRichard Tran Mills Vec column_means; 443*34b254c5SRichard Tran Mills PetscInt m, N, istart, i, kspits; 444*34b254c5SRichard Tran Mills 445*34b254c5SRichard Tran Mills PetscFunctionBegin; 446*34b254c5SRichard Tran Mills if (linear->use_ksp) PetscCall(PetscRegressorLinearGetKSP(regressor, &linear->ksp)); 447*34b254c5SRichard Tran Mills ksp = linear->ksp; 448*34b254c5SRichard Tran Mills 449*34b254c5SRichard Tran Mills /* Solve the least-squares problem (previously set up in PetscRegressorSetUp_Linear()) without finding the intercept. */ 450*34b254c5SRichard Tran Mills if (linear->use_ksp) { 451*34b254c5SRichard Tran Mills PetscCall(KSPSolve(ksp, linear->rhs, linear->coefficients)); 452*34b254c5SRichard Tran Mills PetscCall(KSPGetIterationNumber(ksp, &kspits)); 453*34b254c5SRichard Tran Mills linear->ksp_its += kspits; 454*34b254c5SRichard Tran Mills linear->ksp_tot_its += kspits; 455*34b254c5SRichard Tran Mills } else { 456*34b254c5SRichard Tran Mills PetscCall(TaoSolve(regressor->tao)); 457*34b254c5SRichard Tran Mills } 458*34b254c5SRichard Tran Mills 459*34b254c5SRichard Tran Mills /* Calculate the intercept. */ 460*34b254c5SRichard Tran Mills if (linear->fit_intercept) { 461*34b254c5SRichard Tran Mills PetscCall(MatGetSize(regressor->training, NULL, &N)); 462*34b254c5SRichard Tran Mills PetscCall(PetscMalloc1(N, &column_means_global)); 463*34b254c5SRichard Tran Mills PetscCall(VecMean(regressor->target, &target_mean)); 464*34b254c5SRichard Tran Mills /* We need the means of all columns of regressor->training, placed into a Vec compatible with linear->coefficients. 465*34b254c5SRichard Tran Mills * Note the potential scalability issue: MatGetColumnMeans() computes means of ALL colummns. */ 466*34b254c5SRichard Tran Mills PetscCall(MatGetColumnMeans(regressor->training, column_means_global)); 467*34b254c5SRichard Tran Mills /* TODO: Calculation of the Vec and matrix column means should probably go into the SetUp phase, and also be placed 468*34b254c5SRichard Tran Mills * into a routine that is callable from outside of PetscRegressorFit_Linear(), because we'll want to do the same 469*34b254c5SRichard Tran Mills * thing for other models, such as ridge and LASSO regression, and should avoid code duplication. 470*34b254c5SRichard Tran Mills * What we are calling 'target_mean' and 'column_means' should be stashed in the base linear regressor struct, 471*34b254c5SRichard Tran Mills * and perhaps renamed to make it clear they are offsets that should be applied (though the current naming 472*34b254c5SRichard Tran Mills * makes sense since it makes it clear where these come from.) */ 473*34b254c5SRichard Tran Mills PetscCall(VecDuplicate(linear->coefficients, &column_means)); 474*34b254c5SRichard Tran Mills PetscCall(VecGetLocalSize(column_means, &m)); 475*34b254c5SRichard Tran Mills PetscCall(VecGetOwnershipRange(column_means, &istart, NULL)); 476*34b254c5SRichard Tran Mills PetscCall(VecGetArrayWrite(column_means, &column_means_local)); 477*34b254c5SRichard Tran Mills for (i = 0; i < m; i++) column_means_local[i] = column_means_global[istart + i]; 478*34b254c5SRichard Tran Mills PetscCall(VecRestoreArrayWrite(column_means, &column_means_local)); 479*34b254c5SRichard Tran Mills PetscCall(VecDot(column_means, linear->coefficients, &column_means_dot_coefficients)); 480*34b254c5SRichard Tran Mills PetscCall(VecDestroy(&column_means)); 481*34b254c5SRichard Tran Mills PetscCall(PetscFree(column_means_global)); 482*34b254c5SRichard Tran Mills linear->intercept = target_mean - column_means_dot_coefficients; 483*34b254c5SRichard Tran Mills } else { 484*34b254c5SRichard Tran Mills linear->intercept = 0.0; 485*34b254c5SRichard Tran Mills } 486*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 487*34b254c5SRichard Tran Mills } 488*34b254c5SRichard Tran Mills 489*34b254c5SRichard Tran Mills static PetscErrorCode PetscRegressorPredict_Linear(PetscRegressor regressor, Mat X, Vec y) 490*34b254c5SRichard Tran Mills { 491*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear = (PetscRegressor_Linear *)regressor->data; 492*34b254c5SRichard Tran Mills 493*34b254c5SRichard Tran Mills PetscFunctionBegin; 494*34b254c5SRichard Tran Mills PetscCall(MatMult(X, linear->coefficients, y)); 495*34b254c5SRichard Tran Mills PetscCall(VecShift(y, linear->intercept)); 496*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 497*34b254c5SRichard Tran Mills } 498*34b254c5SRichard Tran Mills 499*34b254c5SRichard Tran Mills /*MC 500*34b254c5SRichard Tran Mills PETSCREGRESSORLINEAR - Linear regression model (ordinary least squares or regularized variants) 501*34b254c5SRichard Tran Mills 502*34b254c5SRichard Tran Mills Options Database: 503*34b254c5SRichard Tran Mills + -regressor_linear_fit_intercept - Calculate the intercept for the linear model 504*34b254c5SRichard Tran Mills - -regressor_linear_use_ksp - Use `KSP` instead of `Tao` for linear model fitting (non-regularized variants only) 505*34b254c5SRichard Tran Mills 506*34b254c5SRichard Tran Mills Level: beginner 507*34b254c5SRichard Tran Mills 508*34b254c5SRichard Tran Mills Note: 509*34b254c5SRichard Tran Mills This is the default regressor in `PetscRegressor`. 510*34b254c5SRichard Tran Mills 511*34b254c5SRichard Tran Mills .seealso: `PetscRegressorCreate()`, `PetscRegressor`, `PetscRegressorSetType()` 512*34b254c5SRichard Tran Mills M*/ 513*34b254c5SRichard Tran Mills PETSC_EXTERN PetscErrorCode PetscRegressorCreate_Linear(PetscRegressor regressor) 514*34b254c5SRichard Tran Mills { 515*34b254c5SRichard Tran Mills PetscRegressor_Linear *linear; 516*34b254c5SRichard Tran Mills 517*34b254c5SRichard Tran Mills PetscFunctionBegin; 518*34b254c5SRichard Tran Mills PetscCall(PetscNew(&linear)); 519*34b254c5SRichard Tran Mills regressor->data = (void *)linear; 520*34b254c5SRichard Tran Mills 521*34b254c5SRichard Tran Mills regressor->ops->setup = PetscRegressorSetUp_Linear; 522*34b254c5SRichard Tran Mills regressor->ops->reset = PetscRegressorReset_Linear; 523*34b254c5SRichard Tran Mills regressor->ops->destroy = PetscRegressorDestroy_Linear; 524*34b254c5SRichard Tran Mills regressor->ops->setfromoptions = PetscRegressorSetFromOptions_Linear; 525*34b254c5SRichard Tran Mills regressor->ops->view = PetscRegressorView_Linear; 526*34b254c5SRichard Tran Mills regressor->ops->fit = PetscRegressorFit_Linear; 527*34b254c5SRichard Tran Mills regressor->ops->predict = PetscRegressorPredict_Linear; 528*34b254c5SRichard Tran Mills 529*34b254c5SRichard Tran Mills linear->intercept = 0.0; 530*34b254c5SRichard Tran Mills linear->fit_intercept = PETSC_TRUE; /* Default to calculating the intercept. */ 531*34b254c5SRichard Tran Mills linear->use_ksp = PETSC_FALSE; /* Do not default to using KSP for solving the model-fitting problem (use TAO instead). */ 532*34b254c5SRichard Tran Mills linear->type = REGRESSOR_LINEAR_OLS; 533*34b254c5SRichard Tran Mills /* Above, manually set the default linear regressor type. 534*34b254c5SRichard Tran Mills We don't use PetscRegressorLinearSetType() here, because that expects the SetUp event to already have happened. */ 535*34b254c5SRichard Tran Mills 536*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearSetFitIntercept_C", PetscRegressorLinearSetFitIntercept_Linear)); 537*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearSetUseKSP_C", PetscRegressorLinearSetUseKSP_Linear)); 538*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearGetKSP_C", PetscRegressorLinearGetKSP_Linear)); 539*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearGetCoefficients_C", PetscRegressorLinearGetCoefficients_Linear)); 540*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearGetIntercept_C", PetscRegressorLinearGetIntercept_Linear)); 541*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearSetType_C", PetscRegressorLinearSetType_Linear)); 542*34b254c5SRichard Tran Mills PetscCall(PetscObjectComposeFunction((PetscObject)regressor, "PetscRegressorLinearGetType_C", PetscRegressorLinearGetType_Linear)); 543*34b254c5SRichard Tran Mills PetscFunctionReturn(PETSC_SUCCESS); 544*34b254c5SRichard Tran Mills } 545