1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: snesj.c,v 1.54 1998/12/09 16:06:57 balay Exp curfman $"; 3 #endif 4 5 #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 6 7 #undef __FUNC__ 8 #define __FUNC__ "SNESDefaultComputeJacobian" 9 /*@C 10 SNESDefaultComputeJacobian - Computes the Jacobian using finite differences. 11 12 Collective on SNES 13 14 Input Parameters: 15 + x1 - compute Jacobian at this point 16 - ctx - application's function context, as set with SNESSetFunction() 17 18 Output Parameters: 19 + J - Jacobian matrix (not altered in this routine) 20 . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 21 - flag - flag indicating whether the matrix sparsity structure has changed 22 23 Options Database Key: 24 . -snes_fd - Activates SNESDefaultComputeJacobian() 25 26 Notes: 27 This routine is slow and expensive, and is not currently optimized 28 to take advantage of sparsity in the problem. Although 29 SNESDefaultComputeJacobian() is not recommended for general use 30 in large-scale applications, It can be useful in checking the 31 correctness of a user-provided Jacobian. 32 33 An alternative routine that uses coloring to explot matrix sparsity is 34 SNESDefaultComputeJacobianWithColoring(). 35 36 Level: intermediate 37 38 .keywords: SNES, finite differences, Jacobian 39 40 .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianWithColoring() 41 @*/ 42 int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B, 43 MatStructure *flag,void *ctx) 44 { 45 Vec j1a,j2a,x2; 46 int i,ierr,N,start,end,j; 47 Scalar dx, mone = -1.0,*y,scale,*xx,wscale; 48 double amax, epsilon = 1.e-8; /* assumes double precision */ 49 double dx_min = 1.e-16, dx_par = 1.e-1; 50 MPI_Comm comm; 51 int (*eval_fct)(SNES,Vec,Vec)=0; 52 53 PetscFunctionBegin; 54 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction; 55 else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient; 56 else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class"); 57 58 PetscObjectGetComm((PetscObject)x1,&comm); 59 MatZeroEntries(*B); 60 if (!snes->nvwork) { 61 ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr); 62 snes->nvwork = 3; 63 PLogObjectParents(snes,3,snes->vwork); 64 } 65 j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 66 67 ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 68 ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 69 ierr = eval_fct(snes,x1,j1a); CHKERRQ(ierr); 70 71 /* Compute Jacobian approximation, 1 column at a time. 72 x1 = current iterate, j1a = F(x1) 73 x2 = perturbed iterate, j2a = F(x2) 74 */ 75 for ( i=0; i<N; i++ ) { 76 ierr = VecCopy(x1,x2); CHKERRQ(ierr); 77 if ( i>= start && i<end) { 78 ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 79 dx = xx[i-start]; 80 ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 81 #if !defined(USE_PETSC_COMPLEX) 82 if (dx < dx_min && dx >= 0.0) dx = dx_par; 83 else if (dx < 0.0 && dx > -dx_min) dx = -dx_par; 84 #else 85 if (PetscAbsScalar(dx) < dx_min && PetscReal(dx) >= 0.0) dx = dx_par; 86 else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par; 87 #endif 88 dx *= epsilon; 89 wscale = 1.0/dx; 90 ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES); CHKERRQ(ierr); 91 } else { 92 wscale = 0.0; 93 } 94 ierr = eval_fct(snes,x2,j2a); CHKERRQ(ierr); 95 ierr = VecAXPY(&mone,j1a,j2a); CHKERRQ(ierr); 96 /* Communicate scale to all processors */ 97 #if !defined(USE_PETSC_COMPLEX) 98 ierr = MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 99 #else 100 ierr = MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 101 #endif 102 ierr = VecScale(&scale,j2a);CHKERRQ(ierr); 103 ierr = VecGetArray(j2a,&y);CHKERRQ(ierr); 104 ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14; 105 for ( j=start; j<end; j++ ) { 106 if (PetscAbsScalar(y[j-start]) > amax) { 107 ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES); CHKERRQ(ierr); 108 } 109 } 110 ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr); 111 } 112 ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 113 ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 114 *flag = DIFFERENT_NONZERO_PATTERN; 115 PetscFunctionReturn(0); 116 } 117 118 #undef __FUNC__ 119 #define __FUNC__ "SNESDefaultComputeHessian" 120 /*@C 121 SNESDefaultComputeHessian - Computes the Hessian using finite differences. 122 123 Collective on SNES 124 125 Input Parameters: 126 + x1 - compute Hessian at this point 127 - ctx - application's gradient context, as set with SNESSetGradient() 128 129 Output Parameters: 130 + J - Hessian matrix (not altered in this routine) 131 . B - newly computed Hessian matrix to use with preconditioner (generally the same as J) 132 - flag - flag indicating whether the matrix sparsity structure has changed 133 134 Options Database Key: 135 $ -snes_fd - Activates SNESDefaultComputeHessian() 136 137 Notes: 138 This routine is slow and expensive, and is not currently optimized 139 to take advantage of sparsity in the problem. Although 140 SNESDefaultComputeHessian() is not recommended for general use 141 in large-scale applications, It can be useful in checking the 142 correctness of a user-provided Hessian. 143 144 .keywords: SNES, finite differences, Hessian 145 146 .seealso: SNESSetHessian() 147 @*/ 148 int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B, 149 MatStructure *flag,void *ctx) 150 { 151 int ierr; 152 153 PetscFunctionBegin; 154 ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr); 155 PetscFunctionReturn(0); 156 } 157