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