1 /* 2 Defines the basic matrix operations for the BAIJ (compressed row) 3 matrix storage format. 4 */ 5 #include "src/mat/impls/baij/seq/baij.h" 6 #include "src/inline/spops.h" 7 #include "petscsys.h" /*I "petscmat.h" I*/ 8 9 #include "src/inline/ilu.h" 10 11 #undef __FUNCT__ 12 #define __FUNCT__ "MatInvertBlockDiagonal_SeqBAIJ" 13 PetscErrorCode MatInvertBlockDiagonal_SeqBAIJ(Mat A) 14 { 15 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*) A->data; 16 PetscErrorCode ierr; 17 PetscInt *diag_offset,i,bs = A->bs,mbs = a->mbs; 18 PetscScalar *v = a->a,*odiag,*diag,*mdiag; 19 20 PetscFunctionBegin; 21 if (a->idiagvalid) PetscFunctionReturn(0); 22 ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 23 diag_offset = a->diag; 24 if (!a->idiag) { 25 ierr = PetscMalloc(2*bs*bs*mbs*sizeof(PetscScalar),&a->idiag);CHKERRQ(ierr); 26 } 27 diag = a->idiag; 28 mdiag = a->idiag+bs*bs*mbs; 29 /* factor and invert each block */ 30 switch (bs){ 31 case 2: 32 for (i=0; i<mbs; i++) { 33 odiag = v + 4*diag_offset[i]; 34 diag[0] = odiag[0]; diag[1] = odiag[1]; diag[2] = odiag[2]; diag[3] = odiag[3]; 35 mdiag[0] = odiag[0]; mdiag[1] = odiag[1]; mdiag[2] = odiag[2]; mdiag[3] = odiag[3]; 36 ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr); 37 diag += 4; 38 mdiag += 4; 39 } 40 break; 41 case 3: 42 for (i=0; i<mbs; i++) { 43 odiag = v + 9*diag_offset[i]; 44 diag[0] = odiag[0]; diag[1] = odiag[1]; diag[2] = odiag[2]; diag[3] = odiag[3]; 45 diag[4] = odiag[4]; diag[5] = odiag[5]; diag[6] = odiag[6]; diag[7] = odiag[7]; 46 diag[8] = odiag[8]; 47 mdiag[0] = odiag[0]; mdiag[1] = odiag[1]; mdiag[2] = odiag[2]; mdiag[3] = odiag[3]; 48 mdiag[4] = odiag[4]; mdiag[5] = odiag[5]; mdiag[6] = odiag[6]; mdiag[7] = odiag[7]; 49 mdiag[8] = odiag[8]; 50 ierr = Kernel_A_gets_inverse_A_3(diag);CHKERRQ(ierr); 51 diag += 9; 52 mdiag += 9; 53 } 54 break; 55 case 4: 56 for (i=0; i<mbs; i++) { 57 odiag = v + 16*diag_offset[i]; 58 ierr = PetscMemcpy(diag,odiag,16*sizeof(PetscScalar));CHKERRQ(ierr); 59 ierr = PetscMemcpy(mdiag,odiag,16*sizeof(PetscScalar));CHKERRQ(ierr); 60 ierr = Kernel_A_gets_inverse_A_4(diag);CHKERRQ(ierr); 61 diag += 16; 62 mdiag += 16; 63 } 64 break; 65 case 5: 66 for (i=0; i<mbs; i++) { 67 odiag = v + 25*diag_offset[i]; 68 ierr = PetscMemcpy(diag,odiag,25*sizeof(PetscScalar));CHKERRQ(ierr); 69 ierr = PetscMemcpy(mdiag,odiag,25*sizeof(PetscScalar));CHKERRQ(ierr); 70 ierr = Kernel_A_gets_inverse_A_5(diag);CHKERRQ(ierr); 71 diag += 25; 72 mdiag += 25; 73 } 74 break; 75 default: 76 SETERRQ1(PETSC_ERR_SUP,"not supported for block size %D",bs); 77 } 78 a->idiagvalid = PETSC_TRUE; 79 PetscFunctionReturn(0); 80 } 81 82 #undef __FUNCT__ 83 #define __FUNCT__ "MatPBRelax_SeqBAIJ_2" 84 PetscErrorCode MatPBRelax_SeqBAIJ_2(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 85 { 86 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 87 PetscScalar *x,x1,x2,s1,s2; 88 const PetscScalar *v,*aa = a->a, *b, *idiag,*mdiag; 89 PetscErrorCode ierr; 90 PetscInt m = a->mbs,i,i2,nz,idx; 91 const PetscInt *diag,*ai = a->i,*aj = a->j,*vi; 92 93 PetscFunctionBegin; 94 its = its*lits; 95 if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 96 if (fshift) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for diagonal shift"); 97 if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for non-trivial relaxation factor"); 98 if ((flag & SOR_EISENSTAT) ||(flag & SOR_APPLY_UPPER) || (flag & SOR_APPLY_LOWER) ) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for Eisenstat trick"); 99 if (its > 1) SETERRQ(PETSC_ERR_SUP,"Sorry, no support yet for multiple point block SOR iterations"); 100 101 if (!a->idiagvalid){ierr = MatInvertBlockDiagonal_SeqBAIJ(A);CHKERRQ(ierr);} 102 103 diag = a->diag; 104 idiag = a->idiag; 105 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 106 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 107 108 if (flag & SOR_ZERO_INITIAL_GUESS) { 109 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 110 x[0] = b[0]*idiag[0] + b[1]*idiag[2]; 111 x[1] = b[0]*idiag[1] + b[1]*idiag[3]; 112 i2 = 2; 113 idiag += 4; 114 for (i=1; i<m; i++) { 115 v = aa + 4*ai[i]; 116 vi = aj + ai[i]; 117 nz = diag[i] - ai[i]; 118 s1 = b[i2]; s2 = b[i2+1]; 119 while (nz--) { 120 idx = 2*(*vi++); 121 x1 = x[idx]; x2 = x[1+idx]; 122 s1 -= v[0]*x1 + v[2]*x2; 123 s2 -= v[1]*x1 + v[3]*x2; 124 v += 4; 125 } 126 x[i2] = idiag[0]*s1 + idiag[2]*s2; 127 x[i2+1] = idiag[1]*s1 + idiag[3]*s2; 128 idiag += 4; 129 i2 += 2; 130 } 131 /* for logging purposes assume number of nonzero in lower half is 1/2 of total */ 132 PetscLogFlops(4*(a->nz)); 133 } 134 if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 135 (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 136 i2 = 0; 137 mdiag = a->idiag+4*a->mbs; 138 for (i=0; i<m; i++) { 139 x1 = x[i2]; x2 = x[i2+1]; 140 x[i2] = mdiag[0]*x1 + mdiag[2]*x2; 141 x[i2+1] = mdiag[1]*x1 + mdiag[3]*x2; 142 mdiag += 4; 143 i2 += 2; 144 } 145 PetscLogFlops(6*m); 146 } else if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 147 ierr = PetscMemcpy(x,b,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 148 } 149 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 150 idiag = a->idiag+4*a->mbs - 4; 151 i2 = 2*m - 2; 152 x1 = x[i2]; x2 = x[i2+1]; 153 x[i2] = idiag[0]*x1 + idiag[2]*x2; 154 x[i2+1] = idiag[1]*x1 + idiag[3]*x2; 155 idiag -= 4; 156 i2 -= 2; 157 for (i=m-2; i>=0; i--) { 158 v = aa + 4*(diag[i]+1); 159 vi = aj + diag[i] + 1; 160 nz = ai[i+1] - diag[i] - 1; 161 s1 = x[i2]; s2 = x[i2+1]; 162 while (nz--) { 163 idx = 2*(*vi++); 164 x1 = x[idx]; x2 = x[1+idx]; 165 s1 -= v[0]*x1 + v[2]*x2; 166 s2 -= v[1]*x1 + v[3]*x2; 167 v += 4; 168 } 169 x[i2] = idiag[0]*s1 + idiag[2]*s2; 170 x[i2+1] = idiag[1]*s1 + idiag[3]*s2; 171 idiag -= 4; 172 i2 -= 2; 173 } 174 PetscLogFlops(4*(a->nz)); 175 } 176 } else { 177 SETERRQ(PETSC_ERR_SUP,"Only supports point block SOR with zero initial guess"); 178 } 179 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 180 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 181 PetscFunctionReturn(0); 182 } 183 184 #undef __FUNCT__ 185 #define __FUNCT__ "MatPBRelax_SeqBAIJ_3" 186 PetscErrorCode MatPBRelax_SeqBAIJ_3(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 187 { 188 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 189 PetscScalar *x,x1,x2,x3,s1,s2,s3; 190 const PetscScalar *v,*aa = a->a, *b, *idiag,*mdiag; 191 PetscErrorCode ierr; 192 PetscInt m = a->mbs,i,i2,nz,idx; 193 const PetscInt *diag,*ai = a->i,*aj = a->j,*vi; 194 195 PetscFunctionBegin; 196 its = its*lits; 197 if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 198 if (fshift) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for diagonal shift"); 199 if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for non-trivial relaxation factor"); 200 if ((flag & SOR_EISENSTAT) ||(flag & SOR_APPLY_UPPER) || (flag & SOR_APPLY_LOWER) ) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for Eisenstat trick"); 201 if (its > 1) SETERRQ(PETSC_ERR_SUP,"Sorry, no support yet for multiple point block SOR iterations"); 202 203 if (!a->idiagvalid){ierr = MatInvertBlockDiagonal_SeqBAIJ(A);CHKERRQ(ierr);} 204 205 diag = a->diag; 206 idiag = a->idiag; 207 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 208 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 209 210 if (flag & SOR_ZERO_INITIAL_GUESS) { 211 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 212 x[0] = b[0]*idiag[0] + b[1]*idiag[3] + b[2]*idiag[6]; 213 x[1] = b[0]*idiag[1] + b[1]*idiag[4] + b[2]*idiag[7]; 214 x[2] = b[0]*idiag[2] + b[1]*idiag[5] + b[2]*idiag[8]; 215 i2 = 3; 216 idiag += 9; 217 for (i=1; i<m; i++) { 218 v = aa + 9*ai[i]; 219 vi = aj + ai[i]; 220 nz = diag[i] - ai[i]; 221 s1 = b[i2]; s2 = b[i2+1]; s3 = b[i2+2]; 222 while (nz--) { 223 idx = 3*(*vi++); 224 x1 = x[idx]; x2 = x[1+idx];x3 = x[2+idx]; 225 s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3; 226 s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3; 227 s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3; 228 v += 9; 229 } 230 x[i2] = idiag[0]*s1 + idiag[3]*s2 + idiag[6]*s3; 231 x[i2+1] = idiag[1]*s1 + idiag[4]*s2 + idiag[7]*s3; 232 x[i2+2] = idiag[2]*s1 + idiag[5]*s2 + idiag[8]*s3; 233 idiag += 9; 234 i2 += 3; 235 } 236 /* for logging purposes assume number of nonzero in lower half is 1/2 of total */ 237 PetscLogFlops(9*(a->nz)); 238 } 239 if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 240 (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 241 i2 = 0; 242 mdiag = a->idiag+9*a->mbs; 243 for (i=0; i<m; i++) { 244 x1 = x[i2]; x2 = x[i2+1]; x3 = x[i2+2]; 245 x[i2] = mdiag[0]*x1 + mdiag[3]*x2 + mdiag[6]*x3; 246 x[i2+1] = mdiag[1]*x1 + mdiag[4]*x2 + mdiag[7]*x3; 247 x[i2+2] = mdiag[2]*x1 + mdiag[5]*x2 + mdiag[8]*x3; 248 mdiag += 9; 249 i2 += 3; 250 } 251 PetscLogFlops(15*m); 252 } else if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 253 ierr = PetscMemcpy(x,b,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 254 } 255 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 256 idiag = a->idiag+9*a->mbs - 9; 257 i2 = 3*m - 3; 258 x1 = x[i2]; x2 = x[i2+1]; x3 = x[i2+2]; 259 x[i2] = idiag[0]*x1 + idiag[3]*x2 + idiag[6]*x3; 260 x[i2+1] = idiag[1]*x1 + idiag[4]*x2 + idiag[7]*x3; 261 x[i2+2] = idiag[2]*x1 + idiag[5]*x2 + idiag[8]*x3; 262 idiag -= 9; 263 i2 -= 3; 264 for (i=m-2; i>=0; i--) { 265 v = aa + 9*(diag[i]+1); 266 vi = aj + diag[i] + 1; 267 nz = ai[i+1] - diag[i] - 1; 268 s1 = x[i2]; s2 = x[i2+1]; s3 = x[i2+2]; 269 while (nz--) { 270 idx = 3*(*vi++); 271 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; 272 s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3; 273 s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3; 274 s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3; 275 v += 9; 276 } 277 x[i2] = idiag[0]*s1 + idiag[3]*s2 + idiag[6]*s3; 278 x[i2+1] = idiag[1]*s1 + idiag[4]*s2 + idiag[7]*s3; 279 x[i2+2] = idiag[2]*s1 + idiag[5]*s2 + idiag[8]*s3; 280 idiag -= 9; 281 i2 -= 3; 282 } 283 PetscLogFlops(9*(a->nz)); 284 } 285 } else { 286 SETERRQ(PETSC_ERR_SUP,"Only supports point block SOR with zero initial guess"); 287 } 288 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 289 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 290 PetscFunctionReturn(0); 291 } 292 293 #undef __FUNCT__ 294 #define __FUNCT__ "MatPBRelax_SeqBAIJ_4" 295 PetscErrorCode MatPBRelax_SeqBAIJ_4(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 296 { 297 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 298 PetscScalar *x,x1,x2,x3,x4,s1,s2,s3,s4; 299 const PetscScalar *v,*aa = a->a, *b, *idiag,*mdiag; 300 PetscErrorCode ierr; 301 PetscInt m = a->mbs,i,i2,nz,idx; 302 const PetscInt *diag,*ai = a->i,*aj = a->j,*vi; 303 304 PetscFunctionBegin; 305 its = its*lits; 306 if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 307 if (fshift) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for diagonal shift"); 308 if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for non-trivial relaxation factor"); 309 if ((flag & SOR_EISENSTAT) ||(flag & SOR_APPLY_UPPER) || (flag & SOR_APPLY_LOWER) ) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for Eisenstat trick"); 310 if (its > 1) SETERRQ(PETSC_ERR_SUP,"Sorry, no support yet for multiple point block SOR iterations"); 311 312 if (!a->idiagvalid){ierr = MatInvertBlockDiagonal_SeqBAIJ(A);CHKERRQ(ierr);} 313 314 diag = a->diag; 315 idiag = a->idiag; 316 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 317 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 318 319 if (flag & SOR_ZERO_INITIAL_GUESS) { 320 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 321 x[0] = b[0]*idiag[0] + b[1]*idiag[4] + b[2]*idiag[8] + b[3]*idiag[12]; 322 x[1] = b[0]*idiag[1] + b[1]*idiag[5] + b[2]*idiag[9] + b[3]*idiag[13]; 323 x[2] = b[0]*idiag[2] + b[1]*idiag[6] + b[2]*idiag[10] + b[3]*idiag[14]; 324 x[3] = b[0]*idiag[3] + b[1]*idiag[7] + b[2]*idiag[11] + b[3]*idiag[15]; 325 i2 = 4; 326 idiag += 16; 327 for (i=1; i<m; i++) { 328 v = aa + 16*ai[i]; 329 vi = aj + ai[i]; 330 nz = diag[i] - ai[i]; 331 s1 = b[i2]; s2 = b[i2+1]; s3 = b[i2+2]; s4 = b[i2+3]; 332 while (nz--) { 333 idx = 4*(*vi++); 334 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; x4 = x[3+idx]; 335 s1 -= v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 336 s2 -= v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 337 s3 -= v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 338 s4 -= v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 339 v += 16; 340 } 341 x[i2] = idiag[0]*s1 + idiag[4]*s2 + idiag[8]*s3 + idiag[12]*s4; 342 x[i2+1] = idiag[1]*s1 + idiag[5]*s2 + idiag[9]*s3 + idiag[13]*s4; 343 x[i2+2] = idiag[2]*s1 + idiag[6]*s2 + idiag[10]*s3 + idiag[14]*s4; 344 x[i2+3] = idiag[3]*s1 + idiag[7]*s2 + idiag[11]*s3 + idiag[15]*s4; 345 idiag += 16; 346 i2 += 4; 347 } 348 /* for logging purposes assume number of nonzero in lower half is 1/2 of total */ 349 PetscLogFlops(16*(a->nz)); 350 } 351 if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 352 (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 353 i2 = 0; 354 mdiag = a->idiag+16*a->mbs; 355 for (i=0; i<m; i++) { 356 x1 = x[i2]; x2 = x[i2+1]; x3 = x[i2+2]; x4 = x[i2+3]; 357 x[i2] = mdiag[0]*x1 + mdiag[4]*x2 + mdiag[8]*x3 + mdiag[12]*x4; 358 x[i2+1] = mdiag[1]*x1 + mdiag[5]*x2 + mdiag[9]*x3 + mdiag[13]*x4; 359 x[i2+2] = mdiag[2]*x1 + mdiag[6]*x2 + mdiag[10]*x3 + mdiag[14]*x4; 360 x[i2+3] = mdiag[3]*x1 + mdiag[7]*x2 + mdiag[11]*x3 + mdiag[15]*x4; 361 mdiag += 16; 362 i2 += 4; 363 } 364 PetscLogFlops(28*m); 365 } else if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 366 ierr = PetscMemcpy(x,b,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 367 } 368 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 369 idiag = a->idiag+16*a->mbs - 16; 370 i2 = 4*m - 4; 371 x1 = x[i2]; x2 = x[i2+1]; x3 = x[i2+2]; x4 = x[i2+3]; 372 x[i2] = idiag[0]*x1 + idiag[4]*x2 + idiag[8]*x3 + idiag[12]*x4; 373 x[i2+1] = idiag[1]*x1 + idiag[5]*x2 + idiag[9]*x3 + idiag[13]*x4; 374 x[i2+2] = idiag[2]*x1 + idiag[6]*x2 + idiag[10]*x3 + idiag[14]*x4; 375 x[i2+3] = idiag[3]*x1 + idiag[7]*x2 + idiag[11]*x3 + idiag[15]*x4; 376 idiag -= 16; 377 i2 -= 4; 378 for (i=m-2; i>=0; i--) { 379 v = aa + 16*(diag[i]+1); 380 vi = aj + diag[i] + 1; 381 nz = ai[i+1] - diag[i] - 1; 382 s1 = x[i2]; s2 = x[i2+1]; s3 = x[i2+2]; s4 = x[i2+3]; 383 while (nz--) { 384 idx = 4*(*vi++); 385 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; x4 = x[3+idx]; 386 s1 -= v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 387 s2 -= v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 388 s3 -= v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 389 s4 -= v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 390 v += 16; 391 } 392 x[i2] = idiag[0]*s1 + idiag[4]*s2 + idiag[8]*s3 + idiag[12]*s4; 393 x[i2+1] = idiag[1]*s1 + idiag[5]*s2 + idiag[9]*s3 + idiag[13]*s4; 394 x[i2+2] = idiag[2]*s1 + idiag[6]*s2 + idiag[10]*s3 + idiag[14]*s4; 395 x[i2+3] = idiag[3]*s1 + idiag[7]*s2 + idiag[11]*s3 + idiag[15]*s4; 396 idiag -= 16; 397 i2 -= 4; 398 } 399 PetscLogFlops(16*(a->nz)); 400 } 401 } else { 402 SETERRQ(PETSC_ERR_SUP,"Only supports point block SOR with zero initial guess"); 403 } 404 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 405 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 406 PetscFunctionReturn(0); 407 } 408 409 #undef __FUNCT__ 410 #define __FUNCT__ "MatPBRelax_SeqBAIJ_5" 411 PetscErrorCode MatPBRelax_SeqBAIJ_5(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 412 { 413 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 414 PetscScalar *x,x1,x2,x3,x4,x5,s1,s2,s3,s4,s5; 415 const PetscScalar *v,*aa = a->a, *b, *idiag,*mdiag; 416 PetscErrorCode ierr; 417 PetscInt m = a->mbs,i,i2,nz,idx; 418 const PetscInt *diag,*ai = a->i,*aj = a->j,*vi; 419 420 PetscFunctionBegin; 421 its = its*lits; 422 if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 423 if (fshift) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for diagonal shift"); 424 if (omega != 1.0) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for non-trivial relaxation factor"); 425 if ((flag & SOR_EISENSTAT) ||(flag & SOR_APPLY_UPPER) || (flag & SOR_APPLY_LOWER) ) SETERRQ(PETSC_ERR_SUP,"Sorry, no support for Eisenstat trick"); 426 if (its > 1) SETERRQ(PETSC_ERR_SUP,"Sorry, no support yet for multiple point block SOR iterations"); 427 428 if (!a->idiagvalid){ierr = MatInvertBlockDiagonal_SeqBAIJ(A);CHKERRQ(ierr);} 429 430 diag = a->diag; 431 idiag = a->idiag; 432 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 433 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 434 435 if (flag & SOR_ZERO_INITIAL_GUESS) { 436 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 437 x[0] = b[0]*idiag[0] + b[1]*idiag[5] + b[2]*idiag[10] + b[3]*idiag[15] + b[4]*idiag[20]; 438 x[1] = b[0]*idiag[1] + b[1]*idiag[6] + b[2]*idiag[11] + b[3]*idiag[16] + b[4]*idiag[21]; 439 x[2] = b[0]*idiag[2] + b[1]*idiag[7] + b[2]*idiag[12] + b[3]*idiag[17] + b[4]*idiag[22]; 440 x[3] = b[0]*idiag[3] + b[1]*idiag[8] + b[2]*idiag[13] + b[3]*idiag[18] + b[4]*idiag[23]; 441 x[4] = b[0]*idiag[4] + b[1]*idiag[9] + b[2]*idiag[14] + b[3]*idiag[19] + b[4]*idiag[24]; 442 i2 = 5; 443 idiag += 25; 444 for (i=1; i<m; i++) { 445 v = aa + 25*ai[i]; 446 vi = aj + ai[i]; 447 nz = diag[i] - ai[i]; 448 s1 = b[i2]; s2 = b[i2+1]; s3 = b[i2+2]; s4 = b[i2+3]; s5 = b[i2+4]; 449 while (nz--) { 450 idx = 5*(*vi++); 451 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; x4 = x[3+idx]; x5 = x[4+idx]; 452 s1 -= v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 453 s2 -= v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 454 s3 -= v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 455 s4 -= v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 456 s5 -= v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 457 v += 25; 458 } 459 x[i2] = idiag[0]*s1 + idiag[5]*s2 + idiag[10]*s3 + idiag[15]*s4 + idiag[20]*s5; 460 x[i2+1] = idiag[1]*s1 + idiag[6]*s2 + idiag[11]*s3 + idiag[16]*s4 + idiag[21]*s5; 461 x[i2+2] = idiag[2]*s1 + idiag[7]*s2 + idiag[12]*s3 + idiag[17]*s4 + idiag[22]*s5; 462 x[i2+3] = idiag[3]*s1 + idiag[8]*s2 + idiag[13]*s3 + idiag[18]*s4 + idiag[23]*s5; 463 x[i2+4] = idiag[4]*s1 + idiag[9]*s2 + idiag[14]*s3 + idiag[19]*s4 + idiag[24]*s5; 464 idiag += 25; 465 i2 += 5; 466 } 467 /* for logging purposes assume number of nonzero in lower half is 1/2 of total */ 468 PetscLogFlops(25*(a->nz)); 469 } 470 if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 471 (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 472 i2 = 0; 473 mdiag = a->idiag+25*a->mbs; 474 for (i=0; i<m; i++) { 475 x1 = x[i2]; x2 = x[i2+1]; x3 = x[i2+2]; x4 = x[i2+3]; x5 = x[i2+4]; 476 x[i2] = mdiag[0]*x1 + mdiag[5]*x2 + mdiag[10]*x3 + mdiag[15]*x4 + mdiag[20]*x5; 477 x[i2+1] = mdiag[1]*x1 + mdiag[6]*x2 + mdiag[11]*x3 + mdiag[16]*x4 + mdiag[21]*x5; 478 x[i2+2] = mdiag[2]*x1 + mdiag[7]*x2 + mdiag[12]*x3 + mdiag[17]*x4 + mdiag[22]*x5; 479 x[i2+3] = mdiag[3]*x1 + mdiag[8]*x2 + mdiag[13]*x3 + mdiag[18]*x4 + mdiag[23]*x5; 480 x[i2+4] = mdiag[4]*x1 + mdiag[9]*x2 + mdiag[14]*x3 + mdiag[19]*x4 + mdiag[24]*x5; 481 mdiag += 25; 482 i2 += 5; 483 } 484 PetscLogFlops(45*m); 485 } else if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 486 ierr = PetscMemcpy(x,b,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 487 } 488 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 489 idiag = a->idiag+25*a->mbs - 25; 490 i2 = 5*m - 5; 491 x1 = x[i2]; x2 = x[i2+1]; x3 = x[i2+2]; x4 = x[i2+3]; x5 = x[i2+4]; 492 x[i2] = idiag[0]*x1 + idiag[5]*x2 + idiag[10]*x3 + idiag[15]*x4 + idiag[20]*x5; 493 x[i2+1] = idiag[1]*x1 + idiag[6]*x2 + idiag[11]*x3 + idiag[16]*x4 + idiag[21]*x5; 494 x[i2+2] = idiag[2]*x1 + idiag[7]*x2 + idiag[12]*x3 + idiag[17]*x4 + idiag[22]*x5; 495 x[i2+3] = idiag[3]*x1 + idiag[8]*x2 + idiag[13]*x3 + idiag[18]*x4 + idiag[23]*x5; 496 x[i2+4] = idiag[4]*x1 + idiag[9]*x2 + idiag[14]*x3 + idiag[19]*x4 + idiag[24]*x5; 497 idiag -= 25; 498 i2 -= 5; 499 for (i=m-2; i>=0; i--) { 500 v = aa + 25*(diag[i]+1); 501 vi = aj + diag[i] + 1; 502 nz = ai[i+1] - diag[i] - 1; 503 s1 = x[i2]; s2 = x[i2+1]; s3 = x[i2+2]; s4 = x[i2+3]; s5 = x[i2+4]; 504 while (nz--) { 505 idx = 5*(*vi++); 506 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; x4 = x[3+idx]; x5 = x[4+idx]; 507 s1 -= v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 508 s2 -= v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 509 s3 -= v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 510 s4 -= v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 511 s5 -= v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 512 v += 25; 513 } 514 x[i2] = idiag[0]*s1 + idiag[5]*s2 + idiag[10]*s3 + idiag[15]*s4 + idiag[20]*s5; 515 x[i2+1] = idiag[1]*s1 + idiag[6]*s2 + idiag[11]*s3 + idiag[16]*s4 + idiag[21]*s5; 516 x[i2+2] = idiag[2]*s1 + idiag[7]*s2 + idiag[12]*s3 + idiag[17]*s4 + idiag[22]*s5; 517 x[i2+3] = idiag[3]*s1 + idiag[8]*s2 + idiag[13]*s3 + idiag[18]*s4 + idiag[23]*s5; 518 x[i2+4] = idiag[4]*s1 + idiag[9]*s2 + idiag[14]*s3 + idiag[19]*s4 + idiag[24]*s5; 519 idiag -= 25; 520 i2 -= 5; 521 } 522 PetscLogFlops(25*(a->nz)); 523 } 524 } else { 525 SETERRQ(PETSC_ERR_SUP,"Only supports point block SOR with zero initial guess"); 526 } 527 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 528 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 529 PetscFunctionReturn(0); 530 } 531 532 /* 533 Special version for Fun3d sequential benchmark 534 */ 535 #if defined(PETSC_HAVE_FORTRAN_CAPS) 536 #define matsetvaluesblocked4_ MATSETVALUESBLOCKED4 537 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 538 #define matsetvaluesblocked4_ matsetvaluesblocked4 539 #endif 540 541 EXTERN_C_BEGIN 542 #undef __FUNCT__ 543 #define __FUNCT__ "matsetvaluesblocked4_" 544 void matsetvaluesblocked4_(Mat *AA,PetscInt *mm,const PetscInt im[],PetscInt *nn,const PetscInt in[],const PetscScalar v[]) 545 { 546 Mat A = *AA; 547 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 548 PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,N,m = *mm,n = *nn; 549 PetscInt *ai=a->i,*ailen=a->ilen; 550 PetscInt *aj=a->j,stepval; 551 const PetscScalar *value = v; 552 MatScalar *ap,*aa = a->a,*bap; 553 554 PetscFunctionBegin; 555 stepval = (n-1)*4; 556 for (k=0; k<m; k++) { /* loop over added rows */ 557 row = im[k]; 558 rp = aj + ai[row]; 559 ap = aa + 16*ai[row]; 560 nrow = ailen[row]; 561 low = 0; 562 for (l=0; l<n; l++) { /* loop over added columns */ 563 col = in[l]; 564 value = v + k*(stepval+4)*4 + l*4; 565 low = 0; high = nrow; 566 while (high-low > 7) { 567 t = (low+high)/2; 568 if (rp[t] > col) high = t; 569 else low = t; 570 } 571 for (i=low; i<high; i++) { 572 if (rp[i] > col) break; 573 if (rp[i] == col) { 574 bap = ap + 16*i; 575 for (ii=0; ii<4; ii++,value+=stepval) { 576 for (jj=ii; jj<16; jj+=4) { 577 bap[jj] += *value++; 578 } 579 } 580 goto noinsert2; 581 } 582 } 583 N = nrow++ - 1; 584 /* shift up all the later entries in this row */ 585 for (ii=N; ii>=i; ii--) { 586 rp[ii+1] = rp[ii]; 587 PetscMemcpy(ap+16*(ii+1),ap+16*(ii),16*sizeof(MatScalar)); 588 } 589 if (N >= i) { 590 PetscMemzero(ap+16*i,16*sizeof(MatScalar)); 591 } 592 rp[i] = col; 593 bap = ap + 16*i; 594 for (ii=0; ii<4; ii++,value+=stepval) { 595 for (jj=ii; jj<16; jj+=4) { 596 bap[jj] = *value++; 597 } 598 } 599 noinsert2:; 600 low = i; 601 } 602 ailen[row] = nrow; 603 } 604 } 605 EXTERN_C_END 606 607 #if defined(PETSC_HAVE_FORTRAN_CAPS) 608 #define matsetvalues4_ MATSETVALUES4 609 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 610 #define matsetvalues4_ matsetvalues4 611 #endif 612 613 EXTERN_C_BEGIN 614 #undef __FUNCT__ 615 #define __FUNCT__ "MatSetValues4_" 616 void matsetvalues4_(Mat *AA,PetscInt *mm,PetscInt *im,PetscInt *nn,PetscInt *in,PetscScalar *v) 617 { 618 Mat A = *AA; 619 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 620 PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,N,n = *nn,m = *mm; 621 PetscInt *ai=a->i,*ailen=a->ilen; 622 PetscInt *aj=a->j,brow,bcol; 623 PetscInt ridx,cidx; 624 MatScalar *ap,value,*aa=a->a,*bap; 625 626 PetscFunctionBegin; 627 for (k=0; k<m; k++) { /* loop over added rows */ 628 row = im[k]; brow = row/4; 629 rp = aj + ai[brow]; 630 ap = aa + 16*ai[brow]; 631 nrow = ailen[brow]; 632 low = 0; 633 for (l=0; l<n; l++) { /* loop over added columns */ 634 col = in[l]; bcol = col/4; 635 ridx = row % 4; cidx = col % 4; 636 value = v[l + k*n]; 637 low = 0; high = nrow; 638 while (high-low > 7) { 639 t = (low+high)/2; 640 if (rp[t] > bcol) high = t; 641 else low = t; 642 } 643 for (i=low; i<high; i++) { 644 if (rp[i] > bcol) break; 645 if (rp[i] == bcol) { 646 bap = ap + 16*i + 4*cidx + ridx; 647 *bap += value; 648 goto noinsert1; 649 } 650 } 651 N = nrow++ - 1; 652 /* shift up all the later entries in this row */ 653 for (ii=N; ii>=i; ii--) { 654 rp[ii+1] = rp[ii]; 655 PetscMemcpy(ap+16*(ii+1),ap+16*(ii),16*sizeof(MatScalar)); 656 } 657 if (N>=i) { 658 PetscMemzero(ap+16*i,16*sizeof(MatScalar)); 659 } 660 rp[i] = bcol; 661 ap[16*i + 4*cidx + ridx] = value; 662 noinsert1:; 663 low = i; 664 } 665 ailen[brow] = nrow; 666 } 667 } 668 EXTERN_C_END 669 670 /* UGLY, ugly, ugly 671 When MatScalar == PetscScalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does 672 not exist. Otherwise ..._MatScalar() takes matrix dlements in single precision and 673 inserts them into the single precision data structure. The function MatSetValuesBlocked_SeqBAIJ() 674 converts the entries into single precision and then calls ..._MatScalar() to put them 675 into the single precision data structures. 676 */ 677 #if defined(PETSC_USE_MAT_SINGLE) 678 EXTERN PetscErrorCode MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const MatScalar[],InsertMode); 679 #else 680 #define MatSetValuesBlocked_SeqBAIJ_MatScalar MatSetValuesBlocked_SeqBAIJ 681 #endif 682 683 #define CHUNKSIZE 10 684 685 /* 686 Checks for missing diagonals 687 */ 688 #undef __FUNCT__ 689 #define __FUNCT__ "MatMissingDiagonal_SeqBAIJ" 690 PetscErrorCode MatMissingDiagonal_SeqBAIJ(Mat A) 691 { 692 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 693 PetscErrorCode ierr; 694 PetscInt *diag,*jj = a->j,i; 695 696 PetscFunctionBegin; 697 ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 698 diag = a->diag; 699 for (i=0; i<a->mbs; i++) { 700 if (jj[diag[i]] != i) { 701 SETERRQ1(PETSC_ERR_PLIB,"Matrix is missing diagonal number %D",i); 702 } 703 } 704 PetscFunctionReturn(0); 705 } 706 707 #undef __FUNCT__ 708 #define __FUNCT__ "MatMarkDiagonal_SeqBAIJ" 709 PetscErrorCode MatMarkDiagonal_SeqBAIJ(Mat A) 710 { 711 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 712 PetscErrorCode ierr; 713 PetscInt i,j,*diag,m = a->mbs; 714 715 PetscFunctionBegin; 716 if (a->diag) PetscFunctionReturn(0); 717 718 ierr = PetscMalloc((m+1)*sizeof(PetscInt),&diag);CHKERRQ(ierr); 719 PetscLogObjectMemory(A,(m+1)*sizeof(PetscInt)); 720 for (i=0; i<m; i++) { 721 diag[i] = a->i[i+1]; 722 for (j=a->i[i]; j<a->i[i+1]; j++) { 723 if (a->j[j] == i) { 724 diag[i] = j; 725 break; 726 } 727 } 728 } 729 a->diag = diag; 730 PetscFunctionReturn(0); 731 } 732 733 734 EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 735 736 #undef __FUNCT__ 737 #define __FUNCT__ "MatGetRowIJ_SeqBAIJ" 738 static PetscErrorCode MatGetRowIJ_SeqBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 739 { 740 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 741 PetscErrorCode ierr; 742 PetscInt n = a->mbs,i; 743 744 PetscFunctionBegin; 745 *nn = n; 746 if (!ia) PetscFunctionReturn(0); 747 if (symmetric) { 748 ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 749 } else if (oshift == 1) { 750 /* temporarily add 1 to i and j indices */ 751 PetscInt nz = a->i[n]; 752 for (i=0; i<nz; i++) a->j[i]++; 753 for (i=0; i<n+1; i++) a->i[i]++; 754 *ia = a->i; *ja = a->j; 755 } else { 756 *ia = a->i; *ja = a->j; 757 } 758 759 PetscFunctionReturn(0); 760 } 761 762 #undef __FUNCT__ 763 #define __FUNCT__ "MatRestoreRowIJ_SeqBAIJ" 764 static PetscErrorCode MatRestoreRowIJ_SeqBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 765 { 766 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 767 PetscErrorCode ierr; 768 PetscInt i,n = a->mbs; 769 770 PetscFunctionBegin; 771 if (!ia) PetscFunctionReturn(0); 772 if (symmetric) { 773 ierr = PetscFree(*ia);CHKERRQ(ierr); 774 ierr = PetscFree(*ja);CHKERRQ(ierr); 775 } else if (oshift == 1) { 776 PetscInt nz = a->i[n]-1; 777 for (i=0; i<nz; i++) a->j[i]--; 778 for (i=0; i<n+1; i++) a->i[i]--; 779 } 780 PetscFunctionReturn(0); 781 } 782 783 #undef __FUNCT__ 784 #define __FUNCT__ "MatDestroy_SeqBAIJ" 785 PetscErrorCode MatDestroy_SeqBAIJ(Mat A) 786 { 787 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 788 PetscErrorCode ierr; 789 790 PetscFunctionBegin; 791 #if defined(PETSC_USE_LOG) 792 PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->m,A->n,a->nz); 793 #endif 794 ierr = PetscFree(a->a);CHKERRQ(ierr); 795 if (!a->singlemalloc) { 796 ierr = PetscFree(a->i);CHKERRQ(ierr); 797 ierr = PetscFree(a->j);CHKERRQ(ierr); 798 } 799 if (a->row) { 800 ierr = ISDestroy(a->row);CHKERRQ(ierr); 801 } 802 if (a->col) { 803 ierr = ISDestroy(a->col);CHKERRQ(ierr); 804 } 805 if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 806 if (a->idiag) {ierr = PetscFree(a->idiag);CHKERRQ(ierr);} 807 if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 808 if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 809 if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 810 if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 811 if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 812 if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 813 #if defined(PETSC_USE_MAT_SINGLE) 814 if (a->setvaluescopy) {ierr = PetscFree(a->setvaluescopy);CHKERRQ(ierr);} 815 #endif 816 if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);} 817 if (a->compressedrow.use){ierr = PetscFree(a->compressedrow.i);} 818 819 ierr = PetscFree(a);CHKERRQ(ierr); 820 821 ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr); 822 ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr); 823 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr); 824 ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 825 ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqbaij_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr); 826 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 827 PetscFunctionReturn(0); 828 } 829 830 #undef __FUNCT__ 831 #define __FUNCT__ "MatSetOption_SeqBAIJ" 832 PetscErrorCode MatSetOption_SeqBAIJ(Mat A,MatOption op) 833 { 834 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 835 836 PetscFunctionBegin; 837 switch (op) { 838 case MAT_ROW_ORIENTED: 839 a->roworiented = PETSC_TRUE; 840 break; 841 case MAT_COLUMN_ORIENTED: 842 a->roworiented = PETSC_FALSE; 843 break; 844 case MAT_COLUMNS_SORTED: 845 a->sorted = PETSC_TRUE; 846 break; 847 case MAT_COLUMNS_UNSORTED: 848 a->sorted = PETSC_FALSE; 849 break; 850 case MAT_KEEP_ZEROED_ROWS: 851 a->keepzeroedrows = PETSC_TRUE; 852 break; 853 case MAT_NO_NEW_NONZERO_LOCATIONS: 854 a->nonew = 1; 855 break; 856 case MAT_NEW_NONZERO_LOCATION_ERR: 857 a->nonew = -1; 858 break; 859 case MAT_NEW_NONZERO_ALLOCATION_ERR: 860 a->nonew = -2; 861 break; 862 case MAT_YES_NEW_NONZERO_LOCATIONS: 863 a->nonew = 0; 864 break; 865 case MAT_ROWS_SORTED: 866 case MAT_ROWS_UNSORTED: 867 case MAT_YES_NEW_DIAGONALS: 868 case MAT_IGNORE_OFF_PROC_ENTRIES: 869 case MAT_USE_HASH_TABLE: 870 PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 871 break; 872 case MAT_NO_NEW_DIAGONALS: 873 SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 874 case MAT_SYMMETRIC: 875 case MAT_STRUCTURALLY_SYMMETRIC: 876 case MAT_NOT_SYMMETRIC: 877 case MAT_NOT_STRUCTURALLY_SYMMETRIC: 878 case MAT_HERMITIAN: 879 case MAT_NOT_HERMITIAN: 880 case MAT_SYMMETRY_ETERNAL: 881 case MAT_NOT_SYMMETRY_ETERNAL: 882 break; 883 default: 884 SETERRQ(PETSC_ERR_SUP,"unknown option"); 885 } 886 PetscFunctionReturn(0); 887 } 888 889 #undef __FUNCT__ 890 #define __FUNCT__ "MatGetRow_SeqBAIJ" 891 PetscErrorCode MatGetRow_SeqBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 892 { 893 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 894 PetscErrorCode ierr; 895 PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 896 MatScalar *aa,*aa_i; 897 PetscScalar *v_i; 898 899 PetscFunctionBegin; 900 bs = A->bs; 901 ai = a->i; 902 aj = a->j; 903 aa = a->a; 904 bs2 = a->bs2; 905 906 if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %D out of range", row); 907 908 bn = row/bs; /* Block number */ 909 bp = row % bs; /* Block Position */ 910 M = ai[bn+1] - ai[bn]; 911 *nz = bs*M; 912 913 if (v) { 914 *v = 0; 915 if (*nz) { 916 ierr = PetscMalloc((*nz)*sizeof(PetscScalar),v);CHKERRQ(ierr); 917 for (i=0; i<M; i++) { /* for each block in the block row */ 918 v_i = *v + i*bs; 919 aa_i = aa + bs2*(ai[bn] + i); 920 for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 921 } 922 } 923 } 924 925 if (idx) { 926 *idx = 0; 927 if (*nz) { 928 ierr = PetscMalloc((*nz)*sizeof(PetscInt),idx);CHKERRQ(ierr); 929 for (i=0; i<M; i++) { /* for each block in the block row */ 930 idx_i = *idx + i*bs; 931 itmp = bs*aj[ai[bn] + i]; 932 for (j=0; j<bs; j++) {idx_i[j] = itmp++;} 933 } 934 } 935 } 936 PetscFunctionReturn(0); 937 } 938 939 #undef __FUNCT__ 940 #define __FUNCT__ "MatRestoreRow_SeqBAIJ" 941 PetscErrorCode MatRestoreRow_SeqBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 942 { 943 PetscErrorCode ierr; 944 945 PetscFunctionBegin; 946 if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 947 if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 948 PetscFunctionReturn(0); 949 } 950 951 #undef __FUNCT__ 952 #define __FUNCT__ "MatTranspose_SeqBAIJ" 953 PetscErrorCode MatTranspose_SeqBAIJ(Mat A,Mat *B) 954 { 955 Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 956 Mat C; 957 PetscErrorCode ierr; 958 PetscInt i,j,k,*aj=a->j,*ai=a->i,bs=A->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 959 PetscInt *rows,*cols,bs2=a->bs2; 960 PetscScalar *array; 961 962 PetscFunctionBegin; 963 if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Square matrix only for in-place"); 964 ierr = PetscMalloc((1+nbs)*sizeof(PetscInt),&col);CHKERRQ(ierr); 965 ierr = PetscMemzero(col,(1+nbs)*sizeof(PetscInt));CHKERRQ(ierr); 966 967 #if defined(PETSC_USE_MAT_SINGLE) 968 ierr = PetscMalloc(a->bs2*a->nz*sizeof(PetscScalar),&array);CHKERRQ(ierr); 969 for (i=0; i<a->bs2*a->nz; i++) array[i] = (PetscScalar)a->a[i]; 970 #else 971 array = a->a; 972 #endif 973 974 for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1; 975 ierr = MatCreate(A->comm,A->n,A->m,A->n,A->m,&C);CHKERRQ(ierr); 976 ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 977 ierr = MatSeqBAIJSetPreallocation(C,bs,PETSC_NULL,col);CHKERRQ(ierr); 978 ierr = PetscFree(col);CHKERRQ(ierr); 979 ierr = PetscMalloc(2*bs*sizeof(PetscInt),&rows);CHKERRQ(ierr); 980 cols = rows + bs; 981 for (i=0; i<mbs; i++) { 982 cols[0] = i*bs; 983 for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1; 984 len = ai[i+1] - ai[i]; 985 for (j=0; j<len; j++) { 986 rows[0] = (*aj++)*bs; 987 for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1; 988 ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 989 array += bs2; 990 } 991 } 992 ierr = PetscFree(rows);CHKERRQ(ierr); 993 #if defined(PETSC_USE_MAT_SINGLE) 994 ierr = PetscFree(array); 995 #endif 996 997 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 998 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 999 1000 if (B) { 1001 *B = C; 1002 } else { 1003 ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 1004 } 1005 PetscFunctionReturn(0); 1006 } 1007 1008 #undef __FUNCT__ 1009 #define __FUNCT__ "MatView_SeqBAIJ_Binary" 1010 static PetscErrorCode MatView_SeqBAIJ_Binary(Mat A,PetscViewer viewer) 1011 { 1012 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1013 PetscErrorCode ierr; 1014 PetscInt i,*col_lens,bs = A->bs,count,*jj,j,k,l,bs2=a->bs2; 1015 int fd; 1016 PetscScalar *aa; 1017 FILE *file; 1018 1019 PetscFunctionBegin; 1020 ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 1021 ierr = PetscMalloc((4+A->m)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 1022 col_lens[0] = MAT_FILE_COOKIE; 1023 1024 col_lens[1] = A->m; 1025 col_lens[2] = A->n; 1026 col_lens[3] = a->nz*bs2; 1027 1028 /* store lengths of each row and write (including header) to file */ 1029 count = 0; 1030 for (i=0; i<a->mbs; i++) { 1031 for (j=0; j<bs; j++) { 1032 col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 1033 } 1034 } 1035 ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 1036 ierr = PetscFree(col_lens);CHKERRQ(ierr); 1037 1038 /* store column indices (zero start index) */ 1039 ierr = PetscMalloc((a->nz+1)*bs2*sizeof(PetscInt),&jj);CHKERRQ(ierr); 1040 count = 0; 1041 for (i=0; i<a->mbs; i++) { 1042 for (j=0; j<bs; j++) { 1043 for (k=a->i[i]; k<a->i[i+1]; k++) { 1044 for (l=0; l<bs; l++) { 1045 jj[count++] = bs*a->j[k] + l; 1046 } 1047 } 1048 } 1049 } 1050 ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 1051 ierr = PetscFree(jj);CHKERRQ(ierr); 1052 1053 /* store nonzero values */ 1054 ierr = PetscMalloc((a->nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 1055 count = 0; 1056 for (i=0; i<a->mbs; i++) { 1057 for (j=0; j<bs; j++) { 1058 for (k=a->i[i]; k<a->i[i+1]; k++) { 1059 for (l=0; l<bs; l++) { 1060 aa[count++] = a->a[bs2*k + l*bs + j]; 1061 } 1062 } 1063 } 1064 } 1065 ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 1066 ierr = PetscFree(aa);CHKERRQ(ierr); 1067 1068 ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 1069 if (file) { 1070 fprintf(file,"-matload_block_size %d\n",(int)A->bs); 1071 } 1072 PetscFunctionReturn(0); 1073 } 1074 1075 #undef __FUNCT__ 1076 #define __FUNCT__ "MatView_SeqBAIJ_ASCII" 1077 static PetscErrorCode MatView_SeqBAIJ_ASCII(Mat A,PetscViewer viewer) 1078 { 1079 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1080 PetscErrorCode ierr; 1081 PetscInt i,j,bs = A->bs,k,l,bs2=a->bs2; 1082 PetscViewerFormat format; 1083 1084 PetscFunctionBegin; 1085 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1086 if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1087 ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 1088 } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 1089 Mat aij; 1090 ierr = MatConvert(A,MATSEQAIJ,&aij);CHKERRQ(ierr); 1091 ierr = MatView(aij,viewer);CHKERRQ(ierr); 1092 ierr = MatDestroy(aij);CHKERRQ(ierr); 1093 } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 1094 PetscFunctionReturn(0); 1095 } else if (format == PETSC_VIEWER_ASCII_COMMON) { 1096 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 1097 for (i=0; i<a->mbs; i++) { 1098 for (j=0; j<bs; j++) { 1099 ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 1100 for (k=a->i[i]; k<a->i[i+1]; k++) { 1101 for (l=0; l<bs; l++) { 1102 #if defined(PETSC_USE_COMPLEX) 1103 if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 1104 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %gi) ",bs*a->j[k]+l, 1105 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 1106 } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 1107 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %gi) ",bs*a->j[k]+l, 1108 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 1109 } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 1110 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 1111 } 1112 #else 1113 if (a->a[bs2*k + l*bs + j] != 0.0) { 1114 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 1115 } 1116 #endif 1117 } 1118 } 1119 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1120 } 1121 } 1122 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 1123 } else { 1124 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 1125 for (i=0; i<a->mbs; i++) { 1126 for (j=0; j<bs; j++) { 1127 ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 1128 for (k=a->i[i]; k<a->i[i+1]; k++) { 1129 for (l=0; l<bs; l++) { 1130 #if defined(PETSC_USE_COMPLEX) 1131 if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 1132 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 1133 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 1134 } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 1135 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l, 1136 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 1137 } else { 1138 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 1139 } 1140 #else 1141 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 1142 #endif 1143 } 1144 } 1145 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1146 } 1147 } 1148 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 1149 } 1150 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1151 PetscFunctionReturn(0); 1152 } 1153 1154 #undef __FUNCT__ 1155 #define __FUNCT__ "MatView_SeqBAIJ_Draw_Zoom" 1156 static PetscErrorCode MatView_SeqBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 1157 { 1158 Mat A = (Mat) Aa; 1159 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 1160 PetscErrorCode ierr; 1161 PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->bs,bs2=a->bs2; 1162 PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 1163 MatScalar *aa; 1164 PetscViewer viewer; 1165 1166 PetscFunctionBegin; 1167 1168 /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/ 1169 ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 1170 1171 ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1172 1173 /* loop over matrix elements drawing boxes */ 1174 color = PETSC_DRAW_BLUE; 1175 for (i=0,row=0; i<mbs; i++,row+=bs) { 1176 for (j=a->i[i]; j<a->i[i+1]; j++) { 1177 y_l = A->m - row - 1.0; y_r = y_l + 1.0; 1178 x_l = a->j[j]*bs; x_r = x_l + 1.0; 1179 aa = a->a + j*bs2; 1180 for (k=0; k<bs; k++) { 1181 for (l=0; l<bs; l++) { 1182 if (PetscRealPart(*aa++) >= 0.) continue; 1183 ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 1184 } 1185 } 1186 } 1187 } 1188 color = PETSC_DRAW_CYAN; 1189 for (i=0,row=0; i<mbs; i++,row+=bs) { 1190 for (j=a->i[i]; j<a->i[i+1]; j++) { 1191 y_l = A->m - row - 1.0; y_r = y_l + 1.0; 1192 x_l = a->j[j]*bs; x_r = x_l + 1.0; 1193 aa = a->a + j*bs2; 1194 for (k=0; k<bs; k++) { 1195 for (l=0; l<bs; l++) { 1196 if (PetscRealPart(*aa++) != 0.) continue; 1197 ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 1198 } 1199 } 1200 } 1201 } 1202 1203 color = PETSC_DRAW_RED; 1204 for (i=0,row=0; i<mbs; i++,row+=bs) { 1205 for (j=a->i[i]; j<a->i[i+1]; j++) { 1206 y_l = A->m - row - 1.0; y_r = y_l + 1.0; 1207 x_l = a->j[j]*bs; x_r = x_l + 1.0; 1208 aa = a->a + j*bs2; 1209 for (k=0; k<bs; k++) { 1210 for (l=0; l<bs; l++) { 1211 if (PetscRealPart(*aa++) <= 0.) continue; 1212 ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 1213 } 1214 } 1215 } 1216 } 1217 PetscFunctionReturn(0); 1218 } 1219 1220 #undef __FUNCT__ 1221 #define __FUNCT__ "MatView_SeqBAIJ_Draw" 1222 static PetscErrorCode MatView_SeqBAIJ_Draw(Mat A,PetscViewer viewer) 1223 { 1224 PetscErrorCode ierr; 1225 PetscReal xl,yl,xr,yr,w,h; 1226 PetscDraw draw; 1227 PetscTruth isnull; 1228 1229 PetscFunctionBegin; 1230 1231 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1232 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 1233 1234 ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1235 xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 1236 xr += w; yr += h; xl = -w; yl = -h; 1237 ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1238 ierr = PetscDrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 1239 ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 1240 PetscFunctionReturn(0); 1241 } 1242 1243 #undef __FUNCT__ 1244 #define __FUNCT__ "MatView_SeqBAIJ" 1245 PetscErrorCode MatView_SeqBAIJ(Mat A,PetscViewer viewer) 1246 { 1247 PetscErrorCode ierr; 1248 PetscTruth iascii,isbinary,isdraw; 1249 1250 PetscFunctionBegin; 1251 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1252 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 1253 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 1254 if (iascii){ 1255 ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 1256 } else if (isbinary) { 1257 ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 1258 } else if (isdraw) { 1259 ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 1260 } else { 1261 Mat B; 1262 ierr = MatConvert(A,MATSEQAIJ,&B);CHKERRQ(ierr); 1263 ierr = MatView(B,viewer);CHKERRQ(ierr); 1264 ierr = MatDestroy(B);CHKERRQ(ierr); 1265 } 1266 PetscFunctionReturn(0); 1267 } 1268 1269 1270 #undef __FUNCT__ 1271 #define __FUNCT__ "MatGetValues_SeqBAIJ" 1272 PetscErrorCode MatGetValues_SeqBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 1273 { 1274 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1275 PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 1276 PetscInt *ai = a->i,*ailen = a->ilen; 1277 PetscInt brow,bcol,ridx,cidx,bs=A->bs,bs2=a->bs2; 1278 MatScalar *ap,*aa = a->a,zero = 0.0; 1279 1280 PetscFunctionBegin; 1281 for (k=0; k<m; k++) { /* loop over rows */ 1282 row = im[k]; brow = row/bs; 1283 if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 1284 if (row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %D too large", row); 1285 rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 1286 nrow = ailen[brow]; 1287 for (l=0; l<n; l++) { /* loop over columns */ 1288 if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 1289 if (in[l] >= A->n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column %D too large", in[l]); 1290 col = in[l] ; 1291 bcol = col/bs; 1292 cidx = col%bs; 1293 ridx = row%bs; 1294 high = nrow; 1295 low = 0; /* assume unsorted */ 1296 while (high-low > 5) { 1297 t = (low+high)/2; 1298 if (rp[t] > bcol) high = t; 1299 else low = t; 1300 } 1301 for (i=low; i<high; i++) { 1302 if (rp[i] > bcol) break; 1303 if (rp[i] == bcol) { 1304 *v++ = ap[bs2*i+bs*cidx+ridx]; 1305 goto finished; 1306 } 1307 } 1308 *v++ = zero; 1309 finished:; 1310 } 1311 } 1312 PetscFunctionReturn(0); 1313 } 1314 1315 #if defined(PETSC_USE_MAT_SINGLE) 1316 #undef __FUNCT__ 1317 #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ" 1318 PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv) 1319 { 1320 Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data; 1321 PetscErrorCode ierr; 1322 PetscInt i,N = m*n*b->bs2; 1323 MatScalar *vsingle; 1324 1325 PetscFunctionBegin; 1326 if (N > b->setvalueslen) { 1327 if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);} 1328 ierr = PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);CHKERRQ(ierr); 1329 b->setvalueslen = N; 1330 } 1331 vsingle = b->setvaluescopy; 1332 for (i=0; i<N; i++) { 1333 vsingle[i] = v[i]; 1334 } 1335 ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr); 1336 PetscFunctionReturn(0); 1337 } 1338 #endif 1339 1340 1341 #undef __FUNCT__ 1342 #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ" 1343 PetscErrorCode MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const MatScalar v[],InsertMode is) 1344 { 1345 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1346 PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 1347 PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 1348 PetscErrorCode ierr; 1349 PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->bs,stepval; 1350 PetscTruth roworiented=a->roworiented; 1351 const MatScalar *value = v; 1352 MatScalar *ap,*aa = a->a,*bap; 1353 1354 PetscFunctionBegin; 1355 if (roworiented) { 1356 stepval = (n-1)*bs; 1357 } else { 1358 stepval = (m-1)*bs; 1359 } 1360 for (k=0; k<m; k++) { /* loop over added rows */ 1361 row = im[k]; 1362 if (row < 0) continue; 1363 #if defined(PETSC_USE_DEBUG) 1364 if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 1365 #endif 1366 rp = aj + ai[row]; 1367 ap = aa + bs2*ai[row]; 1368 rmax = imax[row]; 1369 nrow = ailen[row]; 1370 low = 0; 1371 for (l=0; l<n; l++) { /* loop over added columns */ 1372 if (in[l] < 0) continue; 1373 #if defined(PETSC_USE_DEBUG) 1374 if (in[l] >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],a->nbs-1); 1375 #endif 1376 col = in[l]; 1377 if (roworiented) { 1378 value = v + k*(stepval+bs)*bs + l*bs; 1379 } else { 1380 value = v + l*(stepval+bs)*bs + k*bs; 1381 } 1382 if (!sorted) low = 0; high = nrow; 1383 while (high-low > 7) { 1384 t = (low+high)/2; 1385 if (rp[t] > col) high = t; 1386 else low = t; 1387 } 1388 for (i=low; i<high; i++) { 1389 if (rp[i] > col) break; 1390 if (rp[i] == col) { 1391 bap = ap + bs2*i; 1392 if (roworiented) { 1393 if (is == ADD_VALUES) { 1394 for (ii=0; ii<bs; ii++,value+=stepval) { 1395 for (jj=ii; jj<bs2; jj+=bs) { 1396 bap[jj] += *value++; 1397 } 1398 } 1399 } else { 1400 for (ii=0; ii<bs; ii++,value+=stepval) { 1401 for (jj=ii; jj<bs2; jj+=bs) { 1402 bap[jj] = *value++; 1403 } 1404 } 1405 } 1406 } else { 1407 if (is == ADD_VALUES) { 1408 for (ii=0; ii<bs; ii++,value+=stepval) { 1409 for (jj=0; jj<bs; jj++) { 1410 *bap++ += *value++; 1411 } 1412 } 1413 } else { 1414 for (ii=0; ii<bs; ii++,value+=stepval) { 1415 for (jj=0; jj<bs; jj++) { 1416 *bap++ = *value++; 1417 } 1418 } 1419 } 1420 } 1421 goto noinsert2; 1422 } 1423 } 1424 if (nonew == 1) goto noinsert2; 1425 else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 1426 if (nrow >= rmax) { 1427 /* there is no extra room in row, therefore enlarge */ 1428 PetscInt new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 1429 MatScalar *new_a; 1430 1431 if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 1432 1433 /* malloc new storage space */ 1434 len = new_nz*(sizeof(PetscInt)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(PetscInt); 1435 ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 1436 new_j = (PetscInt*)(new_a + bs2*new_nz); 1437 new_i = new_j + new_nz; 1438 1439 /* copy over old data into new slots */ 1440 for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 1441 for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1442 ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(PetscInt));CHKERRQ(ierr); 1443 len = (new_nz - CHUNKSIZE - ai[row] - nrow); 1444 ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(PetscInt));CHKERRQ(ierr); 1445 ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1446 ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1447 ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 1448 /* free up old matrix storage */ 1449 ierr = PetscFree(a->a);CHKERRQ(ierr); 1450 if (!a->singlemalloc) { 1451 ierr = PetscFree(a->i);CHKERRQ(ierr); 1452 ierr = PetscFree(a->j);CHKERRQ(ierr); 1453 } 1454 aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 1455 a->singlemalloc = PETSC_TRUE; 1456 1457 rp = aj + ai[row]; ap = aa + bs2*ai[row]; 1458 rmax = imax[row] = imax[row] + CHUNKSIZE; 1459 PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(PetscInt) + bs2*sizeof(MatScalar))); 1460 a->maxnz += bs2*CHUNKSIZE; 1461 a->reallocs++; 1462 a->nz++; 1463 } 1464 N = nrow++ - 1; 1465 /* shift up all the later entries in this row */ 1466 for (ii=N; ii>=i; ii--) { 1467 rp[ii+1] = rp[ii]; 1468 ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 1469 } 1470 if (N >= i) { 1471 ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1472 } 1473 rp[i] = col; 1474 bap = ap + bs2*i; 1475 if (roworiented) { 1476 for (ii=0; ii<bs; ii++,value+=stepval) { 1477 for (jj=ii; jj<bs2; jj+=bs) { 1478 bap[jj] = *value++; 1479 } 1480 } 1481 } else { 1482 for (ii=0; ii<bs; ii++,value+=stepval) { 1483 for (jj=0; jj<bs; jj++) { 1484 *bap++ = *value++; 1485 } 1486 } 1487 } 1488 noinsert2:; 1489 low = i; 1490 } 1491 ailen[row] = nrow; 1492 } 1493 PetscFunctionReturn(0); 1494 } 1495 1496 #undef __FUNCT__ 1497 #define __FUNCT__ "MatAssemblyEnd_SeqBAIJ" 1498 PetscErrorCode MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 1499 { 1500 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1501 PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 1502 PetscInt m = A->m,*ip,N,*ailen = a->ilen; 1503 PetscErrorCode ierr; 1504 PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 1505 MatScalar *aa = a->a,*ap; 1506 PetscReal ratio=0.6; 1507 1508 PetscFunctionBegin; 1509 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 1510 1511 if (m) rmax = ailen[0]; 1512 for (i=1; i<mbs; i++) { 1513 /* move each row back by the amount of empty slots (fshift) before it*/ 1514 fshift += imax[i-1] - ailen[i-1]; 1515 rmax = PetscMax(rmax,ailen[i]); 1516 if (fshift) { 1517 ip = aj + ai[i]; ap = aa + bs2*ai[i]; 1518 N = ailen[i]; 1519 for (j=0; j<N; j++) { 1520 ip[j-fshift] = ip[j]; 1521 ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1522 } 1523 } 1524 ai[i] = ai[i-1] + ailen[i-1]; 1525 } 1526 if (mbs) { 1527 fshift += imax[mbs-1] - ailen[mbs-1]; 1528 ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 1529 } 1530 /* reset ilen and imax for each row */ 1531 for (i=0; i<mbs; i++) { 1532 ailen[i] = imax[i] = ai[i+1] - ai[i]; 1533 } 1534 a->nz = ai[mbs]; 1535 1536 /* diagonals may have moved, so kill the diagonal pointers */ 1537 a->idiagvalid = PETSC_FALSE; 1538 if (fshift && a->diag) { 1539 ierr = PetscFree(a->diag);CHKERRQ(ierr); 1540 PetscLogObjectMemory(A,-(mbs+1)*sizeof(PetscInt)); 1541 a->diag = 0; 1542 } 1543 PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",m,A->n,A->bs,fshift*bs2,a->nz*bs2); 1544 PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %D\n",a->reallocs); 1545 PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %D\n",rmax); 1546 a->reallocs = 0; 1547 A->info.nz_unneeded = (PetscReal)fshift*bs2; 1548 1549 /* check for zero rows. If found a large number of zero rows, use CompressedRow functions */ 1550 if (a->compressedrow.use){ 1551 ierr = Mat_CheckCompressedRow(A,&a->compressedrow,a->i,mbs,ratio);CHKERRQ(ierr); 1552 } 1553 1554 A->same_nonzero = PETSC_TRUE; 1555 PetscFunctionReturn(0); 1556 } 1557 1558 /* 1559 This function returns an array of flags which indicate the locations of contiguous 1560 blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 1561 then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 1562 Assume: sizes should be long enough to hold all the values. 1563 */ 1564 #undef __FUNCT__ 1565 #define __FUNCT__ "MatZeroRows_SeqBAIJ_Check_Blocks" 1566 static PetscErrorCode MatZeroRows_SeqBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 1567 { 1568 PetscInt i,j,k,row; 1569 PetscTruth flg; 1570 1571 PetscFunctionBegin; 1572 for (i=0,j=0; i<n; j++) { 1573 row = idx[i]; 1574 if (row%bs!=0) { /* Not the begining of a block */ 1575 sizes[j] = 1; 1576 i++; 1577 } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 1578 sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 1579 i++; 1580 } else { /* Begining of the block, so check if the complete block exists */ 1581 flg = PETSC_TRUE; 1582 for (k=1; k<bs; k++) { 1583 if (row+k != idx[i+k]) { /* break in the block */ 1584 flg = PETSC_FALSE; 1585 break; 1586 } 1587 } 1588 if (flg) { /* No break in the bs */ 1589 sizes[j] = bs; 1590 i+= bs; 1591 } else { 1592 sizes[j] = 1; 1593 i++; 1594 } 1595 } 1596 } 1597 *bs_max = j; 1598 PetscFunctionReturn(0); 1599 } 1600 1601 #undef __FUNCT__ 1602 #define __FUNCT__ "MatZeroRows_SeqBAIJ" 1603 PetscErrorCode MatZeroRows_SeqBAIJ(Mat A,IS is,const PetscScalar *diag) 1604 { 1605 Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 1606 PetscErrorCode ierr; 1607 PetscInt i,j,k,count,is_n,*is_idx,*rows; 1608 PetscInt bs=A->bs,bs2=baij->bs2,*sizes,row,bs_max; 1609 PetscScalar zero = 0.0; 1610 MatScalar *aa; 1611 1612 PetscFunctionBegin; 1613 /* Make a copy of the IS and sort it */ 1614 ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr); 1615 ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 1616 1617 /* allocate memory for rows,sizes */ 1618 ierr = PetscMalloc((3*is_n+1)*sizeof(PetscInt),&rows);CHKERRQ(ierr); 1619 sizes = rows + is_n; 1620 1621 /* copy IS values to rows, and sort them */ 1622 for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 1623 ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 1624 if (baij->keepzeroedrows) { 1625 for (i=0; i<is_n; i++) { sizes[i] = 1; } 1626 bs_max = is_n; 1627 A->same_nonzero = PETSC_TRUE; 1628 } else { 1629 ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 1630 A->same_nonzero = PETSC_FALSE; 1631 } 1632 ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 1633 1634 for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 1635 row = rows[j]; 1636 if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range",row); 1637 count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 1638 aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 1639 if (sizes[i] == bs && !baij->keepzeroedrows) { 1640 if (diag) { 1641 if (baij->ilen[row/bs] > 0) { 1642 baij->ilen[row/bs] = 1; 1643 baij->j[baij->i[row/bs]] = row/bs; 1644 ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 1645 } 1646 /* Now insert all the diagonal values for this bs */ 1647 for (k=0; k<bs; k++) { 1648 ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 1649 } 1650 } else { /* (!diag) */ 1651 baij->ilen[row/bs] = 0; 1652 } /* end (!diag) */ 1653 } else { /* (sizes[i] != bs) */ 1654 #if defined (PETSC_USE_DEBUG) 1655 if (sizes[i] != 1) SETERRQ(PETSC_ERR_PLIB,"Internal Error. Value should be 1"); 1656 #endif 1657 for (k=0; k<count; k++) { 1658 aa[0] = zero; 1659 aa += bs; 1660 } 1661 if (diag) { 1662 ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 1663 } 1664 } 1665 } 1666 1667 ierr = PetscFree(rows);CHKERRQ(ierr); 1668 ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1669 PetscFunctionReturn(0); 1670 } 1671 1672 #undef __FUNCT__ 1673 #define __FUNCT__ "MatSetValues_SeqBAIJ" 1674 PetscErrorCode MatSetValues_SeqBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 1675 { 1676 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1677 PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 1678 PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 1679 PetscInt *aj=a->j,nonew=a->nonew,bs=A->bs,brow,bcol; 1680 PetscErrorCode ierr; 1681 PetscInt ridx,cidx,bs2=a->bs2; 1682 PetscTruth roworiented=a->roworiented; 1683 MatScalar *ap,value,*aa=a->a,*bap; 1684 1685 PetscFunctionBegin; 1686 for (k=0; k<m; k++) { /* loop over added rows */ 1687 row = im[k]; brow = row/bs; 1688 if (row < 0) continue; 1689 #if defined(PETSC_USE_DEBUG) 1690 if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 1691 #endif 1692 rp = aj + ai[brow]; 1693 ap = aa + bs2*ai[brow]; 1694 rmax = imax[brow]; 1695 nrow = ailen[brow]; 1696 low = 0; 1697 for (l=0; l<n; l++) { /* loop over added columns */ 1698 if (in[l] < 0) continue; 1699 #if defined(PETSC_USE_DEBUG) 1700 if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->n-1); 1701 #endif 1702 col = in[l]; bcol = col/bs; 1703 ridx = row % bs; cidx = col % bs; 1704 if (roworiented) { 1705 value = v[l + k*n]; 1706 } else { 1707 value = v[k + l*m]; 1708 } 1709 if (!sorted) low = 0; high = nrow; 1710 while (high-low > 7) { 1711 t = (low+high)/2; 1712 if (rp[t] > bcol) high = t; 1713 else low = t; 1714 } 1715 for (i=low; i<high; i++) { 1716 if (rp[i] > bcol) break; 1717 if (rp[i] == bcol) { 1718 bap = ap + bs2*i + bs*cidx + ridx; 1719 if (is == ADD_VALUES) *bap += value; 1720 else *bap = value; 1721 goto noinsert1; 1722 } 1723 } 1724 if (nonew == 1) goto noinsert1; 1725 else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 1726 if (nrow >= rmax) { 1727 /* there is no extra room in row, therefore enlarge */ 1728 PetscInt new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 1729 MatScalar *new_a; 1730 1731 if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 1732 1733 /* Malloc new storage space */ 1734 len = new_nz*(sizeof(PetscInt)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(PetscInt); 1735 ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 1736 new_j = (PetscInt*)(new_a + bs2*new_nz); 1737 new_i = new_j + new_nz; 1738 1739 /* copy over old data into new slots */ 1740 for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 1741 for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1742 ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(PetscInt));CHKERRQ(ierr); 1743 len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1744 ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(PetscInt));CHKERRQ(ierr); 1745 ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1746 ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1747 ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 1748 /* free up old matrix storage */ 1749 ierr = PetscFree(a->a);CHKERRQ(ierr); 1750 if (!a->singlemalloc) { 1751 ierr = PetscFree(a->i);CHKERRQ(ierr); 1752 ierr = PetscFree(a->j);CHKERRQ(ierr); 1753 } 1754 aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 1755 a->singlemalloc = PETSC_TRUE; 1756 1757 rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 1758 rmax = imax[brow] = imax[brow] + CHUNKSIZE; 1759 PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(PetscInt) + bs2*sizeof(MatScalar))); 1760 a->maxnz += bs2*CHUNKSIZE; 1761 a->reallocs++; 1762 a->nz++; 1763 } 1764 N = nrow++ - 1; 1765 /* shift up all the later entries in this row */ 1766 for (ii=N; ii>=i; ii--) { 1767 rp[ii+1] = rp[ii]; 1768 ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 1769 } 1770 if (N>=i) { 1771 ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1772 } 1773 rp[i] = bcol; 1774 ap[bs2*i + bs*cidx + ridx] = value; 1775 noinsert1:; 1776 low = i; 1777 } 1778 ailen[brow] = nrow; 1779 } 1780 A->same_nonzero = PETSC_FALSE; 1781 PetscFunctionReturn(0); 1782 } 1783 1784 1785 #undef __FUNCT__ 1786 #define __FUNCT__ "MatILUFactor_SeqBAIJ" 1787 PetscErrorCode MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatFactorInfo *info) 1788 { 1789 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 1790 Mat outA; 1791 PetscErrorCode ierr; 1792 PetscTruth row_identity,col_identity; 1793 1794 PetscFunctionBegin; 1795 if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 1796 ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1797 ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1798 if (!row_identity || !col_identity) { 1799 SETERRQ(PETSC_ERR_ARG_WRONG,"Row and column permutations must be identity for in-place ILU"); 1800 } 1801 1802 outA = inA; 1803 inA->factor = FACTOR_LU; 1804 1805 if (!a->diag) { 1806 ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr); 1807 } 1808 1809 a->row = row; 1810 a->col = col; 1811 ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1812 ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 1813 1814 /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1815 ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr); 1816 PetscLogObjectParent(inA,a->icol); 1817 1818 /* 1819 Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 1820 for ILU(0) factorization with natural ordering 1821 */ 1822 if (inA->bs < 8) { 1823 ierr = MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(inA);CHKERRQ(ierr); 1824 } else { 1825 if (!a->solve_work) { 1826 ierr = PetscMalloc((inA->m+inA->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 1827 PetscLogObjectMemory(inA,(inA->m+inA->bs)*sizeof(PetscScalar)); 1828 } 1829 } 1830 1831 ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1832 1833 PetscFunctionReturn(0); 1834 } 1835 #undef __FUNCT__ 1836 #define __FUNCT__ "MatPrintHelp_SeqBAIJ" 1837 PetscErrorCode MatPrintHelp_SeqBAIJ(Mat A) 1838 { 1839 static PetscTruth called = PETSC_FALSE; 1840 MPI_Comm comm = A->comm; 1841 PetscErrorCode ierr; 1842 1843 PetscFunctionBegin; 1844 if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 1845 ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1846 ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 1847 PetscFunctionReturn(0); 1848 } 1849 1850 EXTERN_C_BEGIN 1851 #undef __FUNCT__ 1852 #define __FUNCT__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 1853 PetscErrorCode MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,PetscInt *indices) 1854 { 1855 Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 1856 PetscInt i,nz,nbs; 1857 1858 PetscFunctionBegin; 1859 nz = baij->maxnz/baij->bs2; 1860 nbs = baij->nbs; 1861 for (i=0; i<nz; i++) { 1862 baij->j[i] = indices[i]; 1863 } 1864 baij->nz = nz; 1865 for (i=0; i<nbs; i++) { 1866 baij->ilen[i] = baij->imax[i]; 1867 } 1868 1869 PetscFunctionReturn(0); 1870 } 1871 EXTERN_C_END 1872 1873 #undef __FUNCT__ 1874 #define __FUNCT__ "MatSeqBAIJSetColumnIndices" 1875 /*@ 1876 MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 1877 in the matrix. 1878 1879 Input Parameters: 1880 + mat - the SeqBAIJ matrix 1881 - indices - the column indices 1882 1883 Level: advanced 1884 1885 Notes: 1886 This can be called if you have precomputed the nonzero structure of the 1887 matrix and want to provide it to the matrix object to improve the performance 1888 of the MatSetValues() operation. 1889 1890 You MUST have set the correct numbers of nonzeros per row in the call to 1891 MatCreateSeqBAIJ(). 1892 1893 MUST be called before any calls to MatSetValues(); 1894 1895 @*/ 1896 PetscErrorCode MatSeqBAIJSetColumnIndices(Mat mat,PetscInt *indices) 1897 { 1898 PetscErrorCode ierr,(*f)(Mat,PetscInt *); 1899 1900 PetscFunctionBegin; 1901 PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 1902 PetscValidPointer(indices,2); 1903 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 1904 if (f) { 1905 ierr = (*f)(mat,indices);CHKERRQ(ierr); 1906 } else { 1907 SETERRQ(PETSC_ERR_ARG_WRONG,"Wrong type of matrix to set column indices"); 1908 } 1909 PetscFunctionReturn(0); 1910 } 1911 1912 #undef __FUNCT__ 1913 #define __FUNCT__ "MatGetRowMax_SeqBAIJ" 1914 PetscErrorCode MatGetRowMax_SeqBAIJ(Mat A,Vec v) 1915 { 1916 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1917 PetscErrorCode ierr; 1918 PetscInt i,j,n,row,bs,*ai,*aj,mbs; 1919 PetscReal atmp; 1920 PetscScalar *x,zero = 0.0; 1921 MatScalar *aa; 1922 PetscInt ncols,brow,krow,kcol; 1923 1924 PetscFunctionBegin; 1925 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1926 bs = A->bs; 1927 aa = a->a; 1928 ai = a->i; 1929 aj = a->j; 1930 mbs = a->mbs; 1931 1932 ierr = VecSet(&zero,v);CHKERRQ(ierr); 1933 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1934 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1935 if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1936 for (i=0; i<mbs; i++) { 1937 ncols = ai[1] - ai[0]; ai++; 1938 brow = bs*i; 1939 for (j=0; j<ncols; j++){ 1940 /* bcol = bs*(*aj); */ 1941 for (kcol=0; kcol<bs; kcol++){ 1942 for (krow=0; krow<bs; krow++){ 1943 atmp = PetscAbsScalar(*aa); aa++; 1944 row = brow + krow; /* row index */ 1945 /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */ 1946 if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp; 1947 } 1948 } 1949 aj++; 1950 } 1951 } 1952 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1953 PetscFunctionReturn(0); 1954 } 1955 1956 #undef __FUNCT__ 1957 #define __FUNCT__ "MatSetUpPreallocation_SeqBAIJ" 1958 PetscErrorCode MatSetUpPreallocation_SeqBAIJ(Mat A) 1959 { 1960 PetscErrorCode ierr; 1961 1962 PetscFunctionBegin; 1963 ierr = MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1964 PetscFunctionReturn(0); 1965 } 1966 1967 #undef __FUNCT__ 1968 #define __FUNCT__ "MatGetArray_SeqBAIJ" 1969 PetscErrorCode MatGetArray_SeqBAIJ(Mat A,PetscScalar *array[]) 1970 { 1971 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1972 PetscFunctionBegin; 1973 *array = a->a; 1974 PetscFunctionReturn(0); 1975 } 1976 1977 #undef __FUNCT__ 1978 #define __FUNCT__ "MatRestoreArray_SeqBAIJ" 1979 PetscErrorCode MatRestoreArray_SeqBAIJ(Mat A,PetscScalar *array[]) 1980 { 1981 PetscFunctionBegin; 1982 PetscFunctionReturn(0); 1983 } 1984 1985 #include "petscblaslapack.h" 1986 #undef __FUNCT__ 1987 #define __FUNCT__ "MatAXPY_SeqBAIJ" 1988 PetscErrorCode MatAXPY_SeqBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str) 1989 { 1990 Mat_SeqBAIJ *x = (Mat_SeqBAIJ *)X->data,*y = (Mat_SeqBAIJ *)Y->data; 1991 PetscErrorCode ierr; 1992 PetscInt i,bs=Y->bs,j,bs2; 1993 PetscBLASInt one=1,bnz = (PetscBLASInt)x->nz; 1994 1995 PetscFunctionBegin; 1996 if (str == SAME_NONZERO_PATTERN) { 1997 BLASaxpy_(&bnz,(PetscScalar*)a,x->a,&one,y->a,&one); 1998 } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1999 if (y->xtoy && y->XtoY != X) { 2000 ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 2001 ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 2002 } 2003 if (!y->xtoy) { /* get xtoy */ 2004 ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 2005 y->XtoY = X; 2006 } 2007 bs2 = bs*bs; 2008 for (i=0; i<x->nz; i++) { 2009 j = 0; 2010 while (j < bs2){ 2011 y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]); 2012 j++; 2013 } 2014 } 2015 PetscLogInfo(0,"MatAXPY_SeqBAIJ: ratio of nnz(X)/nnz(Y): %D/%D = %g\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz)); 2016 } else { 2017 ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 2018 } 2019 PetscFunctionReturn(0); 2020 } 2021 2022 /* -------------------------------------------------------------------*/ 2023 static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 2024 MatGetRow_SeqBAIJ, 2025 MatRestoreRow_SeqBAIJ, 2026 MatMult_SeqBAIJ_N, 2027 /* 4*/ MatMultAdd_SeqBAIJ_N, 2028 MatMultTranspose_SeqBAIJ, 2029 MatMultTransposeAdd_SeqBAIJ, 2030 MatSolve_SeqBAIJ_N, 2031 0, 2032 0, 2033 /*10*/ 0, 2034 MatLUFactor_SeqBAIJ, 2035 0, 2036 0, 2037 MatTranspose_SeqBAIJ, 2038 /*15*/ MatGetInfo_SeqBAIJ, 2039 MatEqual_SeqBAIJ, 2040 MatGetDiagonal_SeqBAIJ, 2041 MatDiagonalScale_SeqBAIJ, 2042 MatNorm_SeqBAIJ, 2043 /*20*/ 0, 2044 MatAssemblyEnd_SeqBAIJ, 2045 0, 2046 MatSetOption_SeqBAIJ, 2047 MatZeroEntries_SeqBAIJ, 2048 /*25*/ MatZeroRows_SeqBAIJ, 2049 MatLUFactorSymbolic_SeqBAIJ, 2050 MatLUFactorNumeric_SeqBAIJ_N, 2051 MatCholeskyFactorSymbolic_SeqBAIJ, 2052 MatCholeskyFactorNumeric_SeqBAIJ_N, 2053 /*30*/ MatSetUpPreallocation_SeqBAIJ, 2054 MatILUFactorSymbolic_SeqBAIJ, 2055 MatICCFactorSymbolic_SeqBAIJ, 2056 MatGetArray_SeqBAIJ, 2057 MatRestoreArray_SeqBAIJ, 2058 /*35*/ MatDuplicate_SeqBAIJ, 2059 0, 2060 0, 2061 MatILUFactor_SeqBAIJ, 2062 0, 2063 /*40*/ MatAXPY_SeqBAIJ, 2064 MatGetSubMatrices_SeqBAIJ, 2065 MatIncreaseOverlap_SeqBAIJ, 2066 MatGetValues_SeqBAIJ, 2067 0, 2068 /*45*/ MatPrintHelp_SeqBAIJ, 2069 MatScale_SeqBAIJ, 2070 0, 2071 0, 2072 0, 2073 /*50*/ 0, 2074 MatGetRowIJ_SeqBAIJ, 2075 MatRestoreRowIJ_SeqBAIJ, 2076 0, 2077 0, 2078 /*55*/ 0, 2079 0, 2080 0, 2081 0, 2082 MatSetValuesBlocked_SeqBAIJ, 2083 /*60*/ MatGetSubMatrix_SeqBAIJ, 2084 MatDestroy_SeqBAIJ, 2085 MatView_SeqBAIJ, 2086 MatGetPetscMaps_Petsc, 2087 0, 2088 /*65*/ 0, 2089 0, 2090 0, 2091 0, 2092 0, 2093 /*70*/ MatGetRowMax_SeqBAIJ, 2094 MatConvert_Basic, 2095 0, 2096 0, 2097 0, 2098 /*75*/ 0, 2099 0, 2100 0, 2101 0, 2102 0, 2103 /*80*/ 0, 2104 0, 2105 0, 2106 0, 2107 MatLoad_SeqBAIJ, 2108 /*85*/ 0, 2109 0, 2110 0, 2111 0, 2112 0, 2113 /*90*/ 0, 2114 0, 2115 0, 2116 0, 2117 0, 2118 /*95*/ 0, 2119 0, 2120 0, 2121 0}; 2122 2123 EXTERN_C_BEGIN 2124 #undef __FUNCT__ 2125 #define __FUNCT__ "MatStoreValues_SeqBAIJ" 2126 PetscErrorCode MatStoreValues_SeqBAIJ(Mat mat) 2127 { 2128 Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 2129 PetscInt nz = aij->i[mat->m]*mat->bs*aij->bs2; 2130 PetscErrorCode ierr; 2131 2132 PetscFunctionBegin; 2133 if (aij->nonew != 1) { 2134 SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 2135 } 2136 2137 /* allocate space for values if not already there */ 2138 if (!aij->saved_values) { 2139 ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 2140 } 2141 2142 /* copy values over */ 2143 ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 2144 PetscFunctionReturn(0); 2145 } 2146 EXTERN_C_END 2147 2148 EXTERN_C_BEGIN 2149 #undef __FUNCT__ 2150 #define __FUNCT__ "MatRetrieveValues_SeqBAIJ" 2151 PetscErrorCode MatRetrieveValues_SeqBAIJ(Mat mat) 2152 { 2153 Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 2154 PetscErrorCode ierr; 2155 PetscInt nz = aij->i[mat->m]*mat->bs*aij->bs2; 2156 2157 PetscFunctionBegin; 2158 if (aij->nonew != 1) { 2159 SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 2160 } 2161 if (!aij->saved_values) { 2162 SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 2163 } 2164 2165 /* copy values over */ 2166 ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 2167 PetscFunctionReturn(0); 2168 } 2169 EXTERN_C_END 2170 2171 EXTERN_C_BEGIN 2172 extern PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat,const MatType,Mat*); 2173 extern PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat,const MatType,Mat*); 2174 EXTERN_C_END 2175 2176 EXTERN_C_BEGIN 2177 #undef __FUNCT__ 2178 #define __FUNCT__ "MatSeqBAIJSetPreallocation_SeqBAIJ" 2179 PetscErrorCode MatSeqBAIJSetPreallocation_SeqBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 2180 { 2181 Mat_SeqBAIJ *b; 2182 PetscErrorCode ierr; 2183 PetscInt i,len,mbs,nbs,bs2,newbs = bs; 2184 PetscTruth flg; 2185 2186 PetscFunctionBegin; 2187 2188 B->preallocated = PETSC_TRUE; 2189 ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr); 2190 if (nnz && newbs != bs) { 2191 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot change blocksize from command line if setting nnz"); 2192 } 2193 bs = newbs; 2194 2195 mbs = B->m/bs; 2196 nbs = B->n/bs; 2197 bs2 = bs*bs; 2198 2199 if (mbs*bs!=B->m || nbs*bs!=B->n) { 2200 SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %D, cols %D must be divisible by blocksize %D",B->m,B->n,bs); 2201 } 2202 2203 if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 2204 if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 2205 if (nnz) { 2206 for (i=0; i<mbs; i++) { 2207 if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 2208 if (nnz[i] > nbs) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %D value %D rowlength %D",i,nnz[i],nbs); 2209 } 2210 } 2211 2212 b = (Mat_SeqBAIJ*)B->data; 2213 ierr = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 2214 B->ops->solve = MatSolve_SeqBAIJ_Update; 2215 B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_Update; 2216 if (!flg) { 2217 switch (bs) { 2218 case 1: 2219 B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 2220 B->ops->mult = MatMult_SeqBAIJ_1; 2221 B->ops->multadd = MatMultAdd_SeqBAIJ_1; 2222 break; 2223 case 2: 2224 B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 2225 B->ops->mult = MatMult_SeqBAIJ_2; 2226 B->ops->multadd = MatMultAdd_SeqBAIJ_2; 2227 B->ops->pbrelax = MatPBRelax_SeqBAIJ_2; 2228 break; 2229 case 3: 2230 B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 2231 B->ops->mult = MatMult_SeqBAIJ_3; 2232 B->ops->multadd = MatMultAdd_SeqBAIJ_3; 2233 B->ops->pbrelax = MatPBRelax_SeqBAIJ_3; 2234 break; 2235 case 4: 2236 B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 2237 B->ops->mult = MatMult_SeqBAIJ_4; 2238 B->ops->multadd = MatMultAdd_SeqBAIJ_4; 2239 B->ops->pbrelax = MatPBRelax_SeqBAIJ_4; 2240 break; 2241 case 5: 2242 B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 2243 B->ops->mult = MatMult_SeqBAIJ_5; 2244 B->ops->multadd = MatMultAdd_SeqBAIJ_5; 2245 B->ops->pbrelax = MatPBRelax_SeqBAIJ_5; 2246 break; 2247 case 6: 2248 B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 2249 B->ops->mult = MatMult_SeqBAIJ_6; 2250 B->ops->multadd = MatMultAdd_SeqBAIJ_6; 2251 break; 2252 case 7: 2253 B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7; 2254 B->ops->mult = MatMult_SeqBAIJ_7; 2255 B->ops->multadd = MatMultAdd_SeqBAIJ_7; 2256 break; 2257 default: 2258 B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N; 2259 B->ops->mult = MatMult_SeqBAIJ_N; 2260 B->ops->multadd = MatMultAdd_SeqBAIJ_N; 2261 break; 2262 } 2263 } 2264 B->bs = bs; 2265 b->mbs = mbs; 2266 b->nbs = nbs; 2267 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&b->imax);CHKERRQ(ierr); 2268 if (!nnz) { 2269 if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 2270 else if (nz <= 0) nz = 1; 2271 for (i=0; i<mbs; i++) b->imax[i] = nz; 2272 nz = nz*mbs; 2273 } else { 2274 nz = 0; 2275 for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 2276 } 2277 2278 /* allocate the matrix space */ 2279 len = nz*sizeof(PetscInt) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(PetscInt); 2280 ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 2281 ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 2282 b->j = (PetscInt*)(b->a + nz*bs2); 2283 ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 2284 b->i = b->j + nz; 2285 b->singlemalloc = PETSC_TRUE; 2286 2287 b->i[0] = 0; 2288 for (i=1; i<mbs+1; i++) { 2289 b->i[i] = b->i[i-1] + b->imax[i-1]; 2290 } 2291 2292 /* b->ilen will count nonzeros in each block row so far. */ 2293 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&b->ilen);CHKERRQ(ierr); 2294 PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 2295 for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 2296 2297 B->bs = bs; 2298 b->bs2 = bs2; 2299 b->mbs = mbs; 2300 b->nz = 0; 2301 b->maxnz = nz*bs2; 2302 B->info.nz_unneeded = (PetscReal)b->maxnz; 2303 PetscFunctionReturn(0); 2304 } 2305 EXTERN_C_END 2306 2307 /*MC 2308 MATSEQBAIJ - MATSEQBAIJ = "seqbaij" - A matrix type to be used for sequential block sparse matrices, based on 2309 block sparse compressed row format. 2310 2311 Options Database Keys: 2312 . -mat_type seqbaij - sets the matrix type to "seqbaij" during a call to MatSetFromOptions() 2313 2314 Level: beginner 2315 2316 .seealso: MatCreateSeqBAIJ 2317 M*/ 2318 2319 EXTERN_C_BEGIN 2320 #undef __FUNCT__ 2321 #define __FUNCT__ "MatCreate_SeqBAIJ" 2322 PetscErrorCode MatCreate_SeqBAIJ(Mat B) 2323 { 2324 PetscErrorCode ierr; 2325 PetscMPIInt size; 2326 Mat_SeqBAIJ *b; 2327 2328 PetscFunctionBegin; 2329 ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 2330 if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 2331 2332 B->m = B->M = PetscMax(B->m,B->M); 2333 B->n = B->N = PetscMax(B->n,B->N); 2334 ierr = PetscNew(Mat_SeqBAIJ,&b);CHKERRQ(ierr); 2335 B->data = (void*)b; 2336 ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 2337 B->factor = 0; 2338 B->lupivotthreshold = 1.0; 2339 B->mapping = 0; 2340 b->row = 0; 2341 b->col = 0; 2342 b->icol = 0; 2343 b->reallocs = 0; 2344 b->saved_values = 0; 2345 #if defined(PETSC_USE_MAT_SINGLE) 2346 b->setvalueslen = 0; 2347 b->setvaluescopy = PETSC_NULL; 2348 #endif 2349 2350 ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 2351 ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 2352 2353 b->sorted = PETSC_FALSE; 2354 b->roworiented = PETSC_TRUE; 2355 b->nonew = 0; 2356 b->diag = 0; 2357 b->solve_work = 0; 2358 b->mult_work = 0; 2359 B->spptr = 0; 2360 B->info.nz_unneeded = (PetscReal)b->maxnz; 2361 b->keepzeroedrows = PETSC_FALSE; 2362 b->xtoy = 0; 2363 b->XtoY = 0; 2364 b->compressedrow.use = PETSC_FALSE; 2365 b->compressedrow.nrows = 0; 2366 b->compressedrow.i = PETSC_NULL; 2367 b->compressedrow.rindex = PETSC_NULL; 2368 b->compressedrow.checked = PETSC_FALSE; 2369 B->same_nonzero = PETSC_FALSE; 2370 2371 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 2372 "MatStoreValues_SeqBAIJ", 2373 MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 2374 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 2375 "MatRetrieveValues_SeqBAIJ", 2376 MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 2377 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 2378 "MatSeqBAIJSetColumnIndices_SeqBAIJ", 2379 MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 2380 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqbaij_seqaij_C", 2381 "MatConvert_SeqBAIJ_SeqAIJ", 2382 MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr); 2383 #if !defined(PETSC_USE_64BIT_INT) 2384 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqbaij_seqsbaij_C", 2385 "MatConvert_SeqBAIJ_SeqSBAIJ", 2386 MatConvert_SeqBAIJ_SeqSBAIJ);CHKERRQ(ierr); 2387 #endif 2388 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetPreallocation_C", 2389 "MatSeqBAIJSetPreallocation_SeqBAIJ", 2390 MatSeqBAIJSetPreallocation_SeqBAIJ);CHKERRQ(ierr); 2391 PetscFunctionReturn(0); 2392 } 2393 EXTERN_C_END 2394 2395 #undef __FUNCT__ 2396 #define __FUNCT__ "MatDuplicate_SeqBAIJ" 2397 PetscErrorCode MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 2398 { 2399 Mat C; 2400 Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data; 2401 PetscErrorCode ierr; 2402 PetscInt i,len,mbs = a->mbs,nz = a->nz,bs2 = a->bs2; 2403 2404 PetscFunctionBegin; 2405 if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 2406 2407 *B = 0; 2408 ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 2409 ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 2410 ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 2411 c = (Mat_SeqBAIJ*)C->data; 2412 2413 C->M = A->M; 2414 C->N = A->N; 2415 C->bs = A->bs; 2416 c->bs2 = a->bs2; 2417 c->mbs = a->mbs; 2418 c->nbs = a->nbs; 2419 2420 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->imax);CHKERRQ(ierr); 2421 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->ilen);CHKERRQ(ierr); 2422 for (i=0; i<mbs; i++) { 2423 c->imax[i] = a->imax[i]; 2424 c->ilen[i] = a->ilen[i]; 2425 } 2426 2427 /* allocate the matrix space */ 2428 c->singlemalloc = PETSC_TRUE; 2429 len = (mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)); 2430 ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 2431 c->j = (PetscInt*)(c->a + nz*bs2); 2432 c->i = c->j + nz; 2433 ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 2434 if (mbs > 0) { 2435 ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 2436 if (cpvalues == MAT_COPY_VALUES) { 2437 ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 2438 } else { 2439 ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 2440 } 2441 } 2442 2443 PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 2444 c->sorted = a->sorted; 2445 c->roworiented = a->roworiented; 2446 c->nonew = a->nonew; 2447 2448 if (a->diag) { 2449 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 2450 PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt)); 2451 for (i=0; i<mbs; i++) { 2452 c->diag[i] = a->diag[i]; 2453 } 2454 } else c->diag = 0; 2455 c->nz = a->nz; 2456 c->maxnz = a->maxnz; 2457 c->solve_work = 0; 2458 c->mult_work = 0; 2459 C->preallocated = PETSC_TRUE; 2460 C->assembled = PETSC_TRUE; 2461 2462 c->compressedrow.use = a->compressedrow.use; 2463 c->compressedrow.nrows = a->compressedrow.nrows; 2464 c->compressedrow.checked = a->compressedrow.checked; 2465 if ( a->compressedrow.checked && a->compressedrow.use){ 2466 i = a->compressedrow.nrows; 2467 ierr = PetscMalloc((2*i+1)*sizeof(PetscInt),&c->compressedrow.i);CHKERRQ(ierr); 2468 c->compressedrow.rindex = c->compressedrow.i + i + 1; 2469 ierr = PetscMemcpy(c->compressedrow.i,a->compressedrow.i,(i+1)*sizeof(PetscInt));CHKERRQ(ierr); 2470 ierr = PetscMemcpy(c->compressedrow.rindex,a->compressedrow.rindex,i*sizeof(PetscInt));CHKERRQ(ierr); 2471 } else { 2472 c->compressedrow.use = PETSC_FALSE; 2473 c->compressedrow.i = PETSC_NULL; 2474 c->compressedrow.rindex = PETSC_NULL; 2475 } 2476 C->same_nonzero = A->same_nonzero; 2477 *B = C; 2478 ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 2479 PetscFunctionReturn(0); 2480 } 2481 2482 #undef __FUNCT__ 2483 #define __FUNCT__ "MatLoad_SeqBAIJ" 2484 PetscErrorCode MatLoad_SeqBAIJ(PetscViewer viewer,const MatType type,Mat *A) 2485 { 2486 Mat_SeqBAIJ *a; 2487 Mat B; 2488 PetscErrorCode ierr; 2489 PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 2490 PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 2491 PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 2492 PetscInt *masked,nmask,tmp,bs2,ishift; 2493 PetscMPIInt size; 2494 int fd; 2495 PetscScalar *aa; 2496 MPI_Comm comm = ((PetscObject)viewer)->comm; 2497 2498 PetscFunctionBegin; 2499 ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 2500 bs2 = bs*bs; 2501 2502 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 2503 if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 2504 ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 2505 ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 2506 if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 2507 M = header[1]; N = header[2]; nz = header[3]; 2508 2509 if (header[3] < 0) { 2510 SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ"); 2511 } 2512 2513 if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 2514 2515 /* 2516 This code adds extra rows to make sure the number of rows is 2517 divisible by the blocksize 2518 */ 2519 mbs = M/bs; 2520 extra_rows = bs - M + bs*(mbs); 2521 if (extra_rows == bs) extra_rows = 0; 2522 else mbs++; 2523 if (extra_rows) { 2524 PetscLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 2525 } 2526 2527 /* read in row lengths */ 2528 ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 2529 ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 2530 for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 2531 2532 /* read in column indices */ 2533 ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 2534 ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 2535 for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 2536 2537 /* loop over row lengths determining block row lengths */ 2538 ierr = PetscMalloc(mbs*sizeof(PetscInt),&browlengths);CHKERRQ(ierr); 2539 ierr = PetscMemzero(browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 2540 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 2541 ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 2542 masked = mask + mbs; 2543 rowcount = 0; nzcount = 0; 2544 for (i=0; i<mbs; i++) { 2545 nmask = 0; 2546 for (j=0; j<bs; j++) { 2547 kmax = rowlengths[rowcount]; 2548 for (k=0; k<kmax; k++) { 2549 tmp = jj[nzcount++]/bs; 2550 if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 2551 } 2552 rowcount++; 2553 } 2554 browlengths[i] += nmask; 2555 /* zero out the mask elements we set */ 2556 for (j=0; j<nmask; j++) mask[masked[j]] = 0; 2557 } 2558 2559 /* create our matrix */ 2560 ierr = MatCreate(comm,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows,&B); 2561 ierr = MatSetType(B,type);CHKERRQ(ierr); 2562 ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); 2563 a = (Mat_SeqBAIJ*)B->data; 2564 2565 /* set matrix "i" values */ 2566 a->i[0] = 0; 2567 for (i=1; i<= mbs; i++) { 2568 a->i[i] = a->i[i-1] + browlengths[i-1]; 2569 a->ilen[i-1] = browlengths[i-1]; 2570 } 2571 a->nz = 0; 2572 for (i=0; i<mbs; i++) a->nz += browlengths[i]; 2573 2574 /* read in nonzero values */ 2575 ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 2576 ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 2577 for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 2578 2579 /* set "a" and "j" values into matrix */ 2580 nzcount = 0; jcount = 0; 2581 for (i=0; i<mbs; i++) { 2582 nzcountb = nzcount; 2583 nmask = 0; 2584 for (j=0; j<bs; j++) { 2585 kmax = rowlengths[i*bs+j]; 2586 for (k=0; k<kmax; k++) { 2587 tmp = jj[nzcount++]/bs; 2588 if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 2589 } 2590 } 2591 /* sort the masked values */ 2592 ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 2593 2594 /* set "j" values into matrix */ 2595 maskcount = 1; 2596 for (j=0; j<nmask; j++) { 2597 a->j[jcount++] = masked[j]; 2598 mask[masked[j]] = maskcount++; 2599 } 2600 /* set "a" values into matrix */ 2601 ishift = bs2*a->i[i]; 2602 for (j=0; j<bs; j++) { 2603 kmax = rowlengths[i*bs+j]; 2604 for (k=0; k<kmax; k++) { 2605 tmp = jj[nzcountb]/bs ; 2606 block = mask[tmp] - 1; 2607 point = jj[nzcountb] - bs*tmp; 2608 idx = ishift + bs2*block + j + bs*point; 2609 a->a[idx] = (MatScalar)aa[nzcountb++]; 2610 } 2611 } 2612 /* zero out the mask elements we set */ 2613 for (j=0; j<nmask; j++) mask[masked[j]] = 0; 2614 } 2615 if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 2616 2617 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 2618 ierr = PetscFree(browlengths);CHKERRQ(ierr); 2619 ierr = PetscFree(aa);CHKERRQ(ierr); 2620 ierr = PetscFree(jj);CHKERRQ(ierr); 2621 ierr = PetscFree(mask);CHKERRQ(ierr); 2622 2623 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2624 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2625 ierr = MatView_Private(B);CHKERRQ(ierr); 2626 2627 *A = B; 2628 PetscFunctionReturn(0); 2629 } 2630 2631 #undef __FUNCT__ 2632 #define __FUNCT__ "MatCreateSeqBAIJ" 2633 /*@C 2634 MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 2635 compressed row) format. For good matrix assembly performance the 2636 user should preallocate the matrix storage by setting the parameter nz 2637 (or the array nnz). By setting these parameters accurately, performance 2638 during matrix assembly can be increased by more than a factor of 50. 2639 2640 Collective on MPI_Comm 2641 2642 Input Parameters: 2643 + comm - MPI communicator, set to PETSC_COMM_SELF 2644 . bs - size of block 2645 . m - number of rows 2646 . n - number of columns 2647 . nz - number of nonzero blocks per block row (same for all rows) 2648 - nnz - array containing the number of nonzero blocks in the various block rows 2649 (possibly different for each block row) or PETSC_NULL 2650 2651 Output Parameter: 2652 . A - the matrix 2653 2654 Options Database Keys: 2655 . -mat_no_unroll - uses code that does not unroll the loops in the 2656 block calculations (much slower) 2657 . -mat_block_size - size of the blocks to use 2658 2659 Level: intermediate 2660 2661 Notes: 2662 If the nnz parameter is given then the nz parameter is ignored 2663 2664 A nonzero block is any block that as 1 or more nonzeros in it 2665 2666 The block AIJ format is fully compatible with standard Fortran 77 2667 storage. That is, the stored row and column indices can begin at 2668 either one (as in Fortran) or zero. See the users' manual for details. 2669 2670 Specify the preallocated storage with either nz or nnz (not both). 2671 Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 2672 allocation. For additional details, see the users manual chapter on 2673 matrices. 2674 2675 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 2676 @*/ 2677 PetscErrorCode MatCreateSeqBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 2678 { 2679 PetscErrorCode ierr; 2680 2681 PetscFunctionBegin; 2682 ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 2683 ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr); 2684 ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 2685 PetscFunctionReturn(0); 2686 } 2687 2688 #undef __FUNCT__ 2689 #define __FUNCT__ "MatSeqBAIJSetPreallocation" 2690 /*@C 2691 MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros 2692 per row in the matrix. For good matrix assembly performance the 2693 user should preallocate the matrix storage by setting the parameter nz 2694 (or the array nnz). By setting these parameters accurately, performance 2695 during matrix assembly can be increased by more than a factor of 50. 2696 2697 Collective on MPI_Comm 2698 2699 Input Parameters: 2700 + A - the matrix 2701 . bs - size of block 2702 . nz - number of block nonzeros per block row (same for all rows) 2703 - nnz - array containing the number of block nonzeros in the various block rows 2704 (possibly different for each block row) or PETSC_NULL 2705 2706 Options Database Keys: 2707 . -mat_no_unroll - uses code that does not unroll the loops in the 2708 block calculations (much slower) 2709 . -mat_block_size - size of the blocks to use 2710 2711 Level: intermediate 2712 2713 Notes: 2714 If the nnz parameter is given then the nz parameter is ignored 2715 2716 The block AIJ format is fully compatible with standard Fortran 77 2717 storage. That is, the stored row and column indices can begin at 2718 either one (as in Fortran) or zero. See the users' manual for details. 2719 2720 Specify the preallocated storage with either nz or nnz (not both). 2721 Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 2722 allocation. For additional details, see the users manual chapter on 2723 matrices. 2724 2725 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 2726 @*/ 2727 PetscErrorCode MatSeqBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 2728 { 2729 PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 2730 2731 PetscFunctionBegin; 2732 ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 2733 if (f) { 2734 ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 2735 } 2736 PetscFunctionReturn(0); 2737 } 2738 2739