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