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; 120d1c53f1SHong Zhang PetscScalar dx=0.0,*y,*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; 200d1c53f1SHong Zhang const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; 210d1c53f1SHong Zhang PetscInt bs=J->rmap->bs; 220d1c53f1SHong Zhang 230d1c53f1SHong Zhang PetscFunctionBegin; 240d1c53f1SHong Zhang /* create vscale for storing dx */ 250d1c53f1SHong Zhang if (!vscale) { 260d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 276c145f34SHong Zhang /* contain the "diagonal" on processor scalings followed by the off processor* - garray must be non-bloced! */ 286c145f34SHong Zhang Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)J->data; 296c145f34SHong Zhang PetscInt *garray; 306c145f34SHong Zhang ierr = PetscMalloc(baij->B->cmap->n*sizeof(PetscInt),&garray);CHKERRQ(ierr); 316c145f34SHong Zhang for (i=0; i<baij->B->cmap->n/bs; i++) { 326c145f34SHong Zhang for (j=0; j<bs; j++) { 336c145f34SHong Zhang garray[i*bs+j] = bs*baij->garray[i]+j; 346c145f34SHong Zhang } 356c145f34SHong Zhang } 366c145f34SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,baij->B->cmap->n,garray,&vscale);CHKERRQ(ierr); 376c145f34SHong Zhang ierr = PetscFree(garray);CHKERRQ(ierr); 380d1c53f1SHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 390d1c53f1SHong Zhang ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr); 400d1c53f1SHong Zhang } 410d1c53f1SHong Zhang coloring->vscale = vscale; 420d1c53f1SHong Zhang } 430d1c53f1SHong Zhang 440d1c53f1SHong Zhang /* (1) Set w1 = F(x1) */ 450d1c53f1SHong Zhang if (!coloring->fset) { 460d1c53f1SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 470d1c53f1SHong Zhang ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 480d1c53f1SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 490d1c53f1SHong Zhang } else { 500d1c53f1SHong Zhang coloring->fset = PETSC_FALSE; 510d1c53f1SHong Zhang } 520d1c53f1SHong Zhang 530d1c53f1SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 540d1c53f1SHong Zhang ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 550d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 560d1c53f1SHong Zhang /* vscale = dx is a constant scalar */ 570d1c53f1SHong Zhang ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 580d1c53f1SHong Zhang dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); 590d1c53f1SHong Zhang } else { 600d1c53f1SHong Zhang ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 610d1c53f1SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 620d1c53f1SHong Zhang for (col=0; col<nxloc; col++) { 630d1c53f1SHong Zhang dx = xx[col]; 640d1c53f1SHong Zhang if (PetscAbsScalar(dx) < umin) { 650d1c53f1SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 660d1c53f1SHong Zhang else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 670d1c53f1SHong Zhang } 680d1c53f1SHong Zhang dx *= epsilon; 690d1c53f1SHong Zhang vscale_array[col] = 1.0/dx; 700d1c53f1SHong Zhang } 710d1c53f1SHong Zhang ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 720d1c53f1SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 730d1c53f1SHong Zhang } 74*684f2004SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 750d1c53f1SHong Zhang ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 760d1c53f1SHong Zhang ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 770d1c53f1SHong Zhang } 780d1c53f1SHong Zhang 790d1c53f1SHong Zhang /* (3) Loop over each color */ 800d1c53f1SHong Zhang if (!coloring->w3) { 810d1c53f1SHong Zhang ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 820d1c53f1SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 830d1c53f1SHong Zhang } 840d1c53f1SHong Zhang w3 = coloring->w3; 850d1c53f1SHong Zhang 860d1c53f1SHong Zhang ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 870d1c53f1SHong Zhang if (vscale) { 880d1c53f1SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 890d1c53f1SHong Zhang } 900d1c53f1SHong Zhang nz = 0; 910d1c53f1SHong Zhang for (k=0; k<ncolors; k++) { 920d1c53f1SHong Zhang coloring->currentcolor = k; 930d1c53f1SHong Zhang 940d1c53f1SHong Zhang /* 950d1c53f1SHong Zhang (3-1) Loop over each column associated with color 960d1c53f1SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 970d1c53f1SHong Zhang */ 980d1c53f1SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 990d1c53f1SHong Zhang dy_i = dy; 1000d1c53f1SHong Zhang for (i=0; i<bs; i++) { /* Loop over a block of columns */ 1010d1c53f1SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 1020d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 1030d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 1040d1c53f1SHong Zhang for (l=0; l<ncolumns[k]; l++) { 1050d1c53f1SHong Zhang col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 1060d1c53f1SHong Zhang w3_array[col] += 1.0/dx; 107*684f2004SHong Zhang if (i) w3_array[col-1] -= 1.0/dx; /* resume original w3[col-1] */ 1080d1c53f1SHong Zhang } 1090d1c53f1SHong Zhang } else { /* htype == 'ds' */ 1100d1c53f1SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 1110d1c53f1SHong Zhang for (l=0; l<ncolumns[k]; l++) { 112f8c2866eSHong Zhang col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 1130d1c53f1SHong Zhang w3_array[col] += 1.0/vscale_array[col]; 114*684f2004SHong Zhang if (i) w3_array[col-1] -= 1.0/vscale_array[col-1]; /* resume original w3[col-1] */ 1150d1c53f1SHong Zhang } 1160d1c53f1SHong Zhang vscale_array += cstart; 1170d1c53f1SHong Zhang } 1180d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 1190d1c53f1SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 1200d1c53f1SHong Zhang 1210d1c53f1SHong Zhang /* 1220d1c53f1SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 1230d1c53f1SHong Zhang w2 = F(x1 + dx) - F(x1) 1240d1c53f1SHong Zhang */ 1250d1c53f1SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 1260d1c53f1SHong Zhang ierr = VecPlaceArray(w2,dy_i);CHKERRQ(ierr); /* place w2 to the array dy_i */ 1270d1c53f1SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 1280d1c53f1SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 1290d1c53f1SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 1300d1c53f1SHong Zhang ierr = VecResetArray(w2);CHKERRQ(ierr); 1310d1c53f1SHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1320d1c53f1SHong Zhang } 1330d1c53f1SHong Zhang 1340d1c53f1SHong Zhang /* 1350d1c53f1SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 1360d1c53f1SHong Zhang */ 1370d1c53f1SHong Zhang nrows_k = nrows[k]; 1380d1c53f1SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 1390d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 1400d1c53f1SHong Zhang for (l=0; l<nrows_k; l++) { 1410d1c53f1SHong Zhang row = bs*Jentry[nz].row; /* local row index */ 142*684f2004SHong Zhang valaddr = Jentry[nz++].valaddr; 1430d1c53f1SHong Zhang spidx = 0; 1440d1c53f1SHong Zhang dy_i = dy; 1450d1c53f1SHong Zhang for (i=0; i<bs; i++) { /* column of the block */ 1460d1c53f1SHong Zhang for (j=0; j<bs; j++) { /* row of the block */ 1470d1c53f1SHong Zhang valaddr[spidx++] = dy_i[row+j]*dx; 1480d1c53f1SHong Zhang } 149f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1500d1c53f1SHong Zhang } 1510d1c53f1SHong Zhang } 1520d1c53f1SHong Zhang } else { /* htype == 'ds' */ 1530d1c53f1SHong Zhang for (l=0; l<nrows_k; l++) { 1540d1c53f1SHong Zhang row = bs*Jentry[nz].row; /* local row index */ 1550d1c53f1SHong Zhang col = bs*Jentry[nz].col; /* local column index */ 156*684f2004SHong Zhang valaddr = Jentry[nz++].valaddr; 157f8c2866eSHong Zhang spidx = 0; 158f8c2866eSHong Zhang dy_i = dy; 159f8c2866eSHong Zhang for (i=0; i<bs; i++) { /* column of the block */ 160f8c2866eSHong Zhang for (j=0; j<bs; j++) { /* row of the block */ 161f8c2866eSHong Zhang valaddr[spidx++] = dy_i[row+j]*vscale_array[col+i]; 162f8c2866eSHong Zhang } 163f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 164f8c2866eSHong Zhang } 1650d1c53f1SHong Zhang } 1660d1c53f1SHong Zhang } 1670d1c53f1SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 1680d1c53f1SHong Zhang } 1690d1c53f1SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1700d1c53f1SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1710d1c53f1SHong Zhang if (vscale) { 1720d1c53f1SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 1730d1c53f1SHong Zhang } 1740d1c53f1SHong Zhang 1750d1c53f1SHong Zhang coloring->currentcolor = -1; 1760d1c53f1SHong Zhang PetscFunctionReturn(0); 1770d1c53f1SHong Zhang } 178a64fbb32SBarry Smith 1794a2ae208SSatish Balay #undef __FUNCT__ 180040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringApply_AIJ" 181040ebd07SHong Zhang PetscErrorCode MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 182fcd7ac73SHong Zhang { 183fcd7ac73SHong Zhang PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; 184fcd7ac73SHong Zhang PetscErrorCode ierr; 185a2f2d239SHong Zhang PetscInt k,cstart,cend,l,row,col,nz; 1869e917edbSHong Zhang PetscScalar dx=0.0,*y,*xx,*w3_array; 187fcd7ac73SHong Zhang PetscScalar *vscale_array; 188fcd7ac73SHong Zhang PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; 1898bc97078SHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; 190fcd7ac73SHong Zhang void *fctx=coloring->fctx; 191c53567a0SHong Zhang PetscInt ctype=coloring->ctype,nxloc,nrows_k; 192573f477fSHong Zhang MatEntry *Jentry=coloring->matentry; 1938bc97078SHong Zhang const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; 194fcd7ac73SHong Zhang 195fcd7ac73SHong Zhang PetscFunctionBegin; 1969e917edbSHong Zhang /* create vscale for storing dx */ 1979e917edbSHong Zhang if (!vscale) { 1989e917edbSHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 1996c145f34SHong Zhang Mat_MPIAIJ *aij=(Mat_MPIAIJ*)J->data; 2009e917edbSHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr); 2019e917edbSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 2028bc97078SHong Zhang ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr); 2039e917edbSHong Zhang } 2048bc97078SHong Zhang coloring->vscale = vscale; 205fcd7ac73SHong Zhang } 206fcd7ac73SHong Zhang 2078bc97078SHong Zhang /* (1) Set w1 = F(x1) */ 208fcd7ac73SHong Zhang if (!coloring->fset) { 209fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 210f6af9589SHong Zhang ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 211fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 212fcd7ac73SHong Zhang } else { 213fcd7ac73SHong Zhang coloring->fset = PETSC_FALSE; 214fcd7ac73SHong Zhang } 215fcd7ac73SHong Zhang 2168bc97078SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 217f6af9589SHong Zhang if (coloring->htype[0] == 'w') { 218*684f2004SHong Zhang /* vscale = 1./dx is a constant scalar */ 219f6af9589SHong Zhang ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 220c53567a0SHong Zhang dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); 22170e7395fSHong Zhang } else { 22274d3cef9SHong Zhang ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 223f6af9589SHong Zhang ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 2248bc97078SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 22574d3cef9SHong Zhang for (col=0; col<nxloc; col++) { 226fcd7ac73SHong Zhang dx = xx[col]; 22774d3cef9SHong Zhang if (PetscAbsScalar(dx) < umin) { 22874d3cef9SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 22974d3cef9SHong Zhang else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 23074d3cef9SHong Zhang } 231fcd7ac73SHong Zhang dx *= epsilon; 23274d3cef9SHong Zhang vscale_array[col] = 1.0/dx; 233f6af9589SHong Zhang } 23474d3cef9SHong Zhang ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 2358bc97078SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 23670e7395fSHong Zhang } 237*684f2004SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 2388bc97078SHong Zhang ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2398bc97078SHong Zhang ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 240fcd7ac73SHong Zhang } 241fcd7ac73SHong Zhang 2428bc97078SHong Zhang /* (3) Loop over each color */ 2438bc97078SHong Zhang if (!coloring->w3) { 2448bc97078SHong Zhang ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 2458bc97078SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 2468bc97078SHong Zhang } 2478bc97078SHong Zhang w3 = coloring->w3; 248fcd7ac73SHong Zhang 2498bc97078SHong Zhang ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 2509e917edbSHong Zhang if (vscale) { 2518bc97078SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 2529e917edbSHong Zhang } 2538bc97078SHong Zhang nz = 0; 2548bc97078SHong Zhang for (k=0; k<ncolors; k++) { 2558bc97078SHong Zhang coloring->currentcolor = k; 2569e917edbSHong Zhang 257fcd7ac73SHong Zhang /* 2588bc97078SHong Zhang (3-1) Loop over each column associated with color 259f6af9589SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 260fcd7ac73SHong Zhang */ 261f6af9589SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 262f6af9589SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 263a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 264b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 2658bc97078SHong Zhang for (l=0; l<ncolumns[k]; l++) { 266f6af9589SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 2679e917edbSHong Zhang w3_array[col] += 1.0/dx; 2689e917edbSHong Zhang } 269b039d6c4SHong Zhang } else { /* htype == 'ds' */ 270a2f2d239SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 271b039d6c4SHong Zhang for (l=0; l<ncolumns[k]; l++) { 272b039d6c4SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 273a2f2d239SHong Zhang w3_array[col] += 1.0/vscale_array[col]; 27470e7395fSHong Zhang } 275a2f2d239SHong Zhang vscale_array += cstart; 276fcd7ac73SHong Zhang } 277a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 278fcd7ac73SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 2799e917edbSHong Zhang 280fcd7ac73SHong Zhang /* 2818bc97078SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 282fcd7ac73SHong Zhang w2 = F(x1 + dx) - F(x1) 283fcd7ac73SHong Zhang */ 284fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 285fcd7ac73SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 286fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 287fcd7ac73SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 2889e917edbSHong Zhang 2898bc97078SHong Zhang /* 2908bc97078SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 2918bc97078SHong Zhang */ 292c53567a0SHong Zhang nrows_k = nrows[k]; 293fcd7ac73SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 294b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 295c53567a0SHong Zhang for (l=0; l<nrows_k; l++) { 296573f477fSHong Zhang row = Jentry[nz].row; /* local row index */ 297b039d6c4SHong Zhang *(Jentry[nz++].valaddr) = y[row]*dx; 2989e917edbSHong Zhang } 299b039d6c4SHong Zhang } else { /* htype == 'ds' */ 300c53567a0SHong Zhang for (l=0; l<nrows_k; l++) { 301b039d6c4SHong Zhang row = Jentry[nz].row; /* local row index */ 302b039d6c4SHong Zhang *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col]; 303573f477fSHong Zhang nz++; 304fcd7ac73SHong Zhang } 305b039d6c4SHong Zhang } 306fcd7ac73SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 3078bc97078SHong Zhang } 308f5aae955SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 309f5aae955SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3109e917edbSHong Zhang if (vscale) { 3118bc97078SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 3129e917edbSHong Zhang } 313fcd7ac73SHong Zhang 314fcd7ac73SHong Zhang coloring->currentcolor = -1; 315fcd7ac73SHong Zhang PetscFunctionReturn(0); 316fcd7ac73SHong Zhang } 317fcd7ac73SHong Zhang 318fcd7ac73SHong Zhang #undef __FUNCT__ 319040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIXAIJ" 320040ebd07SHong Zhang PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 321fcd7ac73SHong Zhang { 322fcd7ac73SHong Zhang PetscErrorCode ierr; 323fcd7ac73SHong Zhang PetscMPIInt size,*ncolsonproc,*disp,nn; 324e3225e9dSHong Zhang PetscInt i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb; 32572c15787SHong Zhang const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL; 326fcd7ac73SHong Zhang PetscInt nis=iscoloring->n,nctot,*cols; 327fcd7ac73SHong Zhang IS *isa; 328fcd7ac73SHong Zhang ISLocalToGlobalMapping map=mat->cmap->mapping; 329e3225e9dSHong Zhang PetscInt ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx; 3300d1c53f1SHong Zhang Mat A,B; 331e3225e9dSHong Zhang PetscScalar *A_val,*B_val,**valaddrhit; 332573f477fSHong Zhang MatEntry *Jentry; 3330d1c53f1SHong Zhang PetscBool isBAIJ; 3340d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE) 3350d1c53f1SHong Zhang PetscTable colmap=NULL; 3360d1c53f1SHong Zhang #else 3370d1c53f1SHong Zhang PetscInt *colmap=NULL; /* local col number of off-diag col */ 3380d1c53f1SHong Zhang #endif 339fcd7ac73SHong Zhang 340fcd7ac73SHong Zhang PetscFunctionBegin; 34172c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 34272c15787SHong 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"); 34372c15787SHong Zhang ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr); 34472c15787SHong Zhang } 345fcd7ac73SHong Zhang 3460d1c53f1SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 3470d1c53f1SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr); 348e3225e9dSHong Zhang if (isBAIJ) { 3490d1c53f1SHong Zhang Mat_MPIBAIJ *aij=(Mat_MPIBAIJ*)mat->data; 3506c145f34SHong Zhang Mat_SeqBAIJ *spA,*spB; 3516c145f34SHong Zhang A = aij->A; spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a; 3526c145f34SHong Zhang B = aij->B; spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a; 3530d1c53f1SHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 3540d1c53f1SHong Zhang if (!aij->colmap) { 3550d1c53f1SHong Zhang ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr); 3560d1c53f1SHong Zhang colmap = aij->colmap; 3570d1c53f1SHong Zhang } 3580d1c53f1SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 3590d1c53f1SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 360e3225e9dSHong Zhang } else { 361e3225e9dSHong Zhang Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data; 362e3225e9dSHong Zhang Mat_SeqAIJ *spA,*spB; 363e3225e9dSHong Zhang A = aij->A; spA = (Mat_SeqAIJ*)A->data; A_val = spA->a; 364e3225e9dSHong Zhang B = aij->B; spB = (Mat_SeqAIJ*)B->data; B_val = spB->a; 365e3225e9dSHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 366e3225e9dSHong Zhang if (!aij->colmap) { 367e3225e9dSHong Zhang /* Allow access to data structures of local part of matrix 368e3225e9dSHong Zhang - creates aij->colmap which maps global column number to local number in part B */ 369e3225e9dSHong Zhang ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 370e3225e9dSHong Zhang colmap = aij->colmap; 371e3225e9dSHong Zhang } 372e3225e9dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 373e3225e9dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 374e3225e9dSHong Zhang 375e3225e9dSHong Zhang bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */ 3760d1c53f1SHong Zhang } 3770d1c53f1SHong Zhang 3780d1c53f1SHong Zhang m = mat->rmap->n/bs; 3790d1c53f1SHong Zhang cstart = mat->cmap->rstart/bs; 3800d1c53f1SHong Zhang cend = mat->cmap->rend/bs; 3810d1c53f1SHong Zhang c->M = mat->rmap->N/bs; /* set the global rows and columns and local rows */ 3820d1c53f1SHong Zhang c->N = mat->cmap->N/bs; 383d3825b63SHong Zhang c->m = m; 3840d1c53f1SHong Zhang c->rstart = mat->rmap->rstart/bs; 385fcd7ac73SHong Zhang 386fcd7ac73SHong Zhang c->ncolors = nis; 387fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 388fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 389fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 3908bc97078SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 391fcd7ac73SHong Zhang 392573f477fSHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 3938bc97078SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 3948bc97078SHong Zhang c->matentry = Jentry; 395a774a6f1SHong Zhang 396d3825b63SHong Zhang ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 397fcd7ac73SHong Zhang nz = 0; 39872c15787SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 399d3825b63SHong Zhang for (i=0; i<nis; i++) { /* for each local color */ 400fcd7ac73SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 401fcd7ac73SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 402fcd7ac73SHong Zhang 403fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 404fcd7ac73SHong Zhang if (n) { 405fcd7ac73SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 406fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 407fcd7ac73SHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 408fcd7ac73SHong Zhang } else { 409fcd7ac73SHong Zhang c->columns[i] = 0; 410fcd7ac73SHong Zhang } 411fcd7ac73SHong Zhang 412fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 413d3825b63SHong Zhang /* Determine nctot, the total (parallel) number of columns of this color */ 414fcd7ac73SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 415fcd7ac73SHong Zhang ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 416fcd7ac73SHong Zhang 417d3825b63SHong Zhang /* ncolsonproc[j]: local ncolumns on proc[j] of this color */ 418fcd7ac73SHong Zhang ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 419fcd7ac73SHong Zhang ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 420fcd7ac73SHong Zhang nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 421fcd7ac73SHong Zhang if (!nctot) { 422fcd7ac73SHong Zhang ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 423fcd7ac73SHong Zhang } 424fcd7ac73SHong Zhang 425fcd7ac73SHong Zhang disp[0] = 0; 426fcd7ac73SHong Zhang for (j=1; j<size; j++) { 427fcd7ac73SHong Zhang disp[j] = disp[j-1] + ncolsonproc[j-1]; 428fcd7ac73SHong Zhang } 429fcd7ac73SHong Zhang 430d3825b63SHong Zhang /* Get cols, the complete list of columns for this color on each process */ 431fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 432fcd7ac73SHong Zhang ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 433fcd7ac73SHong Zhang ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 434fcd7ac73SHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 435fcd7ac73SHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 436fcd7ac73SHong Zhang nctot = n; 437fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 438fcd7ac73SHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 439fcd7ac73SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 440fcd7ac73SHong Zhang 4411b97d346SHong Zhang /* Mark all rows affect by these columns */ 442d3825b63SHong Zhang ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr); 4430d1c53f1SHong Zhang bs2 = bs*bs; 4444b2e90caSHong Zhang nrows_i = 0; 4451b97d346SHong Zhang for (j=0; j<nctot; j++) { /* loop over columns*/ 446fcd7ac73SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 447fcd7ac73SHong Zhang col = ltog[cols[j]]; 448fcd7ac73SHong Zhang } else { 449fcd7ac73SHong Zhang col = cols[j]; 450fcd7ac73SHong Zhang } 451e3225e9dSHong Zhang if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */ 452d3825b63SHong Zhang row = A_cj + A_ci[col-cstart]; 453d3825b63SHong Zhang nrows = A_ci[col-cstart+1] - A_ci[col-cstart]; 4544b2e90caSHong Zhang nrows_i += nrows; 455fcd7ac73SHong Zhang /* loop over columns of A marking them in rowhit */ 456d3825b63SHong Zhang for (k=0; k<nrows; k++) { 457d3825b63SHong Zhang /* set valaddrhit for part A */ 458e3225e9dSHong Zhang spidx = bs2*spidxA[A_ci[col-cstart] + k]; 459e3225e9dSHong Zhang valaddrhit[*row] = &A_val[spidx]; 460a774a6f1SHong Zhang rowhit[*row++] = col - cstart + 1; /* local column index */ 461fcd7ac73SHong Zhang } 462e3225e9dSHong Zhang } else { /* column is in B, off-diagonal block of mat */ 463fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE) 4640d1c53f1SHong Zhang ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr); 465fcd7ac73SHong Zhang colb--; 466fcd7ac73SHong Zhang #else 4670d1c53f1SHong Zhang colb = colmap[col] - 1; /* local column index */ 468fcd7ac73SHong Zhang #endif 469fcd7ac73SHong Zhang if (colb == -1) { 470d3825b63SHong Zhang nrows = 0; 471fcd7ac73SHong Zhang } else { 4720d1c53f1SHong Zhang colb = colb/bs; 473d3825b63SHong Zhang row = B_cj + B_ci[colb]; 474d3825b63SHong Zhang nrows = B_ci[colb+1] - B_ci[colb]; 475fcd7ac73SHong Zhang } 4764b2e90caSHong Zhang nrows_i += nrows; 477fcd7ac73SHong Zhang /* loop over columns of B marking them in rowhit */ 478d3825b63SHong Zhang for (k=0; k<nrows; k++) { 479d3825b63SHong Zhang /* set valaddrhit for part B */ 480e3225e9dSHong Zhang spidx = bs2*spidxB[B_ci[colb] + k]; 481e3225e9dSHong Zhang valaddrhit[*row] = &B_val[spidx]; 48270e7395fSHong Zhang rowhit[*row++] = colb + 1 + cend - cstart; /* local column index */ 483fcd7ac73SHong Zhang } 484fcd7ac73SHong Zhang } 485fcd7ac73SHong Zhang } 4864b2e90caSHong Zhang c->nrows[i] = nrows_i; 4878bc97078SHong Zhang 488d3825b63SHong Zhang for (j=0; j<m; j++) { 489fcd7ac73SHong Zhang if (rowhit[j]) { 490573f477fSHong Zhang Jentry[nz].row = j; /* local row index */ 491573f477fSHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 492573f477fSHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 493573f477fSHong Zhang nz++; 494fcd7ac73SHong Zhang } 495fcd7ac73SHong Zhang } 496fcd7ac73SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 497fcd7ac73SHong Zhang } 4988bc97078SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 499fcd7ac73SHong Zhang 500d3825b63SHong Zhang ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 5010d1c53f1SHong Zhang if (isBAIJ) { 5020d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 5030d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 504e3225e9dSHong Zhang ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 5050d1c53f1SHong Zhang } else { 5060d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 5070d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 5080d1c53f1SHong Zhang } 5090d1c53f1SHong Zhang 51072c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 51172c15787SHong Zhang ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr); 51272c15787SHong Zhang } 513fcd7ac73SHong Zhang PetscFunctionReturn(0); 514fcd7ac73SHong Zhang } 515