1 #define PETSCSNES_DLL 2 3 #include "private/matimpl.h" /*I "petscmat.h" I*/ 4 #include "private/snesimpl.h" /*I "petscsnes.h" I*/ 5 6 #undef __FUNCT__ 7 #define __FUNCT__ "SNESDefaultComputeJacobianColor" 8 /*@C 9 SNESDefaultComputeJacobianColor - Computes the Jacobian using 10 finite differences and coloring to exploit matrix sparsity. 11 12 Collective on SNES 13 14 Input Parameters: 15 + snes - nonlinear solver object 16 . x1 - location at which to evaluate Jacobian 17 - ctx - coloring context, where ctx must have type MatFDColoring, 18 as created via MatFDColoringCreate() 19 20 Output Parameters: 21 + J - Jacobian matrix (not altered in this routine) 22 . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 23 - flag - flag indicating whether the matrix sparsity structure has changed 24 25 Level: intermediate 26 27 .keywords: SNES, finite differences, Jacobian, coloring, sparse 28 29 .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESDefaultComputeJacobian() 30 TSDefaultComputeJacobianColor(), MatFDColoringCreate(), 31 MatFDColoringSetFunction() 32 33 @*/ 34 35 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultComputeJacobianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 36 { 37 MatFDColoring color = (MatFDColoring) ctx; 38 PetscErrorCode ierr; 39 Vec f; 40 PetscErrorCode (*ff)(void),(*fd)(void); 41 42 PetscFunctionBegin; 43 *flag = SAME_NONZERO_PATTERN; 44 ierr = SNESGetFunction(snes,&f,(PetscErrorCode (**)(SNES,Vec,Vec,void*))&ff,0);CHKERRQ(ierr); 45 ierr = MatFDColoringGetFunction(color,&fd,PETSC_NULL);CHKERRQ(ierr); 46 if (fd == ff) { /* reuse function value computed in SNES */ 47 ierr = MatFDColoringSetF(color,f);CHKERRQ(ierr); 48 } 49 ierr = MatFDColoringApply(*B,color,x1,flag,snes);CHKERRQ(ierr); 50 if (*J != *B) { 51 ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 52 ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 53 } 54 PetscFunctionReturn(0); 55 } 56 57 58 59