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 { 11a7e14dcfSSatish Balay PetscErrorCode ierr; 12441846f8SBarry Smith Tao tao = (Tao)ctx; 135ca45b2bSBarry Smith 14a7e14dcfSSatish Balay PetscFunctionBegin; 15f4c1ad5cSStefano 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 30f899ff85SJose E. Roman Output Parameter: 31a7e14dcfSSatish Balay . G - Gradient Vector 32a7e14dcfSSatish Balay 33a7e14dcfSSatish Balay Options Database Key: 34f4c1ad5cSStefano Zampini + -tao_fd_gradient - activates TaoDefaultComputeGradient() 35f4c1ad5cSStefano Zampini - -tao_fd_delta <delta> - change in X used to calculate finite differences 36a7e14dcfSSatish Balay 37a7e14dcfSSatish Balay Level: advanced 38a7e14dcfSSatish Balay 39f4c1ad5cSStefano 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 42f4c1ad5cSStefano 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 @*/ 50f4c1ad5cSStefano Zampini PetscErrorCode TaoDefaultComputeGradient(Tao tao,Vec Xin,Vec G,void *dummy) 51a7e14dcfSSatish Balay { 52f4c1ad5cSStefano Zampini Vec X; 53f4c1ad5cSStefano 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); 62f4c1ad5cSStefano Zampini ierr = VecDuplicate(Xin,&X);CHKERRQ(ierr); 63f4c1ad5cSStefano Zampini ierr = VecCopy(Xin,X);CHKERRQ(ierr); 64a7e14dcfSSatish Balay ierr = VecGetSize(X,&N);CHKERRQ(ierr); 65a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(X,&low,&high);CHKERRQ(ierr); 66f4c1ad5cSStefano 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++) { 69f4c1ad5cSStefano Zampini ierr = VecSetValue(X,i,-h,ADD_VALUES);CHKERRQ(ierr); 70f4c1ad5cSStefano Zampini ierr = VecAssemblyBegin(X);CHKERRQ(ierr); 71f4c1ad5cSStefano Zampini ierr = VecAssemblyEnd(X);CHKERRQ(ierr); 72d66999acSBarry Smith ierr = TaoComputeObjective(tao,X,&f);CHKERRQ(ierr); 73f4c1ad5cSStefano Zampini ierr = VecSetValue(X,i,2.0*h,ADD_VALUES);CHKERRQ(ierr); 74f4c1ad5cSStefano Zampini ierr = VecAssemblyBegin(X);CHKERRQ(ierr); 75f4c1ad5cSStefano Zampini ierr = VecAssemblyEnd(X);CHKERRQ(ierr); 76a7e14dcfSSatish Balay ierr = TaoComputeObjective(tao,X,&f2);CHKERRQ(ierr); 77f4c1ad5cSStefano Zampini ierr = VecSetValue(X,i,-h,ADD_VALUES);CHKERRQ(ierr); 78f4c1ad5cSStefano Zampini ierr = VecAssemblyBegin(X);CHKERRQ(ierr); 79f4c1ad5cSStefano 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); 85f4c1ad5cSStefano 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: 104f4c1ad5cSStefano 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 SNES snes; 121f4c1ad5cSStefano Zampini DM dm; 122a7e14dcfSSatish Balay 123a7e14dcfSSatish Balay PetscFunctionBegin; 124a7e14dcfSSatish Balay ierr = PetscInfo(tao,"TAO Using finite differences w/o coloring to compute Hessian matrix\n");CHKERRQ(ierr); 125f4c1ad5cSStefano Zampini ierr = SNESCreate(PetscObjectComm((PetscObject)H),&snes);CHKERRQ(ierr); 126dc822f48SStefano Zampini ierr = SNESSetFunction(snes,NULL,Fsnes,tao);CHKERRQ(ierr); 127f4c1ad5cSStefano Zampini ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 128f4c1ad5cSStefano Zampini ierr = DMShellSetGlobalVector(dm,V);CHKERRQ(ierr); 129f4c1ad5cSStefano Zampini ierr = SNESSetUp(snes);CHKERRQ(ierr); 130f4c1ad5cSStefano Zampini if (H) { 131f4c1ad5cSStefano Zampini PetscInt n,N; 132f4c1ad5cSStefano Zampini 133f4c1ad5cSStefano Zampini ierr = VecGetSize(V,&N);CHKERRQ(ierr); 134f4c1ad5cSStefano Zampini ierr = VecGetLocalSize(V,&n);CHKERRQ(ierr); 135f4c1ad5cSStefano Zampini ierr = MatSetSizes(H,n,n,N,N);CHKERRQ(ierr); 136f4c1ad5cSStefano Zampini ierr = MatSetUp(H);CHKERRQ(ierr); 137f4c1ad5cSStefano Zampini } 138f4c1ad5cSStefano Zampini if (B && B != H) { 139f4c1ad5cSStefano Zampini PetscInt n,N; 140f4c1ad5cSStefano Zampini 141f4c1ad5cSStefano Zampini ierr = VecGetSize(V,&N);CHKERRQ(ierr); 142f4c1ad5cSStefano Zampini ierr = VecGetLocalSize(V,&n);CHKERRQ(ierr); 143f4c1ad5cSStefano Zampini ierr = MatSetSizes(B,n,n,N,N);CHKERRQ(ierr); 144f4c1ad5cSStefano Zampini ierr = MatSetUp(B);CHKERRQ(ierr); 145f4c1ad5cSStefano Zampini } 146f4c1ad5cSStefano Zampini ierr = SNESComputeJacobianDefault(snes,V,H,B,NULL);CHKERRQ(ierr); 147a7e14dcfSSatish Balay ierr = SNESDestroy(&snes);CHKERRQ(ierr); 148a7e14dcfSSatish Balay PetscFunctionReturn(0); 149a7e14dcfSSatish Balay } 150a7e14dcfSSatish Balay 151a7e14dcfSSatish Balay /*@C 152a7e14dcfSSatish Balay TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences. 153a7e14dcfSSatish Balay 154441846f8SBarry Smith Collective on Tao 155a7e14dcfSSatish Balay 156a7e14dcfSSatish Balay Input Parameters: 157441846f8SBarry Smith + tao - the Tao context 158a7e14dcfSSatish Balay . V - compute Hessian at this point 159a7e14dcfSSatish Balay - ctx - the PetscColoring object (must be of type MatFDColoring) 160a7e14dcfSSatish Balay 161a7e14dcfSSatish Balay Output Parameters: 162a7e14dcfSSatish Balay + H - Hessian matrix (not altered in this routine) 16356744491SBarry Smith - B - newly computed Hessian matrix to use with preconditioner (generally the same as H) 164a7e14dcfSSatish Balay 165a7e14dcfSSatish Balay Level: advanced 166a7e14dcfSSatish Balay 167a7e14dcfSSatish Balay .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessian(),SNESComputeJacobianDefaultColor(), TaoSetGradientRoutine() 168a7e14dcfSSatish Balay @*/ 169ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessianColor(Tao tao,Vec V,Mat H,Mat B,void *ctx) 1705ca45b2bSBarry Smith { 171a7e14dcfSSatish Balay PetscErrorCode ierr; 172a7e14dcfSSatish Balay MatFDColoring coloring = (MatFDColoring)ctx; 173a7e14dcfSSatish Balay 174a7e14dcfSSatish Balay PetscFunctionBegin; 175064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,5); 176a7e14dcfSSatish Balay ierr = PetscInfo(tao,"TAO computing matrix using finite differences Hessian and coloring\n");CHKERRQ(ierr); 177d1e9a80fSBarry Smith ierr = MatFDColoringApply(B,coloring,V,ctx);CHKERRQ(ierr); 17894ab13aaSBarry Smith if (H != B) { 17994ab13aaSBarry Smith ierr = MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 18094ab13aaSBarry Smith ierr = MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 181a7e14dcfSSatish Balay } 182a7e14dcfSSatish Balay PetscFunctionReturn(0); 183a7e14dcfSSatish Balay } 184a7e14dcfSSatish Balay 185f4c1ad5cSStefano Zampini PetscErrorCode TaoDefaultComputeHessianMFFD(Tao tao,Vec X,Mat H,Mat B,void *ctx) 186f4c1ad5cSStefano Zampini { 187f4c1ad5cSStefano Zampini PetscInt n,N; 188c5e9d94cSAlp Dener PetscBool assembled; 189f4c1ad5cSStefano Zampini PetscErrorCode ierr; 190a7e14dcfSSatish Balay 191f4c1ad5cSStefano Zampini PetscFunctionBegin; 192f4c1ad5cSStefano Zampini if (B && B != H) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Preconditioning Hessian matrix"); 193c5e9d94cSAlp Dener ierr = MatAssembled(H, &assembled);CHKERRQ(ierr); 194c5e9d94cSAlp Dener if (!assembled) { 195f4c1ad5cSStefano Zampini ierr = VecGetSize(X,&N);CHKERRQ(ierr); 196f4c1ad5cSStefano Zampini ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 197f4c1ad5cSStefano Zampini ierr = MatSetSizes(H,n,n,N,N);CHKERRQ(ierr); 198f4c1ad5cSStefano Zampini ierr = MatSetType(H,MATMFFD);CHKERRQ(ierr); 199f4c1ad5cSStefano Zampini ierr = MatSetUp(H);CHKERRQ(ierr); 200f4c1ad5cSStefano Zampini ierr = MatMFFDSetFunction(H,(PetscErrorCode (*)(void*,Vec,Vec))TaoComputeGradient,tao);CHKERRQ(ierr); 201c5e9d94cSAlp Dener } 202*7faa0519SAlp Dener ierr = MatMFFDSetBase(H,X,NULL);CHKERRQ(ierr); 203f4c1ad5cSStefano Zampini ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 204f4c1ad5cSStefano Zampini ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 205f4c1ad5cSStefano Zampini PetscFunctionReturn(0); 206f4c1ad5cSStefano Zampini } 207