Actual source code: ex110.c
petsc-3.6.1 2015-08-06
1: static char help[] = "Testing MatCreateMPIAIJWithSplitArrays().\n\n";
3: #include <petscmat.h>
4: #include <../src/mat/impls/aij/mpi/mpiaij.h>
8: int main(int argc,char **argv)
9: {
10: Mat A,B;
11: PetscInt i,j,column;
12: PetscInt *di,*dj,*oi,*oj;
13: PetscScalar *oa,*da,value;
14: PetscRandom rctx;
16: PetscBool equal;
17: Mat_SeqAIJ *daij,*oaij;
18: Mat_MPIAIJ *Ampiaij;
19: PetscMPIInt size,rank;
21: PetscInitialize(&argc,&argv,(char*)0,help);
22: MPI_Comm_size(PETSC_COMM_WORLD,&size);
23: if (size == 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Must run with 2 or more processes");
24: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
26: /* Create a mpiaij matrix for checking */
27: MatCreateAIJ(PETSC_COMM_WORLD,5,5,PETSC_DECIDE,PETSC_DECIDE,0,NULL,0,NULL,&A);
28: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
29: MatSetUp(A);
30: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
31: PetscRandomSetFromOptions(rctx);
33: for (i=5*rank; i<5*rank+5; i++) {
34: for (j=0; j<5*size; j++) {
35: PetscRandomGetValue(rctx,&value);
36: column = (PetscInt) (5*size*PetscRealPart(value));
37: PetscRandomGetValue(rctx,&value);
38: MatSetValues(A,1,&i,1,&column,&value,INSERT_VALUES);
39: }
40: }
41: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
42: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
44: Ampiaij = (Mat_MPIAIJ*) A->data;
45: daij = (Mat_SeqAIJ*) Ampiaij->A->data;
46: oaij = (Mat_SeqAIJ*) Ampiaij->B->data;
48: di = daij->i;
49: dj = daij->j;
50: da = daij->a;
52: oi = oaij->i;
53: oa = oaij->a;
54: PetscMalloc1(oi[5],&oj);
55: PetscMemcpy(oj,oaij->j,oi[5]*sizeof(PetscInt));
56: /* modify the column entries in the non-diagonal portion back to global numbering */
57: for (i=0; i<oi[5]; i++) {
58: oj[i] = Ampiaij->garray[oj[i]];
59: }
61: MatCreateMPIAIJWithSplitArrays(PETSC_COMM_WORLD,5,5,PETSC_DETERMINE,PETSC_DETERMINE,di,dj,da,oi,oj,oa,&B);
62: MatSetUp(B);
63: MatEqual(A,B,&equal);
65: if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Likely a bug in MatCreateMPIAIJWithSplitArrays()");
67: /* Free spaces */
68: PetscRandomDestroy(&rctx);
69: MatDestroy(&A);
70: MatDestroy(&B);
71: PetscFree(oj);
72: PetscFinalize();
73: return(0);
74: }