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