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> 7*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 8*d71ae5a4SJacob Faibussowitsch { 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 PetscInt i, j, m = 3, n = 2, its; 15c4762a1bSJed Brown PetscMPIInt size, rank; 16c4762a1bSJed Brown PetscBool mat_nonsymmetric; 17c4762a1bSJed Brown PetscInt its_max; 18c4762a1bSJed Brown MatFactorInfo factinfo; 19c4762a1bSJed Brown IS perm, iperm; 20c4762a1bSJed Brown 21327415f7SBarry Smith PetscFunctionBeginUser; 229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 259566063dSJacob Faibussowitsch PetscCallMPI(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 */ 319566063dSJacob Faibussowitsch PetscCall(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 */ 399566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 409566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 419566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 429566063dSJacob Faibussowitsch PetscCall(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++) { 529371c9d4SSatish Balay v = -1.0; 539371c9d4SSatish Balay i = I / n; 549371c9d4SSatish Balay j = I - i * n; 559371c9d4SSatish Balay if (i > 0) { 569371c9d4SSatish Balay J = I - n; 579371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &I, 1, &J, &v, ADD_VALUES)); 589371c9d4SSatish Balay } 599371c9d4SSatish Balay if (i < m - 1) { 609371c9d4SSatish Balay J = I + n; 619371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &I, 1, &J, &v, ADD_VALUES)); 629371c9d4SSatish Balay } 639371c9d4SSatish Balay if (j > 0) { 649371c9d4SSatish Balay J = I - 1; 659371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &I, 1, &J, &v, ADD_VALUES)); 669371c9d4SSatish Balay } 679371c9d4SSatish Balay if (j < n - 1) { 689371c9d4SSatish Balay J = I + 1; 699371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &I, 1, &J, &v, ADD_VALUES)); 709371c9d4SSatish Balay } 719371c9d4SSatish Balay v = 4.0; 729371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &I, 1, &I, &v, ADD_VALUES)); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* 76c4762a1bSJed Brown Make the matrix nonsymmetric if desired 77c4762a1bSJed Brown */ 78c4762a1bSJed Brown if (mat_nonsymmetric) { 79c4762a1bSJed Brown for (I = Istart; I < Iend; I++) { 809371c9d4SSatish Balay v = -1.5; 819371c9d4SSatish Balay i = I / n; 829371c9d4SSatish Balay if (i > 1) { 839371c9d4SSatish Balay J = I - n - 1; 849371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &I, 1, &J, &v, ADD_VALUES)); 859371c9d4SSatish Balay } 86c4762a1bSJed Brown } 87c4762a1bSJed Brown } else { 889566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_SYMMETRIC, PETSC_TRUE)); 899566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* 93c4762a1bSJed Brown Assemble matrix, using the 2-step process: 94c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 95c4762a1bSJed Brown Computations can be done while messages are in transition 96c4762a1bSJed Brown by placing code between these two statements. 97c4762a1bSJed Brown */ 989566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 999566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown its_max = 1000; 102c4762a1bSJed Brown /* 103c4762a1bSJed Brown Create parallel vectors. 104c4762a1bSJed Brown - When using VecSetSizes(), we specify only the vector's global 105c4762a1bSJed Brown dimension; the parallel partitioning is determined at runtime. 106c4762a1bSJed Brown - Note: We form 1 vector from scratch and then duplicate as needed. 107c4762a1bSJed Brown */ 1089566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &u)); 1099566063dSJacob Faibussowitsch PetscCall(VecSetSizes(u, PETSC_DECIDE, m * n)); 1109566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(u)); 1119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &b)); 1129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b, &x)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* 115c4762a1bSJed Brown Currently, all parallel PETSc vectors are partitioned by 116c4762a1bSJed Brown contiguous chunks across the processors. Determine which 117c4762a1bSJed Brown range of entries are locally owned. 118c4762a1bSJed Brown */ 1199566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low, &high)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* 122c4762a1bSJed Brown Set elements within the exact solution vector in parallel. 123c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 124c4762a1bSJed Brown locally (but any non-local entries will be sent to the 125c4762a1bSJed Brown appropriate processor during vector assembly). 126c4762a1bSJed Brown - Always specify global locations of vector entries. 127c4762a1bSJed Brown */ 1289566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &ldim)); 129c4762a1bSJed Brown for (i = 0; i < ldim; i++) { 130c4762a1bSJed Brown iglobal = i + low; 131c4762a1bSJed Brown v = (PetscScalar)(i + 100 * rank); 1329566063dSJacob Faibussowitsch PetscCall(VecSetValues(u, 1, &iglobal, &v, INSERT_VALUES)); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* 136c4762a1bSJed Brown Assemble vector, using the 2-step process: 137c4762a1bSJed Brown VecAssemblyBegin(), VecAssemblyEnd() 138c4762a1bSJed Brown Computations can be done while messages are in transition, 139c4762a1bSJed Brown by placing code between these two statements. 140c4762a1bSJed Brown */ 1419566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(u)); 1429566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(u)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* Compute right-hand-side vector */ 1459566063dSJacob Faibussowitsch PetscCall(MatMult(C, u, b)); 146c4762a1bSJed Brown 1479566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGNATURAL, &perm, &iperm)); 148c4762a1bSJed Brown its_max = 2000; 149c4762a1bSJed Brown for (i = 0; i < its_max; i++) { 1509566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_LU, &F)); 1519566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, C, perm, iperm, &factinfo)); 15248a46eb9SPierre Jolivet for (j = 0; j < 1; j++) PetscCall(MatLUFactorNumeric(F, C, &factinfo)); 1539566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, x)); 1549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 155c4762a1bSJed Brown } 1569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iperm)); 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* Check the error */ 1609566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, none, u)); 1619566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm)); 1629566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %t\n", (double)norm)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* Free work space. */ 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 1669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1699566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 170b122ec5aSJacob Faibussowitsch return 0; 171c4762a1bSJed Brown } 172