Actual source code: ex110.c
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: Mat A,B;
10: PetscInt i,j,column;
11: PetscInt *di,*dj,*oi,*oj;
12: PetscScalar *oa,*da,value;
13: PetscRandom rctx;
15: PetscTruth equal;
16: Mat_SeqAIJ *daij,*oaij;
17: Mat_MPIAIJ *Ampiaij;
18: PetscMPIInt size,rank;
20: PetscInitialize(&argc,&argv,(char *)0,help);
21: MPI_Comm_size(PETSC_COMM_WORLD,&size);
22: if (size == 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must run with 2 or more processes");
23: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
25: /* Create a mpiaij matrix for checking */
26: MatCreateMPIAIJ(PETSC_COMM_WORLD,5,5,PETSC_DECIDE,PETSC_DECIDE,0,PETSC_NULL,0,PETSC_NULL,&A);
27: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
28: PetscRandomSetFromOptions(rctx);
29:
30: for (i=5*rank; i<5*rank+5; i++) {
31: for (j=0; j<5*size; j++){
32: PetscRandomGetValue(rctx,&value);
33: column = (PetscInt) (5*size*PetscRealPart(value));
34: PetscRandomGetValue(rctx,&value);
35: MatSetValues(A,1,&i,1,&column,&value,INSERT_VALUES);
36: }
37: }
38: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
39: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
41: Ampiaij = (Mat_MPIAIJ*) A->data;
42: daij = (Mat_SeqAIJ*) Ampiaij->A->data;
43: oaij = (Mat_SeqAIJ*) Ampiaij->B->data;
44:
45: di = daij->i;
46: dj = daij->j;
47: da = daij->a;
49: oi = oaij->i;
50: oa = oaij->a;
51: PetscMalloc(oi[5]*sizeof(PetscInt),&oj);
52: PetscMemcpy(oj,oaij->j,oi[5]*sizeof(PetscInt));
53: /* modify the column entries in the non-diagonal portion back to global numbering */
54: for (i=0; i<oi[5]; i++) {
55: oj[i] = Ampiaij->garray[oj[i]];
56: }
58: MatCreateMPIAIJWithSplitArrays(PETSC_COMM_WORLD,5,5,PETSC_DETERMINE,PETSC_DETERMINE,di,dj,da,oi,oj,oa,&B);
60: MatEqual(A,B,&equal);
62: if (!equal) SETERRQ(PETSC_ERR_PLIB,"Likely a bug in MatCreateMPIAIJWithSplitArrays()");
64: /* Free spaces */
65: PetscRandomDestroy(rctx);
66: MatDestroy(A);
67: MatDestroy(B);
68: PetscFree(oj);
69: PetscFinalize();
70: return(0);
71: }