Actual source code: ex55.c

  1: static const char help[]="Example demonstrating PCCOMPOSITE where one of the inner PCs uses a different operator\n\
  2: \n";

  4: /*T
  5:    Concepts: KSP^using nested solves
  6:    Concepts: PC^using composite PCs
  7:    Processors: n
  8: T*/
  9: #include <petscksp.h>

 11: int main(int argc, char **argv)
 12: {
 13:   PetscInt       n=10,i,col[3];
 14:   Vec            x,b;
 15:   Mat            A,B;
 16:   KSP            ksp;
 17:   PC             pc,subpc;
 18:   PetscScalar    value[3];

 20:   PetscInitialize(&argc,&argv,(char*)0,help);

 22:   /* Create a diagonal matrix with a given distribution of diagonal elements */
 23:   MatCreate(PETSC_COMM_WORLD,&A);
 24:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 25:   MatSetFromOptions(A);
 26:   MatSetUp(A);
 27:    /*
 28:      Assemble matrix
 29:   */
 30:   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 31:   for (i=1; i<n-1; i++) {
 32:     col[0] = i-1; col[1] = i; col[2] = i+1;
 33:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 34:   }
 35:   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
 36:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 37:   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 38:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 39:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 40:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 42:   MatCreateVecs(A,&x,&b);

 44:   /* Create a KSP object */
 45:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 46:   KSPSetOperators(ksp,A,A);

 48:   /* Set up a composite preconditioner */
 49:   KSPGetPC(ksp,&pc);
 50:   PCSetType(pc,PCCOMPOSITE); /* default composite with single Identity PC */
 51:   PCCompositeSetType(pc,PC_COMPOSITE_ADDITIVE);
 52:   PCCompositeAddPCType(pc,PCLU);
 53:   PCCompositeGetPC(pc,0,&subpc);
 54:   /*  B is set to the diagonal of A; this demonstrates that setting the operator for a subpc changes the preconditioning */
 55:   MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&B);
 56:   MatGetDiagonal(A,b);
 57:   MatDiagonalSet(B,b,ADD_VALUES);
 58:   PCSetOperators(subpc,B,B);
 59:   PCCompositeAddPCType(pc,PCNONE);

 61:   KSPSetFromOptions(ksp);
 62:   KSPSolve(ksp,b,x);

 64:   KSPDestroy(&ksp);
 65:   MatDestroy(&A);
 66:   MatDestroy(&B);
 67:   VecDestroy(&x);
 68:   VecDestroy(&b);
 69:   PetscFinalize();
 70:   return 0;
 71: }

 73: /*TEST

 75:    test:
 76:      args: -ksp_monitor -pc_composite_type multiplicative

 78: TEST*/