Actual source code: ex3.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2: static char help[] = "Tests relaxation for dense matrices.\n\n";

  4:  #include <petscmat.h>

  6: int main(int argc,char **args)
  7: {
  8:   Mat            C;
  9:   Vec            u,x,b,e;
 10:   PetscInt       i,n = 10,midx[3];
 12:   PetscScalar    v[3];
 13:   PetscReal      omega = 1.0,norm;

 15:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 16:   PetscOptionsGetReal(NULL,NULL,"-omega",&omega,NULL);
 17:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);

 19:   MatCreate(PETSC_COMM_SELF,&C);
 20:   MatSetSizes(C,n,n,n,n);
 21:   MatSetType(C,MATSEQDENSE);
 22:   MatSetUp(C);
 23:   VecCreateSeq(PETSC_COMM_SELF,n,&b);
 24:   VecCreateSeq(PETSC_COMM_SELF,n,&x);
 25:   VecCreateSeq(PETSC_COMM_SELF,n,&u);
 26:   VecCreateSeq(PETSC_COMM_SELF,n,&e);
 27:   VecSet(u,1.0);
 28:   VecSet(x,0.0);

 30:   v[0] = -1.; v[1] = 2.; v[2] = -1.;
 31:   for (i=1; i<n-1; i++) {
 32:     midx[0] = i-1; midx[1] = i; midx[2] = i+1;
 33:     MatSetValues(C,1,&i,3,midx,v,INSERT_VALUES);
 34:   }
 35:   i    = 0; midx[0] = 0; midx[1] = 1;
 36:   v[0] = 2.0; v[1] = -1.;
 37:   MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES);
 38:   i    = n-1; midx[0] = n-2; midx[1] = n-1;
 39:   v[0] = -1.0; v[1] = 2.;
 40:   MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES);

 42:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 43:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 45:   MatMult(C,u,b);

 47:   for (i=0; i<n; i++) {
 48:     MatSOR(C,b,omega,SOR_FORWARD_SWEEP,0.0,1,1,x);
 49:     VecWAXPY(e,-1.0,x,u);
 50:     VecNorm(e,NORM_2,&norm);
 51:     PetscPrintf(PETSC_COMM_SELF,"2-norm of error %g\n",(double)norm);
 52:   }
 53:   MatDestroy(&C);
 54:   VecDestroy(&x);
 55:   VecDestroy(&b);
 56:   VecDestroy(&u);
 57:   VecDestroy(&e);
 58:   PetscFinalize();
 59:   return ierr;
 60: }

 62: /*TEST

 64:    test:

 66: TEST*/