static char help[] = "Test using MatPreallocator to preallocate two matrices\n\n"; /*T Concepts: Mat^matrix preallocation Processors: n T*/ #include PetscErrorCode ex1_nonsquare_bs1(void) { Mat A,A_duplicate,preallocator; PetscInt M,N,m,n,bs; PetscErrorCode ierr; PetscFunctionBegin; M = 10; N = 8; ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); ierr = MatSetBlockSize(A,1);CHKERRQ(ierr); ierr = MatSetFromOptions(A);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_WORLD,&A_duplicate);CHKERRQ(ierr); ierr = MatSetType(A_duplicate,MATAIJ);CHKERRQ(ierr); ierr = MatSetSizes(A_duplicate,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); ierr = MatSetBlockSize(A_duplicate,2);CHKERRQ(ierr); ierr = MatSetFromOptions(A_duplicate);CHKERRQ(ierr); ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); ierr = MatCreate(PetscObjectComm((PetscObject)A),&preallocator);CHKERRQ(ierr); ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); ierr = MatSetSizes(preallocator,m,n,M,N);CHKERRQ(ierr); ierr = MatSetBlockSize(preallocator,bs);CHKERRQ(ierr); ierr = MatSetUp(preallocator);CHKERRQ(ierr); ierr = MatSetValue(preallocator,3,3,0,INSERT_VALUES);CHKERRQ(ierr); ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A);CHKERRQ(ierr); ierr = MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A_duplicate);CHKERRQ(ierr); ierr = MatDestroy(&preallocator);CHKERRQ(ierr); ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = MatSetValue(A,3,3,0.3,INSERT_VALUES);CHKERRQ(ierr); ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = MatDestroy(&A);CHKERRQ(ierr); ierr = MatView(A_duplicate,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = MatSetValue(A_duplicate,3,3,0.3,INSERT_VALUES);CHKERRQ(ierr); ierr = MatAssemblyBegin(A_duplicate,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(A_duplicate,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatView(A_duplicate,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = MatDestroy(&A_duplicate);CHKERRQ(ierr); PetscFunctionReturn(0); } int main(int argc, char **args) { PetscErrorCode ierr; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; ierr = ex1_nonsquare_bs1();CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; } /*TEST test: nsize: {{1 2}separate output} TEST*/