1ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/taoimpl.h> 3aaa7dc30SBarry Smith #include <petscsnes.h> 4f4c1ad5cSStefano 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 { 11441846f8SBarry Smith Tao tao = (Tao)ctx; 125ca45b2bSBarry Smith 13a7e14dcfSSatish Balay PetscFunctionBegin; 14f4c1ad5cSStefano Zampini PetscValidHeaderSpecific(tao,TAO_CLASSID,4); 159566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(tao,X,G)); 16a7e14dcfSSatish Balay PetscFunctionReturn(0); 17a7e14dcfSSatish Balay } 18a7e14dcfSSatish Balay 19a7e14dcfSSatish Balay /*@C 20a7e14dcfSSatish Balay TaoDefaultComputeGradient - computes the gradient using finite differences. 21a7e14dcfSSatish Balay 22*65ba42b6SBarry Smith Collective on tao 23a7e14dcfSSatish Balay 24a7e14dcfSSatish Balay Input Parameters: 25441846f8SBarry Smith + tao - the Tao context 26a7e14dcfSSatish Balay . X - compute gradient at this point 27a7e14dcfSSatish Balay - dummy - not used 28a7e14dcfSSatish Balay 29f899ff85SJose E. Roman Output Parameter: 30a7e14dcfSSatish Balay . G - Gradient Vector 31a7e14dcfSSatish Balay 32a7e14dcfSSatish Balay Options Database Key: 33f4c1ad5cSStefano Zampini + -tao_fd_gradient - activates TaoDefaultComputeGradient() 34f4c1ad5cSStefano Zampini - -tao_fd_delta <delta> - change in X used to calculate finite differences 35a7e14dcfSSatish Balay 36a7e14dcfSSatish Balay Level: advanced 37a7e14dcfSSatish Balay 38f4c1ad5cSStefano Zampini Notes: 39*65ba42b6SBarry Smith This routine is slow and expensive, and is not optimized 40a7e14dcfSSatish Balay to take advantage of sparsity in the problem. Although 41*65ba42b6SBarry Smith not recommended for general use 42*65ba42b6SBarry Smith in large-scale applications, it can be useful in checking the 4358417fe7SBarry Smith correctness of a user-provided gradient. Use the tao method TAOTEST 44a7e14dcfSSatish Balay to get an indication of whether your gradient is correct. 45*65ba42b6SBarry Smith This finite difference gradient evaluation can be set using the routine `TaoSetGradient()` or by using the command line option -tao_fd_gradient 46a7e14dcfSSatish Balay 47*65ba42b6SBarry Smith .seealso: `Tao`, `TaoSetGradient()` 48a7e14dcfSSatish Balay @*/ 49f4c1ad5cSStefano Zampini PetscErrorCode TaoDefaultComputeGradient(Tao tao,Vec Xin,Vec G,void *dummy) 50a7e14dcfSSatish Balay { 51f4c1ad5cSStefano Zampini Vec X; 52f4c1ad5cSStefano Zampini PetscScalar *g; 53a7e14dcfSSatish Balay PetscReal f, f2; 54a7e14dcfSSatish Balay PetscInt low,high,N,i; 55a7e14dcfSSatish Balay PetscBool flg; 56d66999acSBarry Smith PetscReal h=.5*PETSC_SQRT_MACHINE_EPSILON; 575ca45b2bSBarry Smith 58a7e14dcfSSatish Balay PetscFunctionBegin; 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)tao)->options,((PetscObject)tao)->prefix,"-tao_fd_delta",&h,&flg)); 609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Xin,&X)); 619566063dSJacob Faibussowitsch PetscCall(VecCopy(Xin,X)); 629566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&N)); 639566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X,&low,&high)); 649566063dSJacob Faibussowitsch PetscCall(VecSetOption(X,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE)); 659566063dSJacob Faibussowitsch PetscCall(VecGetArray(G,&g)); 66a7e14dcfSSatish Balay for (i=0;i<N;i++) { 679566063dSJacob Faibussowitsch PetscCall(VecSetValue(X,i,-h,ADD_VALUES)); 689566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 699566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 709566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao,X,&f)); 719566063dSJacob Faibussowitsch PetscCall(VecSetValue(X,i,2.0*h,ADD_VALUES)); 729566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 739566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 749566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao,X,&f2)); 759566063dSJacob Faibussowitsch PetscCall(VecSetValue(X,i,-h,ADD_VALUES)); 769566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 779566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 78a7e14dcfSSatish Balay if (i>=low && i<high) { 79d66999acSBarry Smith g[i-low]=(f2-f)/(2.0*h); 80a7e14dcfSSatish Balay } 81a7e14dcfSSatish Balay } 829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G,&g)); 839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 84a7e14dcfSSatish Balay PetscFunctionReturn(0); 85a7e14dcfSSatish Balay } 86a7e14dcfSSatish Balay 87a7e14dcfSSatish Balay /*@C 88a7e14dcfSSatish Balay TaoDefaultComputeHessian - Computes the Hessian using finite differences. 89a7e14dcfSSatish Balay 90*65ba42b6SBarry Smith Collective on tao 91a7e14dcfSSatish Balay 92a7e14dcfSSatish Balay Input Parameters: 93441846f8SBarry Smith + tao - the Tao context 94a7e14dcfSSatish Balay . V - compute Hessian at this point 95a7e14dcfSSatish Balay - dummy - not used 96a7e14dcfSSatish Balay 97a7e14dcfSSatish Balay Output Parameters: 98a7e14dcfSSatish Balay + H - Hessian matrix (not altered in this routine) 9956744491SBarry Smith - B - newly computed Hessian matrix to use with preconditioner (generally the same as H) 100a7e14dcfSSatish Balay 101a7e14dcfSSatish Balay Options Database Key: 102f4c1ad5cSStefano Zampini . -tao_fd_hessian - activates TaoDefaultComputeHessian() 103a7e14dcfSSatish Balay 104a7e14dcfSSatish Balay Level: advanced 105a7e14dcfSSatish Balay 106a7e14dcfSSatish Balay Notes: 107*65ba42b6SBarry Smith This routine is slow and expensive, and is not optimized 108a7e14dcfSSatish Balay to take advantage of sparsity in the problem. Although 109*65ba42b6SBarry Smith it is not recommended for general use 110a7e14dcfSSatish Balay in large-scale applications, It can be useful in checking the 111a7e14dcfSSatish Balay correctness of a user-provided Hessian. 112a7e14dcfSSatish Balay 113*65ba42b6SBarry Smith .seealso: `Tao`, `TaoSetHessian()`, `TaoDefaultComputeHessianColor()`, `SNESComputeJacobianDefault()`, `TaoSetGradient()`, `TaoDefaultComputeGradient()` 114a7e14dcfSSatish Balay @*/ 115ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessian(Tao tao,Vec V,Mat H,Mat B,void *dummy) 116d1e9a80fSBarry Smith { 117a7e14dcfSSatish Balay SNES snes; 118f4c1ad5cSStefano Zampini DM dm; 119a7e14dcfSSatish Balay 120a7e14dcfSSatish Balay PetscFunctionBegin; 1219566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"TAO Using finite differences w/o coloring to compute Hessian matrix\n")); 1229566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)H),&snes)); 1239566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,NULL,Fsnes,tao)); 1249566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 1259566063dSJacob Faibussowitsch PetscCall(DMShellSetGlobalVector(dm,V)); 1269566063dSJacob Faibussowitsch PetscCall(SNESSetUp(snes)); 127f4c1ad5cSStefano Zampini if (H) { 128f4c1ad5cSStefano Zampini PetscInt n,N; 129f4c1ad5cSStefano Zampini 1309566063dSJacob Faibussowitsch PetscCall(VecGetSize(V,&N)); 1319566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(V,&n)); 1329566063dSJacob Faibussowitsch PetscCall(MatSetSizes(H,n,n,N,N)); 1339566063dSJacob Faibussowitsch PetscCall(MatSetUp(H)); 134f4c1ad5cSStefano Zampini } 135f4c1ad5cSStefano Zampini if (B && B != H) { 136f4c1ad5cSStefano Zampini PetscInt n,N; 137f4c1ad5cSStefano Zampini 1389566063dSJacob Faibussowitsch PetscCall(VecGetSize(V,&N)); 1399566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(V,&n)); 1409566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,n,n,N,N)); 1419566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 142f4c1ad5cSStefano Zampini } 1439566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobianDefault(snes,V,H,B,NULL)); 1449566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 145a7e14dcfSSatish Balay PetscFunctionReturn(0); 146a7e14dcfSSatish Balay } 147a7e14dcfSSatish Balay 148a7e14dcfSSatish Balay /*@C 149a7e14dcfSSatish Balay TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences. 150a7e14dcfSSatish Balay 151441846f8SBarry Smith Collective on Tao 152a7e14dcfSSatish Balay 153a7e14dcfSSatish Balay Input Parameters: 154441846f8SBarry Smith + tao - the Tao context 155a7e14dcfSSatish Balay . V - compute Hessian at this point 156*65ba42b6SBarry Smith - ctx - the color object of type `MatFDColoring` 157a7e14dcfSSatish Balay 158a7e14dcfSSatish Balay Output Parameters: 159a7e14dcfSSatish Balay + H - Hessian matrix (not altered in this routine) 16056744491SBarry Smith - B - newly computed Hessian matrix to use with preconditioner (generally the same as H) 161a7e14dcfSSatish Balay 162a7e14dcfSSatish Balay Level: advanced 163a7e14dcfSSatish Balay 164*65ba42b6SBarry Smith .seealso: `Tao`, `MatColoring`, `TaoSetHessian()`, `TaoDefaultComputeHessian()`, `SNESComputeJacobianDefaultColor()`, `TaoSetGradient()` 165a7e14dcfSSatish Balay @*/ 166ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessianColor(Tao tao,Vec V,Mat H,Mat B,void *ctx) 1675ca45b2bSBarry Smith { 168a7e14dcfSSatish Balay MatFDColoring coloring = (MatFDColoring)ctx; 169a7e14dcfSSatish Balay 170a7e14dcfSSatish Balay PetscFunctionBegin; 171064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,5); 1729566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"TAO computing matrix using finite differences Hessian and coloring\n")); 1739566063dSJacob Faibussowitsch PetscCall(MatFDColoringApply(B,coloring,V,ctx)); 17494ab13aaSBarry Smith if (H != B) { 1759566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY)); 1769566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY)); 177a7e14dcfSSatish Balay } 178a7e14dcfSSatish Balay PetscFunctionReturn(0); 179a7e14dcfSSatish Balay } 180a7e14dcfSSatish Balay 181f4c1ad5cSStefano Zampini PetscErrorCode TaoDefaultComputeHessianMFFD(Tao tao,Vec X,Mat H,Mat B,void *ctx) 182f4c1ad5cSStefano Zampini { 183f4c1ad5cSStefano Zampini PetscInt n,N; 184c5e9d94cSAlp Dener PetscBool assembled; 185a7e14dcfSSatish Balay 186f4c1ad5cSStefano Zampini PetscFunctionBegin; 1873c859ba3SBarry Smith PetscCheck(!B || B == H,PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Preconditioning Hessian matrix"); 1889566063dSJacob Faibussowitsch PetscCall(MatAssembled(H, &assembled)); 189c5e9d94cSAlp Dener if (!assembled) { 1909566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&N)); 1919566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X,&n)); 1929566063dSJacob Faibussowitsch PetscCall(MatSetSizes(H,n,n,N,N)); 1939566063dSJacob Faibussowitsch PetscCall(MatSetType(H,MATMFFD)); 1949566063dSJacob Faibussowitsch PetscCall(MatSetUp(H)); 1959566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(H,(PetscErrorCode (*)(void*,Vec,Vec))TaoComputeGradient,tao)); 196c5e9d94cSAlp Dener } 1979566063dSJacob Faibussowitsch PetscCall(MatMFFDSetBase(H,X,NULL)); 1989566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY)); 1999566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY)); 200f4c1ad5cSStefano Zampini PetscFunctionReturn(0); 201f4c1ad5cSStefano Zampini } 202