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