Actual source code: ex98.c
petsc-3.13.6 2020-09-29
2: static char help[] = "Tests MatMPIAIJSetPreallocationCSR()\n\n";
4: /*T
5: Concepts: partitioning
6: Processors: 4
7: T*/
9: /*
10: Include "petscmat.h" so that we can use matrices. Note that this file
11: automatically includes:
12: petscsys.h - base PETSc routines petscvec.h - vectors
13: petscmat.h - matrices
14: petscis.h - index sets
15: petscviewer.h - viewers
16: */
17: #include <petscmat.h>
19: int main(int argc,char **args)
20: {
21: Mat A;
22: PetscInt *ia,*ja;
24: PetscMPIInt rank,size;
26: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
27: MPI_Comm_size(PETSC_COMM_WORLD,&size);
28: if (size != 4) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors");
29: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
31: PetscMalloc1(5,&ia);
32: PetscMalloc1(16,&ja);
33: if (!rank) {
34: ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6;
35: ja[8] = 2; ja[9] = 7;
36: ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
37: } else if (rank == 1) {
38: ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2;
39: ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11;
40: ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
41: } else if (rank == 2) {
42: ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6;
43: ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15;
44: ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
45: } else {
46: ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15;
47: ja[8] = 11; ja[9] = 14;
48: ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
49: }
51: MatCreate(PETSC_COMM_WORLD,&A);
52: MatSetSizes(A,4,4,16,16);
53: MatSetType(A,MATMPIAIJ);
54: MatMPIAIJSetPreallocationCSR(A,ia,ja,NULL);
55: PetscFree(ia);
56: PetscFree(ja);
57: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
58: MatDestroy(&A);
59: PetscFinalize();
60: return ierr;
61: }
63: /*TEST
65: test:
66: nsize: 4
68: TEST*/