xref: /petsc/src/tao/leastsquares/tutorials/chwirut1.c (revision 9fe8556826ed0dc8a9fe15a3f2d3eda474dc9338)
1 /*
2    Include "petsctao.h" so that we can use TAO solvers.  Note that this
3    file automatically includes libraries such as:
4      petsc.h       - base PETSc routines   petscvec.h - vectors
5      petscsys.h    - system routines        petscmat.h - matrices
6      petscis.h     - index sets            petscksp.h - Krylov subspace methods
7      petscviewer.h - viewers               petscpc.h  - preconditioners
8 
9 */
10 
11 #include <petsctao.h>
12 
13 /*
14 Description:   These data are the result of a NIST study involving
15                ultrasonic calibration.  The response variable is
16                ultrasonic response, and the predictor variable is
17                metal distance.
18 
19 Reference:     Chwirut, D., NIST (197?).
20                Ultrasonic Reference Block Study.
21 */
22 
23 static char help[]="Finds the nonlinear least-squares solution to the model \n\
24             y = exp[-b1*x]/(b2+b3*x)  +  e \n";
25 
26 #define NOBSERVATIONS 214
27 #define NPARAMETERS 3
28 
29 /* User-defined application context */
30 typedef struct {
31   /* Working space */
32   PetscReal t[NOBSERVATIONS];   /* array of independent variables of observation */
33   PetscReal y[NOBSERVATIONS];   /* array of dependent variables */
34   PetscReal j[NOBSERVATIONS][NPARAMETERS]; /* dense jacobian matrix array*/
35   PetscInt idm[NOBSERVATIONS];  /* Matrix indices for jacobian */
36   PetscInt idn[NPARAMETERS];
37 } AppCtx;
38 
39 /* User provided Routines */
40 PetscErrorCode InitializeData(AppCtx *user);
41 PetscErrorCode FormStartingPoint(Vec);
42 PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *);
43 PetscErrorCode EvaluateJacobian(Tao, Vec, Mat, Mat, void *);
44 
45 /*--------------------------------------------------------------------*/
46 int main(int argc,char **argv)
47 {
48   Vec            x, f;               /* solution, function */
49   Mat            J;                  /* Jacobian matrix */
50   Tao            tao;                /* Tao solver context */
51   PetscInt       i;               /* iteration information */
52   PetscReal      hist[100],resid[100];
53   PetscInt       lits[100];
54   AppCtx         user;               /* user-defined work context */
55 
56   PetscFunctionBeginUser;
57   PetscCall(PetscInitialize(&argc,&argv,(char *)0,help));
58   /* Allocate vectors */
59   PetscCall(VecCreateSeq(MPI_COMM_SELF,NPARAMETERS,&x));
60   PetscCall(VecCreateSeq(MPI_COMM_SELF,NOBSERVATIONS,&f));
61 
62   /* Create the Jacobian matrix. */
63   PetscCall(MatCreateSeqDense(MPI_COMM_SELF,NOBSERVATIONS,NPARAMETERS,NULL,&J));
64 
65   for (i=0;i<NOBSERVATIONS;i++) user.idm[i] = i;
66 
67   for (i=0;i<NPARAMETERS;i++) user.idn[i] = i;
68 
69   /* Create TAO solver and set desired solution method */
70   PetscCall(TaoCreate(PETSC_COMM_SELF,&tao));
71   PetscCall(TaoSetType(tao,TAOPOUNDERS));
72 
73  /* Set the function and Jacobian routines. */
74   PetscCall(InitializeData(&user));
75   PetscCall(FormStartingPoint(x));
76   PetscCall(TaoSetSolution(tao,x));
77   PetscCall(TaoSetResidualRoutine(tao,f,EvaluateFunction,(void*)&user));
78   PetscCall(TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void*)&user));
79 
80   /* Check for any TAO command line arguments */
81   PetscCall(TaoSetFromOptions(tao));
82 
83   PetscCall(TaoSetConvergenceHistory(tao,hist,resid,0,lits,100,PETSC_TRUE));
84   /* Perform the Solve */
85   PetscCall(TaoSolve(tao));
86 
87   /* View the vector; then destroy it.  */
88   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_SELF));
89 
90   /* Free TAO data structures */
91   PetscCall(TaoDestroy(&tao));
92 
93    /* Free PETSc data structures */
94   PetscCall(VecDestroy(&x));
95   PetscCall(VecDestroy(&f));
96   PetscCall(MatDestroy(&J));
97 
98   PetscCall(PetscFinalize());
99   return 0;
100 }
101 
102 /*--------------------------------------------------------------------*/
103 PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr)
104 {
105   AppCtx         *user = (AppCtx *)ptr;
106   PetscInt       i;
107   const PetscReal *x;
108   PetscReal      *y=user->y,*f,*t=user->t;
109 
110   PetscFunctionBegin;
111   PetscCall(VecGetArrayRead(X,&x));
112   PetscCall(VecGetArray(F,&f));
113 
114   for (i=0;i<NOBSERVATIONS;i++) {
115     f[i] = y[i] - PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]);
116   }
117   PetscCall(VecRestoreArrayRead(X,&x));
118   PetscCall(VecRestoreArray(F,&f));
119   PetscLogFlops(6*NOBSERVATIONS);
120   PetscFunctionReturn(0);
121 }
122 
123 /*------------------------------------------------------------*/
124 /* J[i][j] = df[i]/dt[j] */
125 PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr)
126 {
127   AppCtx         *user = (AppCtx *)ptr;
128   PetscInt       i;
129   const PetscReal *x;
130   PetscReal      *t=user->t;
131   PetscReal      base;
132 
133   PetscFunctionBegin;
134   PetscCall(VecGetArrayRead(X,&x));
135   for (i=0;i<NOBSERVATIONS;i++) {
136     base = PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]);
137 
138     user->j[i][0] = t[i]*base;
139     user->j[i][1] = base/(x[1] + x[2]*t[i]);
140     user->j[i][2] = base*t[i]/(x[1] + x[2]*t[i]);
141   }
142 
143   /* Assemble the matrix */
144   PetscCall(MatSetValues(J,NOBSERVATIONS,user->idm, NPARAMETERS, user->idn,(PetscReal *)user->j,INSERT_VALUES));
145   PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
146   PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
147 
148   PetscCall(VecRestoreArrayRead(X,&x));
149   PetscLogFlops(NOBSERVATIONS * 13);
150   PetscFunctionReturn(0);
151 }
152 
153 /* ------------------------------------------------------------ */
154 PetscErrorCode FormStartingPoint(Vec X)
155 {
156   PetscReal      *x;
157 
158   PetscFunctionBegin;
159   PetscCall(VecGetArray(X,&x));
160   x[0] = 0.15;
161   x[1] = 0.008;
162   x[2] = 0.010;
163   PetscCall(VecRestoreArray(X,&x));
164   PetscFunctionReturn(0);
165 }
166 
167 /* ---------------------------------------------------------------------- */
168 PetscErrorCode InitializeData(AppCtx *user)
169 {
170   PetscReal *t=user->t,*y=user->y;
171   PetscInt  i=0;
172 
173   PetscFunctionBegin;
174   y[i] =   92.9000;   t[i++] =  0.5000;
175   y[i] =    78.7000;  t[i++] =   0.6250;
176   y[i] =    64.2000;  t[i++] =   0.7500;
177   y[i] =    64.9000;  t[i++] =   0.8750;
178   y[i] =    57.1000;  t[i++] =   1.0000;
179   y[i] =    43.3000;  t[i++] =   1.2500;
180   y[i] =    31.1000;   t[i++] =  1.7500;
181   y[i] =    23.6000;   t[i++] =  2.2500;
182   y[i] =    31.0500;   t[i++] =  1.7500;
183   y[i] =    23.7750;   t[i++] =  2.2500;
184   y[i] =    17.7375;   t[i++] =  2.7500;
185   y[i] =    13.8000;   t[i++] =  3.2500;
186   y[i] =    11.5875;   t[i++] =  3.7500;
187   y[i] =     9.4125;   t[i++] =  4.2500;
188   y[i] =     7.7250;   t[i++] =  4.7500;
189   y[i] =     7.3500;   t[i++] =  5.2500;
190   y[i] =     8.0250;   t[i++] =  5.7500;
191   y[i] =    90.6000;   t[i++] =  0.5000;
192   y[i] =    76.9000;   t[i++] =  0.6250;
193   y[i] =    71.6000;   t[i++] = 0.7500;
194   y[i] =    63.6000;   t[i++] =  0.8750;
195   y[i] =    54.0000;   t[i++] =  1.0000;
196   y[i] =    39.2000;   t[i++] =  1.2500;
197   y[i] =    29.3000;   t[i++] = 1.7500;
198   y[i] =    21.4000;   t[i++] =  2.2500;
199   y[i] =    29.1750;   t[i++] =  1.7500;
200   y[i] =    22.1250;   t[i++] =  2.2500;
201   y[i] =    17.5125;   t[i++] =  2.7500;
202   y[i] =    14.2500;   t[i++] =  3.2500;
203   y[i] =     9.4500;   t[i++] =  3.7500;
204   y[i] =     9.1500;   t[i++] =  4.2500;
205   y[i] =     7.9125;   t[i++] =  4.7500;
206   y[i] =     8.4750;   t[i++] =  5.2500;
207   y[i] =     6.1125;   t[i++] =  5.7500;
208   y[i] =    80.0000;   t[i++] =  0.5000;
209   y[i] =    79.0000;   t[i++] =  0.6250;
210   y[i] =    63.8000;   t[i++] =  0.7500;
211   y[i] =    57.2000;   t[i++] =  0.8750;
212   y[i] =    53.2000;   t[i++] =  1.0000;
213   y[i] =   42.5000;   t[i++] =  1.2500;
214   y[i] =   26.8000;   t[i++] =  1.7500;
215   y[i] =    20.4000;   t[i++] =  2.2500;
216   y[i] =    26.8500;  t[i++] =   1.7500;
217   y[i] =    21.0000;  t[i++] =   2.2500;
218   y[i] =    16.4625;  t[i++] =   2.7500;
219   y[i] =    12.5250;  t[i++] =   3.2500;
220   y[i] =    10.5375;  t[i++] =   3.7500;
221   y[i] =     8.5875;  t[i++] =   4.2500;
222   y[i] =     7.1250;  t[i++] =   4.7500;
223   y[i] =     6.1125;  t[i++] =   5.2500;
224   y[i] =     5.9625;  t[i++] =   5.7500;
225   y[i] =    74.1000;  t[i++] =   0.5000;
226   y[i] =    67.3000;  t[i++] =   0.6250;
227   y[i] =    60.8000;  t[i++] =   0.7500;
228   y[i] =    55.5000;  t[i++] =   0.8750;
229   y[i] =    50.3000;  t[i++] =   1.0000;
230   y[i] =    41.0000;  t[i++] =   1.2500;
231   y[i] =    29.4000;  t[i++] =   1.7500;
232   y[i] =    20.4000;  t[i++] =   2.2500;
233   y[i] =    29.3625;  t[i++] =   1.7500;
234   y[i] =    21.1500;  t[i++] =   2.2500;
235   y[i] =    16.7625;  t[i++] =   2.7500;
236   y[i] =    13.2000;  t[i++] =   3.2500;
237   y[i] =    10.8750;  t[i++] =   3.7500;
238   y[i] =     8.1750;  t[i++] =   4.2500;
239   y[i] =     7.3500;  t[i++] =   4.7500;
240   y[i] =     5.9625;  t[i++] =  5.2500;
241   y[i] =     5.6250;  t[i++] =   5.7500;
242   y[i] =    81.5000;  t[i++] =    .5000;
243   y[i] =    62.4000;  t[i++] =    .7500;
244   y[i] =    32.5000;  t[i++] =   1.5000;
245   y[i] =    12.4100;  t[i++] =   3.0000;
246   y[i] =    13.1200;  t[i++] =   3.0000;
247   y[i] =    15.5600;  t[i++] =   3.0000;
248   y[i] =     5.6300;  t[i++] =   6.0000;
249   y[i] =    78.0000;   t[i++] =   .5000;
250   y[i] =    59.9000;  t[i++] =    .7500;
251   y[i] =    33.2000;  t[i++] =   1.5000;
252   y[i] =    13.8400;  t[i++] =   3.0000;
253   y[i] =    12.7500;  t[i++] =   3.0000;
254   y[i] =    14.6200;  t[i++] =   3.0000;
255   y[i] =     3.9400;  t[i++] =   6.0000;
256   y[i] =    76.8000;  t[i++] =    .5000;
257   y[i] =    61.0000;  t[i++] =    .7500;
258   y[i] =    32.9000;  t[i++] =   1.5000;
259   y[i] =   13.8700;   t[i++] = 3.0000;
260   y[i] =    11.8100;  t[i++] =   3.0000;
261   y[i] =    13.3100;  t[i++] =   3.0000;
262   y[i] =     5.4400;  t[i++] =   6.0000;
263   y[i] =    78.0000;  t[i++] =    .5000;
264   y[i] =    63.5000;  t[i++] =    .7500;
265   y[i] =    33.8000;  t[i++] =   1.5000;
266   y[i] =    12.5600;  t[i++] =   3.0000;
267   y[i] =     5.6300;  t[i++] =   6.0000;
268   y[i] =    12.7500;  t[i++] =   3.0000;
269   y[i] =    13.1200;  t[i++] =   3.0000;
270   y[i] =     5.4400;  t[i++] =   6.0000;
271   y[i] =    76.8000;  t[i++] =    .5000;
272   y[i] =    60.0000;  t[i++] =    .7500;
273   y[i] =    47.8000;  t[i++] =   1.0000;
274   y[i] =    32.0000;  t[i++] =   1.5000;
275   y[i] =    22.2000;  t[i++] =   2.0000;
276   y[i] =    22.5700;  t[i++] =   2.0000;
277   y[i] =    18.8200;  t[i++] =   2.5000;
278   y[i] =    13.9500;  t[i++] =   3.0000;
279   y[i] =    11.2500;  t[i++] =   4.0000;
280   y[i] =     9.0000;  t[i++] =   5.0000;
281   y[i] =     6.6700;  t[i++] =   6.0000;
282   y[i] =    75.8000;  t[i++] =    .5000;
283   y[i] =    62.0000;  t[i++] =    .7500;
284   y[i] =    48.8000;  t[i++] =   1.0000;
285   y[i] =    35.2000;  t[i++] =   1.5000;
286   y[i] =    20.0000;  t[i++] =   2.0000;
287   y[i] =    20.3200;  t[i++] =   2.0000;
288   y[i] =    19.3100;  t[i++] =   2.5000;
289   y[i] =    12.7500;  t[i++] =   3.0000;
290   y[i] =    10.4200;  t[i++] =   4.0000;
291   y[i] =     7.3100;  t[i++] =   5.0000;
292   y[i] =     7.4200;  t[i++] =   6.0000;
293   y[i] =    70.5000;  t[i++] =    .5000;
294   y[i] =    59.5000;  t[i++] =    .7500;
295   y[i] =    48.5000;  t[i++] =   1.0000;
296   y[i] =    35.8000;  t[i++] =   1.5000;
297   y[i] =    21.0000;  t[i++] =   2.0000;
298   y[i] =    21.6700;  t[i++] =   2.0000;
299   y[i] =    21.0000;  t[i++] =   2.5000;
300   y[i] =    15.6400;  t[i++] =   3.0000;
301   y[i] =     8.1700;  t[i++] =   4.0000;
302   y[i] =     8.5500;  t[i++] =   5.0000;
303   y[i] =    10.1200;  t[i++] =   6.0000;
304   y[i] =    78.0000;  t[i++] =    .5000;
305   y[i] =    66.0000;  t[i++] =    .6250;
306   y[i] =    62.0000;  t[i++] =    .7500;
307   y[i] =    58.0000;  t[i++] =    .8750;
308   y[i] =    47.7000;  t[i++] =   1.0000;
309   y[i] =    37.8000;  t[i++] =   1.2500;
310   y[i] =    20.2000;  t[i++] =   2.2500;
311   y[i] =    21.0700;  t[i++] =   2.2500;
312   y[i] =    13.8700;  t[i++] =   2.7500;
313   y[i] =     9.6700;  t[i++] =   3.2500;
314   y[i] =     7.7600;  t[i++] =   3.7500;
315   y[i] =    5.4400;   t[i++] =  4.2500;
316   y[i] =    4.8700;   t[i++] =  4.7500;
317   y[i] =     4.0100;  t[i++] =   5.2500;
318   y[i] =     3.7500;  t[i++] =   5.7500;
319   y[i] =    24.1900;  t[i++] =   3.0000;
320   y[i] =    25.7600;  t[i++] =   3.0000;
321   y[i] =    18.0700;  t[i++] =   3.0000;
322   y[i] =    11.8100;  t[i++] =   3.0000;
323   y[i] =    12.0700;  t[i++] =   3.0000;
324   y[i] =    16.1200;  t[i++] =   3.0000;
325   y[i] =    70.8000;  t[i++] =    .5000;
326   y[i] =    54.7000;  t[i++] =    .7500;
327   y[i] =    48.0000;  t[i++] =   1.0000;
328   y[i] =    39.8000;  t[i++] =   1.5000;
329   y[i] =    29.8000;  t[i++] =   2.0000;
330   y[i] =    23.7000;  t[i++] =   2.5000;
331   y[i] =    29.6200;  t[i++] =   2.0000;
332   y[i] =    23.8100;  t[i++] =   2.5000;
333   y[i] =    17.7000;  t[i++] =   3.0000;
334   y[i] =    11.5500;  t[i++] =   4.0000;
335   y[i] =    12.0700;  t[i++] =   5.0000;
336   y[i] =     8.7400;  t[i++] =   6.0000;
337   y[i] =    80.7000;  t[i++] =    .5000;
338   y[i] =    61.3000;  t[i++] =    .7500;
339   y[i] =    47.5000;  t[i++] =   1.0000;
340    y[i] =   29.0000;  t[i++] =   1.5000;
341    y[i] =   24.0000;  t[i++] =   2.0000;
342   y[i] =    17.7000;  t[i++] =   2.5000;
343   y[i] =    24.5600;  t[i++] =   2.0000;
344   y[i] =    18.6700;  t[i++] =   2.5000;
345    y[i] =   16.2400;  t[i++] =   3.0000;
346   y[i] =     8.7400;  t[i++] =   4.0000;
347   y[i] =     7.8700;  t[i++] =   5.0000;
348   y[i] =     8.5100;  t[i++] =   6.0000;
349   y[i] =    66.7000;  t[i++] =    .5000;
350   y[i] =    59.2000;  t[i++] =    .7500;
351   y[i] =    40.8000;  t[i++] =   1.0000;
352   y[i] =    30.7000;  t[i++] =   1.5000;
353   y[i] =    25.7000;  t[i++] =   2.0000;
354   y[i] =    16.3000;  t[i++] =   2.5000;
355   y[i] =    25.9900;  t[i++] =   2.0000;
356   y[i] =    16.9500;  t[i++] =   2.5000;
357   y[i] =    13.3500;  t[i++] =   3.0000;
358   y[i] =     8.6200;  t[i++] =   4.0000;
359   y[i] =     7.2000;  t[i++] =   5.0000;
360   y[i] =     6.6400;  t[i++] =   6.0000;
361   y[i] =    13.6900;  t[i++] =   3.0000;
362   y[i] =    81.0000;  t[i++] =    .5000;
363   y[i] =    64.5000;  t[i++] =    .7500;
364   y[i] =    35.5000;  t[i++] =   1.5000;
365    y[i] =   13.3100;  t[i++] =   3.0000;
366   y[i] =     4.8700;  t[i++] =   6.0000;
367   y[i] =    12.9400;  t[i++] =   3.0000;
368   y[i] =     5.0600;  t[i++] =   6.0000;
369   y[i] =    15.1900;  t[i++] =   3.0000;
370   y[i] =    14.6200;  t[i++] =   3.0000;
371   y[i] =    15.6400;  t[i++] =   3.0000;
372   y[i] =    25.5000;  t[i++] =   1.7500;
373   y[i] =    25.9500;  t[i++] =   1.7500;
374   y[i] =    81.7000;  t[i++] =    .5000;
375   y[i] =    61.6000;  t[i++] =    .7500;
376   y[i] =    29.8000;  t[i++] =   1.7500;
377   y[i] =    29.8100;  t[i++] =   1.7500;
378   y[i] =    17.1700;  t[i++] =   2.7500;
379   y[i] =    10.3900;  t[i++] =   3.7500;
380   y[i] =    28.4000;  t[i++] =   1.7500;
381   y[i] =    28.6900;  t[i++] =   1.7500;
382   y[i] =    81.3000;  t[i++] =    .5000;
383   y[i] =    60.9000;  t[i++] =    .7500;
384   y[i] =    16.6500;  t[i++] =   2.7500;
385   y[i] =    10.0500;  t[i++] =   3.7500;
386   y[i] =    28.9000;  t[i++] =   1.7500;
387   y[i] =    28.9500;  t[i++] =   1.7500;
388   PetscFunctionReturn(0);
389 }
390 
391 /*TEST
392 
393    build:
394       requires: !complex !single
395 
396    test:
397       args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5
398 
399    test:
400       suffix: 2
401       args: -tao_smonitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l2prox -tao_brgn_regularizer_weight 1e-4 -tao_gatol 1.e-5
402 
403    test:
404       suffix: 3
405       args: -tao_smonitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l1dict -tao_brgn_regularizer_weight 1e-4 -tao_brgn_l1_smooth_epsilon 1e-6 -tao_gatol 1.e-5
406 
407    test:
408       suffix: 4
409       args: -tao_smonitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type lm -tao_gatol 1.e-5 -tao_brgn_subsolver_tao_type bnls
410 
411 TEST*/
412