1 #define PETSCSNES_DLL 2 3 #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 4 #include "src/snes/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 Options Database Keys: 26 . -mat_fd_coloring_freq <freq> - Activates SNESDefaultComputeJacobianColor() 27 28 Level: intermediate 29 30 .keywords: SNES, finite differences, Jacobian, coloring, sparse 31 32 .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESDefaultComputeJacobian() 33 TSDefaultComputeJacobianColor(), MatFDColoringCreate(), 34 MatFDColoringSetFunction() 35 36 @*/ 37 38 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultComputeJacobianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 39 { 40 MatFDColoring color = (MatFDColoring) ctx; 41 PetscErrorCode ierr; 42 PetscInt freq,it; 43 Vec f; 44 PetscErrorCode (*ff)(void),(*fd)(void); 45 46 PetscFunctionBegin; 47 ierr = MatFDColoringGetFrequency(color,&freq);CHKERRQ(ierr); 48 ierr = SNESGetIterationNumber(snes,&it);CHKERRQ(ierr); 49 50 if ((freq > 1) && ((it % freq))) { 51 ierr = PetscInfo2(color,"Skipping Jacobian recomputation, it %D, freq %D\n",it,freq);CHKERRQ(ierr); 52 *flag = SAME_PRECONDITIONER; 53 } else { 54 ierr = PetscInfo2(color,"Computing Jacobian, it %D, freq %D\n",it,freq);CHKERRQ(ierr); 55 *flag = SAME_NONZERO_PATTERN; 56 ierr = SNESGetFunction(snes,&f,(PetscErrorCode (**)(SNES,Vec,Vec,void*))&ff,0);CHKERRQ(ierr); 57 ierr = MatFDColoringGetFunction(color,&fd,PETSC_NULL);CHKERRQ(ierr); 58 if (fd == ff) { /* reuse function value computed in SNES */ 59 ierr = MatFDColoringSetF(color,f);CHKERRQ(ierr); 60 } 61 ierr = MatFDColoringApply(*B,color,x1,flag,snes);CHKERRQ(ierr); 62 } 63 if (*J != *B) { 64 ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 65 ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 66 } 67 PetscFunctionReturn(0); 68 } 69 70 71 72