191627157SBarry Smith 291627157SBarry Smith #ifndef lint 3*639f9d9dSBarry Smith static char vcid[] = "$Id: snesj2.c,v 1.1 1996/10/09 19:38:25 bsmith Exp bsmith $"; 491627157SBarry Smith #endif 591627157SBarry Smith 6*639f9d9dSBarry Smith #include "src/mat/matimpl.h" /*I "mat.h" I*/ 7*639f9d9dSBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 891627157SBarry Smith 991627157SBarry Smith 1091627157SBarry Smith /*@C 1191627157SBarry Smith SNESDefaultComputeJacobianWithColoring 1291627157SBarry Smith 1391627157SBarry Smith Input Parameters: 1491627157SBarry Smith . snes - nonlinear solver object 1591627157SBarry Smith . x1 - location at which to evaluate Jacobian 1691627157SBarry Smith . ctx - MatFDColoring contex 1791627157SBarry Smith 1891627157SBarry Smith Output Parameters: 1991627157SBarry Smith . J - Jacobian matrix 2091627157SBarry Smith . B - Jacobian preconditioner 2191627157SBarry Smith . flag - flag indicating if the matrix nonzero structure has changed 2291627157SBarry Smith 2391627157SBarry Smith .keywords: SNES, finite differences, Jacobian 2491627157SBarry Smith 2591627157SBarry Smith .seealso: SNESSetJacobian(), SNESTestJacobian() 2691627157SBarry Smith @*/ 2791627157SBarry Smith int SNESDefaultComputeJacobianWithColoring(SNES snes,Vec x1,Mat *JJ,Mat *B,MatStructure *flag,void *ctx) 2891627157SBarry Smith { 2991627157SBarry Smith MatFDColoring color = (MatFDColoring) ctx; 30*639f9d9dSBarry Smith Mat J = *B; 3191627157SBarry Smith Vec jj1,jj2,x2; 32*639f9d9dSBarry Smith int k, ierr,N,start,end,l,row,col,srow; 3391627157SBarry Smith Scalar dx, mone = -1.0,*y,*scale = color->scale,*xx,*wscale = color->wscale; 34*639f9d9dSBarry Smith double epsilon = color->error_rel, umin = color->umin; 3591627157SBarry Smith MPI_Comm comm = color->comm; 3691627157SBarry Smith 3791627157SBarry Smith ierr = MatZeroEntries(J); CHKERRQ(ierr); 3891627157SBarry Smith if (!snes->nvwork) { 3991627157SBarry Smith ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr); 4091627157SBarry Smith snes->nvwork = 3; 4191627157SBarry Smith PLogObjectParents(snes,3,snes->vwork); 4291627157SBarry Smith } 4391627157SBarry Smith jj1 = snes->vwork[0]; jj2 = snes->vwork[1]; x2 = snes->vwork[2]; 4491627157SBarry Smith 4591627157SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 4691627157SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 4791627157SBarry Smith ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 4891627157SBarry Smith ierr = SNESComputeFunction(snes,x1,jj1); CHKERRQ(ierr); 4991627157SBarry Smith 5091627157SBarry Smith PetscMemzero(wscale,N*sizeof(Scalar)); 5191627157SBarry Smith /* 5291627157SBarry Smith Loop over each color 5391627157SBarry Smith */ 54*639f9d9dSBarry Smith 5591627157SBarry Smith for (k=0; k<color->ncolors; k++) { 5691627157SBarry Smith ierr = VecCopy(x1,x2); CHKERRQ(ierr); 5791627157SBarry Smith /* 5891627157SBarry Smith Loop over each column associated with color adding the 5991627157SBarry Smith perturbation to the vector x2. 6091627157SBarry Smith */ 6191627157SBarry Smith for (l=0; l<color->ncolumns[k]; l++) { 6291627157SBarry Smith col = color->columns[k][l]; /* column of the matrix we are probing for */ 6391627157SBarry Smith dx = xx[col-start]; 6491627157SBarry Smith #if !defined(PETSC_COMPLEX) 65*639f9d9dSBarry Smith if (dx < umin && dx >= 0.0) dx = .1; 66*639f9d9dSBarry Smith else if (dx < 0.0 && dx > -umin) dx = -.1; 6791627157SBarry Smith #else 68*639f9d9dSBarry Smith if (abs(dx) < umin && real(dx) >= 0.0) dx = .1; 69*639f9d9dSBarry Smith else if (real(dx) < 0.0 && abs(dx) < umin) dx = -.1; 7091627157SBarry Smith #endif 7191627157SBarry Smith dx *= epsilon; 7291627157SBarry Smith wscale[col] = 1.0/dx; 7391627157SBarry Smith VecSetValues(x2,1,&col,&dx,ADD_VALUES); 7491627157SBarry Smith } 7591627157SBarry Smith VecRestoreArray(x1,&xx); 7691627157SBarry Smith /* 7791627157SBarry Smith Evaluate function at x1 + dx (here dx is a vector, of perturbations) 7891627157SBarry Smith */ 7991627157SBarry Smith ierr = SNESComputeFunction(snes,x2,jj2); CHKERRQ(ierr); 8091627157SBarry Smith ierr = VecAXPY(&mone,jj1,jj2); CHKERRQ(ierr); 8191627157SBarry Smith /* Communicate scale to all processors */ 8291627157SBarry Smith #if !defined(PETSC_COMPLEX) 8391627157SBarry Smith MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm); 8491627157SBarry Smith #else 8591627157SBarry Smith MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm); 8691627157SBarry Smith #endif 8791627157SBarry Smith /* 8891627157SBarry Smith Loop over rows of vector putting results into Jacobian matrix 8991627157SBarry Smith */ 9091627157SBarry Smith VecGetArray(jj2,&y); 9191627157SBarry Smith for (l=0; l<color->nrows[k]; l++) { 9291627157SBarry Smith row = color->rows[k][l]; 9391627157SBarry Smith col = color->columnsforrow[k][l]; 94*639f9d9dSBarry Smith y[row] *= scale[col]; 95*639f9d9dSBarry Smith srow = row + start; 96*639f9d9dSBarry Smith ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 9791627157SBarry Smith } 9891627157SBarry Smith VecRestoreArray(jj2,&y); 9991627157SBarry Smith } 10091627157SBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 10191627157SBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 10291627157SBarry Smith *flag = SAME_NONZERO_PATTERN; 10391627157SBarry Smith return 0; 10491627157SBarry Smith } 10591627157SBarry Smith 10691627157SBarry Smith 10791627157SBarry Smith 108