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*327415f7SBarry Smith PetscFunctionBeginUser; 159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-omega",&omega,NULL)); 179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 18c4762a1bSJed Brown 199566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&C)); 209566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,n,n,n,n)); 219566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATSEQDENSE)); 229566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 239566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&b)); 249566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&x)); 259566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&u)); 269566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&e)); 279566063dSJacob Faibussowitsch PetscCall(VecSet(u,1.0)); 289566063dSJacob Faibussowitsch PetscCall(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; 339566063dSJacob Faibussowitsch PetscCall(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.; 379566063dSJacob Faibussowitsch PetscCall(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.; 409566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,2,midx,v,INSERT_VALUES)); 41c4762a1bSJed Brown 429566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 439566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(MatMult(C,u,b)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown for (i=0; i<n; i++) { 489566063dSJacob Faibussowitsch PetscCall(MatSOR(C,b,omega,SOR_FORWARD_SWEEP,0.0,1,1,x)); 499566063dSJacob Faibussowitsch PetscCall(VecWAXPY(e,-1.0,x,u)); 509566063dSJacob Faibussowitsch PetscCall(VecNorm(e,NORM_2,&norm)); 519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"2-norm of error %g\n",(double)norm)); 52c4762a1bSJed Brown } 539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 589566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 59b122ec5aSJacob Faibussowitsch return 0; 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62c4762a1bSJed Brown /*TEST 63c4762a1bSJed Brown 64c4762a1bSJed Brown test: 65c4762a1bSJed Brown 66c4762a1bSJed Brown TEST*/ 67