xref: /petsc/src/mat/tests/ex110.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Testing MatCreateMPIAIJWithSplitArrays().\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat             A, B;
8c4762a1bSJed Brown   PetscInt        i, j, column, *ooj;
9c4762a1bSJed Brown   PetscInt       *di, *dj, *oi, *oj, nd;
10c4762a1bSJed Brown   const PetscInt *garray;
11c4762a1bSJed Brown   PetscScalar    *oa, *da;
12c4762a1bSJed Brown   PetscScalar     value;
13c4762a1bSJed Brown   PetscRandom     rctx;
14c4762a1bSJed Brown   PetscBool       equal, done;
15c4762a1bSJed Brown   Mat             AA, AB;
16c4762a1bSJed Brown   PetscMPIInt     size, rank;
17c4762a1bSJed Brown 
18327415f7SBarry Smith   PetscFunctionBeginUser;
19c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
21be096a46SBarry Smith   PetscCheck(size > 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run with 2 or more processes");
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* Create a mpiaij matrix for checking */
259566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, 5, 5, PETSC_DECIDE, PETSC_DECIDE, 0, NULL, 0, NULL, &A));
269566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
279566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
289566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx));
299566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   for (i = 5 * rank; i < 5 * rank + 5; i++) {
32c4762a1bSJed Brown     for (j = 0; j < 5 * size; j++) {
339566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rctx, &value));
34c4762a1bSJed Brown       column = (PetscInt)(5 * size * PetscRealPart(value));
359566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rctx, &value));
369566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 1, &column, &value, INSERT_VALUES));
37c4762a1bSJed Brown     }
38c4762a1bSJed Brown   }
399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
41c4762a1bSJed Brown 
429566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJGetSeqAIJ(A, &AA, &AB, &garray));
439566063dSJacob Faibussowitsch   PetscCall(MatGetRowIJ(AA, 0, PETSC_FALSE, PETSC_FALSE, &nd, (const PetscInt **)&di, (const PetscInt **)&dj, &done));
449566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArray(AA, &da));
459566063dSJacob Faibussowitsch   PetscCall(MatGetRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nd, (const PetscInt **)&oi, (const PetscInt **)&oj, &done));
469566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArray(AB, &oa));
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(oi[5], &ooj));
499566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(ooj, oj, oi[5]));
50c4762a1bSJed Brown   /* modify the column entries in the non-diagonal portion back to global numbering */
51ad540459SPierre Jolivet   for (i = 0; i < oi[5]; i++) ooj[i] = garray[ooj[i]];
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIAIJWithSplitArrays(PETSC_COMM_WORLD, 5, 5, PETSC_DETERMINE, PETSC_DETERMINE, di, dj, da, oi, ooj, oa, &B));
549566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
559566063dSJacob Faibussowitsch   PetscCall(MatEqual(A, B, &equal));
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowIJ(AA, 0, PETSC_FALSE, PETSC_FALSE, &nd, (const PetscInt **)&di, (const PetscInt **)&dj, &done));
589566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArray(AA, &da));
599566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nd, (const PetscInt **)&oi, (const PetscInt **)&oj, &done));
609566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArray(AB, &oa));
61c4762a1bSJed Brown 
6228b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Likely a bug in MatCreateMPIAIJWithSplitArrays()");
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Free spaces */
659566063dSJacob Faibussowitsch   PetscCall(PetscFree(ooj));
669566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
699566063dSJacob Faibussowitsch   PetscCall(PetscFree(oj));
709566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
71b122ec5aSJacob Faibussowitsch   return 0;
72c4762a1bSJed Brown }
73c4762a1bSJed Brown 
74c4762a1bSJed Brown /*TEST
75c4762a1bSJed Brown 
76c4762a1bSJed Brown    test:
77c4762a1bSJed Brown       nsize: 3
78*3886731fSPierre Jolivet       output_file: output/empty.out
79c4762a1bSJed Brown 
80c4762a1bSJed Brown TEST*/
81