1c4762a1bSJed Brown static char help[] = "Example of using MatPreallocator\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /*T 4c4762a1bSJed Brown Concepts: Mat^matrix preallocation 5c4762a1bSJed Brown Processors: n 6c4762a1bSJed Brown T*/ 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include <petscmat.h> 9c4762a1bSJed Brown 10c4762a1bSJed Brown PetscErrorCode ex1_nonsquare_bs1(void) 11c4762a1bSJed Brown { 12c4762a1bSJed Brown Mat A,preallocator; 13c4762a1bSJed Brown PetscInt M,N,m,n,bs; 14c4762a1bSJed Brown PetscErrorCode ierr; 15c4762a1bSJed Brown 16c4762a1bSJed Brown /* 17c4762a1bSJed Brown Create the Jacobian matrix 18c4762a1bSJed Brown */ 19*362febeeSStefano Zampini PetscFunctionBegin; 20c4762a1bSJed Brown M = 10; 21c4762a1bSJed Brown N = 8; 22c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 23c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 24c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = MatSetBlockSize(A,1);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* 29c4762a1bSJed Brown Get the sizes of the jacobian. 30c4762a1bSJed Brown */ 31c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown Create a preallocator matrix with sizes (local and global) matching the jacobian A. 36c4762a1bSJed Brown */ 37c4762a1bSJed Brown ierr = MatCreate(PetscObjectComm((PetscObject)A),&preallocator);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = MatSetSizes(preallocator,m,n,M,N);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = MatSetBlockSize(preallocator,bs);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = MatSetUp(preallocator);CHKERRQ(ierr); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* 44c4762a1bSJed Brown Insert non-zero pattern (e.g. perform a sweep over the grid). 45c4762a1bSJed Brown You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue(). 46c4762a1bSJed Brown */ 47c4762a1bSJed Brown { 48c4762a1bSJed Brown PetscInt ii,jj; 49c4762a1bSJed Brown PetscScalar vv = 0.0; 50c4762a1bSJed Brown 51c4762a1bSJed Brown ii = 3; jj = 3; 52c4762a1bSJed Brown ierr = MatSetValues(preallocator,1,&ii,1,&jj,&vv,INSERT_VALUES);CHKERRQ(ierr); 53c4762a1bSJed Brown 54c4762a1bSJed Brown ii = 7; jj = 4; 55c4762a1bSJed Brown ierr = MatSetValues(preallocator,1,&ii,1,&jj,&vv,INSERT_VALUES);CHKERRQ(ierr); 56c4762a1bSJed Brown 57c4762a1bSJed Brown ii = 9; jj = 7; 58c4762a1bSJed Brown ierr = MatSetValue(preallocator,ii,jj,vv,INSERT_VALUES);CHKERRQ(ierr); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* 64c4762a1bSJed Brown Push the non-zero pattern defined within preallocator into A. 65c4762a1bSJed Brown Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE). 66c4762a1bSJed Brown The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically. 67c4762a1bSJed Brown */ 68c4762a1bSJed Brown ierr = MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A);CHKERRQ(ierr); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* 71c4762a1bSJed Brown We no longer require the preallocator object so destroy it. 72c4762a1bSJed Brown */ 73c4762a1bSJed Brown ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 74c4762a1bSJed Brown 75c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* 78c4762a1bSJed Brown Insert non-zero values into A. 79c4762a1bSJed Brown */ 80c4762a1bSJed Brown { 81c4762a1bSJed Brown PetscInt ii,jj; 82c4762a1bSJed Brown PetscScalar vv; 83c4762a1bSJed Brown 84c4762a1bSJed Brown ii = 3; jj = 3; vv = 0.3; 85c4762a1bSJed Brown ierr = MatSetValue(A,ii,jj,vv,INSERT_VALUES);CHKERRQ(ierr); 86c4762a1bSJed Brown 87c4762a1bSJed Brown ii = 7; jj = 4; vv = 3.3; 88c4762a1bSJed Brown ierr = MatSetValues(A,1,&ii,1,&jj,&vv,INSERT_VALUES);CHKERRQ(ierr); 89c4762a1bSJed Brown 90c4762a1bSJed Brown ii = 9; jj = 7; vv = 4.3; 91c4762a1bSJed Brown ierr = MatSetValues(A,1,&ii,1,&jj,&vv,INSERT_VALUES);CHKERRQ(ierr); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 95c4762a1bSJed Brown 96c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 97c4762a1bSJed Brown 98c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 99c4762a1bSJed Brown PetscFunctionReturn(0); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown PetscErrorCode ex2_square_bsvariable(void) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown Mat A,preallocator; 105c4762a1bSJed Brown PetscInt M,N,m,n,bs = 1; 106c4762a1bSJed Brown PetscErrorCode ierr; 107c4762a1bSJed Brown 108*362febeeSStefano Zampini PetscFunctionBegin; 109c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-block_size",&bs,NULL);CHKERRQ(ierr); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* 112c4762a1bSJed Brown Create the Jacobian matrix. 113c4762a1bSJed Brown */ 114c4762a1bSJed Brown M = 10 * bs; 115c4762a1bSJed Brown N = 10 * bs; 116c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* 123c4762a1bSJed Brown Get the sizes of the jacobian. 124c4762a1bSJed Brown */ 125c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* 129c4762a1bSJed Brown Create a preallocator matrix with dimensions matching the jacobian A. 130c4762a1bSJed Brown */ 131c4762a1bSJed Brown ierr = MatCreate(PetscObjectComm((PetscObject)A),&preallocator);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = MatSetSizes(preallocator,m,n,M,N);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = MatSetBlockSize(preallocator,bs);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = MatSetUp(preallocator);CHKERRQ(ierr); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* 138c4762a1bSJed Brown Insert non-zero pattern (e.g. perform a sweep over the grid). 139c4762a1bSJed Brown You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue(). 140c4762a1bSJed Brown */ 141c4762a1bSJed Brown { 142c4762a1bSJed Brown PetscInt ii,jj; 143c4762a1bSJed Brown PetscScalar *vv; 144c4762a1bSJed Brown 145c4762a1bSJed Brown ierr = PetscCalloc1(bs*bs,&vv);CHKERRQ(ierr); 146c4762a1bSJed Brown 147c4762a1bSJed Brown ii = 0; jj = 9; 148c4762a1bSJed Brown ierr = MatSetValue(preallocator,ii,jj,vv[0],INSERT_VALUES);CHKERRQ(ierr); 149c4762a1bSJed Brown 150c4762a1bSJed Brown ii = 0; jj = 0; 151c4762a1bSJed Brown ierr = MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES);CHKERRQ(ierr); 152c4762a1bSJed Brown 153c4762a1bSJed Brown ii = 2; jj = 4; 154c4762a1bSJed Brown ierr = MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES);CHKERRQ(ierr); 155c4762a1bSJed Brown 156c4762a1bSJed Brown ii = 4; jj = 2; 157c4762a1bSJed Brown ierr = MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES);CHKERRQ(ierr); 158c4762a1bSJed Brown 159c4762a1bSJed Brown ii = 4; jj = 4; 160c4762a1bSJed Brown ierr = MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES);CHKERRQ(ierr); 161c4762a1bSJed Brown 162c4762a1bSJed Brown ii = 9; jj = 9; 163c4762a1bSJed Brown ierr = MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES);CHKERRQ(ierr); 164c4762a1bSJed Brown 165c4762a1bSJed Brown ierr = PetscFree(vv);CHKERRQ(ierr); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* 171c4762a1bSJed Brown Push non-zero pattern defined from preallocator into A. 172c4762a1bSJed Brown Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE). 173c4762a1bSJed Brown The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically. 174c4762a1bSJed Brown */ 175c4762a1bSJed Brown ierr = MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A);CHKERRQ(ierr); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* 178c4762a1bSJed Brown We no longer require the preallocator object so destroy it. 179c4762a1bSJed Brown */ 180c4762a1bSJed Brown ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 181c4762a1bSJed Brown 182c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 183c4762a1bSJed Brown 184c4762a1bSJed Brown { 185c4762a1bSJed Brown PetscInt ii,jj; 186c4762a1bSJed Brown PetscScalar *vv; 187c4762a1bSJed Brown 188c4762a1bSJed Brown ierr = PetscCalloc1(bs*bs,&vv);CHKERRQ(ierr); 189c4762a1bSJed Brown for (ii=0; ii<bs*bs; ii++) { 190c4762a1bSJed Brown vv[ii] = (PetscReal)(ii+1); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown ii = 0; jj = 9; 194c4762a1bSJed Brown ierr = MatSetValue(A,ii,jj,33.3,INSERT_VALUES);CHKERRQ(ierr); 195c4762a1bSJed Brown 196c4762a1bSJed Brown ii = 0; jj = 0; 197c4762a1bSJed Brown ierr = MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES);CHKERRQ(ierr); 198c4762a1bSJed Brown 199c4762a1bSJed Brown ii = 2; jj = 4; 200c4762a1bSJed Brown ierr = MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES);CHKERRQ(ierr); 201c4762a1bSJed Brown 202c4762a1bSJed Brown ii = 4; jj = 2; 203c4762a1bSJed Brown ierr = MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES);CHKERRQ(ierr); 204c4762a1bSJed Brown 205c4762a1bSJed Brown ii = 4; jj = 4; 206c4762a1bSJed Brown ierr = MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES);CHKERRQ(ierr); 207c4762a1bSJed Brown 208c4762a1bSJed Brown ii = 9; jj = 9; 209c4762a1bSJed Brown ierr = MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES);CHKERRQ(ierr); 210c4762a1bSJed Brown 211c4762a1bSJed Brown ierr = PetscFree(vv);CHKERRQ(ierr); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 215c4762a1bSJed Brown 216c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 217c4762a1bSJed Brown 218c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 219c4762a1bSJed Brown PetscFunctionReturn(0); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222c4762a1bSJed Brown int main(int argc, char **args) 223c4762a1bSJed Brown { 224c4762a1bSJed Brown PetscErrorCode ierr; 225c4762a1bSJed Brown PetscInt testid = 0; 226c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 227c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-test_id",&testid,NULL);CHKERRQ(ierr); 228c4762a1bSJed Brown switch (testid) { 229c4762a1bSJed Brown case 0: 230c4762a1bSJed Brown ierr = ex1_nonsquare_bs1();CHKERRQ(ierr); 231c4762a1bSJed Brown break; 232c4762a1bSJed Brown case 1: 233c4762a1bSJed Brown ierr = ex2_square_bsvariable();CHKERRQ(ierr); 234c4762a1bSJed Brown break; 235c4762a1bSJed Brown default: 236c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Invalid value for -test_id. Must be {0,1}"); 237c4762a1bSJed Brown break; 238c4762a1bSJed Brown } 239c4762a1bSJed Brown ierr = PetscFinalize(); 240c4762a1bSJed Brown return ierr; 241c4762a1bSJed Brown } 242c4762a1bSJed Brown 243c4762a1bSJed Brown /*TEST 244c4762a1bSJed Brown 245c4762a1bSJed Brown test: 246c4762a1bSJed Brown suffix: t0_a_aij 247c4762a1bSJed Brown nsize: 1 248c4762a1bSJed Brown args: -test_id 0 -mat_type aij 249c4762a1bSJed Brown 250c4762a1bSJed Brown test: 251c4762a1bSJed Brown suffix: t0_b_aij 252c4762a1bSJed Brown nsize: 6 253c4762a1bSJed Brown args: -test_id 0 -mat_type aij 254c4762a1bSJed Brown 255c4762a1bSJed Brown test: 256c4762a1bSJed Brown suffix: t1_a_aij 257c4762a1bSJed Brown nsize: 1 258c4762a1bSJed Brown args: -test_id 1 -mat_type aij 259c4762a1bSJed Brown 260c4762a1bSJed Brown test: 261c4762a1bSJed Brown suffix: t2_a_baij 262c4762a1bSJed Brown nsize: 1 263c4762a1bSJed Brown args: -test_id 1 -mat_type baij 264c4762a1bSJed Brown 265c4762a1bSJed Brown test: 266c4762a1bSJed Brown suffix: t3_a_sbaij 267c4762a1bSJed Brown nsize: 1 268c4762a1bSJed Brown args: -test_id 1 -mat_type sbaij 269c4762a1bSJed Brown 270c4762a1bSJed Brown test: 271c4762a1bSJed Brown suffix: t4_a_aij_bs3 272c4762a1bSJed Brown nsize: 1 273c4762a1bSJed Brown args: -test_id 1 -mat_type aij -block_size 3 274c4762a1bSJed Brown 275c4762a1bSJed Brown test: 276c4762a1bSJed Brown suffix: t5_a_baij_bs3 277c4762a1bSJed Brown nsize: 1 278c4762a1bSJed Brown args: -test_id 1 -mat_type baij -block_size 3 279c4762a1bSJed Brown 280c4762a1bSJed Brown test: 281c4762a1bSJed Brown suffix: t6_a_sbaij_bs3 282c4762a1bSJed Brown nsize: 1 283c4762a1bSJed Brown args: -test_id 1 -mat_type sbaij -block_size 3 284c4762a1bSJed Brown 285c4762a1bSJed Brown test: 286c4762a1bSJed Brown suffix: t4_b_aij_bs3 287c4762a1bSJed Brown nsize: 6 288c4762a1bSJed Brown args: -test_id 1 -mat_type aij -block_size 3 289c4762a1bSJed Brown 290c4762a1bSJed Brown test: 291c4762a1bSJed Brown suffix: t5_b_baij_bs3 292c4762a1bSJed Brown nsize: 6 293c4762a1bSJed Brown args: -test_id 1 -mat_type baij -block_size 3 294c4762a1bSJed Brown filter: grep -v Mat_ 295c4762a1bSJed Brown 296c4762a1bSJed Brown test: 297c4762a1bSJed Brown suffix: t6_b_sbaij_bs3 298c4762a1bSJed Brown nsize: 6 299c4762a1bSJed Brown args: -test_id 1 -mat_type sbaij -block_size 3 300c4762a1bSJed Brown filter: grep -v Mat_ 301c4762a1bSJed Brown 302c4762a1bSJed Brown TEST*/ 303