1a64fbb32SBarry Smith 2c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 30d1c53f1SHong Zhang #include <../src/mat/impls/baij/mpi/mpibaij.h> 4af0996ceSBarry Smith #include <petsc/private/isimpl.h> 50d1c53f1SHong Zhang 60d1c53f1SHong Zhang #undef __FUNCT__ 7040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringApply_BAIJ" 8d1e9a80fSBarry Smith PetscErrorCode MatFDColoringApply_BAIJ(Mat J,MatFDColoring coloring,Vec x1,void *sctx) 90d1c53f1SHong Zhang { 100d1c53f1SHong Zhang PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; 110d1c53f1SHong Zhang PetscErrorCode ierr; 120d1c53f1SHong Zhang PetscInt k,cstart,cend,l,row,col,nz,spidx,i,j; 135edff71fSBarry Smith PetscScalar dx=0.0,*w3_array,*dy_i,*dy=coloring->dy; 140d1c53f1SHong Zhang PetscScalar *vscale_array; 155edff71fSBarry Smith const PetscScalar *xx; 160d1c53f1SHong Zhang PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; 170d1c53f1SHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; 180d1c53f1SHong Zhang void *fctx=coloring->fctx; 190d1c53f1SHong Zhang PetscInt ctype=coloring->ctype,nxloc,nrows_k; 200d1c53f1SHong Zhang PetscScalar *valaddr; 210d1c53f1SHong Zhang MatEntry *Jentry=coloring->matentry; 220df34763SHong Zhang MatEntry2 *Jentry2=coloring->matentry2; 230d1c53f1SHong Zhang const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; 240d1c53f1SHong Zhang PetscInt bs=J->rmap->bs; 250d1c53f1SHong Zhang 260d1c53f1SHong Zhang PetscFunctionBegin; 270d1c53f1SHong Zhang /* (1) Set w1 = F(x1) */ 280d1c53f1SHong Zhang if (!coloring->fset) { 290d1c53f1SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 300d1c53f1SHong Zhang ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 310d1c53f1SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 320d1c53f1SHong Zhang } else { 330d1c53f1SHong Zhang coloring->fset = PETSC_FALSE; 340d1c53f1SHong Zhang } 350d1c53f1SHong Zhang 360d1c53f1SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 370d1c53f1SHong Zhang ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 380d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 390d1c53f1SHong Zhang /* vscale = dx is a constant scalar */ 400d1c53f1SHong Zhang ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 410d1c53f1SHong Zhang dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); 420d1c53f1SHong Zhang } else { 435edff71fSBarry Smith ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr); 440d1c53f1SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 450d1c53f1SHong Zhang for (col=0; col<nxloc; col++) { 460d1c53f1SHong Zhang dx = xx[col]; 470d1c53f1SHong Zhang if (PetscAbsScalar(dx) < umin) { 480d1c53f1SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 490d1c53f1SHong Zhang else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 500d1c53f1SHong Zhang } 510d1c53f1SHong Zhang dx *= epsilon; 520d1c53f1SHong Zhang vscale_array[col] = 1.0/dx; 530d1c53f1SHong Zhang } 545edff71fSBarry Smith ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr); 550d1c53f1SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 560d1c53f1SHong Zhang } 57684f2004SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 580d1c53f1SHong Zhang ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 590d1c53f1SHong Zhang ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 600d1c53f1SHong Zhang } 610d1c53f1SHong Zhang 620d1c53f1SHong Zhang /* (3) Loop over each color */ 630d1c53f1SHong Zhang if (!coloring->w3) { 640d1c53f1SHong Zhang ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 650d1c53f1SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 660d1c53f1SHong Zhang } 670d1c53f1SHong Zhang w3 = coloring->w3; 680d1c53f1SHong Zhang 690d1c53f1SHong Zhang ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 700d1c53f1SHong Zhang if (vscale) { 710d1c53f1SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 720d1c53f1SHong Zhang } 730d1c53f1SHong Zhang nz = 0; 740d1c53f1SHong Zhang for (k=0; k<ncolors; k++) { 750d1c53f1SHong Zhang coloring->currentcolor = k; 760d1c53f1SHong Zhang 770d1c53f1SHong Zhang /* 780d1c53f1SHong Zhang (3-1) Loop over each column associated with color 790d1c53f1SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 800d1c53f1SHong Zhang */ 810d1c53f1SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 820d1c53f1SHong Zhang dy_i = dy; 830d1c53f1SHong Zhang for (i=0; i<bs; i++) { /* Loop over a block of columns */ 840d1c53f1SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 850d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 860d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 870d1c53f1SHong Zhang for (l=0; l<ncolumns[k]; l++) { 880d1c53f1SHong Zhang col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 890d1c53f1SHong Zhang w3_array[col] += 1.0/dx; 90684f2004SHong Zhang if (i) w3_array[col-1] -= 1.0/dx; /* resume original w3[col-1] */ 910d1c53f1SHong Zhang } 920d1c53f1SHong Zhang } else { /* htype == 'ds' */ 930d1c53f1SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 940d1c53f1SHong Zhang for (l=0; l<ncolumns[k]; l++) { 95f8c2866eSHong Zhang col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 960d1c53f1SHong Zhang w3_array[col] += 1.0/vscale_array[col]; 97684f2004SHong Zhang if (i) w3_array[col-1] -= 1.0/vscale_array[col-1]; /* resume original w3[col-1] */ 980d1c53f1SHong Zhang } 990d1c53f1SHong Zhang vscale_array += cstart; 1000d1c53f1SHong Zhang } 1010d1c53f1SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 1020d1c53f1SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 1030d1c53f1SHong Zhang 1040d1c53f1SHong Zhang /* 1050d1c53f1SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 1060d1c53f1SHong Zhang w2 = F(x1 + dx) - F(x1) 1070d1c53f1SHong Zhang */ 1080d1c53f1SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 1090d1c53f1SHong Zhang ierr = VecPlaceArray(w2,dy_i);CHKERRQ(ierr); /* place w2 to the array dy_i */ 1100d1c53f1SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 1110d1c53f1SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 1120d1c53f1SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 1130d1c53f1SHong Zhang ierr = VecResetArray(w2);CHKERRQ(ierr); 1140d1c53f1SHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1150d1c53f1SHong Zhang } 1160d1c53f1SHong Zhang 1170d1c53f1SHong Zhang /* 1180d1c53f1SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 1190d1c53f1SHong Zhang */ 1200d1c53f1SHong Zhang nrows_k = nrows[k]; 1210d1c53f1SHong Zhang if (coloring->htype[0] == 'w') { 1220d1c53f1SHong Zhang for (l=0; l<nrows_k; l++) { 1230df34763SHong Zhang row = bs*Jentry2[nz].row; /* local row index */ 1240df34763SHong Zhang valaddr = Jentry2[nz++].valaddr; 1250d1c53f1SHong Zhang spidx = 0; 1260d1c53f1SHong Zhang dy_i = dy; 1270d1c53f1SHong Zhang for (i=0; i<bs; i++) { /* column of the block */ 1280d1c53f1SHong Zhang for (j=0; j<bs; j++) { /* row of the block */ 1290d1c53f1SHong Zhang valaddr[spidx++] = dy_i[row+j]*dx; 1300d1c53f1SHong Zhang } 131f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 1320d1c53f1SHong Zhang } 1330d1c53f1SHong Zhang } 1340d1c53f1SHong Zhang } else { /* htype == 'ds' */ 1350d1c53f1SHong Zhang for (l=0; l<nrows_k; l++) { 1360d1c53f1SHong Zhang row = bs*Jentry[nz].row; /* local row index */ 1370d1c53f1SHong Zhang col = bs*Jentry[nz].col; /* local column index */ 138684f2004SHong Zhang valaddr = Jentry[nz++].valaddr; 139f8c2866eSHong Zhang spidx = 0; 140f8c2866eSHong Zhang dy_i = dy; 141f8c2866eSHong Zhang for (i=0; i<bs; i++) { /* column of the block */ 142f8c2866eSHong Zhang for (j=0; j<bs; j++) { /* row of the block */ 143f8c2866eSHong Zhang valaddr[spidx++] = dy_i[row+j]*vscale_array[col+i]; 144f8c2866eSHong Zhang } 145f8c2866eSHong Zhang dy_i += nxloc; /* points to dy+i*nxloc */ 146f8c2866eSHong Zhang } 1470d1c53f1SHong Zhang } 1480d1c53f1SHong Zhang } 1490d1c53f1SHong Zhang } 1500d1c53f1SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1510d1c53f1SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1520d1c53f1SHong Zhang if (vscale) { 1530d1c53f1SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 1540d1c53f1SHong Zhang } 1550d1c53f1SHong Zhang 1560d1c53f1SHong Zhang coloring->currentcolor = -1; 1570d1c53f1SHong Zhang PetscFunctionReturn(0); 1580d1c53f1SHong Zhang } 159a64fbb32SBarry Smith 1604a2ae208SSatish Balay #undef __FUNCT__ 161040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringApply_AIJ" 162d1e9a80fSBarry Smith PetscErrorCode MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,void *sctx) 163fcd7ac73SHong Zhang { 164fcd7ac73SHong Zhang PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; 165fcd7ac73SHong Zhang PetscErrorCode ierr; 166a2f2d239SHong Zhang PetscInt k,cstart,cend,l,row,col,nz; 1675edff71fSBarry Smith PetscScalar dx=0.0,*y,*w3_array; 1685edff71fSBarry Smith const PetscScalar *xx; 169fcd7ac73SHong Zhang PetscScalar *vscale_array; 170fcd7ac73SHong Zhang PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; 1718bc97078SHong Zhang Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; 172fcd7ac73SHong Zhang void *fctx=coloring->fctx; 173c53567a0SHong Zhang PetscInt ctype=coloring->ctype,nxloc,nrows_k; 174573f477fSHong Zhang MatEntry *Jentry=coloring->matentry; 1750df34763SHong Zhang MatEntry2 *Jentry2=coloring->matentry2; 1768bc97078SHong Zhang const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; 177fcd7ac73SHong Zhang 178fcd7ac73SHong Zhang PetscFunctionBegin; 1798bc97078SHong Zhang /* (1) Set w1 = F(x1) */ 180fcd7ac73SHong Zhang if (!coloring->fset) { 181fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 182f6af9589SHong Zhang ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 183fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 184fcd7ac73SHong Zhang } else { 185fcd7ac73SHong Zhang coloring->fset = PETSC_FALSE; 186fcd7ac73SHong Zhang } 187fcd7ac73SHong Zhang 1888bc97078SHong Zhang /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ 189f6af9589SHong Zhang if (coloring->htype[0] == 'w') { 190684f2004SHong Zhang /* vscale = 1./dx is a constant scalar */ 191f6af9589SHong Zhang ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 192c53567a0SHong Zhang dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); 19370e7395fSHong Zhang } else { 19474d3cef9SHong Zhang ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 1955edff71fSBarry Smith ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr); 1968bc97078SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 19774d3cef9SHong Zhang for (col=0; col<nxloc; col++) { 198fcd7ac73SHong Zhang dx = xx[col]; 19974d3cef9SHong Zhang if (PetscAbsScalar(dx) < umin) { 20074d3cef9SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 20174d3cef9SHong Zhang else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 20274d3cef9SHong Zhang } 203fcd7ac73SHong Zhang dx *= epsilon; 20474d3cef9SHong Zhang vscale_array[col] = 1.0/dx; 205f6af9589SHong Zhang } 2065edff71fSBarry Smith ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr); 2078bc97078SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 20870e7395fSHong Zhang } 209684f2004SHong Zhang if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { 2108bc97078SHong Zhang ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2118bc97078SHong Zhang ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 212fcd7ac73SHong Zhang } 213fcd7ac73SHong Zhang 2148bc97078SHong Zhang /* (3) Loop over each color */ 2158bc97078SHong Zhang if (!coloring->w3) { 2168bc97078SHong Zhang ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 2178bc97078SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 2188bc97078SHong Zhang } 2198bc97078SHong Zhang w3 = coloring->w3; 220fcd7ac73SHong Zhang 2218bc97078SHong Zhang ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ 2229e917edbSHong Zhang if (vscale) { 2238bc97078SHong Zhang ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); 2249e917edbSHong Zhang } 2258bc97078SHong Zhang nz = 0; 226e2c857f8SHong Zhang 227b3e1f37bSHong Zhang if (coloring->bcols > 1) { /* use blocked insertion of Jentry */ 228b3e1f37bSHong Zhang PetscInt i,m=J->rmap->n,nbcols,bcols=coloring->bcols; 229e2c857f8SHong Zhang PetscScalar *dy=coloring->dy,*dy_k; 230e2c857f8SHong Zhang 231e2c857f8SHong Zhang nbcols = 0; 232e2c857f8SHong Zhang for (k=0; k<ncolors; k+=bcols) { 233e2c857f8SHong Zhang coloring->currentcolor = k; 234e2c857f8SHong Zhang 235e2c857f8SHong Zhang /* 236e2c857f8SHong Zhang (3-1) Loop over each column associated with color 237e2c857f8SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 238e2c857f8SHong Zhang */ 239e2c857f8SHong Zhang 240e2c857f8SHong Zhang dy_k = dy; 241e2c857f8SHong Zhang if (k + bcols > ncolors) bcols = ncolors - k; 242e2c857f8SHong Zhang for (i=0; i<bcols; i++) { 243e2c857f8SHong Zhang 244e2c857f8SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 245e2c857f8SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 246e2c857f8SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 247e2c857f8SHong Zhang if (coloring->htype[0] == 'w') { 248e2c857f8SHong Zhang for (l=0; l<ncolumns[k+i]; l++) { 249e2c857f8SHong Zhang col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */ 250e2c857f8SHong Zhang w3_array[col] += 1.0/dx; 251e2c857f8SHong Zhang } 252e2c857f8SHong Zhang } else { /* htype == 'ds' */ 253e2c857f8SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 254e2c857f8SHong Zhang for (l=0; l<ncolumns[k+i]; l++) { 255e2c857f8SHong Zhang col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */ 256e2c857f8SHong Zhang w3_array[col] += 1.0/vscale_array[col]; 257e2c857f8SHong Zhang } 258e2c857f8SHong Zhang vscale_array += cstart; 259e2c857f8SHong Zhang } 260e2c857f8SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 261e2c857f8SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 262e2c857f8SHong Zhang 263e2c857f8SHong Zhang /* 264e2c857f8SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 265e2c857f8SHong Zhang w2 = F(x1 + dx) - F(x1) 266e2c857f8SHong Zhang */ 267e2c857f8SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 268e2c857f8SHong Zhang ierr = VecPlaceArray(w2,dy_k);CHKERRQ(ierr); /* place w2 to the array dy_i */ 269e2c857f8SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 270e2c857f8SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 271e2c857f8SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 272e2c857f8SHong Zhang ierr = VecResetArray(w2);CHKERRQ(ierr); 273e2c857f8SHong Zhang dy_k += m; /* points to dy+i*nxloc */ 274e2c857f8SHong Zhang } 275e2c857f8SHong Zhang 276e2c857f8SHong Zhang /* 277e2c857f8SHong Zhang (3-3) Loop over block rows of vector, putting results into Jacobian matrix 278e2c857f8SHong Zhang */ 279d880da65SHong Zhang nrows_k = nrows[nbcols++]; 280e2c857f8SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 281d880da65SHong Zhang 282e2c857f8SHong Zhang if (coloring->htype[0] == 'w') { 283e2c857f8SHong Zhang for (l=0; l<nrows_k; l++) { 2840df34763SHong Zhang row = Jentry2[nz].row; /* local row index */ 2850df34763SHong Zhang *(Jentry2[nz++].valaddr) = dy[row]*dx; 286e2c857f8SHong Zhang } 287e2c857f8SHong Zhang } else { /* htype == 'ds' */ 288e2c857f8SHong Zhang for (l=0; l<nrows_k; l++) { 289e2c857f8SHong Zhang row = Jentry[nz].row; /* local row index */ 290d880da65SHong Zhang *(Jentry[nz].valaddr) = dy[row]*vscale_array[Jentry[nz].col]; 291e2c857f8SHong Zhang nz++; 292e2c857f8SHong Zhang } 293e2c857f8SHong Zhang } 294e2c857f8SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 295e2c857f8SHong Zhang } 296b3e1f37bSHong Zhang } else { /* bcols == 1 */ 2978bc97078SHong Zhang for (k=0; k<ncolors; k++) { 2988bc97078SHong Zhang coloring->currentcolor = k; 2999e917edbSHong Zhang 300fcd7ac73SHong Zhang /* 3018bc97078SHong Zhang (3-1) Loop over each column associated with color 302f6af9589SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 303fcd7ac73SHong Zhang */ 304f6af9589SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 305f6af9589SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 306a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ 307b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 3088bc97078SHong Zhang for (l=0; l<ncolumns[k]; l++) { 309f6af9589SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 3109e917edbSHong Zhang w3_array[col] += 1.0/dx; 3119e917edbSHong Zhang } 312b039d6c4SHong Zhang } else { /* htype == 'ds' */ 313a2f2d239SHong Zhang vscale_array -= cstart; /* shift pointer so global index can be used */ 314b039d6c4SHong Zhang for (l=0; l<ncolumns[k]; l++) { 315b039d6c4SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 316a2f2d239SHong Zhang w3_array[col] += 1.0/vscale_array[col]; 31770e7395fSHong Zhang } 318a2f2d239SHong Zhang vscale_array += cstart; 319fcd7ac73SHong Zhang } 320a2f2d239SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; 321fcd7ac73SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 3229e917edbSHong Zhang 323fcd7ac73SHong Zhang /* 3248bc97078SHong Zhang (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 325fcd7ac73SHong Zhang w2 = F(x1 + dx) - F(x1) 326fcd7ac73SHong Zhang */ 327fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 328fcd7ac73SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 329fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 330fcd7ac73SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 3319e917edbSHong Zhang 3328bc97078SHong Zhang /* 3338bc97078SHong Zhang (3-3) Loop over rows of vector, putting results into Jacobian matrix 3348bc97078SHong Zhang */ 335c53567a0SHong Zhang nrows_k = nrows[k]; 336fcd7ac73SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 337b039d6c4SHong Zhang if (coloring->htype[0] == 'w') { 338c53567a0SHong Zhang for (l=0; l<nrows_k; l++) { 3390df34763SHong Zhang row = Jentry2[nz].row; /* local row index */ 3400df34763SHong Zhang *(Jentry2[nz++].valaddr) = y[row]*dx; 3419e917edbSHong Zhang } 342b039d6c4SHong Zhang } else { /* htype == 'ds' */ 343c53567a0SHong Zhang for (l=0; l<nrows_k; l++) { 344b039d6c4SHong Zhang row = Jentry[nz].row; /* local row index */ 345b039d6c4SHong Zhang *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col]; 346573f477fSHong Zhang nz++; 347fcd7ac73SHong Zhang } 348b039d6c4SHong Zhang } 349fcd7ac73SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 3508bc97078SHong Zhang } 351b3e1f37bSHong Zhang } 352b3e1f37bSHong Zhang 353f5aae955SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 354f5aae955SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3559e917edbSHong Zhang if (vscale) { 3568bc97078SHong Zhang ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); 3579e917edbSHong Zhang } 358fcd7ac73SHong Zhang coloring->currentcolor = -1; 359fcd7ac73SHong Zhang PetscFunctionReturn(0); 360fcd7ac73SHong Zhang } 361fcd7ac73SHong Zhang 362fcd7ac73SHong Zhang #undef __FUNCT__ 363f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetUp_MPIXAIJ" 364f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 365fcd7ac73SHong Zhang { 366fcd7ac73SHong Zhang PetscErrorCode ierr; 367fcd7ac73SHong Zhang PetscMPIInt size,*ncolsonproc,*disp,nn; 368e3225e9dSHong Zhang PetscInt i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb; 36972c15787SHong Zhang const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL; 370fcd7ac73SHong Zhang PetscInt nis=iscoloring->n,nctot,*cols; 371fcd7ac73SHong Zhang IS *isa; 372fcd7ac73SHong Zhang ISLocalToGlobalMapping map=mat->cmap->mapping; 373a8971b87SHong Zhang PetscInt ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx; 3740d1c53f1SHong Zhang Mat A,B; 375e3225e9dSHong Zhang PetscScalar *A_val,*B_val,**valaddrhit; 376573f477fSHong Zhang MatEntry *Jentry; 3770df34763SHong Zhang MatEntry2 *Jentry2; 3780d1c53f1SHong Zhang PetscBool isBAIJ; 379531e53bdSHong Zhang PetscInt bcols=c->bcols; 3800d1c53f1SHong Zhang #if defined(PETSC_USE_CTABLE) 3810d1c53f1SHong Zhang PetscTable colmap=NULL; 3820d1c53f1SHong Zhang #else 3830d1c53f1SHong Zhang PetscInt *colmap=NULL; /* local col number of off-diag col */ 3840d1c53f1SHong Zhang #endif 385fcd7ac73SHong Zhang 386fcd7ac73SHong Zhang PetscFunctionBegin; 38772c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 38872c15787SHong 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"); 38972c15787SHong Zhang ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr); 39072c15787SHong Zhang } 391fcd7ac73SHong Zhang 3920d1c53f1SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 3930d1c53f1SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr); 394e3225e9dSHong Zhang if (isBAIJ) { 395195687c5SHong Zhang Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)mat->data; 3966c145f34SHong Zhang Mat_SeqBAIJ *spA,*spB; 397195687c5SHong Zhang A = baij->A; spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a; 398195687c5SHong Zhang B = baij->B; spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a; 3990d1c53f1SHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 400195687c5SHong Zhang if (!baij->colmap) { 4010d1c53f1SHong Zhang ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr); 402195687c5SHong Zhang colmap = baij->colmap; 4030d1c53f1SHong Zhang } 4040d1c53f1SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 4050d1c53f1SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 406195687c5SHong Zhang 407195687c5SHong Zhang if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */ 408195687c5SHong Zhang PetscInt *garray; 409785e854fSJed Brown ierr = PetscMalloc1(B->cmap->n,&garray);CHKERRQ(ierr); 410195687c5SHong Zhang for (i=0; i<baij->B->cmap->n/bs; i++) { 411195687c5SHong Zhang for (j=0; j<bs; j++) { 412195687c5SHong Zhang garray[i*bs+j] = bs*baij->garray[i]+j; 413195687c5SHong Zhang } 414195687c5SHong Zhang } 415195687c5SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,garray,&c->vscale);CHKERRQ(ierr); 416195687c5SHong Zhang ierr = PetscFree(garray);CHKERRQ(ierr); 417195687c5SHong Zhang } 418e3225e9dSHong Zhang } else { 419e3225e9dSHong Zhang Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data; 420e3225e9dSHong Zhang Mat_SeqAIJ *spA,*spB; 421e3225e9dSHong Zhang A = aij->A; spA = (Mat_SeqAIJ*)A->data; A_val = spA->a; 422e3225e9dSHong Zhang B = aij->B; spB = (Mat_SeqAIJ*)B->data; B_val = spB->a; 423e3225e9dSHong Zhang nz = spA->nz + spB->nz; /* total nonzero entries of mat */ 424e3225e9dSHong Zhang if (!aij->colmap) { 425e3225e9dSHong Zhang /* Allow access to data structures of local part of matrix 426e3225e9dSHong Zhang - creates aij->colmap which maps global column number to local number in part B */ 427e3225e9dSHong Zhang ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 428e3225e9dSHong Zhang colmap = aij->colmap; 429e3225e9dSHong Zhang } 430e3225e9dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 431e3225e9dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 432e3225e9dSHong Zhang 433e3225e9dSHong Zhang bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */ 434195687c5SHong Zhang 435195687c5SHong Zhang if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */ 436195687c5SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 437195687c5SHong Zhang } 4380d1c53f1SHong Zhang } 4390d1c53f1SHong Zhang 4400d1c53f1SHong Zhang m = mat->rmap->n/bs; 4410d1c53f1SHong Zhang cstart = mat->cmap->rstart/bs; 4420d1c53f1SHong Zhang cend = mat->cmap->rend/bs; 443fcd7ac73SHong Zhang 444785e854fSJed Brown ierr = PetscMalloc1(nis,&c->ncolumns);CHKERRQ(ierr); 445785e854fSJed Brown ierr = PetscMalloc1(nis,&c->columns);CHKERRQ(ierr); 446785e854fSJed Brown ierr = PetscMalloc1(nis,&c->nrows);CHKERRQ(ierr); 4478bc97078SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 448fcd7ac73SHong Zhang 4490df34763SHong Zhang if (c->htype[0] == 'd') { 450785e854fSJed Brown ierr = PetscMalloc1(nz,&Jentry);CHKERRQ(ierr); 4518bc97078SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 4528bc97078SHong Zhang c->matentry = Jentry; 4530df34763SHong Zhang } else if (c->htype[0] == 'w') { 454785e854fSJed Brown ierr = PetscMalloc1(nz,&Jentry2);CHKERRQ(ierr); 4550df34763SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry2));CHKERRQ(ierr); 4560df34763SHong Zhang c->matentry2 = Jentry2; 4570df34763SHong Zhang } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"htype is not supported"); 458a774a6f1SHong Zhang 459dcca6d9dSJed Brown ierr = PetscMalloc2(m+1,&rowhit,m+1,&valaddrhit);CHKERRQ(ierr); 460fcd7ac73SHong Zhang nz = 0; 46172c15787SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 462d3825b63SHong Zhang for (i=0; i<nis; i++) { /* for each local color */ 463fcd7ac73SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 464fcd7ac73SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 465fcd7ac73SHong Zhang 466fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 467fcd7ac73SHong Zhang if (n) { 468785e854fSJed Brown ierr = PetscMalloc1(n,&c->columns[i]);CHKERRQ(ierr); 469fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 470fcd7ac73SHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 471fcd7ac73SHong Zhang } else { 472fcd7ac73SHong Zhang c->columns[i] = 0; 473fcd7ac73SHong Zhang } 474fcd7ac73SHong Zhang 475fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 476d3825b63SHong Zhang /* Determine nctot, the total (parallel) number of columns of this color */ 477fcd7ac73SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 478dcca6d9dSJed Brown ierr = PetscMalloc2(size,&ncolsonproc,size,&disp);CHKERRQ(ierr); 479fcd7ac73SHong Zhang 480d3825b63SHong Zhang /* ncolsonproc[j]: local ncolumns on proc[j] of this color */ 481fcd7ac73SHong Zhang ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 482fcd7ac73SHong Zhang ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 483fcd7ac73SHong Zhang nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 484fcd7ac73SHong Zhang if (!nctot) { 485fcd7ac73SHong Zhang ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 486fcd7ac73SHong Zhang } 487fcd7ac73SHong Zhang 488fcd7ac73SHong Zhang disp[0] = 0; 489fcd7ac73SHong Zhang for (j=1; j<size; j++) { 490fcd7ac73SHong Zhang disp[j] = disp[j-1] + ncolsonproc[j-1]; 491fcd7ac73SHong Zhang } 492fcd7ac73SHong Zhang 493d3825b63SHong Zhang /* Get cols, the complete list of columns for this color on each process */ 494854ce69bSBarry Smith ierr = PetscMalloc1(nctot+1,&cols);CHKERRQ(ierr); 495fcd7ac73SHong Zhang ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 496fcd7ac73SHong Zhang ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 497fcd7ac73SHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 498fcd7ac73SHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 499fcd7ac73SHong Zhang nctot = n; 500854ce69bSBarry Smith ierr = PetscMalloc1(nctot+1,&cols);CHKERRQ(ierr); 501fcd7ac73SHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 502fcd7ac73SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 503fcd7ac73SHong Zhang 5041b97d346SHong Zhang /* Mark all rows affect by these columns */ 505d3825b63SHong Zhang ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr); 5060d1c53f1SHong Zhang bs2 = bs*bs; 5074b2e90caSHong Zhang nrows_i = 0; 5081b97d346SHong Zhang for (j=0; j<nctot; j++) { /* loop over columns*/ 509fcd7ac73SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 510fcd7ac73SHong Zhang col = ltog[cols[j]]; 511fcd7ac73SHong Zhang } else { 512fcd7ac73SHong Zhang col = cols[j]; 513fcd7ac73SHong Zhang } 514e3225e9dSHong Zhang if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */ 515d3825b63SHong Zhang row = A_cj + A_ci[col-cstart]; 516d3825b63SHong Zhang nrows = A_ci[col-cstart+1] - A_ci[col-cstart]; 5174b2e90caSHong Zhang nrows_i += nrows; 518fcd7ac73SHong Zhang /* loop over columns of A marking them in rowhit */ 519d3825b63SHong Zhang for (k=0; k<nrows; k++) { 520d3825b63SHong Zhang /* set valaddrhit for part A */ 521e3225e9dSHong Zhang spidx = bs2*spidxA[A_ci[col-cstart] + k]; 522e3225e9dSHong Zhang valaddrhit[*row] = &A_val[spidx]; 523a774a6f1SHong Zhang rowhit[*row++] = col - cstart + 1; /* local column index */ 524fcd7ac73SHong Zhang } 525e3225e9dSHong Zhang } else { /* column is in B, off-diagonal block of mat */ 526fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE) 5270d1c53f1SHong Zhang ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr); 528fcd7ac73SHong Zhang colb--; 529fcd7ac73SHong Zhang #else 5300d1c53f1SHong Zhang colb = colmap[col] - 1; /* local column index */ 531fcd7ac73SHong Zhang #endif 532fcd7ac73SHong Zhang if (colb == -1) { 533d3825b63SHong Zhang nrows = 0; 534fcd7ac73SHong Zhang } else { 5350d1c53f1SHong Zhang colb = colb/bs; 536d3825b63SHong Zhang row = B_cj + B_ci[colb]; 537d3825b63SHong Zhang nrows = B_ci[colb+1] - B_ci[colb]; 538fcd7ac73SHong Zhang } 5394b2e90caSHong Zhang nrows_i += nrows; 540fcd7ac73SHong Zhang /* loop over columns of B marking them in rowhit */ 541d3825b63SHong Zhang for (k=0; k<nrows; k++) { 542d3825b63SHong Zhang /* set valaddrhit for part B */ 543e3225e9dSHong Zhang spidx = bs2*spidxB[B_ci[colb] + k]; 544e3225e9dSHong Zhang valaddrhit[*row] = &B_val[spidx]; 54570e7395fSHong Zhang rowhit[*row++] = colb + 1 + cend - cstart; /* local column index */ 546fcd7ac73SHong Zhang } 547fcd7ac73SHong Zhang } 548fcd7ac73SHong Zhang } 5494b2e90caSHong Zhang c->nrows[i] = nrows_i; 5508bc97078SHong Zhang 5510df34763SHong Zhang if (c->htype[0] == 'd') { 552d3825b63SHong Zhang for (j=0; j<m; j++) { 553fcd7ac73SHong Zhang if (rowhit[j]) { 554573f477fSHong Zhang Jentry[nz].row = j; /* local row index */ 555573f477fSHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 556573f477fSHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 557573f477fSHong Zhang nz++; 558fcd7ac73SHong Zhang } 559fcd7ac73SHong Zhang } 5600df34763SHong Zhang } else { /* c->htype == 'wp' */ 5610df34763SHong Zhang for (j=0; j<m; j++) { 5620df34763SHong Zhang if (rowhit[j]) { 5630df34763SHong Zhang Jentry2[nz].row = j; /* local row index */ 5640df34763SHong Zhang Jentry2[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 5650df34763SHong Zhang nz++; 5660df34763SHong Zhang } 5670df34763SHong Zhang } 5680df34763SHong Zhang } 569fcd7ac73SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 570fcd7ac73SHong Zhang } 571fcd7ac73SHong Zhang 572a8971b87SHong Zhang if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */ 573a8971b87SHong Zhang ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr); 574531e53bdSHong Zhang } 575531e53bdSHong Zhang 5760d1c53f1SHong Zhang if (isBAIJ) { 5770d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 5780d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 579785e854fSJed Brown ierr = PetscMalloc1(bs*mat->rmap->n,&c->dy);CHKERRQ(ierr); 5800d1c53f1SHong Zhang } else { 5810d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 5820d1c53f1SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 5830d1c53f1SHong Zhang } 5840d1c53f1SHong Zhang 585a8971b87SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 586a8971b87SHong Zhang ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 587a8971b87SHong Zhang 58872c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 58972c15787SHong Zhang ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr); 59072c15787SHong Zhang } 591b3e1f37bSHong Zhang ierr = PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr); 592fcd7ac73SHong Zhang PetscFunctionReturn(0); 593fcd7ac73SHong Zhang } 59493dfae19SHong Zhang 59593dfae19SHong Zhang #undef __FUNCT__ 59693dfae19SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIXAIJ" 59793dfae19SHong Zhang PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 59893dfae19SHong Zhang { 599531e53bdSHong Zhang PetscErrorCode ierr; 600a8971b87SHong Zhang PetscInt bs,nis=iscoloring->n,m=mat->rmap->n; 601531e53bdSHong Zhang PetscBool isBAIJ; 602531e53bdSHong Zhang 60393dfae19SHong Zhang PetscFunctionBegin; 604531e53bdSHong Zhang /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian; 605531e53bdSHong Zhang bcols is chosen s.t. dy-array takes 50% of memory space as mat */ 606531e53bdSHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 607531e53bdSHong Zhang ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr); 608*b8dfa574SBarry Smith if (isBAIJ || m == 0) { 609a8971b87SHong Zhang c->brows = m; 610a8971b87SHong Zhang c->bcols = 1; 611531e53bdSHong Zhang } else { /* mpiaij matrix */ 612531e53bdSHong Zhang /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */ 613531e53bdSHong Zhang Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data; 614531e53bdSHong Zhang Mat_SeqAIJ *spA,*spB; 615531e53bdSHong Zhang Mat A,B; 616a8971b87SHong Zhang PetscInt nz,brows,bcols; 617531e53bdSHong Zhang PetscReal mem; 618531e53bdSHong Zhang 619531e53bdSHong Zhang bs = 1; /* only bs=1 is supported for MPIAIJ matrix */ 620531e53bdSHong Zhang 621531e53bdSHong Zhang A = aij->A; spA = (Mat_SeqAIJ*)A->data; 622531e53bdSHong Zhang B = aij->B; spB = (Mat_SeqAIJ*)B->data; 623531e53bdSHong Zhang nz = spA->nz + spB->nz; /* total local nonzero entries of mat */ 624531e53bdSHong Zhang mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt); 625531e53bdSHong Zhang bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar))); 626531e53bdSHong Zhang brows = 1000/bcols; 627531e53bdSHong Zhang if (bcols > nis) bcols = nis; 628531e53bdSHong Zhang if (brows == 0 || brows > m) brows = m; 629531e53bdSHong Zhang c->brows = brows; 630531e53bdSHong Zhang c->bcols = bcols; 631531e53bdSHong Zhang } 632a8971b87SHong Zhang 633a8971b87SHong Zhang c->M = mat->rmap->N/bs; /* set the global rows and columns and local rows */ 634a8971b87SHong Zhang c->N = mat->cmap->N/bs; 635a8971b87SHong Zhang c->m = mat->rmap->n/bs; 636a8971b87SHong Zhang c->rstart = mat->rmap->rstart/bs; 637a8971b87SHong Zhang c->ncolors = nis; 63893dfae19SHong Zhang PetscFunctionReturn(0); 63993dfae19SHong Zhang } 640