xref: /petsc/src/mat/tests/ex3.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscScalar    v[3];
12c4762a1bSJed Brown   PetscReal      omega = 1.0,norm;
13c4762a1bSJed Brown 
14*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-omega",&omega,NULL));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
17c4762a1bSJed Brown 
185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&C));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,n,n,n,n));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(C,MATSEQDENSE));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&b));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&x));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&u));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&e));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u,1.0));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,0.0));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   v[0] = -1.; v[1] = 2.; v[2] = -1.;
30c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
31c4762a1bSJed Brown     midx[0] = i-1; midx[1] = i; midx[2] = i+1;
325f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(C,1,&i,3,midx,v,INSERT_VALUES));
33c4762a1bSJed Brown   }
34c4762a1bSJed Brown   i    = 0; midx[0] = 0; midx[1] = 1;
35c4762a1bSJed Brown   v[0] = 2.0; v[1] = -1.;
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES));
37c4762a1bSJed Brown   i    = n-1; midx[0] = n-2; midx[1] = n-1;
38c4762a1bSJed Brown   v[0] = -1.0; v[1] = 2.;
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES));
40c4762a1bSJed Brown 
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
43c4762a1bSJed Brown 
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(C,u,b));
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   for (i=0; i<n; i++) {
475f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSOR(C,b,omega,SOR_FORWARD_SWEEP,0.0,1,1,x));
485f80ce2aSJacob Faibussowitsch     CHKERRQ(VecWAXPY(e,-1.0,x,u));
495f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(e,NORM_2,&norm));
505f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"2-norm of error %g\n",(double)norm));
51c4762a1bSJed Brown   }
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&e));
57*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
58*b122ec5aSJacob Faibussowitsch   return 0;
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown /*TEST
62c4762a1bSJed Brown 
63c4762a1bSJed Brown    test:
64c4762a1bSJed Brown 
65c4762a1bSJed Brown TEST*/
66