xref: /petsc/src/ksp/pc/tests/ex4.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry Smith   PetscFunctionBeginUser;
169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
17c4762a1bSJed Brown   /* Create vectors */
189566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&b));
199566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&u));
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /* Create and assemble matrix */
229566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&mat));
23c4762a1bSJed Brown   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
24c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
25c4762a1bSJed Brown     col[0] = i-1; col[1] = i; col[2] = i+1;
269566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES));
27c4762a1bSJed Brown   }
28c4762a1bSJed Brown   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
299566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
30c4762a1bSJed Brown   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
319566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   /* Create PC context and set up data structures */
369566063dSJacob Faibussowitsch   PetscCall(PCCreate(PETSC_COMM_WORLD,&pc));
379566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCSOR));
389566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions(pc));
399566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc,mat,mat));
409566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   value[0] = 1.0;
43c4762a1bSJed Brown   for (i=0; i<n; i++) {
449566063dSJacob Faibussowitsch     PetscCall(VecSet(u,0.0));
459566063dSJacob Faibussowitsch     PetscCall(VecSetValues(u,1,&i,value,INSERT_VALUES));
469566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(u));
479566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(u));
489566063dSJacob Faibussowitsch     PetscCall(PCApply(pc,u,b));
499566063dSJacob Faibussowitsch     PetscCall(VecView(b,PETSC_VIEWER_STDOUT_SELF));
50c4762a1bSJed Brown   }
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Free data structures */
539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
549566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&pc));
559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
579566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
58b122ec5aSJacob Faibussowitsch   return 0;
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown /*TEST
62c4762a1bSJed Brown 
63c4762a1bSJed Brown    test:
64c4762a1bSJed Brown 
65c4762a1bSJed Brown TEST*/
66