1 2 /* 3 Defines the basic matrix operations for the BAIJ (compressed row) 4 matrix storage format. 5 */ 6 #include <../src/mat/impls/baij/seq/baij.h> /*I "petscmat.h" I*/ 7 #include <petscblaslapack.h> 8 #include <petsc/private/kernels/blockinvert.h> 9 #include <petsc/private/kernels/blockmatmult.h> 10 11 #if defined(PETSC_HAVE_HYPRE) 12 PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*); 13 #endif 14 15 #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) 16 PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqBAIJMKL(Mat,MatType,MatReuse,Mat*); 17 #endif 18 PETSC_INTERN PetscErrorCode MatConvert_XAIJ_IS(Mat,MatType,MatReuse,Mat*); 19 20 PetscErrorCode MatGetColumnReductions_SeqBAIJ(Mat A,PetscInt type,PetscReal *reductions) 21 { 22 Mat_SeqBAIJ *a_aij = (Mat_SeqBAIJ*) A->data; 23 PetscInt m,n,i; 24 PetscInt ib,jb,bs = A->rmap->bs; 25 MatScalar *a_val = a_aij->a; 26 27 PetscFunctionBegin; 28 PetscCall(MatGetSize(A,&m,&n)); 29 for (i=0; i<n; i++) reductions[i] = 0.0; 30 if (type == NORM_2) { 31 for (i=a_aij->i[0]; i<a_aij->i[A->rmap->n/bs]; i++) { 32 for (jb=0; jb<bs; jb++) { 33 for (ib=0; ib<bs; ib++) { 34 reductions[A->cmap->rstart + a_aij->j[i] * bs + jb] += PetscAbsScalar(*a_val * *a_val); 35 a_val++; 36 } 37 } 38 } 39 } else if (type == NORM_1) { 40 for (i=a_aij->i[0]; i<a_aij->i[A->rmap->n/bs]; i++) { 41 for (jb=0; jb<bs; jb++) { 42 for (ib=0; ib<bs; ib++) { 43 reductions[A->cmap->rstart + a_aij->j[i] * bs + jb] += PetscAbsScalar(*a_val); 44 a_val++; 45 } 46 } 47 } 48 } else if (type == NORM_INFINITY) { 49 for (i=a_aij->i[0]; i<a_aij->i[A->rmap->n/bs]; i++) { 50 for (jb=0; jb<bs; jb++) { 51 for (ib=0; ib<bs; ib++) { 52 int col = A->cmap->rstart + a_aij->j[i] * bs + jb; 53 reductions[col] = PetscMax(PetscAbsScalar(*a_val), reductions[col]); 54 a_val++; 55 } 56 } 57 } 58 } else if (type == REDUCTION_SUM_REALPART || type == REDUCTION_MEAN_REALPART) { 59 for (i=a_aij->i[0]; i<a_aij->i[A->rmap->n/bs]; i++) { 60 for (jb=0; jb<bs; jb++) { 61 for (ib=0; ib<bs; ib++) { 62 reductions[A->cmap->rstart + a_aij->j[i] * bs + jb] += PetscRealPart(*a_val); 63 a_val++; 64 } 65 } 66 } 67 } else if (type == REDUCTION_SUM_IMAGINARYPART || type == REDUCTION_MEAN_IMAGINARYPART) { 68 for (i=a_aij->i[0]; i<a_aij->i[A->rmap->n/bs]; i++) { 69 for (jb=0; jb<bs; jb++) { 70 for (ib=0; ib<bs; ib++) { 71 reductions[A->cmap->rstart + a_aij->j[i] * bs + jb] += PetscImaginaryPart(*a_val); 72 a_val++; 73 } 74 } 75 } 76 } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown reduction type"); 77 if (type == NORM_2) { 78 for (i=0; i<n; i++) reductions[i] = PetscSqrtReal(reductions[i]); 79 } else if (type == REDUCTION_MEAN_REALPART || type == REDUCTION_MEAN_IMAGINARYPART) { 80 for (i=0; i<n; i++) reductions[i] /= m; 81 } 82 PetscFunctionReturn(0); 83 } 84 85 PetscErrorCode MatInvertBlockDiagonal_SeqBAIJ(Mat A,const PetscScalar **values) 86 { 87 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*) A->data; 88 PetscInt *diag_offset,i,bs = A->rmap->bs,mbs = a->mbs,ipvt[5],bs2 = bs*bs,*v_pivots; 89 MatScalar *v = a->a,*odiag,*diag,work[25],*v_work; 90 PetscReal shift = 0.0; 91 PetscBool allowzeropivot,zeropivotdetected=PETSC_FALSE; 92 93 PetscFunctionBegin; 94 allowzeropivot = PetscNot(A->erroriffailure); 95 96 if (a->idiagvalid) { 97 if (values) *values = a->idiag; 98 PetscFunctionReturn(0); 99 } 100 PetscCall(MatMarkDiagonal_SeqBAIJ(A)); 101 diag_offset = a->diag; 102 if (!a->idiag) { 103 PetscCall(PetscMalloc1(bs2*mbs,&a->idiag)); 104 PetscCall(PetscLogObjectMemory((PetscObject)A,bs2*mbs*sizeof(PetscScalar))); 105 } 106 diag = a->idiag; 107 if (values) *values = a->idiag; 108 /* factor and invert each block */ 109 switch (bs) { 110 case 1: 111 for (i=0; i<mbs; i++) { 112 odiag = v + 1*diag_offset[i]; 113 diag[0] = odiag[0]; 114 115 if (PetscAbsScalar(diag[0] + shift) < PETSC_MACHINE_EPSILON) { 116 if (allowzeropivot) { 117 A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 118 A->factorerror_zeropivot_value = PetscAbsScalar(diag[0]); 119 A->factorerror_zeropivot_row = i; 120 PetscCall(PetscInfo(A,"Zero pivot, row %" PetscInt_FMT "\n",i)); 121 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %" PetscInt_FMT " pivot value %g tolerance %g",i,(double)PetscAbsScalar(diag[0]),(double)PETSC_MACHINE_EPSILON); 122 } 123 124 diag[0] = (PetscScalar)1.0 / (diag[0] + shift); 125 diag += 1; 126 } 127 break; 128 case 2: 129 for (i=0; i<mbs; i++) { 130 odiag = v + 4*diag_offset[i]; 131 diag[0] = odiag[0]; diag[1] = odiag[1]; diag[2] = odiag[2]; diag[3] = odiag[3]; 132 PetscCall(PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected)); 133 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 134 diag += 4; 135 } 136 break; 137 case 3: 138 for (i=0; i<mbs; i++) { 139 odiag = v + 9*diag_offset[i]; 140 diag[0] = odiag[0]; diag[1] = odiag[1]; diag[2] = odiag[2]; diag[3] = odiag[3]; 141 diag[4] = odiag[4]; diag[5] = odiag[5]; diag[6] = odiag[6]; diag[7] = odiag[7]; 142 diag[8] = odiag[8]; 143 PetscCall(PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected)); 144 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 145 diag += 9; 146 } 147 break; 148 case 4: 149 for (i=0; i<mbs; i++) { 150 odiag = v + 16*diag_offset[i]; 151 PetscCall(PetscArraycpy(diag,odiag,16)); 152 PetscCall(PetscKernel_A_gets_inverse_A_4(diag,shift,allowzeropivot,&zeropivotdetected)); 153 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 154 diag += 16; 155 } 156 break; 157 case 5: 158 for (i=0; i<mbs; i++) { 159 odiag = v + 25*diag_offset[i]; 160 PetscCall(PetscArraycpy(diag,odiag,25)); 161 PetscCall(PetscKernel_A_gets_inverse_A_5(diag,ipvt,work,shift,allowzeropivot,&zeropivotdetected)); 162 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 163 diag += 25; 164 } 165 break; 166 case 6: 167 for (i=0; i<mbs; i++) { 168 odiag = v + 36*diag_offset[i]; 169 PetscCall(PetscArraycpy(diag,odiag,36)); 170 PetscCall(PetscKernel_A_gets_inverse_A_6(diag,shift,allowzeropivot,&zeropivotdetected)); 171 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 172 diag += 36; 173 } 174 break; 175 case 7: 176 for (i=0; i<mbs; i++) { 177 odiag = v + 49*diag_offset[i]; 178 PetscCall(PetscArraycpy(diag,odiag,49)); 179 PetscCall(PetscKernel_A_gets_inverse_A_7(diag,shift,allowzeropivot,&zeropivotdetected)); 180 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 181 diag += 49; 182 } 183 break; 184 default: 185 PetscCall(PetscMalloc2(bs,&v_work,bs,&v_pivots)); 186 for (i=0; i<mbs; i++) { 187 odiag = v + bs2*diag_offset[i]; 188 PetscCall(PetscArraycpy(diag,odiag,bs2)); 189 PetscCall(PetscKernel_A_gets_inverse_A(bs,diag,v_pivots,v_work,allowzeropivot,&zeropivotdetected)); 190 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 191 diag += bs2; 192 } 193 PetscCall(PetscFree2(v_work,v_pivots)); 194 } 195 a->idiagvalid = PETSC_TRUE; 196 PetscFunctionReturn(0); 197 } 198 199 PetscErrorCode MatSOR_SeqBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 200 { 201 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 202 PetscScalar *x,*work,*w,*workt,*t; 203 const MatScalar *v,*aa = a->a, *idiag; 204 const PetscScalar *b,*xb; 205 PetscScalar s[7], xw[7]={0}; /* avoid some compilers thinking xw is uninitialized */ 206 PetscInt m = a->mbs,i,i2,nz,bs = A->rmap->bs,bs2 = bs*bs,k,j,idx,it; 207 const PetscInt *diag,*ai = a->i,*aj = a->j,*vi; 208 209 PetscFunctionBegin; 210 its = its*lits; 211 PetscCheck(!(flag & SOR_EISENSTAT),PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 212 PetscCheck(its > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits); 213 PetscCheck(!fshift,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for diagonal shift"); 214 PetscCheck(omega == 1.0,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for non-trivial relaxation factor"); 215 PetscCheck(!(flag & SOR_APPLY_UPPER) && !(flag & SOR_APPLY_LOWER),PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for applying upper or lower triangular parts"); 216 217 if (!a->idiagvalid) PetscCall(MatInvertBlockDiagonal(A,NULL)); 218 219 if (!m) PetscFunctionReturn(0); 220 diag = a->diag; 221 idiag = a->idiag; 222 k = PetscMax(A->rmap->n,A->cmap->n); 223 if (!a->mult_work) { 224 PetscCall(PetscMalloc1(k+1,&a->mult_work)); 225 } 226 if (!a->sor_workt) { 227 PetscCall(PetscMalloc1(k,&a->sor_workt)); 228 } 229 if (!a->sor_work) { 230 PetscCall(PetscMalloc1(bs,&a->sor_work)); 231 } 232 work = a->mult_work; 233 t = a->sor_workt; 234 w = a->sor_work; 235 236 PetscCall(VecGetArray(xx,&x)); 237 PetscCall(VecGetArrayRead(bb,&b)); 238 239 if (flag & SOR_ZERO_INITIAL_GUESS) { 240 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 241 switch (bs) { 242 case 1: 243 PetscKernel_v_gets_A_times_w_1(x,idiag,b); 244 t[0] = b[0]; 245 i2 = 1; 246 idiag += 1; 247 for (i=1; i<m; i++) { 248 v = aa + ai[i]; 249 vi = aj + ai[i]; 250 nz = diag[i] - ai[i]; 251 s[0] = b[i2]; 252 for (j=0; j<nz; j++) { 253 xw[0] = x[vi[j]]; 254 PetscKernel_v_gets_v_minus_A_times_w_1(s,(v+j),xw); 255 } 256 t[i2] = s[0]; 257 PetscKernel_v_gets_A_times_w_1(xw,idiag,s); 258 x[i2] = xw[0]; 259 idiag += 1; 260 i2 += 1; 261 } 262 break; 263 case 2: 264 PetscKernel_v_gets_A_times_w_2(x,idiag,b); 265 t[0] = b[0]; t[1] = b[1]; 266 i2 = 2; 267 idiag += 4; 268 for (i=1; i<m; i++) { 269 v = aa + 4*ai[i]; 270 vi = aj + ai[i]; 271 nz = diag[i] - ai[i]; 272 s[0] = b[i2]; s[1] = b[i2+1]; 273 for (j=0; j<nz; j++) { 274 idx = 2*vi[j]; 275 it = 4*j; 276 xw[0] = x[idx]; xw[1] = x[1+idx]; 277 PetscKernel_v_gets_v_minus_A_times_w_2(s,(v+it),xw); 278 } 279 t[i2] = s[0]; t[i2+1] = s[1]; 280 PetscKernel_v_gets_A_times_w_2(xw,idiag,s); 281 x[i2] = xw[0]; x[i2+1] = xw[1]; 282 idiag += 4; 283 i2 += 2; 284 } 285 break; 286 case 3: 287 PetscKernel_v_gets_A_times_w_3(x,idiag,b); 288 t[0] = b[0]; t[1] = b[1]; t[2] = b[2]; 289 i2 = 3; 290 idiag += 9; 291 for (i=1; i<m; i++) { 292 v = aa + 9*ai[i]; 293 vi = aj + ai[i]; 294 nz = diag[i] - ai[i]; 295 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; 296 while (nz--) { 297 idx = 3*(*vi++); 298 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 299 PetscKernel_v_gets_v_minus_A_times_w_3(s,v,xw); 300 v += 9; 301 } 302 t[i2] = s[0]; t[i2+1] = s[1]; t[i2+2] = s[2]; 303 PetscKernel_v_gets_A_times_w_3(xw,idiag,s); 304 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; 305 idiag += 9; 306 i2 += 3; 307 } 308 break; 309 case 4: 310 PetscKernel_v_gets_A_times_w_4(x,idiag,b); 311 t[0] = b[0]; t[1] = b[1]; t[2] = b[2]; t[3] = b[3]; 312 i2 = 4; 313 idiag += 16; 314 for (i=1; i<m; i++) { 315 v = aa + 16*ai[i]; 316 vi = aj + ai[i]; 317 nz = diag[i] - ai[i]; 318 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; s[3] = b[i2+3]; 319 while (nz--) { 320 idx = 4*(*vi++); 321 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; xw[3] = x[3+idx]; 322 PetscKernel_v_gets_v_minus_A_times_w_4(s,v,xw); 323 v += 16; 324 } 325 t[i2] = s[0]; t[i2+1] = s[1]; t[i2+2] = s[2]; t[i2 + 3] = s[3]; 326 PetscKernel_v_gets_A_times_w_4(xw,idiag,s); 327 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; x[i2+3] = xw[3]; 328 idiag += 16; 329 i2 += 4; 330 } 331 break; 332 case 5: 333 PetscKernel_v_gets_A_times_w_5(x,idiag,b); 334 t[0] = b[0]; t[1] = b[1]; t[2] = b[2]; t[3] = b[3]; t[4] = b[4]; 335 i2 = 5; 336 idiag += 25; 337 for (i=1; i<m; i++) { 338 v = aa + 25*ai[i]; 339 vi = aj + ai[i]; 340 nz = diag[i] - ai[i]; 341 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; s[3] = b[i2+3]; s[4] = b[i2+4]; 342 while (nz--) { 343 idx = 5*(*vi++); 344 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; xw[3] = x[3+idx]; xw[4] = x[4+idx]; 345 PetscKernel_v_gets_v_minus_A_times_w_5(s,v,xw); 346 v += 25; 347 } 348 t[i2] = s[0]; t[i2+1] = s[1]; t[i2+2] = s[2]; t[i2+3] = s[3]; t[i2+4] = s[4]; 349 PetscKernel_v_gets_A_times_w_5(xw,idiag,s); 350 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; x[i2+3] = xw[3]; x[i2+4] = xw[4]; 351 idiag += 25; 352 i2 += 5; 353 } 354 break; 355 case 6: 356 PetscKernel_v_gets_A_times_w_6(x,idiag,b); 357 t[0] = b[0]; t[1] = b[1]; t[2] = b[2]; t[3] = b[3]; t[4] = b[4]; t[5] = b[5]; 358 i2 = 6; 359 idiag += 36; 360 for (i=1; i<m; i++) { 361 v = aa + 36*ai[i]; 362 vi = aj + ai[i]; 363 nz = diag[i] - ai[i]; 364 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; s[3] = b[i2+3]; s[4] = b[i2+4]; s[5] = b[i2+5]; 365 while (nz--) { 366 idx = 6*(*vi++); 367 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 368 xw[3] = x[3+idx]; xw[4] = x[4+idx]; xw[5] = x[5+idx]; 369 PetscKernel_v_gets_v_minus_A_times_w_6(s,v,xw); 370 v += 36; 371 } 372 t[i2] = s[0]; t[i2+1] = s[1]; t[i2+2] = s[2]; 373 t[i2+3] = s[3]; t[i2+4] = s[4]; t[i2+5] = s[5]; 374 PetscKernel_v_gets_A_times_w_6(xw,idiag,s); 375 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; x[i2+3] = xw[3]; x[i2+4] = xw[4]; x[i2+5] = xw[5]; 376 idiag += 36; 377 i2 += 6; 378 } 379 break; 380 case 7: 381 PetscKernel_v_gets_A_times_w_7(x,idiag,b); 382 t[0] = b[0]; t[1] = b[1]; t[2] = b[2]; 383 t[3] = b[3]; t[4] = b[4]; t[5] = b[5]; t[6] = b[6]; 384 i2 = 7; 385 idiag += 49; 386 for (i=1; i<m; i++) { 387 v = aa + 49*ai[i]; 388 vi = aj + ai[i]; 389 nz = diag[i] - ai[i]; 390 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; 391 s[3] = b[i2+3]; s[4] = b[i2+4]; s[5] = b[i2+5]; s[6] = b[i2+6]; 392 while (nz--) { 393 idx = 7*(*vi++); 394 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 395 xw[3] = x[3+idx]; xw[4] = x[4+idx]; xw[5] = x[5+idx]; xw[6] = x[6+idx]; 396 PetscKernel_v_gets_v_minus_A_times_w_7(s,v,xw); 397 v += 49; 398 } 399 t[i2] = s[0]; t[i2+1] = s[1]; t[i2+2] = s[2]; 400 t[i2+3] = s[3]; t[i2+4] = s[4]; t[i2+5] = s[5]; t[i2+6] = s[6]; 401 PetscKernel_v_gets_A_times_w_7(xw,idiag,s); 402 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; 403 x[i2+3] = xw[3]; x[i2+4] = xw[4]; x[i2+5] = xw[5]; x[i2+6] = xw[6]; 404 idiag += 49; 405 i2 += 7; 406 } 407 break; 408 default: 409 PetscKernel_w_gets_Ar_times_v(bs,bs,b,idiag,x); 410 PetscCall(PetscArraycpy(t,b,bs)); 411 i2 = bs; 412 idiag += bs2; 413 for (i=1; i<m; i++) { 414 v = aa + bs2*ai[i]; 415 vi = aj + ai[i]; 416 nz = diag[i] - ai[i]; 417 418 PetscCall(PetscArraycpy(w,b+i2,bs)); 419 /* copy all rows of x that are needed into contiguous space */ 420 workt = work; 421 for (j=0; j<nz; j++) { 422 PetscCall(PetscArraycpy(workt,x + bs*(*vi++),bs)); 423 workt += bs; 424 } 425 PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,v,work); 426 PetscCall(PetscArraycpy(t+i2,w,bs)); 427 PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,x+i2); 428 429 idiag += bs2; 430 i2 += bs; 431 } 432 break; 433 } 434 /* for logging purposes assume number of nonzero in lower half is 1/2 of total */ 435 PetscCall(PetscLogFlops(1.0*bs2*a->nz)); 436 xb = t; 437 } 438 else xb = b; 439 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 440 idiag = a->idiag+bs2*(a->mbs-1); 441 i2 = bs * (m-1); 442 switch (bs) { 443 case 1: 444 s[0] = xb[i2]; 445 PetscKernel_v_gets_A_times_w_1(xw,idiag,s); 446 x[i2] = xw[0]; 447 i2 -= 1; 448 for (i=m-2; i>=0; i--) { 449 v = aa + (diag[i]+1); 450 vi = aj + diag[i] + 1; 451 nz = ai[i+1] - diag[i] - 1; 452 s[0] = xb[i2]; 453 for (j=0; j<nz; j++) { 454 xw[0] = x[vi[j]]; 455 PetscKernel_v_gets_v_minus_A_times_w_1(s,(v+j),xw); 456 } 457 PetscKernel_v_gets_A_times_w_1(xw,idiag,s); 458 x[i2] = xw[0]; 459 idiag -= 1; 460 i2 -= 1; 461 } 462 break; 463 case 2: 464 s[0] = xb[i2]; s[1] = xb[i2+1]; 465 PetscKernel_v_gets_A_times_w_2(xw,idiag,s); 466 x[i2] = xw[0]; x[i2+1] = xw[1]; 467 i2 -= 2; 468 idiag -= 4; 469 for (i=m-2; i>=0; i--) { 470 v = aa + 4*(diag[i] + 1); 471 vi = aj + diag[i] + 1; 472 nz = ai[i+1] - diag[i] - 1; 473 s[0] = xb[i2]; s[1] = xb[i2+1]; 474 for (j=0; j<nz; j++) { 475 idx = 2*vi[j]; 476 it = 4*j; 477 xw[0] = x[idx]; xw[1] = x[1+idx]; 478 PetscKernel_v_gets_v_minus_A_times_w_2(s,(v+it),xw); 479 } 480 PetscKernel_v_gets_A_times_w_2(xw,idiag,s); 481 x[i2] = xw[0]; x[i2+1] = xw[1]; 482 idiag -= 4; 483 i2 -= 2; 484 } 485 break; 486 case 3: 487 s[0] = xb[i2]; s[1] = xb[i2+1]; s[2] = xb[i2+2]; 488 PetscKernel_v_gets_A_times_w_3(xw,idiag,s); 489 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; 490 i2 -= 3; 491 idiag -= 9; 492 for (i=m-2; i>=0; i--) { 493 v = aa + 9*(diag[i]+1); 494 vi = aj + diag[i] + 1; 495 nz = ai[i+1] - diag[i] - 1; 496 s[0] = xb[i2]; s[1] = xb[i2+1]; s[2] = xb[i2+2]; 497 while (nz--) { 498 idx = 3*(*vi++); 499 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 500 PetscKernel_v_gets_v_minus_A_times_w_3(s,v,xw); 501 v += 9; 502 } 503 PetscKernel_v_gets_A_times_w_3(xw,idiag,s); 504 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; 505 idiag -= 9; 506 i2 -= 3; 507 } 508 break; 509 case 4: 510 s[0] = xb[i2]; s[1] = xb[i2+1]; s[2] = xb[i2+2]; s[3] = xb[i2+3]; 511 PetscKernel_v_gets_A_times_w_4(xw,idiag,s); 512 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; x[i2+3] = xw[3]; 513 i2 -= 4; 514 idiag -= 16; 515 for (i=m-2; i>=0; i--) { 516 v = aa + 16*(diag[i]+1); 517 vi = aj + diag[i] + 1; 518 nz = ai[i+1] - diag[i] - 1; 519 s[0] = xb[i2]; s[1] = xb[i2+1]; s[2] = xb[i2+2]; s[3] = xb[i2+3]; 520 while (nz--) { 521 idx = 4*(*vi++); 522 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; xw[3] = x[3+idx]; 523 PetscKernel_v_gets_v_minus_A_times_w_4(s,v,xw); 524 v += 16; 525 } 526 PetscKernel_v_gets_A_times_w_4(xw,idiag,s); 527 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; x[i2+3] = xw[3]; 528 idiag -= 16; 529 i2 -= 4; 530 } 531 break; 532 case 5: 533 s[0] = xb[i2]; s[1] = xb[i2+1]; s[2] = xb[i2+2]; s[3] = xb[i2+3]; s[4] = xb[i2+4]; 534 PetscKernel_v_gets_A_times_w_5(xw,idiag,s); 535 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; x[i2+3] = xw[3]; x[i2+4] = xw[4]; 536 i2 -= 5; 537 idiag -= 25; 538 for (i=m-2; i>=0; i--) { 539 v = aa + 25*(diag[i]+1); 540 vi = aj + diag[i] + 1; 541 nz = ai[i+1] - diag[i] - 1; 542 s[0] = xb[i2]; s[1] = xb[i2+1]; s[2] = xb[i2+2]; s[3] = xb[i2+3]; s[4] = xb[i2+4]; 543 while (nz--) { 544 idx = 5*(*vi++); 545 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; xw[3] = x[3+idx]; xw[4] = x[4+idx]; 546 PetscKernel_v_gets_v_minus_A_times_w_5(s,v,xw); 547 v += 25; 548 } 549 PetscKernel_v_gets_A_times_w_5(xw,idiag,s); 550 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; x[i2+3] = xw[3]; x[i2+4] = xw[4]; 551 idiag -= 25; 552 i2 -= 5; 553 } 554 break; 555 case 6: 556 s[0] = xb[i2]; s[1] = xb[i2+1]; s[2] = xb[i2+2]; s[3] = xb[i2+3]; s[4] = xb[i2+4]; s[5] = xb[i2+5]; 557 PetscKernel_v_gets_A_times_w_6(xw,idiag,s); 558 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; x[i2+3] = xw[3]; x[i2+4] = xw[4]; x[i2+5] = xw[5]; 559 i2 -= 6; 560 idiag -= 36; 561 for (i=m-2; i>=0; i--) { 562 v = aa + 36*(diag[i]+1); 563 vi = aj + diag[i] + 1; 564 nz = ai[i+1] - diag[i] - 1; 565 s[0] = xb[i2]; s[1] = xb[i2+1]; s[2] = xb[i2+2]; s[3] = xb[i2+3]; s[4] = xb[i2+4]; s[5] = xb[i2+5]; 566 while (nz--) { 567 idx = 6*(*vi++); 568 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 569 xw[3] = x[3+idx]; xw[4] = x[4+idx]; xw[5] = x[5+idx]; 570 PetscKernel_v_gets_v_minus_A_times_w_6(s,v,xw); 571 v += 36; 572 } 573 PetscKernel_v_gets_A_times_w_6(xw,idiag,s); 574 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; x[i2+3] = xw[3]; x[i2+4] = xw[4]; x[i2+5] = xw[5]; 575 idiag -= 36; 576 i2 -= 6; 577 } 578 break; 579 case 7: 580 s[0] = xb[i2]; s[1] = xb[i2+1]; s[2] = xb[i2+2]; 581 s[3] = xb[i2+3]; s[4] = xb[i2+4]; s[5] = xb[i2+5]; s[6] = xb[i2+6]; 582 PetscKernel_v_gets_A_times_w_7(x,idiag,b); 583 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; 584 x[i2+3] = xw[3]; x[i2+4] = xw[4]; x[i2+5] = xw[5]; x[i2+6] = xw[6]; 585 i2 -= 7; 586 idiag -= 49; 587 for (i=m-2; i>=0; i--) { 588 v = aa + 49*(diag[i]+1); 589 vi = aj + diag[i] + 1; 590 nz = ai[i+1] - diag[i] - 1; 591 s[0] = xb[i2]; s[1] = xb[i2+1]; s[2] = xb[i2+2]; 592 s[3] = xb[i2+3]; s[4] = xb[i2+4]; s[5] = xb[i2+5]; s[6] = xb[i2+6]; 593 while (nz--) { 594 idx = 7*(*vi++); 595 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 596 xw[3] = x[3+idx]; xw[4] = x[4+idx]; xw[5] = x[5+idx]; xw[6] = x[6+idx]; 597 PetscKernel_v_gets_v_minus_A_times_w_7(s,v,xw); 598 v += 49; 599 } 600 PetscKernel_v_gets_A_times_w_7(xw,idiag,s); 601 x[i2] = xw[0]; x[i2+1] = xw[1]; x[i2+2] = xw[2]; 602 x[i2+3] = xw[3]; x[i2+4] = xw[4]; x[i2+5] = xw[5]; x[i2+6] = xw[6]; 603 idiag -= 49; 604 i2 -= 7; 605 } 606 break; 607 default: 608 PetscCall(PetscArraycpy(w,xb+i2,bs)); 609 PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,x+i2); 610 i2 -= bs; 611 idiag -= bs2; 612 for (i=m-2; i>=0; i--) { 613 v = aa + bs2*(diag[i]+1); 614 vi = aj + diag[i] + 1; 615 nz = ai[i+1] - diag[i] - 1; 616 617 PetscCall(PetscArraycpy(w,xb+i2,bs)); 618 /* copy all rows of x that are needed into contiguous space */ 619 workt = work; 620 for (j=0; j<nz; j++) { 621 PetscCall(PetscArraycpy(workt,x + bs*(*vi++),bs)); 622 workt += bs; 623 } 624 PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,v,work); 625 PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,x+i2); 626 627 idiag -= bs2; 628 i2 -= bs; 629 } 630 break; 631 } 632 PetscCall(PetscLogFlops(1.0*bs2*(a->nz))); 633 } 634 its--; 635 } 636 while (its--) { 637 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 638 idiag = a->idiag; 639 i2 = 0; 640 switch (bs) { 641 case 1: 642 for (i=0; i<m; i++) { 643 v = aa + ai[i]; 644 vi = aj + ai[i]; 645 nz = ai[i+1] - ai[i]; 646 s[0] = b[i2]; 647 for (j=0; j<nz; j++) { 648 xw[0] = x[vi[j]]; 649 PetscKernel_v_gets_v_minus_A_times_w_1(s,(v+j),xw); 650 } 651 PetscKernel_v_gets_A_times_w_1(xw,idiag,s); 652 x[i2] += xw[0]; 653 idiag += 1; 654 i2 += 1; 655 } 656 break; 657 case 2: 658 for (i=0; i<m; i++) { 659 v = aa + 4*ai[i]; 660 vi = aj + ai[i]; 661 nz = ai[i+1] - ai[i]; 662 s[0] = b[i2]; s[1] = b[i2+1]; 663 for (j=0; j<nz; j++) { 664 idx = 2*vi[j]; 665 it = 4*j; 666 xw[0] = x[idx]; xw[1] = x[1+idx]; 667 PetscKernel_v_gets_v_minus_A_times_w_2(s,(v+it),xw); 668 } 669 PetscKernel_v_gets_A_times_w_2(xw,idiag,s); 670 x[i2] += xw[0]; x[i2+1] += xw[1]; 671 idiag += 4; 672 i2 += 2; 673 } 674 break; 675 case 3: 676 for (i=0; i<m; i++) { 677 v = aa + 9*ai[i]; 678 vi = aj + ai[i]; 679 nz = ai[i+1] - ai[i]; 680 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; 681 while (nz--) { 682 idx = 3*(*vi++); 683 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 684 PetscKernel_v_gets_v_minus_A_times_w_3(s,v,xw); 685 v += 9; 686 } 687 PetscKernel_v_gets_A_times_w_3(xw,idiag,s); 688 x[i2] += xw[0]; x[i2+1] += xw[1]; x[i2+2] += xw[2]; 689 idiag += 9; 690 i2 += 3; 691 } 692 break; 693 case 4: 694 for (i=0; i<m; i++) { 695 v = aa + 16*ai[i]; 696 vi = aj + ai[i]; 697 nz = ai[i+1] - ai[i]; 698 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; s[3] = b[i2+3]; 699 while (nz--) { 700 idx = 4*(*vi++); 701 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; xw[3] = x[3+idx]; 702 PetscKernel_v_gets_v_minus_A_times_w_4(s,v,xw); 703 v += 16; 704 } 705 PetscKernel_v_gets_A_times_w_4(xw,idiag,s); 706 x[i2] += xw[0]; x[i2+1] += xw[1]; x[i2+2] += xw[2]; x[i2+3] += xw[3]; 707 idiag += 16; 708 i2 += 4; 709 } 710 break; 711 case 5: 712 for (i=0; i<m; i++) { 713 v = aa + 25*ai[i]; 714 vi = aj + ai[i]; 715 nz = ai[i+1] - ai[i]; 716 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; s[3] = b[i2+3]; s[4] = b[i2+4]; 717 while (nz--) { 718 idx = 5*(*vi++); 719 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; xw[3] = x[3+idx]; xw[4] = x[4+idx]; 720 PetscKernel_v_gets_v_minus_A_times_w_5(s,v,xw); 721 v += 25; 722 } 723 PetscKernel_v_gets_A_times_w_5(xw,idiag,s); 724 x[i2] += xw[0]; x[i2+1] += xw[1]; x[i2+2] += xw[2]; x[i2+3] += xw[3]; x[i2+4] += xw[4]; 725 idiag += 25; 726 i2 += 5; 727 } 728 break; 729 case 6: 730 for (i=0; i<m; i++) { 731 v = aa + 36*ai[i]; 732 vi = aj + ai[i]; 733 nz = ai[i+1] - ai[i]; 734 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; s[3] = b[i2+3]; s[4] = b[i2+4]; s[5] = b[i2+5]; 735 while (nz--) { 736 idx = 6*(*vi++); 737 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 738 xw[3] = x[3+idx]; xw[4] = x[4+idx]; xw[5] = x[5+idx]; 739 PetscKernel_v_gets_v_minus_A_times_w_6(s,v,xw); 740 v += 36; 741 } 742 PetscKernel_v_gets_A_times_w_6(xw,idiag,s); 743 x[i2] += xw[0]; x[i2+1] += xw[1]; x[i2+2] += xw[2]; 744 x[i2+3] += xw[3]; x[i2+4] += xw[4]; x[i2+5] += xw[5]; 745 idiag += 36; 746 i2 += 6; 747 } 748 break; 749 case 7: 750 for (i=0; i<m; i++) { 751 v = aa + 49*ai[i]; 752 vi = aj + ai[i]; 753 nz = ai[i+1] - ai[i]; 754 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; 755 s[3] = b[i2+3]; s[4] = b[i2+4]; s[5] = b[i2+5]; s[6] = b[i2+6]; 756 while (nz--) { 757 idx = 7*(*vi++); 758 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 759 xw[3] = x[3+idx]; xw[4] = x[4+idx]; xw[5] = x[5+idx]; xw[6] = x[6+idx]; 760 PetscKernel_v_gets_v_minus_A_times_w_7(s,v,xw); 761 v += 49; 762 } 763 PetscKernel_v_gets_A_times_w_7(xw,idiag,s); 764 x[i2] += xw[0]; x[i2+1] += xw[1]; x[i2+2] += xw[2]; 765 x[i2+3] += xw[3]; x[i2+4] += xw[4]; x[i2+5] += xw[5]; x[i2+6] += xw[6]; 766 idiag += 49; 767 i2 += 7; 768 } 769 break; 770 default: 771 for (i=0; i<m; i++) { 772 v = aa + bs2*ai[i]; 773 vi = aj + ai[i]; 774 nz = ai[i+1] - ai[i]; 775 776 PetscCall(PetscArraycpy(w,b+i2,bs)); 777 /* copy all rows of x that are needed into contiguous space */ 778 workt = work; 779 for (j=0; j<nz; j++) { 780 PetscCall(PetscArraycpy(workt,x + bs*(*vi++),bs)); 781 workt += bs; 782 } 783 PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,v,work); 784 PetscKernel_w_gets_w_plus_Ar_times_v(bs,bs,w,idiag,x+i2); 785 786 idiag += bs2; 787 i2 += bs; 788 } 789 break; 790 } 791 PetscCall(PetscLogFlops(2.0*bs2*a->nz)); 792 } 793 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 794 idiag = a->idiag+bs2*(a->mbs-1); 795 i2 = bs * (m-1); 796 switch (bs) { 797 case 1: 798 for (i=m-1; i>=0; i--) { 799 v = aa + ai[i]; 800 vi = aj + ai[i]; 801 nz = ai[i+1] - ai[i]; 802 s[0] = b[i2]; 803 for (j=0; j<nz; j++) { 804 xw[0] = x[vi[j]]; 805 PetscKernel_v_gets_v_minus_A_times_w_1(s,(v+j),xw); 806 } 807 PetscKernel_v_gets_A_times_w_1(xw,idiag,s); 808 x[i2] += xw[0]; 809 idiag -= 1; 810 i2 -= 1; 811 } 812 break; 813 case 2: 814 for (i=m-1; i>=0; i--) { 815 v = aa + 4*ai[i]; 816 vi = aj + ai[i]; 817 nz = ai[i+1] - ai[i]; 818 s[0] = b[i2]; s[1] = b[i2+1]; 819 for (j=0; j<nz; j++) { 820 idx = 2*vi[j]; 821 it = 4*j; 822 xw[0] = x[idx]; xw[1] = x[1+idx]; 823 PetscKernel_v_gets_v_minus_A_times_w_2(s,(v+it),xw); 824 } 825 PetscKernel_v_gets_A_times_w_2(xw,idiag,s); 826 x[i2] += xw[0]; x[i2+1] += xw[1]; 827 idiag -= 4; 828 i2 -= 2; 829 } 830 break; 831 case 3: 832 for (i=m-1; i>=0; i--) { 833 v = aa + 9*ai[i]; 834 vi = aj + ai[i]; 835 nz = ai[i+1] - ai[i]; 836 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; 837 while (nz--) { 838 idx = 3*(*vi++); 839 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 840 PetscKernel_v_gets_v_minus_A_times_w_3(s,v,xw); 841 v += 9; 842 } 843 PetscKernel_v_gets_A_times_w_3(xw,idiag,s); 844 x[i2] += xw[0]; x[i2+1] += xw[1]; x[i2+2] += xw[2]; 845 idiag -= 9; 846 i2 -= 3; 847 } 848 break; 849 case 4: 850 for (i=m-1; i>=0; i--) { 851 v = aa + 16*ai[i]; 852 vi = aj + ai[i]; 853 nz = ai[i+1] - ai[i]; 854 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; s[3] = b[i2+3]; 855 while (nz--) { 856 idx = 4*(*vi++); 857 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; xw[3] = x[3+idx]; 858 PetscKernel_v_gets_v_minus_A_times_w_4(s,v,xw); 859 v += 16; 860 } 861 PetscKernel_v_gets_A_times_w_4(xw,idiag,s); 862 x[i2] += xw[0]; x[i2+1] += xw[1]; x[i2+2] += xw[2]; x[i2+3] += xw[3]; 863 idiag -= 16; 864 i2 -= 4; 865 } 866 break; 867 case 5: 868 for (i=m-1; i>=0; i--) { 869 v = aa + 25*ai[i]; 870 vi = aj + ai[i]; 871 nz = ai[i+1] - ai[i]; 872 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; s[3] = b[i2+3]; s[4] = b[i2+4]; 873 while (nz--) { 874 idx = 5*(*vi++); 875 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; xw[3] = x[3+idx]; xw[4] = x[4+idx]; 876 PetscKernel_v_gets_v_minus_A_times_w_5(s,v,xw); 877 v += 25; 878 } 879 PetscKernel_v_gets_A_times_w_5(xw,idiag,s); 880 x[i2] += xw[0]; x[i2+1] += xw[1]; x[i2+2] += xw[2]; x[i2+3] += xw[3]; x[i2+4] += xw[4]; 881 idiag -= 25; 882 i2 -= 5; 883 } 884 break; 885 case 6: 886 for (i=m-1; i>=0; i--) { 887 v = aa + 36*ai[i]; 888 vi = aj + ai[i]; 889 nz = ai[i+1] - ai[i]; 890 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; s[3] = b[i2+3]; s[4] = b[i2+4]; s[5] = b[i2+5]; 891 while (nz--) { 892 idx = 6*(*vi++); 893 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 894 xw[3] = x[3+idx]; xw[4] = x[4+idx]; xw[5] = x[5+idx]; 895 PetscKernel_v_gets_v_minus_A_times_w_6(s,v,xw); 896 v += 36; 897 } 898 PetscKernel_v_gets_A_times_w_6(xw,idiag,s); 899 x[i2] += xw[0]; x[i2+1] += xw[1]; x[i2+2] += xw[2]; 900 x[i2+3] += xw[3]; x[i2+4] += xw[4]; x[i2+5] += xw[5]; 901 idiag -= 36; 902 i2 -= 6; 903 } 904 break; 905 case 7: 906 for (i=m-1; i>=0; i--) { 907 v = aa + 49*ai[i]; 908 vi = aj + ai[i]; 909 nz = ai[i+1] - ai[i]; 910 s[0] = b[i2]; s[1] = b[i2+1]; s[2] = b[i2+2]; 911 s[3] = b[i2+3]; s[4] = b[i2+4]; s[5] = b[i2+5]; s[6] = b[i2+6]; 912 while (nz--) { 913 idx = 7*(*vi++); 914 xw[0] = x[idx]; xw[1] = x[1+idx]; xw[2] = x[2+idx]; 915 xw[3] = x[3+idx]; xw[4] = x[4+idx]; xw[5] = x[5+idx]; xw[6] = x[6+idx]; 916 PetscKernel_v_gets_v_minus_A_times_w_7(s,v,xw); 917 v += 49; 918 } 919 PetscKernel_v_gets_A_times_w_7(xw,idiag,s); 920 x[i2] += xw[0]; x[i2+1] += xw[1]; x[i2+2] += xw[2]; 921 x[i2+3] += xw[3]; x[i2+4] += xw[4]; x[i2+5] += xw[5]; x[i2+6] += xw[6]; 922 idiag -= 49; 923 i2 -= 7; 924 } 925 break; 926 default: 927 for (i=m-1; i>=0; i--) { 928 v = aa + bs2*ai[i]; 929 vi = aj + ai[i]; 930 nz = ai[i+1] - ai[i]; 931 932 PetscCall(PetscArraycpy(w,b+i2,bs)); 933 /* copy all rows of x that are needed into contiguous space */ 934 workt = work; 935 for (j=0; j<nz; j++) { 936 PetscCall(PetscArraycpy(workt,x + bs*(*vi++),bs)); 937 workt += bs; 938 } 939 PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,v,work); 940 PetscKernel_w_gets_w_plus_Ar_times_v(bs,bs,w,idiag,x+i2); 941 942 idiag -= bs2; 943 i2 -= bs; 944 } 945 break; 946 } 947 PetscCall(PetscLogFlops(2.0*bs2*(a->nz))); 948 } 949 } 950 PetscCall(VecRestoreArray(xx,&x)); 951 PetscCall(VecRestoreArrayRead(bb,&b)); 952 PetscFunctionReturn(0); 953 } 954 955 /* 956 Special version for direct calls from Fortran (Used in PETSc-fun3d) 957 */ 958 #if defined(PETSC_HAVE_FORTRAN_CAPS) 959 #define matsetvaluesblocked4_ MATSETVALUESBLOCKED4 960 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 961 #define matsetvaluesblocked4_ matsetvaluesblocked4 962 #endif 963 964 PETSC_EXTERN void matsetvaluesblocked4_(Mat *AA,PetscInt *mm,const PetscInt im[],PetscInt *nn,const PetscInt in[],const PetscScalar v[]) 965 { 966 Mat A = *AA; 967 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 968 PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,N,m = *mm,n = *nn; 969 PetscInt *ai =a->i,*ailen=a->ilen; 970 PetscInt *aj =a->j,stepval,lastcol = -1; 971 const PetscScalar *value = v; 972 MatScalar *ap,*aa = a->a,*bap; 973 974 PetscFunctionBegin; 975 if (A->rmap->bs != 4) SETERRABORT(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Can only be called with a block size of 4"); 976 stepval = (n-1)*4; 977 for (k=0; k<m; k++) { /* loop over added rows */ 978 row = im[k]; 979 rp = aj + ai[row]; 980 ap = aa + 16*ai[row]; 981 nrow = ailen[row]; 982 low = 0; 983 high = nrow; 984 for (l=0; l<n; l++) { /* loop over added columns */ 985 col = in[l]; 986 if (col <= lastcol) low = 0; 987 else high = nrow; 988 lastcol = col; 989 value = v + k*(stepval+4 + l)*4; 990 while (high-low > 7) { 991 t = (low+high)/2; 992 if (rp[t] > col) high = t; 993 else low = t; 994 } 995 for (i=low; i<high; i++) { 996 if (rp[i] > col) break; 997 if (rp[i] == col) { 998 bap = ap + 16*i; 999 for (ii=0; ii<4; ii++,value+=stepval) { 1000 for (jj=ii; jj<16; jj+=4) { 1001 bap[jj] += *value++; 1002 } 1003 } 1004 goto noinsert2; 1005 } 1006 } 1007 N = nrow++ - 1; 1008 high++; /* added new column index thus must search to one higher than before */ 1009 /* shift up all the later entries in this row */ 1010 for (ii=N; ii>=i; ii--) { 1011 rp[ii+1] = rp[ii]; 1012 PetscCallVoid(PetscArraycpy(ap+16*(ii+1),ap+16*(ii),16)); 1013 } 1014 if (N >= i) { 1015 PetscCallVoid(PetscArrayzero(ap+16*i,16)); 1016 } 1017 rp[i] = col; 1018 bap = ap + 16*i; 1019 for (ii=0; ii<4; ii++,value+=stepval) { 1020 for (jj=ii; jj<16; jj+=4) { 1021 bap[jj] = *value++; 1022 } 1023 } 1024 noinsert2:; 1025 low = i; 1026 } 1027 ailen[row] = nrow; 1028 } 1029 PetscFunctionReturnVoid(); 1030 } 1031 1032 #if defined(PETSC_HAVE_FORTRAN_CAPS) 1033 #define matsetvalues4_ MATSETVALUES4 1034 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 1035 #define matsetvalues4_ matsetvalues4 1036 #endif 1037 1038 PETSC_EXTERN void matsetvalues4_(Mat *AA,PetscInt *mm,PetscInt *im,PetscInt *nn,PetscInt *in,PetscScalar *v) 1039 { 1040 Mat A = *AA; 1041 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1042 PetscInt *rp,k,low,high,t,row,nrow,i,col,l,N,n = *nn,m = *mm; 1043 PetscInt *ai=a->i,*ailen=a->ilen; 1044 PetscInt *aj=a->j,brow,bcol; 1045 PetscInt ridx,cidx,lastcol = -1; 1046 MatScalar *ap,value,*aa=a->a,*bap; 1047 1048 PetscFunctionBegin; 1049 for (k=0; k<m; k++) { /* loop over added rows */ 1050 row = im[k]; brow = row/4; 1051 rp = aj + ai[brow]; 1052 ap = aa + 16*ai[brow]; 1053 nrow = ailen[brow]; 1054 low = 0; 1055 high = nrow; 1056 for (l=0; l<n; l++) { /* loop over added columns */ 1057 col = in[l]; bcol = col/4; 1058 ridx = row % 4; cidx = col % 4; 1059 value = v[l + k*n]; 1060 if (col <= lastcol) low = 0; 1061 else high = nrow; 1062 lastcol = col; 1063 while (high-low > 7) { 1064 t = (low+high)/2; 1065 if (rp[t] > bcol) high = t; 1066 else low = t; 1067 } 1068 for (i=low; i<high; i++) { 1069 if (rp[i] > bcol) break; 1070 if (rp[i] == bcol) { 1071 bap = ap + 16*i + 4*cidx + ridx; 1072 *bap += value; 1073 goto noinsert1; 1074 } 1075 } 1076 N = nrow++ - 1; 1077 high++; /* added new column thus must search to one higher than before */ 1078 /* shift up all the later entries in this row */ 1079 PetscCallVoid(PetscArraymove(rp+i+1,rp+i,N-i+1)); 1080 PetscCallVoid(PetscArraymove(ap+16*i+16,ap+16*i,16*(N-i+1))); 1081 PetscCallVoid(PetscArrayzero(ap+16*i,16)); 1082 rp[i] = bcol; 1083 ap[16*i + 4*cidx + ridx] = value; 1084 noinsert1:; 1085 low = i; 1086 } 1087 ailen[brow] = nrow; 1088 } 1089 PetscFunctionReturnVoid(); 1090 } 1091 1092 /* 1093 Checks for missing diagonals 1094 */ 1095 PetscErrorCode MatMissingDiagonal_SeqBAIJ(Mat A,PetscBool *missing,PetscInt *d) 1096 { 1097 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1098 PetscInt *diag,*ii = a->i,i; 1099 1100 PetscFunctionBegin; 1101 PetscCall(MatMarkDiagonal_SeqBAIJ(A)); 1102 *missing = PETSC_FALSE; 1103 if (A->rmap->n > 0 && !ii) { 1104 *missing = PETSC_TRUE; 1105 if (d) *d = 0; 1106 PetscCall(PetscInfo(A,"Matrix has no entries therefore is missing diagonal\n")); 1107 } else { 1108 PetscInt n; 1109 n = PetscMin(a->mbs, a->nbs); 1110 diag = a->diag; 1111 for (i=0; i<n; i++) { 1112 if (diag[i] >= ii[i+1]) { 1113 *missing = PETSC_TRUE; 1114 if (d) *d = i; 1115 PetscCall(PetscInfo(A,"Matrix is missing block diagonal number %" PetscInt_FMT "\n",i)); 1116 break; 1117 } 1118 } 1119 } 1120 PetscFunctionReturn(0); 1121 } 1122 1123 PetscErrorCode MatMarkDiagonal_SeqBAIJ(Mat A) 1124 { 1125 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1126 PetscInt i,j,m = a->mbs; 1127 1128 PetscFunctionBegin; 1129 if (!a->diag) { 1130 PetscCall(PetscMalloc1(m,&a->diag)); 1131 PetscCall(PetscLogObjectMemory((PetscObject)A,m*sizeof(PetscInt))); 1132 a->free_diag = PETSC_TRUE; 1133 } 1134 for (i=0; i<m; i++) { 1135 a->diag[i] = a->i[i+1]; 1136 for (j=a->i[i]; j<a->i[i+1]; j++) { 1137 if (a->j[j] == i) { 1138 a->diag[i] = j; 1139 break; 1140 } 1141 } 1142 } 1143 PetscFunctionReturn(0); 1144 } 1145 1146 static PetscErrorCode MatGetRowIJ_SeqBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 1147 { 1148 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1149 PetscInt i,j,n = a->mbs,nz = a->i[n],*tia,*tja,bs = A->rmap->bs,k,l,cnt; 1150 PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 1151 1152 PetscFunctionBegin; 1153 *nn = n; 1154 if (!ia) PetscFunctionReturn(0); 1155 if (symmetric) { 1156 PetscCall(MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,PETSC_TRUE,0,0,&tia,&tja)); 1157 nz = tia[n]; 1158 } else { 1159 tia = a->i; tja = a->j; 1160 } 1161 1162 if (!blockcompressed && bs > 1) { 1163 (*nn) *= bs; 1164 /* malloc & create the natural set of indices */ 1165 PetscCall(PetscMalloc1((n+1)*bs,ia)); 1166 if (n) { 1167 (*ia)[0] = oshift; 1168 for (j=1; j<bs; j++) { 1169 (*ia)[j] = (tia[1]-tia[0])*bs+(*ia)[j-1]; 1170 } 1171 } 1172 1173 for (i=1; i<n; i++) { 1174 (*ia)[i*bs] = (tia[i]-tia[i-1])*bs + (*ia)[i*bs-1]; 1175 for (j=1; j<bs; j++) { 1176 (*ia)[i*bs+j] = (tia[i+1]-tia[i])*bs + (*ia)[i*bs+j-1]; 1177 } 1178 } 1179 if (n) { 1180 (*ia)[n*bs] = (tia[n]-tia[n-1])*bs + (*ia)[n*bs-1]; 1181 } 1182 1183 if (inja) { 1184 PetscCall(PetscMalloc1(nz*bs*bs,ja)); 1185 cnt = 0; 1186 for (i=0; i<n; i++) { 1187 for (j=0; j<bs; j++) { 1188 for (k=tia[i]; k<tia[i+1]; k++) { 1189 for (l=0; l<bs; l++) { 1190 (*ja)[cnt++] = bs*tja[k] + l; 1191 } 1192 } 1193 } 1194 } 1195 } 1196 1197 if (symmetric) { /* deallocate memory allocated in MatToSymmetricIJ_SeqAIJ() */ 1198 PetscCall(PetscFree(tia)); 1199 PetscCall(PetscFree(tja)); 1200 } 1201 } else if (oshift == 1) { 1202 if (symmetric) { 1203 nz = tia[A->rmap->n/bs]; 1204 /* add 1 to i and j indices */ 1205 for (i=0; i<A->rmap->n/bs+1; i++) tia[i] = tia[i] + 1; 1206 *ia = tia; 1207 if (ja) { 1208 for (i=0; i<nz; i++) tja[i] = tja[i] + 1; 1209 *ja = tja; 1210 } 1211 } else { 1212 nz = a->i[A->rmap->n/bs]; 1213 /* malloc space and add 1 to i and j indices */ 1214 PetscCall(PetscMalloc1(A->rmap->n/bs+1,ia)); 1215 for (i=0; i<A->rmap->n/bs+1; i++) (*ia)[i] = a->i[i] + 1; 1216 if (ja) { 1217 PetscCall(PetscMalloc1(nz,ja)); 1218 for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1; 1219 } 1220 } 1221 } else { 1222 *ia = tia; 1223 if (ja) *ja = tja; 1224 } 1225 PetscFunctionReturn(0); 1226 } 1227 1228 static PetscErrorCode MatRestoreRowIJ_SeqBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscBool *done) 1229 { 1230 PetscFunctionBegin; 1231 if (!ia) PetscFunctionReturn(0); 1232 if ((!blockcompressed && A->rmap->bs > 1) || (symmetric || oshift == 1)) { 1233 PetscCall(PetscFree(*ia)); 1234 if (ja) PetscCall(PetscFree(*ja)); 1235 } 1236 PetscFunctionReturn(0); 1237 } 1238 1239 PetscErrorCode MatDestroy_SeqBAIJ(Mat A) 1240 { 1241 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1242 1243 PetscFunctionBegin; 1244 #if defined(PETSC_USE_LOG) 1245 PetscLogObjectState((PetscObject)A,"Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT,A->rmap->N,A->cmap->n,a->nz); 1246 #endif 1247 PetscCall(MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i)); 1248 PetscCall(ISDestroy(&a->row)); 1249 PetscCall(ISDestroy(&a->col)); 1250 if (a->free_diag) PetscCall(PetscFree(a->diag)); 1251 PetscCall(PetscFree(a->idiag)); 1252 if (a->free_imax_ilen) PetscCall(PetscFree2(a->imax,a->ilen)); 1253 PetscCall(PetscFree(a->solve_work)); 1254 PetscCall(PetscFree(a->mult_work)); 1255 PetscCall(PetscFree(a->sor_workt)); 1256 PetscCall(PetscFree(a->sor_work)); 1257 PetscCall(ISDestroy(&a->icol)); 1258 PetscCall(PetscFree(a->saved_values)); 1259 PetscCall(PetscFree2(a->compressedrow.i,a->compressedrow.rindex)); 1260 1261 PetscCall(MatDestroy(&a->sbaijMat)); 1262 PetscCall(MatDestroy(&a->parent)); 1263 PetscCall(PetscFree(A->data)); 1264 1265 PetscCall(PetscObjectChangeTypeName((PetscObject)A,NULL)); 1266 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSeqBAIJGetArray_C",NULL)); 1267 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSeqBAIJRestoreArray_C",NULL)); 1268 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL)); 1269 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL)); 1270 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSeqBAIJSetColumnIndices_C",NULL)); 1271 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqbaij_seqaij_C",NULL)); 1272 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqbaij_seqsbaij_C",NULL)); 1273 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",NULL)); 1274 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSeqBAIJSetPreallocationCSR_C",NULL)); 1275 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqbaij_seqbstrm_C",NULL)); 1276 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatIsTranspose_C",NULL)); 1277 #if defined(PETSC_HAVE_HYPRE) 1278 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqbaij_hypre_C",NULL)); 1279 #endif 1280 PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqbaij_is_C",NULL)); 1281 PetscFunctionReturn(0); 1282 } 1283 1284 PetscErrorCode MatSetOption_SeqBAIJ(Mat A,MatOption op,PetscBool flg) 1285 { 1286 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1287 1288 PetscFunctionBegin; 1289 switch (op) { 1290 case MAT_ROW_ORIENTED: 1291 a->roworiented = flg; 1292 break; 1293 case MAT_KEEP_NONZERO_PATTERN: 1294 a->keepnonzeropattern = flg; 1295 break; 1296 case MAT_NEW_NONZERO_LOCATIONS: 1297 a->nonew = (flg ? 0 : 1); 1298 break; 1299 case MAT_NEW_NONZERO_LOCATION_ERR: 1300 a->nonew = (flg ? -1 : 0); 1301 break; 1302 case MAT_NEW_NONZERO_ALLOCATION_ERR: 1303 a->nonew = (flg ? -2 : 0); 1304 break; 1305 case MAT_UNUSED_NONZERO_LOCATION_ERR: 1306 a->nounused = (flg ? -1 : 0); 1307 break; 1308 case MAT_FORCE_DIAGONAL_ENTRIES: 1309 case MAT_IGNORE_OFF_PROC_ENTRIES: 1310 case MAT_USE_HASH_TABLE: 1311 case MAT_SORTED_FULL: 1312 PetscCall(PetscInfo(A,"Option %s ignored\n",MatOptions[op])); 1313 break; 1314 case MAT_SPD: 1315 case MAT_SYMMETRIC: 1316 case MAT_STRUCTURALLY_SYMMETRIC: 1317 case MAT_HERMITIAN: 1318 case MAT_SYMMETRY_ETERNAL: 1319 case MAT_SUBMAT_SINGLEIS: 1320 case MAT_STRUCTURE_ONLY: 1321 /* These options are handled directly by MatSetOption() */ 1322 break; 1323 default: 1324 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op); 1325 } 1326 PetscFunctionReturn(0); 1327 } 1328 1329 /* used for both SeqBAIJ and SeqSBAIJ matrices */ 1330 PetscErrorCode MatGetRow_SeqBAIJ_private(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v,PetscInt *ai,PetscInt *aj,PetscScalar *aa) 1331 { 1332 PetscInt itmp,i,j,k,M,bn,bp,*idx_i,bs,bs2; 1333 MatScalar *aa_i; 1334 PetscScalar *v_i; 1335 1336 PetscFunctionBegin; 1337 bs = A->rmap->bs; 1338 bs2 = bs*bs; 1339 PetscCheck(row >= 0 && row < A->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %" PetscInt_FMT " out of range", row); 1340 1341 bn = row/bs; /* Block number */ 1342 bp = row % bs; /* Block Position */ 1343 M = ai[bn+1] - ai[bn]; 1344 *nz = bs*M; 1345 1346 if (v) { 1347 *v = NULL; 1348 if (*nz) { 1349 PetscCall(PetscMalloc1(*nz,v)); 1350 for (i=0; i<M; i++) { /* for each block in the block row */ 1351 v_i = *v + i*bs; 1352 aa_i = aa + bs2*(ai[bn] + i); 1353 for (j=bp,k=0; j<bs2; j+=bs,k++) v_i[k] = aa_i[j]; 1354 } 1355 } 1356 } 1357 1358 if (idx) { 1359 *idx = NULL; 1360 if (*nz) { 1361 PetscCall(PetscMalloc1(*nz,idx)); 1362 for (i=0; i<M; i++) { /* for each block in the block row */ 1363 idx_i = *idx + i*bs; 1364 itmp = bs*aj[ai[bn] + i]; 1365 for (j=0; j<bs; j++) idx_i[j] = itmp++; 1366 } 1367 } 1368 } 1369 PetscFunctionReturn(0); 1370 } 1371 1372 PetscErrorCode MatGetRow_SeqBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 1373 { 1374 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1375 1376 PetscFunctionBegin; 1377 PetscCall(MatGetRow_SeqBAIJ_private(A,row,nz,idx,v,a->i,a->j,a->a)); 1378 PetscFunctionReturn(0); 1379 } 1380 1381 PetscErrorCode MatRestoreRow_SeqBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 1382 { 1383 PetscFunctionBegin; 1384 if (nz) *nz = 0; 1385 if (idx) PetscCall(PetscFree(*idx)); 1386 if (v) PetscCall(PetscFree(*v)); 1387 PetscFunctionReturn(0); 1388 } 1389 1390 PetscErrorCode MatTranspose_SeqBAIJ(Mat A,MatReuse reuse,Mat *B) 1391 { 1392 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*at; 1393 Mat C; 1394 PetscInt i,j,k,*aj=a->j,*ai=a->i,bs=A->rmap->bs,mbs=a->mbs,nbs=a->nbs,*atfill; 1395 PetscInt bs2=a->bs2,*ati,*atj,anzj,kr; 1396 MatScalar *ata,*aa=a->a; 1397 1398 PetscFunctionBegin; 1399 PetscCall(PetscCalloc1(1+nbs,&atfill)); 1400 if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) { 1401 for (i=0; i<ai[mbs]; i++) atfill[aj[i]] += 1; /* count num of non-zeros in row aj[i] */ 1402 1403 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&C)); 1404 PetscCall(MatSetSizes(C,A->cmap->n,A->rmap->N,A->cmap->n,A->rmap->N)); 1405 PetscCall(MatSetType(C,((PetscObject)A)->type_name)); 1406 PetscCall(MatSeqBAIJSetPreallocation(C,bs,0,atfill)); 1407 1408 at = (Mat_SeqBAIJ*)C->data; 1409 ati = at->i; 1410 for (i=0; i<nbs; i++) at->ilen[i] = at->imax[i] = ati[i+1] - ati[i]; 1411 } else { 1412 C = *B; 1413 at = (Mat_SeqBAIJ*)C->data; 1414 ati = at->i; 1415 } 1416 1417 atj = at->j; 1418 ata = at->a; 1419 1420 /* Copy ati into atfill so we have locations of the next free space in atj */ 1421 PetscCall(PetscArraycpy(atfill,ati,nbs)); 1422 1423 /* Walk through A row-wise and mark nonzero entries of A^T. */ 1424 for (i=0; i<mbs; i++) { 1425 anzj = ai[i+1] - ai[i]; 1426 for (j=0; j<anzj; j++) { 1427 atj[atfill[*aj]] = i; 1428 for (kr=0; kr<bs; kr++) { 1429 for (k=0; k<bs; k++) { 1430 ata[bs2*atfill[*aj]+k*bs+kr] = *aa++; 1431 } 1432 } 1433 atfill[*aj++] += 1; 1434 } 1435 } 1436 PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 1437 PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 1438 1439 /* Clean up temporary space and complete requests. */ 1440 PetscCall(PetscFree(atfill)); 1441 1442 if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) { 1443 PetscCall(MatSetBlockSizes(C,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs))); 1444 *B = C; 1445 } else { 1446 PetscCall(MatHeaderMerge(A,&C)); 1447 } 1448 PetscFunctionReturn(0); 1449 } 1450 1451 PetscErrorCode MatIsTranspose_SeqBAIJ(Mat A,Mat B,PetscReal tol,PetscBool *f) 1452 { 1453 Mat Btrans; 1454 1455 PetscFunctionBegin; 1456 *f = PETSC_FALSE; 1457 PetscCall(MatTranspose_SeqBAIJ(A,MAT_INITIAL_MATRIX,&Btrans)); 1458 PetscCall(MatEqual_SeqBAIJ(B,Btrans,f)); 1459 PetscCall(MatDestroy(&Btrans)); 1460 PetscFunctionReturn(0); 1461 } 1462 1463 /* Used for both SeqBAIJ and SeqSBAIJ matrices */ 1464 PetscErrorCode MatView_SeqBAIJ_Binary(Mat mat,PetscViewer viewer) 1465 { 1466 Mat_SeqBAIJ *A = (Mat_SeqBAIJ*)mat->data; 1467 PetscInt header[4],M,N,m,bs,nz,cnt,i,j,k,l; 1468 PetscInt *rowlens,*colidxs; 1469 PetscScalar *matvals; 1470 1471 PetscFunctionBegin; 1472 PetscCall(PetscViewerSetUp(viewer)); 1473 1474 M = mat->rmap->N; 1475 N = mat->cmap->N; 1476 m = mat->rmap->n; 1477 bs = mat->rmap->bs; 1478 nz = bs*bs*A->nz; 1479 1480 /* write matrix header */ 1481 header[0] = MAT_FILE_CLASSID; 1482 header[1] = M; header[2] = N; header[3] = nz; 1483 PetscCall(PetscViewerBinaryWrite(viewer,header,4,PETSC_INT)); 1484 1485 /* store row lengths */ 1486 PetscCall(PetscMalloc1(m,&rowlens)); 1487 for (cnt=0, i=0; i<A->mbs; i++) 1488 for (j=0; j<bs; j++) 1489 rowlens[cnt++] = bs*(A->i[i+1] - A->i[i]); 1490 PetscCall(PetscViewerBinaryWrite(viewer,rowlens,m,PETSC_INT)); 1491 PetscCall(PetscFree(rowlens)); 1492 1493 /* store column indices */ 1494 PetscCall(PetscMalloc1(nz,&colidxs)); 1495 for (cnt=0, i=0; i<A->mbs; i++) 1496 for (k=0; k<bs; k++) 1497 for (j=A->i[i]; j<A->i[i+1]; j++) 1498 for (l=0; l<bs; l++) 1499 colidxs[cnt++] = bs*A->j[j] + l; 1500 PetscCheck(cnt == nz,PETSC_COMM_SELF,PETSC_ERR_LIB,"Internal PETSc error: cnt = %" PetscInt_FMT " nz = %" PetscInt_FMT,cnt,nz); 1501 PetscCall(PetscViewerBinaryWrite(viewer,colidxs,nz,PETSC_INT)); 1502 PetscCall(PetscFree(colidxs)); 1503 1504 /* store nonzero values */ 1505 PetscCall(PetscMalloc1(nz,&matvals)); 1506 for (cnt=0, i=0; i<A->mbs; i++) 1507 for (k=0; k<bs; k++) 1508 for (j=A->i[i]; j<A->i[i+1]; j++) 1509 for (l=0; l<bs; l++) 1510 matvals[cnt++] = A->a[bs*(bs*j + l) + k]; 1511 PetscCheck(cnt == nz,PETSC_COMM_SELF,PETSC_ERR_LIB,"Internal PETSc error: cnt = %" PetscInt_FMT " nz = %" PetscInt_FMT,cnt,nz); 1512 PetscCall(PetscViewerBinaryWrite(viewer,matvals,nz,PETSC_SCALAR)); 1513 PetscCall(PetscFree(matvals)); 1514 1515 /* write block size option to the viewer's .info file */ 1516 PetscCall(MatView_Binary_BlockSizes(mat,viewer)); 1517 PetscFunctionReturn(0); 1518 } 1519 1520 static PetscErrorCode MatView_SeqBAIJ_ASCII_structonly(Mat A,PetscViewer viewer) 1521 { 1522 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1523 PetscInt i,bs = A->rmap->bs,k; 1524 1525 PetscFunctionBegin; 1526 PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); 1527 for (i=0; i<a->mbs; i++) { 1528 PetscCall(PetscViewerASCIIPrintf(viewer,"row %" PetscInt_FMT "-%" PetscInt_FMT ":",i*bs,i*bs+bs-1)); 1529 for (k=a->i[i]; k<a->i[i+1]; k++) { 1530 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT "-%" PetscInt_FMT ") ",bs*a->j[k],bs*a->j[k]+bs-1)); 1531 } 1532 PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); 1533 } 1534 PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); 1535 PetscFunctionReturn(0); 1536 } 1537 1538 static PetscErrorCode MatView_SeqBAIJ_ASCII(Mat A,PetscViewer viewer) 1539 { 1540 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1541 PetscInt i,j,bs = A->rmap->bs,k,l,bs2=a->bs2; 1542 PetscViewerFormat format; 1543 1544 PetscFunctionBegin; 1545 if (A->structure_only) { 1546 PetscCall(MatView_SeqBAIJ_ASCII_structonly(A,viewer)); 1547 PetscFunctionReturn(0); 1548 } 1549 1550 PetscCall(PetscViewerGetFormat(viewer,&format)); 1551 if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1552 PetscCall(PetscViewerASCIIPrintf(viewer," block size is %" PetscInt_FMT "\n",bs)); 1553 } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 1554 const char *matname; 1555 Mat aij; 1556 PetscCall(MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij)); 1557 PetscCall(PetscObjectGetName((PetscObject)A,&matname)); 1558 PetscCall(PetscObjectSetName((PetscObject)aij,matname)); 1559 PetscCall(MatView(aij,viewer)); 1560 PetscCall(MatDestroy(&aij)); 1561 } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 1562 PetscFunctionReturn(0); 1563 } else if (format == PETSC_VIEWER_ASCII_COMMON) { 1564 PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); 1565 for (i=0; i<a->mbs; i++) { 1566 for (j=0; j<bs; j++) { 1567 PetscCall(PetscViewerASCIIPrintf(viewer,"row %" PetscInt_FMT ":",i*bs+j)); 1568 for (k=a->i[i]; k<a->i[i+1]; k++) { 1569 for (l=0; l<bs; l++) { 1570 #if defined(PETSC_USE_COMPLEX) 1571 if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 1572 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g + %gi) ",bs*a->j[k]+l, 1573 (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]))); 1574 } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 1575 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g - %gi) ",bs*a->j[k]+l, 1576 (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]))); 1577 } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 1578 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]))); 1579 } 1580 #else 1581 if (a->a[bs2*k + l*bs + j] != 0.0) { 1582 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j])); 1583 } 1584 #endif 1585 } 1586 } 1587 PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); 1588 } 1589 } 1590 PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); 1591 } else { 1592 PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); 1593 for (i=0; i<a->mbs; i++) { 1594 for (j=0; j<bs; j++) { 1595 PetscCall(PetscViewerASCIIPrintf(viewer,"row %" PetscInt_FMT ":",i*bs+j)); 1596 for (k=a->i[i]; k<a->i[i+1]; k++) { 1597 for (l=0; l<bs; l++) { 1598 #if defined(PETSC_USE_COMPLEX) 1599 if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 1600 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g + %g i) ",bs*a->j[k]+l, 1601 (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]))); 1602 } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 1603 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g - %g i) ",bs*a->j[k]+l, 1604 (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]))); 1605 } else { 1606 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]))); 1607 } 1608 #else 1609 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j])); 1610 #endif 1611 } 1612 } 1613 PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); 1614 } 1615 } 1616 PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); 1617 } 1618 PetscCall(PetscViewerFlush(viewer)); 1619 PetscFunctionReturn(0); 1620 } 1621 1622 #include <petscdraw.h> 1623 static PetscErrorCode MatView_SeqBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 1624 { 1625 Mat A = (Mat) Aa; 1626 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 1627 PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap->bs,bs2=a->bs2; 1628 PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 1629 MatScalar *aa; 1630 PetscViewer viewer; 1631 PetscViewerFormat format; 1632 1633 PetscFunctionBegin; 1634 PetscCall(PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer)); 1635 PetscCall(PetscViewerGetFormat(viewer,&format)); 1636 PetscCall(PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr)); 1637 1638 /* loop over matrix elements drawing boxes */ 1639 1640 if (format != PETSC_VIEWER_DRAW_CONTOUR) { 1641 PetscDrawCollectiveBegin(draw); 1642 /* Blue for negative, Cyan for zero and Red for positive */ 1643 color = PETSC_DRAW_BLUE; 1644 for (i=0,row=0; i<mbs; i++,row+=bs) { 1645 for (j=a->i[i]; j<a->i[i+1]; j++) { 1646 y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0; 1647 x_l = a->j[j]*bs; x_r = x_l + 1.0; 1648 aa = a->a + j*bs2; 1649 for (k=0; k<bs; k++) { 1650 for (l=0; l<bs; l++) { 1651 if (PetscRealPart(*aa++) >= 0.) continue; 1652 PetscCall(PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color)); 1653 } 1654 } 1655 } 1656 } 1657 color = PETSC_DRAW_CYAN; 1658 for (i=0,row=0; i<mbs; i++,row+=bs) { 1659 for (j=a->i[i]; j<a->i[i+1]; j++) { 1660 y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0; 1661 x_l = a->j[j]*bs; x_r = x_l + 1.0; 1662 aa = a->a + j*bs2; 1663 for (k=0; k<bs; k++) { 1664 for (l=0; l<bs; l++) { 1665 if (PetscRealPart(*aa++) != 0.) continue; 1666 PetscCall(PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color)); 1667 } 1668 } 1669 } 1670 } 1671 color = PETSC_DRAW_RED; 1672 for (i=0,row=0; i<mbs; i++,row+=bs) { 1673 for (j=a->i[i]; j<a->i[i+1]; j++) { 1674 y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0; 1675 x_l = a->j[j]*bs; x_r = x_l + 1.0; 1676 aa = a->a + j*bs2; 1677 for (k=0; k<bs; k++) { 1678 for (l=0; l<bs; l++) { 1679 if (PetscRealPart(*aa++) <= 0.) continue; 1680 PetscCall(PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color)); 1681 } 1682 } 1683 } 1684 } 1685 PetscDrawCollectiveEnd(draw); 1686 } else { 1687 /* use contour shading to indicate magnitude of values */ 1688 /* first determine max of all nonzero values */ 1689 PetscReal minv = 0.0, maxv = 0.0; 1690 PetscDraw popup; 1691 1692 for (i=0; i<a->nz*a->bs2; i++) { 1693 if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]); 1694 } 1695 if (minv >= maxv) maxv = minv + PETSC_SMALL; 1696 PetscCall(PetscDrawGetPopup(draw,&popup)); 1697 PetscCall(PetscDrawScalePopup(popup,0.0,maxv)); 1698 1699 PetscDrawCollectiveBegin(draw); 1700 for (i=0,row=0; i<mbs; i++,row+=bs) { 1701 for (j=a->i[i]; j<a->i[i+1]; j++) { 1702 y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0; 1703 x_l = a->j[j]*bs; x_r = x_l + 1.0; 1704 aa = a->a + j*bs2; 1705 for (k=0; k<bs; k++) { 1706 for (l=0; l<bs; l++) { 1707 MatScalar v = *aa++; 1708 color = PetscDrawRealToColor(PetscAbsScalar(v),minv,maxv); 1709 PetscCall(PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color)); 1710 } 1711 } 1712 } 1713 } 1714 PetscDrawCollectiveEnd(draw); 1715 } 1716 PetscFunctionReturn(0); 1717 } 1718 1719 static PetscErrorCode MatView_SeqBAIJ_Draw(Mat A,PetscViewer viewer) 1720 { 1721 PetscReal xl,yl,xr,yr,w,h; 1722 PetscDraw draw; 1723 PetscBool isnull; 1724 1725 PetscFunctionBegin; 1726 PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 1727 PetscCall(PetscDrawIsNull(draw,&isnull)); 1728 if (isnull) PetscFunctionReturn(0); 1729 1730 xr = A->cmap->n; yr = A->rmap->N; h = yr/10.0; w = xr/10.0; 1731 xr += w; yr += h; xl = -w; yl = -h; 1732 PetscCall(PetscDrawSetCoordinates(draw,xl,yl,xr,yr)); 1733 PetscCall(PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer)); 1734 PetscCall(PetscDrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A)); 1735 PetscCall(PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL)); 1736 PetscCall(PetscDrawSave(draw)); 1737 PetscFunctionReturn(0); 1738 } 1739 1740 PetscErrorCode MatView_SeqBAIJ(Mat A,PetscViewer viewer) 1741 { 1742 PetscBool iascii,isbinary,isdraw; 1743 1744 PetscFunctionBegin; 1745 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 1746 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 1747 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 1748 if (iascii) { 1749 PetscCall(MatView_SeqBAIJ_ASCII(A,viewer)); 1750 } else if (isbinary) { 1751 PetscCall(MatView_SeqBAIJ_Binary(A,viewer)); 1752 } else if (isdraw) { 1753 PetscCall(MatView_SeqBAIJ_Draw(A,viewer)); 1754 } else { 1755 Mat B; 1756 PetscCall(MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B)); 1757 PetscCall(MatView(B,viewer)); 1758 PetscCall(MatDestroy(&B)); 1759 } 1760 PetscFunctionReturn(0); 1761 } 1762 1763 PetscErrorCode MatGetValues_SeqBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 1764 { 1765 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1766 PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 1767 PetscInt *ai = a->i,*ailen = a->ilen; 1768 PetscInt brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2; 1769 MatScalar *ap,*aa = a->a; 1770 1771 PetscFunctionBegin; 1772 for (k=0; k<m; k++) { /* loop over rows */ 1773 row = im[k]; brow = row/bs; 1774 if (row < 0) {v += n; continue;} /* negative row */ 1775 PetscCheck(row < A->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %" PetscInt_FMT " too large", row); 1776 rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 1777 nrow = ailen[brow]; 1778 for (l=0; l<n; l++) { /* loop over columns */ 1779 if (in[l] < 0) {v++; continue;} /* negative column */ 1780 PetscCheck(in[l] < A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column %" PetscInt_FMT " too large", in[l]); 1781 col = in[l]; 1782 bcol = col/bs; 1783 cidx = col%bs; 1784 ridx = row%bs; 1785 high = nrow; 1786 low = 0; /* assume unsorted */ 1787 while (high-low > 5) { 1788 t = (low+high)/2; 1789 if (rp[t] > bcol) high = t; 1790 else low = t; 1791 } 1792 for (i=low; i<high; i++) { 1793 if (rp[i] > bcol) break; 1794 if (rp[i] == bcol) { 1795 *v++ = ap[bs2*i+bs*cidx+ridx]; 1796 goto finished; 1797 } 1798 } 1799 *v++ = 0.0; 1800 finished:; 1801 } 1802 } 1803 PetscFunctionReturn(0); 1804 } 1805 1806 PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 1807 { 1808 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1809 PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 1810 PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 1811 PetscInt *aj =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval; 1812 PetscBool roworiented=a->roworiented; 1813 const PetscScalar *value = v; 1814 MatScalar *ap=NULL,*aa = a->a,*bap; 1815 1816 PetscFunctionBegin; 1817 if (roworiented) { 1818 stepval = (n-1)*bs; 1819 } else { 1820 stepval = (m-1)*bs; 1821 } 1822 for (k=0; k<m; k++) { /* loop over added rows */ 1823 row = im[k]; 1824 if (row < 0) continue; 1825 PetscCheck(row < a->mbs,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block row index too large %" PetscInt_FMT " max %" PetscInt_FMT,row,a->mbs-1); 1826 rp = aj + ai[row]; 1827 if (!A->structure_only) ap = aa + bs2*ai[row]; 1828 rmax = imax[row]; 1829 nrow = ailen[row]; 1830 low = 0; 1831 high = nrow; 1832 for (l=0; l<n; l++) { /* loop over added columns */ 1833 if (in[l] < 0) continue; 1834 PetscCheck(in[l] < a->nbs,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block column index too large %" PetscInt_FMT " max %" PetscInt_FMT,in[l],a->nbs-1); 1835 col = in[l]; 1836 if (!A->structure_only) { 1837 if (roworiented) { 1838 value = v + (k*(stepval+bs) + l)*bs; 1839 } else { 1840 value = v + (l*(stepval+bs) + k)*bs; 1841 } 1842 } 1843 if (col <= lastcol) low = 0; 1844 else high = nrow; 1845 lastcol = col; 1846 while (high-low > 7) { 1847 t = (low+high)/2; 1848 if (rp[t] > col) high = t; 1849 else low = t; 1850 } 1851 for (i=low; i<high; i++) { 1852 if (rp[i] > col) break; 1853 if (rp[i] == col) { 1854 if (A->structure_only) goto noinsert2; 1855 bap = ap + bs2*i; 1856 if (roworiented) { 1857 if (is == ADD_VALUES) { 1858 for (ii=0; ii<bs; ii++,value+=stepval) { 1859 for (jj=ii; jj<bs2; jj+=bs) { 1860 bap[jj] += *value++; 1861 } 1862 } 1863 } else { 1864 for (ii=0; ii<bs; ii++,value+=stepval) { 1865 for (jj=ii; jj<bs2; jj+=bs) { 1866 bap[jj] = *value++; 1867 } 1868 } 1869 } 1870 } else { 1871 if (is == ADD_VALUES) { 1872 for (ii=0; ii<bs; ii++,value+=bs+stepval) { 1873 for (jj=0; jj<bs; jj++) { 1874 bap[jj] += value[jj]; 1875 } 1876 bap += bs; 1877 } 1878 } else { 1879 for (ii=0; ii<bs; ii++,value+=bs+stepval) { 1880 for (jj=0; jj<bs; jj++) { 1881 bap[jj] = value[jj]; 1882 } 1883 bap += bs; 1884 } 1885 } 1886 } 1887 goto noinsert2; 1888 } 1889 } 1890 if (nonew == 1) goto noinsert2; 1891 PetscCheck(nonew != -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new blocked index new nonzero block (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col); 1892 if (A->structure_only) { 1893 MatSeqXAIJReallocateAIJ_structure_only(A,a->mbs,bs2,nrow,row,col,rmax,ai,aj,rp,imax,nonew,MatScalar); 1894 } else { 1895 MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 1896 } 1897 N = nrow++ - 1; high++; 1898 /* shift up all the later entries in this row */ 1899 PetscCall(PetscArraymove(rp+i+1,rp+i,N-i+1)); 1900 rp[i] = col; 1901 if (!A->structure_only) { 1902 PetscCall(PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1))); 1903 bap = ap + bs2*i; 1904 if (roworiented) { 1905 for (ii=0; ii<bs; ii++,value+=stepval) { 1906 for (jj=ii; jj<bs2; jj+=bs) { 1907 bap[jj] = *value++; 1908 } 1909 } 1910 } else { 1911 for (ii=0; ii<bs; ii++,value+=stepval) { 1912 for (jj=0; jj<bs; jj++) { 1913 *bap++ = *value++; 1914 } 1915 } 1916 } 1917 } 1918 noinsert2:; 1919 low = i; 1920 } 1921 ailen[row] = nrow; 1922 } 1923 PetscFunctionReturn(0); 1924 } 1925 1926 PetscErrorCode MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 1927 { 1928 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1929 PetscInt fshift = 0,i,*ai = a->i,*aj = a->j,*imax = a->imax; 1930 PetscInt m = A->rmap->N,*ip,N,*ailen = a->ilen; 1931 PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 1932 MatScalar *aa = a->a,*ap; 1933 PetscReal ratio=0.6; 1934 1935 PetscFunctionBegin; 1936 if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 1937 1938 if (m) rmax = ailen[0]; 1939 for (i=1; i<mbs; i++) { 1940 /* move each row back by the amount of empty slots (fshift) before it*/ 1941 fshift += imax[i-1] - ailen[i-1]; 1942 rmax = PetscMax(rmax,ailen[i]); 1943 if (fshift) { 1944 ip = aj + ai[i]; 1945 ap = aa + bs2*ai[i]; 1946 N = ailen[i]; 1947 PetscCall(PetscArraymove(ip-fshift,ip,N)); 1948 if (!A->structure_only) { 1949 PetscCall(PetscArraymove(ap-bs2*fshift,ap,bs2*N)); 1950 } 1951 } 1952 ai[i] = ai[i-1] + ailen[i-1]; 1953 } 1954 if (mbs) { 1955 fshift += imax[mbs-1] - ailen[mbs-1]; 1956 ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 1957 } 1958 1959 /* reset ilen and imax for each row */ 1960 a->nonzerorowcnt = 0; 1961 if (A->structure_only) { 1962 PetscCall(PetscFree2(a->imax,a->ilen)); 1963 } else { /* !A->structure_only */ 1964 for (i=0; i<mbs; i++) { 1965 ailen[i] = imax[i] = ai[i+1] - ai[i]; 1966 a->nonzerorowcnt += ((ai[i+1] - ai[i]) > 0); 1967 } 1968 } 1969 a->nz = ai[mbs]; 1970 1971 /* diagonals may have moved, so kill the diagonal pointers */ 1972 a->idiagvalid = PETSC_FALSE; 1973 if (fshift && a->diag) { 1974 PetscCall(PetscFree(a->diag)); 1975 PetscCall(PetscLogObjectMemory((PetscObject)A,-(mbs+1)*sizeof(PetscInt))); 1976 a->diag = NULL; 1977 } 1978 if (fshift) PetscCheck(a->nounused != -1,PETSC_COMM_SELF,PETSC_ERR_PLIB, "Unused space detected in matrix: %" PetscInt_FMT " X %" PetscInt_FMT " block size %" PetscInt_FMT ", %" PetscInt_FMT " unneeded", m, A->cmap->n, A->rmap->bs, fshift*bs2); 1979 PetscCall(PetscInfo(A,"Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT ", block size %" PetscInt_FMT "; storage space: %" PetscInt_FMT " unneeded, %" PetscInt_FMT " used\n",m,A->cmap->n,A->rmap->bs,fshift*bs2,a->nz*bs2)); 1980 PetscCall(PetscInfo(A,"Number of mallocs during MatSetValues is %" PetscInt_FMT "\n",a->reallocs)); 1981 PetscCall(PetscInfo(A,"Most nonzeros blocks in any row is %" PetscInt_FMT "\n",rmax)); 1982 1983 A->info.mallocs += a->reallocs; 1984 a->reallocs = 0; 1985 A->info.nz_unneeded = (PetscReal)fshift*bs2; 1986 a->rmax = rmax; 1987 1988 if (!A->structure_only) { 1989 PetscCall(MatCheckCompressedRow(A,a->nonzerorowcnt,&a->compressedrow,a->i,mbs,ratio)); 1990 } 1991 PetscFunctionReturn(0); 1992 } 1993 1994 /* 1995 This function returns an array of flags which indicate the locations of contiguous 1996 blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 1997 then the resulting sizes = [3,1,1,3,1] corresponding to sets [(0,1,2),(3),(5),(6,7,8),(9)] 1998 Assume: sizes should be long enough to hold all the values. 1999 */ 2000 static PetscErrorCode MatZeroRows_SeqBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 2001 { 2002 PetscInt i,j,k,row; 2003 PetscBool flg; 2004 2005 PetscFunctionBegin; 2006 for (i=0,j=0; i<n; j++) { 2007 row = idx[i]; 2008 if (row%bs!=0) { /* Not the beginning of a block */ 2009 sizes[j] = 1; 2010 i++; 2011 } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 2012 sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 2013 i++; 2014 } else { /* Begining of the block, so check if the complete block exists */ 2015 flg = PETSC_TRUE; 2016 for (k=1; k<bs; k++) { 2017 if (row+k != idx[i+k]) { /* break in the block */ 2018 flg = PETSC_FALSE; 2019 break; 2020 } 2021 } 2022 if (flg) { /* No break in the bs */ 2023 sizes[j] = bs; 2024 i += bs; 2025 } else { 2026 sizes[j] = 1; 2027 i++; 2028 } 2029 } 2030 } 2031 *bs_max = j; 2032 PetscFunctionReturn(0); 2033 } 2034 2035 PetscErrorCode MatZeroRows_SeqBAIJ(Mat A,PetscInt is_n,const PetscInt is_idx[],PetscScalar diag,Vec x, Vec b) 2036 { 2037 Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 2038 PetscInt i,j,k,count,*rows; 2039 PetscInt bs=A->rmap->bs,bs2=baij->bs2,*sizes,row,bs_max; 2040 PetscScalar zero = 0.0; 2041 MatScalar *aa; 2042 const PetscScalar *xx; 2043 PetscScalar *bb; 2044 2045 PetscFunctionBegin; 2046 /* fix right hand side if needed */ 2047 if (x && b) { 2048 PetscCall(VecGetArrayRead(x,&xx)); 2049 PetscCall(VecGetArray(b,&bb)); 2050 for (i=0; i<is_n; i++) { 2051 bb[is_idx[i]] = diag*xx[is_idx[i]]; 2052 } 2053 PetscCall(VecRestoreArrayRead(x,&xx)); 2054 PetscCall(VecRestoreArray(b,&bb)); 2055 } 2056 2057 /* Make a copy of the IS and sort it */ 2058 /* allocate memory for rows,sizes */ 2059 PetscCall(PetscMalloc2(is_n,&rows,2*is_n,&sizes)); 2060 2061 /* copy IS values to rows, and sort them */ 2062 for (i=0; i<is_n; i++) rows[i] = is_idx[i]; 2063 PetscCall(PetscSortInt(is_n,rows)); 2064 2065 if (baij->keepnonzeropattern) { 2066 for (i=0; i<is_n; i++) sizes[i] = 1; 2067 bs_max = is_n; 2068 } else { 2069 PetscCall(MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max)); 2070 A->nonzerostate++; 2071 } 2072 2073 for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 2074 row = rows[j]; 2075 PetscCheck(row >= 0 && row <= A->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %" PetscInt_FMT " out of range",row); 2076 count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 2077 aa = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs); 2078 if (sizes[i] == bs && !baij->keepnonzeropattern) { 2079 if (diag != (PetscScalar)0.0) { 2080 if (baij->ilen[row/bs] > 0) { 2081 baij->ilen[row/bs] = 1; 2082 baij->j[baij->i[row/bs]] = row/bs; 2083 2084 PetscCall(PetscArrayzero(aa,count*bs)); 2085 } 2086 /* Now insert all the diagonal values for this bs */ 2087 for (k=0; k<bs; k++) { 2088 PetscCall((*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,&diag,INSERT_VALUES)); 2089 } 2090 } else { /* (diag == 0.0) */ 2091 baij->ilen[row/bs] = 0; 2092 } /* end (diag == 0.0) */ 2093 } else { /* (sizes[i] != bs) */ 2094 PetscAssert(sizes[i] == 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal Error. Value should be 1"); 2095 for (k=0; k<count; k++) { 2096 aa[0] = zero; 2097 aa += bs; 2098 } 2099 if (diag != (PetscScalar)0.0) { 2100 PetscCall((*A->ops->setvalues)(A,1,rows+j,1,rows+j,&diag,INSERT_VALUES)); 2101 } 2102 } 2103 } 2104 2105 PetscCall(PetscFree2(rows,sizes)); 2106 PetscCall(MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY)); 2107 PetscFunctionReturn(0); 2108 } 2109 2110 PetscErrorCode MatZeroRowsColumns_SeqBAIJ(Mat A,PetscInt is_n,const PetscInt is_idx[],PetscScalar diag,Vec x, Vec b) 2111 { 2112 Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 2113 PetscInt i,j,k,count; 2114 PetscInt bs =A->rmap->bs,bs2=baij->bs2,row,col; 2115 PetscScalar zero = 0.0; 2116 MatScalar *aa; 2117 const PetscScalar *xx; 2118 PetscScalar *bb; 2119 PetscBool *zeroed,vecs = PETSC_FALSE; 2120 2121 PetscFunctionBegin; 2122 /* fix right hand side if needed */ 2123 if (x && b) { 2124 PetscCall(VecGetArrayRead(x,&xx)); 2125 PetscCall(VecGetArray(b,&bb)); 2126 vecs = PETSC_TRUE; 2127 } 2128 2129 /* zero the columns */ 2130 PetscCall(PetscCalloc1(A->rmap->n,&zeroed)); 2131 for (i=0; i<is_n; i++) { 2132 PetscCheck(is_idx[i] >= 0 && is_idx[i] < A->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %" PetscInt_FMT " out of range",is_idx[i]); 2133 zeroed[is_idx[i]] = PETSC_TRUE; 2134 } 2135 for (i=0; i<A->rmap->N; i++) { 2136 if (!zeroed[i]) { 2137 row = i/bs; 2138 for (j=baij->i[row]; j<baij->i[row+1]; j++) { 2139 for (k=0; k<bs; k++) { 2140 col = bs*baij->j[j] + k; 2141 if (zeroed[col]) { 2142 aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k; 2143 if (vecs) bb[i] -= aa[0]*xx[col]; 2144 aa[0] = 0.0; 2145 } 2146 } 2147 } 2148 } else if (vecs) bb[i] = diag*xx[i]; 2149 } 2150 PetscCall(PetscFree(zeroed)); 2151 if (vecs) { 2152 PetscCall(VecRestoreArrayRead(x,&xx)); 2153 PetscCall(VecRestoreArray(b,&bb)); 2154 } 2155 2156 /* zero the rows */ 2157 for (i=0; i<is_n; i++) { 2158 row = is_idx[i]; 2159 count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 2160 aa = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs); 2161 for (k=0; k<count; k++) { 2162 aa[0] = zero; 2163 aa += bs; 2164 } 2165 if (diag != (PetscScalar)0.0) { 2166 PetscCall((*A->ops->setvalues)(A,1,&row,1,&row,&diag,INSERT_VALUES)); 2167 } 2168 } 2169 PetscCall(MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY)); 2170 PetscFunctionReturn(0); 2171 } 2172 2173 PetscErrorCode MatSetValues_SeqBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 2174 { 2175 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2176 PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 2177 PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 2178 PetscInt *aj =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol; 2179 PetscInt ridx,cidx,bs2=a->bs2; 2180 PetscBool roworiented=a->roworiented; 2181 MatScalar *ap=NULL,value=0.0,*aa=a->a,*bap; 2182 2183 PetscFunctionBegin; 2184 for (k=0; k<m; k++) { /* loop over added rows */ 2185 row = im[k]; 2186 brow = row/bs; 2187 if (row < 0) continue; 2188 PetscCheck(row < A->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT,row,A->rmap->N-1); 2189 rp = aj + ai[brow]; 2190 if (!A->structure_only) ap = aa + bs2*ai[brow]; 2191 rmax = imax[brow]; 2192 nrow = ailen[brow]; 2193 low = 0; 2194 high = nrow; 2195 for (l=0; l<n; l++) { /* loop over added columns */ 2196 if (in[l] < 0) continue; 2197 PetscCheck(in[l] < A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT,in[l],A->cmap->n-1); 2198 col = in[l]; bcol = col/bs; 2199 ridx = row % bs; cidx = col % bs; 2200 if (!A->structure_only) { 2201 if (roworiented) { 2202 value = v[l + k*n]; 2203 } else { 2204 value = v[k + l*m]; 2205 } 2206 } 2207 if (col <= lastcol) low = 0; else high = nrow; 2208 lastcol = col; 2209 while (high-low > 7) { 2210 t = (low+high)/2; 2211 if (rp[t] > bcol) high = t; 2212 else low = t; 2213 } 2214 for (i=low; i<high; i++) { 2215 if (rp[i] > bcol) break; 2216 if (rp[i] == bcol) { 2217 bap = ap + bs2*i + bs*cidx + ridx; 2218 if (!A->structure_only) { 2219 if (is == ADD_VALUES) *bap += value; 2220 else *bap = value; 2221 } 2222 goto noinsert1; 2223 } 2224 } 2225 if (nonew == 1) goto noinsert1; 2226 PetscCheck(nonew != -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col); 2227 if (A->structure_only) { 2228 MatSeqXAIJReallocateAIJ_structure_only(A,a->mbs,bs2,nrow,brow,bcol,rmax,ai,aj,rp,imax,nonew,MatScalar); 2229 } else { 2230 MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 2231 } 2232 N = nrow++ - 1; high++; 2233 /* shift up all the later entries in this row */ 2234 PetscCall(PetscArraymove(rp+i+1,rp+i,N-i+1)); 2235 rp[i] = bcol; 2236 if (!A->structure_only) { 2237 PetscCall(PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1))); 2238 PetscCall(PetscArrayzero(ap+bs2*i,bs2)); 2239 ap[bs2*i + bs*cidx + ridx] = value; 2240 } 2241 a->nz++; 2242 A->nonzerostate++; 2243 noinsert1:; 2244 low = i; 2245 } 2246 ailen[brow] = nrow; 2247 } 2248 PetscFunctionReturn(0); 2249 } 2250 2251 PetscErrorCode MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,const MatFactorInfo *info) 2252 { 2253 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 2254 Mat outA; 2255 PetscBool row_identity,col_identity; 2256 2257 PetscFunctionBegin; 2258 PetscCheck(info->levels == 0,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 2259 PetscCall(ISIdentity(row,&row_identity)); 2260 PetscCall(ISIdentity(col,&col_identity)); 2261 PetscCheck(row_identity && col_identity,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Row and column permutations must be identity for in-place ILU"); 2262 2263 outA = inA; 2264 inA->factortype = MAT_FACTOR_LU; 2265 PetscCall(PetscFree(inA->solvertype)); 2266 PetscCall(PetscStrallocpy(MATSOLVERPETSC,&inA->solvertype)); 2267 2268 PetscCall(MatMarkDiagonal_SeqBAIJ(inA)); 2269 2270 PetscCall(PetscObjectReference((PetscObject)row)); 2271 PetscCall(ISDestroy(&a->row)); 2272 a->row = row; 2273 PetscCall(PetscObjectReference((PetscObject)col)); 2274 PetscCall(ISDestroy(&a->col)); 2275 a->col = col; 2276 2277 /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 2278 PetscCall(ISDestroy(&a->icol)); 2279 PetscCall(ISInvertPermutation(col,PETSC_DECIDE,&a->icol)); 2280 PetscCall(PetscLogObjectParent((PetscObject)inA,(PetscObject)a->icol)); 2281 2282 PetscCall(MatSeqBAIJSetNumericFactorization_inplace(inA,(PetscBool)(row_identity && col_identity))); 2283 if (!a->solve_work) { 2284 PetscCall(PetscMalloc1(inA->rmap->N+inA->rmap->bs,&a->solve_work)); 2285 PetscCall(PetscLogObjectMemory((PetscObject)inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar))); 2286 } 2287 PetscCall(MatLUFactorNumeric(outA,inA,info)); 2288 PetscFunctionReturn(0); 2289 } 2290 2291 PetscErrorCode MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,PetscInt *indices) 2292 { 2293 Mat_SeqBAIJ *baij = (Mat_SeqBAIJ*)mat->data; 2294 PetscInt i,nz,mbs; 2295 2296 PetscFunctionBegin; 2297 nz = baij->maxnz; 2298 mbs = baij->mbs; 2299 for (i=0; i<nz; i++) { 2300 baij->j[i] = indices[i]; 2301 } 2302 baij->nz = nz; 2303 for (i=0; i<mbs; i++) { 2304 baij->ilen[i] = baij->imax[i]; 2305 } 2306 PetscFunctionReturn(0); 2307 } 2308 2309 /*@ 2310 MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 2311 in the matrix. 2312 2313 Input Parameters: 2314 + mat - the SeqBAIJ matrix 2315 - indices - the column indices 2316 2317 Level: advanced 2318 2319 Notes: 2320 This can be called if you have precomputed the nonzero structure of the 2321 matrix and want to provide it to the matrix object to improve the performance 2322 of the MatSetValues() operation. 2323 2324 You MUST have set the correct numbers of nonzeros per row in the call to 2325 MatCreateSeqBAIJ(), and the columns indices MUST be sorted. 2326 2327 MUST be called before any calls to MatSetValues(); 2328 2329 @*/ 2330 PetscErrorCode MatSeqBAIJSetColumnIndices(Mat mat,PetscInt *indices) 2331 { 2332 PetscFunctionBegin; 2333 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 2334 PetscValidIntPointer(indices,2); 2335 PetscUseMethod(mat,"MatSeqBAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices)); 2336 PetscFunctionReturn(0); 2337 } 2338 2339 PetscErrorCode MatGetRowMaxAbs_SeqBAIJ(Mat A,Vec v,PetscInt idx[]) 2340 { 2341 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2342 PetscInt i,j,n,row,bs,*ai,*aj,mbs; 2343 PetscReal atmp; 2344 PetscScalar *x,zero = 0.0; 2345 MatScalar *aa; 2346 PetscInt ncols,brow,krow,kcol; 2347 2348 PetscFunctionBegin; 2349 /* why is this not a macro???????????????????????????????????????????????????????????????? */ 2350 PetscCheck(!A->factortype,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 2351 bs = A->rmap->bs; 2352 aa = a->a; 2353 ai = a->i; 2354 aj = a->j; 2355 mbs = a->mbs; 2356 2357 PetscCall(VecSet(v,zero)); 2358 PetscCall(VecGetArray(v,&x)); 2359 PetscCall(VecGetLocalSize(v,&n)); 2360 PetscCheck(n == A->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 2361 for (i=0; i<mbs; i++) { 2362 ncols = ai[1] - ai[0]; ai++; 2363 brow = bs*i; 2364 for (j=0; j<ncols; j++) { 2365 for (kcol=0; kcol<bs; kcol++) { 2366 for (krow=0; krow<bs; krow++) { 2367 atmp = PetscAbsScalar(*aa);aa++; 2368 row = brow + krow; /* row index */ 2369 if (PetscAbsScalar(x[row]) < atmp) {x[row] = atmp; if (idx) idx[row] = bs*(*aj) + kcol;} 2370 } 2371 } 2372 aj++; 2373 } 2374 } 2375 PetscCall(VecRestoreArray(v,&x)); 2376 PetscFunctionReturn(0); 2377 } 2378 2379 PetscErrorCode MatCopy_SeqBAIJ(Mat A,Mat B,MatStructure str) 2380 { 2381 PetscFunctionBegin; 2382 /* If the two matrices have the same copy implementation, use fast copy. */ 2383 if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 2384 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2385 Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)B->data; 2386 PetscInt ambs=a->mbs,bmbs=b->mbs,abs=A->rmap->bs,bbs=B->rmap->bs,bs2=abs*abs; 2387 2388 PetscCheck(a->i[ambs] == b->i[bmbs],PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of nonzero blocks in matrices A %" PetscInt_FMT " and B %" PetscInt_FMT " are different",a->i[ambs],b->i[bmbs]); 2389 PetscCheck(abs == bbs,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Block size A %" PetscInt_FMT " and B %" PetscInt_FMT " are different",abs,bbs); 2390 PetscCall(PetscArraycpy(b->a,a->a,bs2*a->i[ambs])); 2391 PetscCall(PetscObjectStateIncrease((PetscObject)B)); 2392 } else { 2393 PetscCall(MatCopy_Basic(A,B,str)); 2394 } 2395 PetscFunctionReturn(0); 2396 } 2397 2398 PetscErrorCode MatSetUp_SeqBAIJ(Mat A) 2399 { 2400 PetscFunctionBegin; 2401 PetscCall(MatSeqBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,NULL)); 2402 PetscFunctionReturn(0); 2403 } 2404 2405 static PetscErrorCode MatSeqBAIJGetArray_SeqBAIJ(Mat A,PetscScalar *array[]) 2406 { 2407 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2408 2409 PetscFunctionBegin; 2410 *array = a->a; 2411 PetscFunctionReturn(0); 2412 } 2413 2414 static PetscErrorCode MatSeqBAIJRestoreArray_SeqBAIJ(Mat A,PetscScalar *array[]) 2415 { 2416 PetscFunctionBegin; 2417 *array = NULL; 2418 PetscFunctionReturn(0); 2419 } 2420 2421 PetscErrorCode MatAXPYGetPreallocation_SeqBAIJ(Mat Y,Mat X,PetscInt *nnz) 2422 { 2423 PetscInt bs = Y->rmap->bs,mbs = Y->rmap->N/bs; 2424 Mat_SeqBAIJ *x = (Mat_SeqBAIJ*)X->data; 2425 Mat_SeqBAIJ *y = (Mat_SeqBAIJ*)Y->data; 2426 2427 PetscFunctionBegin; 2428 /* Set the number of nonzeros in the new matrix */ 2429 PetscCall(MatAXPYGetPreallocation_SeqX_private(mbs,x->i,x->j,y->i,y->j,nnz)); 2430 PetscFunctionReturn(0); 2431 } 2432 2433 PetscErrorCode MatAXPY_SeqBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) 2434 { 2435 Mat_SeqBAIJ *x = (Mat_SeqBAIJ*)X->data,*y = (Mat_SeqBAIJ*)Y->data; 2436 PetscInt bs=Y->rmap->bs,bs2=bs*bs; 2437 PetscBLASInt one=1; 2438 2439 PetscFunctionBegin; 2440 if (str == UNKNOWN_NONZERO_PATTERN || (PetscDefined(USE_DEBUG) && str == SAME_NONZERO_PATTERN)) { 2441 PetscBool e = x->nz == y->nz && x->mbs == y->mbs && bs == X->rmap->bs ? PETSC_TRUE : PETSC_FALSE; 2442 if (e) { 2443 PetscCall(PetscArraycmp(x->i,y->i,x->mbs+1,&e)); 2444 if (e) { 2445 PetscCall(PetscArraycmp(x->j,y->j,x->i[x->mbs],&e)); 2446 if (e) str = SAME_NONZERO_PATTERN; 2447 } 2448 } 2449 if (!e) PetscCheck(str != SAME_NONZERO_PATTERN,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MatStructure is not SAME_NONZERO_PATTERN"); 2450 } 2451 if (str == SAME_NONZERO_PATTERN) { 2452 PetscScalar alpha = a; 2453 PetscBLASInt bnz; 2454 PetscCall(PetscBLASIntCast(x->nz*bs2,&bnz)); 2455 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one)); 2456 PetscCall(PetscObjectStateIncrease((PetscObject)Y)); 2457 } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 2458 PetscCall(MatAXPY_Basic(Y,a,X,str)); 2459 } else { 2460 Mat B; 2461 PetscInt *nnz; 2462 PetscCheck(bs == X->rmap->bs,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size"); 2463 PetscCall(PetscMalloc1(Y->rmap->N,&nnz)); 2464 PetscCall(MatCreate(PetscObjectComm((PetscObject)Y),&B)); 2465 PetscCall(PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name)); 2466 PetscCall(MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N)); 2467 PetscCall(MatSetBlockSizesFromMats(B,Y,Y)); 2468 PetscCall(MatSetType(B,(MatType) ((PetscObject)Y)->type_name)); 2469 PetscCall(MatAXPYGetPreallocation_SeqBAIJ(Y,X,nnz)); 2470 PetscCall(MatSeqBAIJSetPreallocation(B,bs,0,nnz)); 2471 PetscCall(MatAXPY_BasicWithPreallocation(B,Y,a,X,str)); 2472 PetscCall(MatHeaderMerge(Y,&B)); 2473 PetscCall(PetscFree(nnz)); 2474 } 2475 PetscFunctionReturn(0); 2476 } 2477 2478 PETSC_INTERN PetscErrorCode MatConjugate_SeqBAIJ(Mat A) 2479 { 2480 #if defined(PETSC_USE_COMPLEX) 2481 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2482 PetscInt i,nz = a->bs2*a->i[a->mbs]; 2483 MatScalar *aa = a->a; 2484 2485 PetscFunctionBegin; 2486 for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]); 2487 #else 2488 PetscFunctionBegin; 2489 #endif 2490 PetscFunctionReturn(0); 2491 } 2492 2493 PetscErrorCode MatRealPart_SeqBAIJ(Mat A) 2494 { 2495 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2496 PetscInt i,nz = a->bs2*a->i[a->mbs]; 2497 MatScalar *aa = a->a; 2498 2499 PetscFunctionBegin; 2500 for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 2501 PetscFunctionReturn(0); 2502 } 2503 2504 PetscErrorCode MatImaginaryPart_SeqBAIJ(Mat A) 2505 { 2506 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2507 PetscInt i,nz = a->bs2*a->i[a->mbs]; 2508 MatScalar *aa = a->a; 2509 2510 PetscFunctionBegin; 2511 for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 2512 PetscFunctionReturn(0); 2513 } 2514 2515 /* 2516 Code almost identical to MatGetColumnIJ_SeqAIJ() should share common code 2517 */ 2518 PetscErrorCode MatGetColumnIJ_SeqBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscBool *done) 2519 { 2520 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2521 PetscInt bs = A->rmap->bs,i,*collengths,*cia,*cja,n = A->cmap->n/bs,m = A->rmap->n/bs; 2522 PetscInt nz = a->i[m],row,*jj,mr,col; 2523 2524 PetscFunctionBegin; 2525 *nn = n; 2526 if (!ia) PetscFunctionReturn(0); 2527 PetscCheck(!symmetric,PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for BAIJ matrices"); 2528 else { 2529 PetscCall(PetscCalloc1(n,&collengths)); 2530 PetscCall(PetscMalloc1(n+1,&cia)); 2531 PetscCall(PetscMalloc1(nz,&cja)); 2532 jj = a->j; 2533 for (i=0; i<nz; i++) { 2534 collengths[jj[i]]++; 2535 } 2536 cia[0] = oshift; 2537 for (i=0; i<n; i++) { 2538 cia[i+1] = cia[i] + collengths[i]; 2539 } 2540 PetscCall(PetscArrayzero(collengths,n)); 2541 jj = a->j; 2542 for (row=0; row<m; row++) { 2543 mr = a->i[row+1] - a->i[row]; 2544 for (i=0; i<mr; i++) { 2545 col = *jj++; 2546 2547 cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 2548 } 2549 } 2550 PetscCall(PetscFree(collengths)); 2551 *ia = cia; *ja = cja; 2552 } 2553 PetscFunctionReturn(0); 2554 } 2555 2556 PetscErrorCode MatRestoreColumnIJ_SeqBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool *done) 2557 { 2558 PetscFunctionBegin; 2559 if (!ia) PetscFunctionReturn(0); 2560 PetscCall(PetscFree(*ia)); 2561 PetscCall(PetscFree(*ja)); 2562 PetscFunctionReturn(0); 2563 } 2564 2565 /* 2566 MatGetColumnIJ_SeqBAIJ_Color() and MatRestoreColumnIJ_SeqBAIJ_Color() are customized from 2567 MatGetColumnIJ_SeqBAIJ() and MatRestoreColumnIJ_SeqBAIJ() by adding an output 2568 spidx[], index of a->a, to be used in MatTransposeColoringCreate() and MatFDColoringCreate() 2569 */ 2570 PetscErrorCode MatGetColumnIJ_SeqBAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 2571 { 2572 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2573 PetscInt i,*collengths,*cia,*cja,n=a->nbs,m=a->mbs; 2574 PetscInt nz = a->i[m],row,*jj,mr,col; 2575 PetscInt *cspidx; 2576 2577 PetscFunctionBegin; 2578 *nn = n; 2579 if (!ia) PetscFunctionReturn(0); 2580 2581 PetscCall(PetscCalloc1(n,&collengths)); 2582 PetscCall(PetscMalloc1(n+1,&cia)); 2583 PetscCall(PetscMalloc1(nz,&cja)); 2584 PetscCall(PetscMalloc1(nz,&cspidx)); 2585 jj = a->j; 2586 for (i=0; i<nz; i++) { 2587 collengths[jj[i]]++; 2588 } 2589 cia[0] = oshift; 2590 for (i=0; i<n; i++) { 2591 cia[i+1] = cia[i] + collengths[i]; 2592 } 2593 PetscCall(PetscArrayzero(collengths,n)); 2594 jj = a->j; 2595 for (row=0; row<m; row++) { 2596 mr = a->i[row+1] - a->i[row]; 2597 for (i=0; i<mr; i++) { 2598 col = *jj++; 2599 cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 2600 cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 2601 } 2602 } 2603 PetscCall(PetscFree(collengths)); 2604 *ia = cia; 2605 *ja = cja; 2606 *spidx = cspidx; 2607 PetscFunctionReturn(0); 2608 } 2609 2610 PetscErrorCode MatRestoreColumnIJ_SeqBAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 2611 { 2612 PetscFunctionBegin; 2613 PetscCall(MatRestoreColumnIJ_SeqBAIJ(A,oshift,symmetric,inodecompressed,n,ia,ja,done)); 2614 PetscCall(PetscFree(*spidx)); 2615 PetscFunctionReturn(0); 2616 } 2617 2618 PetscErrorCode MatShift_SeqBAIJ(Mat Y,PetscScalar a) 2619 { 2620 Mat_SeqBAIJ *aij = (Mat_SeqBAIJ*)Y->data; 2621 2622 PetscFunctionBegin; 2623 if (!Y->preallocated || !aij->nz) { 2624 PetscCall(MatSeqBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL)); 2625 } 2626 PetscCall(MatShift_Basic(Y,a)); 2627 PetscFunctionReturn(0); 2628 } 2629 2630 /* -------------------------------------------------------------------*/ 2631 static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 2632 MatGetRow_SeqBAIJ, 2633 MatRestoreRow_SeqBAIJ, 2634 MatMult_SeqBAIJ_N, 2635 /* 4*/ MatMultAdd_SeqBAIJ_N, 2636 MatMultTranspose_SeqBAIJ, 2637 MatMultTransposeAdd_SeqBAIJ, 2638 NULL, 2639 NULL, 2640 NULL, 2641 /* 10*/ NULL, 2642 MatLUFactor_SeqBAIJ, 2643 NULL, 2644 NULL, 2645 MatTranspose_SeqBAIJ, 2646 /* 15*/ MatGetInfo_SeqBAIJ, 2647 MatEqual_SeqBAIJ, 2648 MatGetDiagonal_SeqBAIJ, 2649 MatDiagonalScale_SeqBAIJ, 2650 MatNorm_SeqBAIJ, 2651 /* 20*/ NULL, 2652 MatAssemblyEnd_SeqBAIJ, 2653 MatSetOption_SeqBAIJ, 2654 MatZeroEntries_SeqBAIJ, 2655 /* 24*/ MatZeroRows_SeqBAIJ, 2656 NULL, 2657 NULL, 2658 NULL, 2659 NULL, 2660 /* 29*/ MatSetUp_SeqBAIJ, 2661 NULL, 2662 NULL, 2663 NULL, 2664 NULL, 2665 /* 34*/ MatDuplicate_SeqBAIJ, 2666 NULL, 2667 NULL, 2668 MatILUFactor_SeqBAIJ, 2669 NULL, 2670 /* 39*/ MatAXPY_SeqBAIJ, 2671 MatCreateSubMatrices_SeqBAIJ, 2672 MatIncreaseOverlap_SeqBAIJ, 2673 MatGetValues_SeqBAIJ, 2674 MatCopy_SeqBAIJ, 2675 /* 44*/ NULL, 2676 MatScale_SeqBAIJ, 2677 MatShift_SeqBAIJ, 2678 NULL, 2679 MatZeroRowsColumns_SeqBAIJ, 2680 /* 49*/ NULL, 2681 MatGetRowIJ_SeqBAIJ, 2682 MatRestoreRowIJ_SeqBAIJ, 2683 MatGetColumnIJ_SeqBAIJ, 2684 MatRestoreColumnIJ_SeqBAIJ, 2685 /* 54*/ MatFDColoringCreate_SeqXAIJ, 2686 NULL, 2687 NULL, 2688 NULL, 2689 MatSetValuesBlocked_SeqBAIJ, 2690 /* 59*/ MatCreateSubMatrix_SeqBAIJ, 2691 MatDestroy_SeqBAIJ, 2692 MatView_SeqBAIJ, 2693 NULL, 2694 NULL, 2695 /* 64*/ NULL, 2696 NULL, 2697 NULL, 2698 NULL, 2699 NULL, 2700 /* 69*/ MatGetRowMaxAbs_SeqBAIJ, 2701 NULL, 2702 MatConvert_Basic, 2703 NULL, 2704 NULL, 2705 /* 74*/ NULL, 2706 MatFDColoringApply_BAIJ, 2707 NULL, 2708 NULL, 2709 NULL, 2710 /* 79*/ NULL, 2711 NULL, 2712 NULL, 2713 NULL, 2714 MatLoad_SeqBAIJ, 2715 /* 84*/ NULL, 2716 NULL, 2717 NULL, 2718 NULL, 2719 NULL, 2720 /* 89*/ NULL, 2721 NULL, 2722 NULL, 2723 NULL, 2724 NULL, 2725 /* 94*/ NULL, 2726 NULL, 2727 NULL, 2728 NULL, 2729 NULL, 2730 /* 99*/ NULL, 2731 NULL, 2732 NULL, 2733 MatConjugate_SeqBAIJ, 2734 NULL, 2735 /*104*/ NULL, 2736 MatRealPart_SeqBAIJ, 2737 MatImaginaryPart_SeqBAIJ, 2738 NULL, 2739 NULL, 2740 /*109*/ NULL, 2741 NULL, 2742 NULL, 2743 NULL, 2744 MatMissingDiagonal_SeqBAIJ, 2745 /*114*/ NULL, 2746 NULL, 2747 NULL, 2748 NULL, 2749 NULL, 2750 /*119*/ NULL, 2751 NULL, 2752 MatMultHermitianTranspose_SeqBAIJ, 2753 MatMultHermitianTransposeAdd_SeqBAIJ, 2754 NULL, 2755 /*124*/ NULL, 2756 MatGetColumnReductions_SeqBAIJ, 2757 MatInvertBlockDiagonal_SeqBAIJ, 2758 NULL, 2759 NULL, 2760 /*129*/ NULL, 2761 NULL, 2762 NULL, 2763 NULL, 2764 NULL, 2765 /*134*/ NULL, 2766 NULL, 2767 NULL, 2768 NULL, 2769 NULL, 2770 /*139*/ MatSetBlockSizes_Default, 2771 NULL, 2772 NULL, 2773 MatFDColoringSetUp_SeqXAIJ, 2774 NULL, 2775 /*144*/MatCreateMPIMatConcatenateSeqMat_SeqBAIJ, 2776 MatDestroySubMatrices_SeqBAIJ, 2777 NULL, 2778 NULL 2779 }; 2780 2781 PetscErrorCode MatStoreValues_SeqBAIJ(Mat mat) 2782 { 2783 Mat_SeqBAIJ *aij = (Mat_SeqBAIJ*)mat->data; 2784 PetscInt nz = aij->i[aij->mbs]*aij->bs2; 2785 2786 PetscFunctionBegin; 2787 PetscCheck(aij->nonew == 1,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 2788 2789 /* allocate space for values if not already there */ 2790 if (!aij->saved_values) { 2791 PetscCall(PetscMalloc1(nz+1,&aij->saved_values)); 2792 PetscCall(PetscLogObjectMemory((PetscObject)mat,(nz+1)*sizeof(PetscScalar))); 2793 } 2794 2795 /* copy values over */ 2796 PetscCall(PetscArraycpy(aij->saved_values,aij->a,nz)); 2797 PetscFunctionReturn(0); 2798 } 2799 2800 PetscErrorCode MatRetrieveValues_SeqBAIJ(Mat mat) 2801 { 2802 Mat_SeqBAIJ *aij = (Mat_SeqBAIJ*)mat->data; 2803 PetscInt nz = aij->i[aij->mbs]*aij->bs2; 2804 2805 PetscFunctionBegin; 2806 PetscCheck(aij->nonew == 1,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 2807 PetscCheck(aij->saved_values,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 2808 2809 /* copy values over */ 2810 PetscCall(PetscArraycpy(aij->a,aij->saved_values,nz)); 2811 PetscFunctionReturn(0); 2812 } 2813 2814 PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*); 2815 PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat, MatType,MatReuse,Mat*); 2816 2817 PetscErrorCode MatSeqBAIJSetPreallocation_SeqBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 2818 { 2819 Mat_SeqBAIJ *b; 2820 PetscInt i,mbs,nbs,bs2; 2821 PetscBool flg = PETSC_FALSE,skipallocation = PETSC_FALSE,realalloc = PETSC_FALSE; 2822 2823 PetscFunctionBegin; 2824 if (nz >= 0 || nnz) realalloc = PETSC_TRUE; 2825 if (nz == MAT_SKIP_ALLOCATION) { 2826 skipallocation = PETSC_TRUE; 2827 nz = 0; 2828 } 2829 2830 PetscCall(MatSetBlockSize(B,PetscAbs(bs))); 2831 PetscCall(PetscLayoutSetUp(B->rmap)); 2832 PetscCall(PetscLayoutSetUp(B->cmap)); 2833 PetscCall(PetscLayoutGetBlockSize(B->rmap,&bs)); 2834 2835 B->preallocated = PETSC_TRUE; 2836 2837 mbs = B->rmap->n/bs; 2838 nbs = B->cmap->n/bs; 2839 bs2 = bs*bs; 2840 2841 PetscCheck(mbs*bs==B->rmap->n && nbs*bs==B->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows %" PetscInt_FMT ", cols %" PetscInt_FMT " must be divisible by blocksize %" PetscInt_FMT,B->rmap->N,B->cmap->n,bs); 2842 2843 if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 2844 PetscCheck(nz >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %" PetscInt_FMT,nz); 2845 if (nnz) { 2846 for (i=0; i<mbs; i++) { 2847 PetscCheck(nnz[i] >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %" PetscInt_FMT " value %" PetscInt_FMT,i,nnz[i]); 2848 PetscCheck(nnz[i] <= nbs,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %" PetscInt_FMT " value %" PetscInt_FMT " rowlength %" PetscInt_FMT,i,nnz[i],nbs); 2849 } 2850 } 2851 2852 b = (Mat_SeqBAIJ*)B->data; 2853 PetscOptionsBegin(PetscObjectComm((PetscObject)B),NULL,"Optimize options for SEQBAIJ matrix 2 ","Mat"); 2854 PetscCall(PetscOptionsBool("-mat_no_unroll","Do not optimize for block size (slow)",NULL,flg,&flg,NULL)); 2855 PetscOptionsEnd(); 2856 2857 if (!flg) { 2858 switch (bs) { 2859 case 1: 2860 B->ops->mult = MatMult_SeqBAIJ_1; 2861 B->ops->multadd = MatMultAdd_SeqBAIJ_1; 2862 break; 2863 case 2: 2864 B->ops->mult = MatMult_SeqBAIJ_2; 2865 B->ops->multadd = MatMultAdd_SeqBAIJ_2; 2866 break; 2867 case 3: 2868 B->ops->mult = MatMult_SeqBAIJ_3; 2869 B->ops->multadd = MatMultAdd_SeqBAIJ_3; 2870 break; 2871 case 4: 2872 B->ops->mult = MatMult_SeqBAIJ_4; 2873 B->ops->multadd = MatMultAdd_SeqBAIJ_4; 2874 break; 2875 case 5: 2876 B->ops->mult = MatMult_SeqBAIJ_5; 2877 B->ops->multadd = MatMultAdd_SeqBAIJ_5; 2878 break; 2879 case 6: 2880 B->ops->mult = MatMult_SeqBAIJ_6; 2881 B->ops->multadd = MatMultAdd_SeqBAIJ_6; 2882 break; 2883 case 7: 2884 B->ops->mult = MatMult_SeqBAIJ_7; 2885 B->ops->multadd = MatMultAdd_SeqBAIJ_7; 2886 break; 2887 case 9: 2888 { 2889 PetscInt version = 1; 2890 PetscCall(PetscOptionsGetInt(NULL,((PetscObject)B)->prefix,"-mat_baij_mult_version",&version,NULL)); 2891 switch (version) { 2892 #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX2__) && defined(__FMA__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES) 2893 case 1: 2894 B->ops->mult = MatMult_SeqBAIJ_9_AVX2; 2895 B->ops->multadd = MatMultAdd_SeqBAIJ_9_AVX2; 2896 PetscCall(PetscInfo((PetscObject)B,"Using AVX2 for MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",bs)); 2897 break; 2898 #endif 2899 default: 2900 B->ops->mult = MatMult_SeqBAIJ_N; 2901 B->ops->multadd = MatMultAdd_SeqBAIJ_N; 2902 PetscCall(PetscInfo((PetscObject)B,"Using BLAS for MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",bs)); 2903 break; 2904 } 2905 break; 2906 } 2907 case 11: 2908 B->ops->mult = MatMult_SeqBAIJ_11; 2909 B->ops->multadd = MatMultAdd_SeqBAIJ_11; 2910 break; 2911 case 12: 2912 { 2913 PetscInt version = 1; 2914 PetscCall(PetscOptionsGetInt(NULL,((PetscObject)B)->prefix,"-mat_baij_mult_version",&version,NULL)); 2915 switch (version) { 2916 case 1: 2917 B->ops->mult = MatMult_SeqBAIJ_12_ver1; 2918 B->ops->multadd = MatMultAdd_SeqBAIJ_12_ver1; 2919 PetscCall(PetscInfo((PetscObject)B,"Using version %" PetscInt_FMT " of MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",version,bs)); 2920 break; 2921 case 2: 2922 B->ops->mult = MatMult_SeqBAIJ_12_ver2; 2923 B->ops->multadd = MatMultAdd_SeqBAIJ_12_ver2; 2924 PetscCall(PetscInfo((PetscObject)B,"Using version %" PetscInt_FMT " of MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",version,bs)); 2925 break; 2926 #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX2__) && defined(__FMA__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES) 2927 case 3: 2928 B->ops->mult = MatMult_SeqBAIJ_12_AVX2; 2929 B->ops->multadd = MatMultAdd_SeqBAIJ_12_ver1; 2930 PetscCall(PetscInfo((PetscObject)B,"Using AVX2 for MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",bs)); 2931 break; 2932 #endif 2933 default: 2934 B->ops->mult = MatMult_SeqBAIJ_N; 2935 B->ops->multadd = MatMultAdd_SeqBAIJ_N; 2936 PetscCall(PetscInfo((PetscObject)B,"Using BLAS for MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",bs)); 2937 break; 2938 } 2939 break; 2940 } 2941 case 15: 2942 { 2943 PetscInt version = 1; 2944 PetscCall(PetscOptionsGetInt(NULL,((PetscObject)B)->prefix,"-mat_baij_mult_version",&version,NULL)); 2945 switch (version) { 2946 case 1: 2947 B->ops->mult = MatMult_SeqBAIJ_15_ver1; 2948 PetscCall(PetscInfo((PetscObject)B,"Using version %" PetscInt_FMT " of MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",version,bs)); 2949 break; 2950 case 2: 2951 B->ops->mult = MatMult_SeqBAIJ_15_ver2; 2952 PetscCall(PetscInfo((PetscObject)B,"Using version %" PetscInt_FMT " of MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",version,bs)); 2953 break; 2954 case 3: 2955 B->ops->mult = MatMult_SeqBAIJ_15_ver3; 2956 PetscCall(PetscInfo((PetscObject)B,"Using version %" PetscInt_FMT " of MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",version,bs)); 2957 break; 2958 case 4: 2959 B->ops->mult = MatMult_SeqBAIJ_15_ver4; 2960 PetscCall(PetscInfo((PetscObject)B,"Using version %" PetscInt_FMT " of MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",version,bs)); 2961 break; 2962 default: 2963 B->ops->mult = MatMult_SeqBAIJ_N; 2964 PetscCall(PetscInfo((PetscObject)B,"Using BLAS for MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",bs)); 2965 break; 2966 } 2967 B->ops->multadd = MatMultAdd_SeqBAIJ_N; 2968 break; 2969 } 2970 default: 2971 B->ops->mult = MatMult_SeqBAIJ_N; 2972 B->ops->multadd = MatMultAdd_SeqBAIJ_N; 2973 PetscCall(PetscInfo((PetscObject)B,"Using BLAS for MatMult for BAIJ for blocksize %" PetscInt_FMT "\n",bs)); 2974 break; 2975 } 2976 } 2977 B->ops->sor = MatSOR_SeqBAIJ; 2978 b->mbs = mbs; 2979 b->nbs = nbs; 2980 if (!skipallocation) { 2981 if (!b->imax) { 2982 PetscCall(PetscMalloc2(mbs,&b->imax,mbs,&b->ilen)); 2983 PetscCall(PetscLogObjectMemory((PetscObject)B,2*mbs*sizeof(PetscInt))); 2984 2985 b->free_imax_ilen = PETSC_TRUE; 2986 } 2987 /* b->ilen will count nonzeros in each block row so far. */ 2988 for (i=0; i<mbs; i++) b->ilen[i] = 0; 2989 if (!nnz) { 2990 if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 2991 else if (nz < 0) nz = 1; 2992 nz = PetscMin(nz,nbs); 2993 for (i=0; i<mbs; i++) b->imax[i] = nz; 2994 PetscCall(PetscIntMultError(nz,mbs,&nz)); 2995 } else { 2996 PetscInt64 nz64 = 0; 2997 for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz64 += nnz[i];} 2998 PetscCall(PetscIntCast(nz64,&nz)); 2999 } 3000 3001 /* allocate the matrix space */ 3002 PetscCall(MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i)); 3003 if (B->structure_only) { 3004 PetscCall(PetscMalloc1(nz,&b->j)); 3005 PetscCall(PetscMalloc1(B->rmap->N+1,&b->i)); 3006 PetscCall(PetscLogObjectMemory((PetscObject)B,(B->rmap->N+1)*sizeof(PetscInt)+nz*sizeof(PetscInt))); 3007 } else { 3008 PetscInt nzbs2 = 0; 3009 PetscCall(PetscIntMultError(nz,bs2,&nzbs2)); 3010 PetscCall(PetscMalloc3(nzbs2,&b->a,nz,&b->j,B->rmap->N+1,&b->i)); 3011 PetscCall(PetscLogObjectMemory((PetscObject)B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)))); 3012 PetscCall(PetscArrayzero(b->a,nz*bs2)); 3013 } 3014 PetscCall(PetscArrayzero(b->j,nz)); 3015 3016 if (B->structure_only) { 3017 b->singlemalloc = PETSC_FALSE; 3018 b->free_a = PETSC_FALSE; 3019 } else { 3020 b->singlemalloc = PETSC_TRUE; 3021 b->free_a = PETSC_TRUE; 3022 } 3023 b->free_ij = PETSC_TRUE; 3024 3025 b->i[0] = 0; 3026 for (i=1; i<mbs+1; i++) { 3027 b->i[i] = b->i[i-1] + b->imax[i-1]; 3028 } 3029 3030 } else { 3031 b->free_a = PETSC_FALSE; 3032 b->free_ij = PETSC_FALSE; 3033 } 3034 3035 b->bs2 = bs2; 3036 b->mbs = mbs; 3037 b->nz = 0; 3038 b->maxnz = nz; 3039 B->info.nz_unneeded = (PetscReal)b->maxnz*bs2; 3040 B->was_assembled = PETSC_FALSE; 3041 B->assembled = PETSC_FALSE; 3042 if (realalloc) PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 3043 PetscFunctionReturn(0); 3044 } 3045 3046 PetscErrorCode MatSeqBAIJSetPreallocationCSR_SeqBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[]) 3047 { 3048 PetscInt i,m,nz,nz_max=0,*nnz; 3049 PetscScalar *values=NULL; 3050 PetscBool roworiented = ((Mat_SeqBAIJ*)B->data)->roworiented; 3051 3052 PetscFunctionBegin; 3053 PetscCheck(bs >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %" PetscInt_FMT,bs); 3054 PetscCall(PetscLayoutSetBlockSize(B->rmap,bs)); 3055 PetscCall(PetscLayoutSetBlockSize(B->cmap,bs)); 3056 PetscCall(PetscLayoutSetUp(B->rmap)); 3057 PetscCall(PetscLayoutSetUp(B->cmap)); 3058 PetscCall(PetscLayoutGetBlockSize(B->rmap,&bs)); 3059 m = B->rmap->n/bs; 3060 3061 PetscCheck(ii[0] == 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "ii[0] must be 0 but it is %" PetscInt_FMT,ii[0]); 3062 PetscCall(PetscMalloc1(m+1, &nnz)); 3063 for (i=0; i<m; i++) { 3064 nz = ii[i+1]- ii[i]; 3065 PetscCheck(nz >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Local row %" PetscInt_FMT " has a negative number of columns %" PetscInt_FMT,i,nz); 3066 nz_max = PetscMax(nz_max, nz); 3067 nnz[i] = nz; 3068 } 3069 PetscCall(MatSeqBAIJSetPreallocation(B,bs,0,nnz)); 3070 PetscCall(PetscFree(nnz)); 3071 3072 values = (PetscScalar*)V; 3073 if (!values) { 3074 PetscCall(PetscCalloc1(bs*bs*(nz_max+1),&values)); 3075 } 3076 for (i=0; i<m; i++) { 3077 PetscInt ncols = ii[i+1] - ii[i]; 3078 const PetscInt *icols = jj + ii[i]; 3079 if (bs == 1 || !roworiented) { 3080 const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0); 3081 PetscCall(MatSetValuesBlocked_SeqBAIJ(B,1,&i,ncols,icols,svals,INSERT_VALUES)); 3082 } else { 3083 PetscInt j; 3084 for (j=0; j<ncols; j++) { 3085 const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0); 3086 PetscCall(MatSetValuesBlocked_SeqBAIJ(B,1,&i,1,&icols[j],svals,INSERT_VALUES)); 3087 } 3088 } 3089 } 3090 if (!V) PetscCall(PetscFree(values)); 3091 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 3092 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 3093 PetscCall(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 3094 PetscFunctionReturn(0); 3095 } 3096 3097 /*@C 3098 MatSeqBAIJGetArray - gives access to the array where the data for a MATSEQBAIJ matrix is stored 3099 3100 Not Collective 3101 3102 Input Parameter: 3103 . mat - a MATSEQBAIJ matrix 3104 3105 Output Parameter: 3106 . array - pointer to the data 3107 3108 Level: intermediate 3109 3110 .seealso: MatSeqBAIJRestoreArray(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray() 3111 @*/ 3112 PetscErrorCode MatSeqBAIJGetArray(Mat A,PetscScalar **array) 3113 { 3114 PetscFunctionBegin; 3115 PetscUseMethod(A,"MatSeqBAIJGetArray_C",(Mat,PetscScalar**),(A,array)); 3116 PetscFunctionReturn(0); 3117 } 3118 3119 /*@C 3120 MatSeqBAIJRestoreArray - returns access to the array where the data for a MATSEQBAIJ matrix is stored obtained by MatSeqBAIJGetArray() 3121 3122 Not Collective 3123 3124 Input Parameters: 3125 + mat - a MATSEQBAIJ matrix 3126 - array - pointer to the data 3127 3128 Level: intermediate 3129 3130 .seealso: MatSeqBAIJGetArray(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray() 3131 @*/ 3132 PetscErrorCode MatSeqBAIJRestoreArray(Mat A,PetscScalar **array) 3133 { 3134 PetscFunctionBegin; 3135 PetscUseMethod(A,"MatSeqBAIJRestoreArray_C",(Mat,PetscScalar**),(A,array)); 3136 PetscFunctionReturn(0); 3137 } 3138 3139 /*MC 3140 MATSEQBAIJ - MATSEQBAIJ = "seqbaij" - A matrix type to be used for sequential block sparse matrices, based on 3141 block sparse compressed row format. 3142 3143 Options Database Keys: 3144 + -mat_type seqbaij - sets the matrix type to "seqbaij" during a call to MatSetFromOptions() 3145 - -mat_baij_mult_version version - indicate the version of the matrix-vector product to use (0 often indicates using BLAS) 3146 3147 Level: beginner 3148 3149 Notes: 3150 MatSetOptions(,MAT_STRUCTURE_ONLY,PETSC_TRUE) may be called for this matrix type. In this no 3151 space is allocated for the nonzero entries and any entries passed with MatSetValues() are ignored 3152 3153 Run with -info to see what version of the matrix-vector product is being used 3154 3155 .seealso: MatCreateSeqBAIJ() 3156 M*/ 3157 3158 PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqBSTRM(Mat, MatType,MatReuse,Mat*); 3159 3160 PETSC_EXTERN PetscErrorCode MatCreate_SeqBAIJ(Mat B) 3161 { 3162 PetscMPIInt size; 3163 Mat_SeqBAIJ *b; 3164 3165 PetscFunctionBegin; 3166 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B),&size)); 3167 PetscCheck(size == 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 3168 3169 PetscCall(PetscNewLog(B,&b)); 3170 B->data = (void*)b; 3171 PetscCall(PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps))); 3172 3173 b->row = NULL; 3174 b->col = NULL; 3175 b->icol = NULL; 3176 b->reallocs = 0; 3177 b->saved_values = NULL; 3178 3179 b->roworiented = PETSC_TRUE; 3180 b->nonew = 0; 3181 b->diag = NULL; 3182 B->spptr = NULL; 3183 B->info.nz_unneeded = (PetscReal)b->maxnz*b->bs2; 3184 b->keepnonzeropattern = PETSC_FALSE; 3185 3186 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJGetArray_C",MatSeqBAIJGetArray_SeqBAIJ)); 3187 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJRestoreArray_C",MatSeqBAIJRestoreArray_SeqBAIJ)); 3188 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqBAIJ)); 3189 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqBAIJ)); 3190 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",MatSeqBAIJSetColumnIndices_SeqBAIJ)); 3191 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqbaij_seqaij_C",MatConvert_SeqBAIJ_SeqAIJ)); 3192 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqbaij_seqsbaij_C",MatConvert_SeqBAIJ_SeqSBAIJ)); 3193 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetPreallocation_C",MatSeqBAIJSetPreallocation_SeqBAIJ)); 3194 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetPreallocationCSR_C",MatSeqBAIJSetPreallocationCSR_SeqBAIJ)); 3195 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatIsTranspose_C",MatIsTranspose_SeqBAIJ)); 3196 #if defined(PETSC_HAVE_HYPRE) 3197 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqbaij_hypre_C",MatConvert_AIJ_HYPRE)); 3198 #endif 3199 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqbaij_is_C",MatConvert_XAIJ_IS)); 3200 PetscCall(PetscObjectChangeTypeName((PetscObject)B,MATSEQBAIJ)); 3201 PetscFunctionReturn(0); 3202 } 3203 3204 PetscErrorCode MatDuplicateNoCreate_SeqBAIJ(Mat C,Mat A,MatDuplicateOption cpvalues,PetscBool mallocmatspace) 3205 { 3206 Mat_SeqBAIJ *c = (Mat_SeqBAIJ*)C->data,*a = (Mat_SeqBAIJ*)A->data; 3207 PetscInt i,mbs = a->mbs,nz = a->nz,bs2 = a->bs2; 3208 3209 PetscFunctionBegin; 3210 PetscCheck(a->i[mbs] == nz,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt matrix"); 3211 3212 if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 3213 c->imax = a->imax; 3214 c->ilen = a->ilen; 3215 c->free_imax_ilen = PETSC_FALSE; 3216 } else { 3217 PetscCall(PetscMalloc2(mbs,&c->imax,mbs,&c->ilen)); 3218 PetscCall(PetscLogObjectMemory((PetscObject)C,2*mbs*sizeof(PetscInt))); 3219 for (i=0; i<mbs; i++) { 3220 c->imax[i] = a->imax[i]; 3221 c->ilen[i] = a->ilen[i]; 3222 } 3223 c->free_imax_ilen = PETSC_TRUE; 3224 } 3225 3226 /* allocate the matrix space */ 3227 if (mallocmatspace) { 3228 if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 3229 PetscCall(PetscCalloc1(bs2*nz,&c->a)); 3230 PetscCall(PetscLogObjectMemory((PetscObject)C,a->i[mbs]*bs2*sizeof(PetscScalar))); 3231 3232 c->i = a->i; 3233 c->j = a->j; 3234 c->singlemalloc = PETSC_FALSE; 3235 c->free_a = PETSC_TRUE; 3236 c->free_ij = PETSC_FALSE; 3237 c->parent = A; 3238 C->preallocated = PETSC_TRUE; 3239 C->assembled = PETSC_TRUE; 3240 3241 PetscCall(PetscObjectReference((PetscObject)A)); 3242 PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 3243 PetscCall(MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 3244 } else { 3245 PetscCall(PetscMalloc3(bs2*nz,&c->a,nz,&c->j,mbs+1,&c->i)); 3246 PetscCall(PetscLogObjectMemory((PetscObject)C,a->i[mbs]*(bs2*sizeof(PetscScalar)+sizeof(PetscInt))+(mbs+1)*sizeof(PetscInt))); 3247 3248 c->singlemalloc = PETSC_TRUE; 3249 c->free_a = PETSC_TRUE; 3250 c->free_ij = PETSC_TRUE; 3251 3252 PetscCall(PetscArraycpy(c->i,a->i,mbs+1)); 3253 if (mbs > 0) { 3254 PetscCall(PetscArraycpy(c->j,a->j,nz)); 3255 if (cpvalues == MAT_COPY_VALUES) { 3256 PetscCall(PetscArraycpy(c->a,a->a,bs2*nz)); 3257 } else { 3258 PetscCall(PetscArrayzero(c->a,bs2*nz)); 3259 } 3260 } 3261 C->preallocated = PETSC_TRUE; 3262 C->assembled = PETSC_TRUE; 3263 } 3264 } 3265 3266 c->roworiented = a->roworiented; 3267 c->nonew = a->nonew; 3268 3269 PetscCall(PetscLayoutReference(A->rmap,&C->rmap)); 3270 PetscCall(PetscLayoutReference(A->cmap,&C->cmap)); 3271 3272 c->bs2 = a->bs2; 3273 c->mbs = a->mbs; 3274 c->nbs = a->nbs; 3275 3276 if (a->diag) { 3277 if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 3278 c->diag = a->diag; 3279 c->free_diag = PETSC_FALSE; 3280 } else { 3281 PetscCall(PetscMalloc1(mbs+1,&c->diag)); 3282 PetscCall(PetscLogObjectMemory((PetscObject)C,(mbs+1)*sizeof(PetscInt))); 3283 for (i=0; i<mbs; i++) c->diag[i] = a->diag[i]; 3284 c->free_diag = PETSC_TRUE; 3285 } 3286 } else c->diag = NULL; 3287 3288 c->nz = a->nz; 3289 c->maxnz = a->nz; /* Since we allocate exactly the right amount */ 3290 c->solve_work = NULL; 3291 c->mult_work = NULL; 3292 c->sor_workt = NULL; 3293 c->sor_work = NULL; 3294 3295 c->compressedrow.use = a->compressedrow.use; 3296 c->compressedrow.nrows = a->compressedrow.nrows; 3297 if (a->compressedrow.use) { 3298 i = a->compressedrow.nrows; 3299 PetscCall(PetscMalloc2(i+1,&c->compressedrow.i,i+1,&c->compressedrow.rindex)); 3300 PetscCall(PetscLogObjectMemory((PetscObject)C,(2*i+1)*sizeof(PetscInt))); 3301 PetscCall(PetscArraycpy(c->compressedrow.i,a->compressedrow.i,i+1)); 3302 PetscCall(PetscArraycpy(c->compressedrow.rindex,a->compressedrow.rindex,i)); 3303 } else { 3304 c->compressedrow.use = PETSC_FALSE; 3305 c->compressedrow.i = NULL; 3306 c->compressedrow.rindex = NULL; 3307 } 3308 C->nonzerostate = A->nonzerostate; 3309 3310 PetscCall(PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist)); 3311 PetscFunctionReturn(0); 3312 } 3313 3314 PetscErrorCode MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 3315 { 3316 PetscFunctionBegin; 3317 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),B)); 3318 PetscCall(MatSetSizes(*B,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n)); 3319 PetscCall(MatSetType(*B,MATSEQBAIJ)); 3320 PetscCall(MatDuplicateNoCreate_SeqBAIJ(*B,A,cpvalues,PETSC_TRUE)); 3321 PetscFunctionReturn(0); 3322 } 3323 3324 /* Used for both SeqBAIJ and SeqSBAIJ matrices */ 3325 PetscErrorCode MatLoad_SeqBAIJ_Binary(Mat mat,PetscViewer viewer) 3326 { 3327 PetscInt header[4],M,N,nz,bs,m,n,mbs,nbs,rows,cols,sum,i,j,k; 3328 PetscInt *rowidxs,*colidxs; 3329 PetscScalar *matvals; 3330 3331 PetscFunctionBegin; 3332 PetscCall(PetscViewerSetUp(viewer)); 3333 3334 /* read matrix header */ 3335 PetscCall(PetscViewerBinaryRead(viewer,header,4,NULL,PETSC_INT)); 3336 PetscCheck(header[0] == MAT_FILE_CLASSID,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not a matrix object in file"); 3337 M = header[1]; N = header[2]; nz = header[3]; 3338 PetscCheck(M >= 0,PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Matrix row size (%" PetscInt_FMT ") in file is negative",M); 3339 PetscCheck(N >= 0,PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Matrix column size (%" PetscInt_FMT ") in file is negative",N); 3340 PetscCheck(nz >= 0,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format on disk, cannot load as SeqBAIJ"); 3341 3342 /* set block sizes from the viewer's .info file */ 3343 PetscCall(MatLoad_Binary_BlockSizes(mat,viewer)); 3344 /* set local and global sizes if not set already */ 3345 if (mat->rmap->n < 0) mat->rmap->n = M; 3346 if (mat->cmap->n < 0) mat->cmap->n = N; 3347 if (mat->rmap->N < 0) mat->rmap->N = M; 3348 if (mat->cmap->N < 0) mat->cmap->N = N; 3349 PetscCall(PetscLayoutSetUp(mat->rmap)); 3350 PetscCall(PetscLayoutSetUp(mat->cmap)); 3351 3352 /* check if the matrix sizes are correct */ 3353 PetscCall(MatGetSize(mat,&rows,&cols)); 3354 PetscCheck(M == rows && N == cols,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different sizes (%" PetscInt_FMT ", %" PetscInt_FMT ") than the input matrix (%" PetscInt_FMT ", %" PetscInt_FMT ")",M,N,rows,cols); 3355 PetscCall(MatGetBlockSize(mat,&bs)); 3356 PetscCall(MatGetLocalSize(mat,&m,&n)); 3357 mbs = m/bs; nbs = n/bs; 3358 3359 /* read in row lengths, column indices and nonzero values */ 3360 PetscCall(PetscMalloc1(m+1,&rowidxs)); 3361 PetscCall(PetscViewerBinaryRead(viewer,rowidxs+1,m,NULL,PETSC_INT)); 3362 rowidxs[0] = 0; for (i=0; i<m; i++) rowidxs[i+1] += rowidxs[i]; 3363 sum = rowidxs[m]; 3364 PetscCheck(sum == nz,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Inconsistent matrix data in file: nonzeros = %" PetscInt_FMT ", sum-row-lengths = %" PetscInt_FMT,nz,sum); 3365 3366 /* read in column indices and nonzero values */ 3367 PetscCall(PetscMalloc2(rowidxs[m],&colidxs,nz,&matvals)); 3368 PetscCall(PetscViewerBinaryRead(viewer,colidxs,rowidxs[m],NULL,PETSC_INT)); 3369 PetscCall(PetscViewerBinaryRead(viewer,matvals,rowidxs[m],NULL,PETSC_SCALAR)); 3370 3371 { /* preallocate matrix storage */ 3372 PetscBT bt; /* helper bit set to count nonzeros */ 3373 PetscInt *nnz; 3374 PetscBool sbaij; 3375 3376 PetscCall(PetscBTCreate(nbs,&bt)); 3377 PetscCall(PetscCalloc1(mbs,&nnz)); 3378 PetscCall(PetscObjectTypeCompare((PetscObject)mat,MATSEQSBAIJ,&sbaij)); 3379 for (i=0; i<mbs; i++) { 3380 PetscCall(PetscBTMemzero(nbs,bt)); 3381 for (k=0; k<bs; k++) { 3382 PetscInt row = bs*i + k; 3383 for (j=rowidxs[row]; j<rowidxs[row+1]; j++) { 3384 PetscInt col = colidxs[j]; 3385 if (!sbaij || col >= row) 3386 if (!PetscBTLookupSet(bt,col/bs)) nnz[i]++; 3387 } 3388 } 3389 } 3390 PetscCall(PetscBTDestroy(&bt)); 3391 PetscCall(MatSeqBAIJSetPreallocation(mat,bs,0,nnz)); 3392 PetscCall(MatSeqSBAIJSetPreallocation(mat,bs,0,nnz)); 3393 PetscCall(PetscFree(nnz)); 3394 } 3395 3396 /* store matrix values */ 3397 for (i=0; i<m; i++) { 3398 PetscInt row = i, s = rowidxs[i], e = rowidxs[i+1]; 3399 PetscCall((*mat->ops->setvalues)(mat,1,&row,e-s,colidxs+s,matvals+s,INSERT_VALUES)); 3400 } 3401 3402 PetscCall(PetscFree(rowidxs)); 3403 PetscCall(PetscFree2(colidxs,matvals)); 3404 PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 3405 PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 3406 PetscFunctionReturn(0); 3407 } 3408 3409 PetscErrorCode MatLoad_SeqBAIJ(Mat mat,PetscViewer viewer) 3410 { 3411 PetscBool isbinary; 3412 3413 PetscFunctionBegin; 3414 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 3415 PetscCheck(isbinary,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)mat)->type_name); 3416 PetscCall(MatLoad_SeqBAIJ_Binary(mat,viewer)); 3417 PetscFunctionReturn(0); 3418 } 3419 3420 /*@C 3421 MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 3422 compressed row) format. For good matrix assembly performance the 3423 user should preallocate the matrix storage by setting the parameter nz 3424 (or the array nnz). By setting these parameters accurately, performance 3425 during matrix assembly can be increased by more than a factor of 50. 3426 3427 Collective 3428 3429 Input Parameters: 3430 + comm - MPI communicator, set to PETSC_COMM_SELF 3431 . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row 3432 blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs() 3433 . m - number of rows 3434 . n - number of columns 3435 . nz - number of nonzero blocks per block row (same for all rows) 3436 - nnz - array containing the number of nonzero blocks in the various block rows 3437 (possibly different for each block row) or NULL 3438 3439 Output Parameter: 3440 . A - the matrix 3441 3442 It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(), 3443 MatXXXXSetPreallocation() paradigm instead of this routine directly. 3444 [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation] 3445 3446 Options Database Keys: 3447 + -mat_no_unroll - uses code that does not unroll the loops in the 3448 block calculations (much slower) 3449 - -mat_block_size - size of the blocks to use 3450 3451 Level: intermediate 3452 3453 Notes: 3454 The number of rows and columns must be divisible by blocksize. 3455 3456 If the nnz parameter is given then the nz parameter is ignored 3457 3458 A nonzero block is any block that as 1 or more nonzeros in it 3459 3460 The block AIJ format is fully compatible with standard Fortran 77 3461 storage. That is, the stored row and column indices can begin at 3462 either one (as in Fortran) or zero. See the users' manual for details. 3463 3464 Specify the preallocated storage with either nz or nnz (not both). 3465 Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory 3466 allocation. See Users-Manual: ch_mat for details. 3467 matrices. 3468 3469 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateBAIJ() 3470 @*/ 3471 PetscErrorCode MatCreateSeqBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 3472 { 3473 PetscFunctionBegin; 3474 PetscCall(MatCreate(comm,A)); 3475 PetscCall(MatSetSizes(*A,m,n,m,n)); 3476 PetscCall(MatSetType(*A,MATSEQBAIJ)); 3477 PetscCall(MatSeqBAIJSetPreallocation(*A,bs,nz,(PetscInt*)nnz)); 3478 PetscFunctionReturn(0); 3479 } 3480 3481 /*@C 3482 MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros 3483 per row in the matrix. For good matrix assembly performance the 3484 user should preallocate the matrix storage by setting the parameter nz 3485 (or the array nnz). By setting these parameters accurately, performance 3486 during matrix assembly can be increased by more than a factor of 50. 3487 3488 Collective 3489 3490 Input Parameters: 3491 + B - the matrix 3492 . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row 3493 blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs() 3494 . nz - number of block nonzeros per block row (same for all rows) 3495 - nnz - array containing the number of block nonzeros in the various block rows 3496 (possibly different for each block row) or NULL 3497 3498 Options Database Keys: 3499 + -mat_no_unroll - uses code that does not unroll the loops in the 3500 block calculations (much slower) 3501 - -mat_block_size - size of the blocks to use 3502 3503 Level: intermediate 3504 3505 Notes: 3506 If the nnz parameter is given then the nz parameter is ignored 3507 3508 You can call MatGetInfo() to get information on how effective the preallocation was; 3509 for example the fields mallocs,nz_allocated,nz_used,nz_unneeded; 3510 You can also run with the option -info and look for messages with the string 3511 malloc in them to see if additional memory allocation was needed. 3512 3513 The block AIJ format is fully compatible with standard Fortran 77 3514 storage. That is, the stored row and column indices can begin at 3515 either one (as in Fortran) or zero. See the users' manual for details. 3516 3517 Specify the preallocated storage with either nz or nnz (not both). 3518 Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory 3519 allocation. See Users-Manual: ch_mat for details. 3520 3521 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateBAIJ(), MatGetInfo() 3522 @*/ 3523 PetscErrorCode MatSeqBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 3524 { 3525 PetscFunctionBegin; 3526 PetscValidHeaderSpecific(B,MAT_CLASSID,1); 3527 PetscValidType(B,1); 3528 PetscValidLogicalCollectiveInt(B,bs,2); 3529 PetscTryMethod(B,"MatSeqBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz)); 3530 PetscFunctionReturn(0); 3531 } 3532 3533 /*@C 3534 MatSeqBAIJSetPreallocationCSR - Creates a sparse parallel matrix in BAIJ format using the given nonzero structure and (optional) numerical values 3535 3536 Collective 3537 3538 Input Parameters: 3539 + B - the matrix 3540 . i - the indices into j for the start of each local row (starts with zero) 3541 . j - the column indices for each local row (starts with zero) these must be sorted for each row 3542 - v - optional values in the matrix 3543 3544 Level: advanced 3545 3546 Notes: 3547 The order of the entries in values is specified by the MatOption MAT_ROW_ORIENTED. For example, C programs 3548 may want to use the default MAT_ROW_ORIENTED=PETSC_TRUE and use an array v[nnz][bs][bs] where the second index is 3549 over rows within a block and the last index is over columns within a block row. Fortran programs will likely set 3550 MAT_ROW_ORIENTED=PETSC_FALSE and use a Fortran array v(bs,bs,nnz) in which the first index is over rows within a 3551 block column and the second index is over columns within a block. 3552 3553 Though this routine has Preallocation() in the name it also sets the exact nonzero locations of the matrix entries and usually the numerical values as well 3554 3555 .seealso: MatCreate(), MatCreateSeqBAIJ(), MatSetValues(), MatSeqBAIJSetPreallocation(), MATSEQBAIJ 3556 @*/ 3557 PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[]) 3558 { 3559 PetscFunctionBegin; 3560 PetscValidHeaderSpecific(B,MAT_CLASSID,1); 3561 PetscValidType(B,1); 3562 PetscValidLogicalCollectiveInt(B,bs,2); 3563 PetscTryMethod(B,"MatSeqBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v)); 3564 PetscFunctionReturn(0); 3565 } 3566 3567 /*@ 3568 MatCreateSeqBAIJWithArrays - Creates an sequential BAIJ matrix using matrix elements provided by the user. 3569 3570 Collective 3571 3572 Input Parameters: 3573 + comm - must be an MPI communicator of size 1 3574 . bs - size of block 3575 . m - number of rows 3576 . n - number of columns 3577 . i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row block row of the matrix 3578 . j - column indices 3579 - a - matrix values 3580 3581 Output Parameter: 3582 . mat - the matrix 3583 3584 Level: advanced 3585 3586 Notes: 3587 The i, j, and a arrays are not copied by this routine, the user must free these arrays 3588 once the matrix is destroyed 3589 3590 You cannot set new nonzero locations into this matrix, that will generate an error. 3591 3592 The i and j indices are 0 based 3593 3594 When block size is greater than 1 the matrix values must be stored using the BAIJ storage format (see the BAIJ code to determine this). 3595 3596 The order of the entries in values is the same as the block compressed sparse row storage format; that is, it is 3597 the same as a three dimensional array in Fortran values(bs,bs,nnz) that contains the first column of the first 3598 block, followed by the second column of the first block etc etc. That is, the blocks are contiguous in memory 3599 with column-major ordering within blocks. 3600 3601 .seealso: MatCreate(), MatCreateBAIJ(), MatCreateSeqBAIJ() 3602 3603 @*/ 3604 PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt i[],PetscInt j[],PetscScalar a[],Mat *mat) 3605 { 3606 PetscInt ii; 3607 Mat_SeqBAIJ *baij; 3608 3609 PetscFunctionBegin; 3610 PetscCheck(bs == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"block size %" PetscInt_FMT " > 1 is not supported yet",bs); 3611 if (m > 0) PetscCheck(i[0] == 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); 3612 3613 PetscCall(MatCreate(comm,mat)); 3614 PetscCall(MatSetSizes(*mat,m,n,m,n)); 3615 PetscCall(MatSetType(*mat,MATSEQBAIJ)); 3616 PetscCall(MatSeqBAIJSetPreallocation(*mat,bs,MAT_SKIP_ALLOCATION,NULL)); 3617 baij = (Mat_SeqBAIJ*)(*mat)->data; 3618 PetscCall(PetscMalloc2(m,&baij->imax,m,&baij->ilen)); 3619 PetscCall(PetscLogObjectMemory((PetscObject)*mat,2*m*sizeof(PetscInt))); 3620 3621 baij->i = i; 3622 baij->j = j; 3623 baij->a = a; 3624 3625 baij->singlemalloc = PETSC_FALSE; 3626 baij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 3627 baij->free_a = PETSC_FALSE; 3628 baij->free_ij = PETSC_FALSE; 3629 3630 for (ii=0; ii<m; ii++) { 3631 baij->ilen[ii] = baij->imax[ii] = i[ii+1] - i[ii]; 3632 PetscCheck(i[ii+1] - i[ii] >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %" PetscInt_FMT " length = %" PetscInt_FMT,ii,i[ii+1] - i[ii]); 3633 } 3634 if (PetscDefined(USE_DEBUG)) { 3635 for (ii=0; ii<baij->i[m]; ii++) { 3636 PetscCheck(j[ii] >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %" PetscInt_FMT " index = %" PetscInt_FMT,ii,j[ii]); 3637 PetscCheck(j[ii] <= n - 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %" PetscInt_FMT " index = %" PetscInt_FMT,ii,j[ii]); 3638 } 3639 } 3640 3641 PetscCall(MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY)); 3642 PetscCall(MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY)); 3643 PetscFunctionReturn(0); 3644 } 3645 3646 PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat) 3647 { 3648 PetscMPIInt size; 3649 3650 PetscFunctionBegin; 3651 PetscCallMPI(MPI_Comm_size(comm,&size)); 3652 if (size == 1 && scall == MAT_REUSE_MATRIX) { 3653 PetscCall(MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN)); 3654 } else { 3655 PetscCall(MatCreateMPIMatConcatenateSeqMat_MPIBAIJ(comm,inmat,n,scall,outmat)); 3656 } 3657 PetscFunctionReturn(0); 3658 } 3659