xref: /petsc/src/tao/interface/taosolver_fg.c (revision ae75128c749e3b81916ddbb0b06a8d54f9604792)
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);
4909baa881SHong Zhang   if (!test) PetscFunctionReturn(0);
5009baa881SHong Zhang 
5109baa881SHong Zhang   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
5209baa881SHong Zhang   ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr);
5309baa881SHong Zhang   ierr = PetscViewerASCIIGetTab(viewer, &tabs);CHKERRQ(ierr);
5409baa881SHong Zhang   ierr = PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel);CHKERRQ(ierr);
5509baa881SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ---------- Testing Gradient -------------\n");CHKERRQ(ierr);
5609baa881SHong Zhang   if (!complete_print && !directionsprinted) {
5709baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Run with -tao_test_gradient_view and optionally -tao_test_gradient <threshold> to show difference\n");CHKERRQ(ierr);
5809baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    of hand-coded and finite difference gradient entries greater than <threshold>.\n");CHKERRQ(ierr);
5909baa881SHong Zhang   }
6009baa881SHong Zhang   if (!directionsprinted) {
6109baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Testing hand-coded Gradient, if (for double precision runs) ||G - Gfd||_F/||G||_F is\n");CHKERRQ(ierr);
62*ae75128cSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"    O(1.e-8), the hand-coded Gradient is probably correct.\n");CHKERRQ(ierr);
6309baa881SHong Zhang     directionsprinted = PETSC_TRUE;
6409baa881SHong Zhang   }
65913eda9aSHong Zhang   if (complete_print) {
66913eda9aSHong Zhang     ierr = PetscViewerPushFormat(mviewer,format);CHKERRQ(ierr);
67913eda9aSHong Zhang   }
6809baa881SHong Zhang 
6909baa881SHong Zhang   ierr = VecDuplicate(x,&g2);CHKERRQ(ierr);
7009baa881SHong Zhang   ierr = VecDuplicate(x,&g3);CHKERRQ(ierr);
7109baa881SHong Zhang 
7209baa881SHong Zhang   /* Compute finite difference gradient, assume the gradient is already computed by TaoComputeGradient() and put into g1 */
7309baa881SHong Zhang   ierr = TaoDefaultComputeGradient(tao,x,g2,NULL);CHKERRQ(ierr);
7409baa881SHong Zhang 
7509baa881SHong Zhang   ierr = VecNorm(g2,NORM_2,&fdnorm);CHKERRQ(ierr);
7609baa881SHong Zhang   ierr = VecNorm(g1,NORM_2,&hcnorm);CHKERRQ(ierr);
7709baa881SHong Zhang   ierr = VecNorm(g2,NORM_INFINITY,&fdmax);CHKERRQ(ierr);
7809baa881SHong Zhang   ierr = VecNorm(g1,NORM_INFINITY,&hcmax);CHKERRQ(ierr);
7909baa881SHong Zhang   ierr = VecDot(g1,g2,&dot);CHKERRQ(ierr);
8009baa881SHong Zhang   ierr = VecCopy(g1,g3);CHKERRQ(ierr);
8109baa881SHong Zhang   ierr = VecAXPY(g3,-1.0,g2);CHKERRQ(ierr);
8209baa881SHong Zhang   ierr = VecNorm(g3,NORM_2,&diffnorm);CHKERRQ(ierr);
8309baa881SHong Zhang   ierr = VecNorm(g3,NORM_INFINITY,&diffmax);CHKERRQ(ierr);
8409baa881SHong 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);
8509baa881SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  2-norm ||G - Gfd||/||G|| = %g, ||G - Gfd|| = %g\n",(double)(diffnorm/PetscMax(hcnorm,fdnorm)),(double)diffnorm);CHKERRQ(ierr);
8609baa881SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  max-norm ||G - Gfd||/||G|| = %g, ||G - Gfd|| = %g\n",(double)(diffmax/PetscMax(hcmax,fdmax)),(double)diffmax);CHKERRQ(ierr);
8709baa881SHong Zhang 
8809baa881SHong Zhang   if (complete_print) {
8909baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Hand-coded gradient ----------\n");CHKERRQ(ierr);
90913eda9aSHong Zhang     ierr = VecView(g1,mviewer);CHKERRQ(ierr);
9109baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Finite difference gradient ----------\n");CHKERRQ(ierr);
92913eda9aSHong Zhang     ierr = VecView(g2,mviewer);CHKERRQ(ierr);
9309baa881SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Hand-coded minus finite-difference gradient ----------\n");CHKERRQ(ierr);
94913eda9aSHong Zhang     ierr = VecView(g3,mviewer);CHKERRQ(ierr);
9509baa881SHong Zhang   }
9609baa881SHong Zhang   ierr = VecDestroy(&g2);CHKERRQ(ierr);
9709baa881SHong Zhang   ierr = VecDestroy(&g3);CHKERRQ(ierr);
98913eda9aSHong Zhang 
99913eda9aSHong Zhang   if (complete_print) {
100913eda9aSHong Zhang     ierr = PetscViewerPopFormat(mviewer);CHKERRQ(ierr);
101913eda9aSHong Zhang   }
102913eda9aSHong Zhang   ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
10309baa881SHong Zhang   PetscFunctionReturn(0);
10409baa881SHong Zhang }
10509baa881SHong Zhang 
106a7e14dcfSSatish Balay /*@
107a7e14dcfSSatish Balay   TaoComputeGradient - Computes the gradient of the objective function
108a7e14dcfSSatish Balay 
109441846f8SBarry Smith   Collective on Tao
110a7e14dcfSSatish Balay 
111a7e14dcfSSatish Balay   Input Parameters:
112441846f8SBarry Smith + tao - the Tao context
113a7e14dcfSSatish Balay - X - input vector
114a7e14dcfSSatish Balay 
115a7e14dcfSSatish Balay   Output Parameter:
116a7e14dcfSSatish Balay . G - gradient vector
117a7e14dcfSSatish Balay 
11809baa881SHong Zhang   Options Database Keys:
11909baa881SHong Zhang +    -tao_test_gradient - compare the user provided gradient with one compute via finite differences to check for errors
120dfe02fe6SHong 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
12109baa881SHong Zhang 
12295452b02SPatrick Sanan   Notes:
12395452b02SPatrick Sanan     TaoComputeGradient() is typically used within minimization implementations,
124a7e14dcfSSatish Balay   so most users would not generally call this routine themselves.
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay   Level: advanced
127a7e14dcfSSatish Balay 
128a7e14dcfSSatish Balay .seealso: TaoComputeObjective(), TaoComputeObjectiveAndGradient(), TaoSetGradientRoutine()
129a7e14dcfSSatish Balay @*/
130441846f8SBarry Smith PetscErrorCode TaoComputeGradient(Tao tao, Vec X, Vec G)
131a7e14dcfSSatish Balay {
132a7e14dcfSSatish Balay   PetscErrorCode ierr;
133a7e14dcfSSatish Balay   PetscReal      dummy;
13487f595a5SBarry Smith 
135a7e14dcfSSatish Balay   PetscFunctionBegin;
136441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
137a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
138a7e14dcfSSatish Balay   PetscValidHeaderSpecific(G,VEC_CLASSID,2);
139a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
140a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,G,3);
1418860a134SJunchao Zhang   ierr = VecLockReadPush(X);CHKERRQ(ierr);
142a7e14dcfSSatish Balay   if (tao->ops->computegradient) {
1430ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
144441846f8SBarry Smith     PetscStackPush("Tao user gradient evaluation routine");
145a7e14dcfSSatish Balay     ierr = (*tao->ops->computegradient)(tao,X,G,tao->user_gradP);CHKERRQ(ierr);
146a7e14dcfSSatish Balay     PetscStackPop;
1470ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
148a7e14dcfSSatish Balay     tao->ngrads++;
149a7e14dcfSSatish Balay   } else if (tao->ops->computeobjectiveandgradient) {
1500ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjGradEval,tao,X,G,NULL);CHKERRQ(ierr);
151a7e14dcfSSatish Balay     PetscStackPush("Tao user objective/gradient evaluation routine");
152a7e14dcfSSatish Balay     ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,&dummy,G,tao->user_objgradP);CHKERRQ(ierr);
153a7e14dcfSSatish Balay     PetscStackPop;
1540ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjGradEval,tao,X,G,NULL);CHKERRQ(ierr);
155a7e14dcfSSatish Balay     tao->nfuncgrads++;
15687f595a5SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetGradientRoutine() has not been called");
1578860a134SJunchao Zhang   ierr = VecLockReadPop(X);CHKERRQ(ierr);
15809baa881SHong Zhang 
159412cdd55SHong Zhang   ierr = TaoTestGradient(tao,X,G);CHKERRQ(ierr);
160a7e14dcfSSatish Balay   PetscFunctionReturn(0);
161a7e14dcfSSatish Balay }
162a7e14dcfSSatish Balay 
163a7e14dcfSSatish Balay /*@
164a7e14dcfSSatish Balay   TaoComputeObjective - Computes the objective function value at a given point
165a7e14dcfSSatish Balay 
166441846f8SBarry Smith   Collective on Tao
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay   Input Parameters:
169441846f8SBarry Smith + tao - the Tao context
170a7e14dcfSSatish Balay - X - input vector
171a7e14dcfSSatish Balay 
172a7e14dcfSSatish Balay   Output Parameter:
173a7e14dcfSSatish Balay . f - Objective value at X
174a7e14dcfSSatish Balay 
17595452b02SPatrick Sanan   Notes:
17695452b02SPatrick Sanan     TaoComputeObjective() is typically used within minimization implementations,
177a7e14dcfSSatish Balay   so most users would not generally call this routine themselves.
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay   Level: advanced
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
182a7e14dcfSSatish Balay @*/
183441846f8SBarry Smith PetscErrorCode TaoComputeObjective(Tao tao, Vec X, PetscReal *f)
184a7e14dcfSSatish Balay {
185a7e14dcfSSatish Balay   PetscErrorCode ierr;
186a7e14dcfSSatish Balay   Vec            temp;
18787f595a5SBarry Smith 
188a7e14dcfSSatish Balay   PetscFunctionBegin;
189441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
190a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
191a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
1928860a134SJunchao Zhang   ierr = VecLockReadPush(X);CHKERRQ(ierr);
193a7e14dcfSSatish Balay   if (tao->ops->computeobjective) {
1940ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
195441846f8SBarry Smith     PetscStackPush("Tao user objective evaluation routine");
196a7e14dcfSSatish Balay     ierr = (*tao->ops->computeobjective)(tao,X,f,tao->user_objP);CHKERRQ(ierr);
197a7e14dcfSSatish Balay     PetscStackPop;
1980ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
199a7e14dcfSSatish Balay     tao->nfuncs++;
200a7e14dcfSSatish Balay   } else if (tao->ops->computeobjectiveandgradient) {
201955c1f14SBarry Smith     ierr = PetscInfo(tao,"Duplicating variable vector in order to call func/grad routine\n");CHKERRQ(ierr);
202a7e14dcfSSatish Balay     ierr = VecDuplicate(X,&temp);CHKERRQ(ierr);
2030ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjGradEval,tao,X,NULL,NULL);CHKERRQ(ierr);
204441846f8SBarry Smith     PetscStackPush("Tao user objective/gradient evaluation routine");
205a7e14dcfSSatish Balay     ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,f,temp,tao->user_objgradP);CHKERRQ(ierr);
206a7e14dcfSSatish Balay     PetscStackPop;
2070ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjGradEval,tao,X,NULL,NULL);CHKERRQ(ierr);
208a7e14dcfSSatish Balay     ierr = VecDestroy(&temp);CHKERRQ(ierr);
209a7e14dcfSSatish Balay     tao->nfuncgrads++;
21087f595a5SBarry Smith   }  else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjectiveRoutine() has not been called");
211e356b196STodd Munson   ierr = PetscInfo1(tao,"TAO Function evaluation: %20.19e\n",(double)(*f));CHKERRQ(ierr);
2128860a134SJunchao Zhang   ierr = VecLockReadPop(X);CHKERRQ(ierr);
213a7e14dcfSSatish Balay   PetscFunctionReturn(0);
214a7e14dcfSSatish Balay }
215a7e14dcfSSatish Balay 
216a7e14dcfSSatish Balay /*@
217a7e14dcfSSatish Balay   TaoComputeObjectiveAndGradient - Computes the objective function value at a given point
218a7e14dcfSSatish Balay 
219441846f8SBarry Smith   Collective on Tao
220a7e14dcfSSatish Balay 
221a7e14dcfSSatish Balay   Input Parameters:
222441846f8SBarry Smith + tao - the Tao context
223a7e14dcfSSatish Balay - X - input vector
224a7e14dcfSSatish Balay 
225a7e14dcfSSatish Balay   Output Parameter:
226a7e14dcfSSatish Balay + f - Objective value at X
227a7e14dcfSSatish Balay - g - Gradient vector at X
228a7e14dcfSSatish Balay 
22995452b02SPatrick Sanan   Notes:
23095452b02SPatrick Sanan     TaoComputeObjectiveAndGradient() is typically used within minimization implementations,
231a7e14dcfSSatish Balay   so most users would not generally call this routine themselves.
232a7e14dcfSSatish Balay 
233a7e14dcfSSatish Balay   Level: advanced
234a7e14dcfSSatish Balay 
235a7e14dcfSSatish Balay .seealso: TaoComputeGradient(), TaoComputeObjectiveAndGradient(), TaoSetObjectiveRoutine()
236a7e14dcfSSatish Balay @*/
237441846f8SBarry Smith PetscErrorCode TaoComputeObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec G)
238a7e14dcfSSatish Balay {
239a7e14dcfSSatish Balay   PetscErrorCode ierr;
24087f595a5SBarry Smith 
241a7e14dcfSSatish Balay   PetscFunctionBegin;
242441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
243a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
244a7e14dcfSSatish Balay   PetscValidHeaderSpecific(G,VEC_CLASSID,4);
245a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
246a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,G,4);
2478860a134SJunchao Zhang   ierr = VecLockReadPush(X);CHKERRQ(ierr);
248a7e14dcfSSatish Balay   if (tao->ops->computeobjectiveandgradient) {
2490ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjGradEval,tao,X,G,NULL);CHKERRQ(ierr);
250f4c1ad5cSStefano Zampini     if (tao->ops->computegradient == TaoDefaultComputeGradient) {
251f4c1ad5cSStefano Zampini       ierr = TaoComputeObjective(tao,X,f);CHKERRQ(ierr);
252f4c1ad5cSStefano Zampini       ierr = TaoDefaultComputeGradient(tao,X,G,NULL);CHKERRQ(ierr);
253f4c1ad5cSStefano Zampini     } else {
254441846f8SBarry Smith       PetscStackPush("Tao user objective/gradient evaluation routine");
255a7e14dcfSSatish Balay       ierr = (*tao->ops->computeobjectiveandgradient)(tao,X,f,G,tao->user_objgradP);CHKERRQ(ierr);
2560cbffdbaSBarry Smith       PetscStackPop;
257a7e14dcfSSatish Balay     }
2580ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjGradEval,tao,X,G,NULL);CHKERRQ(ierr);
259a7e14dcfSSatish Balay     tao->nfuncgrads++;
260a7e14dcfSSatish Balay   } else if (tao->ops->computeobjective && tao->ops->computegradient) {
2610ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
262441846f8SBarry Smith     PetscStackPush("Tao user objective evaluation routine");
263a7e14dcfSSatish Balay     ierr = (*tao->ops->computeobjective)(tao,X,f,tao->user_objP);CHKERRQ(ierr);
264a7e14dcfSSatish Balay     PetscStackPop;
2650ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
266a7e14dcfSSatish Balay     tao->nfuncs++;
2670ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
268441846f8SBarry Smith     PetscStackPush("Tao user gradient evaluation routine");
269a7e14dcfSSatish Balay     ierr = (*tao->ops->computegradient)(tao,X,G,tao->user_gradP);CHKERRQ(ierr);
270a7e14dcfSSatish Balay     PetscStackPop;
2710ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_GradientEval,tao,X,G,NULL);CHKERRQ(ierr);
272a7e14dcfSSatish Balay     tao->ngrads++;
27387f595a5SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetObjectiveRoutine() or TaoSetGradientRoutine() not set");
274e356b196STodd Munson   ierr = PetscInfo1(tao,"TAO Function evaluation: %20.19e\n",(double)(*f));CHKERRQ(ierr);
2758860a134SJunchao Zhang   ierr = VecLockReadPop(X);CHKERRQ(ierr);
2761657496cSHong Zhang 
277412cdd55SHong Zhang   ierr = TaoTestGradient(tao,X,G);CHKERRQ(ierr);
278a7e14dcfSSatish Balay   PetscFunctionReturn(0);
279a7e14dcfSSatish Balay }
280a7e14dcfSSatish Balay 
281a7e14dcfSSatish Balay /*@C
282a7e14dcfSSatish Balay   TaoSetObjectiveRoutine - Sets the function evaluation routine for minimization
283a7e14dcfSSatish Balay 
284441846f8SBarry Smith   Logically collective on Tao
285a7e14dcfSSatish Balay 
286a7e14dcfSSatish Balay   Input Parameter:
287441846f8SBarry Smith + tao - the Tao context
288a7e14dcfSSatish Balay . func - the objective function
289a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the function evaluation
2906c23d075SBarry Smith         routine (may be NULL)
291a7e14dcfSSatish Balay 
292a7e14dcfSSatish Balay   Calling sequence of func:
293441846f8SBarry Smith $      func (Tao tao, Vec x, PetscReal *f, void *ctx);
294a7e14dcfSSatish Balay 
295a7e14dcfSSatish Balay + x - input vector
296a7e14dcfSSatish Balay . f - function value
297a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
298a7e14dcfSSatish Balay 
299a7e14dcfSSatish Balay   Level: beginner
300a7e14dcfSSatish Balay 
301a7e14dcfSSatish Balay .seealso: TaoSetGradientRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
302a7e14dcfSSatish Balay @*/
303441846f8SBarry Smith PetscErrorCode TaoSetObjectiveRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*,void*),void *ctx)
304a7e14dcfSSatish Balay {
305a7e14dcfSSatish Balay   PetscFunctionBegin;
306441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
307a7e14dcfSSatish Balay   tao->user_objP = ctx;
308a7e14dcfSSatish Balay   tao->ops->computeobjective = func;
309a7e14dcfSSatish Balay   PetscFunctionReturn(0);
310a7e14dcfSSatish Balay }
311a7e14dcfSSatish Balay 
312a7e14dcfSSatish Balay /*@C
3134a48860cSAlp Dener   TaoSetResidualRoutine - Sets the residual evaluation routine for least-square applications
314a7e14dcfSSatish Balay 
315441846f8SBarry Smith   Logically collective on Tao
316a7e14dcfSSatish Balay 
317a7e14dcfSSatish Balay   Input Parameter:
318441846f8SBarry Smith + tao - the Tao context
3194a48860cSAlp Dener . func - the residual evaluation routine
320a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the function evaluation
3216c23d075SBarry Smith         routine (may be NULL)
322a7e14dcfSSatish Balay 
323a7e14dcfSSatish Balay   Calling sequence of func:
324441846f8SBarry Smith $      func (Tao tao, Vec x, Vec f, void *ctx);
325a7e14dcfSSatish Balay 
326a7e14dcfSSatish Balay + x - input vector
327a7e14dcfSSatish Balay . f - function value vector
328a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
329a7e14dcfSSatish Balay 
330a7e14dcfSSatish Balay   Level: beginner
331a7e14dcfSSatish Balay 
332a7e14dcfSSatish Balay .seealso: TaoSetObjectiveRoutine(), TaoSetJacobianRoutine()
333a7e14dcfSSatish Balay @*/
3344a48860cSAlp Dener PetscErrorCode TaoSetResidualRoutine(Tao tao, Vec res, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
335a7e14dcfSSatish Balay {
336737f463aSAlp Dener   PetscErrorCode ierr;
337737f463aSAlp Dener 
338a7e14dcfSSatish Balay   PetscFunctionBegin;
339441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
3404a48860cSAlp Dener   PetscValidHeaderSpecific(res,VEC_CLASSID,2);
3414ffbe8acSAlp Dener   ierr = PetscObjectReference((PetscObject)res);CHKERRQ(ierr);
342737f463aSAlp Dener   if (tao->ls_res) {
343737f463aSAlp Dener     ierr = VecDestroy(&tao->ls_res);CHKERRQ(ierr);
344737f463aSAlp Dener   }
3454a48860cSAlp Dener   tao->ls_res = res;
3464ffbe8acSAlp Dener   tao->user_lsresP = ctx;
3474a48860cSAlp Dener   tao->ops->computeresidual = func;
348737f463aSAlp Dener 
349737f463aSAlp Dener   PetscFunctionReturn(0);
350737f463aSAlp Dener }
351737f463aSAlp Dener 
352737f463aSAlp Dener /*@
353737f463aSAlp 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.
354737f463aSAlp Dener 
355737f463aSAlp Dener   Collective on Tao
356737f463aSAlp Dener 
357737f463aSAlp Dener   Input Parameters:
358737f463aSAlp Dener + tao - the Tao context
359737f463aSAlp Dener . sigma_v - vector of weights (diagonal terms only)
360737f463aSAlp Dener . n       - the number of weights (if using off-diagonal)
361737f463aSAlp Dener . rows    - index list of rows for sigma_w
362737f463aSAlp Dener . cols    - index list of columns for sigma_w
363737f463aSAlp Dener - vals - array of weights
364737f463aSAlp Dener 
365737f463aSAlp Dener 
366737f463aSAlp Dener 
367737f463aSAlp Dener   Note: Either sigma_v or sigma_w (or both) should be NULL
368737f463aSAlp Dener 
369737f463aSAlp Dener   Level: intermediate
370737f463aSAlp Dener 
371737f463aSAlp Dener .seealso: TaoSetResidualRoutine()
372737f463aSAlp Dener @*/
373737f463aSAlp Dener PetscErrorCode TaoSetResidualWeights(Tao tao, Vec sigma_v, PetscInt n, PetscInt *rows, PetscInt *cols, PetscReal *vals)
374737f463aSAlp Dener {
375737f463aSAlp Dener   PetscErrorCode ierr;
376737f463aSAlp Dener   PetscInt       i;
377737f463aSAlp Dener   PetscFunctionBegin;
378737f463aSAlp Dener   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
379737f463aSAlp Dener   if (sigma_v) {
3804ffbe8acSAlp Dener     PetscValidHeaderSpecific(sigma_v,VEC_CLASSID,2);
381737f463aSAlp Dener     ierr = PetscObjectReference((PetscObject)sigma_v);CHKERRQ(ierr);
382737f463aSAlp Dener   }
3834ffbe8acSAlp Dener   if (tao->res_weights_v) {
3844ffbe8acSAlp Dener     ierr = VecDestroy(&tao->res_weights_v);CHKERRQ(ierr);
3854ffbe8acSAlp Dener   }
3864ffbe8acSAlp Dener   tao->res_weights_v=sigma_v;
387737f463aSAlp Dener   if (vals) {
388737f463aSAlp Dener     if (tao->res_weights_n) {
389737f463aSAlp Dener       ierr = PetscFree(tao->res_weights_rows);CHKERRQ(ierr);
390737f463aSAlp Dener       ierr = PetscFree(tao->res_weights_cols);CHKERRQ(ierr);
391737f463aSAlp Dener       ierr = PetscFree(tao->res_weights_w);CHKERRQ(ierr);
392737f463aSAlp Dener     }
393737f463aSAlp Dener     ierr = PetscMalloc1(n,&tao->res_weights_rows);CHKERRQ(ierr);
394737f463aSAlp Dener     ierr = PetscMalloc1(n,&tao->res_weights_cols);CHKERRQ(ierr);
395737f463aSAlp Dener     ierr = PetscMalloc1(n,&tao->res_weights_w);CHKERRQ(ierr);
396737f463aSAlp Dener     tao->res_weights_n=n;
397737f463aSAlp Dener     for (i=0;i<n;i++) {
398737f463aSAlp Dener       tao->res_weights_rows[i]=rows[i];
399737f463aSAlp Dener       tao->res_weights_cols[i]=cols[i];
400737f463aSAlp Dener       tao->res_weights_w[i]=vals[i];
401737f463aSAlp Dener     }
402737f463aSAlp Dener   } else {
403737f463aSAlp Dener     tao->res_weights_n=0;
404737f463aSAlp Dener     tao->res_weights_rows=0;
405737f463aSAlp Dener     tao->res_weights_cols=0;
406737f463aSAlp Dener   }
407a7e14dcfSSatish Balay   PetscFunctionReturn(0);
408a7e14dcfSSatish Balay }
409a7e14dcfSSatish Balay 
4108b7a9b22SJason Sarich /*@
4114a48860cSAlp Dener   TaoComputeResidual - Computes a least-squares residual vector at a given point
412a7e14dcfSSatish Balay 
413441846f8SBarry Smith   Collective on Tao
414a7e14dcfSSatish Balay 
415a7e14dcfSSatish Balay   Input Parameters:
416441846f8SBarry Smith + tao - the Tao context
417a7e14dcfSSatish Balay - X - input vector
418a7e14dcfSSatish Balay 
419a7e14dcfSSatish Balay   Output Parameter:
420a7e14dcfSSatish Balay . f - Objective vector at X
421a7e14dcfSSatish Balay 
42295452b02SPatrick Sanan   Notes:
4234a48860cSAlp Dener     TaoComputeResidual() is typically used within minimization implementations,
424a7e14dcfSSatish Balay   so most users would not generally call this routine themselves.
425a7e14dcfSSatish Balay 
426a7e14dcfSSatish Balay   Level: advanced
427a7e14dcfSSatish Balay 
4284a48860cSAlp Dener .seealso: TaoSetResidualRoutine()
429a7e14dcfSSatish Balay @*/
4304a48860cSAlp Dener PetscErrorCode TaoComputeResidual(Tao tao, Vec X, Vec F)
431a7e14dcfSSatish Balay {
432a7e14dcfSSatish Balay   PetscErrorCode ierr;
43387f595a5SBarry Smith 
434a7e14dcfSSatish Balay   PetscFunctionBegin;
435441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
436a7e14dcfSSatish Balay   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
437a7e14dcfSSatish Balay   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
438a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,X,2);
439a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,F,3);
4404a48860cSAlp Dener   if (tao->ops->computeresidual) {
4410ebee16dSLisandro Dalcin     ierr = PetscLogEventBegin(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
4424a48860cSAlp Dener     PetscStackPush("Tao user least-squares residual evaluation routine");
4434a48860cSAlp Dener     ierr = (*tao->ops->computeresidual)(tao,X,F,tao->user_lsresP);CHKERRQ(ierr);
444a7e14dcfSSatish Balay     PetscStackPop;
4450ebee16dSLisandro Dalcin     ierr = PetscLogEventEnd(TAO_ObjectiveEval,tao,X,NULL,NULL);CHKERRQ(ierr);
446a7e14dcfSSatish Balay     tao->nfuncs++;
4474a48860cSAlp Dener   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"TaoSetResidualRoutine() has not been called");
4484a48860cSAlp Dener   ierr = PetscInfo(tao,"TAO least-squares residual evaluation.\n");CHKERRQ(ierr);
449a7e14dcfSSatish Balay   PetscFunctionReturn(0);
450a7e14dcfSSatish Balay }
451a7e14dcfSSatish Balay 
452a7e14dcfSSatish Balay /*@C
453a7e14dcfSSatish Balay   TaoSetGradientRoutine - Sets the gradient evaluation routine for minimization
454a7e14dcfSSatish Balay 
455441846f8SBarry Smith   Logically collective on Tao
456a7e14dcfSSatish Balay 
457a7e14dcfSSatish Balay   Input Parameter:
458441846f8SBarry Smith + tao - the Tao context
459a7e14dcfSSatish Balay . func - the gradient function
460a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the gradient evaluation
4616c23d075SBarry Smith         routine (may be NULL)
462a7e14dcfSSatish Balay 
463a7e14dcfSSatish Balay   Calling sequence of func:
464441846f8SBarry Smith $      func (Tao tao, Vec x, Vec g, void *ctx);
465a7e14dcfSSatish Balay 
466a7e14dcfSSatish Balay + x - input vector
467a7e14dcfSSatish Balay . g - gradient value (output)
468a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
469a7e14dcfSSatish Balay 
470a7e14dcfSSatish Balay   Level: beginner
471a7e14dcfSSatish Balay 
472a7e14dcfSSatish Balay .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
473a7e14dcfSSatish Balay @*/
474441846f8SBarry Smith PetscErrorCode TaoSetGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, Vec, void*),void *ctx)
475a7e14dcfSSatish Balay {
476a7e14dcfSSatish Balay   PetscFunctionBegin;
477441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
478a7e14dcfSSatish Balay   tao->user_gradP = ctx;
479a7e14dcfSSatish Balay   tao->ops->computegradient = func;
480a7e14dcfSSatish Balay   PetscFunctionReturn(0);
481a7e14dcfSSatish Balay }
482a7e14dcfSSatish Balay 
483a7e14dcfSSatish Balay /*@C
484a7e14dcfSSatish Balay   TaoSetObjectiveAndGradientRoutine - Sets a combined objective function and gradient evaluation routine for minimization
485a7e14dcfSSatish Balay 
486441846f8SBarry Smith   Logically collective on Tao
487a7e14dcfSSatish Balay 
488a7e14dcfSSatish Balay   Input Parameter:
489441846f8SBarry Smith + tao - the Tao context
490a7e14dcfSSatish Balay . func - the gradient function
491a7e14dcfSSatish Balay - ctx - [optional] user-defined context for private data for the gradient evaluation
4926c23d075SBarry Smith         routine (may be NULL)
493a7e14dcfSSatish Balay 
494a7e14dcfSSatish Balay   Calling sequence of func:
49517477c02SJason Sarich $      func (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx);
496a7e14dcfSSatish Balay 
497a7e14dcfSSatish Balay + x - input vector
49817477c02SJason Sarich . f - objective value (output)
499a7e14dcfSSatish Balay . g - gradient value (output)
500a7e14dcfSSatish Balay - ctx - [optional] user-defined function context
501a7e14dcfSSatish Balay 
502a7e14dcfSSatish Balay   Level: beginner
503a7e14dcfSSatish Balay 
504a7e14dcfSSatish Balay .seealso: TaoSetObjectiveRoutine(), TaoSetHessianRoutine() TaoSetObjectiveAndGradientRoutine()
505a7e14dcfSSatish Balay @*/
506441846f8SBarry Smith PetscErrorCode TaoSetObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal *, Vec, void*), void *ctx)
507a7e14dcfSSatish Balay {
508a7e14dcfSSatish Balay   PetscFunctionBegin;
509441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
510a7e14dcfSSatish Balay   tao->user_objgradP = ctx;
511a7e14dcfSSatish Balay   tao->ops->computeobjectiveandgradient = func;
512a7e14dcfSSatish Balay   PetscFunctionReturn(0);
513a7e14dcfSSatish Balay }
514a7e14dcfSSatish Balay 
515a7e14dcfSSatish Balay /*@
516a7e14dcfSSatish Balay   TaoIsObjectiveDefined -- Checks to see if the user has
517a7e14dcfSSatish Balay   declared an objective-only routine.  Useful for determining when
518a7e14dcfSSatish Balay   it is appropriate to call TaoComputeObjective() or
519a7e14dcfSSatish Balay   TaoComputeObjectiveAndGradient()
520a7e14dcfSSatish Balay 
521441846f8SBarry Smith   Collective on Tao
522a7e14dcfSSatish Balay 
523a7e14dcfSSatish Balay   Input Parameter:
524441846f8SBarry Smith + tao - the Tao context
525a7e14dcfSSatish Balay - ctx - PETSC_TRUE if objective function routine is set by user,
526a7e14dcfSSatish Balay         PETSC_FALSE otherwise
527a7e14dcfSSatish Balay   Level: developer
528a7e14dcfSSatish Balay 
529a7e14dcfSSatish Balay .seealso: TaoSetObjectiveRoutine(), TaoIsGradientDefined(), TaoIsObjectiveAndGradientDefined()
530a7e14dcfSSatish Balay @*/
531441846f8SBarry Smith PetscErrorCode TaoIsObjectiveDefined(Tao tao, PetscBool *flg)
532a7e14dcfSSatish Balay {
533a7e14dcfSSatish Balay   PetscFunctionBegin;
534441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
53545cf516eSBarry Smith   if (tao->ops->computeobjective == 0) *flg = PETSC_FALSE;
53645cf516eSBarry Smith   else *flg = PETSC_TRUE;
537a7e14dcfSSatish Balay   PetscFunctionReturn(0);
538a7e14dcfSSatish Balay }
539a7e14dcfSSatish Balay 
540a7e14dcfSSatish Balay /*@
541a7e14dcfSSatish Balay   TaoIsGradientDefined -- Checks to see if the user has
542a7e14dcfSSatish Balay   declared an objective-only routine.  Useful for determining when
543a7e14dcfSSatish Balay   it is appropriate to call TaoComputeGradient() or
544a7e14dcfSSatish Balay   TaoComputeGradientAndGradient()
545a7e14dcfSSatish Balay 
546a7e14dcfSSatish Balay   Not Collective
547a7e14dcfSSatish Balay 
548a7e14dcfSSatish Balay   Input Parameter:
549441846f8SBarry Smith + tao - the Tao context
550a7e14dcfSSatish Balay - ctx - PETSC_TRUE if gradient routine is set by user, PETSC_FALSE otherwise
551a7e14dcfSSatish Balay   Level: developer
552a7e14dcfSSatish Balay 
553a7e14dcfSSatish Balay .seealso: TaoSetGradientRoutine(), TaoIsObjectiveDefined(), TaoIsObjectiveAndGradientDefined()
554a7e14dcfSSatish Balay @*/
555441846f8SBarry Smith PetscErrorCode TaoIsGradientDefined(Tao tao, PetscBool *flg)
556a7e14dcfSSatish Balay {
557a7e14dcfSSatish Balay   PetscFunctionBegin;
558441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
55945cf516eSBarry Smith   if (tao->ops->computegradient == 0) *flg = PETSC_FALSE;
56045cf516eSBarry Smith   else *flg = PETSC_TRUE;
561a7e14dcfSSatish Balay   PetscFunctionReturn(0);
562a7e14dcfSSatish Balay }
563a7e14dcfSSatish Balay 
564a7e14dcfSSatish Balay /*@
565a7e14dcfSSatish Balay   TaoIsObjectiveAndGradientDefined -- Checks to see if the user has
566a7e14dcfSSatish Balay   declared a joint objective/gradient routine.  Useful for determining when
567a7e14dcfSSatish Balay   it is appropriate to call TaoComputeObjective() or
568a7e14dcfSSatish Balay   TaoComputeObjectiveAndGradient()
569a7e14dcfSSatish Balay 
570a7e14dcfSSatish Balay   Not Collective
571a7e14dcfSSatish Balay 
572a7e14dcfSSatish Balay   Input Parameter:
573441846f8SBarry Smith + tao - the Tao context
574a7e14dcfSSatish Balay - ctx - PETSC_TRUE if objective/gradient routine is set by user, PETSC_FALSE otherwise
575a7e14dcfSSatish Balay   Level: developer
576a7e14dcfSSatish Balay 
577a7e14dcfSSatish Balay .seealso: TaoSetObjectiveAndGradientRoutine(), TaoIsObjectiveDefined(), TaoIsGradientDefined()
578a7e14dcfSSatish Balay @*/
579441846f8SBarry Smith PetscErrorCode TaoIsObjectiveAndGradientDefined(Tao tao, PetscBool *flg)
580a7e14dcfSSatish Balay {
581a7e14dcfSSatish Balay   PetscFunctionBegin;
582441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
58345cf516eSBarry Smith   if (tao->ops->computeobjectiveandgradient == 0) *flg = PETSC_FALSE;
58445cf516eSBarry Smith   else *flg = PETSC_TRUE;
585a7e14dcfSSatish Balay   PetscFunctionReturn(0);
586a7e14dcfSSatish Balay }
587