1b9ac7092SAlp Dener /* 2b9ac7092SAlp Dener This provides a thin wrapper around LMVM matrices in order to use their MatSolve 3b9ac7092SAlp Dener methods as preconditioner applications in KSP solves. 4b9ac7092SAlp Dener */ 5b9ac7092SAlp Dener 6b9ac7092SAlp Dener #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 7b9ac7092SAlp Dener 8b9ac7092SAlp Dener typedef struct { 9b9ac7092SAlp Dener Vec xwork, ywork; 10b9ac7092SAlp Dener IS inactive; 11b9ac7092SAlp Dener Mat B; 12*b2d8c577SAlp Dener PetscBool allocated; 13b9ac7092SAlp Dener } PC_LMVM; 14b9ac7092SAlp Dener 15b9ac7092SAlp Dener /*@ 16b9ac7092SAlp Dener PCLMVMSetMatLMVM - Replaces the LMVM matrix inside the preconditioner with 17b9ac7092SAlp Dener the one provided by the user. 18b9ac7092SAlp Dener 19b9ac7092SAlp Dener Input Parameters: 20b9ac7092SAlp Dener + pc - An LMVM preconditioner 21b9ac7092SAlp Dener - B - An LMVM-type matrix (MATLDFP, MATLBFGS, MATLSR1, MATLBRDN, MATLMBRDN, MATLSBRDN) 22b9ac7092SAlp Dener 23b9ac7092SAlp Dener Level: intermediate 24b9ac7092SAlp Dener @*/ 25b9ac7092SAlp Dener PetscErrorCode PCLMVMSetMatLMVM(PC pc, Mat B) 26b9ac7092SAlp Dener { 27b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM*)pc->data; 28b9ac7092SAlp Dener PetscErrorCode ierr; 29b9ac7092SAlp Dener PetscBool same; 30b9ac7092SAlp Dener 31b9ac7092SAlp Dener PetscFunctionBegin; 32b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 33b9ac7092SAlp Dener PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 34b9ac7092SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same);CHKERRQ(ierr); 35b9ac7092SAlp Dener if (!same) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 36b9ac7092SAlp Dener ierr = PetscObjectBaseTypeCompare((PetscObject)B, MATLMVM, &same);CHKERRQ(ierr); 37b9ac7092SAlp Dener if (!same) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Matrix must be an LMVM-type."); 38b9ac7092SAlp Dener ierr = MatDestroy(&ctx->B);CHKERRQ(ierr); 39b9ac7092SAlp Dener ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 40b9ac7092SAlp Dener ctx->B = B; 41b9ac7092SAlp Dener PetscFunctionReturn(0); 42b9ac7092SAlp Dener } 43b9ac7092SAlp Dener 44b9ac7092SAlp Dener /*@ 45b9ac7092SAlp Dener PCLMVMGetMatLMVM - Returns a pointer to the underlying LMVM matrix. 46b9ac7092SAlp Dener 47b9ac7092SAlp Dener Input Parameters: 48b9ac7092SAlp Dener . pc - An LMVM preconditioner 49b9ac7092SAlp Dener 50b9ac7092SAlp Dener Output Parameters: 51b9ac7092SAlp Dener . B - LMVM matrix inside the preconditioner 52b9ac7092SAlp Dener 53b9ac7092SAlp Dener Level: intermediate 54b9ac7092SAlp Dener @*/ 55b9ac7092SAlp Dener PetscErrorCode PCLMVMGetMatLMVM(PC pc, Mat *B) 56b9ac7092SAlp Dener { 57b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM*)pc->data; 58b9ac7092SAlp Dener PetscErrorCode ierr; 59b9ac7092SAlp Dener PetscBool same; 60b9ac7092SAlp Dener 61b9ac7092SAlp Dener PetscFunctionBegin; 62b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 63b9ac7092SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same);CHKERRQ(ierr); 64b9ac7092SAlp Dener if (!same) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 65b9ac7092SAlp Dener *B = ctx->B; 66b9ac7092SAlp Dener PetscFunctionReturn(0); 67b9ac7092SAlp Dener } 68b9ac7092SAlp Dener 69b9ac7092SAlp Dener /*@ 70b9ac7092SAlp Dener PCLMVMSetIS - Sets the index sets that reduce the PC application. 71b9ac7092SAlp Dener 72b9ac7092SAlp Dener Input Parameters: 73b9ac7092SAlp Dener + pc - An LMVM preconditioner 74b9ac7092SAlp Dener - inactive - Index set defining the variables removed from the problem 75b9ac7092SAlp Dener 76b9ac7092SAlp Dener Level: intermediate 77b9ac7092SAlp Dener 78b9ac7092SAlp Dener .seealso: MatLMVMUpdate() 79b9ac7092SAlp Dener @*/ 80b9ac7092SAlp Dener PetscErrorCode PCLMVMSetIS(PC pc, IS inactive) 81b9ac7092SAlp Dener { 82b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM*)pc->data; 83b9ac7092SAlp Dener PetscErrorCode ierr; 84b9ac7092SAlp Dener PetscBool same; 85b9ac7092SAlp Dener 86b9ac7092SAlp Dener PetscFunctionBegin; 87b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 88*b2d8c577SAlp Dener PetscValidHeaderSpecific(inactive, IS_CLASSID, 2); 89b9ac7092SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same);CHKERRQ(ierr); 90b9ac7092SAlp Dener if (!same) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 91b9ac7092SAlp Dener ierr = PCLMVMClearIS(pc); 92b9ac7092SAlp Dener ierr = PetscObjectReference((PetscObject)inactive);CHKERRQ(ierr); 93b9ac7092SAlp Dener ctx->inactive = inactive; 94b9ac7092SAlp Dener PetscFunctionReturn(0); 95b9ac7092SAlp Dener } 96b9ac7092SAlp Dener 97b9ac7092SAlp Dener /*@ 98b9ac7092SAlp Dener PCLMVMClearIS - Removes the inactive variable index set. 99b9ac7092SAlp Dener 100b9ac7092SAlp Dener Input Parameters: 101b9ac7092SAlp Dener . pc - An LMVM preconditioner 102b9ac7092SAlp Dener 103b9ac7092SAlp Dener Level: intermediate 104b9ac7092SAlp Dener 105b9ac7092SAlp Dener .seealso: MatLMVMUpdate() 106b9ac7092SAlp Dener @*/ 107b9ac7092SAlp Dener PetscErrorCode PCLMVMClearIS(PC pc) 108b9ac7092SAlp Dener { 109b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM*)pc->data; 110b9ac7092SAlp Dener PetscErrorCode ierr; 111b9ac7092SAlp Dener PetscBool same; 112b9ac7092SAlp Dener 113b9ac7092SAlp Dener PetscFunctionBegin; 114b9ac7092SAlp Dener PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 115b9ac7092SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &same);CHKERRQ(ierr); 116b9ac7092SAlp Dener if (!same) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC must be a PCLMVM type."); 117b9ac7092SAlp Dener if (ctx->inactive) { 118b9ac7092SAlp Dener ierr = ISDestroy(&ctx->inactive);CHKERRQ(ierr); 119b9ac7092SAlp Dener } 120b9ac7092SAlp Dener PetscFunctionReturn(0); 121b9ac7092SAlp Dener } 122b9ac7092SAlp Dener 123b9ac7092SAlp Dener static PetscErrorCode PCApply_LMVM(PC pc,Vec x,Vec y) 124b9ac7092SAlp Dener { 125b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM*)pc->data; 126b9ac7092SAlp Dener PetscErrorCode ierr; 127b9ac7092SAlp Dener Vec xsub, ysub; 128b9ac7092SAlp Dener 129b9ac7092SAlp Dener PetscFunctionBegin; 130b9ac7092SAlp Dener if (ctx->inactive) { 131b9ac7092SAlp Dener ierr = VecZeroEntries(ctx->xwork);CHKERRQ(ierr); 132b9ac7092SAlp Dener ierr = VecGetSubVector(ctx->xwork, ctx->inactive, &xsub);CHKERRQ(ierr); 133b9ac7092SAlp Dener ierr = VecCopy(x, xsub);CHKERRQ(ierr); 134b9ac7092SAlp Dener ierr = VecRestoreSubVector(ctx->xwork, ctx->inactive, &xsub);CHKERRQ(ierr); 135b9ac7092SAlp Dener } else { 136b9ac7092SAlp Dener ierr = VecCopy(x, ctx->xwork);CHKERRQ(ierr); 137b9ac7092SAlp Dener } 138b9ac7092SAlp Dener ierr = MatSolve(ctx->B, ctx->xwork, ctx->ywork);CHKERRQ(ierr); 139b9ac7092SAlp Dener if (ctx->inactive) { 140b9ac7092SAlp Dener ierr = VecGetSubVector(ctx->ywork, ctx->inactive, &ysub);CHKERRQ(ierr); 141b9ac7092SAlp Dener ierr = VecCopy(ysub, y);CHKERRQ(ierr); 142b9ac7092SAlp Dener ierr = VecRestoreSubVector(ctx->ywork, ctx->inactive, &ysub);CHKERRQ(ierr); 143b9ac7092SAlp Dener } else { 144b9ac7092SAlp Dener ierr = VecCopy(ctx->ywork, y);CHKERRQ(ierr); 145b9ac7092SAlp Dener } 146b9ac7092SAlp Dener PetscFunctionReturn(0); 147b9ac7092SAlp Dener } 148b9ac7092SAlp Dener 149b9ac7092SAlp Dener static PetscErrorCode PCReset_LMVM(PC pc) 150b9ac7092SAlp Dener { 151b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM*)pc->data; 152b9ac7092SAlp Dener PetscErrorCode ierr; 153b9ac7092SAlp Dener 154b9ac7092SAlp Dener PetscFunctionBegin; 155b9ac7092SAlp Dener if (ctx->xwork) { 156b9ac7092SAlp Dener ierr = VecDestroy(&ctx->xwork);CHKERRQ(ierr); 157b9ac7092SAlp Dener } 158b9ac7092SAlp Dener if (ctx->ywork) { 159b9ac7092SAlp Dener ierr = VecDestroy(&ctx->ywork);CHKERRQ(ierr); 160b9ac7092SAlp Dener } 161b9ac7092SAlp Dener PetscFunctionReturn(0); 162b9ac7092SAlp Dener } 163b9ac7092SAlp Dener 164b9ac7092SAlp Dener static PetscErrorCode PCSetUp_LMVM(PC pc) 165b9ac7092SAlp Dener { 166b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM*)pc->data; 167b9ac7092SAlp Dener PetscErrorCode ierr; 168b9ac7092SAlp Dener PetscInt n, N; 169*b2d8c577SAlp Dener PetscBool allocated; 170b9ac7092SAlp Dener 171b9ac7092SAlp Dener PetscFunctionBegin; 172*b2d8c577SAlp Dener ierr = MatLMVMIsAllocated(ctx->B, &allocated);CHKERRQ(ierr); 173*b2d8c577SAlp Dener if (!allocated) { 174b9ac7092SAlp Dener ierr = MatCreateVecs(pc->mat, &ctx->xwork, &ctx->ywork);CHKERRQ(ierr); 175b9ac7092SAlp Dener ierr = VecGetLocalSize(ctx->xwork, &n);CHKERRQ(ierr); 176b9ac7092SAlp Dener ierr = VecGetSize(ctx->xwork, &N);CHKERRQ(ierr); 177b9ac7092SAlp Dener ierr = MatSetSizes(ctx->B, n, n, N, N);CHKERRQ(ierr); 178b9ac7092SAlp Dener ierr = MatLMVMAllocate(ctx->B, ctx->xwork, ctx->ywork);CHKERRQ(ierr); 179*b2d8c577SAlp Dener } else { 180*b2d8c577SAlp Dener ierr = MatCreateVecs(ctx->B, &ctx->xwork, &ctx->ywork);CHKERRQ(ierr); 181*b2d8c577SAlp Dener } 182b9ac7092SAlp Dener PetscFunctionReturn(0); 183b9ac7092SAlp Dener } 184b9ac7092SAlp Dener 185b9ac7092SAlp Dener static PetscErrorCode PCView_LMVM(PC pc,PetscViewer viewer) 186b9ac7092SAlp Dener { 187b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM*)pc->data; 188b9ac7092SAlp Dener PetscErrorCode ierr; 189b9ac7092SAlp Dener PetscBool iascii; 190b9ac7092SAlp Dener 191b9ac7092SAlp Dener PetscFunctionBegin; 192b9ac7092SAlp Dener ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 193b9ac7092SAlp Dener if (iascii) { 194b9ac7092SAlp Dener ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 195b9ac7092SAlp Dener ierr = MatView(ctx->B, viewer);CHKERRQ(ierr); 196b9ac7092SAlp Dener ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 197b9ac7092SAlp Dener } 198b9ac7092SAlp Dener PetscFunctionReturn(0); 199b9ac7092SAlp Dener } 200b9ac7092SAlp Dener 201b9ac7092SAlp Dener static PetscErrorCode PCSetFromOptions_LMVM(PetscOptionItems* PetscOptionsObject, PC pc) 202b9ac7092SAlp Dener { 203b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM*)pc->data; 204b9ac7092SAlp Dener PetscErrorCode ierr; 205b9ac7092SAlp Dener 206b9ac7092SAlp Dener PetscFunctionBegin; 207b9ac7092SAlp Dener ierr = MatSetFromOptions(ctx->B);CHKERRQ(ierr); 208b9ac7092SAlp Dener PetscFunctionReturn(0); 209b9ac7092SAlp Dener } 210b9ac7092SAlp Dener 211b9ac7092SAlp Dener static PetscErrorCode PCDestroy_LMVM(PC pc) 212b9ac7092SAlp Dener { 213b9ac7092SAlp Dener PC_LMVM *ctx = (PC_LMVM*)pc->data; 214b9ac7092SAlp Dener PetscErrorCode ierr; 215b9ac7092SAlp Dener 216b9ac7092SAlp Dener PetscFunctionBegin; 217b9ac7092SAlp Dener if (ctx->inactive) { 218b9ac7092SAlp Dener ierr = ISDestroy(&ctx->inactive);CHKERRQ(ierr); 219b9ac7092SAlp Dener } 220b9ac7092SAlp Dener if (pc->setupcalled) { 221b9ac7092SAlp Dener ierr = VecDestroy(&ctx->xwork);CHKERRQ(ierr); 222b9ac7092SAlp Dener ierr = VecDestroy(&ctx->ywork);CHKERRQ(ierr); 223b9ac7092SAlp Dener } 224b9ac7092SAlp Dener ierr = MatDestroy(&ctx->B);CHKERRQ(ierr); 225b9ac7092SAlp Dener ierr = PetscFree(pc->data);CHKERRQ(ierr); 226b9ac7092SAlp Dener PetscFunctionReturn(0); 227b9ac7092SAlp Dener } 228b9ac7092SAlp Dener 229b9ac7092SAlp Dener /*MC 230b9ac7092SAlp Dener PCLMVM - Creates a preconditioner around an LMVM matrix. Options for the 231b9ac7092SAlp Dener underlying LMVM matrix can be access with the "-pc_lmvm_" prefix. 232b9ac7092SAlp Dener 233b9ac7092SAlp Dener .seealso: PCCreate(), PCSetType(), PCType (for list of available types), 234b9ac7092SAlp Dener PC, MATLMVM, PCLMVMUpdate(), PCLMVMSetMatLMVM(), PCLMVMGetMatLMVM() 235b9ac7092SAlp Dener M*/ 236b9ac7092SAlp Dener PETSC_EXTERN PetscErrorCode PCCreate_LMVM(PC pc) 237b9ac7092SAlp Dener { 238b9ac7092SAlp Dener PetscErrorCode ierr; 239b9ac7092SAlp Dener PC_LMVM *ctx; 240b9ac7092SAlp Dener 241b9ac7092SAlp Dener PetscFunctionBegin; 242b9ac7092SAlp Dener ierr = PetscNewLog(pc,&ctx);CHKERRQ(ierr); 243b9ac7092SAlp Dener pc->data = (void*)ctx; 244b9ac7092SAlp Dener 245b9ac7092SAlp Dener pc->ops->reset = PCReset_LMVM; 246b9ac7092SAlp Dener pc->ops->setup = PCSetUp_LMVM; 247b9ac7092SAlp Dener pc->ops->destroy = PCDestroy_LMVM; 248b9ac7092SAlp Dener pc->ops->view = PCView_LMVM; 249b9ac7092SAlp Dener pc->ops->apply = PCApply_LMVM; 250b9ac7092SAlp Dener pc->ops->setfromoptions = PCSetFromOptions_LMVM; 251b9ac7092SAlp Dener pc->ops->applysymmetricleft = 0; 252b9ac7092SAlp Dener pc->ops->applysymmetricright = 0; 253b9ac7092SAlp Dener pc->ops->applytranspose = 0; 254b9ac7092SAlp Dener pc->ops->applyrichardson = 0; 255b9ac7092SAlp Dener pc->ops->presolve = 0; 256b9ac7092SAlp Dener pc->ops->postsolve = 0; 257b9ac7092SAlp Dener 258b9ac7092SAlp Dener ierr = PCSetReusePreconditioner(pc, PETSC_TRUE); 259b9ac7092SAlp Dener 260b9ac7092SAlp Dener ierr = MatCreate(PetscObjectComm((PetscObject)pc), &ctx->B);CHKERRQ(ierr); 261b9ac7092SAlp Dener ierr = MatSetType(ctx->B, MATLBFGS);CHKERRQ(ierr); 262b9ac7092SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)ctx->B, (PetscObject)pc, 1);CHKERRQ(ierr); 263b9ac7092SAlp Dener ierr = MatSetOptionsPrefix(ctx->B, "pc_lmvm_");CHKERRQ(ierr); 264b9ac7092SAlp Dener PetscFunctionReturn(0); 265b9ac7092SAlp Dener } 266b9ac7092SAlp Dener 267b9ac7092SAlp Dener 268b9ac7092SAlp Dener 269b9ac7092SAlp Dener 270b9ac7092SAlp Dener 271b9ac7092SAlp Dener 272