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 15c4762a1bSJed Brown /* 16c4762a1bSJed Brown Create the Jacobian matrix 17c4762a1bSJed Brown */ 18362febeeSStefano Zampini PetscFunctionBegin; 19c4762a1bSJed Brown M = 10; 20c4762a1bSJed Brown N = 8; 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(A,1)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* 28c4762a1bSJed Brown Get the sizes of the jacobian. 29c4762a1bSJed Brown */ 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(A,&bs)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* 34c4762a1bSJed Brown Create a preallocator matrix with sizes (local and global) matching the jacobian A. 35c4762a1bSJed Brown */ 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),&preallocator)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(preallocator,MATPREALLOCATOR)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(preallocator,m,n,M,N)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(preallocator,bs)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(preallocator)); 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* 43c4762a1bSJed Brown Insert non-zero pattern (e.g. perform a sweep over the grid). 44c4762a1bSJed Brown You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue(). 45c4762a1bSJed Brown */ 46c4762a1bSJed Brown { 47c4762a1bSJed Brown PetscInt ii,jj; 48c4762a1bSJed Brown PetscScalar vv = 0.0; 49c4762a1bSJed Brown 50c4762a1bSJed Brown ii = 3; jj = 3; 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(preallocator,1,&ii,1,&jj,&vv,INSERT_VALUES)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown ii = 7; jj = 4; 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(preallocator,1,&ii,1,&jj,&vv,INSERT_VALUES)); 55c4762a1bSJed Brown 56c4762a1bSJed Brown ii = 9; jj = 7; 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(preallocator,ii,jj,vv,INSERT_VALUES)); 58c4762a1bSJed Brown } 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY)); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* 63c4762a1bSJed Brown Push the non-zero pattern defined within preallocator into A. 64c4762a1bSJed Brown Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE). 65c4762a1bSJed Brown The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically. 66c4762a1bSJed Brown */ 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* 70c4762a1bSJed Brown We no longer require the preallocator object so destroy it. 71c4762a1bSJed Brown */ 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&preallocator)); 73c4762a1bSJed Brown 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* 77c4762a1bSJed Brown Insert non-zero values into A. 78c4762a1bSJed Brown */ 79c4762a1bSJed Brown { 80c4762a1bSJed Brown PetscInt ii,jj; 81c4762a1bSJed Brown PetscScalar vv; 82c4762a1bSJed Brown 83c4762a1bSJed Brown ii = 3; jj = 3; vv = 0.3; 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,ii,jj,vv,INSERT_VALUES)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown ii = 7; jj = 4; vv = 3.3; 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&ii,1,&jj,&vv,INSERT_VALUES)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown ii = 9; jj = 7; vv = 4.3; 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&ii,1,&jj,&vv,INSERT_VALUES)); 91c4762a1bSJed Brown } 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 94c4762a1bSJed Brown 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 96c4762a1bSJed Brown 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 98c4762a1bSJed Brown PetscFunctionReturn(0); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown PetscErrorCode ex2_square_bsvariable(void) 102c4762a1bSJed Brown { 103c4762a1bSJed Brown Mat A,preallocator; 104c4762a1bSJed Brown PetscInt M,N,m,n,bs = 1; 105c4762a1bSJed Brown 106362febeeSStefano Zampini PetscFunctionBegin; 1075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-block_size",&bs,NULL)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* 110c4762a1bSJed Brown Create the Jacobian matrix. 111c4762a1bSJed Brown */ 112c4762a1bSJed Brown M = 10 * bs; 113c4762a1bSJed Brown N = 10 * bs; 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(A,bs)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* 121c4762a1bSJed Brown Get the sizes of the jacobian. 122c4762a1bSJed Brown */ 1235f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(A,&bs)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* 127c4762a1bSJed Brown Create a preallocator matrix with dimensions matching the jacobian A. 128c4762a1bSJed Brown */ 1295f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),&preallocator)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(preallocator,MATPREALLOCATOR)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(preallocator,m,n,M,N)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(preallocator,bs)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(preallocator)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* 136c4762a1bSJed Brown Insert non-zero pattern (e.g. perform a sweep over the grid). 137c4762a1bSJed Brown You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue(). 138c4762a1bSJed Brown */ 139c4762a1bSJed Brown { 140c4762a1bSJed Brown PetscInt ii,jj; 141c4762a1bSJed Brown PetscScalar *vv; 142c4762a1bSJed Brown 1435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(bs*bs,&vv)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown ii = 0; jj = 9; 1465f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(preallocator,ii,jj,vv[0],INSERT_VALUES)); 147c4762a1bSJed Brown 148c4762a1bSJed Brown ii = 0; jj = 0; 1495f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown ii = 2; jj = 4; 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown ii = 4; jj = 2; 1555f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES)); 156c4762a1bSJed Brown 157c4762a1bSJed Brown ii = 4; jj = 4; 1585f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES)); 159c4762a1bSJed Brown 160c4762a1bSJed Brown ii = 9; jj = 9; 1615f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES)); 162c4762a1bSJed Brown 1635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(vv)); 164c4762a1bSJed Brown } 1655f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY)); 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* 169c4762a1bSJed Brown Push non-zero pattern defined from preallocator into A. 170c4762a1bSJed Brown Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE). 171c4762a1bSJed Brown The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically. 172c4762a1bSJed Brown */ 1735f80ce2aSJacob Faibussowitsch CHKERRQ(MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A)); 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* 176c4762a1bSJed Brown We no longer require the preallocator object so destroy it. 177c4762a1bSJed Brown */ 1785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&preallocator)); 179c4762a1bSJed Brown 1805f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 181c4762a1bSJed Brown 182c4762a1bSJed Brown { 183c4762a1bSJed Brown PetscInt ii,jj; 184c4762a1bSJed Brown PetscScalar *vv; 185c4762a1bSJed Brown 1865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(bs*bs,&vv)); 187c4762a1bSJed Brown for (ii=0; ii<bs*bs; ii++) { 188c4762a1bSJed Brown vv[ii] = (PetscReal)(ii+1); 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown ii = 0; jj = 9; 1925f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,ii,jj,33.3,INSERT_VALUES)); 193c4762a1bSJed Brown 194c4762a1bSJed Brown ii = 0; jj = 0; 1955f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES)); 196c4762a1bSJed Brown 197c4762a1bSJed Brown ii = 2; jj = 4; 1985f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES)); 199c4762a1bSJed Brown 200c4762a1bSJed Brown ii = 4; jj = 2; 2015f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown ii = 4; jj = 4; 2045f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES)); 205c4762a1bSJed Brown 206c4762a1bSJed Brown ii = 9; jj = 9; 2075f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES)); 208c4762a1bSJed Brown 2095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(vv)); 210c4762a1bSJed Brown } 2115f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 213c4762a1bSJed Brown 2145f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 215c4762a1bSJed Brown 2165f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 217c4762a1bSJed Brown PetscFunctionReturn(0); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220c4762a1bSJed Brown int main(int argc, char **args) 221c4762a1bSJed Brown { 222c4762a1bSJed Brown PetscInt testid = 0; 223*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-test_id",&testid,NULL)); 225c4762a1bSJed Brown switch (testid) { 226c4762a1bSJed Brown case 0: 2275f80ce2aSJacob Faibussowitsch CHKERRQ(ex1_nonsquare_bs1()); 228c4762a1bSJed Brown break; 229c4762a1bSJed Brown case 1: 2305f80ce2aSJacob Faibussowitsch CHKERRQ(ex2_square_bsvariable()); 231c4762a1bSJed Brown break; 232c4762a1bSJed Brown default: 233c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Invalid value for -test_id. Must be {0,1}"); 234c4762a1bSJed Brown } 235*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 236*b122ec5aSJacob Faibussowitsch return 0; 237c4762a1bSJed Brown } 238c4762a1bSJed Brown 239c4762a1bSJed Brown /*TEST 240c4762a1bSJed Brown 241c4762a1bSJed Brown test: 242c4762a1bSJed Brown suffix: t0_a_aij 243c4762a1bSJed Brown nsize: 1 244c4762a1bSJed Brown args: -test_id 0 -mat_type aij 245c4762a1bSJed Brown 246c4762a1bSJed Brown test: 247c4762a1bSJed Brown suffix: t0_b_aij 248c4762a1bSJed Brown nsize: 6 249c4762a1bSJed Brown args: -test_id 0 -mat_type aij 250c4762a1bSJed Brown 251c4762a1bSJed Brown test: 252c4762a1bSJed Brown suffix: t1_a_aij 253c4762a1bSJed Brown nsize: 1 254c4762a1bSJed Brown args: -test_id 1 -mat_type aij 255c4762a1bSJed Brown 256c4762a1bSJed Brown test: 257c4762a1bSJed Brown suffix: t2_a_baij 258c4762a1bSJed Brown nsize: 1 259c4762a1bSJed Brown args: -test_id 1 -mat_type baij 260c4762a1bSJed Brown 261c4762a1bSJed Brown test: 262c4762a1bSJed Brown suffix: t3_a_sbaij 263c4762a1bSJed Brown nsize: 1 264c4762a1bSJed Brown args: -test_id 1 -mat_type sbaij 265c4762a1bSJed Brown 266c4762a1bSJed Brown test: 267c4762a1bSJed Brown suffix: t4_a_aij_bs3 268c4762a1bSJed Brown nsize: 1 269c4762a1bSJed Brown args: -test_id 1 -mat_type aij -block_size 3 270c4762a1bSJed Brown 271c4762a1bSJed Brown test: 272c4762a1bSJed Brown suffix: t5_a_baij_bs3 273c4762a1bSJed Brown nsize: 1 274c4762a1bSJed Brown args: -test_id 1 -mat_type baij -block_size 3 275c4762a1bSJed Brown 276c4762a1bSJed Brown test: 277c4762a1bSJed Brown suffix: t6_a_sbaij_bs3 278c4762a1bSJed Brown nsize: 1 279c4762a1bSJed Brown args: -test_id 1 -mat_type sbaij -block_size 3 280c4762a1bSJed Brown 281c4762a1bSJed Brown test: 282c4762a1bSJed Brown suffix: t4_b_aij_bs3 283c4762a1bSJed Brown nsize: 6 284c4762a1bSJed Brown args: -test_id 1 -mat_type aij -block_size 3 285c4762a1bSJed Brown 286c4762a1bSJed Brown test: 287c4762a1bSJed Brown suffix: t5_b_baij_bs3 288c4762a1bSJed Brown nsize: 6 289c4762a1bSJed Brown args: -test_id 1 -mat_type baij -block_size 3 290c4762a1bSJed Brown filter: grep -v Mat_ 291c4762a1bSJed Brown 292c4762a1bSJed Brown test: 293c4762a1bSJed Brown suffix: t6_b_sbaij_bs3 294c4762a1bSJed Brown nsize: 6 295c4762a1bSJed Brown args: -test_id 1 -mat_type sbaij -block_size 3 296c4762a1bSJed Brown filter: grep -v Mat_ 297c4762a1bSJed Brown 298c4762a1bSJed Brown TEST*/ 299