xref: /petsc/src/mat/tests/ex3.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests relaxation for dense matrices.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            C;
9c4762a1bSJed Brown   Vec            u,x,b,e;
10c4762a1bSJed Brown   PetscInt       i,n = 10,midx[3];
11c4762a1bSJed Brown   PetscErrorCode ierr;
12c4762a1bSJed Brown   PetscScalar    v[3];
13c4762a1bSJed Brown   PetscReal      omega = 1.0,norm;
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
16*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-omega",&omega,NULL));
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
18c4762a1bSJed Brown 
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&C));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,n,n,n,n));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(C,MATSEQDENSE));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&b));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&x));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&u));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&e));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u,1.0));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,0.0));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   v[0] = -1.; v[1] = 2.; v[2] = -1.;
31c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
32c4762a1bSJed Brown     midx[0] = i-1; midx[1] = i; midx[2] = i+1;
33*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(C,1,&i,3,midx,v,INSERT_VALUES));
34c4762a1bSJed Brown   }
35c4762a1bSJed Brown   i    = 0; midx[0] = 0; midx[1] = 1;
36c4762a1bSJed Brown   v[0] = 2.0; v[1] = -1.;
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES));
38c4762a1bSJed Brown   i    = n-1; midx[0] = n-2; midx[1] = n-1;
39c4762a1bSJed Brown   v[0] = -1.0; v[1] = 2.;
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES));
41c4762a1bSJed Brown 
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
44c4762a1bSJed Brown 
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(C,u,b));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   for (i=0; i<n; i++) {
48*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSOR(C,b,omega,SOR_FORWARD_SWEEP,0.0,1,1,x));
49*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecWAXPY(e,-1.0,x,u));
50*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(e,NORM_2,&norm));
51*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"2-norm of error %g\n",(double)norm));
52c4762a1bSJed Brown   }
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&e));
58c4762a1bSJed Brown   ierr = PetscFinalize();
59c4762a1bSJed Brown   return ierr;
60c4762a1bSJed Brown }
61c4762a1bSJed Brown 
62c4762a1bSJed Brown /*TEST
63c4762a1bSJed Brown 
64c4762a1bSJed Brown    test:
65c4762a1bSJed Brown 
66c4762a1bSJed Brown TEST*/
67