1*c0c276a7Ssdargavi static char help[] = "Testing MatCreateMPIAIJWithSeqAIJ().\n\n"; 2*c0c276a7Ssdargavi 3*c0c276a7Ssdargavi #include <petscmat.h> 4*c0c276a7Ssdargavi 5*c0c276a7Ssdargavi int main(int argc, char **argv) 6*c0c276a7Ssdargavi { 7*c0c276a7Ssdargavi Mat A, B; 8*c0c276a7Ssdargavi PetscInt i, j, column, M, N, m, n; 9*c0c276a7Ssdargavi PetscInt *oi, *oj, nd; 10*c0c276a7Ssdargavi const PetscInt *garray; 11*c0c276a7Ssdargavi PetscInt *garray_h; 12*c0c276a7Ssdargavi PetscScalar value; 13*c0c276a7Ssdargavi PetscScalar *oa; 14*c0c276a7Ssdargavi PetscRandom rctx; 15*c0c276a7Ssdargavi PetscBool equal, done; 16*c0c276a7Ssdargavi Mat AA, AB; 17*c0c276a7Ssdargavi PetscMPIInt size, rank; 18*c0c276a7Ssdargavi MatType mat_type; 19*c0c276a7Ssdargavi 20*c0c276a7Ssdargavi PetscFunctionBeginUser; 21*c0c276a7Ssdargavi PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 22*c0c276a7Ssdargavi PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 23*c0c276a7Ssdargavi PetscCheck(size > 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run with 2 or more processes"); 24*c0c276a7Ssdargavi PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 25*c0c276a7Ssdargavi 26*c0c276a7Ssdargavi /* Create a mpiaij matrix for checking */ 27*c0c276a7Ssdargavi PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, 5, 5, PETSC_DECIDE, PETSC_DECIDE, 0, NULL, 0, NULL, &A)); 28*c0c276a7Ssdargavi PetscCall(MatSetFromOptions(A)); 29*c0c276a7Ssdargavi PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE)); 30*c0c276a7Ssdargavi PetscCall(MatSetUp(A)); 31*c0c276a7Ssdargavi PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 32*c0c276a7Ssdargavi PetscCall(PetscRandomSetFromOptions(rctx)); 33*c0c276a7Ssdargavi 34*c0c276a7Ssdargavi for (i = 5 * rank; i < 5 * rank + 5; i++) { 35*c0c276a7Ssdargavi for (j = 0; j < 5 * size; j++) { 36*c0c276a7Ssdargavi PetscCall(PetscRandomGetValue(rctx, &value)); 37*c0c276a7Ssdargavi column = (PetscInt)(5 * size * PetscRealPart(value)); 38*c0c276a7Ssdargavi PetscCall(PetscRandomGetValue(rctx, &value)); 39*c0c276a7Ssdargavi PetscCall(MatSetValues(A, 1, &i, 1, &column, &value, INSERT_VALUES)); 40*c0c276a7Ssdargavi } 41*c0c276a7Ssdargavi } 42*c0c276a7Ssdargavi PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 43*c0c276a7Ssdargavi PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 44*c0c276a7Ssdargavi PetscCall(MatGetSize(A, &M, &N)); 45*c0c276a7Ssdargavi PetscCall(MatGetLocalSize(A, &m, &n)); 46*c0c276a7Ssdargavi 47*c0c276a7Ssdargavi PetscCall(MatMPIAIJGetSeqAIJ(A, &AA, &AB, &garray)); 48*c0c276a7Ssdargavi 49*c0c276a7Ssdargavi Mat output_mat_local, output_mat_nonlocal, output_mat_local_copy, output_mat_nonlocal_copy; 50*c0c276a7Ssdargavi 51*c0c276a7Ssdargavi PetscCall(MatConvert(AA, MATSAME, MAT_INITIAL_MATRIX, &output_mat_local)); 52*c0c276a7Ssdargavi PetscCall(MatConvert(AB, MATSAME, MAT_INITIAL_MATRIX, &output_mat_nonlocal)); 53*c0c276a7Ssdargavi PetscCall(MatConvert(AA, MATSAME, MAT_INITIAL_MATRIX, &output_mat_local_copy)); 54*c0c276a7Ssdargavi 55*c0c276a7Ssdargavi // The garray passed in has to be on the host, but it can be created 56*c0c276a7Ssdargavi // on device and copied to the host 57*c0c276a7Ssdargavi // We're just going to copy the existing host values here 58*c0c276a7Ssdargavi PetscInt nonlocalCols; 59*c0c276a7Ssdargavi PetscCall(MatGetLocalSize(AB, NULL, &nonlocalCols)); 60*c0c276a7Ssdargavi PetscCall(PetscMalloc1(nonlocalCols, &garray_h)); 61*c0c276a7Ssdargavi for (int i = 0; i < nonlocalCols; i++) { garray_h[i] = garray[i]; } 62*c0c276a7Ssdargavi 63*c0c276a7Ssdargavi // Build our MPI matrix 64*c0c276a7Ssdargavi // If we provide garray and output_mat_nonlocal with local indices and the compactified size 65*c0c276a7Ssdargavi // it doesn't compactify 66*c0c276a7Ssdargavi PetscCall(MatCreateMPIAIJWithSeqAIJ(PETSC_COMM_WORLD, M, N, output_mat_local, output_mat_nonlocal, garray_h, &B)); 67*c0c276a7Ssdargavi 68*c0c276a7Ssdargavi PetscCall(MatEqual(A, B, &equal)); 69*c0c276a7Ssdargavi PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Likely a bug in MatCreateMPIAIJWithSeqAIJ()"); 70*c0c276a7Ssdargavi PetscCall(MatDestroy(&B)); 71*c0c276a7Ssdargavi 72*c0c276a7Ssdargavi // ~~~~~~~~~~~~~~~~~ 73*c0c276a7Ssdargavi // Test MatCreateMPIAIJWithSeqAIJ with compactification 74*c0c276a7Ssdargavi // This is just for testing - would be silly to do this in practice 75*c0c276a7Ssdargavi // ~~~~~~~~~~~~~~~~~ 76*c0c276a7Ssdargavi garray_h = NULL; 77*c0c276a7Ssdargavi PetscCall(MatGetRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nd, (const PetscInt **)&oi, (const PetscInt **)&oj, &done)); 78*c0c276a7Ssdargavi PetscCall(MatSeqAIJGetArray(AB, &oa)); 79*c0c276a7Ssdargavi 80*c0c276a7Ssdargavi // Create a version of AB of size N with global indices 81*c0c276a7Ssdargavi PetscCall(MatGetType(AB, &mat_type)); 82*c0c276a7Ssdargavi PetscCall(MatCreate(PETSC_COMM_SELF, &output_mat_nonlocal_copy)); 83*c0c276a7Ssdargavi PetscCall(MatSetSizes(output_mat_nonlocal_copy, m, N, m, N)); 84*c0c276a7Ssdargavi PetscCall(MatSetType(output_mat_nonlocal_copy, mat_type)); 85*c0c276a7Ssdargavi PetscCall(MatSeqAIJSetPreallocation(output_mat_nonlocal_copy, oi[5], NULL)); 86*c0c276a7Ssdargavi 87*c0c276a7Ssdargavi // Fill the matrix 88*c0c276a7Ssdargavi for (int i = 0; i < 5; i++) { 89*c0c276a7Ssdargavi for (int j = 0; j < oi[i + 1] - oi[i]; j++) { PetscCall(MatSetValue(output_mat_nonlocal_copy, i, garray[oj[oi[i] + j]], oa[oi[i] + j], INSERT_VALUES)); } 90*c0c276a7Ssdargavi } 91*c0c276a7Ssdargavi PetscCall(MatAssemblyBegin(output_mat_nonlocal_copy, MAT_FINAL_ASSEMBLY)); 92*c0c276a7Ssdargavi PetscCall(MatAssemblyEnd(output_mat_nonlocal_copy, MAT_FINAL_ASSEMBLY)); 93*c0c276a7Ssdargavi 94*c0c276a7Ssdargavi PetscCall(MatRestoreRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nd, (const PetscInt **)&oi, (const PetscInt **)&oj, &done)); 95*c0c276a7Ssdargavi PetscCall(MatSeqAIJRestoreArray(AB, &oa)); 96*c0c276a7Ssdargavi 97*c0c276a7Ssdargavi // Build our MPI matrix 98*c0c276a7Ssdargavi // If we don't provide garray and output_mat_local_copy with global indices and size N 99*c0c276a7Ssdargavi // it will do compactification 100*c0c276a7Ssdargavi PetscCall(MatCreateMPIAIJWithSeqAIJ(PETSC_COMM_WORLD, M, N, output_mat_local_copy, output_mat_nonlocal_copy, garray_h, &B)); 101*c0c276a7Ssdargavi 102*c0c276a7Ssdargavi PetscCall(MatEqual(A, B, &equal)); 103*c0c276a7Ssdargavi PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Likely a bug in MatCreateMPIAIJWithSeqAIJ()"); 104*c0c276a7Ssdargavi PetscCall(MatDestroy(&B)); 105*c0c276a7Ssdargavi 106*c0c276a7Ssdargavi // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 107*c0c276a7Ssdargavi 108*c0c276a7Ssdargavi /* Free spaces */ 109*c0c276a7Ssdargavi PetscCall(PetscRandomDestroy(&rctx)); 110*c0c276a7Ssdargavi PetscCall(MatDestroy(&A)); 111*c0c276a7Ssdargavi PetscCall(PetscFinalize()); 112*c0c276a7Ssdargavi return 0; 113*c0c276a7Ssdargavi } 114*c0c276a7Ssdargavi 115*c0c276a7Ssdargavi /*TEST 116*c0c276a7Ssdargavi 117*c0c276a7Ssdargavi test: 118*c0c276a7Ssdargavi nsize: 2 119*c0c276a7Ssdargavi args: -mat_type aij 120*c0c276a7Ssdargavi output_file: output/empty.out 121*c0c276a7Ssdargavi 122*c0c276a7Ssdargavi TEST*/ 123