xref: /petsc/src/ksp/pc/tests/ex4.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Demonstrates the use of fast Richardson for SOR. And tests\n\
3c4762a1bSJed Brown the MatSOR() routines.\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscpc.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            mat;
10c4762a1bSJed Brown   Vec            b,u;
11c4762a1bSJed Brown   PC             pc;
12c4762a1bSJed Brown   PetscInt       n = 5,i,col[3];
13c4762a1bSJed Brown   PetscScalar    value[3];
14c4762a1bSJed Brown 
15*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
16c4762a1bSJed Brown   /* Create vectors */
175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&b));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&u));
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   /* Create and assemble matrix */
215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&mat));
22c4762a1bSJed Brown   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
23c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
24c4762a1bSJed Brown     col[0] = i-1; col[1] = i; col[2] = i+1;
255f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES));
26c4762a1bSJed Brown   }
27c4762a1bSJed Brown   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
29c4762a1bSJed Brown   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /* Create PC context and set up data structures */
355f80ce2aSJacob Faibussowitsch   CHKERRQ(PCCreate(PETSC_COMM_WORLD,&pc));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(pc,PCSOR));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetFromOptions(pc));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetOperators(pc,mat,mat));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetUp(pc));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   value[0] = 1.0;
42c4762a1bSJed Brown   for (i=0; i<n; i++) {
435f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(u,0.0));
445f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(u,1,&i,value,INSERT_VALUES));
455f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAssemblyBegin(u));
465f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAssemblyEnd(u));
475f80ce2aSJacob Faibussowitsch     CHKERRQ(PCApply(pc,u,b));
485f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(b,PETSC_VIEWER_STDOUT_SELF));
49c4762a1bSJed Brown   }
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Free data structures */
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&mat));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(PCDestroy(&pc));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
56*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
57*b122ec5aSJacob Faibussowitsch   return 0;
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown /*TEST
61c4762a1bSJed Brown 
62c4762a1bSJed Brown    test:
63c4762a1bSJed Brown 
64c4762a1bSJed Brown TEST*/
65