1 2 #include <../src/mat/impls/aij/mpi/mpiaij.h> 3 #include <../src/mat/impls/baij/mpi/mpibaij.h> 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "MatFDColoringApply_BAIJ" 7 PetscErrorCode MatFDColoringApply_BAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 8 { 9 PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; 10 PetscErrorCode ierr; 11 PetscInt k,cstart,cend,l,row,col,nz,spidx,i,j; 12 PetscScalar dx=0.0,*xx,*w3_array,*dy_i,*dy=coloring->dy; 13 PetscScalar *vscale_array; 14 PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; 15 Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; 16 void *fctx=coloring->fctx; 17 PetscInt ctype=coloring->ctype,nxloc,nrows_k; 18 PetscScalar *valaddr; 19 MatEntry *Jentry=coloring->matentry; 20 const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; 21 PetscInt bs=J->rmap->bs; 22 23 PetscFunctionBegin; 24 /* create vscale for storing dx */ 25 if (!vscale) { 26 if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 27 /* contain the "diagonal" on processor scalings followed by the off processor* - garray must be non-bloced! */ 28 Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)J->data; 29 PetscInt *garray; 30 ierr = PetscMalloc(baij->B->cmap->n*sizeof(PetscInt),&garray);CHKERRQ(ierr); 31 for (i=0; i<baij->B->cmap->n/bs; i++) { 32 for (j=0; j<bs; j++) { 33 garray[i*bs+j] = bs*baij->garray[i]+j; 34 } 35 } 36 ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,baij->B->cmap->n,garray,&vscale);CHKERRQ(ierr); 37 ierr = PetscFree(garray);CHKERRQ(ierr); 38 } else if (ctype == IS_COLORING_GHOSTED) { 39 ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr); 40 } 41 coloring->vscale = vscale; 42 } 43 44 /* (1) Set w1 = F(x1) */ 45 if (!coloring->fset) { 46 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 47 ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 48 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 49 } else { 50 coloring->fset = PETSC_FALSE; 51 } 52 53 /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 54 ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 55 if (coloring->htype[0] == 'w') { 56 /* vscale = dx is a constant scalar */ 57 ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 58 dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); 59 } else { 60 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 61 ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 62 for (col=0; col<nxloc; col++) { 63 dx = xx[col]; 64 if (PetscAbsScalar(dx) < umin) { 65 if (PetscRealPart(dx) >= 0.0) dx = umin; 66 else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 67 } 68 dx *= epsilon; 69 vscale_array[col] = 1.0/dx; 70 } 71 ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 72 ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 73 } 74 if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 75 ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 76 ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 77 } 78 79 /* (3) Loop over each color */ 80 if (!coloring->w3) { 81 ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 82 ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 83 } 84 w3 = coloring->w3; 85 86 ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 87 if (vscale) { 88 ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 89 } 90 nz = 0; 91 for (k=0; k<ncolors; k++) { 92 coloring->currentcolor = k; 93 94 /* 95 (3-1) Loop over each column associated with color 96 adding the perturbation to the vector w3 = x1 + dx. 97 */ 98 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 99 dy_i = dy; 100 for (i=0; i<bs; i++) { /* Loop over a block of columns */ 101 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 102 if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 103 if (coloring->htype[0] == 'w') { 104 for (l=0; l<ncolumns[k]; l++) { 105 col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 106 w3_array[col] += 1.0/dx; 107 if (i) w3_array[col-1] -= 1.0/dx; /* resume original w3[col-1] */ 108 } 109 } else { /* htype == 'ds' */ 110 vscale_array -= cstart; /* shift pointer so global index can be used */ 111 for (l=0; l<ncolumns[k]; l++) { 112 col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 113 w3_array[col] += 1.0/vscale_array[col]; 114 if (i) w3_array[col-1] -= 1.0/vscale_array[col-1]; /* resume original w3[col-1] */ 115 } 116 vscale_array += cstart; 117 } 118 if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 119 ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 120 121 /* 122 (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 123 w2 = F(x1 + dx) - F(x1) 124 */ 125 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 126 ierr = VecPlaceArray(w2,dy_i);CHKERRQ(ierr); /* place w2 to the array dy_i */ 127 ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 128 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 129 ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 130 ierr = VecResetArray(w2);CHKERRQ(ierr); 131 dy_i += nxloc; /* points to dy+i*nxloc */ 132 } 133 134 /* 135 (3-3) Loop over rows of vector, putting results into Jacobian matrix 136 */ 137 nrows_k = nrows[k]; 138 if (coloring->htype[0] == 'w') { 139 for (l=0; l<nrows_k; l++) { 140 row = bs*Jentry[nz].row; /* local row index */ 141 valaddr = Jentry[nz++].valaddr; 142 spidx = 0; 143 dy_i = dy; 144 for (i=0; i<bs; i++) { /* column of the block */ 145 for (j=0; j<bs; j++) { /* row of the block */ 146 valaddr[spidx++] = dy_i[row+j]*dx; 147 } 148 dy_i += nxloc; /* points to dy+i*nxloc */ 149 } 150 } 151 } else { /* htype == 'ds' */ 152 for (l=0; l<nrows_k; l++) { 153 row = bs*Jentry[nz].row; /* local row index */ 154 col = bs*Jentry[nz].col; /* local column index */ 155 valaddr = Jentry[nz++].valaddr; 156 spidx = 0; 157 dy_i = dy; 158 for (i=0; i<bs; i++) { /* column of the block */ 159 for (j=0; j<bs; j++) { /* row of the block */ 160 valaddr[spidx++] = dy_i[row+j]*vscale_array[col+i]; 161 } 162 dy_i += nxloc; /* points to dy+i*nxloc */ 163 } 164 } 165 } 166 } 167 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169 if (vscale) { 170 ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 171 } 172 173 coloring->currentcolor = -1; 174 PetscFunctionReturn(0); 175 } 176 177 #undef __FUNCT__ 178 #define __FUNCT__ "MatFDColoringApply_AIJ" 179 PetscErrorCode MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 180 { 181 PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; 182 PetscErrorCode ierr; 183 PetscInt k,cstart,cend,l,row,col,nz; 184 PetscScalar dx=0.0,*y,*xx,*w3_array; 185 PetscScalar *vscale_array; 186 PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; 187 Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; 188 void *fctx=coloring->fctx; 189 PetscInt ctype=coloring->ctype,nxloc,nrows_k; 190 MatEntry *Jentry=coloring->matentry; 191 const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; 192 193 PetscFunctionBegin; 194 /* create vscale for storing dx */ 195 if (!vscale) { 196 if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 197 Mat_MPIAIJ *aij=(Mat_MPIAIJ*)J->data; 198 ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr); 199 } else if (ctype == IS_COLORING_GHOSTED) { 200 ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr); 201 } 202 coloring->vscale = vscale; 203 } 204 205 /* (1) Set w1 = F(x1) */ 206 if (!coloring->fset) { 207 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 208 ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 209 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 210 } else { 211 coloring->fset = PETSC_FALSE; 212 } 213 214 /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 215 if (coloring->htype[0] == 'w') { 216 /* vscale = 1./dx is a constant scalar */ 217 ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 218 dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); 219 } else { 220 ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 221 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 222 ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 223 for (col=0; col<nxloc; col++) { 224 dx = xx[col]; 225 if (PetscAbsScalar(dx) < umin) { 226 if (PetscRealPart(dx) >= 0.0) dx = umin; 227 else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 228 } 229 dx *= epsilon; 230 vscale_array[col] = 1.0/dx; 231 } 232 ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 233 ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 234 } 235 if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 236 ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 237 ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 238 } 239 240 /* (3) Loop over each color */ 241 if (!coloring->w3) { 242 ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 243 ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 244 } 245 w3 = coloring->w3; 246 247 ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 248 if (vscale) { 249 ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 250 } 251 nz = 0; 252 253 if (coloring->bcols > 1) { /* use blocked insertion of Jentry */ 254 PetscInt i,m=J->rmap->n,nbcols,bcols=coloring->bcols; 255 PetscScalar *dy=coloring->dy,*dy_k; 256 257 nbcols = 0; 258 for (k=0; k<ncolors; k+=bcols) { 259 coloring->currentcolor = k; 260 261 /* 262 (3-1) Loop over each column associated with color 263 adding the perturbation to the vector w3 = x1 + dx. 264 */ 265 266 dy_k = dy; 267 if (k + bcols > ncolors) bcols = ncolors - k; 268 for (i=0; i<bcols; i++) { 269 270 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 271 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 272 if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 273 if (coloring->htype[0] == 'w') { 274 for (l=0; l<ncolumns[k+i]; l++) { 275 col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */ 276 w3_array[col] += 1.0/dx; 277 } 278 } else { /* htype == 'ds' */ 279 vscale_array -= cstart; /* shift pointer so global index can be used */ 280 for (l=0; l<ncolumns[k+i]; l++) { 281 col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */ 282 w3_array[col] += 1.0/vscale_array[col]; 283 } 284 vscale_array += cstart; 285 } 286 if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 287 ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 288 289 /* 290 (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 291 w2 = F(x1 + dx) - F(x1) 292 */ 293 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 294 ierr = VecPlaceArray(w2,dy_k);CHKERRQ(ierr); /* place w2 to the array dy_i */ 295 ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 296 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 297 ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 298 ierr = VecResetArray(w2);CHKERRQ(ierr); 299 dy_k += m; /* points to dy+i*nxloc */ 300 301 } 302 303 /* 304 (3-3) Loop over block rows of vector, putting results into Jacobian matrix 305 */ 306 nrows_k = nrows[nbcols++]; 307 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 308 309 if (coloring->htype[0] == 'w') { 310 for (l=0; l<nrows_k; l++) { 311 row = Jentry[nz].row; /* local row index */ 312 *(Jentry[nz++].valaddr) = dy[row]*dx; 313 } 314 } else { /* htype == 'ds' */ 315 for (l=0; l<nrows_k; l++) { 316 row = Jentry[nz].row; /* local row index */ 317 *(Jentry[nz].valaddr) = dy[row]*vscale_array[Jentry[nz].col]; 318 nz++; 319 } 320 } 321 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 322 } 323 } else { /* bcols == 1 */ 324 for (k=0; k<ncolors; k++) { 325 coloring->currentcolor = k; 326 327 /* 328 (3-1) Loop over each column associated with color 329 adding the perturbation to the vector w3 = x1 + dx. 330 */ 331 ierr = VecCopy(x1,w3);CHKERRQ(ierr); 332 ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 333 if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 334 if (coloring->htype[0] == 'w') { 335 for (l=0; l<ncolumns[k]; l++) { 336 col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 337 w3_array[col] += 1.0/dx; 338 } 339 } else { /* htype == 'ds' */ 340 vscale_array -= cstart; /* shift pointer so global index can be used */ 341 for (l=0; l<ncolumns[k]; l++) { 342 col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 343 w3_array[col] += 1.0/vscale_array[col]; 344 } 345 vscale_array += cstart; 346 } 347 if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 348 ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 349 350 /* 351 (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 352 w2 = F(x1 + dx) - F(x1) 353 */ 354 ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 355 ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 356 ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 357 ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 358 359 /* 360 (3-3) Loop over rows of vector, putting results into Jacobian matrix 361 */ 362 nrows_k = nrows[k]; 363 ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 364 if (coloring->htype[0] == 'w') { 365 for (l=0; l<nrows_k; l++) { 366 row = Jentry[nz].row; /* local row index */ 367 *(Jentry[nz++].valaddr) = y[row]*dx; 368 } 369 } else { /* htype == 'ds' */ 370 for (l=0; l<nrows_k; l++) { 371 row = Jentry[nz].row; /* local row index */ 372 *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col]; 373 nz++; 374 } 375 } 376 ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 377 } 378 } 379 380 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 381 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 382 if (vscale) { 383 ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 384 } 385 coloring->currentcolor = -1; 386 PetscFunctionReturn(0); 387 } 388 389 #undef __FUNCT__ 390 #define __FUNCT__ "MatFDColoringSetUp_MPIXAIJ" 391 PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 392 { 393 PetscErrorCode ierr; 394 PetscMPIInt size,*ncolsonproc,*disp,nn; 395 PetscInt i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb; 396 const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL; 397 PetscInt nis=iscoloring->n,nctot,*cols; 398 IS *isa; 399 ISLocalToGlobalMapping map=mat->cmap->mapping; 400 PetscInt ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx; 401 Mat A,B; 402 PetscScalar *A_val,*B_val,**valaddrhit; 403 MatEntry *Jentry; 404 PetscBool isBAIJ; 405 PetscInt bcols=c->bcols; 406 #if defined(PETSC_USE_CTABLE) 407 PetscTable colmap=NULL; 408 #else 409 PetscInt *colmap=NULL; /* local col number of off-diag col */ 410 #endif 411 412 PetscFunctionBegin; 413 if (ctype == IS_COLORING_GHOSTED) { 414 if (!map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); 415 ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr); 416 } 417 418 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 419 ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr); 420 if (isBAIJ) { 421 Mat_MPIBAIJ *aij=(Mat_MPIBAIJ*)mat->data; 422 Mat_SeqBAIJ *spA,*spB; 423 A = aij->A; spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a; 424 B = aij->B; spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a; 425 nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 426 if (!aij->colmap) { 427 ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr); 428 colmap = aij->colmap; 429 } 430 ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 431 ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 432 } else { 433 Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data; 434 Mat_SeqAIJ *spA,*spB; 435 A = aij->A; spA = (Mat_SeqAIJ*)A->data; A_val = spA->a; 436 B = aij->B; spB = (Mat_SeqAIJ*)B->data; B_val = spB->a; 437 nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 438 if (!aij->colmap) { 439 /* Allow access to data structures of local part of matrix 440 - creates aij->colmap which maps global column number to local number in part B */ 441 ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 442 colmap = aij->colmap; 443 } 444 ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 445 ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 446 447 bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */ 448 } 449 450 m = mat->rmap->n/bs; 451 cstart = mat->cmap->rstart/bs; 452 cend = mat->cmap->rend/bs; 453 454 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 455 ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 456 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 457 ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 458 459 ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 460 ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 461 c->matentry = Jentry; 462 463 ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 464 nz = 0; 465 ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 466 for (i=0; i<nis; i++) { /* for each local color */ 467 ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 468 ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 469 470 c->ncolumns[i] = n; /* local number of columns of this color on this process */ 471 if (n) { 472 ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 473 ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 474 ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 475 } else { 476 c->columns[i] = 0; 477 } 478 479 if (ctype == IS_COLORING_GLOBAL) { 480 /* Determine nctot, the total (parallel) number of columns of this color */ 481 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 482 ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 483 484 /* ncolsonproc[j]: local ncolumns on proc[j] of this color */ 485 ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 486 ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 487 nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 488 if (!nctot) { 489 ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 490 } 491 492 disp[0] = 0; 493 for (j=1; j<size; j++) { 494 disp[j] = disp[j-1] + ncolsonproc[j-1]; 495 } 496 497 /* Get cols, the complete list of columns for this color on each process */ 498 ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 499 ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 500 ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 501 } else if (ctype == IS_COLORING_GHOSTED) { 502 /* Determine local number of columns of this color on this process, including ghost points */ 503 nctot = n; 504 ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 505 ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 506 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 507 508 /* Mark all rows affect by these columns */ 509 ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr); 510 bs2 = bs*bs; 511 nrows_i = 0; 512 for (j=0; j<nctot; j++) { /* loop over columns*/ 513 if (ctype == IS_COLORING_GHOSTED) { 514 col = ltog[cols[j]]; 515 } else { 516 col = cols[j]; 517 } 518 if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */ 519 row = A_cj + A_ci[col-cstart]; 520 nrows = A_ci[col-cstart+1] - A_ci[col-cstart]; 521 nrows_i += nrows; 522 /* loop over columns of A marking them in rowhit */ 523 for (k=0; k<nrows; k++) { 524 /* set valaddrhit for part A */ 525 spidx = bs2*spidxA[A_ci[col-cstart] + k]; 526 valaddrhit[*row] = &A_val[spidx]; 527 rowhit[*row++] = col - cstart + 1; /* local column index */ 528 } 529 } else { /* column is in B, off-diagonal block of mat */ 530 #if defined(PETSC_USE_CTABLE) 531 ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr); 532 colb--; 533 #else 534 colb = colmap[col] - 1; /* local column index */ 535 #endif 536 if (colb == -1) { 537 nrows = 0; 538 } else { 539 colb = colb/bs; 540 row = B_cj + B_ci[colb]; 541 nrows = B_ci[colb+1] - B_ci[colb]; 542 } 543 nrows_i += nrows; 544 /* loop over columns of B marking them in rowhit */ 545 for (k=0; k<nrows; k++) { 546 /* set valaddrhit for part B */ 547 spidx = bs2*spidxB[B_ci[colb] + k]; 548 valaddrhit[*row] = &B_val[spidx]; 549 rowhit[*row++] = colb + 1 + cend - cstart; /* local column index */ 550 } 551 } 552 } 553 c->nrows[i] = nrows_i; 554 555 for (j=0; j<m; j++) { 556 if (rowhit[j]) { 557 Jentry[nz].row = j; /* local row index */ 558 Jentry[nz].col = rowhit[j] - 1; /* local column index */ 559 Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 560 nz++; 561 } 562 } 563 ierr = PetscFree(cols);CHKERRQ(ierr); 564 } 565 566 if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */ 567 ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr); 568 } 569 570 if (isBAIJ) { 571 ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 572 ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 573 ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 574 } else { 575 ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 576 ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 577 } 578 579 ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 580 ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 581 582 if (ctype == IS_COLORING_GHOSTED) { 583 ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr); 584 } 585 #if defined(PETSC_USE_INFO) 586 ierr = PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr); 587 #endif 588 PetscFunctionReturn(0); 589 } 590 591 #undef __FUNCT__ 592 #define __FUNCT__ "MatFDColoringCreate_MPIXAIJ" 593 PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 594 { 595 PetscErrorCode ierr; 596 PetscInt bs,nis=iscoloring->n,m=mat->rmap->n; 597 PetscBool isBAIJ; 598 599 PetscFunctionBegin; 600 /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian; 601 bcols is chosen s.t. dy-array takes 50% of memory space as mat */ 602 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 603 ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr); 604 if (isBAIJ) { 605 c->brows = m; 606 c->bcols = 1; 607 } else { /* mpiaij matrix */ 608 /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */ 609 Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data; 610 Mat_SeqAIJ *spA,*spB; 611 Mat A,B; 612 PetscInt nz,brows,bcols; 613 PetscReal mem; 614 615 bs = 1; /* only bs=1 is supported for MPIAIJ matrix */ 616 617 A = aij->A; spA = (Mat_SeqAIJ*)A->data; 618 B = aij->B; spB = (Mat_SeqAIJ*)B->data; 619 nz = spA->nz + spB->nz; /* total local nonzero entries of mat */ 620 mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt); 621 bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar))); 622 brows = 1000/bcols; 623 if (bcols > nis) bcols = nis; 624 if (brows == 0 || brows > m) brows = m; 625 c->brows = brows; 626 c->bcols = bcols; 627 } 628 629 c->M = mat->rmap->N/bs; /* set the global rows and columns and local rows */ 630 c->N = mat->cmap->N/bs; 631 c->m = mat->rmap->n/bs; 632 c->rstart = mat->rmap->rstart/bs; 633 c->ncolors = nis; 634 PetscFunctionReturn(0); 635 } 636