1 #include <petsctao.h> /*I "petsctao.h" I*/ 2 #include <petsc/private/taoimpl.h> 3 #include <petscsnes.h> 4 5 /* 6 For finited difference computations of the Hessian, we use PETSc's SNESComputeJacobianDefault 7 */ 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "Fsnes" 11 static PetscErrorCode Fsnes(SNES snes ,Vec X,Vec G,void*ctx) 12 { 13 PetscErrorCode ierr; 14 Tao tao = (Tao)ctx; 15 16 PetscFunctionBegin; 17 PetscValidHeaderSpecific(ctx,TAO_CLASSID,4); 18 ierr=TaoComputeGradient(tao,X,G);CHKERRQ(ierr); 19 PetscFunctionReturn(0); 20 } 21 22 #undef __FUNCT__ 23 #define __FUNCT__ "TaoDefaultComputeGradient" 24 /*@C 25 TaoDefaultComputeGradient - computes the gradient using finite differences. 26 27 Collective on Tao 28 29 Input Parameters: 30 + tao - the Tao context 31 . X - compute gradient at this point 32 - dummy - not used 33 34 Output Parameters: 35 . G - Gradient Vector 36 37 Options Database Key: 38 + -tao_fd_gradient - Activates TaoDefaultComputeGradient() 39 - -tao_fd_delta <delta> - change in x used to calculate finite differences 40 41 Level: advanced 42 43 Note: 44 This routine is slow and expensive, and is not currently optimized 45 to take advantage of sparsity in the problem. Although 46 TaoAppDefaultComputeGradient is not recommended for general use 47 in large-scale applications, It can be useful in checking the 48 correctness of a user-provided gradient. Use the tao method TAOTEST 49 to get an indication of whether your gradient is correct. 50 51 52 Note: 53 This finite difference gradient evaluation can be set using the routine TaoSetGradientRoutine() or by using the command line option -tao_fd_gradient 54 55 .seealso: TaoSetGradientRoutine() 56 57 @*/ 58 PetscErrorCode TaoDefaultComputeGradient(Tao tao,Vec X,Vec G,void *dummy) 59 { 60 PetscScalar *x,*g; 61 PetscReal f, f2; 62 PetscErrorCode ierr; 63 PetscInt low,high,N,i; 64 PetscBool flg; 65 PetscReal h=.5*PETSC_SQRT_MACHINE_EPSILON; 66 67 PetscFunctionBegin; 68 ierr = PetscOptionsGetReal(((PetscObject)tao)->options,((PetscObject)tao)->prefix,"-tao_fd_delta",&h,&flg);CHKERRQ(ierr); 69 ierr = VecGetSize(X,&N);CHKERRQ(ierr); 70 ierr = VecGetOwnershipRange(X,&low,&high);CHKERRQ(ierr); 71 ierr = VecGetArray(G,&g);CHKERRQ(ierr); 72 for (i=0;i<N;i++) { 73 if (i>=low && i<high) { 74 ierr = VecGetArray(X,&x);CHKERRQ(ierr); 75 x[i-low] -= h; 76 ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 77 } 78 79 ierr = TaoComputeObjective(tao, X,&f);CHKERRQ(ierr); 80 81 if (i>=low && i<high) { 82 ierr = VecGetArray(X,&x);CHKERRQ(ierr); 83 x[i-low] += 2*h; 84 ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 85 } 86 87 ierr = TaoComputeObjective(tao,X,&f2);CHKERRQ(ierr); 88 89 if (i>=low && i<high) { 90 ierr = VecGetArray(X,&x);CHKERRQ(ierr); 91 x[i-low] -= h; 92 ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 93 } 94 if (i>=low && i<high) { 95 g[i-low]=(f2-f)/(2.0*h); 96 } 97 } 98 ierr = VecRestoreArray(G,&g);CHKERRQ(ierr); 99 PetscFunctionReturn(0); 100 } 101 102 #undef __FUNCT__ 103 #define __FUNCT__ "TaoDefaultComputeHessian" 104 /*@C 105 TaoDefaultComputeHessian - Computes the Hessian using finite differences. 106 107 Collective on Tao 108 109 Input Parameters: 110 + tao - the Tao context 111 . V - compute Hessian at this point 112 - dummy - not used 113 114 Output Parameters: 115 + H - Hessian matrix (not altered in this routine) 116 - B - newly computed Hessian matrix to use with preconditioner (generally the same as H) 117 118 Options Database Key: 119 + -tao_fd - Activates TaoDefaultComputeHessian() 120 - -tao_view_hessian - view the hessian after each evaluation using PETSC_VIEWER_STDOUT_WORLD 121 122 Level: advanced 123 124 Notes: 125 This routine is slow and expensive, and is not currently optimized 126 to take advantage of sparsity in the problem. Although 127 TaoDefaultComputeHessian() is not recommended for general use 128 in large-scale applications, It can be useful in checking the 129 correctness of a user-provided Hessian. 130 131 .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessianColor(), SNESComputeJacobianDefault(), TaoSetGradientRoutine(), TaoDefaultComputeGradient() 132 133 @*/ 134 PetscErrorCode TaoDefaultComputeHessian(Tao tao,Vec V,Mat H,Mat B,void *dummy) 135 { 136 PetscErrorCode ierr; 137 MPI_Comm comm; 138 Vec G; 139 SNES snes; 140 141 PetscFunctionBegin; 142 PetscValidHeaderSpecific(V,VEC_CLASSID,2); 143 ierr = VecDuplicate(V,&G);CHKERRQ(ierr); 144 145 ierr = PetscInfo(tao,"TAO Using finite differences w/o coloring to compute Hessian matrix\n");CHKERRQ(ierr); 146 147 ierr = TaoComputeGradient(tao,V,G);CHKERRQ(ierr); 148 149 ierr = PetscObjectGetComm((PetscObject)H,&comm);CHKERRQ(ierr); 150 ierr = SNESCreate(comm,&snes);CHKERRQ(ierr); 151 152 ierr = SNESSetFunction(snes,G,Fsnes,tao);CHKERRQ(ierr); 153 ierr = SNESComputeJacobianDefault(snes,V,H,B,tao);CHKERRQ(ierr); 154 ierr = SNESDestroy(&snes);CHKERRQ(ierr); 155 ierr = VecDestroy(&G);CHKERRQ(ierr); 156 PetscFunctionReturn(0); 157 } 158 159 #undef __FUNCT__ 160 #define __FUNCT__ "TaoDefaultComputeHessianColor" 161 /*@C 162 TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences. 163 164 Collective on Tao 165 166 Input Parameters: 167 + tao - the Tao context 168 . V - compute Hessian at this point 169 - ctx - the PetscColoring object (must be of type MatFDColoring) 170 171 Output Parameters: 172 + H - Hessian matrix (not altered in this routine) 173 - B - newly computed Hessian matrix to use with preconditioner (generally the same as H) 174 175 Level: advanced 176 177 178 .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessian(),SNESComputeJacobianDefaultColor(), TaoSetGradientRoutine() 179 180 @*/ 181 PetscErrorCode TaoDefaultComputeHessianColor(Tao tao, Vec V, Mat H,Mat B,void *ctx) 182 { 183 PetscErrorCode ierr; 184 MatFDColoring coloring = (MatFDColoring)ctx; 185 186 PetscFunctionBegin; 187 PetscValidHeaderSpecific(ctx,MAT_FDCOLORING_CLASSID,6); 188 ierr=PetscInfo(tao,"TAO computing matrix using finite differences Hessian and coloring\n");CHKERRQ(ierr); 189 ierr = MatFDColoringApply(B,coloring,V,ctx);CHKERRQ(ierr); 190 if (H != B) { 191 ierr = MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 192 ierr = MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 193 } 194 PetscFunctionReturn(0); 195 } 196 197 198