1 2 static char help[] = "Solves a simple data assimilation problem with one dimensional Burger's equation using TSAdjoint\n\n"; 3 4 /* 5 6 Not yet tested in parallel 7 8 */ 9 10 /* ------------------------------------------------------------------------ 11 12 This program uses the one-dimensional Burger's equation 13 u_t = mu*u_xx - u u_x, 14 on the domain 0 <= x <= 1, with periodic boundary conditions 15 16 to demonstrate solving a data assimilation problem of finding the initial conditions 17 to produce a given solution at a fixed time. 18 19 The operators are discretized with the spectral element method 20 21 See the paper PDE-CONSTRAINED OPTIMIZATION WITH SPECTRAL ELEMENTS USING PETSC AND TAO 22 by OANA MARIN, EMIL CONSTANTINESCU, AND BARRY SMITH for details on the exact solution 23 used 24 25 ------------------------------------------------------------------------- */ 26 27 #include <petsctao.h> 28 #include <petscts.h> 29 #include <petscdt.h> 30 #include <petscdraw.h> 31 #include <petscdmda.h> 32 33 /* 34 User-defined application context - contains data needed by the 35 application-provided call-back routines. 36 */ 37 38 typedef struct { 39 PetscInt n; /* number of nodes */ 40 PetscReal *nodes; /* GLL nodes */ 41 PetscReal *weights; /* GLL weights */ 42 } PetscGLL; 43 44 typedef struct { 45 PetscInt N; /* grid points per elements*/ 46 PetscInt E; /* number of elements */ 47 PetscReal tol_L2, tol_max; /* error norms */ 48 PetscInt steps; /* number of timesteps */ 49 PetscReal Tend; /* endtime */ 50 PetscReal mu; /* viscosity */ 51 PetscReal L; /* total length of domain */ 52 PetscReal Le; 53 PetscReal Tadj; 54 } PetscParam; 55 56 typedef struct { 57 Vec obj; /* desired end state */ 58 Vec grid; /* total grid */ 59 Vec grad; 60 Vec ic; 61 Vec curr_sol; 62 Vec true_solution; /* actual initial conditions for the final solution */ 63 } PetscData; 64 65 typedef struct { 66 Vec grid; /* total grid */ 67 Vec mass; /* mass matrix for total integration */ 68 Mat stiff; /* stifness matrix */ 69 Mat keptstiff; 70 Mat grad; 71 PetscGLL gll; 72 } PetscSEMOperators; 73 74 typedef struct { 75 DM da; /* distributed array data structure */ 76 PetscSEMOperators SEMop; 77 PetscParam param; 78 PetscData dat; 79 TS ts; 80 PetscReal initial_dt; 81 } AppCtx; 82 83 /* 84 User-defined routines 85 */ 86 extern PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 87 extern PetscErrorCode RHSMatrixLaplaciangllDM(TS, PetscReal, Vec, Mat, Mat, void *); 88 extern PetscErrorCode RHSMatrixAdvectiongllDM(TS, PetscReal, Vec, Mat, Mat, void *); 89 extern PetscErrorCode InitialConditions(Vec, AppCtx *); 90 extern PetscErrorCode TrueSolution(Vec, AppCtx *); 91 extern PetscErrorCode ComputeObjective(PetscReal, Vec, AppCtx *); 92 extern PetscErrorCode MonitorError(Tao, void *); 93 extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *); 94 extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *); 95 96 int main(int argc, char **argv) { 97 AppCtx appctx; /* user-defined application context */ 98 Tao tao; 99 Vec u; /* approximate solution vector */ 100 PetscInt i, xs, xm, ind, j, lenglob; 101 PetscReal x, *wrk_ptr1, *wrk_ptr2; 102 MatNullSpace nsp; 103 PetscMPIInt size; 104 105 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 106 Initialize program and set problem parameters 107 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 108 PetscFunctionBegin; 109 110 PetscFunctionBeginUser; 111 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 112 113 /*initialize parameters */ 114 appctx.param.N = 10; /* order of the spectral element */ 115 appctx.param.E = 10; /* number of elements */ 116 appctx.param.L = 4.0; /* length of the domain */ 117 appctx.param.mu = 0.01; /* diffusion coefficient */ 118 appctx.initial_dt = 5e-3; 119 appctx.param.steps = PETSC_MAX_INT; 120 appctx.param.Tend = 4; 121 122 PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &appctx.param.N, NULL)); 123 PetscCall(PetscOptionsGetInt(NULL, NULL, "-E", &appctx.param.E, NULL)); 124 PetscCall(PetscOptionsGetReal(NULL, NULL, "-Tend", &appctx.param.Tend, NULL)); 125 PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &appctx.param.mu, NULL)); 126 appctx.param.Le = appctx.param.L / appctx.param.E; 127 128 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 129 PetscCheck((appctx.param.E % size) == 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Number of elements must be divisible by number of processes"); 130 131 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 132 Create GLL data structures 133 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 134 PetscCall(PetscMalloc2(appctx.param.N, &appctx.SEMop.gll.nodes, appctx.param.N, &appctx.SEMop.gll.weights)); 135 PetscCall(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N, PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights)); 136 appctx.SEMop.gll.n = appctx.param.N; 137 lenglob = appctx.param.E * (appctx.param.N - 1); 138 139 /* 140 Create distributed array (DMDA) to manage parallel grid and vectors 141 and to set up the ghost point communication pattern. There are E*(Nl-1)+1 142 total grid values spread equally among all the processors, except first and last 143 */ 144 145 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, lenglob, 1, 1, NULL, &appctx.da)); 146 PetscCall(DMSetFromOptions(appctx.da)); 147 PetscCall(DMSetUp(appctx.da)); 148 149 /* 150 Extract global and local vectors from DMDA; we use these to store the 151 approximate solution. Then duplicate these for remaining vectors that 152 have the same types. 153 */ 154 155 PetscCall(DMCreateGlobalVector(appctx.da, &u)); 156 PetscCall(VecDuplicate(u, &appctx.dat.ic)); 157 PetscCall(VecDuplicate(u, &appctx.dat.true_solution)); 158 PetscCall(VecDuplicate(u, &appctx.dat.obj)); 159 PetscCall(VecDuplicate(u, &appctx.SEMop.grid)); 160 PetscCall(VecDuplicate(u, &appctx.SEMop.mass)); 161 PetscCall(VecDuplicate(u, &appctx.dat.curr_sol)); 162 163 PetscCall(DMDAGetCorners(appctx.da, &xs, NULL, NULL, &xm, NULL, NULL)); 164 PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1)); 165 PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2)); 166 167 /* Compute function over the locally owned part of the grid */ 168 169 xs = xs / (appctx.param.N - 1); 170 xm = xm / (appctx.param.N - 1); 171 172 /* 173 Build total grid and mass over entire mesh (multi-elemental) 174 */ 175 176 for (i = xs; i < xs + xm; i++) { 177 for (j = 0; j < appctx.param.N - 1; j++) { 178 x = (appctx.param.Le / 2.0) * (appctx.SEMop.gll.nodes[j] + 1.0) + appctx.param.Le * i; 179 ind = i * (appctx.param.N - 1) + j; 180 wrk_ptr1[ind] = x; 181 wrk_ptr2[ind] = .5 * appctx.param.Le * appctx.SEMop.gll.weights[j]; 182 if (j == 0) wrk_ptr2[ind] += .5 * appctx.param.Le * appctx.SEMop.gll.weights[j]; 183 } 184 } 185 PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1)); 186 PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2)); 187 188 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 189 Create matrix data structure; set matrix evaluation routine. 190 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 191 PetscCall(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE)); 192 PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.stiff)); 193 PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.grad)); 194 /* 195 For linear problems with a time-dependent f(u,t) in the equation 196 u_t = f(u,t), the user provides the discretized right-hand-side 197 as a time-dependent matrix. 198 */ 199 PetscCall(RHSMatrixLaplaciangllDM(appctx.ts, 0.0, u, appctx.SEMop.stiff, appctx.SEMop.stiff, &appctx)); 200 PetscCall(RHSMatrixAdvectiongllDM(appctx.ts, 0.0, u, appctx.SEMop.grad, appctx.SEMop.grad, &appctx)); 201 /* 202 For linear problems with a time-dependent f(u,t) in the equation 203 u_t = f(u,t), the user provides the discretized right-hand-side 204 as a time-dependent matrix. 205 */ 206 207 PetscCall(MatDuplicate(appctx.SEMop.stiff, MAT_COPY_VALUES, &appctx.SEMop.keptstiff)); 208 209 /* attach the null space to the matrix, this probably is not needed but does no harm */ 210 PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp)); 211 PetscCall(MatSetNullSpace(appctx.SEMop.stiff, nsp)); 212 PetscCall(MatSetNullSpace(appctx.SEMop.keptstiff, nsp)); 213 PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.stiff, NULL)); 214 PetscCall(MatNullSpaceDestroy(&nsp)); 215 /* attach the null space to the matrix, this probably is not needed but does no harm */ 216 PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp)); 217 PetscCall(MatSetNullSpace(appctx.SEMop.grad, nsp)); 218 PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.grad, NULL)); 219 PetscCall(MatNullSpaceDestroy(&nsp)); 220 221 /* Create the TS solver that solves the ODE and its adjoint; set its options */ 222 PetscCall(TSCreate(PETSC_COMM_WORLD, &appctx.ts)); 223 PetscCall(TSSetProblemType(appctx.ts, TS_NONLINEAR)); 224 PetscCall(TSSetType(appctx.ts, TSRK)); 225 PetscCall(TSSetDM(appctx.ts, appctx.da)); 226 PetscCall(TSSetTime(appctx.ts, 0.0)); 227 PetscCall(TSSetTimeStep(appctx.ts, appctx.initial_dt)); 228 PetscCall(TSSetMaxSteps(appctx.ts, appctx.param.steps)); 229 PetscCall(TSSetMaxTime(appctx.ts, appctx.param.Tend)); 230 PetscCall(TSSetExactFinalTime(appctx.ts, TS_EXACTFINALTIME_MATCHSTEP)); 231 PetscCall(TSSetTolerances(appctx.ts, 1e-7, NULL, 1e-7, NULL)); 232 PetscCall(TSSetFromOptions(appctx.ts)); 233 /* Need to save initial timestep user may have set with -ts_dt so it can be reset for each new TSSolve() */ 234 PetscCall(TSGetTimeStep(appctx.ts, &appctx.initial_dt)); 235 PetscCall(TSSetRHSFunction(appctx.ts, NULL, RHSFunction, &appctx)); 236 PetscCall(TSSetRHSJacobian(appctx.ts, appctx.SEMop.stiff, appctx.SEMop.stiff, RHSJacobian, &appctx)); 237 238 /* Set Objective and Initial conditions for the problem and compute Objective function (evolution of true_solution to final time */ 239 PetscCall(InitialConditions(appctx.dat.ic, &appctx)); 240 PetscCall(TrueSolution(appctx.dat.true_solution, &appctx)); 241 PetscCall(ComputeObjective(appctx.param.Tend, appctx.dat.obj, &appctx)); 242 243 PetscCall(TSSetSaveTrajectory(appctx.ts)); 244 PetscCall(TSSetFromOptions(appctx.ts)); 245 246 /* Create TAO solver and set desired solution method */ 247 PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 248 PetscCall(TaoSetMonitor(tao, MonitorError, &appctx, NULL)); 249 PetscCall(TaoSetType(tao, TAOBQNLS)); 250 PetscCall(TaoSetSolution(tao, appctx.dat.ic)); 251 /* Set routine for function and gradient evaluation */ 252 PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&appctx)); 253 /* Check for any TAO command line options */ 254 PetscCall(TaoSetTolerances(tao, 1e-8, PETSC_DEFAULT, PETSC_DEFAULT)); 255 PetscCall(TaoSetFromOptions(tao)); 256 PetscCall(TaoSolve(tao)); 257 258 PetscCall(TaoDestroy(&tao)); 259 PetscCall(MatDestroy(&appctx.SEMop.stiff)); 260 PetscCall(MatDestroy(&appctx.SEMop.keptstiff)); 261 PetscCall(MatDestroy(&appctx.SEMop.grad)); 262 PetscCall(VecDestroy(&u)); 263 PetscCall(VecDestroy(&appctx.dat.ic)); 264 PetscCall(VecDestroy(&appctx.dat.true_solution)); 265 PetscCall(VecDestroy(&appctx.dat.obj)); 266 PetscCall(VecDestroy(&appctx.SEMop.grid)); 267 PetscCall(VecDestroy(&appctx.SEMop.mass)); 268 PetscCall(VecDestroy(&appctx.dat.curr_sol)); 269 PetscCall(PetscFree2(appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights)); 270 PetscCall(DMDestroy(&appctx.da)); 271 PetscCall(TSDestroy(&appctx.ts)); 272 273 /* 274 Always call PetscFinalize() before exiting a program. This routine 275 - finalizes the PETSc libraries as well as MPI 276 - provides summary and diagnostic information if certain runtime 277 options are chosen (e.g., -log_summary). 278 */ 279 PetscCall(PetscFinalize()); 280 return 0; 281 } 282 283 /* --------------------------------------------------------------------- */ 284 /* 285 InitialConditions - Computes the initial conditions for the Tao optimization solve (these are also initial conditions for the first TSSolve() 286 287 The routine TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function 288 289 Input Parameter: 290 u - uninitialized solution vector (global) 291 appctx - user-defined application context 292 293 Output Parameter: 294 u - vector with solution at initial time (global) 295 */ 296 PetscErrorCode InitialConditions(Vec u, AppCtx *appctx) { 297 PetscScalar *s; 298 const PetscScalar *xg; 299 PetscInt i, xs, xn; 300 301 PetscFunctionBegin; 302 PetscCall(DMDAVecGetArray(appctx->da, u, &s)); 303 PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg)); 304 PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 305 for (i = xs; i < xs + xn; i++) { s[i] = 2.0 * appctx->param.mu * PETSC_PI * PetscSinScalar(PETSC_PI * xg[i]) / (2.0 + PetscCosScalar(PETSC_PI * xg[i])) + 0.25 * PetscExpReal(-4.0 * PetscPowReal(xg[i] - 2.0, 2.0)); } 306 PetscCall(DMDAVecRestoreArray(appctx->da, u, &s)); 307 PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg)); 308 PetscFunctionReturn(0); 309 } 310 311 /* 312 TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function. 313 314 InitialConditions() computes the initial conditions for the beginning of the Tao iterations 315 316 Input Parameter: 317 u - uninitialized solution vector (global) 318 appctx - user-defined application context 319 320 Output Parameter: 321 u - vector with solution at initial time (global) 322 */ 323 PetscErrorCode TrueSolution(Vec u, AppCtx *appctx) { 324 PetscScalar *s; 325 const PetscScalar *xg; 326 PetscInt i, xs, xn; 327 328 PetscFunctionBegin; 329 PetscCall(DMDAVecGetArray(appctx->da, u, &s)); 330 PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg)); 331 PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 332 for (i = xs; i < xs + xn; i++) { s[i] = 2.0 * appctx->param.mu * PETSC_PI * PetscSinScalar(PETSC_PI * xg[i]) / (2.0 + PetscCosScalar(PETSC_PI * xg[i])); } 333 PetscCall(DMDAVecRestoreArray(appctx->da, u, &s)); 334 PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg)); 335 PetscFunctionReturn(0); 336 } 337 /* --------------------------------------------------------------------- */ 338 /* 339 Sets the desired profile for the final end time 340 341 Input Parameters: 342 t - final time 343 obj - vector storing the desired profile 344 appctx - user-defined application context 345 346 */ 347 PetscErrorCode ComputeObjective(PetscReal t, Vec obj, AppCtx *appctx) { 348 PetscScalar *s; 349 const PetscScalar *xg; 350 PetscInt i, xs, xn; 351 352 PetscFunctionBegin; 353 PetscCall(DMDAVecGetArray(appctx->da, obj, &s)); 354 PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg)); 355 PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 356 for (i = xs; i < xs + xn; i++) { 357 s[i] = 2.0 * appctx->param.mu * PETSC_PI * PetscSinScalar(PETSC_PI * xg[i]) * PetscExpScalar(-PETSC_PI * PETSC_PI * t * appctx->param.mu) / (2.0 + PetscExpScalar(-PETSC_PI * PETSC_PI * t * appctx->param.mu) * PetscCosScalar(PETSC_PI * xg[i])); 358 } 359 PetscCall(DMDAVecRestoreArray(appctx->da, obj, &s)); 360 PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg)); 361 PetscFunctionReturn(0); 362 } 363 364 PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx) { 365 AppCtx *appctx = (AppCtx *)ctx; 366 367 PetscFunctionBegin; 368 PetscCall(MatMult(appctx->SEMop.grad, globalin, globalout)); /* grad u */ 369 PetscCall(VecPointwiseMult(globalout, globalin, globalout)); /* u grad u */ 370 PetscCall(VecScale(globalout, -1.0)); 371 PetscCall(MatMultAdd(appctx->SEMop.keptstiff, globalin, globalout, globalout)); 372 PetscFunctionReturn(0); 373 } 374 375 /* 376 377 K is the discretiziation of the Laplacian 378 G is the discretization of the gradient 379 380 Computes Jacobian of K u + diag(u) G u which is given by 381 K + diag(u)G + diag(Gu) 382 */ 383 PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec globalin, Mat A, Mat B, void *ctx) { 384 AppCtx *appctx = (AppCtx *)ctx; 385 Vec Gglobalin; 386 387 PetscFunctionBegin; 388 /* A = diag(u) G */ 389 390 PetscCall(MatCopy(appctx->SEMop.grad, A, SAME_NONZERO_PATTERN)); 391 PetscCall(MatDiagonalScale(A, globalin, NULL)); 392 393 /* A = A + diag(Gu) */ 394 PetscCall(VecDuplicate(globalin, &Gglobalin)); 395 PetscCall(MatMult(appctx->SEMop.grad, globalin, Gglobalin)); 396 PetscCall(MatDiagonalSet(A, Gglobalin, ADD_VALUES)); 397 PetscCall(VecDestroy(&Gglobalin)); 398 399 /* A = K - A */ 400 PetscCall(MatScale(A, -1.0)); 401 PetscCall(MatAXPY(A, 1.0, appctx->SEMop.keptstiff, SAME_NONZERO_PATTERN)); 402 PetscFunctionReturn(0); 403 } 404 405 /* --------------------------------------------------------------------- */ 406 407 /* 408 RHSMatrixLaplacian - User-provided routine to compute the right-hand-side 409 matrix for the heat equation. 410 411 Input Parameters: 412 ts - the TS context 413 t - current time (ignored) 414 X - current solution (ignored) 415 dummy - optional user-defined context, as set by TSetRHSJacobian() 416 417 Output Parameters: 418 AA - Jacobian matrix 419 BB - optionally different matrix from which the preconditioner is built 420 str - flag indicating matrix structure 421 422 */ 423 PetscErrorCode RHSMatrixLaplaciangllDM(TS ts, PetscReal t, Vec X, Mat A, Mat BB, void *ctx) { 424 PetscReal **temp; 425 PetscReal vv; 426 AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 427 PetscInt i, xs, xn, l, j; 428 PetscInt *rowsDM; 429 430 PetscFunctionBegin; 431 /* 432 Creates the element stiffness matrix for the given gll 433 */ 434 PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 435 /* workaround for clang analyzer warning: Division by zero */ 436 PetscCheck(appctx->param.N > 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Spectral element order should be > 1"); 437 438 /* scale by the size of the element */ 439 for (i = 0; i < appctx->param.N; i++) { 440 vv = -appctx->param.mu * 2.0 / appctx->param.Le; 441 for (j = 0; j < appctx->param.N; j++) temp[i][j] = temp[i][j] * vv; 442 } 443 444 PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 445 PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 446 447 xs = xs / (appctx->param.N - 1); 448 xn = xn / (appctx->param.N - 1); 449 450 PetscCall(PetscMalloc1(appctx->param.N, &rowsDM)); 451 /* 452 loop over local elements 453 */ 454 for (j = xs; j < xs + xn; j++) { 455 for (l = 0; l < appctx->param.N; l++) { rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l; } 456 PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES)); 457 } 458 PetscCall(PetscFree(rowsDM)); 459 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 460 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 461 PetscCall(VecReciprocal(appctx->SEMop.mass)); 462 PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0)); 463 PetscCall(VecReciprocal(appctx->SEMop.mass)); 464 465 PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 466 PetscFunctionReturn(0); 467 } 468 469 /* 470 RHSMatrixAdvection - User-provided routine to compute the right-hand-side 471 matrix for the Advection equation. 472 473 Input Parameters: 474 ts - the TS context 475 t - current time 476 global_in - global input vector 477 dummy - optional user-defined context, as set by TSetRHSJacobian() 478 479 Output Parameters: 480 AA - Jacobian matrix 481 BB - optionally different preconditioning matrix 482 str - flag indicating matrix structure 483 484 */ 485 PetscErrorCode RHSMatrixAdvectiongllDM(TS ts, PetscReal t, Vec X, Mat A, Mat BB, void *ctx) { 486 PetscReal **temp; 487 AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 488 PetscInt xs, xn, l, j; 489 PetscInt *rowsDM; 490 491 PetscFunctionBegin; 492 /* 493 Creates the advection matrix for the given gll 494 */ 495 PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 496 PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 497 498 PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 499 500 xs = xs / (appctx->param.N - 1); 501 xn = xn / (appctx->param.N - 1); 502 503 PetscCall(PetscMalloc1(appctx->param.N, &rowsDM)); 504 for (j = xs; j < xs + xn; j++) { 505 for (l = 0; l < appctx->param.N; l++) { rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l; } 506 PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES)); 507 } 508 PetscCall(PetscFree(rowsDM)); 509 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 510 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 511 512 PetscCall(VecReciprocal(appctx->SEMop.mass)); 513 PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0)); 514 PetscCall(VecReciprocal(appctx->SEMop.mass)); 515 PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 516 PetscFunctionReturn(0); 517 } 518 /* ------------------------------------------------------------------ */ 519 /* 520 FormFunctionGradient - Evaluates the function and corresponding gradient. 521 522 Input Parameters: 523 tao - the Tao context 524 IC - the input vector 525 ctx - optional user-defined context, as set when calling TaoSetObjectiveAndGradient() 526 527 Output Parameters: 528 f - the newly evaluated function 529 G - the newly evaluated gradient 530 531 Notes: 532 533 The forward equation is 534 M u_t = F(U) 535 which is converted to 536 u_t = M^{-1} F(u) 537 in the user code since TS has no direct way of providing a mass matrix. The Jacobian of this is 538 M^{-1} J 539 where J is the Jacobian of F. Now the adjoint equation is 540 M v_t = J^T v 541 but TSAdjoint does not solve this since it can only solve the transposed system for the 542 Jacobian the user provided. Hence TSAdjoint solves 543 w_t = J^T M^{-1} w (where w = M v) 544 since there is no way to indicate the mass matrix as a separate entity to TS. Thus one 545 must be careful in initializing the "adjoint equation" and using the result. This is 546 why 547 G = -2 M(u(T) - u_d) 548 below (instead of -2(u(T) - u_d) and why the result is 549 G = G/appctx->SEMop.mass (that is G = M^{-1}w) 550 below (instead of just the result of the "adjoint solve"). 551 552 */ 553 PetscErrorCode FormFunctionGradient(Tao tao, Vec IC, PetscReal *f, Vec G, void *ctx) { 554 AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 555 Vec temp; 556 PetscInt its; 557 PetscReal ff, gnorm, cnorm, xdiff, errex; 558 TaoConvergedReason reason; 559 560 PetscFunctionBegin; 561 PetscCall(TSSetTime(appctx->ts, 0.0)); 562 PetscCall(TSSetStepNumber(appctx->ts, 0)); 563 PetscCall(TSSetTimeStep(appctx->ts, appctx->initial_dt)); 564 PetscCall(VecCopy(IC, appctx->dat.curr_sol)); 565 566 PetscCall(TSSolve(appctx->ts, appctx->dat.curr_sol)); 567 568 PetscCall(VecWAXPY(G, -1.0, appctx->dat.curr_sol, appctx->dat.obj)); 569 570 /* 571 Compute the L2-norm of the objective function, cost function is f 572 */ 573 PetscCall(VecDuplicate(G, &temp)); 574 PetscCall(VecPointwiseMult(temp, G, G)); 575 PetscCall(VecDot(temp, appctx->SEMop.mass, f)); 576 577 /* local error evaluation */ 578 PetscCall(VecWAXPY(temp, -1.0, appctx->dat.ic, appctx->dat.true_solution)); 579 PetscCall(VecPointwiseMult(temp, temp, temp)); 580 /* for error evaluation */ 581 PetscCall(VecDot(temp, appctx->SEMop.mass, &errex)); 582 PetscCall(VecDestroy(&temp)); 583 errex = PetscSqrtReal(errex); 584 585 /* 586 Compute initial conditions for the adjoint integration. See Notes above 587 */ 588 589 PetscCall(VecScale(G, -2.0)); 590 PetscCall(VecPointwiseMult(G, G, appctx->SEMop.mass)); 591 PetscCall(TSSetCostGradients(appctx->ts, 1, &G, NULL)); 592 PetscCall(TSAdjointSolve(appctx->ts)); 593 PetscCall(VecPointwiseDivide(G, G, appctx->SEMop.mass)); 594 595 PetscCall(TaoGetSolutionStatus(tao, &its, &ff, &gnorm, &cnorm, &xdiff, &reason)); 596 PetscFunctionReturn(0); 597 } 598 599 PetscErrorCode MonitorError(Tao tao, void *ctx) { 600 AppCtx *appctx = (AppCtx *)ctx; 601 Vec temp; 602 PetscReal nrm; 603 604 PetscFunctionBegin; 605 PetscCall(VecDuplicate(appctx->dat.ic, &temp)); 606 PetscCall(VecWAXPY(temp, -1.0, appctx->dat.ic, appctx->dat.true_solution)); 607 PetscCall(VecPointwiseMult(temp, temp, temp)); 608 PetscCall(VecDot(temp, appctx->SEMop.mass, &nrm)); 609 PetscCall(VecDestroy(&temp)); 610 nrm = PetscSqrtReal(nrm); 611 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error for initial conditions %g\n", (double)nrm)); 612 PetscFunctionReturn(0); 613 } 614 615 /*TEST 616 617 build: 618 requires: !complex 619 620 test: 621 args: -tao_max_it 5 -tao_gatol 1.e-4 622 requires: !single 623 624 test: 625 suffix: 2 626 nsize: 2 627 args: -tao_max_it 5 -tao_gatol 1.e-4 628 requires: !single 629 630 TEST*/ 631