1 #define PETSCSNES_DLL 2 3 #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "SNESDefaultComputeJacobian" 7 /*@C 8 SNESDefaultComputeJacobian - Computes the Jacobian using finite differences. 9 10 Collective on SNES 11 12 Input Parameters: 13 + x1 - compute Jacobian at this point 14 - ctx - application's function context, as set with SNESSetFunction() 15 16 Output Parameters: 17 + J - Jacobian matrix (not altered in this routine) 18 . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 19 - flag - flag indicating whether the matrix sparsity structure has changed 20 21 Options Database Key: 22 + -snes_fd - Activates SNESDefaultComputeJacobian() 23 . -snes_test_err - Square root of function error tolerance, default square root of machine 24 epsilon (1.e-8 in double, 3.e-4 in single) 25 - -mat_fd_type - Either wp or ds (see MATSNESMF_WP or MATSNESMF_DS) 26 27 Notes: 28 This routine is slow and expensive, and is not currently optimized 29 to take advantage of sparsity in the problem. Although 30 SNESDefaultComputeJacobian() is not recommended for general use 31 in large-scale applications, It can be useful in checking the 32 correctness of a user-provided Jacobian. 33 34 An alternative routine that uses coloring to exploit matrix sparsity is 35 SNESDefaultComputeJacobianColor(). 36 37 Level: intermediate 38 39 .keywords: SNES, finite differences, Jacobian 40 41 .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor(), MatCreateSNESMF() 42 @*/ 43 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 44 { 45 Vec j1a,j2a,x2; 46 PetscErrorCode ierr; 47 PetscInt i,N,start,end,j,value; 48 PetscScalar dx,*y,scale,*xx,wscale; 49 PetscReal amax,epsilon = PETSC_SQRT_MACHINE_EPSILON; 50 PetscReal dx_min = 1.e-16,dx_par = 1.e-1,unorm; 51 MPI_Comm comm; 52 PetscErrorCode (*eval_fct)(SNES,Vec,Vec)=0; 53 PetscTruth assembled,use_wp = PETSC_TRUE,flg; 54 const char *list[2] = {"ds","wp"}; 55 56 PetscFunctionBegin; 57 ierr = PetscOptionsGetReal(snes->prefix,"-snes_test_err",&epsilon,0);CHKERRQ(ierr); 58 eval_fct = SNESComputeFunction; 59 60 ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr); 61 ierr = MatAssembled(*B,&assembled);CHKERRQ(ierr); 62 if (assembled) { 63 ierr = MatZeroEntries(*B);CHKERRQ(ierr); 64 } 65 if (!snes->nvwork) { 66 ierr = VecDuplicateVecs(x1,3,&snes->vwork);CHKERRQ(ierr); 67 snes->nvwork = 3; 68 ierr = PetscLogObjectParents(snes,3,snes->vwork);CHKERRQ(ierr); 69 } 70 j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 71 72 ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 73 ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 74 ierr = (*eval_fct)(snes,x1,j1a);CHKERRQ(ierr); 75 76 ierr = PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESDefaultComputeJacobian",list,2,"wp",&value,&flg);CHKERRQ(ierr); 77 if (flg && !value) { 78 use_wp = PETSC_FALSE; 79 } 80 if (use_wp) { 81 ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 82 } 83 /* Compute Jacobian approximation, 1 column at a time. 84 x1 = current iterate, j1a = F(x1) 85 x2 = perturbed iterate, j2a = F(x2) 86 */ 87 for (i=0; i<N; i++) { 88 ierr = VecCopy(x1,x2);CHKERRQ(ierr); 89 if (i>= start && i<end) { 90 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 91 if (use_wp) { 92 dx = 1.0 + unorm; 93 } else { 94 dx = xx[i-start]; 95 } 96 ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 97 #if !defined(PETSC_USE_COMPLEX) 98 if (dx < dx_min && dx >= 0.0) dx = dx_par; 99 else if (dx < 0.0 && dx > -dx_min) dx = -dx_par; 100 #else 101 if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par; 102 else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par; 103 #endif 104 dx *= epsilon; 105 wscale = 1.0/dx; 106 ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr); 107 } else { 108 wscale = 0.0; 109 } 110 ierr = (*eval_fct)(snes,x2,j2a);CHKERRQ(ierr); 111 ierr = VecAXPY(j2a,-1.0,j1a);CHKERRQ(ierr); 112 /* Communicate scale to all processors */ 113 ierr = MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr); 114 ierr = VecScale(j2a,scale);CHKERRQ(ierr); 115 ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14; 116 ierr = VecGetArray(j2a,&y);CHKERRQ(ierr); 117 for (j=start; j<end; j++) { 118 if (PetscAbsScalar(y[j-start]) > amax) { 119 ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr); 120 } 121 } 122 ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr); 123 } 124 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 125 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126 if (*B != *J) { 127 ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 128 ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 129 } 130 *flag = DIFFERENT_NONZERO_PATTERN; 131 PetscFunctionReturn(0); 132 } 133 134 135