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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 15*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-omega",&omega,NULL)); 16*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 17c4762a1bSJed Brown 18*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&C)); 19*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,n,n,n,n)); 20*9566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATSEQDENSE)); 21*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 22*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&b)); 23*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&x)); 24*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&u)); 25*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&e)); 26*9566063dSJacob Faibussowitsch PetscCall(VecSet(u,1.0)); 27*9566063dSJacob Faibussowitsch PetscCall(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; 32*9566063dSJacob Faibussowitsch PetscCall(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.; 36*9566063dSJacob Faibussowitsch PetscCall(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.; 39*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES)); 40c4762a1bSJed Brown 41*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 42*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 43c4762a1bSJed Brown 44*9566063dSJacob Faibussowitsch PetscCall(MatMult(C,u,b)); 45c4762a1bSJed Brown 46c4762a1bSJed Brown for (i=0; i<n; i++) { 47*9566063dSJacob Faibussowitsch PetscCall(MatSOR(C,b,omega,SOR_FORWARD_SWEEP,0.0,1,1,x)); 48*9566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e,-1.0,x,u)); 49*9566063dSJacob Faibussowitsch PetscCall(VecNorm(e,NORM_2,&norm)); 50*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"2-norm of error %g\n",(double)norm)); 51c4762a1bSJed Brown } 52*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 53*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 54*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 55*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 56*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 57*9566063dSJacob 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