1 2 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 3 #include <petscdm.h> 4 5 /*@C 6 SNESComputeJacobianDefault - Computes the Jacobian using finite differences. 7 8 Collective on SNES 9 10 Input Parameters: 11 + snes - the SNES context 12 . x1 - compute Jacobian at this point 13 - ctx - application's function context, as set with SNESSetFunction() 14 15 Output Parameters: 16 + J - Jacobian matrix (not altered in this routine) 17 - B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 18 19 Options Database Key: 20 + -snes_fd - Activates SNESComputeJacobianDefault() 21 . -snes_test_err - Square root of function error tolerance, default square root of machine 22 epsilon (1.e-8 in double, 3.e-4 in single) 23 - -mat_fd_type - Either wp or ds (see MATMFFD_WP or MATMFFD_DS) 24 25 Notes: 26 This routine is slow and expensive, and is not currently optimized 27 to take advantage of sparsity in the problem. Although 28 SNESComputeJacobianDefault() is not recommended for general use 29 in large-scale applications, It can be useful in checking the 30 correctness of a user-provided Jacobian. 31 32 An alternative routine that uses coloring to exploit matrix sparsity is 33 SNESComputeJacobianDefaultColor(). 34 35 This routine ignores the maximum number of function evaluations set with SNESSetTolerances() and the function 36 evaluations it performs are not counted in what is returned by of SNESGetNumberFunctionEvals(). 37 38 Level: intermediate 39 40 .seealso: `SNESSetJacobian()`, `SNESComputeJacobianDefaultColor()`, `MatCreateSNESMF()` 41 @*/ 42 PetscErrorCode SNESComputeJacobianDefault(SNES snes,Vec x1,Mat J,Mat B,void *ctx) 43 { 44 Vec j1a,j2a,x2; 45 PetscInt i,N,start,end,j,value,root,max_funcs = snes->max_funcs; 46 PetscScalar dx,*y,wscale; 47 const PetscScalar *xx; 48 PetscReal amax,epsilon = PETSC_SQRT_MACHINE_EPSILON; 49 PetscReal dx_min = 1.e-16,dx_par = 1.e-1,unorm; 50 MPI_Comm comm; 51 PetscBool assembled,use_wp = PETSC_TRUE,flg; 52 const char *list[2] = {"ds","wp"}; 53 PetscMPIInt size; 54 const PetscInt *ranges; 55 DM dm; 56 DMSNES dms; 57 58 PetscFunctionBegin; 59 snes->max_funcs = PETSC_MAX_INT; 60 /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */ 61 PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 62 PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_err",&epsilon,NULL)); 63 64 PetscCall(PetscObjectGetComm((PetscObject)x1,&comm)); 65 PetscCallMPI(MPI_Comm_size(comm,&size)); 66 PetscCall(MatAssembled(B,&assembled)); 67 if (assembled) { 68 PetscCall(MatZeroEntries(B)); 69 } 70 if (!snes->nvwork) { 71 if (snes->dm) { 72 PetscCall(DMGetGlobalVector(snes->dm,&j1a)); 73 PetscCall(DMGetGlobalVector(snes->dm,&j2a)); 74 PetscCall(DMGetGlobalVector(snes->dm,&x2)); 75 } else { 76 snes->nvwork = 3; 77 PetscCall(VecDuplicateVecs(x1,snes->nvwork,&snes->vwork)); 78 PetscCall(PetscLogObjectParents(snes,snes->nvwork,snes->vwork)); 79 j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 80 } 81 } 82 83 PetscCall(VecGetSize(x1,&N)); 84 PetscCall(VecGetOwnershipRange(x1,&start,&end)); 85 PetscCall(SNESGetDM(snes,&dm)); 86 PetscCall(DMGetDMSNES(dm,&dms)); 87 if (dms->ops->computemffunction) { 88 PetscCall(SNESComputeMFFunction(snes,x1,j1a)); 89 } else { 90 PetscCall(SNESComputeFunction(snes,x1,j1a)); 91 } 92 93 PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing options","SNES"); 94 PetscCall(PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESComputeJacobianDefault",list,2,"wp",&value,&flg)); 95 PetscOptionsEnd(); 96 if (flg && !value) use_wp = PETSC_FALSE; 97 98 if (use_wp) { 99 PetscCall(VecNorm(x1,NORM_2,&unorm)); 100 } 101 /* Compute Jacobian approximation, 1 column at a time. 102 x1 = current iterate, j1a = F(x1) 103 x2 = perturbed iterate, j2a = F(x2) 104 */ 105 for (i=0; i<N; i++) { 106 PetscCall(VecCopy(x1,x2)); 107 if (i>= start && i<end) { 108 PetscCall(VecGetArrayRead(x1,&xx)); 109 if (use_wp) dx = PetscSqrtReal(1.0 + unorm); 110 else dx = xx[i-start]; 111 PetscCall(VecRestoreArrayRead(x1,&xx)); 112 if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par; 113 dx *= epsilon; 114 wscale = 1.0/dx; 115 PetscCall(VecSetValues(x2,1,&i,&dx,ADD_VALUES)); 116 } else { 117 wscale = 0.0; 118 } 119 PetscCall(VecAssemblyBegin(x2)); 120 PetscCall(VecAssemblyEnd(x2)); 121 if (dms->ops->computemffunction) { 122 PetscCall(SNESComputeMFFunction(snes,x2,j2a)); 123 } else { 124 PetscCall(SNESComputeFunction(snes,x2,j2a)); 125 } 126 PetscCall(VecAXPY(j2a,-1.0,j1a)); 127 /* Communicate scale=1/dx_i to all processors */ 128 PetscCall(VecGetOwnershipRanges(x1,&ranges)); 129 root = size; 130 for (j=size-1; j>-1; j--) { 131 root--; 132 if (i>=ranges[j]) break; 133 } 134 PetscCallMPI(MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm)); 135 PetscCall(VecScale(j2a,wscale)); 136 PetscCall(VecNorm(j2a,NORM_INFINITY,&amax)); amax *= 1.e-14; 137 PetscCall(VecGetArray(j2a,&y)); 138 for (j=start; j<end; j++) { 139 if (PetscAbsScalar(y[j-start]) > amax || j == i) { 140 PetscCall(MatSetValues(B,1,&j,1,&i,y+j-start,INSERT_VALUES)); 141 } 142 } 143 PetscCall(VecRestoreArray(j2a,&y)); 144 } 145 if (snes->dm) { 146 PetscCall(DMRestoreGlobalVector(snes->dm,&j1a)); 147 PetscCall(DMRestoreGlobalVector(snes->dm,&j2a)); 148 PetscCall(DMRestoreGlobalVector(snes->dm,&x2)); 149 } 150 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 151 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 152 if (B != J) { 153 PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 154 PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 155 } 156 snes->max_funcs = max_funcs; 157 snes->nfuncs -= N; 158 PetscFunctionReturn(0); 159 } 160