static char help[] = "Tests PC and KSP on a tridiagonal matrix. Note that most\n\ users should employ the KSP interface instead of using PC directly.\n\n"; #include int main(int argc,char **args) { Mat mat; /* matrix */ Vec b,ustar,u; /* vectors (RHS, exact solution, approx solution) */ PC pc; /* PC context */ KSP ksp; /* KSP context */ PetscErrorCode ierr; PetscInt n = 10,i,its,col[3]; PetscScalar value[3]; PCType pcname; KSPType kspname; PetscReal norm,tol=1000.*PETSC_MACHINE_EPSILON; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; /* Create and initialize vectors */ ierr = VecCreateSeq(PETSC_COMM_SELF,n,&b);CHKERRQ(ierr); ierr = VecCreateSeq(PETSC_COMM_SELF,n,&ustar);CHKERRQ(ierr); ierr = VecCreateSeq(PETSC_COMM_SELF,n,&u);CHKERRQ(ierr); ierr = VecSet(ustar,1.0);CHKERRQ(ierr); ierr = VecSet(u,0.0);CHKERRQ(ierr); /* Create and assemble matrix */ ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&mat);CHKERRQ(ierr); value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; for (i=1; i tol) { ierr = PetscPrintf(PETSC_COMM_SELF,"2 norm of error %g Number of iterations %D\n",(double)norm,its);CHKERRQ(ierr); } /* Free data structures */ ierr = KSPDestroy(&ksp);CHKERRQ(ierr); ierr = VecDestroy(&u);CHKERRQ(ierr); ierr = VecDestroy(&ustar);CHKERRQ(ierr); ierr = VecDestroy(&b);CHKERRQ(ierr); ierr = MatDestroy(&mat);CHKERRQ(ierr); ierr = PCDestroy(&pc);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; } /*TEST test: args: -ksp_type cg -ksp_monitor_short TEST*/