1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*3a40ed3dSBarry Smith static char vcid[] = "$Id: snesj2.c,v 1.11 1997/10/12 21:45:52 bsmith Exp bsmith $"; 391627157SBarry Smith #endif 491627157SBarry Smith 5639f9d9dSBarry Smith #include "src/mat/matimpl.h" /*I "mat.h" I*/ 6639f9d9dSBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 791627157SBarry Smith 85615d1e5SSatish Balay #undef __FUNC__ 95615d1e5SSatish Balay #define __FUNC__ "SNESDefaultComputeJacobianWithColoring" 1091627157SBarry Smith /*@C 11b4fc646aSLois Curfman McInnes SNESDefaultComputeJacobianWithColoring - Computes the Jacobian using 12b4fc646aSLois Curfman McInnes finite differences and coloring to exploit matrix sparsity. 1391627157SBarry Smith 1491627157SBarry Smith Input Parameters: 1591627157SBarry Smith . snes - nonlinear solver object 1691627157SBarry Smith . x1 - location at which to evaluate Jacobian 17ff0cfa39SBarry Smith . ctx - coloring context, where 18ff0cfa39SBarry Smith $ ctx must have type MatFDColoring, 19ff0cfa39SBarry Smith $ as created via MatFDColoringCreate() 2091627157SBarry Smith 2191627157SBarry Smith Output Parameters: 22b4fc646aSLois Curfman McInnes . J - Jacobian matrix (not altered in this routine) 23b4fc646aSLois Curfman McInnes . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 24b4fc646aSLois Curfman McInnes . flag - flag indicating whether the matrix sparsity structure has changed 2591627157SBarry Smith 26dff777c9SBarry Smith Options Database Keys: 27dff777c9SBarry Smith $ -mat_fd_coloring_freq <freq> 28dff777c9SBarry Smith 29b4fc646aSLois Curfman McInnes .keywords: SNES, finite differences, Jacobian, coloring, sparse 3091627157SBarry Smith 31b4fc646aSLois Curfman McInnes .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESDefaultComputeJacobian() 3291627157SBarry Smith @*/ 33b4fc646aSLois Curfman McInnes int SNESDefaultComputeJacobianWithColoring(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 3491627157SBarry Smith { 3591627157SBarry Smith MatFDColoring color = (MatFDColoring) ctx; 36dff777c9SBarry Smith int ierr,freq,it; 37dff777c9SBarry Smith 38*3a40ed3dSBarry Smith PetscFunctionBegin; 39dff777c9SBarry Smith ierr = MatFDColoringGetFrequency(color,&freq);CHKERRQ(ierr); 40dff777c9SBarry Smith ierr = SNESGetIterationNumber(snes,&it); CHKERRQ(ierr); 41dff777c9SBarry Smith 42dff777c9SBarry Smith if ((freq > 1) && ((it % freq) != 1)) { 43dff777c9SBarry Smith PLogInfo(color,"SNESDefaultComputeJacobianWithColoring:Skipping Jacobian, it %d, freq %d\n",it,freq); 44dff777c9SBarry Smith *flag = SAME_PRECONDITIONER; 45*3a40ed3dSBarry Smith PetscFunctionReturn(0); 46dff777c9SBarry Smith } else { 47dff777c9SBarry Smith PLogInfo(color,"SNESDefaultComputeJacobianWithColoring:Computing Jacobian, it %d, freq %d\n",it,freq); 48dff777c9SBarry Smith *flag = SAME_NONZERO_PATTERN; 49dff777c9SBarry Smith } 5088389382SBarry Smith 51005c665bSBarry Smith ierr = MatFDColoringApply(*B,color,x1,flag,snes); CHKERRQ(ierr); 52*3a40ed3dSBarry Smith PetscFunctionReturn(0); 5391627157SBarry Smith } 5491627157SBarry Smith 5591627157SBarry Smith 5691627157SBarry Smith 57