1c4762a1bSJed Brown /*
2c4762a1bSJed Brown Include "petsctao.h" so that we can use TAO solvers. Note that this
3c4762a1bSJed Brown file automatically includes libraries such as:
4c4762a1bSJed Brown petsc.h - base PETSc routines petscvec.h - vectors
5a5b23f4aSJose E. Roman petscsys.h - system routines petscmat.h - matrices
6c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods
7c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners
8c4762a1bSJed Brown
9c4762a1bSJed Brown */
10c4762a1bSJed Brown
11c4762a1bSJed Brown #include <petsctao.h>
12c4762a1bSJed Brown
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown Description: These data are the result of a NIST study involving
15c4762a1bSJed Brown ultrasonic calibration. The response variable is
16c4762a1bSJed Brown ultrasonic response, and the predictor variable is
17c4762a1bSJed Brown metal distance.
18c4762a1bSJed Brown
19c4762a1bSJed Brown Reference: Chwirut, D., NIST (197?).
20c4762a1bSJed Brown Ultrasonic Reference Block Study.
21c4762a1bSJed Brown */
22c4762a1bSJed Brown
23c4762a1bSJed Brown static char help[] = "Finds the nonlinear least-squares solution to the model \n\
24c4762a1bSJed Brown y = exp[-b1*x]/(b2+b3*x) + e \n";
25c4762a1bSJed Brown
26c4762a1bSJed Brown #define NOBSERVATIONS 214
27c4762a1bSJed Brown #define NPARAMETERS 3
28c4762a1bSJed Brown
29c4762a1bSJed Brown #define DIE_TAG 2000
30c4762a1bSJed Brown #define IDLE_TAG 1000
31c4762a1bSJed Brown
32c4762a1bSJed Brown /* User-defined application context */
33c4762a1bSJed Brown typedef struct {
34c4762a1bSJed Brown /* Working space */
35c4762a1bSJed Brown PetscReal t[NOBSERVATIONS]; /* array of independent variables of observation */
36c4762a1bSJed Brown PetscReal y[NOBSERVATIONS]; /* array of dependent variables */
37c4762a1bSJed Brown PetscMPIInt size, rank;
38c4762a1bSJed Brown } AppCtx;
39c4762a1bSJed Brown
40c4762a1bSJed Brown /* User provided Routines */
41c4762a1bSJed Brown PetscErrorCode InitializeData(AppCtx *user);
42c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec);
43c4762a1bSJed Brown PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *);
44c4762a1bSJed Brown PetscErrorCode TaskWorker(AppCtx *user);
45c4762a1bSJed Brown PetscErrorCode StopWorkers(AppCtx *user);
46c4762a1bSJed Brown PetscErrorCode RunSimulation(PetscReal *x, PetscInt i, PetscReal *f, AppCtx *user);
47c4762a1bSJed Brown
48c4762a1bSJed Brown /*--------------------------------------------------------------------*/
main(int argc,char ** argv)49d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
50d71ae5a4SJacob Faibussowitsch {
51c4762a1bSJed Brown Vec x, f; /* solution, function */
52c4762a1bSJed Brown Tao tao; /* Tao solver context */
53c4762a1bSJed Brown AppCtx user; /* user-defined work context */
54c4762a1bSJed Brown
55c4762a1bSJed Brown /* Initialize TAO and PETSc */
56327415f7SBarry Smith PetscFunctionBeginUser;
57*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
58c4762a1bSJed Brown MPI_Comm_size(MPI_COMM_WORLD, &user.size);
59c4762a1bSJed Brown MPI_Comm_rank(MPI_COMM_WORLD, &user.rank);
609566063dSJacob Faibussowitsch PetscCall(InitializeData(&user));
61c4762a1bSJed Brown
62c4762a1bSJed Brown /* Run optimization on rank 0 */
63c4762a1bSJed Brown if (user.rank == 0) {
64c4762a1bSJed Brown /* Allocate vectors */
659566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, NPARAMETERS, &x));
669566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, NOBSERVATIONS, &f));
67c4762a1bSJed Brown
68c4762a1bSJed Brown /* TAO code begins here */
69c4762a1bSJed Brown
70c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
719566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
729566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOPOUNDERS));
73c4762a1bSJed Brown
74c4762a1bSJed Brown /* Set the function and Jacobian routines. */
759566063dSJacob Faibussowitsch PetscCall(FormStartingPoint(x));
769566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
779566063dSJacob Faibussowitsch PetscCall(TaoSetResidualRoutine(tao, f, EvaluateFunction, (void *)&user));
78c4762a1bSJed Brown
79c4762a1bSJed Brown /* Check for any TAO command line arguments */
809566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
81c4762a1bSJed Brown
82c4762a1bSJed Brown /* Perform the Solve */
839566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
84c4762a1bSJed Brown
85c4762a1bSJed Brown /* Free TAO data structures */
869566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
87c4762a1bSJed Brown
88c4762a1bSJed Brown /* Free PETSc data structures */
899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&f));
913ba16761SJacob Faibussowitsch PetscCall(StopWorkers(&user));
92c4762a1bSJed Brown } else {
933ba16761SJacob Faibussowitsch PetscCall(TaskWorker(&user));
94c4762a1bSJed Brown }
959566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
96b122ec5aSJacob Faibussowitsch return 0;
97c4762a1bSJed Brown }
98c4762a1bSJed Brown
99c4762a1bSJed Brown /*--------------------------------------------------------------------*/
EvaluateFunction(Tao tao,Vec X,Vec F,void * ptr)100d71ae5a4SJacob Faibussowitsch PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr)
101d71ae5a4SJacob Faibussowitsch {
102c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
103c4762a1bSJed Brown PetscInt i;
104c4762a1bSJed Brown PetscReal *x, *f;
105c4762a1bSJed Brown
106c4762a1bSJed Brown PetscFunctionBegin;
1079566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x));
1089566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
109c4762a1bSJed Brown if (user->size == 1) {
110c4762a1bSJed Brown /* Single processor */
11148a46eb9SPierre Jolivet for (i = 0; i < NOBSERVATIONS; i++) PetscCall(RunSimulation(x, i, &f[i], user));
112c4762a1bSJed Brown } else {
1139dddd249SSatish Balay /* Multiprocessor main */
114c4762a1bSJed Brown PetscMPIInt tag;
115c4762a1bSJed Brown PetscInt finishedtasks, next_task, checkedin;
116c4762a1bSJed Brown PetscReal f_i = 0.0;
117c4762a1bSJed Brown MPI_Status status;
118c4762a1bSJed Brown
119c4762a1bSJed Brown next_task = 0;
120c4762a1bSJed Brown finishedtasks = 0;
121c4762a1bSJed Brown checkedin = 0;
122c4762a1bSJed Brown
123c4762a1bSJed Brown while (finishedtasks < NOBSERVATIONS || checkedin < user->size - 1) {
1249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&f_i, 1, MPIU_REAL, MPI_ANY_SOURCE, MPI_ANY_TAG, PETSC_COMM_WORLD, &status));
125c4762a1bSJed Brown if (status.MPI_TAG == IDLE_TAG) {
126c4762a1bSJed Brown checkedin++;
127c4762a1bSJed Brown } else {
128c4762a1bSJed Brown tag = status.MPI_TAG;
129c4762a1bSJed Brown f[tag] = (PetscReal)f_i;
130c4762a1bSJed Brown finishedtasks++;
131c4762a1bSJed Brown }
132c4762a1bSJed Brown
133c4762a1bSJed Brown if (next_task < NOBSERVATIONS) {
1346497c311SBarry Smith PetscCallMPI(MPI_Send(x, NPARAMETERS, MPIU_REAL, status.MPI_SOURCE, (PetscMPIInt)next_task, PETSC_COMM_WORLD));
135c4762a1bSJed Brown next_task++;
136c4762a1bSJed Brown
137c4762a1bSJed Brown } else {
138c4762a1bSJed Brown /* Send idle message */
1399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(x, NPARAMETERS, MPIU_REAL, status.MPI_SOURCE, IDLE_TAG, PETSC_COMM_WORLD));
140c4762a1bSJed Brown }
141c4762a1bSJed Brown }
142c4762a1bSJed Brown }
1439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x));
1449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
1453ba16761SJacob Faibussowitsch PetscCall(PetscLogFlops(6 * NOBSERVATIONS));
1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
147c4762a1bSJed Brown }
148c4762a1bSJed Brown
149c4762a1bSJed Brown /* ------------------------------------------------------------ */
FormStartingPoint(Vec X)150d71ae5a4SJacob Faibussowitsch PetscErrorCode FormStartingPoint(Vec X)
151d71ae5a4SJacob Faibussowitsch {
152c4762a1bSJed Brown PetscReal *x;
153c4762a1bSJed Brown
154c4762a1bSJed Brown PetscFunctionBegin;
1559566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x));
156c4762a1bSJed Brown x[0] = 0.15;
157c4762a1bSJed Brown x[1] = 0.008;
158c4762a1bSJed Brown x[2] = 0.010;
1599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x));
1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
161c4762a1bSJed Brown }
162c4762a1bSJed Brown
163c4762a1bSJed Brown /* ---------------------------------------------------------------------- */
InitializeData(AppCtx * user)164d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeData(AppCtx *user)
165d71ae5a4SJacob Faibussowitsch {
166c4762a1bSJed Brown PetscReal *t = user->t, *y = user->y;
167c4762a1bSJed Brown PetscInt i = 0;
168c4762a1bSJed Brown
169c4762a1bSJed Brown PetscFunctionBegin;
1709371c9d4SSatish Balay y[i] = 92.9000;
1719371c9d4SSatish Balay t[i++] = 0.5000;
1729371c9d4SSatish Balay y[i] = 78.7000;
1739371c9d4SSatish Balay t[i++] = 0.6250;
1749371c9d4SSatish Balay y[i] = 64.2000;
1759371c9d4SSatish Balay t[i++] = 0.7500;
1769371c9d4SSatish Balay y[i] = 64.9000;
1779371c9d4SSatish Balay t[i++] = 0.8750;
1789371c9d4SSatish Balay y[i] = 57.1000;
1799371c9d4SSatish Balay t[i++] = 1.0000;
1809371c9d4SSatish Balay y[i] = 43.3000;
1819371c9d4SSatish Balay t[i++] = 1.2500;
1829371c9d4SSatish Balay y[i] = 31.1000;
1839371c9d4SSatish Balay t[i++] = 1.7500;
1849371c9d4SSatish Balay y[i] = 23.6000;
1859371c9d4SSatish Balay t[i++] = 2.2500;
1869371c9d4SSatish Balay y[i] = 31.0500;
1879371c9d4SSatish Balay t[i++] = 1.7500;
1889371c9d4SSatish Balay y[i] = 23.7750;
1899371c9d4SSatish Balay t[i++] = 2.2500;
1909371c9d4SSatish Balay y[i] = 17.7375;
1919371c9d4SSatish Balay t[i++] = 2.7500;
1929371c9d4SSatish Balay y[i] = 13.8000;
1939371c9d4SSatish Balay t[i++] = 3.2500;
1949371c9d4SSatish Balay y[i] = 11.5875;
1959371c9d4SSatish Balay t[i++] = 3.7500;
1969371c9d4SSatish Balay y[i] = 9.4125;
1979371c9d4SSatish Balay t[i++] = 4.2500;
1989371c9d4SSatish Balay y[i] = 7.7250;
1999371c9d4SSatish Balay t[i++] = 4.7500;
2009371c9d4SSatish Balay y[i] = 7.3500;
2019371c9d4SSatish Balay t[i++] = 5.2500;
2029371c9d4SSatish Balay y[i] = 8.0250;
2039371c9d4SSatish Balay t[i++] = 5.7500;
2049371c9d4SSatish Balay y[i] = 90.6000;
2059371c9d4SSatish Balay t[i++] = 0.5000;
2069371c9d4SSatish Balay y[i] = 76.9000;
2079371c9d4SSatish Balay t[i++] = 0.6250;
2089371c9d4SSatish Balay y[i] = 71.6000;
2099371c9d4SSatish Balay t[i++] = 0.7500;
2109371c9d4SSatish Balay y[i] = 63.6000;
2119371c9d4SSatish Balay t[i++] = 0.8750;
2129371c9d4SSatish Balay y[i] = 54.0000;
2139371c9d4SSatish Balay t[i++] = 1.0000;
2149371c9d4SSatish Balay y[i] = 39.2000;
2159371c9d4SSatish Balay t[i++] = 1.2500;
2169371c9d4SSatish Balay y[i] = 29.3000;
2179371c9d4SSatish Balay t[i++] = 1.7500;
2189371c9d4SSatish Balay y[i] = 21.4000;
2199371c9d4SSatish Balay t[i++] = 2.2500;
2209371c9d4SSatish Balay y[i] = 29.1750;
2219371c9d4SSatish Balay t[i++] = 1.7500;
2229371c9d4SSatish Balay y[i] = 22.1250;
2239371c9d4SSatish Balay t[i++] = 2.2500;
2249371c9d4SSatish Balay y[i] = 17.5125;
2259371c9d4SSatish Balay t[i++] = 2.7500;
2269371c9d4SSatish Balay y[i] = 14.2500;
2279371c9d4SSatish Balay t[i++] = 3.2500;
2289371c9d4SSatish Balay y[i] = 9.4500;
2299371c9d4SSatish Balay t[i++] = 3.7500;
2309371c9d4SSatish Balay y[i] = 9.1500;
2319371c9d4SSatish Balay t[i++] = 4.2500;
2329371c9d4SSatish Balay y[i] = 7.9125;
2339371c9d4SSatish Balay t[i++] = 4.7500;
2349371c9d4SSatish Balay y[i] = 8.4750;
2359371c9d4SSatish Balay t[i++] = 5.2500;
2369371c9d4SSatish Balay y[i] = 6.1125;
2379371c9d4SSatish Balay t[i++] = 5.7500;
2389371c9d4SSatish Balay y[i] = 80.0000;
2399371c9d4SSatish Balay t[i++] = 0.5000;
2409371c9d4SSatish Balay y[i] = 79.0000;
2419371c9d4SSatish Balay t[i++] = 0.6250;
2429371c9d4SSatish Balay y[i] = 63.8000;
2439371c9d4SSatish Balay t[i++] = 0.7500;
2449371c9d4SSatish Balay y[i] = 57.2000;
2459371c9d4SSatish Balay t[i++] = 0.8750;
2469371c9d4SSatish Balay y[i] = 53.2000;
2479371c9d4SSatish Balay t[i++] = 1.0000;
2489371c9d4SSatish Balay y[i] = 42.5000;
2499371c9d4SSatish Balay t[i++] = 1.2500;
2509371c9d4SSatish Balay y[i] = 26.8000;
2519371c9d4SSatish Balay t[i++] = 1.7500;
2529371c9d4SSatish Balay y[i] = 20.4000;
2539371c9d4SSatish Balay t[i++] = 2.2500;
2549371c9d4SSatish Balay y[i] = 26.8500;
2559371c9d4SSatish Balay t[i++] = 1.7500;
2569371c9d4SSatish Balay y[i] = 21.0000;
2579371c9d4SSatish Balay t[i++] = 2.2500;
2589371c9d4SSatish Balay y[i] = 16.4625;
2599371c9d4SSatish Balay t[i++] = 2.7500;
2609371c9d4SSatish Balay y[i] = 12.5250;
2619371c9d4SSatish Balay t[i++] = 3.2500;
2629371c9d4SSatish Balay y[i] = 10.5375;
2639371c9d4SSatish Balay t[i++] = 3.7500;
2649371c9d4SSatish Balay y[i] = 8.5875;
2659371c9d4SSatish Balay t[i++] = 4.2500;
2669371c9d4SSatish Balay y[i] = 7.1250;
2679371c9d4SSatish Balay t[i++] = 4.7500;
2689371c9d4SSatish Balay y[i] = 6.1125;
2699371c9d4SSatish Balay t[i++] = 5.2500;
2709371c9d4SSatish Balay y[i] = 5.9625;
2719371c9d4SSatish Balay t[i++] = 5.7500;
2729371c9d4SSatish Balay y[i] = 74.1000;
2739371c9d4SSatish Balay t[i++] = 0.5000;
2749371c9d4SSatish Balay y[i] = 67.3000;
2759371c9d4SSatish Balay t[i++] = 0.6250;
2769371c9d4SSatish Balay y[i] = 60.8000;
2779371c9d4SSatish Balay t[i++] = 0.7500;
2789371c9d4SSatish Balay y[i] = 55.5000;
2799371c9d4SSatish Balay t[i++] = 0.8750;
2809371c9d4SSatish Balay y[i] = 50.3000;
2819371c9d4SSatish Balay t[i++] = 1.0000;
2829371c9d4SSatish Balay y[i] = 41.0000;
2839371c9d4SSatish Balay t[i++] = 1.2500;
2849371c9d4SSatish Balay y[i] = 29.4000;
2859371c9d4SSatish Balay t[i++] = 1.7500;
2869371c9d4SSatish Balay y[i] = 20.4000;
2879371c9d4SSatish Balay t[i++] = 2.2500;
2889371c9d4SSatish Balay y[i] = 29.3625;
2899371c9d4SSatish Balay t[i++] = 1.7500;
2909371c9d4SSatish Balay y[i] = 21.1500;
2919371c9d4SSatish Balay t[i++] = 2.2500;
2929371c9d4SSatish Balay y[i] = 16.7625;
2939371c9d4SSatish Balay t[i++] = 2.7500;
2949371c9d4SSatish Balay y[i] = 13.2000;
2959371c9d4SSatish Balay t[i++] = 3.2500;
2969371c9d4SSatish Balay y[i] = 10.8750;
2979371c9d4SSatish Balay t[i++] = 3.7500;
2989371c9d4SSatish Balay y[i] = 8.1750;
2999371c9d4SSatish Balay t[i++] = 4.2500;
3009371c9d4SSatish Balay y[i] = 7.3500;
3019371c9d4SSatish Balay t[i++] = 4.7500;
3029371c9d4SSatish Balay y[i] = 5.9625;
3039371c9d4SSatish Balay t[i++] = 5.2500;
3049371c9d4SSatish Balay y[i] = 5.6250;
3059371c9d4SSatish Balay t[i++] = 5.7500;
3069371c9d4SSatish Balay y[i] = 81.5000;
3079371c9d4SSatish Balay t[i++] = .5000;
3089371c9d4SSatish Balay y[i] = 62.4000;
3099371c9d4SSatish Balay t[i++] = .7500;
3109371c9d4SSatish Balay y[i] = 32.5000;
3119371c9d4SSatish Balay t[i++] = 1.5000;
3129371c9d4SSatish Balay y[i] = 12.4100;
3139371c9d4SSatish Balay t[i++] = 3.0000;
3149371c9d4SSatish Balay y[i] = 13.1200;
3159371c9d4SSatish Balay t[i++] = 3.0000;
3169371c9d4SSatish Balay y[i] = 15.5600;
3179371c9d4SSatish Balay t[i++] = 3.0000;
3189371c9d4SSatish Balay y[i] = 5.6300;
3199371c9d4SSatish Balay t[i++] = 6.0000;
3209371c9d4SSatish Balay y[i] = 78.0000;
3219371c9d4SSatish Balay t[i++] = .5000;
3229371c9d4SSatish Balay y[i] = 59.9000;
3239371c9d4SSatish Balay t[i++] = .7500;
3249371c9d4SSatish Balay y[i] = 33.2000;
3259371c9d4SSatish Balay t[i++] = 1.5000;
3269371c9d4SSatish Balay y[i] = 13.8400;
3279371c9d4SSatish Balay t[i++] = 3.0000;
3289371c9d4SSatish Balay y[i] = 12.7500;
3299371c9d4SSatish Balay t[i++] = 3.0000;
3309371c9d4SSatish Balay y[i] = 14.6200;
3319371c9d4SSatish Balay t[i++] = 3.0000;
3329371c9d4SSatish Balay y[i] = 3.9400;
3339371c9d4SSatish Balay t[i++] = 6.0000;
3349371c9d4SSatish Balay y[i] = 76.8000;
3359371c9d4SSatish Balay t[i++] = .5000;
3369371c9d4SSatish Balay y[i] = 61.0000;
3379371c9d4SSatish Balay t[i++] = .7500;
3389371c9d4SSatish Balay y[i] = 32.9000;
3399371c9d4SSatish Balay t[i++] = 1.5000;
3409371c9d4SSatish Balay y[i] = 13.8700;
3419371c9d4SSatish Balay t[i++] = 3.0000;
3429371c9d4SSatish Balay y[i] = 11.8100;
3439371c9d4SSatish Balay t[i++] = 3.0000;
3449371c9d4SSatish Balay y[i] = 13.3100;
3459371c9d4SSatish Balay t[i++] = 3.0000;
3469371c9d4SSatish Balay y[i] = 5.4400;
3479371c9d4SSatish Balay t[i++] = 6.0000;
3489371c9d4SSatish Balay y[i] = 78.0000;
3499371c9d4SSatish Balay t[i++] = .5000;
3509371c9d4SSatish Balay y[i] = 63.5000;
3519371c9d4SSatish Balay t[i++] = .7500;
3529371c9d4SSatish Balay y[i] = 33.8000;
3539371c9d4SSatish Balay t[i++] = 1.5000;
3549371c9d4SSatish Balay y[i] = 12.5600;
3559371c9d4SSatish Balay t[i++] = 3.0000;
3569371c9d4SSatish Balay y[i] = 5.6300;
3579371c9d4SSatish Balay t[i++] = 6.0000;
3589371c9d4SSatish Balay y[i] = 12.7500;
3599371c9d4SSatish Balay t[i++] = 3.0000;
3609371c9d4SSatish Balay y[i] = 13.1200;
3619371c9d4SSatish Balay t[i++] = 3.0000;
3629371c9d4SSatish Balay y[i] = 5.4400;
3639371c9d4SSatish Balay t[i++] = 6.0000;
3649371c9d4SSatish Balay y[i] = 76.8000;
3659371c9d4SSatish Balay t[i++] = .5000;
3669371c9d4SSatish Balay y[i] = 60.0000;
3679371c9d4SSatish Balay t[i++] = .7500;
3689371c9d4SSatish Balay y[i] = 47.8000;
3699371c9d4SSatish Balay t[i++] = 1.0000;
3709371c9d4SSatish Balay y[i] = 32.0000;
3719371c9d4SSatish Balay t[i++] = 1.5000;
3729371c9d4SSatish Balay y[i] = 22.2000;
3739371c9d4SSatish Balay t[i++] = 2.0000;
3749371c9d4SSatish Balay y[i] = 22.5700;
3759371c9d4SSatish Balay t[i++] = 2.0000;
3769371c9d4SSatish Balay y[i] = 18.8200;
3779371c9d4SSatish Balay t[i++] = 2.5000;
3789371c9d4SSatish Balay y[i] = 13.9500;
3799371c9d4SSatish Balay t[i++] = 3.0000;
3809371c9d4SSatish Balay y[i] = 11.2500;
3819371c9d4SSatish Balay t[i++] = 4.0000;
3829371c9d4SSatish Balay y[i] = 9.0000;
3839371c9d4SSatish Balay t[i++] = 5.0000;
3849371c9d4SSatish Balay y[i] = 6.6700;
3859371c9d4SSatish Balay t[i++] = 6.0000;
3869371c9d4SSatish Balay y[i] = 75.8000;
3879371c9d4SSatish Balay t[i++] = .5000;
3889371c9d4SSatish Balay y[i] = 62.0000;
3899371c9d4SSatish Balay t[i++] = .7500;
3909371c9d4SSatish Balay y[i] = 48.8000;
3919371c9d4SSatish Balay t[i++] = 1.0000;
3929371c9d4SSatish Balay y[i] = 35.2000;
3939371c9d4SSatish Balay t[i++] = 1.5000;
3949371c9d4SSatish Balay y[i] = 20.0000;
3959371c9d4SSatish Balay t[i++] = 2.0000;
3969371c9d4SSatish Balay y[i] = 20.3200;
3979371c9d4SSatish Balay t[i++] = 2.0000;
3989371c9d4SSatish Balay y[i] = 19.3100;
3999371c9d4SSatish Balay t[i++] = 2.5000;
4009371c9d4SSatish Balay y[i] = 12.7500;
4019371c9d4SSatish Balay t[i++] = 3.0000;
4029371c9d4SSatish Balay y[i] = 10.4200;
4039371c9d4SSatish Balay t[i++] = 4.0000;
4049371c9d4SSatish Balay y[i] = 7.3100;
4059371c9d4SSatish Balay t[i++] = 5.0000;
4069371c9d4SSatish Balay y[i] = 7.4200;
4079371c9d4SSatish Balay t[i++] = 6.0000;
4089371c9d4SSatish Balay y[i] = 70.5000;
4099371c9d4SSatish Balay t[i++] = .5000;
4109371c9d4SSatish Balay y[i] = 59.5000;
4119371c9d4SSatish Balay t[i++] = .7500;
4129371c9d4SSatish Balay y[i] = 48.5000;
4139371c9d4SSatish Balay t[i++] = 1.0000;
4149371c9d4SSatish Balay y[i] = 35.8000;
4159371c9d4SSatish Balay t[i++] = 1.5000;
4169371c9d4SSatish Balay y[i] = 21.0000;
4179371c9d4SSatish Balay t[i++] = 2.0000;
4189371c9d4SSatish Balay y[i] = 21.6700;
4199371c9d4SSatish Balay t[i++] = 2.0000;
4209371c9d4SSatish Balay y[i] = 21.0000;
4219371c9d4SSatish Balay t[i++] = 2.5000;
4229371c9d4SSatish Balay y[i] = 15.6400;
4239371c9d4SSatish Balay t[i++] = 3.0000;
4249371c9d4SSatish Balay y[i] = 8.1700;
4259371c9d4SSatish Balay t[i++] = 4.0000;
4269371c9d4SSatish Balay y[i] = 8.5500;
4279371c9d4SSatish Balay t[i++] = 5.0000;
4289371c9d4SSatish Balay y[i] = 10.1200;
4299371c9d4SSatish Balay t[i++] = 6.0000;
4309371c9d4SSatish Balay y[i] = 78.0000;
4319371c9d4SSatish Balay t[i++] = .5000;
4329371c9d4SSatish Balay y[i] = 66.0000;
4339371c9d4SSatish Balay t[i++] = .6250;
4349371c9d4SSatish Balay y[i] = 62.0000;
4359371c9d4SSatish Balay t[i++] = .7500;
4369371c9d4SSatish Balay y[i] = 58.0000;
4379371c9d4SSatish Balay t[i++] = .8750;
4389371c9d4SSatish Balay y[i] = 47.7000;
4399371c9d4SSatish Balay t[i++] = 1.0000;
4409371c9d4SSatish Balay y[i] = 37.8000;
4419371c9d4SSatish Balay t[i++] = 1.2500;
4429371c9d4SSatish Balay y[i] = 20.2000;
4439371c9d4SSatish Balay t[i++] = 2.2500;
4449371c9d4SSatish Balay y[i] = 21.0700;
4459371c9d4SSatish Balay t[i++] = 2.2500;
4469371c9d4SSatish Balay y[i] = 13.8700;
4479371c9d4SSatish Balay t[i++] = 2.7500;
4489371c9d4SSatish Balay y[i] = 9.6700;
4499371c9d4SSatish Balay t[i++] = 3.2500;
4509371c9d4SSatish Balay y[i] = 7.7600;
4519371c9d4SSatish Balay t[i++] = 3.7500;
4529371c9d4SSatish Balay y[i] = 5.4400;
4539371c9d4SSatish Balay t[i++] = 4.2500;
4549371c9d4SSatish Balay y[i] = 4.8700;
4559371c9d4SSatish Balay t[i++] = 4.7500;
4569371c9d4SSatish Balay y[i] = 4.0100;
4579371c9d4SSatish Balay t[i++] = 5.2500;
4589371c9d4SSatish Balay y[i] = 3.7500;
4599371c9d4SSatish Balay t[i++] = 5.7500;
4609371c9d4SSatish Balay y[i] = 24.1900;
4619371c9d4SSatish Balay t[i++] = 3.0000;
4629371c9d4SSatish Balay y[i] = 25.7600;
4639371c9d4SSatish Balay t[i++] = 3.0000;
4649371c9d4SSatish Balay y[i] = 18.0700;
4659371c9d4SSatish Balay t[i++] = 3.0000;
4669371c9d4SSatish Balay y[i] = 11.8100;
4679371c9d4SSatish Balay t[i++] = 3.0000;
4689371c9d4SSatish Balay y[i] = 12.0700;
4699371c9d4SSatish Balay t[i++] = 3.0000;
4709371c9d4SSatish Balay y[i] = 16.1200;
4719371c9d4SSatish Balay t[i++] = 3.0000;
4729371c9d4SSatish Balay y[i] = 70.8000;
4739371c9d4SSatish Balay t[i++] = .5000;
4749371c9d4SSatish Balay y[i] = 54.7000;
4759371c9d4SSatish Balay t[i++] = .7500;
4769371c9d4SSatish Balay y[i] = 48.0000;
4779371c9d4SSatish Balay t[i++] = 1.0000;
4789371c9d4SSatish Balay y[i] = 39.8000;
4799371c9d4SSatish Balay t[i++] = 1.5000;
4809371c9d4SSatish Balay y[i] = 29.8000;
4819371c9d4SSatish Balay t[i++] = 2.0000;
4829371c9d4SSatish Balay y[i] = 23.7000;
4839371c9d4SSatish Balay t[i++] = 2.5000;
4849371c9d4SSatish Balay y[i] = 29.6200;
4859371c9d4SSatish Balay t[i++] = 2.0000;
4869371c9d4SSatish Balay y[i] = 23.8100;
4879371c9d4SSatish Balay t[i++] = 2.5000;
4889371c9d4SSatish Balay y[i] = 17.7000;
4899371c9d4SSatish Balay t[i++] = 3.0000;
4909371c9d4SSatish Balay y[i] = 11.5500;
4919371c9d4SSatish Balay t[i++] = 4.0000;
4929371c9d4SSatish Balay y[i] = 12.0700;
4939371c9d4SSatish Balay t[i++] = 5.0000;
4949371c9d4SSatish Balay y[i] = 8.7400;
4959371c9d4SSatish Balay t[i++] = 6.0000;
4969371c9d4SSatish Balay y[i] = 80.7000;
4979371c9d4SSatish Balay t[i++] = .5000;
4989371c9d4SSatish Balay y[i] = 61.3000;
4999371c9d4SSatish Balay t[i++] = .7500;
5009371c9d4SSatish Balay y[i] = 47.5000;
5019371c9d4SSatish Balay t[i++] = 1.0000;
5029371c9d4SSatish Balay y[i] = 29.0000;
5039371c9d4SSatish Balay t[i++] = 1.5000;
5049371c9d4SSatish Balay y[i] = 24.0000;
5059371c9d4SSatish Balay t[i++] = 2.0000;
5069371c9d4SSatish Balay y[i] = 17.7000;
5079371c9d4SSatish Balay t[i++] = 2.5000;
5089371c9d4SSatish Balay y[i] = 24.5600;
5099371c9d4SSatish Balay t[i++] = 2.0000;
5109371c9d4SSatish Balay y[i] = 18.6700;
5119371c9d4SSatish Balay t[i++] = 2.5000;
5129371c9d4SSatish Balay y[i] = 16.2400;
5139371c9d4SSatish Balay t[i++] = 3.0000;
5149371c9d4SSatish Balay y[i] = 8.7400;
5159371c9d4SSatish Balay t[i++] = 4.0000;
5169371c9d4SSatish Balay y[i] = 7.8700;
5179371c9d4SSatish Balay t[i++] = 5.0000;
5189371c9d4SSatish Balay y[i] = 8.5100;
5199371c9d4SSatish Balay t[i++] = 6.0000;
5209371c9d4SSatish Balay y[i] = 66.7000;
5219371c9d4SSatish Balay t[i++] = .5000;
5229371c9d4SSatish Balay y[i] = 59.2000;
5239371c9d4SSatish Balay t[i++] = .7500;
5249371c9d4SSatish Balay y[i] = 40.8000;
5259371c9d4SSatish Balay t[i++] = 1.0000;
5269371c9d4SSatish Balay y[i] = 30.7000;
5279371c9d4SSatish Balay t[i++] = 1.5000;
5289371c9d4SSatish Balay y[i] = 25.7000;
5299371c9d4SSatish Balay t[i++] = 2.0000;
5309371c9d4SSatish Balay y[i] = 16.3000;
5319371c9d4SSatish Balay t[i++] = 2.5000;
5329371c9d4SSatish Balay y[i] = 25.9900;
5339371c9d4SSatish Balay t[i++] = 2.0000;
5349371c9d4SSatish Balay y[i] = 16.9500;
5359371c9d4SSatish Balay t[i++] = 2.5000;
5369371c9d4SSatish Balay y[i] = 13.3500;
5379371c9d4SSatish Balay t[i++] = 3.0000;
5389371c9d4SSatish Balay y[i] = 8.6200;
5399371c9d4SSatish Balay t[i++] = 4.0000;
5409371c9d4SSatish Balay y[i] = 7.2000;
5419371c9d4SSatish Balay t[i++] = 5.0000;
5429371c9d4SSatish Balay y[i] = 6.6400;
5439371c9d4SSatish Balay t[i++] = 6.0000;
5449371c9d4SSatish Balay y[i] = 13.6900;
5459371c9d4SSatish Balay t[i++] = 3.0000;
5469371c9d4SSatish Balay y[i] = 81.0000;
5479371c9d4SSatish Balay t[i++] = .5000;
5489371c9d4SSatish Balay y[i] = 64.5000;
5499371c9d4SSatish Balay t[i++] = .7500;
5509371c9d4SSatish Balay y[i] = 35.5000;
5519371c9d4SSatish Balay t[i++] = 1.5000;
5529371c9d4SSatish Balay y[i] = 13.3100;
5539371c9d4SSatish Balay t[i++] = 3.0000;
5549371c9d4SSatish Balay y[i] = 4.8700;
5559371c9d4SSatish Balay t[i++] = 6.0000;
5569371c9d4SSatish Balay y[i] = 12.9400;
5579371c9d4SSatish Balay t[i++] = 3.0000;
5589371c9d4SSatish Balay y[i] = 5.0600;
5599371c9d4SSatish Balay t[i++] = 6.0000;
5609371c9d4SSatish Balay y[i] = 15.1900;
5619371c9d4SSatish Balay t[i++] = 3.0000;
5629371c9d4SSatish Balay y[i] = 14.6200;
5639371c9d4SSatish Balay t[i++] = 3.0000;
5649371c9d4SSatish Balay y[i] = 15.6400;
5659371c9d4SSatish Balay t[i++] = 3.0000;
5669371c9d4SSatish Balay y[i] = 25.5000;
5679371c9d4SSatish Balay t[i++] = 1.7500;
5689371c9d4SSatish Balay y[i] = 25.9500;
5699371c9d4SSatish Balay t[i++] = 1.7500;
5709371c9d4SSatish Balay y[i] = 81.7000;
5719371c9d4SSatish Balay t[i++] = .5000;
5729371c9d4SSatish Balay y[i] = 61.6000;
5739371c9d4SSatish Balay t[i++] = .7500;
5749371c9d4SSatish Balay y[i] = 29.8000;
5759371c9d4SSatish Balay t[i++] = 1.7500;
5769371c9d4SSatish Balay y[i] = 29.8100;
5779371c9d4SSatish Balay t[i++] = 1.7500;
5789371c9d4SSatish Balay y[i] = 17.1700;
5799371c9d4SSatish Balay t[i++] = 2.7500;
5809371c9d4SSatish Balay y[i] = 10.3900;
5819371c9d4SSatish Balay t[i++] = 3.7500;
5829371c9d4SSatish Balay y[i] = 28.4000;
5839371c9d4SSatish Balay t[i++] = 1.7500;
5849371c9d4SSatish Balay y[i] = 28.6900;
5859371c9d4SSatish Balay t[i++] = 1.7500;
5869371c9d4SSatish Balay y[i] = 81.3000;
5879371c9d4SSatish Balay t[i++] = .5000;
5889371c9d4SSatish Balay y[i] = 60.9000;
5899371c9d4SSatish Balay t[i++] = .7500;
5909371c9d4SSatish Balay y[i] = 16.6500;
5919371c9d4SSatish Balay t[i++] = 2.7500;
5929371c9d4SSatish Balay y[i] = 10.0500;
5939371c9d4SSatish Balay t[i++] = 3.7500;
5949371c9d4SSatish Balay y[i] = 28.9000;
5959371c9d4SSatish Balay t[i++] = 1.7500;
5969371c9d4SSatish Balay y[i] = 28.9500;
5979371c9d4SSatish Balay t[i++] = 1.7500;
5983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
599c4762a1bSJed Brown }
600c4762a1bSJed Brown
TaskWorker(AppCtx * user)601d71ae5a4SJacob Faibussowitsch PetscErrorCode TaskWorker(AppCtx *user)
602d71ae5a4SJacob Faibussowitsch {
603c4762a1bSJed Brown PetscReal x[NPARAMETERS], f = 0.0;
604c4762a1bSJed Brown PetscMPIInt tag = IDLE_TAG;
605c4762a1bSJed Brown PetscInt index;
606c4762a1bSJed Brown MPI_Status status;
607c4762a1bSJed Brown
608c4762a1bSJed Brown PetscFunctionBegin;
6099dddd249SSatish Balay /* Send check-in message to rank-0 */
610c4762a1bSJed Brown
6119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&f, 1, MPIU_REAL, 0, IDLE_TAG, PETSC_COMM_WORLD));
612c4762a1bSJed Brown while (tag != DIE_TAG) {
6139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(x, NPARAMETERS, MPIU_REAL, 0, MPI_ANY_TAG, PETSC_COMM_WORLD, &status));
614c4762a1bSJed Brown tag = status.MPI_TAG;
615c4762a1bSJed Brown if (tag == IDLE_TAG) {
6169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&f, 1, MPIU_REAL, 0, IDLE_TAG, PETSC_COMM_WORLD));
617c4762a1bSJed Brown } else if (tag != DIE_TAG) {
618c4762a1bSJed Brown index = (PetscInt)tag;
6199566063dSJacob Faibussowitsch PetscCall(RunSimulation(x, index, &f, user));
6209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&f, 1, MPIU_REAL, 0, tag, PETSC_COMM_WORLD));
621c4762a1bSJed Brown }
622c4762a1bSJed Brown }
6233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
624c4762a1bSJed Brown }
625c4762a1bSJed Brown
RunSimulation(PetscReal * x,PetscInt i,PetscReal * f,AppCtx * user)626d71ae5a4SJacob Faibussowitsch PetscErrorCode RunSimulation(PetscReal *x, PetscInt i, PetscReal *f, AppCtx *user)
627d71ae5a4SJacob Faibussowitsch {
628c4762a1bSJed Brown PetscReal *t = user->t;
629c4762a1bSJed Brown PetscReal *y = user->y;
6303ba16761SJacob Faibussowitsch
6313ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
632c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE)
633e1dfdf8eSBarry Smith *f = y[i] - exp(-x[0] * t[i]) / (x[1] + x[2] * t[i]); /* expf() for single-precision breaks this example on Freebsd, Valgrind errors on Linux */
634c4762a1bSJed Brown #else
635c4762a1bSJed Brown *f = y[i] - PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]);
636c4762a1bSJed Brown #endif
6373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
638c4762a1bSJed Brown }
639c4762a1bSJed Brown
StopWorkers(AppCtx * user)640d71ae5a4SJacob Faibussowitsch PetscErrorCode StopWorkers(AppCtx *user)
641d71ae5a4SJacob Faibussowitsch {
642c4762a1bSJed Brown PetscInt checkedin;
643c4762a1bSJed Brown MPI_Status status;
644c4762a1bSJed Brown PetscReal f, x[NPARAMETERS];
645c4762a1bSJed Brown
6463ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
647c4762a1bSJed Brown checkedin = 0;
648c4762a1bSJed Brown while (checkedin < user->size - 1) {
6499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&f, 1, MPIU_REAL, MPI_ANY_SOURCE, MPI_ANY_TAG, PETSC_COMM_WORLD, &status));
650c4762a1bSJed Brown checkedin++;
6519566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(x, NPARAMETERS));
6529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(x, NPARAMETERS, MPIU_REAL, status.MPI_SOURCE, DIE_TAG, PETSC_COMM_WORLD));
653c4762a1bSJed Brown }
6543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
655c4762a1bSJed Brown }
656c4762a1bSJed Brown
657c4762a1bSJed Brown /*TEST
658c4762a1bSJed Brown
659c4762a1bSJed Brown build:
660c4762a1bSJed Brown requires: !complex
661c4762a1bSJed Brown
662c4762a1bSJed Brown test:
663c4762a1bSJed Brown nsize: 3
664c4762a1bSJed Brown requires: !single
66510978b7dSBarry Smith args: -tao_monitor_short -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5
666c4762a1bSJed Brown
667c4762a1bSJed Brown TEST*/
668