1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Test repeated LU factorizations. Used for checking memory leak\n\ 3c4762a1bSJed Brown -m <size> : problem size\n\ 4c4762a1bSJed Brown -mat_nonsym : use nonsymmetric matrix (default is symmetric)\n\n"; 5c4762a1bSJed Brown 6c4762a1bSJed Brown #include <petscmat.h> 7c4762a1bSJed Brown int main(int argc,char **args) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown Mat C,F; /* matrix */ 10c4762a1bSJed Brown Vec x,u,b; /* approx solution, RHS, exact solution */ 11c4762a1bSJed Brown PetscReal norm; /* norm of solution error */ 12c4762a1bSJed Brown PetscScalar v,none = -1.0; 13c4762a1bSJed Brown PetscInt I,J,ldim,low,high,iglobal,Istart,Iend; 14c4762a1bSJed Brown PetscErrorCode ierr; 15c4762a1bSJed Brown PetscInt i,j,m = 3,n = 2,its; 16c4762a1bSJed Brown PetscMPIInt size,rank; 17c4762a1bSJed Brown PetscBool mat_nonsymmetric; 18c4762a1bSJed Brown PetscInt its_max; 19c4762a1bSJed Brown MatFactorInfo factinfo; 20c4762a1bSJed Brown IS perm,iperm; 21c4762a1bSJed Brown 22c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 24*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 25*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 26c4762a1bSJed Brown n = 2*size; 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* 29c4762a1bSJed Brown Set flag if we are doing a nonsymmetric problem; the default is symmetric. 30c4762a1bSJed Brown */ 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mat_nonsym",&mat_nonsymmetric)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* 34c4762a1bSJed Brown Create parallel matrix, specifying only its global dimensions. 35c4762a1bSJed Brown When using MatCreate(), the matrix format can be specified at 36c4762a1bSJed Brown runtime. Also, the parallel partitioning of the matrix is 37c4762a1bSJed Brown determined by PETSc at runtime. 38c4762a1bSJed Brown */ 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(C,&Istart,&Iend)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* 45c4762a1bSJed Brown Set matrix entries matrix in parallel. 46c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 47c4762a1bSJed Brown locally (but any non-local elements will be sent to the 48c4762a1bSJed Brown appropriate processor during matrix assembly). 49c4762a1bSJed Brown - Always specify global row and columns of matrix entries. 50c4762a1bSJed Brown */ 51c4762a1bSJed Brown for (I=Istart; I<Iend; I++) { 52c4762a1bSJed Brown v = -1.0; i = I/n; j = I - i*n; 53*5f80ce2aSJacob Faibussowitsch if (i>0) {J = I - n; CHKERRQ(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));} 54*5f80ce2aSJacob Faibussowitsch if (i<m-1) {J = I + n; CHKERRQ(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));} 55*5f80ce2aSJacob Faibussowitsch if (j>0) {J = I - 1; CHKERRQ(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));} 56*5f80ce2aSJacob Faibussowitsch if (j<n-1) {J = I + 1; CHKERRQ(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));} 57*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES)); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* 61c4762a1bSJed Brown Make the matrix nonsymmetric if desired 62c4762a1bSJed Brown */ 63c4762a1bSJed Brown if (mat_nonsymmetric) { 64c4762a1bSJed Brown for (I=Istart; I<Iend; I++) { 65c4762a1bSJed Brown v = -1.5; i = I/n; 66*5f80ce2aSJacob Faibussowitsch if (i>1) {J = I-n-1; CHKERRQ(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));} 67c4762a1bSJed Brown } 68c4762a1bSJed Brown } else { 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE)); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* 74c4762a1bSJed Brown Assemble matrix, using the 2-step process: 75c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 76c4762a1bSJed Brown Computations can be done while messages are in transition 77c4762a1bSJed Brown by placing code between these two statements. 78c4762a1bSJed Brown */ 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown its_max=1000; 83c4762a1bSJed Brown /* 84c4762a1bSJed Brown Create parallel vectors. 85c4762a1bSJed Brown - When using VecSetSizes(), we specify only the vector's global 86c4762a1bSJed Brown dimension; the parallel partitioning is determined at runtime. 87c4762a1bSJed Brown - Note: We form 1 vector from scratch and then duplicate as needed. 88c4762a1bSJed Brown */ 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&u)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(u,PETSC_DECIDE,m*n)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(u)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&b)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(b,&x)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* 96c4762a1bSJed Brown Currently, all parallel PETSc vectors are partitioned by 97c4762a1bSJed Brown contiguous chunks across the processors. Determine which 98c4762a1bSJed Brown range of entries are locally owned. 99c4762a1bSJed Brown */ 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(x,&low,&high)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* 103c4762a1bSJed Brown Set elements within the exact solution vector in parallel. 104c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 105c4762a1bSJed Brown locally (but any non-local entries will be sent to the 106c4762a1bSJed Brown appropriate processor during vector assembly). 107c4762a1bSJed Brown - Always specify global locations of vector entries. 108c4762a1bSJed Brown */ 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(x,&ldim)); 110c4762a1bSJed Brown for (i=0; i<ldim; i++) { 111c4762a1bSJed Brown iglobal = i + low; 112c4762a1bSJed Brown v = (PetscScalar)(i + 100*rank); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(u,1,&iglobal,&v,INSERT_VALUES)); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* 117c4762a1bSJed Brown Assemble vector, using the 2-step process: 118c4762a1bSJed Brown VecAssemblyBegin(), VecAssemblyEnd() 119c4762a1bSJed Brown Computations can be done while messages are in transition, 120c4762a1bSJed Brown by placing code between these two statements. 121c4762a1bSJed Brown */ 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(u)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(u)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* Compute right-hand-side vector */ 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(C,u,b)); 127c4762a1bSJed Brown 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(C,MATORDERINGNATURAL,&perm,&iperm)); 129c4762a1bSJed Brown its_max = 2000; 130c4762a1bSJed Brown for (i=0; i<its_max; i++) { 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&F)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorSymbolic(F,C,perm,iperm,&factinfo)); 133c4762a1bSJed Brown for (j=0; j<1; j++) { 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorNumeric(F,C,&factinfo)); 135c4762a1bSJed Brown } 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,x)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 138c4762a1bSJed Brown } 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iperm)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* Check the error */ 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(x,none,u)); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(x,NORM_2,&norm)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %t\n",(double)norm)); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Free work space. */ 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 152c4762a1bSJed Brown ierr = PetscFinalize(); 153c4762a1bSJed Brown return ierr; 154c4762a1bSJed Brown } 155