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