1ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/taoimpl.h> 3aaa7dc30SBarry Smith #include <petscsnes.h> 4*f4c1ad5cSStefano Zampini #include <petscdmshell.h> 5a7e14dcfSSatish Balay 6a7e14dcfSSatish Balay /* 7a7e14dcfSSatish Balay For finited difference computations of the Hessian, we use PETSc's SNESComputeJacobianDefault 8a7e14dcfSSatish Balay */ 95ca45b2bSBarry Smith static PetscErrorCode Fsnes(SNES snes,Vec X,Vec G,void* ctx) 105ca45b2bSBarry Smith { 11a7e14dcfSSatish Balay PetscErrorCode ierr; 12441846f8SBarry Smith Tao tao = (Tao)ctx; 135ca45b2bSBarry Smith 14a7e14dcfSSatish Balay PetscFunctionBegin; 15*f4c1ad5cSStefano Zampini PetscValidHeaderSpecific(tao,TAO_CLASSID,4); 16a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao,X,G);CHKERRQ(ierr); 17a7e14dcfSSatish Balay PetscFunctionReturn(0); 18a7e14dcfSSatish Balay } 19a7e14dcfSSatish Balay 20a7e14dcfSSatish Balay /*@C 21a7e14dcfSSatish Balay TaoDefaultComputeGradient - computes the gradient using finite differences. 22a7e14dcfSSatish Balay 23441846f8SBarry Smith Collective on Tao 24a7e14dcfSSatish Balay 25a7e14dcfSSatish Balay Input Parameters: 26441846f8SBarry Smith + tao - the Tao context 27a7e14dcfSSatish Balay . X - compute gradient at this point 28a7e14dcfSSatish Balay - dummy - not used 29a7e14dcfSSatish Balay 30a7e14dcfSSatish Balay Output Parameters: 31a7e14dcfSSatish Balay . G - Gradient Vector 32a7e14dcfSSatish Balay 33a7e14dcfSSatish Balay Options Database Key: 34*f4c1ad5cSStefano Zampini + -tao_fd_gradient - activates TaoDefaultComputeGradient() 35*f4c1ad5cSStefano Zampini - -tao_fd_delta <delta> - change in X used to calculate finite differences 36a7e14dcfSSatish Balay 37a7e14dcfSSatish Balay Level: advanced 38a7e14dcfSSatish Balay 39*f4c1ad5cSStefano Zampini Notes: 40a7e14dcfSSatish Balay This routine is slow and expensive, and is not currently optimized 41a7e14dcfSSatish Balay to take advantage of sparsity in the problem. Although 42*f4c1ad5cSStefano Zampini TaoDefaultComputeGradient is not recommended for general use 43a7e14dcfSSatish Balay in large-scale applications, It can be useful in checking the 4458417fe7SBarry Smith correctness of a user-provided gradient. Use the tao method TAOTEST 45a7e14dcfSSatish Balay to get an indication of whether your gradient is correct. 46a7e14dcfSSatish Balay This finite difference gradient evaluation can be set using the routine TaoSetGradientRoutine() or by using the command line option -tao_fd_gradient 47a7e14dcfSSatish Balay 48a7e14dcfSSatish Balay .seealso: TaoSetGradientRoutine() 49a7e14dcfSSatish Balay @*/ 50*f4c1ad5cSStefano Zampini PetscErrorCode TaoDefaultComputeGradient(Tao tao,Vec Xin,Vec G,void *dummy) 51a7e14dcfSSatish Balay { 52*f4c1ad5cSStefano Zampini Vec X; 53*f4c1ad5cSStefano Zampini PetscScalar *g; 54a7e14dcfSSatish Balay PetscReal f, f2; 55a7e14dcfSSatish Balay PetscErrorCode ierr; 56a7e14dcfSSatish Balay PetscInt low,high,N,i; 57a7e14dcfSSatish Balay PetscBool flg; 58d66999acSBarry Smith PetscReal h=.5*PETSC_SQRT_MACHINE_EPSILON; 595ca45b2bSBarry Smith 60a7e14dcfSSatish Balay PetscFunctionBegin; 61c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)tao)->options,((PetscObject)tao)->prefix,"-tao_fd_delta",&h,&flg);CHKERRQ(ierr); 62*f4c1ad5cSStefano Zampini ierr = VecDuplicate(Xin,&X);CHKERRQ(ierr); 63*f4c1ad5cSStefano Zampini ierr = VecCopy(Xin,X);CHKERRQ(ierr); 64a7e14dcfSSatish Balay ierr = VecGetSize(X,&N);CHKERRQ(ierr); 65a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(X,&low,&high);CHKERRQ(ierr); 66*f4c1ad5cSStefano Zampini ierr = VecSetOption(X,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 67a7e14dcfSSatish Balay ierr = VecGetArray(G,&g);CHKERRQ(ierr); 68a7e14dcfSSatish Balay for (i=0;i<N;i++) { 69*f4c1ad5cSStefano Zampini ierr = VecSetValue(X,i,-h,ADD_VALUES);CHKERRQ(ierr); 70*f4c1ad5cSStefano Zampini ierr = VecAssemblyBegin(X);CHKERRQ(ierr); 71*f4c1ad5cSStefano Zampini ierr = VecAssemblyEnd(X);CHKERRQ(ierr); 72d66999acSBarry Smith ierr = TaoComputeObjective(tao,X,&f);CHKERRQ(ierr); 73*f4c1ad5cSStefano Zampini ierr = VecSetValue(X,i,2.0*h,ADD_VALUES);CHKERRQ(ierr); 74*f4c1ad5cSStefano Zampini ierr = VecAssemblyBegin(X);CHKERRQ(ierr); 75*f4c1ad5cSStefano Zampini ierr = VecAssemblyEnd(X);CHKERRQ(ierr); 76a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao,X,&f2);CHKERRQ(ierr); 77*f4c1ad5cSStefano Zampini ierr = VecSetValue(X,i,-h,ADD_VALUES);CHKERRQ(ierr); 78*f4c1ad5cSStefano Zampini ierr = VecAssemblyBegin(X);CHKERRQ(ierr); 79*f4c1ad5cSStefano Zampini ierr = VecAssemblyEnd(X);CHKERRQ(ierr); 80a7e14dcfSSatish Balay if (i>=low && i<high) { 81d66999acSBarry Smith g[i-low]=(f2-f)/(2.0*h); 82a7e14dcfSSatish Balay } 83a7e14dcfSSatish Balay } 84a7e14dcfSSatish Balay ierr = VecRestoreArray(G,&g);CHKERRQ(ierr); 85*f4c1ad5cSStefano Zampini ierr = VecDestroy(&X);CHKERRQ(ierr); 86a7e14dcfSSatish Balay PetscFunctionReturn(0); 87a7e14dcfSSatish Balay } 88a7e14dcfSSatish Balay 89a7e14dcfSSatish Balay /*@C 90a7e14dcfSSatish Balay TaoDefaultComputeHessian - Computes the Hessian using finite differences. 91a7e14dcfSSatish Balay 92441846f8SBarry Smith Collective on Tao 93a7e14dcfSSatish Balay 94a7e14dcfSSatish Balay Input Parameters: 95441846f8SBarry Smith + tao - the Tao context 96a7e14dcfSSatish Balay . V - compute Hessian at this point 97a7e14dcfSSatish Balay - dummy - not used 98a7e14dcfSSatish Balay 99a7e14dcfSSatish Balay Output Parameters: 100a7e14dcfSSatish Balay + H - Hessian matrix (not altered in this routine) 10156744491SBarry Smith - B - newly computed Hessian matrix to use with preconditioner (generally the same as H) 102a7e14dcfSSatish Balay 103a7e14dcfSSatish Balay Options Database Key: 104*f4c1ad5cSStefano Zampini . -tao_fd_hessian - activates TaoDefaultComputeHessian() 105a7e14dcfSSatish Balay 106a7e14dcfSSatish Balay Level: advanced 107a7e14dcfSSatish Balay 108a7e14dcfSSatish Balay Notes: 109a7e14dcfSSatish Balay This routine is slow and expensive, and is not currently optimized 110a7e14dcfSSatish Balay to take advantage of sparsity in the problem. Although 111a7e14dcfSSatish Balay TaoDefaultComputeHessian() is not recommended for general use 112a7e14dcfSSatish Balay in large-scale applications, It can be useful in checking the 113a7e14dcfSSatish Balay correctness of a user-provided Hessian. 114a7e14dcfSSatish Balay 115a7e14dcfSSatish Balay .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessianColor(), SNESComputeJacobianDefault(), TaoSetGradientRoutine(), TaoDefaultComputeGradient() 116a7e14dcfSSatish Balay @*/ 117ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessian(Tao tao,Vec V,Mat H,Mat B,void *dummy) 118d1e9a80fSBarry Smith { 119a7e14dcfSSatish Balay PetscErrorCode ierr; 120a7e14dcfSSatish Balay Vec G; 121a7e14dcfSSatish Balay SNES snes; 122*f4c1ad5cSStefano Zampini DM dm; 123a7e14dcfSSatish Balay 124a7e14dcfSSatish Balay PetscFunctionBegin; 125a7e14dcfSSatish Balay ierr = VecDuplicate(V,&G);CHKERRQ(ierr); 126a7e14dcfSSatish Balay ierr = PetscInfo(tao,"TAO Using finite differences w/o coloring to compute Hessian matrix\n");CHKERRQ(ierr); 127a7e14dcfSSatish Balay ierr = TaoComputeGradient(tao,V,G);CHKERRQ(ierr); 128*f4c1ad5cSStefano Zampini ierr = SNESCreate(PetscObjectComm((PetscObject)H),&snes);CHKERRQ(ierr); 129a7e14dcfSSatish Balay ierr = SNESSetFunction(snes,G,Fsnes,tao);CHKERRQ(ierr); 130*f4c1ad5cSStefano Zampini ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 131*f4c1ad5cSStefano Zampini ierr = DMShellSetGlobalVector(dm,V);CHKERRQ(ierr); 132*f4c1ad5cSStefano Zampini ierr = SNESSetUp(snes);CHKERRQ(ierr); 133*f4c1ad5cSStefano Zampini if (H) { 134*f4c1ad5cSStefano Zampini PetscInt n,N; 135*f4c1ad5cSStefano Zampini 136*f4c1ad5cSStefano Zampini ierr = VecGetSize(V,&N);CHKERRQ(ierr); 137*f4c1ad5cSStefano Zampini ierr = VecGetLocalSize(V,&n);CHKERRQ(ierr); 138*f4c1ad5cSStefano Zampini ierr = MatSetSizes(H,n,n,N,N);CHKERRQ(ierr); 139*f4c1ad5cSStefano Zampini ierr = MatSetUp(H);CHKERRQ(ierr); 140*f4c1ad5cSStefano Zampini } 141*f4c1ad5cSStefano Zampini if (B && B != H) { 142*f4c1ad5cSStefano Zampini PetscInt n,N; 143*f4c1ad5cSStefano Zampini 144*f4c1ad5cSStefano Zampini ierr = VecGetSize(V,&N);CHKERRQ(ierr); 145*f4c1ad5cSStefano Zampini ierr = VecGetLocalSize(V,&n);CHKERRQ(ierr); 146*f4c1ad5cSStefano Zampini ierr = MatSetSizes(B,n,n,N,N);CHKERRQ(ierr); 147*f4c1ad5cSStefano Zampini ierr = MatSetUp(B);CHKERRQ(ierr); 148*f4c1ad5cSStefano Zampini } 149*f4c1ad5cSStefano Zampini ierr = SNESComputeJacobianDefault(snes,V,H,B,NULL);CHKERRQ(ierr); 150a7e14dcfSSatish Balay ierr = SNESDestroy(&snes);CHKERRQ(ierr); 151a7e14dcfSSatish Balay ierr = VecDestroy(&G);CHKERRQ(ierr); 152a7e14dcfSSatish Balay PetscFunctionReturn(0); 153a7e14dcfSSatish Balay } 154a7e14dcfSSatish Balay 155a7e14dcfSSatish Balay /*@C 156a7e14dcfSSatish Balay TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences. 157a7e14dcfSSatish Balay 158441846f8SBarry Smith Collective on Tao 159a7e14dcfSSatish Balay 160a7e14dcfSSatish Balay Input Parameters: 161441846f8SBarry Smith + tao - the Tao context 162a7e14dcfSSatish Balay . V - compute Hessian at this point 163a7e14dcfSSatish Balay - ctx - the PetscColoring object (must be of type MatFDColoring) 164a7e14dcfSSatish Balay 165a7e14dcfSSatish Balay Output Parameters: 166a7e14dcfSSatish Balay + H - Hessian matrix (not altered in this routine) 16756744491SBarry Smith - B - newly computed Hessian matrix to use with preconditioner (generally the same as H) 168a7e14dcfSSatish Balay 169a7e14dcfSSatish Balay Level: advanced 170a7e14dcfSSatish Balay 171a7e14dcfSSatish Balay 172a7e14dcfSSatish Balay .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessian(),SNESComputeJacobianDefaultColor(), TaoSetGradientRoutine() 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 190*f4c1ad5cSStefano Zampini PetscErrorCode TaoDefaultComputeHessianMFFD(Tao tao,Vec X,Mat H,Mat B,void *ctx) 191*f4c1ad5cSStefano Zampini { 192*f4c1ad5cSStefano Zampini PetscInt n,N; 193*f4c1ad5cSStefano Zampini PetscErrorCode ierr; 194a7e14dcfSSatish Balay 195*f4c1ad5cSStefano Zampini PetscFunctionBegin; 196*f4c1ad5cSStefano Zampini if (B && B != H) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Preconditioning Hessian matrix"); 197*f4c1ad5cSStefano Zampini ierr = VecGetSize(X,&N);CHKERRQ(ierr); 198*f4c1ad5cSStefano Zampini ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 199*f4c1ad5cSStefano Zampini ierr = MatSetSizes(H,n,n,N,N);CHKERRQ(ierr); 200*f4c1ad5cSStefano Zampini ierr = MatSetType(H,MATMFFD);CHKERRQ(ierr); 201*f4c1ad5cSStefano Zampini ierr = MatSetUp(H);CHKERRQ(ierr); 202*f4c1ad5cSStefano Zampini ierr = MatMFFDSetBase(H,X,NULL);CHKERRQ(ierr); 203*f4c1ad5cSStefano Zampini ierr = MatMFFDSetFunction(H,(PetscErrorCode (*)(void*,Vec,Vec))TaoComputeGradient,tao);CHKERRQ(ierr); 204*f4c1ad5cSStefano Zampini ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 205*f4c1ad5cSStefano Zampini ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 206*f4c1ad5cSStefano Zampini PetscFunctionReturn(0); 207*f4c1ad5cSStefano Zampini } 208