1b9ac7092SAlp Dener /* 2f7bf01afSAlp Dener This provides a thin wrapper around LMVM matrices in order to use their MatLMVMSolve 3b9ac7092SAlp Dener methods as preconditioner applications in KSP solves. 4b9ac7092SAlp Dener */ 5b9ac7092SAlp Dener 6b9ac7092SAlp Dener #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 7dbf6bb8dSprj- #include <petsc/private/matimpl.h> 8b9ac7092SAlp Dener 9b9ac7092SAlp Dener typedef struct { 10b9ac7092SAlp Dener Vec xwork, ywork; 11b9ac7092SAlp Dener IS inactive; 12b9ac7092SAlp Dener Mat B; 13970ba202SStefano Zampini Vec X; 14970ba202SStefano Zampini PetscObjectState Xstate; 15970ba202SStefano Zampini PetscBool setfromoptionscalled; 16b9ac7092SAlp Dener } PC_LMVM; 17b9ac7092SAlp Dener 18b9ac7092SAlp Dener /*@ 19970ba202SStefano Zampini PCLMVMSetUpdateVec - Set the vector to be used as solution update for the internal LMVM matrix. 20b9ac7092SAlp Dener 21b9ac7092SAlp Dener Input Parameters: 22970ba202SStefano Zampini + pc - The preconditioner 23970ba202SStefano Zampini - X - Solution vector 24b9ac7092SAlp Dener 25b9ac7092SAlp Dener Level: intermediate 26f1580f4eSBarry Smith 27970ba202SStefano Zampini Notes: 28970ba202SStefano Zampini This is only needed if you want the preconditioner to automatically update the internal matrix. 29970ba202SStefano Zampini It is called in some `SNES` implementations to update the preconditioner. 30970ba202SStefano Zampini The right-hand side of the linear system is used as function vector. 31970ba202SStefano Zampini 32970ba202SStefano Zampini .seealso: `MatLMVMUpdate()`, `PCLMVMSetMatLMVM()` 33970ba202SStefano Zampini @*/ 34970ba202SStefano Zampini PetscErrorCode PCLMVMSetUpdateVec(PC pc, Vec X) 35970ba202SStefano Zampini { 36970ba202SStefano Zampini PC_LMVM *ctx; 37970ba202SStefano Zampini PetscBool same; 38970ba202SStefano Zampini 39970ba202SStefano Zampini PetscFunctionBegin; 40970ba202SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 41970ba202SStefano Zampini if (X) PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 42970ba202SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 43970ba202SStefano Zampini if (!same) PetscFunctionReturn(PETSC_SUCCESS); 44970ba202SStefano Zampini ctx = (PC_LMVM *)pc->data; 45970ba202SStefano Zampini PetscCall(PetscObjectReference((PetscObject)X)); 46970ba202SStefano Zampini PetscCall(VecDestroy(&ctx->X)); 47970ba202SStefano Zampini ctx->X = X; 48970ba202SStefano Zampini ctx->Xstate = -1; 49970ba202SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 50970ba202SStefano Zampini } 51970ba202SStefano Zampini 52970ba202SStefano Zampini /*@ 53970ba202SStefano Zampini PCLMVMSetMatLMVM - Replaces the `MATLMVM` matrix inside the preconditioner with the one provided by the user. 54970ba202SStefano Zampini 55970ba202SStefano Zampini Input Parameters: 56970ba202SStefano Zampini + pc - An `PCLMVM` preconditioner 57970ba202SStefano Zampini - B - An `MATLMVM` type matrix 58970ba202SStefano Zampini 59970ba202SStefano Zampini Level: intermediate 60970ba202SStefano Zampini 61970ba202SStefano Zampini .seealso: [](ch_ksp), `PCLMVMGetMatLMVM()` 62b9ac7092SAlp Dener @*/ 63d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLMVMSetMatLMVM(PC pc, Mat B) 64d71ae5a4SJacob Faibussowitsch { 65970ba202SStefano Zampini PC_LMVM *ctx; 66b9ac7092SAlp Dener PetscBool same; 67b9ac7092SAlp Dener 68b9ac7092SAlp Dener PetscFunctionBegin; 69b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 70b9ac7092SAlp Dener PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 72970ba202SStefano Zampini if (!same) PetscFunctionReturn(PETSC_SUCCESS); 739566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATLMVM, &same)); 74970ba202SStefano Zampini PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Matrix must be an MATLMVM."); 75970ba202SStefano Zampini ctx = (PC_LMVM *)pc->data; 769566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)B)); 77970ba202SStefano Zampini PetscCall(MatDestroy(&ctx->B)); 78b9ac7092SAlp Dener ctx->B = B; 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80b9ac7092SAlp Dener } 81b9ac7092SAlp Dener 82b9ac7092SAlp Dener /*@ 83f1580f4eSBarry Smith PCLMVMGetMatLMVM - Returns a pointer to the underlying `MATLMVM` matrix. 84b9ac7092SAlp Dener 85f1580f4eSBarry Smith Input Parameter: 86f1580f4eSBarry Smith . pc - An `PCLMVM` preconditioner 87b9ac7092SAlp Dener 88f1580f4eSBarry Smith Output Parameter: 89970ba202SStefano Zampini . B - `MATLMVM` matrix used by the preconditioner 90b9ac7092SAlp Dener 91b9ac7092SAlp Dener Level: intermediate 92f1580f4eSBarry Smith 93970ba202SStefano Zampini .seealso: [](ch_ksp), `PCLMVMSetMatLMVM()` 94b9ac7092SAlp Dener @*/ 95d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLMVMGetMatLMVM(PC pc, Mat *B) 96d71ae5a4SJacob Faibussowitsch { 97970ba202SStefano Zampini PC_LMVM *ctx; 98b9ac7092SAlp Dener PetscBool same; 99b9ac7092SAlp Dener 100b9ac7092SAlp Dener PetscFunctionBegin; 101b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1029566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 10328b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 104970ba202SStefano Zampini ctx = (PC_LMVM *)pc->data; 105970ba202SStefano Zampini if (!ctx->B) { 106970ba202SStefano Zampini Mat J; 107970ba202SStefano Zampini 108970ba202SStefano Zampini if (pc->useAmat) J = pc->mat; 109970ba202SStefano Zampini else J = pc->pmat; 110970ba202SStefano Zampini PetscCall(PetscObjectBaseTypeCompare((PetscObject)J, MATLMVM, &same)); 111970ba202SStefano Zampini if (same) *B = J; 112970ba202SStefano Zampini else { 113970ba202SStefano Zampini const char *prefix; 114970ba202SStefano Zampini 115970ba202SStefano Zampini PetscCall(PCGetOptionsPrefix(pc, &prefix)); 116970ba202SStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &ctx->B)); 117970ba202SStefano Zampini PetscCall(MatSetOptionsPrefix(ctx->B, prefix)); 118970ba202SStefano Zampini PetscCall(MatAppendOptionsPrefix(ctx->B, "pc_lmvm_")); 119970ba202SStefano Zampini PetscCall(MatSetType(ctx->B, MATLMVMBFGS)); 120970ba202SStefano Zampini PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->B, (PetscObject)pc, 1)); 121b9ac7092SAlp Dener *B = ctx->B; 122970ba202SStefano Zampini } 123970ba202SStefano Zampini } else *B = ctx->B; 1243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125b9ac7092SAlp Dener } 126b9ac7092SAlp Dener 127b9ac7092SAlp Dener /*@ 128f1580f4eSBarry Smith PCLMVMSetIS - Sets the index sets that reduce the `PC` application. 129b9ac7092SAlp Dener 130b9ac7092SAlp Dener Input Parameters: 131f1580f4eSBarry Smith + pc - An `PCLMVM` preconditioner 132b9ac7092SAlp Dener - inactive - Index set defining the variables removed from the problem 133b9ac7092SAlp Dener 134b9ac7092SAlp Dener Level: intermediate 135b9ac7092SAlp Dener 136feefa0e1SJacob Faibussowitsch Developer Notes: 137f1580f4eSBarry Smith Need to explain the purpose of this `IS` 138f1580f4eSBarry Smith 139970ba202SStefano Zampini .seealso: [](ch_ksp), `PCLMVMClearIS()` 140b9ac7092SAlp Dener @*/ 141d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLMVMSetIS(PC pc, IS inactive) 142d71ae5a4SJacob Faibussowitsch { 143970ba202SStefano Zampini PC_LMVM *ctx; 144b9ac7092SAlp Dener PetscBool same; 145b9ac7092SAlp Dener 146b9ac7092SAlp Dener PetscFunctionBegin; 147b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 148b2d8c577SAlp Dener PetscValidHeaderSpecific(inactive, IS_CLASSID, 2); 1499566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 150970ba202SStefano Zampini if (!same) PetscFunctionReturn(PETSC_SUCCESS); 151970ba202SStefano Zampini ctx = (PC_LMVM *)pc->data; 1529566063dSJacob Faibussowitsch PetscCall(PCLMVMClearIS(pc)); 1539566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)inactive)); 154b9ac7092SAlp Dener ctx->inactive = inactive; 1553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 156b9ac7092SAlp Dener } 157b9ac7092SAlp Dener 158b9ac7092SAlp Dener /*@ 159f1580f4eSBarry Smith PCLMVMClearIS - Removes the inactive variable index set from a `PCLMVM` 160b9ac7092SAlp Dener 161f1580f4eSBarry Smith Input Parameter: 162f1580f4eSBarry Smith . pc - An `PCLMVM` preconditioner 163b9ac7092SAlp Dener 164b9ac7092SAlp Dener Level: intermediate 165b9ac7092SAlp Dener 166970ba202SStefano Zampini .seealso: [](ch_ksp), `PCLMVMSetIS()` 167b9ac7092SAlp Dener @*/ 168d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLMVMClearIS(PC pc) 169d71ae5a4SJacob Faibussowitsch { 170970ba202SStefano Zampini PC_LMVM *ctx; 171b9ac7092SAlp Dener PetscBool same; 172b9ac7092SAlp Dener 173b9ac7092SAlp Dener PetscFunctionBegin; 174b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1759566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same)); 176970ba202SStefano Zampini if (!same) PetscFunctionReturn(PETSC_SUCCESS); 177970ba202SStefano Zampini ctx = (PC_LMVM *)pc->data; 178970ba202SStefano Zampini PetscCall(ISDestroy(&ctx->inactive)); 1793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 180b9ac7092SAlp Dener } 181b9ac7092SAlp Dener 182d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_LMVM(PC pc, Vec x, Vec y) 183d71ae5a4SJacob Faibussowitsch { 184b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 185970ba202SStefano Zampini Vec xsub, ysub, Bx = x, By = y; 186970ba202SStefano Zampini Mat B = ctx->B ? ctx->B : (pc->useAmat ? pc->mat : pc->pmat); 187b9ac7092SAlp Dener 188b9ac7092SAlp Dener PetscFunctionBegin; 189b9ac7092SAlp Dener if (ctx->inactive) { 190970ba202SStefano Zampini if (!ctx->xwork) PetscCall(MatCreateVecs(B, &ctx->xwork, &ctx->ywork)); 1919566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(ctx->xwork)); 1929566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ctx->xwork, ctx->inactive, &xsub)); 1939566063dSJacob Faibussowitsch PetscCall(VecCopy(x, xsub)); 1949566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ctx->xwork, ctx->inactive, &xsub)); 195970ba202SStefano Zampini Bx = ctx->xwork; 196970ba202SStefano Zampini By = ctx->ywork; 197b9ac7092SAlp Dener } 198970ba202SStefano Zampini PetscCall(MatSolve(B, Bx, By)); 199b9ac7092SAlp Dener if (ctx->inactive) { 2009566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ctx->ywork, ctx->inactive, &ysub)); 2019566063dSJacob Faibussowitsch PetscCall(VecCopy(ysub, y)); 2029566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ctx->ywork, ctx->inactive, &ysub)); 203b9ac7092SAlp Dener } 2043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 205b9ac7092SAlp Dener } 206b9ac7092SAlp Dener 207d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_LMVM(PC pc) 208d71ae5a4SJacob Faibussowitsch { 209b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 210b9ac7092SAlp Dener 211b9ac7092SAlp Dener PetscFunctionBegin; 212970ba202SStefano Zampini PetscCall(VecDestroy(&ctx->xwork)); 213970ba202SStefano Zampini PetscCall(VecDestroy(&ctx->ywork)); 2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 215b9ac7092SAlp Dener } 216b9ac7092SAlp Dener 217d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_LMVM(PC pc) 218d71ae5a4SJacob Faibussowitsch { 219b9ac7092SAlp Dener PetscInt n, N; 220b2d8c577SAlp Dener PetscBool allocated; 221970ba202SStefano Zampini Mat B; 222970ba202SStefano Zampini PC_LMVM *ctx = (PC_LMVM *)pc->data; 223b9ac7092SAlp Dener 224b9ac7092SAlp Dener PetscFunctionBegin; 225970ba202SStefano Zampini PetscCall(PCLMVMGetMatLMVM(pc, &B)); 226970ba202SStefano Zampini PetscCall(MatLMVMIsAllocated(B, &allocated)); 227b2d8c577SAlp Dener if (!allocated) { 228970ba202SStefano Zampini Vec t1, t2; 229970ba202SStefano Zampini 230970ba202SStefano Zampini PetscCall(MatCreateVecs(pc->mat, &t1, &t2)); 231970ba202SStefano Zampini PetscCall(VecGetLocalSize(t1, &n)); 232970ba202SStefano Zampini PetscCall(VecGetSize(t1, &N)); 233970ba202SStefano Zampini PetscCall(MatSetSizes(B, n, n, N, N)); 234970ba202SStefano Zampini PetscCall(MatLMVMAllocate(B, t1, t2)); 235970ba202SStefano Zampini PetscCall(VecDestroy(&t1)); 236970ba202SStefano Zampini PetscCall(VecDestroy(&t2)); 237b2d8c577SAlp Dener } 238970ba202SStefano Zampini /* Only call SetFromOptions if we internally handle the LMVM matrix */ 239970ba202SStefano Zampini if (B == ctx->B && ctx->setfromoptionscalled) PetscCall(MatSetFromOptions(ctx->B)); 240970ba202SStefano Zampini ctx->setfromoptionscalled = PETSC_FALSE; 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 242b9ac7092SAlp Dener } 243b9ac7092SAlp Dener 244d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_LMVM(PC pc, PetscViewer viewer) 245d71ae5a4SJacob Faibussowitsch { 246b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 247*9f196a02SMartin Diehl PetscBool isascii; 248b9ac7092SAlp Dener 249b9ac7092SAlp Dener PetscFunctionBegin; 250*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 251*9f196a02SMartin Diehl if (isascii && ctx->B && ctx->B->assembled) { 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 2539566063dSJacob Faibussowitsch PetscCall(MatView(ctx->B, viewer)); 2549566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 255b9ac7092SAlp Dener } 2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 257b9ac7092SAlp Dener } 258b9ac7092SAlp Dener 259ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_LMVM(PC pc, PetscOptionItems PetscOptionsObject) 260d71ae5a4SJacob Faibussowitsch { 261b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 262b9ac7092SAlp Dener 263b9ac7092SAlp Dener PetscFunctionBegin; 264970ba202SStefano Zampini /* defer SetFromOptions calls to PCSetUp_LMVM */ 265970ba202SStefano Zampini ctx->setfromoptionscalled = PETSC_TRUE; 266970ba202SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 267970ba202SStefano Zampini } 268970ba202SStefano Zampini 269970ba202SStefano Zampini static PetscErrorCode PCPreSolve_LMVM(PC pc, KSP ksp, Vec F, Vec X) 270970ba202SStefano Zampini { 271970ba202SStefano Zampini PC_LMVM *ctx = (PC_LMVM *)pc->data; 272970ba202SStefano Zampini 273970ba202SStefano Zampini PetscFunctionBegin; 274970ba202SStefano Zampini if (ctx->X && ctx->B) { /* Perform update only if requested. Otherwise we assume the user, e.g. TAO, has already taken care of it */ 275970ba202SStefano Zampini PetscObjectState Xstate; 276970ba202SStefano Zampini 277970ba202SStefano Zampini PetscCall(PetscObjectStateGet((PetscObject)ctx->X, &Xstate)); 278970ba202SStefano Zampini if (ctx->Xstate != Xstate) PetscCall(MatLMVMUpdate(ctx->B, ctx->X, F)); 279970ba202SStefano Zampini ctx->Xstate = Xstate; 280970ba202SStefano Zampini } 2813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 282b9ac7092SAlp Dener } 283b9ac7092SAlp Dener 284d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_LMVM(PC pc) 285d71ae5a4SJacob Faibussowitsch { 286b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM *)pc->data; 287b9ac7092SAlp Dener 288b9ac7092SAlp Dener PetscFunctionBegin; 289970ba202SStefano Zampini PetscCall(ISDestroy(&ctx->inactive)); 2909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->xwork)); 2919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->ywork)); 292970ba202SStefano Zampini PetscCall(VecDestroy(&ctx->X)); 2939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->B)); 2949566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 296b9ac7092SAlp Dener } 297b9ac7092SAlp Dener 298b9ac7092SAlp Dener /*MC 299970ba202SStefano Zampini PCLMVM - A preconditioner constructed from a `MATLMVM` matrix. 3007addb90fSBarry Smith If the matrix used to construct the preconditioner is not of type `MATLMVM`, an internal matrix is used. 3017addb90fSBarry Smith Options for the internal `MATLMVM` matrix can be accessed with the `-pc_lmvm_` prefix. 302970ba202SStefano Zampini Alternatively, the user can pass a suitable matrix with `PCLMVMSetMatLMVM()`. 303b9ac7092SAlp Dener 30406dd6b0eSSatish Balay Level: intermediate 30506dd6b0eSSatish Balay 306970ba202SStefano Zampini .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCLMVMSetUpdateVec()`, `PCLMVMSetMatLMVM()`, `PCLMVMGetMatLMVM()` 307b9ac7092SAlp Dener M*/ 308d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_LMVM(PC pc) 309d71ae5a4SJacob Faibussowitsch { 310b9ac7092SAlp Dener PC_LMVM *ctx; 311b9ac7092SAlp Dener 312b9ac7092SAlp Dener PetscFunctionBegin; 3134dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 314b9ac7092SAlp Dener pc->data = (void *)ctx; 315b9ac7092SAlp Dener 316b9ac7092SAlp Dener pc->ops->reset = PCReset_LMVM; 317b9ac7092SAlp Dener pc->ops->setup = PCSetUp_LMVM; 318b9ac7092SAlp Dener pc->ops->destroy = PCDestroy_LMVM; 319b9ac7092SAlp Dener pc->ops->view = PCView_LMVM; 320b9ac7092SAlp Dener pc->ops->apply = PCApply_LMVM; 321b9ac7092SAlp Dener pc->ops->setfromoptions = PCSetFromOptions_LMVM; 3220a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 3230a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 3240a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 3250a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 326970ba202SStefano Zampini pc->ops->presolve = PCPreSolve_LMVM; 3270a545947SLisandro Dalcin pc->ops->postsolve = NULL; 3283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 329b9ac7092SAlp Dener } 330