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 /*--------------------------------------------------------------------*/ 499371c9d4SSatish Balay int main(int argc, char **argv) { 50c4762a1bSJed Brown Vec x, f; /* solution, function */ 51c4762a1bSJed Brown Tao tao; /* Tao solver context */ 52c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Initialize TAO and PETSc */ 55327415f7SBarry Smith PetscFunctionBeginUser; 569566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 57c4762a1bSJed Brown MPI_Comm_size(MPI_COMM_WORLD, &user.size); 58c4762a1bSJed Brown MPI_Comm_rank(MPI_COMM_WORLD, &user.rank); 599566063dSJacob Faibussowitsch PetscCall(InitializeData(&user)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Run optimization on rank 0 */ 62c4762a1bSJed Brown if (user.rank == 0) { 63c4762a1bSJed Brown /* Allocate vectors */ 649566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, NPARAMETERS, &x)); 659566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, NOBSERVATIONS, &f)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* TAO code begins here */ 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 709566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao)); 719566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOPOUNDERS)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Set the function and Jacobian routines. */ 749566063dSJacob Faibussowitsch PetscCall(FormStartingPoint(x)); 759566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 769566063dSJacob Faibussowitsch PetscCall(TaoSetResidualRoutine(tao, f, EvaluateFunction, (void *)&user)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Check for any TAO command line arguments */ 799566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Perform the Solve */ 829566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Free TAO data structures */ 859566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Free PETSc data structures */ 889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&f)); 90c4762a1bSJed Brown StopWorkers(&user); 91c4762a1bSJed Brown } else { 92c4762a1bSJed Brown TaskWorker(&user); 93c4762a1bSJed Brown } 949566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 95b122ec5aSJacob Faibussowitsch return 0; 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown /*--------------------------------------------------------------------*/ 999371c9d4SSatish Balay PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr) { 100c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 101c4762a1bSJed Brown PetscInt i; 102c4762a1bSJed Brown PetscReal *x, *f; 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscFunctionBegin; 1059566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 1069566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 107c4762a1bSJed Brown if (user->size == 1) { 108c4762a1bSJed Brown /* Single processor */ 109*48a46eb9SPierre Jolivet for (i = 0; i < NOBSERVATIONS; i++) PetscCall(RunSimulation(x, i, &f[i], user)); 110c4762a1bSJed Brown } else { 1119dddd249SSatish Balay /* Multiprocessor main */ 112c4762a1bSJed Brown PetscMPIInt tag; 113c4762a1bSJed Brown PetscInt finishedtasks, next_task, checkedin; 114c4762a1bSJed Brown PetscReal f_i = 0.0; 115c4762a1bSJed Brown MPI_Status status; 116c4762a1bSJed Brown 117c4762a1bSJed Brown next_task = 0; 118c4762a1bSJed Brown finishedtasks = 0; 119c4762a1bSJed Brown checkedin = 0; 120c4762a1bSJed Brown 121c4762a1bSJed Brown while (finishedtasks < NOBSERVATIONS || checkedin < user->size - 1) { 1229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&f_i, 1, MPIU_REAL, MPI_ANY_SOURCE, MPI_ANY_TAG, PETSC_COMM_WORLD, &status)); 123c4762a1bSJed Brown if (status.MPI_TAG == IDLE_TAG) { 124c4762a1bSJed Brown checkedin++; 125c4762a1bSJed Brown } else { 126c4762a1bSJed Brown tag = status.MPI_TAG; 127c4762a1bSJed Brown f[tag] = (PetscReal)f_i; 128c4762a1bSJed Brown finishedtasks++; 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 131c4762a1bSJed Brown if (next_task < NOBSERVATIONS) { 1329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(x, NPARAMETERS, MPIU_REAL, status.MPI_SOURCE, next_task, PETSC_COMM_WORLD)); 133c4762a1bSJed Brown next_task++; 134c4762a1bSJed Brown 135c4762a1bSJed Brown } else { 136c4762a1bSJed Brown /* Send idle message */ 1379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(x, NPARAMETERS, MPIU_REAL, status.MPI_SOURCE, IDLE_TAG, PETSC_COMM_WORLD)); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown } 140c4762a1bSJed Brown } 1419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 1429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 143c4762a1bSJed Brown PetscLogFlops(6 * NOBSERVATIONS); 144c4762a1bSJed Brown PetscFunctionReturn(0); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* ------------------------------------------------------------ */ 1489371c9d4SSatish Balay PetscErrorCode FormStartingPoint(Vec X) { 149c4762a1bSJed Brown PetscReal *x; 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscFunctionBegin; 1529566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 153c4762a1bSJed Brown x[0] = 0.15; 154c4762a1bSJed Brown x[1] = 0.008; 155c4762a1bSJed Brown x[2] = 0.010; 1569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 157c4762a1bSJed Brown PetscFunctionReturn(0); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* ---------------------------------------------------------------------- */ 1619371c9d4SSatish Balay PetscErrorCode InitializeData(AppCtx *user) { 162c4762a1bSJed Brown PetscReal *t = user->t, *y = user->y; 163c4762a1bSJed Brown PetscInt i = 0; 164c4762a1bSJed Brown 165c4762a1bSJed Brown PetscFunctionBegin; 1669371c9d4SSatish Balay y[i] = 92.9000; 1679371c9d4SSatish Balay t[i++] = 0.5000; 1689371c9d4SSatish Balay y[i] = 78.7000; 1699371c9d4SSatish Balay t[i++] = 0.6250; 1709371c9d4SSatish Balay y[i] = 64.2000; 1719371c9d4SSatish Balay t[i++] = 0.7500; 1729371c9d4SSatish Balay y[i] = 64.9000; 1739371c9d4SSatish Balay t[i++] = 0.8750; 1749371c9d4SSatish Balay y[i] = 57.1000; 1759371c9d4SSatish Balay t[i++] = 1.0000; 1769371c9d4SSatish Balay y[i] = 43.3000; 1779371c9d4SSatish Balay t[i++] = 1.2500; 1789371c9d4SSatish Balay y[i] = 31.1000; 1799371c9d4SSatish Balay t[i++] = 1.7500; 1809371c9d4SSatish Balay y[i] = 23.6000; 1819371c9d4SSatish Balay t[i++] = 2.2500; 1829371c9d4SSatish Balay y[i] = 31.0500; 1839371c9d4SSatish Balay t[i++] = 1.7500; 1849371c9d4SSatish Balay y[i] = 23.7750; 1859371c9d4SSatish Balay t[i++] = 2.2500; 1869371c9d4SSatish Balay y[i] = 17.7375; 1879371c9d4SSatish Balay t[i++] = 2.7500; 1889371c9d4SSatish Balay y[i] = 13.8000; 1899371c9d4SSatish Balay t[i++] = 3.2500; 1909371c9d4SSatish Balay y[i] = 11.5875; 1919371c9d4SSatish Balay t[i++] = 3.7500; 1929371c9d4SSatish Balay y[i] = 9.4125; 1939371c9d4SSatish Balay t[i++] = 4.2500; 1949371c9d4SSatish Balay y[i] = 7.7250; 1959371c9d4SSatish Balay t[i++] = 4.7500; 1969371c9d4SSatish Balay y[i] = 7.3500; 1979371c9d4SSatish Balay t[i++] = 5.2500; 1989371c9d4SSatish Balay y[i] = 8.0250; 1999371c9d4SSatish Balay t[i++] = 5.7500; 2009371c9d4SSatish Balay y[i] = 90.6000; 2019371c9d4SSatish Balay t[i++] = 0.5000; 2029371c9d4SSatish Balay y[i] = 76.9000; 2039371c9d4SSatish Balay t[i++] = 0.6250; 2049371c9d4SSatish Balay y[i] = 71.6000; 2059371c9d4SSatish Balay t[i++] = 0.7500; 2069371c9d4SSatish Balay y[i] = 63.6000; 2079371c9d4SSatish Balay t[i++] = 0.8750; 2089371c9d4SSatish Balay y[i] = 54.0000; 2099371c9d4SSatish Balay t[i++] = 1.0000; 2109371c9d4SSatish Balay y[i] = 39.2000; 2119371c9d4SSatish Balay t[i++] = 1.2500; 2129371c9d4SSatish Balay y[i] = 29.3000; 2139371c9d4SSatish Balay t[i++] = 1.7500; 2149371c9d4SSatish Balay y[i] = 21.4000; 2159371c9d4SSatish Balay t[i++] = 2.2500; 2169371c9d4SSatish Balay y[i] = 29.1750; 2179371c9d4SSatish Balay t[i++] = 1.7500; 2189371c9d4SSatish Balay y[i] = 22.1250; 2199371c9d4SSatish Balay t[i++] = 2.2500; 2209371c9d4SSatish Balay y[i] = 17.5125; 2219371c9d4SSatish Balay t[i++] = 2.7500; 2229371c9d4SSatish Balay y[i] = 14.2500; 2239371c9d4SSatish Balay t[i++] = 3.2500; 2249371c9d4SSatish Balay y[i] = 9.4500; 2259371c9d4SSatish Balay t[i++] = 3.7500; 2269371c9d4SSatish Balay y[i] = 9.1500; 2279371c9d4SSatish Balay t[i++] = 4.2500; 2289371c9d4SSatish Balay y[i] = 7.9125; 2299371c9d4SSatish Balay t[i++] = 4.7500; 2309371c9d4SSatish Balay y[i] = 8.4750; 2319371c9d4SSatish Balay t[i++] = 5.2500; 2329371c9d4SSatish Balay y[i] = 6.1125; 2339371c9d4SSatish Balay t[i++] = 5.7500; 2349371c9d4SSatish Balay y[i] = 80.0000; 2359371c9d4SSatish Balay t[i++] = 0.5000; 2369371c9d4SSatish Balay y[i] = 79.0000; 2379371c9d4SSatish Balay t[i++] = 0.6250; 2389371c9d4SSatish Balay y[i] = 63.8000; 2399371c9d4SSatish Balay t[i++] = 0.7500; 2409371c9d4SSatish Balay y[i] = 57.2000; 2419371c9d4SSatish Balay t[i++] = 0.8750; 2429371c9d4SSatish Balay y[i] = 53.2000; 2439371c9d4SSatish Balay t[i++] = 1.0000; 2449371c9d4SSatish Balay y[i] = 42.5000; 2459371c9d4SSatish Balay t[i++] = 1.2500; 2469371c9d4SSatish Balay y[i] = 26.8000; 2479371c9d4SSatish Balay t[i++] = 1.7500; 2489371c9d4SSatish Balay y[i] = 20.4000; 2499371c9d4SSatish Balay t[i++] = 2.2500; 2509371c9d4SSatish Balay y[i] = 26.8500; 2519371c9d4SSatish Balay t[i++] = 1.7500; 2529371c9d4SSatish Balay y[i] = 21.0000; 2539371c9d4SSatish Balay t[i++] = 2.2500; 2549371c9d4SSatish Balay y[i] = 16.4625; 2559371c9d4SSatish Balay t[i++] = 2.7500; 2569371c9d4SSatish Balay y[i] = 12.5250; 2579371c9d4SSatish Balay t[i++] = 3.2500; 2589371c9d4SSatish Balay y[i] = 10.5375; 2599371c9d4SSatish Balay t[i++] = 3.7500; 2609371c9d4SSatish Balay y[i] = 8.5875; 2619371c9d4SSatish Balay t[i++] = 4.2500; 2629371c9d4SSatish Balay y[i] = 7.1250; 2639371c9d4SSatish Balay t[i++] = 4.7500; 2649371c9d4SSatish Balay y[i] = 6.1125; 2659371c9d4SSatish Balay t[i++] = 5.2500; 2669371c9d4SSatish Balay y[i] = 5.9625; 2679371c9d4SSatish Balay t[i++] = 5.7500; 2689371c9d4SSatish Balay y[i] = 74.1000; 2699371c9d4SSatish Balay t[i++] = 0.5000; 2709371c9d4SSatish Balay y[i] = 67.3000; 2719371c9d4SSatish Balay t[i++] = 0.6250; 2729371c9d4SSatish Balay y[i] = 60.8000; 2739371c9d4SSatish Balay t[i++] = 0.7500; 2749371c9d4SSatish Balay y[i] = 55.5000; 2759371c9d4SSatish Balay t[i++] = 0.8750; 2769371c9d4SSatish Balay y[i] = 50.3000; 2779371c9d4SSatish Balay t[i++] = 1.0000; 2789371c9d4SSatish Balay y[i] = 41.0000; 2799371c9d4SSatish Balay t[i++] = 1.2500; 2809371c9d4SSatish Balay y[i] = 29.4000; 2819371c9d4SSatish Balay t[i++] = 1.7500; 2829371c9d4SSatish Balay y[i] = 20.4000; 2839371c9d4SSatish Balay t[i++] = 2.2500; 2849371c9d4SSatish Balay y[i] = 29.3625; 2859371c9d4SSatish Balay t[i++] = 1.7500; 2869371c9d4SSatish Balay y[i] = 21.1500; 2879371c9d4SSatish Balay t[i++] = 2.2500; 2889371c9d4SSatish Balay y[i] = 16.7625; 2899371c9d4SSatish Balay t[i++] = 2.7500; 2909371c9d4SSatish Balay y[i] = 13.2000; 2919371c9d4SSatish Balay t[i++] = 3.2500; 2929371c9d4SSatish Balay y[i] = 10.8750; 2939371c9d4SSatish Balay t[i++] = 3.7500; 2949371c9d4SSatish Balay y[i] = 8.1750; 2959371c9d4SSatish Balay t[i++] = 4.2500; 2969371c9d4SSatish Balay y[i] = 7.3500; 2979371c9d4SSatish Balay t[i++] = 4.7500; 2989371c9d4SSatish Balay y[i] = 5.9625; 2999371c9d4SSatish Balay t[i++] = 5.2500; 3009371c9d4SSatish Balay y[i] = 5.6250; 3019371c9d4SSatish Balay t[i++] = 5.7500; 3029371c9d4SSatish Balay y[i] = 81.5000; 3039371c9d4SSatish Balay t[i++] = .5000; 3049371c9d4SSatish Balay y[i] = 62.4000; 3059371c9d4SSatish Balay t[i++] = .7500; 3069371c9d4SSatish Balay y[i] = 32.5000; 3079371c9d4SSatish Balay t[i++] = 1.5000; 3089371c9d4SSatish Balay y[i] = 12.4100; 3099371c9d4SSatish Balay t[i++] = 3.0000; 3109371c9d4SSatish Balay y[i] = 13.1200; 3119371c9d4SSatish Balay t[i++] = 3.0000; 3129371c9d4SSatish Balay y[i] = 15.5600; 3139371c9d4SSatish Balay t[i++] = 3.0000; 3149371c9d4SSatish Balay y[i] = 5.6300; 3159371c9d4SSatish Balay t[i++] = 6.0000; 3169371c9d4SSatish Balay y[i] = 78.0000; 3179371c9d4SSatish Balay t[i++] = .5000; 3189371c9d4SSatish Balay y[i] = 59.9000; 3199371c9d4SSatish Balay t[i++] = .7500; 3209371c9d4SSatish Balay y[i] = 33.2000; 3219371c9d4SSatish Balay t[i++] = 1.5000; 3229371c9d4SSatish Balay y[i] = 13.8400; 3239371c9d4SSatish Balay t[i++] = 3.0000; 3249371c9d4SSatish Balay y[i] = 12.7500; 3259371c9d4SSatish Balay t[i++] = 3.0000; 3269371c9d4SSatish Balay y[i] = 14.6200; 3279371c9d4SSatish Balay t[i++] = 3.0000; 3289371c9d4SSatish Balay y[i] = 3.9400; 3299371c9d4SSatish Balay t[i++] = 6.0000; 3309371c9d4SSatish Balay y[i] = 76.8000; 3319371c9d4SSatish Balay t[i++] = .5000; 3329371c9d4SSatish Balay y[i] = 61.0000; 3339371c9d4SSatish Balay t[i++] = .7500; 3349371c9d4SSatish Balay y[i] = 32.9000; 3359371c9d4SSatish Balay t[i++] = 1.5000; 3369371c9d4SSatish Balay y[i] = 13.8700; 3379371c9d4SSatish Balay t[i++] = 3.0000; 3389371c9d4SSatish Balay y[i] = 11.8100; 3399371c9d4SSatish Balay t[i++] = 3.0000; 3409371c9d4SSatish Balay y[i] = 13.3100; 3419371c9d4SSatish Balay t[i++] = 3.0000; 3429371c9d4SSatish Balay y[i] = 5.4400; 3439371c9d4SSatish Balay t[i++] = 6.0000; 3449371c9d4SSatish Balay y[i] = 78.0000; 3459371c9d4SSatish Balay t[i++] = .5000; 3469371c9d4SSatish Balay y[i] = 63.5000; 3479371c9d4SSatish Balay t[i++] = .7500; 3489371c9d4SSatish Balay y[i] = 33.8000; 3499371c9d4SSatish Balay t[i++] = 1.5000; 3509371c9d4SSatish Balay y[i] = 12.5600; 3519371c9d4SSatish Balay t[i++] = 3.0000; 3529371c9d4SSatish Balay y[i] = 5.6300; 3539371c9d4SSatish Balay t[i++] = 6.0000; 3549371c9d4SSatish Balay y[i] = 12.7500; 3559371c9d4SSatish Balay t[i++] = 3.0000; 3569371c9d4SSatish Balay y[i] = 13.1200; 3579371c9d4SSatish Balay t[i++] = 3.0000; 3589371c9d4SSatish Balay y[i] = 5.4400; 3599371c9d4SSatish Balay t[i++] = 6.0000; 3609371c9d4SSatish Balay y[i] = 76.8000; 3619371c9d4SSatish Balay t[i++] = .5000; 3629371c9d4SSatish Balay y[i] = 60.0000; 3639371c9d4SSatish Balay t[i++] = .7500; 3649371c9d4SSatish Balay y[i] = 47.8000; 3659371c9d4SSatish Balay t[i++] = 1.0000; 3669371c9d4SSatish Balay y[i] = 32.0000; 3679371c9d4SSatish Balay t[i++] = 1.5000; 3689371c9d4SSatish Balay y[i] = 22.2000; 3699371c9d4SSatish Balay t[i++] = 2.0000; 3709371c9d4SSatish Balay y[i] = 22.5700; 3719371c9d4SSatish Balay t[i++] = 2.0000; 3729371c9d4SSatish Balay y[i] = 18.8200; 3739371c9d4SSatish Balay t[i++] = 2.5000; 3749371c9d4SSatish Balay y[i] = 13.9500; 3759371c9d4SSatish Balay t[i++] = 3.0000; 3769371c9d4SSatish Balay y[i] = 11.2500; 3779371c9d4SSatish Balay t[i++] = 4.0000; 3789371c9d4SSatish Balay y[i] = 9.0000; 3799371c9d4SSatish Balay t[i++] = 5.0000; 3809371c9d4SSatish Balay y[i] = 6.6700; 3819371c9d4SSatish Balay t[i++] = 6.0000; 3829371c9d4SSatish Balay y[i] = 75.8000; 3839371c9d4SSatish Balay t[i++] = .5000; 3849371c9d4SSatish Balay y[i] = 62.0000; 3859371c9d4SSatish Balay t[i++] = .7500; 3869371c9d4SSatish Balay y[i] = 48.8000; 3879371c9d4SSatish Balay t[i++] = 1.0000; 3889371c9d4SSatish Balay y[i] = 35.2000; 3899371c9d4SSatish Balay t[i++] = 1.5000; 3909371c9d4SSatish Balay y[i] = 20.0000; 3919371c9d4SSatish Balay t[i++] = 2.0000; 3929371c9d4SSatish Balay y[i] = 20.3200; 3939371c9d4SSatish Balay t[i++] = 2.0000; 3949371c9d4SSatish Balay y[i] = 19.3100; 3959371c9d4SSatish Balay t[i++] = 2.5000; 3969371c9d4SSatish Balay y[i] = 12.7500; 3979371c9d4SSatish Balay t[i++] = 3.0000; 3989371c9d4SSatish Balay y[i] = 10.4200; 3999371c9d4SSatish Balay t[i++] = 4.0000; 4009371c9d4SSatish Balay y[i] = 7.3100; 4019371c9d4SSatish Balay t[i++] = 5.0000; 4029371c9d4SSatish Balay y[i] = 7.4200; 4039371c9d4SSatish Balay t[i++] = 6.0000; 4049371c9d4SSatish Balay y[i] = 70.5000; 4059371c9d4SSatish Balay t[i++] = .5000; 4069371c9d4SSatish Balay y[i] = 59.5000; 4079371c9d4SSatish Balay t[i++] = .7500; 4089371c9d4SSatish Balay y[i] = 48.5000; 4099371c9d4SSatish Balay t[i++] = 1.0000; 4109371c9d4SSatish Balay y[i] = 35.8000; 4119371c9d4SSatish Balay t[i++] = 1.5000; 4129371c9d4SSatish Balay y[i] = 21.0000; 4139371c9d4SSatish Balay t[i++] = 2.0000; 4149371c9d4SSatish Balay y[i] = 21.6700; 4159371c9d4SSatish Balay t[i++] = 2.0000; 4169371c9d4SSatish Balay y[i] = 21.0000; 4179371c9d4SSatish Balay t[i++] = 2.5000; 4189371c9d4SSatish Balay y[i] = 15.6400; 4199371c9d4SSatish Balay t[i++] = 3.0000; 4209371c9d4SSatish Balay y[i] = 8.1700; 4219371c9d4SSatish Balay t[i++] = 4.0000; 4229371c9d4SSatish Balay y[i] = 8.5500; 4239371c9d4SSatish Balay t[i++] = 5.0000; 4249371c9d4SSatish Balay y[i] = 10.1200; 4259371c9d4SSatish Balay t[i++] = 6.0000; 4269371c9d4SSatish Balay y[i] = 78.0000; 4279371c9d4SSatish Balay t[i++] = .5000; 4289371c9d4SSatish Balay y[i] = 66.0000; 4299371c9d4SSatish Balay t[i++] = .6250; 4309371c9d4SSatish Balay y[i] = 62.0000; 4319371c9d4SSatish Balay t[i++] = .7500; 4329371c9d4SSatish Balay y[i] = 58.0000; 4339371c9d4SSatish Balay t[i++] = .8750; 4349371c9d4SSatish Balay y[i] = 47.7000; 4359371c9d4SSatish Balay t[i++] = 1.0000; 4369371c9d4SSatish Balay y[i] = 37.8000; 4379371c9d4SSatish Balay t[i++] = 1.2500; 4389371c9d4SSatish Balay y[i] = 20.2000; 4399371c9d4SSatish Balay t[i++] = 2.2500; 4409371c9d4SSatish Balay y[i] = 21.0700; 4419371c9d4SSatish Balay t[i++] = 2.2500; 4429371c9d4SSatish Balay y[i] = 13.8700; 4439371c9d4SSatish Balay t[i++] = 2.7500; 4449371c9d4SSatish Balay y[i] = 9.6700; 4459371c9d4SSatish Balay t[i++] = 3.2500; 4469371c9d4SSatish Balay y[i] = 7.7600; 4479371c9d4SSatish Balay t[i++] = 3.7500; 4489371c9d4SSatish Balay y[i] = 5.4400; 4499371c9d4SSatish Balay t[i++] = 4.2500; 4509371c9d4SSatish Balay y[i] = 4.8700; 4519371c9d4SSatish Balay t[i++] = 4.7500; 4529371c9d4SSatish Balay y[i] = 4.0100; 4539371c9d4SSatish Balay t[i++] = 5.2500; 4549371c9d4SSatish Balay y[i] = 3.7500; 4559371c9d4SSatish Balay t[i++] = 5.7500; 4569371c9d4SSatish Balay y[i] = 24.1900; 4579371c9d4SSatish Balay t[i++] = 3.0000; 4589371c9d4SSatish Balay y[i] = 25.7600; 4599371c9d4SSatish Balay t[i++] = 3.0000; 4609371c9d4SSatish Balay y[i] = 18.0700; 4619371c9d4SSatish Balay t[i++] = 3.0000; 4629371c9d4SSatish Balay y[i] = 11.8100; 4639371c9d4SSatish Balay t[i++] = 3.0000; 4649371c9d4SSatish Balay y[i] = 12.0700; 4659371c9d4SSatish Balay t[i++] = 3.0000; 4669371c9d4SSatish Balay y[i] = 16.1200; 4679371c9d4SSatish Balay t[i++] = 3.0000; 4689371c9d4SSatish Balay y[i] = 70.8000; 4699371c9d4SSatish Balay t[i++] = .5000; 4709371c9d4SSatish Balay y[i] = 54.7000; 4719371c9d4SSatish Balay t[i++] = .7500; 4729371c9d4SSatish Balay y[i] = 48.0000; 4739371c9d4SSatish Balay t[i++] = 1.0000; 4749371c9d4SSatish Balay y[i] = 39.8000; 4759371c9d4SSatish Balay t[i++] = 1.5000; 4769371c9d4SSatish Balay y[i] = 29.8000; 4779371c9d4SSatish Balay t[i++] = 2.0000; 4789371c9d4SSatish Balay y[i] = 23.7000; 4799371c9d4SSatish Balay t[i++] = 2.5000; 4809371c9d4SSatish Balay y[i] = 29.6200; 4819371c9d4SSatish Balay t[i++] = 2.0000; 4829371c9d4SSatish Balay y[i] = 23.8100; 4839371c9d4SSatish Balay t[i++] = 2.5000; 4849371c9d4SSatish Balay y[i] = 17.7000; 4859371c9d4SSatish Balay t[i++] = 3.0000; 4869371c9d4SSatish Balay y[i] = 11.5500; 4879371c9d4SSatish Balay t[i++] = 4.0000; 4889371c9d4SSatish Balay y[i] = 12.0700; 4899371c9d4SSatish Balay t[i++] = 5.0000; 4909371c9d4SSatish Balay y[i] = 8.7400; 4919371c9d4SSatish Balay t[i++] = 6.0000; 4929371c9d4SSatish Balay y[i] = 80.7000; 4939371c9d4SSatish Balay t[i++] = .5000; 4949371c9d4SSatish Balay y[i] = 61.3000; 4959371c9d4SSatish Balay t[i++] = .7500; 4969371c9d4SSatish Balay y[i] = 47.5000; 4979371c9d4SSatish Balay t[i++] = 1.0000; 4989371c9d4SSatish Balay y[i] = 29.0000; 4999371c9d4SSatish Balay t[i++] = 1.5000; 5009371c9d4SSatish Balay y[i] = 24.0000; 5019371c9d4SSatish Balay t[i++] = 2.0000; 5029371c9d4SSatish Balay y[i] = 17.7000; 5039371c9d4SSatish Balay t[i++] = 2.5000; 5049371c9d4SSatish Balay y[i] = 24.5600; 5059371c9d4SSatish Balay t[i++] = 2.0000; 5069371c9d4SSatish Balay y[i] = 18.6700; 5079371c9d4SSatish Balay t[i++] = 2.5000; 5089371c9d4SSatish Balay y[i] = 16.2400; 5099371c9d4SSatish Balay t[i++] = 3.0000; 5109371c9d4SSatish Balay y[i] = 8.7400; 5119371c9d4SSatish Balay t[i++] = 4.0000; 5129371c9d4SSatish Balay y[i] = 7.8700; 5139371c9d4SSatish Balay t[i++] = 5.0000; 5149371c9d4SSatish Balay y[i] = 8.5100; 5159371c9d4SSatish Balay t[i++] = 6.0000; 5169371c9d4SSatish Balay y[i] = 66.7000; 5179371c9d4SSatish Balay t[i++] = .5000; 5189371c9d4SSatish Balay y[i] = 59.2000; 5199371c9d4SSatish Balay t[i++] = .7500; 5209371c9d4SSatish Balay y[i] = 40.8000; 5219371c9d4SSatish Balay t[i++] = 1.0000; 5229371c9d4SSatish Balay y[i] = 30.7000; 5239371c9d4SSatish Balay t[i++] = 1.5000; 5249371c9d4SSatish Balay y[i] = 25.7000; 5259371c9d4SSatish Balay t[i++] = 2.0000; 5269371c9d4SSatish Balay y[i] = 16.3000; 5279371c9d4SSatish Balay t[i++] = 2.5000; 5289371c9d4SSatish Balay y[i] = 25.9900; 5299371c9d4SSatish Balay t[i++] = 2.0000; 5309371c9d4SSatish Balay y[i] = 16.9500; 5319371c9d4SSatish Balay t[i++] = 2.5000; 5329371c9d4SSatish Balay y[i] = 13.3500; 5339371c9d4SSatish Balay t[i++] = 3.0000; 5349371c9d4SSatish Balay y[i] = 8.6200; 5359371c9d4SSatish Balay t[i++] = 4.0000; 5369371c9d4SSatish Balay y[i] = 7.2000; 5379371c9d4SSatish Balay t[i++] = 5.0000; 5389371c9d4SSatish Balay y[i] = 6.6400; 5399371c9d4SSatish Balay t[i++] = 6.0000; 5409371c9d4SSatish Balay y[i] = 13.6900; 5419371c9d4SSatish Balay t[i++] = 3.0000; 5429371c9d4SSatish Balay y[i] = 81.0000; 5439371c9d4SSatish Balay t[i++] = .5000; 5449371c9d4SSatish Balay y[i] = 64.5000; 5459371c9d4SSatish Balay t[i++] = .7500; 5469371c9d4SSatish Balay y[i] = 35.5000; 5479371c9d4SSatish Balay t[i++] = 1.5000; 5489371c9d4SSatish Balay y[i] = 13.3100; 5499371c9d4SSatish Balay t[i++] = 3.0000; 5509371c9d4SSatish Balay y[i] = 4.8700; 5519371c9d4SSatish Balay t[i++] = 6.0000; 5529371c9d4SSatish Balay y[i] = 12.9400; 5539371c9d4SSatish Balay t[i++] = 3.0000; 5549371c9d4SSatish Balay y[i] = 5.0600; 5559371c9d4SSatish Balay t[i++] = 6.0000; 5569371c9d4SSatish Balay y[i] = 15.1900; 5579371c9d4SSatish Balay t[i++] = 3.0000; 5589371c9d4SSatish Balay y[i] = 14.6200; 5599371c9d4SSatish Balay t[i++] = 3.0000; 5609371c9d4SSatish Balay y[i] = 15.6400; 5619371c9d4SSatish Balay t[i++] = 3.0000; 5629371c9d4SSatish Balay y[i] = 25.5000; 5639371c9d4SSatish Balay t[i++] = 1.7500; 5649371c9d4SSatish Balay y[i] = 25.9500; 5659371c9d4SSatish Balay t[i++] = 1.7500; 5669371c9d4SSatish Balay y[i] = 81.7000; 5679371c9d4SSatish Balay t[i++] = .5000; 5689371c9d4SSatish Balay y[i] = 61.6000; 5699371c9d4SSatish Balay t[i++] = .7500; 5709371c9d4SSatish Balay y[i] = 29.8000; 5719371c9d4SSatish Balay t[i++] = 1.7500; 5729371c9d4SSatish Balay y[i] = 29.8100; 5739371c9d4SSatish Balay t[i++] = 1.7500; 5749371c9d4SSatish Balay y[i] = 17.1700; 5759371c9d4SSatish Balay t[i++] = 2.7500; 5769371c9d4SSatish Balay y[i] = 10.3900; 5779371c9d4SSatish Balay t[i++] = 3.7500; 5789371c9d4SSatish Balay y[i] = 28.4000; 5799371c9d4SSatish Balay t[i++] = 1.7500; 5809371c9d4SSatish Balay y[i] = 28.6900; 5819371c9d4SSatish Balay t[i++] = 1.7500; 5829371c9d4SSatish Balay y[i] = 81.3000; 5839371c9d4SSatish Balay t[i++] = .5000; 5849371c9d4SSatish Balay y[i] = 60.9000; 5859371c9d4SSatish Balay t[i++] = .7500; 5869371c9d4SSatish Balay y[i] = 16.6500; 5879371c9d4SSatish Balay t[i++] = 2.7500; 5889371c9d4SSatish Balay y[i] = 10.0500; 5899371c9d4SSatish Balay t[i++] = 3.7500; 5909371c9d4SSatish Balay y[i] = 28.9000; 5919371c9d4SSatish Balay t[i++] = 1.7500; 5929371c9d4SSatish Balay y[i] = 28.9500; 5939371c9d4SSatish Balay t[i++] = 1.7500; 594c4762a1bSJed Brown PetscFunctionReturn(0); 595c4762a1bSJed Brown } 596c4762a1bSJed Brown 5979371c9d4SSatish Balay PetscErrorCode TaskWorker(AppCtx *user) { 598c4762a1bSJed Brown PetscReal x[NPARAMETERS], f = 0.0; 599c4762a1bSJed Brown PetscMPIInt tag = IDLE_TAG; 600c4762a1bSJed Brown PetscInt index; 601c4762a1bSJed Brown MPI_Status status; 602c4762a1bSJed Brown 603c4762a1bSJed Brown PetscFunctionBegin; 6049dddd249SSatish Balay /* Send check-in message to rank-0 */ 605c4762a1bSJed Brown 6069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&f, 1, MPIU_REAL, 0, IDLE_TAG, PETSC_COMM_WORLD)); 607c4762a1bSJed Brown while (tag != DIE_TAG) { 6089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(x, NPARAMETERS, MPIU_REAL, 0, MPI_ANY_TAG, PETSC_COMM_WORLD, &status)); 609c4762a1bSJed Brown tag = status.MPI_TAG; 610c4762a1bSJed Brown if (tag == IDLE_TAG) { 6119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&f, 1, MPIU_REAL, 0, IDLE_TAG, PETSC_COMM_WORLD)); 612c4762a1bSJed Brown } else if (tag != DIE_TAG) { 613c4762a1bSJed Brown index = (PetscInt)tag; 6149566063dSJacob Faibussowitsch PetscCall(RunSimulation(x, index, &f, user)); 6159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&f, 1, MPIU_REAL, 0, tag, PETSC_COMM_WORLD)); 616c4762a1bSJed Brown } 617c4762a1bSJed Brown } 618c4762a1bSJed Brown PetscFunctionReturn(0); 619c4762a1bSJed Brown } 620c4762a1bSJed Brown 6219371c9d4SSatish Balay PetscErrorCode RunSimulation(PetscReal *x, PetscInt i, PetscReal *f, AppCtx *user) { 622c4762a1bSJed Brown PetscReal *t = user->t; 623c4762a1bSJed Brown PetscReal *y = user->y; 624c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) 625e1dfdf8eSBarry 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 */ 626c4762a1bSJed Brown #else 627c4762a1bSJed Brown *f = y[i] - PetscExpScalar(-x[0] * t[i]) / (x[1] + x[2] * t[i]); 628c4762a1bSJed Brown #endif 629c4762a1bSJed Brown return (0); 630c4762a1bSJed Brown } 631c4762a1bSJed Brown 6329371c9d4SSatish Balay PetscErrorCode StopWorkers(AppCtx *user) { 633c4762a1bSJed Brown PetscInt checkedin; 634c4762a1bSJed Brown MPI_Status status; 635c4762a1bSJed Brown PetscReal f, x[NPARAMETERS]; 636c4762a1bSJed Brown 637c4762a1bSJed Brown PetscFunctionBegin; 638c4762a1bSJed Brown checkedin = 0; 639c4762a1bSJed Brown while (checkedin < user->size - 1) { 6409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&f, 1, MPIU_REAL, MPI_ANY_SOURCE, MPI_ANY_TAG, PETSC_COMM_WORLD, &status)); 641c4762a1bSJed Brown checkedin++; 6429566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(x, NPARAMETERS)); 6439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(x, NPARAMETERS, MPIU_REAL, status.MPI_SOURCE, DIE_TAG, PETSC_COMM_WORLD)); 644c4762a1bSJed Brown } 645c4762a1bSJed Brown PetscFunctionReturn(0); 646c4762a1bSJed Brown } 647c4762a1bSJed Brown 648c4762a1bSJed Brown /*TEST 649c4762a1bSJed Brown 650c4762a1bSJed Brown build: 651c4762a1bSJed Brown requires: !complex 652c4762a1bSJed Brown 653c4762a1bSJed Brown test: 654c4762a1bSJed Brown nsize: 3 655c4762a1bSJed Brown requires: !single 656c4762a1bSJed Brown args: -tao_smonitor -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-5 657c4762a1bSJed Brown 658c4762a1bSJed Brown TEST*/ 659