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 PetscBool flg=PETSC_FALSE; 180d1c53f1SHong Zhang PetscInt ctype=coloring->ctype,nxloc,nrows_k; 190d1c53f1SHong Zhang PetscScalar *valaddr; 200d1c53f1SHong Zhang MatEntry *Jentry=coloring->matentry; 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 ierr = MatSetUnfactored(J);CHKERRQ(ierr); 260d1c53f1SHong Zhang ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); 270d1c53f1SHong Zhang if (flg) { 280d1c53f1SHong Zhang ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 290d1c53f1SHong Zhang } else { 300d1c53f1SHong Zhang PetscBool assembled; 310d1c53f1SHong Zhang ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 320d1c53f1SHong Zhang if (assembled) { 330d1c53f1SHong Zhang ierr = MatZeroEntries(J);CHKERRQ(ierr); 340d1c53f1SHong Zhang } 350d1c53f1SHong Zhang } 360d1c53f1SHong Zhang 370d1c53f1SHong Zhang /* create vscale for storing dx */ 380d1c53f1SHong Zhang if (!vscale) { 390d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 406c145f34SHong Zhang /* contain the "diagonal" on processor scalings followed by the off processor* - garray must be non-bloced! */ 416c145f34SHong Zhang Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)J->data; 426c145f34SHong Zhang PetscInt *garray; 436c145f34SHong Zhang ierr = PetscMalloc(baij->B->cmap->n*sizeof(PetscInt),&garray);CHKERRQ(ierr); 446c145f34SHong Zhang for (i=0; i<baij->B->cmap->n/bs; i++) { 456c145f34SHong Zhang for (j=0; j<bs; j++) { 466c145f34SHong Zhang garray[i*bs+j] = bs*baij->garray[i]+j; 476c145f34SHong Zhang } 486c145f34SHong Zhang } 496c145f34SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,baij->B->cmap->n,garray,&vscale);CHKERRQ(ierr); 506c145f34SHong Zhang ierr = PetscFree(garray);CHKERRQ(ierr); 510d1c53f1SHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 520d1c53f1SHong Zhang ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr); 530d1c53f1SHong Zhang } 540d1c53f1SHong Zhang coloring->vscale = vscale; 550d1c53f1SHong Zhang } 560d1c53f1SHong Zhang 570d1c53f1SHong Zhang /* (1) Set w1 = F(x1) */ 580d1c53f1SHong Zhang if (!coloring->fset) { 590d1c53f1SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 600d1c53f1SHong Zhang ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 610d1c53f1SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 620d1c53f1SHong Zhang } else { 630d1c53f1SHong Zhang coloring->fset = PETSC_FALSE; 640d1c53f1SHong Zhang } 650d1c53f1SHong Zhang 660d1c53f1SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 670d1c53f1SHong Zhang ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 680d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 690d1c53f1SHong Zhang /* vscale = dx is a constant scalar */ 700d1c53f1SHong Zhang ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 710d1c53f1SHong Zhang dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); 720d1c53f1SHong Zhang } else { 730d1c53f1SHong Zhang ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 740d1c53f1SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 750d1c53f1SHong Zhang for (col=0; col<nxloc; col++) { 760d1c53f1SHong Zhang dx = xx[col]; 770d1c53f1SHong Zhang if (PetscAbsScalar(dx) < umin) { 780d1c53f1SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 790d1c53f1SHong Zhang else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 800d1c53f1SHong Zhang } 810d1c53f1SHong Zhang dx *= epsilon; 820d1c53f1SHong Zhang vscale_array[col] = 1.0/dx; 830d1c53f1SHong Zhang } 840d1c53f1SHong Zhang ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 850d1c53f1SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 860d1c53f1SHong Zhang } 870d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] != 'w') { 880d1c53f1SHong Zhang ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 890d1c53f1SHong Zhang ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 900d1c53f1SHong Zhang } 910d1c53f1SHong Zhang 920d1c53f1SHong Zhang /* (3) Loop over each color */ 930d1c53f1SHong Zhang if (!coloring->w3) { 940d1c53f1SHong Zhang ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 950d1c53f1SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 960d1c53f1SHong Zhang } 970d1c53f1SHong Zhang w3 = coloring->w3; 980d1c53f1SHong Zhang 990d1c53f1SHong Zhang ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 1000d1c53f1SHong Zhang if (vscale) { 1010d1c53f1SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 1020d1c53f1SHong Zhang } 1030d1c53f1SHong Zhang nz = 0; 1040d1c53f1SHong Zhang for (k=0; k<ncolors; k++) { 1050d1c53f1SHong Zhang coloring->currentcolor = k; 1060d1c53f1SHong Zhang 1070d1c53f1SHong Zhang /* 1080d1c53f1SHong Zhang (3-1) Loop over each column associated with color 1090d1c53f1SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 1100d1c53f1SHong Zhang */ 1110d1c53f1SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 1120d1c53f1SHong Zhang dy_i = dy; 1130d1c53f1SHong Zhang for (i=0; i<bs; i++) { /* Loop over a block of columns */ 1140d1c53f1SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 1150d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 1160d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 1170d1c53f1SHong Zhang for (l=0; l<ncolumns[k]; l++) { 1180d1c53f1SHong Zhang col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 1190d1c53f1SHong Zhang w3_array[col] += 1.0/dx; 1200d1c53f1SHong Zhang 1210d1c53f1SHong Zhang if (i) { 1220d1c53f1SHong Zhang w3_array[col-1] -= 1.0/dx; /* resume original w3[col-1] */ 1230d1c53f1SHong Zhang } 1240d1c53f1SHong Zhang 1250d1c53f1SHong Zhang } 1260d1c53f1SHong Zhang } else { /* htype == 'ds' */ 1270d1c53f1SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 1280d1c53f1SHong Zhang for (l=0; l<ncolumns[k]; l++) { 129f8c2866eSHong Zhang col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 1300d1c53f1SHong Zhang w3_array[col] += 1.0/vscale_array[col]; 1310d1c53f1SHong Zhang 1320d1c53f1SHong Zhang if (i) { 1330d1c53f1SHong Zhang w3_array[col-1] -= 1.0/vscale_array[col-1]; /* resume original w3[col-1] */ 1340d1c53f1SHong Zhang } 1350d1c53f1SHong Zhang 1360d1c53f1SHong Zhang } 1370d1c53f1SHong Zhang vscale_array += cstart; 1380d1c53f1SHong Zhang } 1390d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 1400d1c53f1SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 1410d1c53f1SHong Zhang 1420d1c53f1SHong Zhang /* 1430d1c53f1SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 1440d1c53f1SHong Zhang w2 = F(x1 + dx) - F(x1) 1450d1c53f1SHong Zhang */ 1460d1c53f1SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 1470d1c53f1SHong Zhang ierr = VecPlaceArray(w2,dy_i);CHKERRQ(ierr); /* place w2 to the array dy_i */ 1480d1c53f1SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 1490d1c53f1SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 1500d1c53f1SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 1510d1c53f1SHong Zhang ierr = VecResetArray(w2);CHKERRQ(ierr); 1520d1c53f1SHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1530d1c53f1SHong Zhang } 1540d1c53f1SHong Zhang 1550d1c53f1SHong Zhang /* 1560d1c53f1SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 1570d1c53f1SHong Zhang */ 1580d1c53f1SHong Zhang nrows_k = nrows[k]; 1590d1c53f1SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 1600d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 1610d1c53f1SHong Zhang for (l=0; l<nrows_k; l++) { 1620d1c53f1SHong Zhang row = bs*Jentry[nz].row; /* local row index */ 1630d1c53f1SHong Zhang valaddr = Jentry[nz].valaddr; 1640d1c53f1SHong Zhang nz++; 1650d1c53f1SHong Zhang spidx = 0; 1660d1c53f1SHong Zhang dy_i = dy; 1670d1c53f1SHong Zhang for (i=0; i<bs; i++) { /* column of the block */ 1680d1c53f1SHong Zhang for (j=0; j<bs; j++) { /* row of the block */ 1690d1c53f1SHong Zhang valaddr[spidx++] = dy_i[row+j]*dx; 1700d1c53f1SHong Zhang } 171f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1720d1c53f1SHong Zhang } 1730d1c53f1SHong Zhang } 1740d1c53f1SHong Zhang } else { /* htype == 'ds' */ 1750d1c53f1SHong Zhang for (l=0; l<nrows_k; l++) { 1760d1c53f1SHong Zhang row = bs*Jentry[nz].row; /* local row index */ 1770d1c53f1SHong Zhang col = bs*Jentry[nz].col; /* local column index */ 178f8c2866eSHong Zhang valaddr = Jentry[nz].valaddr; 1790d1c53f1SHong Zhang nz++; 180f8c2866eSHong Zhang spidx = 0; 181f8c2866eSHong Zhang dy_i = dy; 182f8c2866eSHong Zhang for (i=0; i<bs; i++) { /* column of the block */ 183f8c2866eSHong Zhang for (j=0; j<bs; j++) { /* row of the block */ 184f8c2866eSHong Zhang valaddr[spidx++] = dy_i[row+j]*vscale_array[col+i]; 185f8c2866eSHong Zhang } 186f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 187f8c2866eSHong Zhang } 1880d1c53f1SHong Zhang } 1890d1c53f1SHong Zhang } 1900d1c53f1SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 1910d1c53f1SHong Zhang } 1920d1c53f1SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1930d1c53f1SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1940d1c53f1SHong Zhang if (vscale) { 1950d1c53f1SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 1960d1c53f1SHong Zhang } 1970d1c53f1SHong Zhang 1980d1c53f1SHong Zhang coloring->currentcolor = -1; 1990d1c53f1SHong Zhang PetscFunctionReturn(0); 2000d1c53f1SHong Zhang } 201a64fbb32SBarry Smith 2024a2ae208SSatish Balay #undef __FUNCT__ 203040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringApply_AIJ" 204040ebd07SHong Zhang PetscErrorCode MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 205fcd7ac73SHong Zhang { 206fcd7ac73SHong Zhang PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; 207fcd7ac73SHong Zhang PetscErrorCode ierr; 208a2f2d239SHong Zhang PetscInt k,cstart,cend,l,row,col,nz; 2099e917edbSHong Zhang PetscScalar dx=0.0,*y,*xx,*w3_array; 210fcd7ac73SHong Zhang PetscScalar *vscale_array; 211fcd7ac73SHong Zhang PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; 2128bc97078SHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; 213fcd7ac73SHong Zhang void *fctx=coloring->fctx; 214fcd7ac73SHong Zhang PetscBool flg=PETSC_FALSE; 215c53567a0SHong Zhang PetscInt ctype=coloring->ctype,nxloc,nrows_k; 216573f477fSHong Zhang MatEntry *Jentry=coloring->matentry; 2178bc97078SHong Zhang const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; 218fcd7ac73SHong Zhang 219fcd7ac73SHong Zhang PetscFunctionBegin; 220fcd7ac73SHong Zhang ierr = MatSetUnfactored(J);CHKERRQ(ierr); 221fcd7ac73SHong Zhang ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); 222fcd7ac73SHong Zhang if (flg) { 223fcd7ac73SHong Zhang ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 224fcd7ac73SHong Zhang } else { 225fcd7ac73SHong Zhang PetscBool assembled; 226fcd7ac73SHong Zhang ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 227fcd7ac73SHong Zhang if (assembled) { 228fcd7ac73SHong Zhang ierr = MatZeroEntries(J);CHKERRQ(ierr); 229fcd7ac73SHong Zhang } 230fcd7ac73SHong Zhang } 231fcd7ac73SHong Zhang 2329e917edbSHong Zhang /* create vscale for storing dx */ 2339e917edbSHong Zhang if (!vscale) { 2349e917edbSHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 2356c145f34SHong Zhang Mat_MPIAIJ *aij=(Mat_MPIAIJ*)J->data; 2369e917edbSHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr); 2379e917edbSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 2388bc97078SHong Zhang ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr); 2399e917edbSHong Zhang } 2408bc97078SHong Zhang coloring->vscale = vscale; 241fcd7ac73SHong Zhang } 242fcd7ac73SHong Zhang 2438bc97078SHong Zhang /* (1) Set w1 = F(x1) */ 244fcd7ac73SHong Zhang if (!coloring->fset) { 245fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 246f6af9589SHong Zhang ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 247fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 248fcd7ac73SHong Zhang } else { 249fcd7ac73SHong Zhang coloring->fset = PETSC_FALSE; 250fcd7ac73SHong Zhang } 251fcd7ac73SHong Zhang 2528bc97078SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 253f6af9589SHong Zhang if (coloring->htype[0] == 'w') { 2549e917edbSHong Zhang /* vscale = dx is a constant scalar */ 255f6af9589SHong Zhang ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 256c53567a0SHong Zhang dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); 25770e7395fSHong Zhang } else { 25874d3cef9SHong Zhang ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 259f6af9589SHong Zhang ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 2608bc97078SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 26174d3cef9SHong Zhang for (col=0; col<nxloc; col++) { 262fcd7ac73SHong Zhang dx = xx[col]; 26374d3cef9SHong Zhang if (PetscAbsScalar(dx) < umin) { 26474d3cef9SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 26574d3cef9SHong Zhang else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 26674d3cef9SHong Zhang } 267fcd7ac73SHong Zhang dx *= epsilon; 26874d3cef9SHong Zhang vscale_array[col] = 1.0/dx; 269f6af9589SHong Zhang } 27074d3cef9SHong Zhang ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 2718bc97078SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 27270e7395fSHong Zhang } 2739e917edbSHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] != 'w') { 2748bc97078SHong Zhang ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2758bc97078SHong Zhang ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 276fcd7ac73SHong Zhang } 277fcd7ac73SHong Zhang 2788bc97078SHong Zhang /* (3) Loop over each color */ 2798bc97078SHong Zhang if (!coloring->w3) { 2808bc97078SHong Zhang ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 2818bc97078SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 2828bc97078SHong Zhang } 2838bc97078SHong Zhang w3 = coloring->w3; 284fcd7ac73SHong Zhang 2858bc97078SHong Zhang ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 2869e917edbSHong Zhang if (vscale) { 2878bc97078SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 2889e917edbSHong Zhang } 2898bc97078SHong Zhang nz = 0; 2908bc97078SHong Zhang for (k=0; k<ncolors; k++) { 2918bc97078SHong Zhang coloring->currentcolor = k; 2929e917edbSHong Zhang 293fcd7ac73SHong Zhang /* 2948bc97078SHong Zhang (3-1) Loop over each column associated with color 295f6af9589SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 296fcd7ac73SHong Zhang */ 297f6af9589SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 298f6af9589SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 299a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 300b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 3018bc97078SHong Zhang for (l=0; l<ncolumns[k]; l++) { 302f6af9589SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 3039e917edbSHong Zhang w3_array[col] += 1.0/dx; 3049e917edbSHong Zhang } 305b039d6c4SHong Zhang } else { /* htype == 'ds' */ 306a2f2d239SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 307b039d6c4SHong Zhang for (l=0; l<ncolumns[k]; l++) { 308b039d6c4SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 309a2f2d239SHong Zhang w3_array[col] += 1.0/vscale_array[col]; 31070e7395fSHong Zhang } 311a2f2d239SHong Zhang vscale_array += cstart; 312fcd7ac73SHong Zhang } 313a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 314fcd7ac73SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 3159e917edbSHong Zhang 316fcd7ac73SHong Zhang /* 3178bc97078SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 318fcd7ac73SHong Zhang w2 = F(x1 + dx) - F(x1) 319fcd7ac73SHong Zhang */ 320fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 321fcd7ac73SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 322fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 323fcd7ac73SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 3249e917edbSHong Zhang 3258bc97078SHong Zhang /* 3268bc97078SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 3278bc97078SHong Zhang */ 328c53567a0SHong Zhang nrows_k = nrows[k]; 329fcd7ac73SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 330b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 331c53567a0SHong Zhang for (l=0; l<nrows_k; l++) { 332573f477fSHong Zhang row = Jentry[nz].row; /* local row index */ 333b039d6c4SHong Zhang *(Jentry[nz++].valaddr) = y[row]*dx; 3349e917edbSHong Zhang } 335b039d6c4SHong Zhang } else { /* htype == 'ds' */ 336c53567a0SHong Zhang for (l=0; l<nrows_k; l++) { 337b039d6c4SHong Zhang row = Jentry[nz].row; /* local row index */ 338b039d6c4SHong Zhang *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col]; 339573f477fSHong Zhang nz++; 340fcd7ac73SHong Zhang } 341b039d6c4SHong Zhang } 342fcd7ac73SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 3438bc97078SHong Zhang } 344f5aae955SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 345f5aae955SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3469e917edbSHong Zhang if (vscale) { 3478bc97078SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 3489e917edbSHong Zhang } 349fcd7ac73SHong Zhang 350fcd7ac73SHong Zhang coloring->currentcolor = -1; 351fcd7ac73SHong Zhang PetscFunctionReturn(0); 352fcd7ac73SHong Zhang } 353fcd7ac73SHong Zhang 354fcd7ac73SHong Zhang #undef __FUNCT__ 355040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIXAIJ" 356040ebd07SHong Zhang PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 357fcd7ac73SHong Zhang { 358fcd7ac73SHong Zhang PetscErrorCode ierr; 359fcd7ac73SHong Zhang PetscMPIInt size,*ncolsonproc,*disp,nn; 360*e3225e9dSHong Zhang PetscInt i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb; 36172c15787SHong Zhang const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL; 362fcd7ac73SHong Zhang PetscInt nis=iscoloring->n,nctot,*cols; 363fcd7ac73SHong Zhang IS *isa; 364fcd7ac73SHong Zhang ISLocalToGlobalMapping map=mat->cmap->mapping; 365*e3225e9dSHong Zhang PetscInt ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx; 3660d1c53f1SHong Zhang Mat A,B; 367*e3225e9dSHong Zhang PetscScalar *A_val,*B_val,**valaddrhit; 368573f477fSHong Zhang MatEntry *Jentry; 3690d1c53f1SHong Zhang PetscBool isBAIJ; 3700d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE) 3710d1c53f1SHong Zhang PetscTable colmap=NULL; 3720d1c53f1SHong Zhang #else 3730d1c53f1SHong Zhang PetscInt *colmap=NULL; /* local col number of off-diag col */ 3740d1c53f1SHong Zhang #endif 375fcd7ac73SHong Zhang 376fcd7ac73SHong Zhang PetscFunctionBegin; 37772c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 37872c15787SHong 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"); 37972c15787SHong Zhang ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr); 38072c15787SHong Zhang } 381fcd7ac73SHong Zhang 3820d1c53f1SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 3830d1c53f1SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr); 384*e3225e9dSHong Zhang if (isBAIJ) { 3850d1c53f1SHong Zhang Mat_MPIBAIJ *aij=(Mat_MPIBAIJ*)mat->data; 3866c145f34SHong Zhang Mat_SeqBAIJ *spA,*spB; 3876c145f34SHong Zhang A = aij->A; spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a; 3886c145f34SHong Zhang B = aij->B; spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a; 3890d1c53f1SHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 3900d1c53f1SHong Zhang if (!aij->colmap) { 3910d1c53f1SHong Zhang ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr); 3920d1c53f1SHong Zhang colmap = aij->colmap; 3930d1c53f1SHong Zhang } 3940d1c53f1SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 3950d1c53f1SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 396*e3225e9dSHong Zhang } else { 397*e3225e9dSHong Zhang Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data; 398*e3225e9dSHong Zhang Mat_SeqAIJ *spA,*spB; 399*e3225e9dSHong Zhang A = aij->A; spA = (Mat_SeqAIJ*)A->data; A_val = spA->a; 400*e3225e9dSHong Zhang B = aij->B; spB = (Mat_SeqAIJ*)B->data; B_val = spB->a; 401*e3225e9dSHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 402*e3225e9dSHong Zhang if (!aij->colmap) { 403*e3225e9dSHong Zhang /* Allow access to data structures of local part of matrix 404*e3225e9dSHong Zhang - creates aij->colmap which maps global column number to local number in part B */ 405*e3225e9dSHong Zhang ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 406*e3225e9dSHong Zhang colmap = aij->colmap; 407*e3225e9dSHong Zhang } 408*e3225e9dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 409*e3225e9dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 410*e3225e9dSHong Zhang 411*e3225e9dSHong Zhang bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */ 4120d1c53f1SHong Zhang } 4130d1c53f1SHong Zhang 4140d1c53f1SHong Zhang m = mat->rmap->n/bs; 4150d1c53f1SHong Zhang cstart = mat->cmap->rstart/bs; 4160d1c53f1SHong Zhang cend = mat->cmap->rend/bs; 4170d1c53f1SHong Zhang c->M = mat->rmap->N/bs; /* set the global rows and columns and local rows */ 4180d1c53f1SHong Zhang c->N = mat->cmap->N/bs; 419d3825b63SHong Zhang c->m = m; 4200d1c53f1SHong Zhang c->rstart = mat->rmap->rstart/bs; 421fcd7ac73SHong Zhang 422fcd7ac73SHong Zhang c->ncolors = nis; 423fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 424fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 425fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 4268bc97078SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 427fcd7ac73SHong Zhang 428573f477fSHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 4298bc97078SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 4308bc97078SHong Zhang c->matentry = Jentry; 431a774a6f1SHong Zhang 432d3825b63SHong Zhang ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 433fcd7ac73SHong Zhang nz = 0; 43472c15787SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 435d3825b63SHong Zhang for (i=0; i<nis; i++) { /* for each local color */ 436fcd7ac73SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 437fcd7ac73SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 438fcd7ac73SHong Zhang 439fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 440fcd7ac73SHong Zhang if (n) { 441fcd7ac73SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 442fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 443fcd7ac73SHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 444fcd7ac73SHong Zhang } else { 445fcd7ac73SHong Zhang c->columns[i] = 0; 446fcd7ac73SHong Zhang } 447fcd7ac73SHong Zhang 448fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 449d3825b63SHong Zhang /* Determine nctot, the total (parallel) number of columns of this color */ 450fcd7ac73SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 451fcd7ac73SHong Zhang ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 452fcd7ac73SHong Zhang 453d3825b63SHong Zhang /* ncolsonproc[j]: local ncolumns on proc[j] of this color */ 454fcd7ac73SHong Zhang ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 455fcd7ac73SHong Zhang ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 456fcd7ac73SHong Zhang nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 457fcd7ac73SHong Zhang if (!nctot) { 458fcd7ac73SHong Zhang ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 459fcd7ac73SHong Zhang } 460fcd7ac73SHong Zhang 461fcd7ac73SHong Zhang disp[0] = 0; 462fcd7ac73SHong Zhang for (j=1; j<size; j++) { 463fcd7ac73SHong Zhang disp[j] = disp[j-1] + ncolsonproc[j-1]; 464fcd7ac73SHong Zhang } 465fcd7ac73SHong Zhang 466d3825b63SHong Zhang /* Get cols, the complete list of columns for this color on each process */ 467fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 468fcd7ac73SHong Zhang ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 469fcd7ac73SHong Zhang ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 470fcd7ac73SHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 471fcd7ac73SHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 472fcd7ac73SHong Zhang nctot = n; 473fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 474fcd7ac73SHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 475fcd7ac73SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 476fcd7ac73SHong Zhang 4771b97d346SHong Zhang /* Mark all rows affect by these columns */ 478d3825b63SHong Zhang ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr); 4790d1c53f1SHong Zhang bs2 = bs*bs; 4804b2e90caSHong Zhang nrows_i = 0; 4811b97d346SHong Zhang for (j=0; j<nctot; j++) { /* loop over columns*/ 482fcd7ac73SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 483fcd7ac73SHong Zhang col = ltog[cols[j]]; 484fcd7ac73SHong Zhang } else { 485fcd7ac73SHong Zhang col = cols[j]; 486fcd7ac73SHong Zhang } 487*e3225e9dSHong Zhang if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */ 488d3825b63SHong Zhang row = A_cj + A_ci[col-cstart]; 489d3825b63SHong Zhang nrows = A_ci[col-cstart+1] - A_ci[col-cstart]; 4904b2e90caSHong Zhang nrows_i += nrows; 491fcd7ac73SHong Zhang /* loop over columns of A marking them in rowhit */ 492d3825b63SHong Zhang for (k=0; k<nrows; k++) { 493d3825b63SHong Zhang /* set valaddrhit for part A */ 494*e3225e9dSHong Zhang spidx = bs2*spidxA[A_ci[col-cstart] + k]; 495*e3225e9dSHong Zhang valaddrhit[*row] = &A_val[spidx]; 496a774a6f1SHong Zhang rowhit[*row++] = col - cstart + 1; /* local column index */ 497fcd7ac73SHong Zhang } 498*e3225e9dSHong Zhang } else { /* column is in B, off-diagonal block of mat */ 499fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE) 5000d1c53f1SHong Zhang ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr); 501fcd7ac73SHong Zhang colb--; 502fcd7ac73SHong Zhang #else 5030d1c53f1SHong Zhang colb = colmap[col] - 1; /* local column index */ 504fcd7ac73SHong Zhang #endif 505fcd7ac73SHong Zhang if (colb == -1) { 506d3825b63SHong Zhang nrows = 0; 507fcd7ac73SHong Zhang } else { 5080d1c53f1SHong Zhang colb = colb/bs; 509d3825b63SHong Zhang row = B_cj + B_ci[colb]; 510d3825b63SHong Zhang nrows = B_ci[colb+1] - B_ci[colb]; 511fcd7ac73SHong Zhang } 5124b2e90caSHong Zhang nrows_i += nrows; 513fcd7ac73SHong Zhang /* loop over columns of B marking them in rowhit */ 514d3825b63SHong Zhang for (k=0; k<nrows; k++) { 515d3825b63SHong Zhang /* set valaddrhit for part B */ 516*e3225e9dSHong Zhang spidx = bs2*spidxB[B_ci[colb] + k]; 517*e3225e9dSHong Zhang valaddrhit[*row] = &B_val[spidx]; 51870e7395fSHong Zhang rowhit[*row++] = colb + 1 + cend - cstart; /* local column index */ 519fcd7ac73SHong Zhang } 520fcd7ac73SHong Zhang } 521fcd7ac73SHong Zhang } 5224b2e90caSHong Zhang c->nrows[i] = nrows_i; 5238bc97078SHong Zhang 524d3825b63SHong Zhang for (j=0; j<m; j++) { 525fcd7ac73SHong Zhang if (rowhit[j]) { 526573f477fSHong Zhang Jentry[nz].row = j; /* local row index */ 527573f477fSHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 528573f477fSHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 529573f477fSHong Zhang nz++; 530fcd7ac73SHong Zhang } 531fcd7ac73SHong Zhang } 532fcd7ac73SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 533fcd7ac73SHong Zhang } 5348bc97078SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 535fcd7ac73SHong Zhang 536d3825b63SHong Zhang ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 5370d1c53f1SHong Zhang if (isBAIJ) { 5380d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 5390d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 540*e3225e9dSHong Zhang ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 5410d1c53f1SHong Zhang } else { 5420d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 5430d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 5440d1c53f1SHong Zhang } 5450d1c53f1SHong Zhang 54672c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 54772c15787SHong Zhang ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr); 54872c15787SHong Zhang } 549fcd7ac73SHong Zhang PetscFunctionReturn(0); 550fcd7ac73SHong Zhang } 551