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