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