xref: /petsc/src/ksp/pc/tests/ex4.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Demonstrates the use of fast Richardson for SOR. And tests\n\
3*c4762a1bSJed Brown the MatSOR() routines.\n\n";
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown #include <petscpc.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown int main(int argc,char **args)
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   Mat            mat;
10*c4762a1bSJed Brown   Vec            b,u;
11*c4762a1bSJed Brown   PC             pc;
12*c4762a1bSJed Brown   PetscErrorCode ierr;
13*c4762a1bSJed Brown   PetscInt       n = 5,i,col[3];
14*c4762a1bSJed Brown   PetscScalar    value[3];
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
17*c4762a1bSJed Brown   /* Create vectors */
18*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&b);CHKERRQ(ierr);
19*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&u);CHKERRQ(ierr);
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown   /* Create and assemble matrix */
22*c4762a1bSJed Brown   ierr     = MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&mat);CHKERRQ(ierr);
23*c4762a1bSJed Brown   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
24*c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
25*c4762a1bSJed Brown     col[0] = i-1; col[1] = i; col[2] = i+1;
26*c4762a1bSJed Brown     ierr   = MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
27*c4762a1bSJed Brown   }
28*c4762a1bSJed Brown   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
29*c4762a1bSJed Brown   ierr = MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
30*c4762a1bSJed Brown   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
31*c4762a1bSJed Brown   ierr = MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown   /* Create PC context and set up data structures */
36*c4762a1bSJed Brown   ierr = PCCreate(PETSC_COMM_WORLD,&pc);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = PCSetType(pc,PCSOR);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = PCSetOperators(pc,mat,mat);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = PCSetUp(pc);CHKERRQ(ierr);
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown   value[0] = 1.0;
43*c4762a1bSJed Brown   for (i=0; i<n; i++) {
44*c4762a1bSJed Brown     ierr = VecSet(u,0.0);CHKERRQ(ierr);
45*c4762a1bSJed Brown     ierr = VecSetValues(u,1,&i,value,INSERT_VALUES);CHKERRQ(ierr);
46*c4762a1bSJed Brown     ierr = VecAssemblyBegin(u);CHKERRQ(ierr);
47*c4762a1bSJed Brown     ierr = VecAssemblyEnd(u);CHKERRQ(ierr);
48*c4762a1bSJed Brown     ierr = PCApply(pc,u,b);CHKERRQ(ierr);
49*c4762a1bSJed Brown     ierr = VecView(b,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
50*c4762a1bSJed Brown   }
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   /* Free data structures */
53*c4762a1bSJed Brown   ierr = MatDestroy(&mat);CHKERRQ(ierr);
54*c4762a1bSJed Brown   ierr = PCDestroy(&pc);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = PetscFinalize();
58*c4762a1bSJed Brown   return ierr;
59*c4762a1bSJed Brown }
60*c4762a1bSJed Brown 
61*c4762a1bSJed Brown /*TEST
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown    test:
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown TEST*/
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown 
68