Actual source code: ex86.c
petsc-3.4.5 2014-06-29
1: static char help[] = "Testing MatCreateMPIAIJConcatenateSeqAIJ().\n\n";
3: #include <petscmat.h>
6: int main(int argc,char **argv)
7: {
9: Mat seqaijmat,mpiaijmat;
10: PetscMPIInt rank;
11: PetscScalar value[3];
12: PetscInt i,col[3],n=10;
14: PetscInitialize(&argc,&argv,(char*)0,help);
15: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
17: /* Create seqaij matrices of size n by (n+rank) */
18: MatCreate(PETSC_COMM_SELF,&seqaijmat);
19: MatSetSizes(seqaijmat,n+rank,PETSC_DECIDE,PETSC_DECIDE,n);
20: MatSetFromOptions(seqaijmat);
21: MatSeqAIJSetPreallocation(seqaijmat,3,NULL);
23: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
24: for (i=1; i<n-1; i++) {
25: col[0] = i-1; col[1] = i; col[2] = i+1;
26: MatSetValues(seqaijmat,1,&i,3,col,value,INSERT_VALUES);
27: }
28: i = n - 1; col[0] = n - 2; col[1] = n - 1;
30: MatSetValues(seqaijmat,1,&i,2,col,value,INSERT_VALUES);
32: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
34: MatSetValues(seqaijmat,1,&i,2,col,value,INSERT_VALUES);
35: MatAssemblyBegin(seqaijmat,MAT_FINAL_ASSEMBLY);
36: MatAssemblyEnd(seqaijmat,MAT_FINAL_ASSEMBLY);
38: /* Concatenate seqaij matrices into a single mpiaij matrix */
39: MatCreateMPIAIJConcatenateSeqAIJ(PETSC_COMM_WORLD,seqaijmat,PETSC_DECIDE,MAT_INITIAL_MATRIX,&mpiaijmat);
41: MatDestroy(&seqaijmat);
42: MatDestroy(&mpiaijmat);
43: PetscFinalize();
44: return 0;
45: }