Actual source code: ex135.c
petsc-3.13.6 2020-09-29
1: static const char help[] = "Test parallel assembly of SBAIJ matrices\n\n";
3: #include <petscmat.h>
5: PetscErrorCode Assemble(MPI_Comm comm,PetscInt n,MatType mtype)
6: {
7: Mat A;
8: PetscInt first,last,i;
10: PetscMPIInt rank,size;
13: MatCreate(PETSC_COMM_WORLD,&A);
14: MatSetSizes(A, PETSC_DECIDE,PETSC_DECIDE,n,n);
15: MatSetType(A,MATMPISBAIJ);
16: MatSetFromOptions(A);
17: MPI_Comm_size(comm,&size);
18: MPI_Comm_rank(comm,&rank);
19: if (rank < size-1) {
20: MatMPISBAIJSetPreallocation(A,1,1,NULL,1,NULL);
21: } else {
22: MatMPISBAIJSetPreallocation(A,1,2,NULL,0,NULL);
23: }
24: MatGetOwnershipRange(A,&first,&last);
25: MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
26: last--;
27: for (i=first; i<=last; i++) {
28: MatSetValue(A,i,i,2.,INSERT_VALUES);
29: if (i != n-1) {MatSetValue(A,i,n-1,-1.,INSERT_VALUES);}
30: }
31: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
32: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
33: MatDestroy(&A);
34: return(0);
35: }
37: int main(int argc,char *argv[])
38: {
40: MPI_Comm comm;
41: PetscInt n = 6;
43: PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
44: comm = PETSC_COMM_WORLD;
45: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
46: Assemble(comm,n,MATMPISBAIJ);
47: PetscFinalize();
48: return ierr;
49: }
52: /*TEST
54: test:
55: nsize: 4
56: args: -n 1000 -mat_view ascii::ascii_info_detail -vecscatter_type sf
57: requires: double !complex !define(PETSC_USE_64BIT_INDICES)
59: TEST*/