1a64fbb32SBarry Smith 2c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 30d1c53f1SHong Zhang #include <../src/mat/impls/baij/mpi/mpibaij.h> 40d1c53f1SHong Zhang 50d1c53f1SHong Zhang #undef __FUNCT__ 6040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringApply_BAIJ" 7040ebd07SHong Zhang PetscErrorCode MatFDColoringApply_BAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 80d1c53f1SHong Zhang { 90d1c53f1SHong Zhang PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; 100d1c53f1SHong Zhang PetscErrorCode ierr; 110d1c53f1SHong Zhang PetscInt k,cstart,cend,l,row,col,nz,spidx,i,j; 12e2c857f8SHong Zhang PetscScalar dx=0.0,*xx,*w3_array,*dy_i,*dy=coloring->dy; 130d1c53f1SHong Zhang PetscScalar *vscale_array; 140d1c53f1SHong Zhang PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; 150d1c53f1SHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; 160d1c53f1SHong Zhang void *fctx=coloring->fctx; 170d1c53f1SHong Zhang PetscInt ctype=coloring->ctype,nxloc,nrows_k; 180d1c53f1SHong Zhang PetscScalar *valaddr; 190d1c53f1SHong Zhang MatEntry *Jentry=coloring->matentry; 200df34763SHong Zhang MatEntry2 *Jentry2=coloring->matentry2; 210d1c53f1SHong Zhang const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; 220d1c53f1SHong Zhang PetscInt bs=J->rmap->bs; 230d1c53f1SHong Zhang 240d1c53f1SHong Zhang PetscFunctionBegin; 250d1c53f1SHong Zhang /* (1) Set w1 = F(x1) */ 260d1c53f1SHong Zhang if (!coloring->fset) { 270d1c53f1SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 280d1c53f1SHong Zhang ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 290d1c53f1SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 300d1c53f1SHong Zhang } else { 310d1c53f1SHong Zhang coloring->fset = PETSC_FALSE; 320d1c53f1SHong Zhang } 330d1c53f1SHong Zhang 340d1c53f1SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 350d1c53f1SHong Zhang ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 360d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 370d1c53f1SHong Zhang /* vscale = dx is a constant scalar */ 380d1c53f1SHong Zhang ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 390d1c53f1SHong Zhang dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); 400d1c53f1SHong Zhang } else { 410d1c53f1SHong Zhang ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 420d1c53f1SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 430d1c53f1SHong Zhang for (col=0; col<nxloc; col++) { 440d1c53f1SHong Zhang dx = xx[col]; 450d1c53f1SHong Zhang if (PetscAbsScalar(dx) < umin) { 460d1c53f1SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 470d1c53f1SHong Zhang else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 480d1c53f1SHong Zhang } 490d1c53f1SHong Zhang dx *= epsilon; 500d1c53f1SHong Zhang vscale_array[col] = 1.0/dx; 510d1c53f1SHong Zhang } 520d1c53f1SHong Zhang ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 530d1c53f1SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 540d1c53f1SHong Zhang } 55684f2004SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 560d1c53f1SHong Zhang ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 570d1c53f1SHong Zhang ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 580d1c53f1SHong Zhang } 590d1c53f1SHong Zhang 600d1c53f1SHong Zhang /* (3) Loop over each color */ 610d1c53f1SHong Zhang if (!coloring->w3) { 620d1c53f1SHong Zhang ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 630d1c53f1SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 640d1c53f1SHong Zhang } 650d1c53f1SHong Zhang w3 = coloring->w3; 660d1c53f1SHong Zhang 670d1c53f1SHong Zhang ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 680d1c53f1SHong Zhang if (vscale) { 690d1c53f1SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 700d1c53f1SHong Zhang } 710d1c53f1SHong Zhang nz = 0; 720d1c53f1SHong Zhang for (k=0; k<ncolors; k++) { 730d1c53f1SHong Zhang coloring->currentcolor = k; 740d1c53f1SHong Zhang 750d1c53f1SHong Zhang /* 760d1c53f1SHong Zhang (3-1) Loop over each column associated with color 770d1c53f1SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 780d1c53f1SHong Zhang */ 790d1c53f1SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 800d1c53f1SHong Zhang dy_i = dy; 810d1c53f1SHong Zhang for (i=0; i<bs; i++) { /* Loop over a block of columns */ 820d1c53f1SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 830d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 840d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 850d1c53f1SHong Zhang for (l=0; l<ncolumns[k]; l++) { 860d1c53f1SHong Zhang col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 870d1c53f1SHong Zhang w3_array[col] += 1.0/dx; 88684f2004SHong Zhang if (i) w3_array[col-1] -= 1.0/dx; /* resume original w3[col-1] */ 890d1c53f1SHong Zhang } 900d1c53f1SHong Zhang } else { /* htype == 'ds' */ 910d1c53f1SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 920d1c53f1SHong Zhang for (l=0; l<ncolumns[k]; l++) { 93f8c2866eSHong Zhang col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 940d1c53f1SHong Zhang w3_array[col] += 1.0/vscale_array[col]; 95684f2004SHong Zhang if (i) w3_array[col-1] -= 1.0/vscale_array[col-1]; /* resume original w3[col-1] */ 960d1c53f1SHong Zhang } 970d1c53f1SHong Zhang vscale_array += cstart; 980d1c53f1SHong Zhang } 990d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 1000d1c53f1SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 1010d1c53f1SHong Zhang 1020d1c53f1SHong Zhang /* 1030d1c53f1SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 1040d1c53f1SHong Zhang w2 = F(x1 + dx) - F(x1) 1050d1c53f1SHong Zhang */ 1060d1c53f1SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 1070d1c53f1SHong Zhang ierr = VecPlaceArray(w2,dy_i);CHKERRQ(ierr); /* place w2 to the array dy_i */ 1080d1c53f1SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 1090d1c53f1SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 1100d1c53f1SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 1110d1c53f1SHong Zhang ierr = VecResetArray(w2);CHKERRQ(ierr); 1120d1c53f1SHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1130d1c53f1SHong Zhang } 1140d1c53f1SHong Zhang 1150d1c53f1SHong Zhang /* 1160d1c53f1SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 1170d1c53f1SHong Zhang */ 1180d1c53f1SHong Zhang nrows_k = nrows[k]; 1190d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 1200d1c53f1SHong Zhang for (l=0; l<nrows_k; l++) { 1210df34763SHong Zhang row = bs*Jentry2[nz].row; /* local row index */ 1220df34763SHong Zhang valaddr = Jentry2[nz++].valaddr; 1230d1c53f1SHong Zhang spidx = 0; 1240d1c53f1SHong Zhang dy_i = dy; 1250d1c53f1SHong Zhang for (i=0; i<bs; i++) { /* column of the block */ 1260d1c53f1SHong Zhang for (j=0; j<bs; j++) { /* row of the block */ 1270d1c53f1SHong Zhang valaddr[spidx++] = dy_i[row+j]*dx; 1280d1c53f1SHong Zhang } 129f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1300d1c53f1SHong Zhang } 1310d1c53f1SHong Zhang } 1320d1c53f1SHong Zhang } else { /* htype == 'ds' */ 1330d1c53f1SHong Zhang for (l=0; l<nrows_k; l++) { 1340d1c53f1SHong Zhang row = bs*Jentry[nz].row; /* local row index */ 1350d1c53f1SHong Zhang col = bs*Jentry[nz].col; /* local column index */ 136684f2004SHong Zhang valaddr = Jentry[nz++].valaddr; 137f8c2866eSHong Zhang spidx = 0; 138f8c2866eSHong Zhang dy_i = dy; 139f8c2866eSHong Zhang for (i=0; i<bs; i++) { /* column of the block */ 140f8c2866eSHong Zhang for (j=0; j<bs; j++) { /* row of the block */ 141f8c2866eSHong Zhang valaddr[spidx++] = dy_i[row+j]*vscale_array[col+i]; 142f8c2866eSHong Zhang } 143f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 144f8c2866eSHong Zhang } 1450d1c53f1SHong Zhang } 1460d1c53f1SHong Zhang } 1470d1c53f1SHong Zhang } 1480d1c53f1SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1490d1c53f1SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1500d1c53f1SHong Zhang if (vscale) { 1510d1c53f1SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 1520d1c53f1SHong Zhang } 1530d1c53f1SHong Zhang 1540d1c53f1SHong Zhang coloring->currentcolor = -1; 1550d1c53f1SHong Zhang PetscFunctionReturn(0); 1560d1c53f1SHong Zhang } 157a64fbb32SBarry Smith 1584a2ae208SSatish Balay #undef __FUNCT__ 159040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringApply_AIJ" 160040ebd07SHong Zhang PetscErrorCode MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 161fcd7ac73SHong Zhang { 162fcd7ac73SHong Zhang PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; 163fcd7ac73SHong Zhang PetscErrorCode ierr; 164a2f2d239SHong Zhang PetscInt k,cstart,cend,l,row,col,nz; 1659e917edbSHong Zhang PetscScalar dx=0.0,*y,*xx,*w3_array; 166fcd7ac73SHong Zhang PetscScalar *vscale_array; 167fcd7ac73SHong Zhang PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; 1688bc97078SHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; 169fcd7ac73SHong Zhang void *fctx=coloring->fctx; 170c53567a0SHong Zhang PetscInt ctype=coloring->ctype,nxloc,nrows_k; 171573f477fSHong Zhang MatEntry *Jentry=coloring->matentry; 1720df34763SHong Zhang MatEntry2 *Jentry2=coloring->matentry2; 1738bc97078SHong Zhang const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; 174fcd7ac73SHong Zhang 175fcd7ac73SHong Zhang PetscFunctionBegin; 1768bc97078SHong Zhang /* (1) Set w1 = F(x1) */ 177fcd7ac73SHong Zhang if (!coloring->fset) { 178fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 179f6af9589SHong Zhang ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 180fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 181fcd7ac73SHong Zhang } else { 182fcd7ac73SHong Zhang coloring->fset = PETSC_FALSE; 183fcd7ac73SHong Zhang } 184fcd7ac73SHong Zhang 1858bc97078SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 186f6af9589SHong Zhang if (coloring->htype[0] == 'w') { 187684f2004SHong Zhang /* vscale = 1./dx is a constant scalar */ 188f6af9589SHong Zhang ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 189c53567a0SHong Zhang dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); 19070e7395fSHong Zhang } else { 19174d3cef9SHong Zhang ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 192f6af9589SHong Zhang ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 1938bc97078SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 19474d3cef9SHong Zhang for (col=0; col<nxloc; col++) { 195fcd7ac73SHong Zhang dx = xx[col]; 19674d3cef9SHong Zhang if (PetscAbsScalar(dx) < umin) { 19774d3cef9SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 19874d3cef9SHong Zhang else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 19974d3cef9SHong Zhang } 200fcd7ac73SHong Zhang dx *= epsilon; 20174d3cef9SHong Zhang vscale_array[col] = 1.0/dx; 202f6af9589SHong Zhang } 20374d3cef9SHong Zhang ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 2048bc97078SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 20570e7395fSHong Zhang } 206684f2004SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 2078bc97078SHong Zhang ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2088bc97078SHong Zhang ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 209fcd7ac73SHong Zhang } 210fcd7ac73SHong Zhang 2118bc97078SHong Zhang /* (3) Loop over each color */ 2128bc97078SHong Zhang if (!coloring->w3) { 2138bc97078SHong Zhang ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 2148bc97078SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 2158bc97078SHong Zhang } 2168bc97078SHong Zhang w3 = coloring->w3; 217fcd7ac73SHong Zhang 2188bc97078SHong Zhang ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 2199e917edbSHong Zhang if (vscale) { 2208bc97078SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 2219e917edbSHong Zhang } 2228bc97078SHong Zhang nz = 0; 223e2c857f8SHong Zhang 224b3e1f37bSHong Zhang if (coloring->bcols > 1) { /* use blocked insertion of Jentry */ 225b3e1f37bSHong Zhang PetscInt i,m=J->rmap->n,nbcols,bcols=coloring->bcols; 226e2c857f8SHong Zhang PetscScalar *dy=coloring->dy,*dy_k; 227e2c857f8SHong Zhang 228e2c857f8SHong Zhang nbcols = 0; 229e2c857f8SHong Zhang for (k=0; k<ncolors; k+=bcols) { 230e2c857f8SHong Zhang coloring->currentcolor = k; 231e2c857f8SHong Zhang 232e2c857f8SHong Zhang /* 233e2c857f8SHong Zhang (3-1) Loop over each column associated with color 234e2c857f8SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 235e2c857f8SHong Zhang */ 236e2c857f8SHong Zhang 237e2c857f8SHong Zhang dy_k = dy; 238e2c857f8SHong Zhang if (k + bcols > ncolors) bcols = ncolors - k; 239e2c857f8SHong Zhang for (i=0; i<bcols; i++) { 240e2c857f8SHong Zhang 241e2c857f8SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 242e2c857f8SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 243e2c857f8SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 244e2c857f8SHong Zhang if (coloring->htype[0] == 'w') { 245e2c857f8SHong Zhang for (l=0; l<ncolumns[k+i]; l++) { 246e2c857f8SHong Zhang col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */ 247e2c857f8SHong Zhang w3_array[col] += 1.0/dx; 248e2c857f8SHong Zhang } 249e2c857f8SHong Zhang } else { /* htype == 'ds' */ 250e2c857f8SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 251e2c857f8SHong Zhang for (l=0; l<ncolumns[k+i]; l++) { 252e2c857f8SHong Zhang col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */ 253e2c857f8SHong Zhang w3_array[col] += 1.0/vscale_array[col]; 254e2c857f8SHong Zhang } 255e2c857f8SHong Zhang vscale_array += cstart; 256e2c857f8SHong Zhang } 257e2c857f8SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 258e2c857f8SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 259e2c857f8SHong Zhang 260e2c857f8SHong Zhang /* 261e2c857f8SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 262e2c857f8SHong Zhang w2 = F(x1 + dx) - F(x1) 263e2c857f8SHong Zhang */ 264e2c857f8SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 265e2c857f8SHong Zhang ierr = VecPlaceArray(w2,dy_k);CHKERRQ(ierr); /* place w2 to the array dy_i */ 266e2c857f8SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 267e2c857f8SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 268e2c857f8SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 269e2c857f8SHong Zhang ierr = VecResetArray(w2);CHKERRQ(ierr); 270e2c857f8SHong Zhang dy_k += m; /* points to dy+i*nxloc */ 271e2c857f8SHong Zhang } 272e2c857f8SHong Zhang 273e2c857f8SHong Zhang /* 274e2c857f8SHong Zhang (3-3) Loop over block rows of vector, putting results into Jacobian matrix 275e2c857f8SHong Zhang */ 276d880da65SHong Zhang nrows_k = nrows[nbcols++]; 277e2c857f8SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 278d880da65SHong Zhang 279e2c857f8SHong Zhang if (coloring->htype[0] == 'w') { 280e2c857f8SHong Zhang for (l=0; l<nrows_k; l++) { 2810df34763SHong Zhang row = Jentry2[nz].row; /* local row index */ 2820df34763SHong Zhang *(Jentry2[nz++].valaddr) = dy[row]*dx; 283e2c857f8SHong Zhang } 284e2c857f8SHong Zhang } else { /* htype == 'ds' */ 285e2c857f8SHong Zhang for (l=0; l<nrows_k; l++) { 286e2c857f8SHong Zhang row = Jentry[nz].row; /* local row index */ 287d880da65SHong Zhang *(Jentry[nz].valaddr) = dy[row]*vscale_array[Jentry[nz].col]; 288e2c857f8SHong Zhang nz++; 289e2c857f8SHong Zhang } 290e2c857f8SHong Zhang } 291e2c857f8SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 292e2c857f8SHong Zhang } 293b3e1f37bSHong Zhang } else { /* bcols == 1 */ 2948bc97078SHong Zhang for (k=0; k<ncolors; k++) { 2958bc97078SHong Zhang coloring->currentcolor = k; 2969e917edbSHong Zhang 297fcd7ac73SHong Zhang /* 2988bc97078SHong Zhang (3-1) Loop over each column associated with color 299f6af9589SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 300fcd7ac73SHong Zhang */ 301f6af9589SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 302f6af9589SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 303a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 304b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 3058bc97078SHong Zhang for (l=0; l<ncolumns[k]; l++) { 306f6af9589SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 3079e917edbSHong Zhang w3_array[col] += 1.0/dx; 3089e917edbSHong Zhang } 309b039d6c4SHong Zhang } else { /* htype == 'ds' */ 310a2f2d239SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 311b039d6c4SHong Zhang for (l=0; l<ncolumns[k]; l++) { 312b039d6c4SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 313a2f2d239SHong Zhang w3_array[col] += 1.0/vscale_array[col]; 31470e7395fSHong Zhang } 315a2f2d239SHong Zhang vscale_array += cstart; 316fcd7ac73SHong Zhang } 317a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 318fcd7ac73SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 3199e917edbSHong Zhang 320fcd7ac73SHong Zhang /* 3218bc97078SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 322fcd7ac73SHong Zhang w2 = F(x1 + dx) - F(x1) 323fcd7ac73SHong Zhang */ 324fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 325fcd7ac73SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 326fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 327fcd7ac73SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 3289e917edbSHong Zhang 3298bc97078SHong Zhang /* 3308bc97078SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 3318bc97078SHong Zhang */ 332c53567a0SHong Zhang nrows_k = nrows[k]; 333fcd7ac73SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 334b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 335c53567a0SHong Zhang for (l=0; l<nrows_k; l++) { 3360df34763SHong Zhang row = Jentry2[nz].row; /* local row index */ 3370df34763SHong Zhang *(Jentry2[nz++].valaddr) = y[row]*dx; 3389e917edbSHong Zhang } 339b039d6c4SHong Zhang } else { /* htype == 'ds' */ 340c53567a0SHong Zhang for (l=0; l<nrows_k; l++) { 341b039d6c4SHong Zhang row = Jentry[nz].row; /* local row index */ 342b039d6c4SHong Zhang *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col]; 343573f477fSHong Zhang nz++; 344fcd7ac73SHong Zhang } 345b039d6c4SHong Zhang } 346fcd7ac73SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 3478bc97078SHong Zhang } 348b3e1f37bSHong Zhang } 349b3e1f37bSHong Zhang 350f5aae955SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 351f5aae955SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3529e917edbSHong Zhang if (vscale) { 3538bc97078SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 3549e917edbSHong Zhang } 355fcd7ac73SHong Zhang coloring->currentcolor = -1; 356fcd7ac73SHong Zhang PetscFunctionReturn(0); 357fcd7ac73SHong Zhang } 358fcd7ac73SHong Zhang 359fcd7ac73SHong Zhang #undef __FUNCT__ 360f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetUp_MPIXAIJ" 361f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 362fcd7ac73SHong Zhang { 363fcd7ac73SHong Zhang PetscErrorCode ierr; 364fcd7ac73SHong Zhang PetscMPIInt size,*ncolsonproc,*disp,nn; 365e3225e9dSHong Zhang PetscInt i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb; 36672c15787SHong Zhang const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL; 367fcd7ac73SHong Zhang PetscInt nis=iscoloring->n,nctot,*cols; 368fcd7ac73SHong Zhang IS *isa; 369fcd7ac73SHong Zhang ISLocalToGlobalMapping map=mat->cmap->mapping; 370a8971b87SHong Zhang PetscInt ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx; 3710d1c53f1SHong Zhang Mat A,B; 372e3225e9dSHong Zhang PetscScalar *A_val,*B_val,**valaddrhit; 373573f477fSHong Zhang MatEntry *Jentry; 3740df34763SHong Zhang MatEntry2 *Jentry2; 3750d1c53f1SHong Zhang PetscBool isBAIJ; 376531e53bdSHong Zhang PetscInt bcols=c->bcols; 3770d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE) 3780d1c53f1SHong Zhang PetscTable colmap=NULL; 3790d1c53f1SHong Zhang #else 3800d1c53f1SHong Zhang PetscInt *colmap=NULL; /* local col number of off-diag col */ 3810d1c53f1SHong Zhang #endif 382fcd7ac73SHong Zhang 383fcd7ac73SHong Zhang PetscFunctionBegin; 38472c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 38572c15787SHong Zhang if (!map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); 38672c15787SHong Zhang ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr); 38772c15787SHong Zhang } 388fcd7ac73SHong Zhang 3890d1c53f1SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 3900d1c53f1SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr); 391e3225e9dSHong Zhang if (isBAIJ) { 392195687c5SHong Zhang Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)mat->data; 3936c145f34SHong Zhang Mat_SeqBAIJ *spA,*spB; 394195687c5SHong Zhang A = baij->A; spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a; 395195687c5SHong Zhang B = baij->B; spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a; 3960d1c53f1SHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 397195687c5SHong Zhang if (!baij->colmap) { 3980d1c53f1SHong Zhang ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr); 399195687c5SHong Zhang colmap = baij->colmap; 4000d1c53f1SHong Zhang } 4010d1c53f1SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 4020d1c53f1SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 403195687c5SHong Zhang 404195687c5SHong Zhang if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */ 405195687c5SHong Zhang PetscInt *garray; 406195687c5SHong Zhang ierr = PetscMalloc(B->cmap->n*sizeof(PetscInt),&garray);CHKERRQ(ierr); 407195687c5SHong Zhang for (i=0; i<baij->B->cmap->n/bs; i++) { 408195687c5SHong Zhang for (j=0; j<bs; j++) { 409195687c5SHong Zhang garray[i*bs+j] = bs*baij->garray[i]+j; 410195687c5SHong Zhang } 411195687c5SHong Zhang } 412195687c5SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,garray,&c->vscale);CHKERRQ(ierr); 413195687c5SHong Zhang ierr = PetscFree(garray);CHKERRQ(ierr); 414195687c5SHong Zhang } 415e3225e9dSHong Zhang } else { 416e3225e9dSHong Zhang Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data; 417e3225e9dSHong Zhang Mat_SeqAIJ *spA,*spB; 418e3225e9dSHong Zhang A = aij->A; spA = (Mat_SeqAIJ*)A->data; A_val = spA->a; 419e3225e9dSHong Zhang B = aij->B; spB = (Mat_SeqAIJ*)B->data; B_val = spB->a; 420e3225e9dSHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 421e3225e9dSHong Zhang if (!aij->colmap) { 422e3225e9dSHong Zhang /* Allow access to data structures of local part of matrix 423e3225e9dSHong Zhang - creates aij->colmap which maps global column number to local number in part B */ 424e3225e9dSHong Zhang ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 425e3225e9dSHong Zhang colmap = aij->colmap; 426e3225e9dSHong Zhang } 427e3225e9dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 428e3225e9dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 429e3225e9dSHong Zhang 430e3225e9dSHong Zhang bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */ 431195687c5SHong Zhang 432195687c5SHong Zhang if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */ 433195687c5SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 434195687c5SHong Zhang } 4350d1c53f1SHong Zhang } 4360d1c53f1SHong Zhang 4370d1c53f1SHong Zhang m = mat->rmap->n/bs; 4380d1c53f1SHong Zhang cstart = mat->cmap->rstart/bs; 4390d1c53f1SHong Zhang cend = mat->cmap->rend/bs; 440fcd7ac73SHong Zhang 441fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 442fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 443fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 4448bc97078SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 445fcd7ac73SHong Zhang 4460df34763SHong Zhang if (c->htype[0] == 'd') { 447573f477fSHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 4488bc97078SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 4498bc97078SHong Zhang c->matentry = Jentry; 4500df34763SHong Zhang } else if (c->htype[0] == 'w') { 4510df34763SHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry2),&Jentry2);CHKERRQ(ierr); 4520df34763SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry2));CHKERRQ(ierr); 4530df34763SHong Zhang c->matentry2 = Jentry2; 4540df34763SHong Zhang } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"htype is not supported"); 455a774a6f1SHong Zhang 456*dcca6d9dSJed Brown ierr = PetscMalloc2(m+1,&rowhit,m+1,&valaddrhit);CHKERRQ(ierr); 457fcd7ac73SHong Zhang nz = 0; 45872c15787SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 459d3825b63SHong Zhang for (i=0; i<nis; i++) { /* for each local color */ 460fcd7ac73SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 461fcd7ac73SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 462fcd7ac73SHong Zhang 463fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 464fcd7ac73SHong Zhang if (n) { 465fcd7ac73SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 466fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 467fcd7ac73SHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 468fcd7ac73SHong Zhang } else { 469fcd7ac73SHong Zhang c->columns[i] = 0; 470fcd7ac73SHong Zhang } 471fcd7ac73SHong Zhang 472fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 473d3825b63SHong Zhang /* Determine nctot, the total (parallel) number of columns of this color */ 474fcd7ac73SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 475*dcca6d9dSJed Brown ierr = PetscMalloc2(size,&ncolsonproc,size,&disp);CHKERRQ(ierr); 476fcd7ac73SHong Zhang 477d3825b63SHong Zhang /* ncolsonproc[j]: local ncolumns on proc[j] of this color */ 478fcd7ac73SHong Zhang ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 479fcd7ac73SHong Zhang ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 480fcd7ac73SHong Zhang nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 481fcd7ac73SHong Zhang if (!nctot) { 482fcd7ac73SHong Zhang ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 483fcd7ac73SHong Zhang } 484fcd7ac73SHong Zhang 485fcd7ac73SHong Zhang disp[0] = 0; 486fcd7ac73SHong Zhang for (j=1; j<size; j++) { 487fcd7ac73SHong Zhang disp[j] = disp[j-1] + ncolsonproc[j-1]; 488fcd7ac73SHong Zhang } 489fcd7ac73SHong Zhang 490d3825b63SHong Zhang /* Get cols, the complete list of columns for this color on each process */ 491fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 492fcd7ac73SHong Zhang ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 493fcd7ac73SHong Zhang ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 494fcd7ac73SHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 495fcd7ac73SHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 496fcd7ac73SHong Zhang nctot = n; 497fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 498fcd7ac73SHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 499fcd7ac73SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 500fcd7ac73SHong Zhang 5011b97d346SHong Zhang /* Mark all rows affect by these columns */ 502d3825b63SHong Zhang ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr); 5030d1c53f1SHong Zhang bs2 = bs*bs; 5044b2e90caSHong Zhang nrows_i = 0; 5051b97d346SHong Zhang for (j=0; j<nctot; j++) { /* loop over columns*/ 506fcd7ac73SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 507fcd7ac73SHong Zhang col = ltog[cols[j]]; 508fcd7ac73SHong Zhang } else { 509fcd7ac73SHong Zhang col = cols[j]; 510fcd7ac73SHong Zhang } 511e3225e9dSHong Zhang if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */ 512d3825b63SHong Zhang row = A_cj + A_ci[col-cstart]; 513d3825b63SHong Zhang nrows = A_ci[col-cstart+1] - A_ci[col-cstart]; 5144b2e90caSHong Zhang nrows_i += nrows; 515fcd7ac73SHong Zhang /* loop over columns of A marking them in rowhit */ 516d3825b63SHong Zhang for (k=0; k<nrows; k++) { 517d3825b63SHong Zhang /* set valaddrhit for part A */ 518e3225e9dSHong Zhang spidx = bs2*spidxA[A_ci[col-cstart] + k]; 519e3225e9dSHong Zhang valaddrhit[*row] = &A_val[spidx]; 520a774a6f1SHong Zhang rowhit[*row++] = col - cstart + 1; /* local column index */ 521fcd7ac73SHong Zhang } 522e3225e9dSHong Zhang } else { /* column is in B, off-diagonal block of mat */ 523fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE) 5240d1c53f1SHong Zhang ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr); 525fcd7ac73SHong Zhang colb--; 526fcd7ac73SHong Zhang #else 5270d1c53f1SHong Zhang colb = colmap[col] - 1; /* local column index */ 528fcd7ac73SHong Zhang #endif 529fcd7ac73SHong Zhang if (colb == -1) { 530d3825b63SHong Zhang nrows = 0; 531fcd7ac73SHong Zhang } else { 5320d1c53f1SHong Zhang colb = colb/bs; 533d3825b63SHong Zhang row = B_cj + B_ci[colb]; 534d3825b63SHong Zhang nrows = B_ci[colb+1] - B_ci[colb]; 535fcd7ac73SHong Zhang } 5364b2e90caSHong Zhang nrows_i += nrows; 537fcd7ac73SHong Zhang /* loop over columns of B marking them in rowhit */ 538d3825b63SHong Zhang for (k=0; k<nrows; k++) { 539d3825b63SHong Zhang /* set valaddrhit for part B */ 540e3225e9dSHong Zhang spidx = bs2*spidxB[B_ci[colb] + k]; 541e3225e9dSHong Zhang valaddrhit[*row] = &B_val[spidx]; 54270e7395fSHong Zhang rowhit[*row++] = colb + 1 + cend - cstart; /* local column index */ 543fcd7ac73SHong Zhang } 544fcd7ac73SHong Zhang } 545fcd7ac73SHong Zhang } 5464b2e90caSHong Zhang c->nrows[i] = nrows_i; 5478bc97078SHong Zhang 5480df34763SHong Zhang if (c->htype[0] == 'd') { 549d3825b63SHong Zhang for (j=0; j<m; j++) { 550fcd7ac73SHong Zhang if (rowhit[j]) { 551573f477fSHong Zhang Jentry[nz].row = j; /* local row index */ 552573f477fSHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 553573f477fSHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 554573f477fSHong Zhang nz++; 555fcd7ac73SHong Zhang } 556fcd7ac73SHong Zhang } 5570df34763SHong Zhang } else { /* c->htype == 'wp' */ 5580df34763SHong Zhang for (j=0; j<m; j++) { 5590df34763SHong Zhang if (rowhit[j]) { 5600df34763SHong Zhang Jentry2[nz].row = j; /* local row index */ 5610df34763SHong Zhang Jentry2[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 5620df34763SHong Zhang nz++; 5630df34763SHong Zhang } 5640df34763SHong Zhang } 5650df34763SHong Zhang } 566fcd7ac73SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 567fcd7ac73SHong Zhang } 568fcd7ac73SHong Zhang 569a8971b87SHong Zhang if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */ 570a8971b87SHong Zhang ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr); 571531e53bdSHong Zhang } 572531e53bdSHong Zhang 5730d1c53f1SHong Zhang if (isBAIJ) { 5740d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 5750d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 576e3225e9dSHong Zhang ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 5770d1c53f1SHong Zhang } else { 5780d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 5790d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 5800d1c53f1SHong Zhang } 5810d1c53f1SHong Zhang 582a8971b87SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 583a8971b87SHong Zhang ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 584a8971b87SHong Zhang 58572c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 58672c15787SHong Zhang ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr); 58772c15787SHong Zhang } 588b3e1f37bSHong Zhang ierr = PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr); 589fcd7ac73SHong Zhang PetscFunctionReturn(0); 590fcd7ac73SHong Zhang } 59193dfae19SHong Zhang 59293dfae19SHong Zhang #undef __FUNCT__ 59393dfae19SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIXAIJ" 59493dfae19SHong Zhang PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 59593dfae19SHong Zhang { 596531e53bdSHong Zhang PetscErrorCode ierr; 597a8971b87SHong Zhang PetscInt bs,nis=iscoloring->n,m=mat->rmap->n; 598531e53bdSHong Zhang PetscBool isBAIJ; 599531e53bdSHong Zhang 60093dfae19SHong Zhang PetscFunctionBegin; 601531e53bdSHong Zhang /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian; 602531e53bdSHong Zhang bcols is chosen s.t. dy-array takes 50% of memory space as mat */ 603531e53bdSHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 604531e53bdSHong Zhang ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr); 605a8971b87SHong Zhang if (isBAIJ) { 606a8971b87SHong Zhang c->brows = m; 607a8971b87SHong Zhang c->bcols = 1; 608531e53bdSHong Zhang } else { /* mpiaij matrix */ 609531e53bdSHong Zhang /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */ 610531e53bdSHong Zhang Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data; 611531e53bdSHong Zhang Mat_SeqAIJ *spA,*spB; 612531e53bdSHong Zhang Mat A,B; 613a8971b87SHong Zhang PetscInt nz,brows,bcols; 614531e53bdSHong Zhang PetscReal mem; 615531e53bdSHong Zhang 616531e53bdSHong Zhang bs = 1; /* only bs=1 is supported for MPIAIJ matrix */ 617531e53bdSHong Zhang 618531e53bdSHong Zhang A = aij->A; spA = (Mat_SeqAIJ*)A->data; 619531e53bdSHong Zhang B = aij->B; spB = (Mat_SeqAIJ*)B->data; 620531e53bdSHong Zhang nz = spA->nz + spB->nz; /* total local nonzero entries of mat */ 621531e53bdSHong Zhang mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt); 622531e53bdSHong Zhang bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar))); 623531e53bdSHong Zhang brows = 1000/bcols; 624531e53bdSHong Zhang if (bcols > nis) bcols = nis; 625531e53bdSHong Zhang if (brows == 0 || brows > m) brows = m; 626531e53bdSHong Zhang c->brows = brows; 627531e53bdSHong Zhang c->bcols = bcols; 628531e53bdSHong Zhang } 629a8971b87SHong Zhang 630a8971b87SHong Zhang c->M = mat->rmap->N/bs; /* set the global rows and columns and local rows */ 631a8971b87SHong Zhang c->N = mat->cmap->N/bs; 632a8971b87SHong Zhang c->m = mat->rmap->n/bs; 633a8971b87SHong Zhang c->rstart = mat->rmap->rstart/bs; 634a8971b87SHong Zhang c->ncolors = nis; 63593dfae19SHong Zhang PetscFunctionReturn(0); 63693dfae19SHong Zhang } 637