1c4762a1bSJed Brown static char help[] = "Solves a simple data assimilation problem with one dimensional advection diffusion equation using TSAdjoint\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown
5c4762a1bSJed Brown Not yet tested in parallel
6c4762a1bSJed Brown
7c4762a1bSJed Brown */
8c4762a1bSJed Brown
9c4762a1bSJed Brown /* ------------------------------------------------------------------------
10c4762a1bSJed Brown
11c4762a1bSJed Brown This program uses the one-dimensional advection-diffusion equation),
12c4762a1bSJed Brown u_t = mu*u_xx - a u_x,
13c4762a1bSJed Brown on the domain 0 <= x <= 1, with periodic boundary conditions
14c4762a1bSJed Brown
15c4762a1bSJed Brown to demonstrate solving a data assimilation problem of finding the initial conditions
16c4762a1bSJed Brown to produce a given solution at a fixed time.
17c4762a1bSJed Brown
18c4762a1bSJed Brown The operators are discretized with the spectral element method
19c4762a1bSJed Brown
20c4762a1bSJed Brown ------------------------------------------------------------------------- */
21c4762a1bSJed Brown
22c4762a1bSJed Brown /*
23c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this file
24c4762a1bSJed Brown automatically includes:
25c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors
26c4762a1bSJed Brown petscmat.h - matrices
27c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods
28c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners
29c4762a1bSJed Brown petscksp.h - linear solvers petscsnes.h - nonlinear solvers
30c4762a1bSJed Brown */
31c4762a1bSJed Brown
32c4762a1bSJed Brown #include <petsctao.h>
33c4762a1bSJed Brown #include <petscts.h>
34c4762a1bSJed Brown #include <petscdt.h>
35c4762a1bSJed Brown #include <petscdraw.h>
36c4762a1bSJed Brown #include <petscdmda.h>
37c4762a1bSJed Brown
38c4762a1bSJed Brown /*
39c4762a1bSJed Brown User-defined application context - contains data needed by the
40c4762a1bSJed Brown application-provided call-back routines.
41c4762a1bSJed Brown */
42c4762a1bSJed Brown
43c4762a1bSJed Brown typedef struct {
44c4762a1bSJed Brown PetscInt n; /* number of nodes */
45c4762a1bSJed Brown PetscReal *nodes; /* GLL nodes */
46c4762a1bSJed Brown PetscReal *weights; /* GLL weights */
47c4762a1bSJed Brown } PetscGLL;
48c4762a1bSJed Brown
49c4762a1bSJed Brown typedef struct {
50c4762a1bSJed Brown PetscInt N; /* grid points per elements*/
51c4762a1bSJed Brown PetscInt E; /* number of elements */
52c4762a1bSJed Brown PetscReal tol_L2, tol_max; /* error norms */
53c4762a1bSJed Brown PetscInt steps; /* number of timesteps */
54c4762a1bSJed Brown PetscReal Tend; /* endtime */
55c4762a1bSJed Brown PetscReal mu; /* viscosity */
56c4762a1bSJed Brown PetscReal a; /* advection speed */
57c4762a1bSJed Brown PetscReal L; /* total length of domain */
58c4762a1bSJed Brown PetscReal Le;
59c4762a1bSJed Brown PetscReal Tadj;
60c4762a1bSJed Brown } PetscParam;
61c4762a1bSJed Brown
62c4762a1bSJed Brown typedef struct {
63c4762a1bSJed Brown Vec reference; /* desired end state */
64c4762a1bSJed Brown Vec grid; /* total grid */
65c4762a1bSJed Brown Vec grad;
66c4762a1bSJed Brown Vec ic;
67c4762a1bSJed Brown Vec curr_sol;
68c4762a1bSJed Brown Vec joe;
69c4762a1bSJed Brown Vec true_solution; /* actual initial conditions for the final solution */
70c4762a1bSJed Brown } PetscData;
71c4762a1bSJed Brown
72c4762a1bSJed Brown typedef struct {
73c4762a1bSJed Brown Vec grid; /* total grid */
74c4762a1bSJed Brown Vec mass; /* mass matrix for total integration */
75943f6d28SPierre Jolivet Mat stiff; /* stiffness matrix */
76c4762a1bSJed Brown Mat advec;
77c4762a1bSJed Brown Mat keptstiff;
78c4762a1bSJed Brown PetscGLL gll;
79c4762a1bSJed Brown } PetscSEMOperators;
80c4762a1bSJed Brown
81c4762a1bSJed Brown typedef struct {
82c4762a1bSJed Brown DM da; /* distributed array data structure */
83c4762a1bSJed Brown PetscSEMOperators SEMop;
84c4762a1bSJed Brown PetscParam param;
85c4762a1bSJed Brown PetscData dat;
86c4762a1bSJed Brown TS ts;
87c4762a1bSJed Brown PetscReal initial_dt;
88c4762a1bSJed Brown PetscReal *solutioncoefficients;
89c4762a1bSJed Brown PetscInt ncoeff;
90c4762a1bSJed Brown } AppCtx;
91c4762a1bSJed Brown
92c4762a1bSJed Brown /*
93c4762a1bSJed Brown User-defined routines
94c4762a1bSJed Brown */
95c4762a1bSJed Brown extern PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
96c4762a1bSJed Brown extern PetscErrorCode RHSLaplacian(TS, PetscReal, Vec, Mat, Mat, void *);
97c4762a1bSJed Brown extern PetscErrorCode RHSAdvection(TS, PetscReal, Vec, Mat, Mat, void *);
98c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *);
99c4762a1bSJed Brown extern PetscErrorCode ComputeReference(TS, PetscReal, Vec, AppCtx *);
100c4762a1bSJed Brown extern PetscErrorCode MonitorError(Tao, void *);
101*2a8381b2SBarry Smith extern PetscErrorCode MonitorDestroy(PetscCtxRt);
102c4762a1bSJed Brown extern PetscErrorCode ComputeSolutionCoefficients(AppCtx *);
103c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
104c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
105c4762a1bSJed Brown
main(int argc,char ** argv)106d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
107d71ae5a4SJacob Faibussowitsch {
108c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */
109c4762a1bSJed Brown Tao tao;
110c4762a1bSJed Brown Vec u; /* approximate solution vector */
111c4762a1bSJed Brown PetscInt i, xs, xm, ind, j, lenglob;
112c4762a1bSJed Brown PetscReal x, *wrk_ptr1, *wrk_ptr2;
113c4762a1bSJed Brown MatNullSpace nsp;
114c4762a1bSJed Brown
115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116c4762a1bSJed Brown Initialize program and set problem parameters
117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118327415f7SBarry Smith PetscFunctionBeginUser;
119c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
120c4762a1bSJed Brown
121c4762a1bSJed Brown /*initialize parameters */
122c4762a1bSJed Brown appctx.param.N = 10; /* order of the spectral element */
123c4762a1bSJed Brown appctx.param.E = 8; /* number of elements */
124c4762a1bSJed Brown appctx.param.L = 1.0; /* length of the domain */
125c4762a1bSJed Brown appctx.param.mu = 0.00001; /* diffusion coefficient */
126c4762a1bSJed Brown appctx.param.a = 0.0; /* advection speed */
127c4762a1bSJed Brown appctx.initial_dt = 1e-4;
1281690c2aeSBarry Smith appctx.param.steps = PETSC_INT_MAX;
129c4762a1bSJed Brown appctx.param.Tend = 0.01;
130c4762a1bSJed Brown appctx.ncoeff = 2;
131c4762a1bSJed Brown
1329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &appctx.param.N, NULL));
1339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-E", &appctx.param.E, NULL));
1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ncoeff", &appctx.ncoeff, NULL));
1359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-Tend", &appctx.param.Tend, NULL));
1369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &appctx.param.mu, NULL));
1379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-a", &appctx.param.a, NULL));
138c4762a1bSJed Brown appctx.param.Le = appctx.param.L / appctx.param.E;
139c4762a1bSJed Brown
140c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141c4762a1bSJed Brown Create GLL data structures
142c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1439566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(appctx.param.N, &appctx.SEMop.gll.nodes, appctx.param.N, &appctx.SEMop.gll.weights));
1449566063dSJacob Faibussowitsch PetscCall(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N, PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights));
145c4762a1bSJed Brown appctx.SEMop.gll.n = appctx.param.N;
146c4762a1bSJed Brown lenglob = appctx.param.E * (appctx.param.N - 1);
147c4762a1bSJed Brown
148c4762a1bSJed Brown /*
149c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
150c4762a1bSJed Brown and to set up the ghost point communication pattern. There are E*(Nl-1)+1
151c4762a1bSJed Brown total grid values spread equally among all the processors, except first and last
152c4762a1bSJed Brown */
153c4762a1bSJed Brown
1549566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, lenglob, 1, 1, NULL, &appctx.da));
1559566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(appctx.da));
1569566063dSJacob Faibussowitsch PetscCall(DMSetUp(appctx.da));
157c4762a1bSJed Brown
158c4762a1bSJed Brown /*
159c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the
160c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that
161c4762a1bSJed Brown have the same types.
162c4762a1bSJed Brown */
163c4762a1bSJed Brown
1649566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(appctx.da, &u));
1659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &appctx.dat.ic));
1669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &appctx.dat.true_solution));
1679566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &appctx.dat.reference));
1689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &appctx.SEMop.grid));
1699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &appctx.SEMop.mass));
1709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &appctx.dat.curr_sol));
1719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &appctx.dat.joe));
172c4762a1bSJed Brown
1739566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx.da, &xs, NULL, NULL, &xm, NULL, NULL));
1749566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1));
1759566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2));
176c4762a1bSJed Brown
177c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
178c4762a1bSJed Brown
179c4762a1bSJed Brown xs = xs / (appctx.param.N - 1);
180c4762a1bSJed Brown xm = xm / (appctx.param.N - 1);
181c4762a1bSJed Brown
182c4762a1bSJed Brown /*
183c4762a1bSJed Brown Build total grid and mass over entire mesh (multi-elemental)
184c4762a1bSJed Brown */
185c4762a1bSJed Brown
186c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
187c4762a1bSJed Brown for (j = 0; j < appctx.param.N - 1; j++) {
188c4762a1bSJed Brown x = (appctx.param.Le / 2.0) * (appctx.SEMop.gll.nodes[j] + 1.0) + appctx.param.Le * i;
189c4762a1bSJed Brown ind = i * (appctx.param.N - 1) + j;
190c4762a1bSJed Brown wrk_ptr1[ind] = x;
191c4762a1bSJed Brown wrk_ptr2[ind] = .5 * appctx.param.Le * appctx.SEMop.gll.weights[j];
192c4762a1bSJed Brown if (j == 0) wrk_ptr2[ind] += .5 * appctx.param.Le * appctx.SEMop.gll.weights[j];
193c4762a1bSJed Brown }
194c4762a1bSJed Brown }
1959566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1));
1969566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2));
197c4762a1bSJed Brown
198c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine.
200c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2019566063dSJacob Faibussowitsch PetscCall(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE));
2029566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.stiff));
2039566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.advec));
204c4762a1bSJed Brown
205c4762a1bSJed Brown /*
206c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation
207dd8e379bSPierre Jolivet u_t = f(u,t), the user provides the discretized right-hand side
208c4762a1bSJed Brown as a time-dependent matrix.
209c4762a1bSJed Brown */
2109566063dSJacob Faibussowitsch PetscCall(RHSLaplacian(appctx.ts, 0.0, u, appctx.SEMop.stiff, appctx.SEMop.stiff, &appctx));
2119566063dSJacob Faibussowitsch PetscCall(RHSAdvection(appctx.ts, 0.0, u, appctx.SEMop.advec, appctx.SEMop.advec, &appctx));
2129566063dSJacob Faibussowitsch PetscCall(MatAXPY(appctx.SEMop.stiff, -1.0, appctx.SEMop.advec, DIFFERENT_NONZERO_PATTERN));
2139566063dSJacob Faibussowitsch PetscCall(MatDuplicate(appctx.SEMop.stiff, MAT_COPY_VALUES, &appctx.SEMop.keptstiff));
214c4762a1bSJed Brown
215c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */
2169566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp));
2179566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.stiff, nsp));
2189566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.stiff, NULL));
2199566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp));
220c4762a1bSJed Brown
221c4762a1bSJed Brown /* Create the TS solver that solves the ODE and its adjoint; set its options */
2229566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &appctx.ts));
2239566063dSJacob Faibussowitsch PetscCall(TSSetSolutionFunction(appctx.ts, (PetscErrorCode (*)(TS, PetscReal, Vec, void *))ComputeReference, &appctx));
2249566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(appctx.ts, TS_LINEAR));
2259566063dSJacob Faibussowitsch PetscCall(TSSetType(appctx.ts, TSRK));
2269566063dSJacob Faibussowitsch PetscCall(TSSetDM(appctx.ts, appctx.da));
2279566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx.ts, 0.0));
2289566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(appctx.ts, appctx.initial_dt));
2299566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(appctx.ts, appctx.param.steps));
2309566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(appctx.ts, appctx.param.Tend));
2319566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(appctx.ts, TS_EXACTFINALTIME_MATCHSTEP));
2329566063dSJacob Faibussowitsch PetscCall(TSSetTolerances(appctx.ts, 1e-7, NULL, 1e-7, NULL));
2339566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx.ts));
234188af4bfSBarry Smith /* Need to save initial timestep user may have set with -ts_time_step so it can be reset for each new TSSolve() */
2359566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(appctx.ts, &appctx.initial_dt));
2369566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(appctx.ts, NULL, TSComputeRHSFunctionLinear, &appctx));
2379566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(appctx.ts, appctx.SEMop.stiff, appctx.SEMop.stiff, TSComputeRHSJacobianConstant, &appctx));
2389566063dSJacob Faibussowitsch /* PetscCall(TSSetRHSFunction(appctx.ts,NULL,RHSFunction,&appctx));
2399566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,RHSJacobian,&appctx)); */
240c4762a1bSJed Brown
241c4762a1bSJed Brown /* Set random initial conditions as initial guess, compute analytic reference solution and analytic (true) initial conditions */
2429566063dSJacob Faibussowitsch PetscCall(ComputeSolutionCoefficients(&appctx));
2439566063dSJacob Faibussowitsch PetscCall(InitialConditions(appctx.dat.ic, &appctx));
2449566063dSJacob Faibussowitsch PetscCall(ComputeReference(appctx.ts, appctx.param.Tend, appctx.dat.reference, &appctx));
2459566063dSJacob Faibussowitsch PetscCall(ComputeReference(appctx.ts, 0.0, appctx.dat.true_solution, &appctx));
246c4762a1bSJed Brown
247f32d6360SSatish Balay /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */
2489566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(appctx.ts));
2499566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx.ts));
250f32d6360SSatish Balay
251c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
2529566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
25310978b7dSBarry Smith PetscCall(TaoMonitorSet(tao, MonitorError, &appctx, MonitorDestroy));
2549566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBQNLS));
2559566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, appctx.dat.ic));
256c4762a1bSJed Brown /* Set routine for function and gradient evaluation */
2579566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&appctx));
258c4762a1bSJed Brown /* Check for any TAO command line options */
259606f75f6SBarry Smith PetscCall(TaoSetTolerances(tao, 1e-8, PETSC_CURRENT, PETSC_CURRENT));
2609566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
2619566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
262c4762a1bSJed Brown
2639566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
2649566063dSJacob Faibussowitsch PetscCall(PetscFree(appctx.solutioncoefficients));
2659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.advec));
2669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.stiff));
2679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.keptstiff));
2689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u));
2699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.ic));
2709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.joe));
2719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.true_solution));
2729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.reference));
2739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.SEMop.grid));
2749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.SEMop.mass));
2759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.curr_sol));
2769566063dSJacob Faibussowitsch PetscCall(PetscFree2(appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights));
2779566063dSJacob Faibussowitsch PetscCall(DMDestroy(&appctx.da));
2789566063dSJacob Faibussowitsch PetscCall(TSDestroy(&appctx.ts));
279c4762a1bSJed Brown
280c4762a1bSJed Brown /*
281c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine
282c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI
283c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime
284d75802c7SJacob Faibussowitsch options are chosen (e.g., -log_view).
285c4762a1bSJed Brown */
2869566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
287b122ec5aSJacob Faibussowitsch return 0;
288c4762a1bSJed Brown }
289c4762a1bSJed Brown
290c4762a1bSJed Brown /*
291c4762a1bSJed Brown Computes the coefficients for the analytic solution to the PDE
292c4762a1bSJed Brown */
ComputeSolutionCoefficients(AppCtx * appctx)293d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeSolutionCoefficients(AppCtx *appctx)
294d71ae5a4SJacob Faibussowitsch {
295c4762a1bSJed Brown PetscRandom rand;
296c4762a1bSJed Brown PetscInt i;
297c4762a1bSJed Brown
298c4762a1bSJed Brown PetscFunctionBegin;
2999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->ncoeff, &appctx->solutioncoefficients));
3009566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
3019566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, .9, 1.0));
30248a46eb9SPierre Jolivet for (i = 0; i < appctx->ncoeff; i++) PetscCall(PetscRandomGetValue(rand, &appctx->solutioncoefficients[i]));
3039566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand));
3043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
305c4762a1bSJed Brown }
306c4762a1bSJed Brown
307c4762a1bSJed Brown /* --------------------------------------------------------------------- */
308c4762a1bSJed Brown /*
309c4762a1bSJed Brown InitialConditions - Computes the (random) initial conditions for the Tao optimization solve (these are also initial conditions for the first TSSolve()
310c4762a1bSJed Brown
311c4762a1bSJed Brown Input Parameter:
312c4762a1bSJed Brown u - uninitialized solution vector (global)
313c4762a1bSJed Brown appctx - user-defined application context
314c4762a1bSJed Brown
315c4762a1bSJed Brown Output Parameter:
316c4762a1bSJed Brown u - vector with solution at initial time (global)
317c4762a1bSJed Brown */
InitialConditions(Vec u,AppCtx * appctx)318d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(Vec u, AppCtx *appctx)
319d71ae5a4SJacob Faibussowitsch {
320c4762a1bSJed Brown PetscScalar *s;
321c4762a1bSJed Brown const PetscScalar *xg;
322c4762a1bSJed Brown PetscInt i, j, lenglob;
323c4762a1bSJed Brown PetscReal sum, val;
324c4762a1bSJed Brown PetscRandom rand;
325c4762a1bSJed Brown
326c4762a1bSJed Brown PetscFunctionBegin;
3279566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
3289566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, .9, 1.0));
3299566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, u, &s));
3309566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
331c4762a1bSJed Brown lenglob = appctx->param.E * (appctx->param.N - 1);
332c4762a1bSJed Brown for (i = 0; i < lenglob; i++) {
333c4762a1bSJed Brown s[i] = 0;
334c4762a1bSJed Brown for (j = 0; j < appctx->ncoeff; j++) {
3359566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &val));
336c4762a1bSJed Brown s[i] += val * PetscSinScalar(2 * (j + 1) * PETSC_PI * xg[i]);
337c4762a1bSJed Brown }
338c4762a1bSJed Brown }
3399566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand));
3409566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, u, &s));
3419566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
342c4762a1bSJed Brown /* make sure initial conditions do not contain the constant functions, since with periodic boundary conditions the constant functions introduce a null space */
3439566063dSJacob Faibussowitsch PetscCall(VecSum(u, &sum));
3449566063dSJacob Faibussowitsch PetscCall(VecShift(u, -sum / lenglob));
3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
346c4762a1bSJed Brown }
347c4762a1bSJed Brown
348c4762a1bSJed Brown /*
349c4762a1bSJed Brown TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function.
350c4762a1bSJed Brown
351a5b23f4aSJose E. Roman InitialConditions() computes the initial conditions for the beginning of the Tao iterations
352c4762a1bSJed Brown
353c4762a1bSJed Brown Input Parameter:
354c4762a1bSJed Brown u - uninitialized solution vector (global)
355c4762a1bSJed Brown appctx - user-defined application context
356c4762a1bSJed Brown
357c4762a1bSJed Brown Output Parameter:
358c4762a1bSJed Brown u - vector with solution at initial time (global)
359c4762a1bSJed Brown */
TrueSolution(Vec u,AppCtx * appctx)360d71ae5a4SJacob Faibussowitsch PetscErrorCode TrueSolution(Vec u, AppCtx *appctx)
361d71ae5a4SJacob Faibussowitsch {
362c4762a1bSJed Brown PetscScalar *s;
363c4762a1bSJed Brown const PetscScalar *xg;
364c4762a1bSJed Brown PetscInt i, j, lenglob;
365c4762a1bSJed Brown PetscReal sum;
366c4762a1bSJed Brown
367c4762a1bSJed Brown PetscFunctionBegin;
3689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, u, &s));
3699566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
370c4762a1bSJed Brown lenglob = appctx->param.E * (appctx->param.N - 1);
371c4762a1bSJed Brown for (i = 0; i < lenglob; i++) {
372c4762a1bSJed Brown s[i] = 0;
373ad540459SPierre Jolivet for (j = 0; j < appctx->ncoeff; j++) s[i] += appctx->solutioncoefficients[j] * PetscSinScalar(2 * (j + 1) * PETSC_PI * xg[i]);
374c4762a1bSJed Brown }
3759566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, u, &s));
3769566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
377c4762a1bSJed Brown /* make sure initial conditions do not contain the constant functions, since with periodic boundary conditions the constant functions introduce a null space */
3789566063dSJacob Faibussowitsch PetscCall(VecSum(u, &sum));
3799566063dSJacob Faibussowitsch PetscCall(VecShift(u, -sum / lenglob));
3803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
381c4762a1bSJed Brown }
382c4762a1bSJed Brown /* --------------------------------------------------------------------- */
383c4762a1bSJed Brown /*
384c4762a1bSJed Brown Sets the desired profile for the final end time
385c4762a1bSJed Brown
386c4762a1bSJed Brown Input Parameters:
387c4762a1bSJed Brown t - final time
388c4762a1bSJed Brown obj - vector storing the desired profile
389c4762a1bSJed Brown appctx - user-defined application context
390c4762a1bSJed Brown
391c4762a1bSJed Brown */
ComputeReference(TS ts,PetscReal t,Vec obj,AppCtx * appctx)392d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeReference(TS ts, PetscReal t, Vec obj, AppCtx *appctx)
393d71ae5a4SJacob Faibussowitsch {
394c4762a1bSJed Brown PetscScalar *s, tc;
395c4762a1bSJed Brown const PetscScalar *xg;
396c4762a1bSJed Brown PetscInt i, j, lenglob;
397c4762a1bSJed Brown
398c4762a1bSJed Brown PetscFunctionBegin;
3999566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, obj, &s));
4009566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
401c4762a1bSJed Brown lenglob = appctx->param.E * (appctx->param.N - 1);
402c4762a1bSJed Brown for (i = 0; i < lenglob; i++) {
403c4762a1bSJed Brown s[i] = 0;
404c4762a1bSJed Brown for (j = 0; j < appctx->ncoeff; j++) {
405c4762a1bSJed Brown tc = -appctx->param.mu * (j + 1) * (j + 1) * 4.0 * PETSC_PI * PETSC_PI * t;
406c4762a1bSJed Brown s[i] += appctx->solutioncoefficients[j] * PetscSinScalar(2 * (j + 1) * PETSC_PI * (xg[i] + appctx->param.a * t)) * PetscExpReal(tc);
407c4762a1bSJed Brown }
408c4762a1bSJed Brown }
4099566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, obj, &s));
4109566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
4113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
412c4762a1bSJed Brown }
413c4762a1bSJed Brown
RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,PetscCtx ctx)414*2a8381b2SBarry Smith PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, PetscCtx ctx)
415d71ae5a4SJacob Faibussowitsch {
416c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx;
417c4762a1bSJed Brown
418c4762a1bSJed Brown PetscFunctionBegin;
4199566063dSJacob Faibussowitsch PetscCall(MatMult(appctx->SEMop.keptstiff, globalin, globalout));
4203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
421c4762a1bSJed Brown }
422c4762a1bSJed Brown
RHSJacobian(TS ts,PetscReal t,Vec globalin,Mat A,Mat B,PetscCtx ctx)423*2a8381b2SBarry Smith PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec globalin, Mat A, Mat B, PetscCtx ctx)
424d71ae5a4SJacob Faibussowitsch {
425c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx;
426c4762a1bSJed Brown
427c4762a1bSJed Brown PetscFunctionBegin;
4289566063dSJacob Faibussowitsch PetscCall(MatCopy(appctx->SEMop.keptstiff, A, DIFFERENT_NONZERO_PATTERN));
4293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
430c4762a1bSJed Brown }
431c4762a1bSJed Brown
432c4762a1bSJed Brown /* --------------------------------------------------------------------- */
433c4762a1bSJed Brown
434c4762a1bSJed Brown /*
435c4762a1bSJed Brown RHSLaplacian - matrix for diffusion
436c4762a1bSJed Brown
437c4762a1bSJed Brown Input Parameters:
438c4762a1bSJed Brown ts - the TS context
439c4762a1bSJed Brown t - current time (ignored)
440c4762a1bSJed Brown X - current solution (ignored)
441c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian()
442c4762a1bSJed Brown
443c4762a1bSJed Brown Output Parameters:
444c4762a1bSJed Brown AA - Jacobian matrix
445c4762a1bSJed Brown BB - optionally different matrix from which the preconditioner is built
446c4762a1bSJed Brown
447c4762a1bSJed Brown Scales by the inverse of the mass matrix (perhaps that should be pulled out)
448c4762a1bSJed Brown
449c4762a1bSJed Brown */
RHSLaplacian(TS ts,PetscReal t,Vec X,Mat A,Mat BB,PetscCtx ctx)450*2a8381b2SBarry Smith PetscErrorCode RHSLaplacian(TS ts, PetscReal t, Vec X, Mat A, Mat BB, PetscCtx ctx)
451d71ae5a4SJacob Faibussowitsch {
452c4762a1bSJed Brown PetscReal **temp;
453c4762a1bSJed Brown PetscReal vv;
454c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
455c4762a1bSJed Brown PetscInt i, xs, xn, l, j;
456c4762a1bSJed Brown PetscInt *rowsDM;
457c4762a1bSJed Brown
458c4762a1bSJed Brown PetscFunctionBegin;
459c4762a1bSJed Brown /*
460c4762a1bSJed Brown Creates the element stiffness matrix for the given gll
461c4762a1bSJed Brown */
4629566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
463c4762a1bSJed Brown
464c4762a1bSJed Brown /* scale by the size of the element */
465c4762a1bSJed Brown for (i = 0; i < appctx->param.N; i++) {
466c4762a1bSJed Brown vv = -appctx->param.mu * 2.0 / appctx->param.Le;
467c4762a1bSJed Brown for (j = 0; j < appctx->param.N; j++) temp[i][j] = temp[i][j] * vv;
468c4762a1bSJed Brown }
469c4762a1bSJed Brown
4709566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
4719566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
472c4762a1bSJed Brown
4733c859ba3SBarry Smith PetscCheck(appctx->param.N - 1 >= 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Polynomial order must be at least 2");
474c4762a1bSJed Brown xs = xs / (appctx->param.N - 1);
475c4762a1bSJed Brown xn = xn / (appctx->param.N - 1);
476c4762a1bSJed Brown
4779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->param.N, &rowsDM));
478c4762a1bSJed Brown /*
479c4762a1bSJed Brown loop over local elements
480c4762a1bSJed Brown */
481c4762a1bSJed Brown for (j = xs; j < xs + xn; j++) {
482ad540459SPierre Jolivet for (l = 0; l < appctx->param.N; l++) rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l;
4839566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES));
484c4762a1bSJed Brown }
4859566063dSJacob Faibussowitsch PetscCall(PetscFree(rowsDM));
4869566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4879566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
4889566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass));
4899566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0));
4909566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass));
491c4762a1bSJed Brown
4929566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
4933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
494c4762a1bSJed Brown }
495c4762a1bSJed Brown
496c4762a1bSJed Brown /*
497c4762a1bSJed Brown Almost identical to Laplacian
498c4762a1bSJed Brown
499c4762a1bSJed Brown Note that the element matrix is NOT scaled by the size of element like the Laplacian term.
500c4762a1bSJed Brown */
RHSAdvection(TS ts,PetscReal t,Vec X,Mat A,Mat BB,PetscCtx ctx)501*2a8381b2SBarry Smith PetscErrorCode RHSAdvection(TS ts, PetscReal t, Vec X, Mat A, Mat BB, PetscCtx ctx)
502d71ae5a4SJacob Faibussowitsch {
503c4762a1bSJed Brown PetscReal **temp;
504c4762a1bSJed Brown PetscReal vv;
505c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
506c4762a1bSJed Brown PetscInt i, xs, xn, l, j;
507c4762a1bSJed Brown PetscInt *rowsDM;
508c4762a1bSJed Brown
509c4762a1bSJed Brown PetscFunctionBegin;
510c4762a1bSJed Brown /*
511c4762a1bSJed Brown Creates the element stiffness matrix for the given gll
512c4762a1bSJed Brown */
5139566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
514c4762a1bSJed Brown
515c4762a1bSJed Brown /* scale by the size of the element */
516c4762a1bSJed Brown for (i = 0; i < appctx->param.N; i++) {
517c4762a1bSJed Brown vv = -appctx->param.a;
518c4762a1bSJed Brown for (j = 0; j < appctx->param.N; j++) temp[i][j] = temp[i][j] * vv;
519c4762a1bSJed Brown }
520c4762a1bSJed Brown
5219566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
5229566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
523c4762a1bSJed Brown
5243c859ba3SBarry Smith PetscCheck(appctx->param.N - 1 >= 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Polynomial order must be at least 2");
525c4762a1bSJed Brown xs = xs / (appctx->param.N - 1);
526c4762a1bSJed Brown xn = xn / (appctx->param.N - 1);
527c4762a1bSJed Brown
5289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->param.N, &rowsDM));
529c4762a1bSJed Brown /*
530c4762a1bSJed Brown loop over local elements
531c4762a1bSJed Brown */
532c4762a1bSJed Brown for (j = xs; j < xs + xn; j++) {
533ad540459SPierre Jolivet for (l = 0; l < appctx->param.N; l++) rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l;
5349566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES));
535c4762a1bSJed Brown }
5369566063dSJacob Faibussowitsch PetscCall(PetscFree(rowsDM));
5379566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5389566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5399566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass));
5409566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0));
5419566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass));
542c4762a1bSJed Brown
5439566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
5443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
545c4762a1bSJed Brown }
546c4762a1bSJed Brown
547c4762a1bSJed Brown /* ------------------------------------------------------------------ */
548c4762a1bSJed Brown /*
549c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient.
550c4762a1bSJed Brown
551c4762a1bSJed Brown Input Parameters:
552c4762a1bSJed Brown tao - the Tao context
553c4762a1bSJed Brown ic - the input vector
554a82e8c82SStefano Zampini ctx - optional user-defined context, as set when calling TaoSetObjectiveAndGradient()
555c4762a1bSJed Brown
556c4762a1bSJed Brown Output Parameters:
557c4762a1bSJed Brown f - the newly evaluated function
558c4762a1bSJed Brown G - the newly evaluated gradient
559c4762a1bSJed Brown
560c4762a1bSJed Brown Notes:
561c4762a1bSJed Brown
562c4762a1bSJed Brown The forward equation is
563c4762a1bSJed Brown M u_t = F(U)
564c4762a1bSJed Brown which is converted to
565c4762a1bSJed Brown u_t = M^{-1} F(u)
566c4762a1bSJed Brown in the user code since TS has no direct way of providing a mass matrix. The Jacobian of this is
567c4762a1bSJed Brown M^{-1} J
568c4762a1bSJed Brown where J is the Jacobian of F. Now the adjoint equation is
569c4762a1bSJed Brown M v_t = J^T v
570c4762a1bSJed Brown but TSAdjoint does not solve this since it can only solve the transposed system for the
571c4762a1bSJed Brown Jacobian the user provided. Hence TSAdjoint solves
572c4762a1bSJed Brown w_t = J^T M^{-1} w (where w = M v)
573a5b23f4aSJose E. Roman since there is no way to indicate the mass matrix as a separate entity to TS. Thus one
574c4762a1bSJed Brown must be careful in initializing the "adjoint equation" and using the result. This is
575c4762a1bSJed Brown why
576c4762a1bSJed Brown G = -2 M(u(T) - u_d)
577c4762a1bSJed Brown below (instead of -2(u(T) - u_d)
578c4762a1bSJed Brown
579c4762a1bSJed Brown */
FormFunctionGradient(Tao tao,Vec ic,PetscReal * f,Vec G,PetscCtx ctx)580*2a8381b2SBarry Smith PetscErrorCode FormFunctionGradient(Tao tao, Vec ic, PetscReal *f, Vec G, PetscCtx ctx)
581d71ae5a4SJacob Faibussowitsch {
582c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
583c4762a1bSJed Brown Vec temp;
584c4762a1bSJed Brown
585c4762a1bSJed Brown PetscFunctionBegin;
5869566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx->ts, 0.0));
5879566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(appctx->ts, 0));
5889566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(appctx->ts, appctx->initial_dt));
5899566063dSJacob Faibussowitsch PetscCall(VecCopy(ic, appctx->dat.curr_sol));
590c4762a1bSJed Brown
5919566063dSJacob Faibussowitsch PetscCall(TSSolve(appctx->ts, appctx->dat.curr_sol));
5929566063dSJacob Faibussowitsch PetscCall(VecCopy(appctx->dat.curr_sol, appctx->dat.joe));
593c4762a1bSJed Brown
594c4762a1bSJed Brown /* Compute the difference between the current ODE solution and target ODE solution */
5959566063dSJacob Faibussowitsch PetscCall(VecWAXPY(G, -1.0, appctx->dat.curr_sol, appctx->dat.reference));
596c4762a1bSJed Brown
597c4762a1bSJed Brown /* Compute the objective/cost function */
5989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(G, &temp));
5999566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(temp, G, G));
6009566063dSJacob Faibussowitsch PetscCall(VecDot(temp, appctx->SEMop.mass, f));
6019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&temp));
602c4762a1bSJed Brown
603c4762a1bSJed Brown /* Compute initial conditions for the adjoint integration. See Notes above */
6049566063dSJacob Faibussowitsch PetscCall(VecScale(G, -2.0));
6059566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(G, G, appctx->SEMop.mass));
6069566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(appctx->ts, 1, &G, NULL));
607c4762a1bSJed Brown
6089566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(appctx->ts));
6099566063dSJacob Faibussowitsch /* PetscCall(VecPointwiseDivide(G,G,appctx->SEMop.mass));*/
6103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
611c4762a1bSJed Brown }
612c4762a1bSJed Brown
MonitorError(Tao tao,PetscCtx ctx)613*2a8381b2SBarry Smith PetscErrorCode MonitorError(Tao tao, PetscCtx ctx)
614d71ae5a4SJacob Faibussowitsch {
615c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx;
616c4762a1bSJed Brown Vec temp, grad;
617c4762a1bSJed Brown PetscReal nrm;
618c4762a1bSJed Brown PetscInt its;
619c4762a1bSJed Brown PetscReal fct, gnorm;
620c4762a1bSJed Brown
621c4762a1bSJed Brown PetscFunctionBegin;
6229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx->dat.ic, &temp));
6239566063dSJacob Faibussowitsch PetscCall(VecWAXPY(temp, -1.0, appctx->dat.ic, appctx->dat.true_solution));
6249566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(temp, temp, temp));
6259566063dSJacob Faibussowitsch PetscCall(VecDot(temp, appctx->SEMop.mass, &nrm));
626c4762a1bSJed Brown nrm = PetscSqrtReal(nrm);
6279566063dSJacob Faibussowitsch PetscCall(TaoGetGradient(tao, &grad, NULL, NULL));
6289566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(temp, temp, temp));
6299566063dSJacob Faibussowitsch PetscCall(VecDot(temp, appctx->SEMop.mass, &gnorm));
630c4762a1bSJed Brown gnorm = PetscSqrtReal(gnorm);
6319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&temp));
6329566063dSJacob Faibussowitsch PetscCall(TaoGetIterationNumber(tao, &its));
6339566063dSJacob Faibussowitsch PetscCall(TaoGetSolutionStatus(tao, NULL, &fct, NULL, NULL, NULL, NULL));
634c4762a1bSJed Brown if (!its) {
6359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%% Iteration Error Objective Gradient-norm\n"));
6369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "history = [\n"));
637c4762a1bSJed Brown }
63863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%3" PetscInt_FMT " %g %g %g\n", its, (double)nrm, (double)fct, (double)gnorm));
6393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
640c4762a1bSJed Brown }
641c4762a1bSJed Brown
MonitorDestroy(PetscCtxRt ctx)642*2a8381b2SBarry Smith PetscErrorCode MonitorDestroy(PetscCtxRt ctx)
643d71ae5a4SJacob Faibussowitsch {
644c4762a1bSJed Brown PetscFunctionBegin;
6459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "];\n"));
6463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
647c4762a1bSJed Brown }
648c4762a1bSJed Brown
649c4762a1bSJed Brown /*TEST
650c4762a1bSJed Brown
651c4762a1bSJed Brown build:
652c4762a1bSJed Brown requires: !complex
653c4762a1bSJed Brown
654c4762a1bSJed Brown test:
655c4762a1bSJed Brown requires: !single
656c4762a1bSJed Brown args: -ts_adapt_dt_max 3.e-3 -E 10 -N 8 -ncoeff 5 -tao_bqnls_mat_lmvm_scale_type none
657c4762a1bSJed Brown
658c4762a1bSJed Brown test:
659c4762a1bSJed Brown suffix: cn
660c4762a1bSJed Brown requires: !single
661188af4bfSBarry Smith args: -ts_type cn -ts_time_step .003 -pc_type lu -E 10 -N 8 -ncoeff 5 -tao_bqnls_mat_lmvm_scale_type none
662c4762a1bSJed Brown
663c4762a1bSJed Brown test:
664c4762a1bSJed Brown suffix: 2
665c4762a1bSJed Brown requires: !single
666c4762a1bSJed Brown args: -ts_adapt_dt_max 3.e-3 -E 10 -N 8 -ncoeff 5 -a .1 -tao_bqnls_mat_lmvm_scale_type none
667c4762a1bSJed Brown
668c4762a1bSJed Brown TEST*/
669