xref: /petsc/src/tao/interface/taosolver_fg.c (revision d8d19677bbccf95218448bee62e6b87f4513e133)
1af0996ceSBarry Smith #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
2a7e14dcfSSatish Balay 
3a7e14dcfSSatish Balay /*@
4a7e14dcfSSatish Balay   TaoSetInitialVector - Sets the initial guess for the solve
5a7e14dcfSSatish Balay 
6441846f8SBarry Smith   Logically collective on Tao
7a7e14dcfSSatish Balay 
8a7e14dcfSSatish Balay   Input Parameters:
9441846f8SBarry Smith + tao - the Tao context
10a7e14dcfSSatish Balay - x0  - the initial guess
11a7e14dcfSSatish Balay 
12a7e14dcfSSatish Balay   Level: beginner
13a7e14dcfSSatish Balay .seealso: TaoCreate(), TaoSolve()
14a7e14dcfSSatish Balay @*/
15a7e14dcfSSatish Balay 
16441846f8SBarry Smith PetscErrorCode TaoSetInitialVector(Tao tao, Vec x0)
1745cf516eSBarry Smith {
18a7e14dcfSSatish Balay   PetscErrorCode ierr;
19a7e14dcfSSatish Balay 
20a7e14dcfSSatish Balay   PetscFunctionBegin;
21441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
22a7e14dcfSSatish Balay   if (x0) {
23a7e14dcfSSatish Balay     PetscValidHeaderSpecific(x0,VEC_CLASSID,2);
24a7e14dcfSSatish Balay     PetscObjectReference((PetscObject)x0);
25a7e14dcfSSatish Balay   }
26a7e14dcfSSatish Balay   ierr = VecDestroy(&tao->solution);CHKERRQ(ierr);
27a7e14dcfSSatish Balay   tao->solution = x0;
28a7e14dcfSSatish Balay   PetscFunctionReturn(0);
29a7e14dcfSSatish Balay }
30a7e14dcfSSatish Balay 
31412cdd55SHong Zhang PetscErrorCode TaoTestGradient(Tao tao,Vec x,Vec g1)
3209baa881SHong Zhang {
33412cdd55SHong Zhang   Vec               g2,g3;
3409baa881SHong Zhang   PetscBool         complete_print = PETSC_FALSE,test = PETSC_FALSE;
3509baa881SHong Zhang   PetscReal         hcnorm,fdnorm,hcmax,fdmax,diffmax,diffnorm;
3609baa881SHong Zhang   PetscScalar       dot;
3709baa881SHong Zhang   MPI_Comm          comm;
38913eda9aSHong Zhang   PetscViewer       viewer,mviewer;
39913eda9aSHong Zhang   PetscViewerFormat format;
4009baa881SHong Zhang   PetscInt          tabs;
4109baa881SHong Zhang   static PetscBool  directionsprinted = PETSC_FALSE;
4209baa881SHong Zhang   PetscErrorCode    ierr;
4309baa881SHong Zhang 
4409baa881SHong Zhang   PetscFunctionBegin;
4509baa881SHong Zhang   ierr = PetscObjectOptionsBegin((PetscObject)tao);CHKERRQ(ierr);
4609baa881SHong Zhang   ierr = PetscOptionsName("-tao_test_gradient","Compare hand-coded and finite difference Gradients","None",&test);CHKERRQ(ierr);
47913eda9aSHong Zhang   ierr = PetscOptionsViewer("-tao_test_gradient_view","View difference between hand-coded and finite difference Gradients element entries","None",&mviewer,&format,&complete_print);CHKERRQ(ierr);
4809baa881SHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
492f4b6201SAlp Dener   if (!test) {
502f4b6201SAlp Dener     if (complete_print) {
512f4b6201SAlp Dener       ierr = PetscViewerDestroy(&mviewer);CHKERRQ(ierr);
522f4b6201SAlp Dener     }
532f4b6201SAlp Dener     PetscFunctionReturn(0);
542f4b6201SAlp Dener   }
5509baa881SHong Zhang 
5609baa881SHong Zhang   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
5709baa881SHong Zhang   ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr);
5809baa881SHong Zhang   ierr = PetscViewerASCIIGetTab(viewer, &tabs);CHKERRQ(ierr);
5909baa881SHong Zhang   ierr = PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel);CHKERRQ(ierr);
6009baa881SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ---------- Testing Gradient -------------\n");CHKERRQ(ierr);
6109baa881SHong Zhang   if (!complete_print && !directionsprinted) {
6209baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Run with -tao_test_gradient_view and optionally -tao_test_gradient <threshold> to show difference\n");CHKERRQ(ierr);
6309baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    of hand-coded and finite difference gradient entries greater than <threshold>.\n");CHKERRQ(ierr);
6409baa881SHong Zhang   }
6509baa881SHong Zhang   if (!directionsprinted) {
6687cb359fSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Testing hand-coded Gradient, if (for double precision runs) ||G - Gfd||/||G|| is\n");CHKERRQ(ierr);
67ae75128cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"    O(1.e-8), the hand-coded Gradient is probably correct.\n");CHKERRQ(ierr);
6809baa881SHong Zhang     directionsprinted = PETSC_TRUE;
6909baa881SHong Zhang   }
70913eda9aSHong Zhang   if (complete_print) {
71913eda9aSHong Zhang     ierr = PetscViewerPushFormat(mviewer,format);CHKERRQ(ierr);
72913eda9aSHong Zhang   }
7309baa881SHong Zhang 
7409baa881SHong Zhang   ierr = VecDuplicate(x,&g2);CHKERRQ(ierr);
7509baa881SHong Zhang   ierr = VecDuplicate(x,&g3);CHKERRQ(ierr);
7609baa881SHong Zhang 
7709baa881SHong Zhang   /* Compute finite difference gradient, assume the gradient is already computed by TaoComputeGradient() and put into g1 */
7809baa881SHong Zhang   ierr = TaoDefaultComputeGradient(tao,x,g2,NULL);CHKERRQ(ierr);
7909baa881SHong Zhang 
8009baa881SHong Zhang   ierr = VecNorm(g2,NORM_2,&fdnorm);CHKERRQ(ierr);
8109baa881SHong Zhang   ierr = VecNorm(g1,NORM_2,&hcnorm);CHKERRQ(ierr);
8209baa881SHong Zhang   ierr = VecNorm(g2,NORM_INFINITY,&fdmax);CHKERRQ(ierr);
8309baa881SHong Zhang   ierr = VecNorm(g1,NORM_INFINITY,&hcmax);CHKERRQ(ierr);
8409baa881SHong Zhang   ierr = VecDot(g1,g2,&dot);CHKERRQ(ierr);
8509baa881SHong Zhang   ierr = VecCopy(g1,g3);CHKERRQ(ierr);
8609baa881SHong Zhang   ierr = VecAXPY(g3,-1.0,g2);CHKERRQ(ierr);
8709baa881SHong Zhang   ierr = VecNorm(g3,NORM_2,&diffnorm);CHKERRQ(ierr);
8809baa881SHong Zhang   ierr = VecNorm(g3,NORM_INFINITY,&diffmax);CHKERRQ(ierr);
8909baa881SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ||Gfd|| %g, ||G|| = %g, angle cosine = (Gfd'G)/||Gfd||||G|| = %g\n", (double)fdnorm, (double)hcnorm, (double)(PetscRealPart(dot)/(fdnorm*hcnorm)));CHKERRQ(ierr);
9009baa881SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  2-norm ||G - Gfd||/||G|| = %g, ||G - Gfd|| = %g\n",(double)(diffnorm/PetscMax(hcnorm,fdnorm)),(double)diffnorm);CHKERRQ(ierr);
9109baa881SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  max-norm ||G - Gfd||/||G|| = %g, ||G - Gfd|| = %g\n",(double)(diffmax/PetscMax(hcmax,fdmax)),(double)diffmax);CHKERRQ(ierr);
9209baa881SHong Zhang 
9309baa881SHong Zhang   if (complete_print) {
9409baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Hand-coded gradient ----------\n");CHKERRQ(ierr);
95913eda9aSHong Zhang     ierr = VecView(g1,mviewer);CHKERRQ(ierr);
9609baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Finite difference gradient ----------\n");CHKERRQ(ierr);
97913eda9aSHong Zhang     ierr = VecView(g2,mviewer);CHKERRQ(ierr);
9809baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Hand-coded minus finite-difference gradient ----------\n");CHKERRQ(ierr);
99913eda9aSHong Zhang     ierr = VecView(g3,mviewer);CHKERRQ(ierr);
10009baa881SHong Zhang   }
10109baa881SHong Zhang   ierr = VecDestroy(&g2);CHKERRQ(ierr);
10209baa881SHong Zhang   ierr = VecDestroy(&g3);CHKERRQ(ierr);
103913eda9aSHong Zhang 
104913eda9aSHong Zhang   if (complete_print) {
105913eda9aSHong Zhang     ierr = PetscViewerPopFormat(mviewer);CHKERRQ(ierr);
1062f4b6201SAlp Dener     ierr = PetscViewerDestroy(&mviewer);CHKERRQ(ierr);
107913eda9aSHong Zhang   }
108913eda9aSHong Zhang   ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
10909baa881SHong Zhang   PetscFunctionReturn(0);
11009baa881SHong Zhang }
11109baa881SHong Zhang 
112a7e14dcfSSatish Balay /*@
113a7e14dcfSSatish Balay   TaoComputeGradient - Computes the gradient of the objective function
114a7e14dcfSSatish Balay 
115441846f8SBarry Smith   Collective on Tao
116a7e14dcfSSatish Balay 
117a7e14dcfSSatish Balay   Input Parameters:
118441846f8SBarry Smith + tao - the Tao context
119a7e14dcfSSatish Balay - X - input vector
120a7e14dcfSSatish Balay 
121a7e14dcfSSatish Balay   Output Parameter:
122a7e14dcfSSatish Balay . G - gradient vector
123a7e14dcfSSatish Balay 
12409baa881SHong Zhang   Options Database Keys:
12509baa881SHong Zhang +    -tao_test_gradient - compare the user provided gradient with one compute via finite differences to check for errors
126dfe02fe6SHong Zhang -    -tao_test_gradient_view - display the user provided gradient, the finite difference gradient and the difference between them to help users detect the location of errors in the user provided gradient
12709baa881SHong Zhang 
12895452b02SPatrick Sanan   Notes:
12995452b02SPatrick Sanan     TaoComputeGradient() is typically used within minimization implementations,
130a7e14dcfSSatish Balay   so most users would not generally call this routine themselves.
131a7e14dcfSSatish Balay 
132a7e14dcfSSatish Balay   Level: advanced
133a7e14dcfSSatish Balay 
134a7e14dcfSSatish Balay .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetGradientRoutine()
135a7e14dcfSSatish Balay @*/
136441846f8SBarry Smith PetscErrorCode TaoComputeGradient(Tao tao, Vec X, Vec G)
137a7e14dcfSSatish Balay {
138a7e14dcfSSatish Balay   PetscErrorCode ierr;
139a7e14dcfSSatish Balay   PetscReal      dummy;
14087f595a5SBarry Smith 
141a7e14dcfSSatish Balay   PetscFunctionBegin;
142441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
143a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
144064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(G,VEC_CLASSID,3);
145a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
146a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,G,3);
1478860a134SJunchao Zhang   ierr = VecLockReadPush(X);CHKERRQ(ierr);
148a7e14dcfSSatish Balay   if (tao->ops->computegradient) {
1490ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
150441846f8SBarry Smith     PetscStackPush("Tao user gradient evaluation routine");
151a7e14dcfSSatish Balay     ierr = (*tao->ops->computegradient)(tao,X,G,tao->user_gradP);CHKERRQ(ierr);
152a7e14dcfSSatish Balay     PetscStackPop;
1530ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
154a7e14dcfSSatish Balay     tao->ngrads++;
155a7e14dcfSSatish Balay   } else if (tao->ops->computeobjectiveandgradient) {
1560ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjGradEval,tao,X,G,NULL);CHKERRQ(ierr);
157a7e14dcfSSatish Balay     PetscStackPush("Tao user objective/gradient evaluation routine");
158a7e14dcfSSatish Balay     ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,&dummy,G,tao->user_objgradP);CHKERRQ(ierr);
159a7e14dcfSSatish Balay     PetscStackPop;
1600ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjGradEval,tao,X,G,NULL);CHKERRQ(ierr);
161a7e14dcfSSatish Balay     tao->nfuncgrads++;
162691b26d3SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetGradientRoutine() has not been called");
1638860a134SJunchao Zhang   ierr = VecLockReadPop(X);CHKERRQ(ierr);
16409baa881SHong Zhang 
165412cdd55SHong Zhang   ierr = TaoTestGradient(tao,X,G);CHKERRQ(ierr);
166a7e14dcfSSatish Balay   PetscFunctionReturn(0);
167a7e14dcfSSatish Balay }
168a7e14dcfSSatish Balay 
169a7e14dcfSSatish Balay /*@
170a7e14dcfSSatish Balay   TaoComputeObjective - Computes the objective function value at a given point
171a7e14dcfSSatish Balay 
172441846f8SBarry Smith   Collective on Tao
173a7e14dcfSSatish Balay 
174a7e14dcfSSatish Balay   Input Parameters:
175441846f8SBarry Smith + tao - the Tao context
176a7e14dcfSSatish Balay - X - input vector
177a7e14dcfSSatish Balay 
178a7e14dcfSSatish Balay   Output Parameter:
179a7e14dcfSSatish Balay . f - Objective value at X
180a7e14dcfSSatish Balay 
18195452b02SPatrick Sanan   Notes:
18295452b02SPatrick Sanan     TaoComputeObjective() is typically used within minimization implementations,
183a7e14dcfSSatish Balay   so most users would not generally call this routine themselves.
184a7e14dcfSSatish Balay 
185a7e14dcfSSatish Balay   Level: advanced
186a7e14dcfSSatish Balay 
187a7e14dcfSSatish Balay .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
188a7e14dcfSSatish Balay @*/
189441846f8SBarry Smith PetscErrorCode TaoComputeObjective(Tao tao, Vec X, PetscReal *f)
190a7e14dcfSSatish Balay {
191a7e14dcfSSatish Balay   PetscErrorCode ierr;
192a7e14dcfSSatish Balay   Vec            temp;
19387f595a5SBarry Smith 
194a7e14dcfSSatish Balay   PetscFunctionBegin;
195441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
196a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
197a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
1988860a134SJunchao Zhang   ierr = VecLockReadPush(X);CHKERRQ(ierr);
199a7e14dcfSSatish Balay   if (tao->ops->computeobjective) {
2000ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
201441846f8SBarry Smith     PetscStackPush("Tao user objective evaluation routine");
202a7e14dcfSSatish Balay     ierr = (*tao->ops->computeobjective)(tao,X,f,tao->user_objP);CHKERRQ(ierr);
203a7e14dcfSSatish Balay     PetscStackPop;
2040ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
205a7e14dcfSSatish Balay     tao->nfuncs++;
206a7e14dcfSSatish Balay   } else if (tao->ops->computeobjectiveandgradient) {
207955c1f14SBarry Smith     ierr = PetscInfo(tao,"Duplicating variable vector in order to call func/grad routine\n");CHKERRQ(ierr);
208a7e14dcfSSatish Balay     ierr = VecDuplicate(X,&temp);CHKERRQ(ierr);
2090ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjGradEval,tao,X,NULL,NULL);CHKERRQ(ierr);
210441846f8SBarry Smith     PetscStackPush("Tao user objective/gradient evaluation routine");
211a7e14dcfSSatish Balay     ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,f,temp,tao->user_objgradP);CHKERRQ(ierr);
212a7e14dcfSSatish Balay     PetscStackPop;
2130ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjGradEval,tao,X,NULL,NULL);CHKERRQ(ierr);
214a7e14dcfSSatish Balay     ierr = VecDestroy(&temp);CHKERRQ(ierr);
215a7e14dcfSSatish Balay     tao->nfuncgrads++;
216691b26d3SBarry Smith   }  else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjectiveRoutine() has not been called");
217e356b196STodd Munson   ierr = PetscInfo1(tao,"TAO Function evaluation: %20.19e\n",(double)(*f));CHKERRQ(ierr);
2188860a134SJunchao Zhang   ierr = VecLockReadPop(X);CHKERRQ(ierr);
219a7e14dcfSSatish Balay   PetscFunctionReturn(0);
220a7e14dcfSSatish Balay }
221a7e14dcfSSatish Balay 
222a7e14dcfSSatish Balay /*@
223a7e14dcfSSatish Balay   TaoComputeObjectiveAndGradient - Computes the objective function value at a given point
224a7e14dcfSSatish Balay 
225441846f8SBarry Smith   Collective on Tao
226a7e14dcfSSatish Balay 
227a7e14dcfSSatish Balay   Input Parameters:
228441846f8SBarry Smith + tao - the Tao context
229a7e14dcfSSatish Balay - X - input vector
230a7e14dcfSSatish Balay 
231*d8d19677SJose E. Roman   Output Parameters:
232a7e14dcfSSatish Balay + f - Objective value at X
233a7e14dcfSSatish Balay - g - Gradient vector at X
234a7e14dcfSSatish Balay 
23595452b02SPatrick Sanan   Notes:
23695452b02SPatrick Sanan     TaoComputeObjectiveAndGradient() is typically used within minimization implementations,
237a7e14dcfSSatish Balay   so most users would not generally call this routine themselves.
238a7e14dcfSSatish Balay 
239a7e14dcfSSatish Balay   Level: advanced
240a7e14dcfSSatish Balay 
241a7e14dcfSSatish Balay .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
242a7e14dcfSSatish Balay @*/
243441846f8SBarry Smith PetscErrorCode TaoComputeObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec G)
244a7e14dcfSSatish Balay {
245a7e14dcfSSatish Balay   PetscErrorCode ierr;
24687f595a5SBarry Smith 
247a7e14dcfSSatish Balay   PetscFunctionBegin;
248441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
249a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
250a7e14dcfSSatish Balay   PetscValidHeaderSpecific(G,VEC_CLASSID,4);
251a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
252a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,G,4);
2538860a134SJunchao Zhang   ierr = VecLockReadPush(X);CHKERRQ(ierr);
254a7e14dcfSSatish Balay   if (tao->ops->computeobjectiveandgradient) {
2550ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjGradEval,tao,X,G,NULL);CHKERRQ(ierr);
256f4c1ad5cSStefano Zampini     if (tao->ops->computegradient == TaoDefaultComputeGradient) {
257f4c1ad5cSStefano Zampini       ierr = TaoComputeObjective(tao,X,f);CHKERRQ(ierr);
258f4c1ad5cSStefano Zampini       ierr = TaoDefaultComputeGradient(tao,X,G,NULL);CHKERRQ(ierr);
259f4c1ad5cSStefano Zampini     } else {
260441846f8SBarry Smith       PetscStackPush("Tao user objective/gradient evaluation routine");
261a7e14dcfSSatish Balay       ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,f,G,tao->user_objgradP);CHKERRQ(ierr);
2620cbffdbaSBarry Smith       PetscStackPop;
263a7e14dcfSSatish Balay     }
2640ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjGradEval,tao,X,G,NULL);CHKERRQ(ierr);
265a7e14dcfSSatish Balay     tao->nfuncgrads++;
266a7e14dcfSSatish Balay   } else if (tao->ops->computeobjective && tao->ops->computegradient) {
2670ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
268441846f8SBarry Smith     PetscStackPush("Tao user objective evaluation routine");
269a7e14dcfSSatish Balay     ierr = (*tao->ops->computeobjective)(tao,X,f,tao->user_objP);CHKERRQ(ierr);
270a7e14dcfSSatish Balay     PetscStackPop;
2710ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
272a7e14dcfSSatish Balay     tao->nfuncs++;
2730ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
274441846f8SBarry Smith     PetscStackPush("Tao user gradient evaluation routine");
275a7e14dcfSSatish Balay     ierr = (*tao->ops->computegradient)(tao,X,G,tao->user_gradP);CHKERRQ(ierr);
276a7e14dcfSSatish Balay     PetscStackPop;
2770ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
278a7e14dcfSSatish Balay     tao->ngrads++;
279691b26d3SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjectiveRoutine() or TaoSetGradientRoutine() not set");
280e356b196STodd Munson   ierr = PetscInfo1(tao,"TAO Function evaluation: %20.19e\n",(double)(*f));CHKERRQ(ierr);
2818860a134SJunchao Zhang   ierr = VecLockReadPop(X);CHKERRQ(ierr);
2821657496cSHong Zhang 
283412cdd55SHong Zhang   ierr = TaoTestGradient(tao,X,G);CHKERRQ(ierr);
284a7e14dcfSSatish Balay   PetscFunctionReturn(0);
285a7e14dcfSSatish Balay }
286a7e14dcfSSatish Balay 
287a7e14dcfSSatish Balay /*@C
288a7e14dcfSSatish Balay   TaoSetObjectiveRoutine - Sets the function evaluation routine for minimization
289a7e14dcfSSatish Balay 
290441846f8SBarry Smith   Logically collective on Tao
291a7e14dcfSSatish Balay 
292*d8d19677SJose E. Roman   Input Parameters:
293441846f8SBarry Smith + tao - the Tao context
294a7e14dcfSSatish Balay . func - the objective function
295a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the function evaluation
2966c23d075SBarry Smith         routine (may be NULL)
297a7e14dcfSSatish Balay 
298a7e14dcfSSatish Balay   Calling sequence of func:
299441846f8SBarry Smith $      func (Tao tao, Vec x, PetscReal *f, void *ctx);
300a7e14dcfSSatish Balay 
301a7e14dcfSSatish Balay + x - input vector
302a7e14dcfSSatish Balay . f - function value
303a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
304a7e14dcfSSatish Balay 
305a7e14dcfSSatish Balay   Level: beginner
306a7e14dcfSSatish Balay 
307a7e14dcfSSatish Balay .seealso: TaoSetGradientRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
308a7e14dcfSSatish Balay @*/
309441846f8SBarry Smith PetscErrorCode TaoSetObjectiveRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*,void*),void *ctx)
310a7e14dcfSSatish Balay {
311a7e14dcfSSatish Balay   PetscFunctionBegin;
312441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
313a7e14dcfSSatish Balay   tao->user_objP = ctx;
314a7e14dcfSSatish Balay   tao->ops->computeobjective = func;
315a7e14dcfSSatish Balay   PetscFunctionReturn(0);
316a7e14dcfSSatish Balay }
317a7e14dcfSSatish Balay 
318a7e14dcfSSatish Balay /*@C
3194a48860cSAlp Dener   TaoSetResidualRoutine - Sets the residual evaluation routine for least-square applications
320a7e14dcfSSatish Balay 
321441846f8SBarry Smith   Logically collective on Tao
322a7e14dcfSSatish Balay 
323*d8d19677SJose E. Roman   Input Parameters:
324441846f8SBarry Smith + tao - the Tao context
3254a48860cSAlp Dener . func - the residual evaluation routine
326a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the function evaluation
3276c23d075SBarry Smith         routine (may be NULL)
328a7e14dcfSSatish Balay 
329a7e14dcfSSatish Balay   Calling sequence of func:
330441846f8SBarry Smith $      func (Tao tao, Vec x, Vec f, void *ctx);
331a7e14dcfSSatish Balay 
332a7e14dcfSSatish Balay + x - input vector
333a7e14dcfSSatish Balay . f - function value vector
334a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
335a7e14dcfSSatish Balay 
336a7e14dcfSSatish Balay   Level: beginner
337a7e14dcfSSatish Balay 
338a7e14dcfSSatish Balay .seealso: TaoSetObjectiveRoutine(), TaoSetJacobianRoutine()
339a7e14dcfSSatish Balay @*/
3404a48860cSAlp Dener PetscErrorCode TaoSetResidualRoutine(Tao tao, Vec res, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
341a7e14dcfSSatish Balay {
342737f463aSAlp Dener   PetscErrorCode ierr;
343737f463aSAlp Dener 
344a7e14dcfSSatish Balay   PetscFunctionBegin;
345441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
3464a48860cSAlp Dener   PetscValidHeaderSpecific(res,VEC_CLASSID,2);
3474ffbe8acSAlp Dener   ierr = PetscObjectReference((PetscObject)res);CHKERRQ(ierr);
348737f463aSAlp Dener   if (tao->ls_res) {
349737f463aSAlp Dener     ierr = VecDestroy(&tao->ls_res);CHKERRQ(ierr);
350737f463aSAlp Dener   }
3514a48860cSAlp Dener   tao->ls_res = res;
3524ffbe8acSAlp Dener   tao->user_lsresP = ctx;
3534a48860cSAlp Dener   tao->ops->computeresidual = func;
354737f463aSAlp Dener 
355737f463aSAlp Dener   PetscFunctionReturn(0);
356737f463aSAlp Dener }
357737f463aSAlp Dener 
358737f463aSAlp Dener /*@
359737f463aSAlp Dener   TaoSetResidualWeights - Give weights for the residual values. A vector can be used if only diagonal terms are used, otherwise a matrix can be give. If this function is not used, or if sigma_v and sigma_w are both NULL, then the default identity matrix will be used for weights.
360737f463aSAlp Dener 
361737f463aSAlp Dener   Collective on Tao
362737f463aSAlp Dener 
363737f463aSAlp Dener   Input Parameters:
364737f463aSAlp Dener + tao - the Tao context
365737f463aSAlp Dener . sigma_v - vector of weights (diagonal terms only)
366737f463aSAlp Dener . n       - the number of weights (if using off-diagonal)
367737f463aSAlp Dener . rows    - index list of rows for sigma_w
368737f463aSAlp Dener . cols    - index list of columns for sigma_w
369737f463aSAlp Dener - vals - array of weights
370737f463aSAlp Dener 
371737f463aSAlp Dener   Note: Either sigma_v or sigma_w (or both) should be NULL
372737f463aSAlp Dener 
373737f463aSAlp Dener   Level: intermediate
374737f463aSAlp Dener 
375737f463aSAlp Dener .seealso: TaoSetResidualRoutine()
376737f463aSAlp Dener @*/
377737f463aSAlp Dener PetscErrorCode TaoSetResidualWeights(Tao tao, Vec sigma_v, PetscInt n, PetscInt *rows, PetscInt *cols, PetscReal *vals)
378737f463aSAlp Dener {
379737f463aSAlp Dener   PetscErrorCode ierr;
380737f463aSAlp Dener   PetscInt       i;
381737f463aSAlp Dener   PetscFunctionBegin;
382737f463aSAlp Dener   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
383737f463aSAlp Dener   if (sigma_v) {
3844ffbe8acSAlp Dener     PetscValidHeaderSpecific(sigma_v,VEC_CLASSID,2);
385737f463aSAlp Dener     ierr = PetscObjectReference((PetscObject)sigma_v);CHKERRQ(ierr);
386737f463aSAlp Dener   }
3874ffbe8acSAlp Dener   if (tao->res_weights_v) {
3884ffbe8acSAlp Dener     ierr = VecDestroy(&tao->res_weights_v);CHKERRQ(ierr);
3894ffbe8acSAlp Dener   }
3904ffbe8acSAlp Dener   tao->res_weights_v=sigma_v;
391737f463aSAlp Dener   if (vals) {
392737f463aSAlp Dener     if (tao->res_weights_n) {
393737f463aSAlp Dener       ierr = PetscFree(tao->res_weights_rows);CHKERRQ(ierr);
394737f463aSAlp Dener       ierr = PetscFree(tao->res_weights_cols);CHKERRQ(ierr);
395737f463aSAlp Dener       ierr = PetscFree(tao->res_weights_w);CHKERRQ(ierr);
396737f463aSAlp Dener     }
397737f463aSAlp Dener     ierr = PetscMalloc1(n,&tao->res_weights_rows);CHKERRQ(ierr);
398737f463aSAlp Dener     ierr = PetscMalloc1(n,&tao->res_weights_cols);CHKERRQ(ierr);
399737f463aSAlp Dener     ierr = PetscMalloc1(n,&tao->res_weights_w);CHKERRQ(ierr);
400737f463aSAlp Dener     tao->res_weights_n=n;
401737f463aSAlp Dener     for (i=0;i<n;i++) {
402737f463aSAlp Dener       tao->res_weights_rows[i]=rows[i];
403737f463aSAlp Dener       tao->res_weights_cols[i]=cols[i];
404737f463aSAlp Dener       tao->res_weights_w[i]=vals[i];
405737f463aSAlp Dener     }
406737f463aSAlp Dener   } else {
407737f463aSAlp Dener     tao->res_weights_n=0;
40883c8fe1dSLisandro Dalcin     tao->res_weights_rows=NULL;
40983c8fe1dSLisandro Dalcin     tao->res_weights_cols=NULL;
410737f463aSAlp Dener   }
411a7e14dcfSSatish Balay   PetscFunctionReturn(0);
412a7e14dcfSSatish Balay }
413a7e14dcfSSatish Balay 
4148b7a9b22SJason Sarich /*@
4154a48860cSAlp Dener   TaoComputeResidual - Computes a least-squares residual vector at a given point
416a7e14dcfSSatish Balay 
417441846f8SBarry Smith   Collective on Tao
418a7e14dcfSSatish Balay 
419a7e14dcfSSatish Balay   Input Parameters:
420441846f8SBarry Smith + tao - the Tao context
421a7e14dcfSSatish Balay - X - input vector
422a7e14dcfSSatish Balay 
423a7e14dcfSSatish Balay   Output Parameter:
424a7e14dcfSSatish Balay . f - Objective vector at X
425a7e14dcfSSatish Balay 
42695452b02SPatrick Sanan   Notes:
4274a48860cSAlp Dener     TaoComputeResidual() is typically used within minimization implementations,
428a7e14dcfSSatish Balay   so most users would not generally call this routine themselves.
429a7e14dcfSSatish Balay 
430a7e14dcfSSatish Balay   Level: advanced
431a7e14dcfSSatish Balay 
4324a48860cSAlp Dener .seealso: TaoSetResidualRoutine()
433a7e14dcfSSatish Balay @*/
4344a48860cSAlp Dener PetscErrorCode TaoComputeResidual(Tao tao, Vec X, Vec F)
435a7e14dcfSSatish Balay {
436a7e14dcfSSatish Balay   PetscErrorCode ierr;
43787f595a5SBarry Smith 
438a7e14dcfSSatish Balay   PetscFunctionBegin;
439441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
440a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
441a7e14dcfSSatish Balay   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
442a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
443a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,F,3);
4444a48860cSAlp Dener   if (tao->ops->computeresidual) {
4450ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
4464a48860cSAlp Dener     PetscStackPush("Tao user least-squares residual evaluation routine");
4474a48860cSAlp Dener     ierr = (*tao->ops->computeresidual)(tao,X,F,tao->user_lsresP);CHKERRQ(ierr);
448a7e14dcfSSatish Balay     PetscStackPop;
4490ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
450a7e14dcfSSatish Balay     tao->nfuncs++;
451691b26d3SBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"TaoSetResidualRoutine() has not been called");
4524a48860cSAlp Dener   ierr = PetscInfo(tao,"TAO least-squares residual evaluation.\n");CHKERRQ(ierr);
453a7e14dcfSSatish Balay   PetscFunctionReturn(0);
454a7e14dcfSSatish Balay }
455a7e14dcfSSatish Balay 
456a7e14dcfSSatish Balay /*@C
457a7e14dcfSSatish Balay   TaoSetGradientRoutine - Sets the gradient evaluation routine for minimization
458a7e14dcfSSatish Balay 
459441846f8SBarry Smith   Logically collective on Tao
460a7e14dcfSSatish Balay 
461*d8d19677SJose E. Roman   Input Parameters:
462441846f8SBarry Smith + tao - the Tao context
463a7e14dcfSSatish Balay . func - the gradient function
464a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the gradient evaluation
4656c23d075SBarry Smith         routine (may be NULL)
466a7e14dcfSSatish Balay 
467a7e14dcfSSatish Balay   Calling sequence of func:
468441846f8SBarry Smith $      func (Tao tao, Vec x, Vec g, void *ctx);
469a7e14dcfSSatish Balay 
470a7e14dcfSSatish Balay + x - input vector
471a7e14dcfSSatish Balay . g - gradient value (output)
472a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
473a7e14dcfSSatish Balay 
474a7e14dcfSSatish Balay   Level: beginner
475a7e14dcfSSatish Balay 
476a7e14dcfSSatish Balay .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
477a7e14dcfSSatish Balay @*/
478441846f8SBarry Smith PetscErrorCode TaoSetGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
479a7e14dcfSSatish Balay {
480a7e14dcfSSatish Balay   PetscFunctionBegin;
481441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
482a7e14dcfSSatish Balay   tao->user_gradP = ctx;
483a7e14dcfSSatish Balay   tao->ops->computegradient = func;
484a7e14dcfSSatish Balay   PetscFunctionReturn(0);
485a7e14dcfSSatish Balay }
486a7e14dcfSSatish Balay 
487a7e14dcfSSatish Balay /*@C
488a7e14dcfSSatish Balay   TaoSetObjectiveAndGradientRoutine - Sets a combined objective function and gradient evaluation routine for minimization
489a7e14dcfSSatish Balay 
490441846f8SBarry Smith   Logically collective on Tao
491a7e14dcfSSatish Balay 
492*d8d19677SJose E. Roman   Input Parameters:
493441846f8SBarry Smith + tao - the Tao context
494a7e14dcfSSatish Balay . func - the gradient function
495a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the gradient evaluation
4966c23d075SBarry Smith         routine (may be NULL)
497a7e14dcfSSatish Balay 
498a7e14dcfSSatish Balay   Calling sequence of func:
49917477c02SJason Sarich $      func (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx);
500a7e14dcfSSatish Balay 
501a7e14dcfSSatish Balay + x - input vector
50217477c02SJason Sarich . f - objective value (output)
503a7e14dcfSSatish Balay . g - gradient value (output)
504a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
505a7e14dcfSSatish Balay 
506a7e14dcfSSatish Balay   Level: beginner
507a7e14dcfSSatish Balay 
508a7e14dcfSSatish Balay .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
509a7e14dcfSSatish Balay @*/
510441846f8SBarry Smith PetscErrorCode TaoSetObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void*), void *ctx)
511a7e14dcfSSatish Balay {
512a7e14dcfSSatish Balay   PetscFunctionBegin;
513441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
514a7e14dcfSSatish Balay   tao->user_objgradP = ctx;
515a7e14dcfSSatish Balay   tao->ops->computeobjectiveandgradient = func;
516a7e14dcfSSatish Balay   PetscFunctionReturn(0);
517a7e14dcfSSatish Balay }
518a7e14dcfSSatish Balay 
519a7e14dcfSSatish Balay /*@
520a7e14dcfSSatish Balay   TaoIsObjectiveDefined -- Checks to see if the user has
521a7e14dcfSSatish Balay   declared an objective-only routine.  Useful for determining when
522a7e14dcfSSatish Balay   it is appropriate to call TaoComputeObjective() or
523a7e14dcfSSatish Balay   TaoComputeObjectiveAndGradient()
524a7e14dcfSSatish Balay 
525441846f8SBarry Smith   Collective on Tao
526a7e14dcfSSatish Balay 
527*d8d19677SJose E. Roman   Input Parameters:
528441846f8SBarry Smith + tao - the Tao context
529a7e14dcfSSatish Balay - ctx - PETSC_TRUE if objective function routine is set by user,
530a7e14dcfSSatish Balay         PETSC_FALSE otherwise
531a7e14dcfSSatish Balay   Level: developer
532a7e14dcfSSatish Balay 
533a7e14dcfSSatish Balay .seealso: TaoSetObjectiveRoutine(), TaoIsGradientDefined(), TaoIsObjectiveAndGradientDefined()
534a7e14dcfSSatish Balay @*/
535441846f8SBarry Smith PetscErrorCode TaoIsObjectiveDefined(Tao tao, PetscBool *flg)
536a7e14dcfSSatish Balay {
537a7e14dcfSSatish Balay   PetscFunctionBegin;
538441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
53983c8fe1dSLisandro Dalcin   if (tao->ops->computeobjective == NULL) *flg = PETSC_FALSE;
54045cf516eSBarry Smith   else *flg = PETSC_TRUE;
541a7e14dcfSSatish Balay   PetscFunctionReturn(0);
542a7e14dcfSSatish Balay }
543a7e14dcfSSatish Balay 
544a7e14dcfSSatish Balay /*@
545a7e14dcfSSatish Balay   TaoIsGradientDefined -- Checks to see if the user has
546a7e14dcfSSatish Balay   declared an objective-only routine.  Useful for determining when
547a7e14dcfSSatish Balay   it is appropriate to call TaoComputeGradient() or
548a7e14dcfSSatish Balay   TaoComputeGradientAndGradient()
549a7e14dcfSSatish Balay 
550a7e14dcfSSatish Balay   Not Collective
551a7e14dcfSSatish Balay 
552*d8d19677SJose E. Roman   Input Parameters:
553441846f8SBarry Smith + tao - the Tao context
554a7e14dcfSSatish Balay - ctx - PETSC_TRUE if gradient routine is set by user, PETSC_FALSE otherwise
555a7e14dcfSSatish Balay   Level: developer
556a7e14dcfSSatish Balay 
557a7e14dcfSSatish Balay .seealso: TaoSetGradientRoutine(), TaoIsObjectiveDefined(), TaoIsObjectiveAndGradientDefined()
558a7e14dcfSSatish Balay @*/
559441846f8SBarry Smith PetscErrorCode TaoIsGradientDefined(Tao tao, PetscBool *flg)
560a7e14dcfSSatish Balay {
561a7e14dcfSSatish Balay   PetscFunctionBegin;
562441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
56383c8fe1dSLisandro Dalcin   if (tao->ops->computegradient == NULL) *flg = PETSC_FALSE;
56445cf516eSBarry Smith   else *flg = PETSC_TRUE;
565a7e14dcfSSatish Balay   PetscFunctionReturn(0);
566a7e14dcfSSatish Balay }
567a7e14dcfSSatish Balay 
568a7e14dcfSSatish Balay /*@
569a7e14dcfSSatish Balay   TaoIsObjectiveAndGradientDefined -- Checks to see if the user has
570a7e14dcfSSatish Balay   declared a joint objective/gradient routine.  Useful for determining when
571a7e14dcfSSatish Balay   it is appropriate to call TaoComputeObjective() or
572a7e14dcfSSatish Balay   TaoComputeObjectiveAndGradient()
573a7e14dcfSSatish Balay 
574a7e14dcfSSatish Balay   Not Collective
575a7e14dcfSSatish Balay 
576*d8d19677SJose E. Roman   Input Parameters:
577441846f8SBarry Smith + tao - the Tao context
578a7e14dcfSSatish Balay - ctx - PETSC_TRUE if objective/gradient routine is set by user, PETSC_FALSE otherwise
579a7e14dcfSSatish Balay   Level: developer
580a7e14dcfSSatish Balay 
581a7e14dcfSSatish Balay .seealso: TaoSetObjectiveAndGradientRoutine(), TaoIsObjectiveDefined(), TaoIsGradientDefined()
582a7e14dcfSSatish Balay @*/
583441846f8SBarry Smith PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao tao, PetscBool *flg)
584a7e14dcfSSatish Balay {
585a7e14dcfSSatish Balay   PetscFunctionBegin;
586441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
58783c8fe1dSLisandro Dalcin   if (tao->ops->computeobjectiveandgradient == NULL) *flg = PETSC_FALSE;
58845cf516eSBarry Smith   else *flg = PETSC_TRUE;
589a7e14dcfSSatish Balay   PetscFunctionReturn(0);
590a7e14dcfSSatish Balay }
591