1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatNorm(), MatLUFactor(), MatSolve() and MatSolveAdd().\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 PetscInt i,j,m = 3,n = 3,Ii,J; 10c4762a1bSJed Brown PetscErrorCode ierr; 11c4762a1bSJed Brown PetscBool flg; 12c4762a1bSJed Brown PetscScalar v; 13c4762a1bSJed Brown IS perm,iperm; 14c4762a1bSJed Brown Vec x,u,b,y; 15c4762a1bSJed Brown PetscReal norm,tol=PETSC_SMALL; 16c4762a1bSJed Brown MatFactorInfo info; 17c4762a1bSJed Brown PetscMPIInt size; 18c4762a1bSJed Brown 19c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 20*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 212c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-symmetric",&flg)); 27c4762a1bSJed Brown if (flg) { /* Treat matrix as symmetric only if we set this flag */ 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE)); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE)); 30c4762a1bSJed Brown } 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Create the matrix for the five point stencil, YET AGAIN */ 33c4762a1bSJed Brown for (i=0; i<m; i++) { 34c4762a1bSJed Brown for (j=0; j<n; j++) { 35c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 36*5f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 37*5f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 38*5f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 39*5f80ce2aSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 40*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown } 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(perm,PETSC_VIEWER_STDOUT_SELF)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m*n,&u)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(u,1.0)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&x)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&b)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&y)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(C,u,b)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(b,y)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(y,2.0)); 56c4762a1bSJed Brown 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_FROBENIUS,&norm)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Frobenius norm of matrix %g\n",(double)norm)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_1,&norm)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"One norm of matrix %g\n",(double)norm)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_INFINITY,&norm)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Infinity norm of matrix %g\n",(double)norm)); 63c4762a1bSJed Brown 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorInfoInitialize(&info)); 65c4762a1bSJed Brown info.fill = 2.0; 66c4762a1bSJed Brown info.dtcol = 0.0; 67c4762a1bSJed Brown info.zeropivot = 1.e-14; 68c4762a1bSJed Brown info.pivotinblocks = 1.0; 69c4762a1bSJed Brown 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactor(C,perm,iperm,&info)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* Test MatSolve */ 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(C,b,x)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(b,PETSC_VIEWER_STDOUT_SELF)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_SELF)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(x,-1.0,u)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(x,NORM_2,&norm)); 78c4762a1bSJed Brown if (norm > tol) { 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatSolve: Norm of error %g\n",(double)norm)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Test MatSolveAdd */ 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolveAdd(C,b,y,x)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(x,-1.0,y)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(x,-1.0,u)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(x,NORM_2,&norm)); 87c4762a1bSJed Brown if (norm > tol) { 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatSolveAdd(): Norm of error %g\n",(double)norm)); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iperm)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 98c4762a1bSJed Brown ierr = PetscFinalize(); 99c4762a1bSJed Brown return ierr; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /*TEST 103c4762a1bSJed Brown 104c4762a1bSJed Brown test: 105c4762a1bSJed Brown 106c4762a1bSJed Brown TEST*/ 107