1c4762a1bSJed Brown static const char help[] = "Test parallel assembly of SBAIJ matrices\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 59371c9d4SSatish Balay PetscErrorCode Assemble(MPI_Comm comm, PetscInt bs, MatType mtype) { 6c4762a1bSJed Brown const PetscInt rc[] = {0, 1, 2, 3}; 79371c9d4SSatish Balay const PetscScalar vals[] = {100, 2, 3, 4, 5, 600, 7, 8, 9, 100, 11, 1200, 13, 14, 15, 1600, 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 2800, 29, 30, 31, 32, 89371c9d4SSatish Balay 33, 34, 35, 36, 37, 38, 39, 40, 41, 42, 43, 44, 45, 46, 47, 48, 49, 50, 51, 52, 53, 54, 55, 56, 57, 58, 49, 60, 61, 62, 63, 64}; 9c4762a1bSJed Brown Mat A; 10c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) || defined(PETSC_HAVE_MKL_CPARDISO) 11c4762a1bSJed Brown Mat F; 1283863137SStefano Zampini MatSolverType stype = MATSOLVERPETSC; 13c4762a1bSJed Brown PetscRandom rdm; 14c4762a1bSJed Brown Vec b, x, y; 15c4762a1bSJed Brown PetscInt i, j; 16560a203cSprj- PetscReal norm2, tol = 100 * PETSC_SMALL; 17c4762a1bSJed Brown PetscBool issbaij; 18c4762a1bSJed Brown #endif 19c4762a1bSJed Brown PetscViewer viewer; 20c4762a1bSJed Brown 21c4762a1bSJed Brown PetscFunctionBegin; 229566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &A)); 239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 4 * bs, 4 * bs)); 249566063dSJacob Faibussowitsch PetscCall(MatSetType(A, mtype)); 259566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A, bs, 2, NULL, 2, NULL)); 269566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(A, bs, 2, NULL, 2, NULL)); 279566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 28c4762a1bSJed Brown /* All processes contribute a global matrix */ 299566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(A, 4, rc, 4, rc, vals, ADD_VALUES)); 309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Matrix %s(%" PetscInt_FMT ")\n", mtype, bs)); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(comm, &viewer)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL)); 359566063dSJacob Faibussowitsch PetscCall(MatView(A, viewer)); 369566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 379566063dSJacob Faibussowitsch PetscCall(MatView(A, viewer)); 38c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) || defined(PETSC_HAVE_MKL_CPARDISO) 399566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mtype, MATMPISBAIJ, &issbaij)); 40*48a46eb9SPierre Jolivet if (!issbaij) PetscCall(MatShift(A, 10)); 419566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 429566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 439566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, &y)); 449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &b)); 45c4762a1bSJed Brown for (j = 0; j < 2; j++) { 46c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 47c4762a1bSJed Brown if (j == 0) stype = MATSOLVERMUMPS; 48c4762a1bSJed Brown #else 49c4762a1bSJed Brown if (j == 0) continue; 50c4762a1bSJed Brown #endif 51c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_CPARDISO) 52c4762a1bSJed Brown if (j == 1) stype = MATSOLVERMKL_CPARDISO; 53c4762a1bSJed Brown #else 54c4762a1bSJed Brown if (j == 1) continue; 55c4762a1bSJed Brown #endif 56c4762a1bSJed Brown if (issbaij) { 579566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, stype, MAT_FACTOR_CHOLESKY, &F)); 589566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL)); 599566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, A, NULL)); 60c4762a1bSJed Brown } else { 619566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, stype, MAT_FACTOR_LU, &F)); 629566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL)); 639566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, A, NULL)); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown for (i = 0; i < 10; i++) { 669566063dSJacob Faibussowitsch PetscCall(VecSetRandom(b, rdm)); 679566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 68c4762a1bSJed Brown /* Check the error */ 699566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, x)); 709566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, -1.0, b)); 719566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm2)); 72*48a46eb9SPierre Jolivet if (norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error:MatSolve(), norm2: %g\n", (double)norm2)); 73c4762a1bSJed Brown } 749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 75c4762a1bSJed Brown } 769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 799566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 80c4762a1bSJed Brown #endif 819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 82c4762a1bSJed Brown PetscFunctionReturn(0); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 859371c9d4SSatish Balay int main(int argc, char *argv[]) { 86c4762a1bSJed Brown MPI_Comm comm; 87c4762a1bSJed Brown PetscMPIInt size; 88c4762a1bSJed Brown 89327415f7SBarry Smith PetscFunctionBeginUser; 909566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 91c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 93be096a46SBarry Smith PetscCheck(size == 2, comm, PETSC_ERR_USER, "This example must be run with exactly two processes"); 949566063dSJacob Faibussowitsch PetscCall(Assemble(comm, 2, MATMPIBAIJ)); 959566063dSJacob Faibussowitsch PetscCall(Assemble(comm, 2, MATMPISBAIJ)); 969566063dSJacob Faibussowitsch PetscCall(Assemble(comm, 1, MATMPIBAIJ)); 979566063dSJacob Faibussowitsch PetscCall(Assemble(comm, 1, MATMPISBAIJ)); 989566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 99b122ec5aSJacob Faibussowitsch return 0; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /*TEST 103c4762a1bSJed Brown 104c4762a1bSJed Brown test: 105c4762a1bSJed Brown nsize: 2 10697929ea7SJunchao Zhang args: -mat_ignore_lower_triangular 107c4762a1bSJed Brown filter: sed -e "s~mem [0-9]*~mem~g" 108c4762a1bSJed Brown 109c4762a1bSJed Brown TEST*/ 110