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; 20ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 21*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 22c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 23c4762a1bSJed Brown ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr); 24c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-symmetric",&flg);CHKERRQ(ierr); 27c4762a1bSJed Brown if (flg) { /* Treat matrix as symmetric only if we set this flag */ 28c4762a1bSJed Brown ierr = MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);CHKERRQ(ierr); 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; 36c4762a1bSJed Brown if (i>0) {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 37c4762a1bSJed Brown if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 38c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 39c4762a1bSJed Brown if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 40c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown } 43c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = ISView(perm,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,m*n,&u);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = VecSet(u,1.0);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = VecDuplicate(u,&x);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = VecDuplicate(u,&b);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = VecDuplicate(u,&y);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = MatMult(C,u,b);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = VecCopy(b,y);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = VecScale(y,2.0);CHKERRQ(ierr); 56c4762a1bSJed Brown 57c4762a1bSJed Brown ierr = MatNorm(C,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Frobenius norm of matrix %g\n",(double)norm);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = MatNorm(C,NORM_1,&norm);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"One norm of matrix %g\n",(double)norm);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = MatNorm(C,NORM_INFINITY,&norm);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Infinity norm of matrix %g\n",(double)norm);CHKERRQ(ierr); 63c4762a1bSJed Brown 64c4762a1bSJed Brown ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 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 70c4762a1bSJed Brown ierr = MatLUFactor(C,perm,iperm,&info);CHKERRQ(ierr); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* Test MatSolve */ 73c4762a1bSJed Brown ierr = MatSolve(C,b,x);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = VecView(b,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = VecAXPY(x,-1.0,u);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr); 78c4762a1bSJed Brown if (norm > tol) { 79c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve: Norm of error %g\n",(double)norm);CHKERRQ(ierr); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Test MatSolveAdd */ 83c4762a1bSJed Brown ierr = MatSolveAdd(C,b,y,x);CHKERRQ(ierr); 84c4762a1bSJed Brown ierr = VecAXPY(x,-1.0,y);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = VecAXPY(x,-1.0,u);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr); 87c4762a1bSJed Brown if (norm > tol) { 88c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolveAdd(): Norm of error %g\n",(double)norm);CHKERRQ(ierr); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown ierr = ISDestroy(&perm);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = ISDestroy(&iperm);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 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