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