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 ierr = PetscOptionsGetInt(NULL,NULL,"-N",&appctx.param.N,NULL);CHKERRQ(ierr); 143 ierr = PetscOptionsGetInt(NULL,NULL,"-E",&appctx.param.E,NULL);CHKERRQ(ierr); 144 ierr = PetscOptionsGetInt(NULL,NULL,"-ncoeff",&appctx.ncoeff,NULL);CHKERRQ(ierr); 145 ierr = PetscOptionsGetReal(NULL,NULL,"-Tend",&appctx.param.Tend,NULL);CHKERRQ(ierr); 146 ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&appctx.param.mu,NULL);CHKERRQ(ierr); 147 ierr = PetscOptionsGetReal(NULL,NULL,"-a",&appctx.param.a,NULL);CHKERRQ(ierr); 148 appctx.param.Le = appctx.param.L/appctx.param.E; 149 150 151 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 152 Create GLL data structures 153 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 154 ierr = PetscMalloc2(appctx.param.N,&appctx.SEMop.gll.nodes,appctx.param.N,&appctx.SEMop.gll.weights);CHKERRQ(ierr); 155 ierr = PetscDTGaussLobattoLegendreQuadrature(appctx.param.N,PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights);CHKERRQ(ierr); 156 appctx.SEMop.gll.n = appctx.param.N; 157 lenglob = appctx.param.E*(appctx.param.N-1); 158 159 /* 160 Create distributed array (DMDA) to manage parallel grid and vectors 161 and to set up the ghost point communication pattern. There are E*(Nl-1)+1 162 total grid values spread equally among all the processors, except first and last 163 */ 164 165 ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,lenglob,1,1,NULL,&appctx.da);CHKERRQ(ierr); 166 ierr = DMSetFromOptions(appctx.da);CHKERRQ(ierr); 167 ierr = DMSetUp(appctx.da);CHKERRQ(ierr); 168 169 /* 170 Extract global and local vectors from DMDA; we use these to store the 171 approximate solution. Then duplicate these for remaining vectors that 172 have the same types. 173 */ 174 175 ierr = DMCreateGlobalVector(appctx.da,&u);CHKERRQ(ierr); 176 ierr = VecDuplicate(u,&appctx.dat.ic);CHKERRQ(ierr); 177 ierr = VecDuplicate(u,&appctx.dat.true_solution);CHKERRQ(ierr); 178 ierr = VecDuplicate(u,&appctx.dat.reference);CHKERRQ(ierr); 179 ierr = VecDuplicate(u,&appctx.SEMop.grid);CHKERRQ(ierr); 180 ierr = VecDuplicate(u,&appctx.SEMop.mass);CHKERRQ(ierr); 181 ierr = VecDuplicate(u,&appctx.dat.curr_sol);CHKERRQ(ierr); 182 ierr = VecDuplicate(u,&appctx.dat.joe);CHKERRQ(ierr); 183 184 ierr = DMDAGetCorners(appctx.da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 185 ierr = DMDAVecGetArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1);CHKERRQ(ierr); 186 ierr = DMDAVecGetArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2);CHKERRQ(ierr); 187 188 /* Compute function over the locally owned part of the grid */ 189 190 xs=xs/(appctx.param.N-1); 191 xm=xm/(appctx.param.N-1); 192 193 /* 194 Build total grid and mass over entire mesh (multi-elemental) 195 */ 196 197 for (i=xs; i<xs+xm; i++) { 198 for (j=0; j<appctx.param.N-1; j++) { 199 x = (appctx.param.Le/2.0)*(appctx.SEMop.gll.nodes[j]+1.0)+appctx.param.Le*i; 200 ind=i*(appctx.param.N-1)+j; 201 wrk_ptr1[ind]=x; 202 wrk_ptr2[ind]=.5*appctx.param.Le*appctx.SEMop.gll.weights[j]; 203 if (j==0) wrk_ptr2[ind]+=.5*appctx.param.Le*appctx.SEMop.gll.weights[j]; 204 } 205 } 206 ierr = DMDAVecRestoreArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1);CHKERRQ(ierr); 207 ierr = DMDAVecRestoreArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2);CHKERRQ(ierr); 208 209 210 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 211 Create matrix data structure; set matrix evaluation routine. 212 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 213 ierr = DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE);CHKERRQ(ierr); 214 ierr = DMCreateMatrix(appctx.da,&appctx.SEMop.stiff);CHKERRQ(ierr); 215 ierr = DMCreateMatrix(appctx.da,&appctx.SEMop.advec);CHKERRQ(ierr); 216 217 /* 218 For linear problems with a time-dependent f(u,t) in the equation 219 u_t = f(u,t), the user provides the discretized right-hand-side 220 as a time-dependent matrix. 221 */ 222 ierr = RHSLaplacian(appctx.ts,0.0,u,appctx.SEMop.stiff,appctx.SEMop.stiff,&appctx);CHKERRQ(ierr); 223 ierr = RHSAdvection(appctx.ts,0.0,u,appctx.SEMop.advec,appctx.SEMop.advec,&appctx);CHKERRQ(ierr); 224 ierr = MatAXPY(appctx.SEMop.stiff,-1.0,appctx.SEMop.advec,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 225 ierr = MatDuplicate(appctx.SEMop.stiff,MAT_COPY_VALUES,&appctx.SEMop.keptstiff);CHKERRQ(ierr); 226 227 /* attach the null space to the matrix, this probably is not needed but does no harm */ 228 ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nsp);CHKERRQ(ierr); 229 ierr = MatSetNullSpace(appctx.SEMop.stiff,nsp);CHKERRQ(ierr); 230 ierr = MatNullSpaceTest(nsp,appctx.SEMop.stiff,NULL);CHKERRQ(ierr); 231 ierr = MatNullSpaceDestroy(&nsp);CHKERRQ(ierr); 232 233 /* Create the TS solver that solves the ODE and its adjoint; set its options */ 234 ierr = TSCreate(PETSC_COMM_WORLD,&appctx.ts);CHKERRQ(ierr); 235 ierr = TSSetSolutionFunction(appctx.ts,(PetscErrorCode (*)(TS,PetscReal,Vec, void *))ComputeReference,&appctx);CHKERRQ(ierr); 236 ierr = TSSetProblemType(appctx.ts,TS_LINEAR);CHKERRQ(ierr); 237 ierr = TSSetType(appctx.ts,TSRK);CHKERRQ(ierr); 238 ierr = TSSetDM(appctx.ts,appctx.da);CHKERRQ(ierr); 239 ierr = TSSetTime(appctx.ts,0.0);CHKERRQ(ierr); 240 ierr = TSSetTimeStep(appctx.ts,appctx.initial_dt);CHKERRQ(ierr); 241 ierr = TSSetMaxSteps(appctx.ts,appctx.param.steps);CHKERRQ(ierr); 242 ierr = TSSetMaxTime(appctx.ts,appctx.param.Tend);CHKERRQ(ierr); 243 ierr = TSSetExactFinalTime(appctx.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 244 ierr = TSSetTolerances(appctx.ts,1e-7,NULL,1e-7,NULL);CHKERRQ(ierr); 245 ierr = TSSetFromOptions(appctx.ts);CHKERRQ(ierr); 246 /* Need to save initial timestep user may have set with -ts_dt so it can be reset for each new TSSolve() */ 247 ierr = TSGetTimeStep(appctx.ts,&appctx.initial_dt);CHKERRQ(ierr); 248 ierr = TSSetRHSFunction(appctx.ts,NULL,TSComputeRHSFunctionLinear,&appctx);CHKERRQ(ierr); 249 ierr = TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,TSComputeRHSJacobianConstant,&appctx);CHKERRQ(ierr); 250 /* ierr = TSSetRHSFunction(appctx.ts,NULL,RHSFunction,&appctx);CHKERRQ(ierr); 251 ierr = TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,RHSJacobian,&appctx);CHKERRQ(ierr); */ 252 ierr = TSSetSaveTrajectory(appctx.ts);CHKERRQ(ierr); 253 254 /* Set random initial conditions as initial guess, compute analytic reference solution and analytic (true) initial conditions */ 255 ierr = ComputeSolutionCoefficients(&appctx);CHKERRQ(ierr); 256 ierr = InitialConditions(appctx.dat.ic,&appctx);CHKERRQ(ierr); 257 ierr = ComputeReference(appctx.ts,appctx.param.Tend,appctx.dat.reference,&appctx);CHKERRQ(ierr); 258 ierr = ComputeReference(appctx.ts,0.0,appctx.dat.true_solution,&appctx);CHKERRQ(ierr); 259 260 /* Create TAO solver and set desired solution method */ 261 ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr); 262 ierr = TaoSetMonitor(tao,MonitorError,&appctx,MonitorDestroy);CHKERRQ(ierr); 263 ierr = TaoSetType(tao,TAOBQNLS);CHKERRQ(ierr); 264 ierr = TaoSetInitialVector(tao,appctx.dat.ic);CHKERRQ(ierr); 265 /* Set routine for function and gradient evaluation */ 266 ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&appctx);CHKERRQ(ierr); 267 /* Check for any TAO command line options */ 268 ierr = TaoSetTolerances(tao,1e-8,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 269 ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 270 ierr = TaoSolve(tao);CHKERRQ(ierr); 271 272 ierr = TaoDestroy(&tao);CHKERRQ(ierr); 273 ierr = PetscFree(appctx.solutioncoefficients);CHKERRQ(ierr); 274 ierr = MatDestroy(&appctx.SEMop.advec);CHKERRQ(ierr); 275 ierr = MatDestroy(&appctx.SEMop.stiff);CHKERRQ(ierr); 276 ierr = MatDestroy(&appctx.SEMop.keptstiff);CHKERRQ(ierr); 277 ierr = VecDestroy(&u);CHKERRQ(ierr); 278 ierr = VecDestroy(&appctx.dat.ic);CHKERRQ(ierr); 279 ierr = VecDestroy(&appctx.dat.joe);CHKERRQ(ierr); 280 ierr = VecDestroy(&appctx.dat.true_solution);CHKERRQ(ierr); 281 ierr = VecDestroy(&appctx.dat.reference);CHKERRQ(ierr); 282 ierr = VecDestroy(&appctx.SEMop.grid);CHKERRQ(ierr); 283 ierr = VecDestroy(&appctx.SEMop.mass);CHKERRQ(ierr); 284 ierr = VecDestroy(&appctx.dat.curr_sol);CHKERRQ(ierr); 285 ierr = PetscFree2(appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights);CHKERRQ(ierr); 286 ierr = DMDestroy(&appctx.da);CHKERRQ(ierr); 287 ierr = TSDestroy(&appctx.ts);CHKERRQ(ierr); 288 289 /* 290 Always call PetscFinalize() before exiting a program. This routine 291 - finalizes the PETSc libraries as well as MPI 292 - provides summary and diagnostic information if certain runtime 293 options are chosen (e.g., -log_summary). 294 */ 295 ierr = PetscFinalize(); 296 return ierr; 297 } 298 299 /* 300 Computes the coefficients for the analytic solution to the PDE 301 */ 302 PetscErrorCode ComputeSolutionCoefficients(AppCtx *appctx) 303 { 304 PetscErrorCode ierr; 305 PetscRandom rand; 306 PetscInt i; 307 308 PetscFunctionBegin; 309 ierr = PetscMalloc1(appctx->ncoeff,&appctx->solutioncoefficients);CHKERRQ(ierr); 310 ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); 311 ierr = PetscRandomSetInterval(rand,.9,1.0);CHKERRQ(ierr); 312 for (i=0; i<appctx->ncoeff; i++) { 313 ierr = PetscRandomGetValue(rand,&appctx->solutioncoefficients[i]);CHKERRQ(ierr); 314 } 315 ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 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 PetscErrorCode ierr; 335 PetscInt i,j,lenglob; 336 PetscReal sum,val; 337 PetscRandom rand; 338 339 PetscFunctionBegin; 340 ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); 341 ierr = PetscRandomSetInterval(rand,.9,1.0);CHKERRQ(ierr); 342 ierr = DMDAVecGetArray(appctx->da,u,&s);CHKERRQ(ierr); 343 ierr = DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg);CHKERRQ(ierr); 344 lenglob = appctx->param.E*(appctx->param.N-1); 345 for (i=0; i<lenglob; i++) { 346 s[i]= 0; 347 for (j=0; j<appctx->ncoeff; j++) { 348 ierr = PetscRandomGetValue(rand,&val);CHKERRQ(ierr); 349 s[i] += val*PetscSinScalar(2*(j+1)*PETSC_PI*xg[i]); 350 } 351 } 352 ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 353 ierr = DMDAVecRestoreArray(appctx->da,u,&s);CHKERRQ(ierr); 354 ierr = DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg);CHKERRQ(ierr); 355 /* make sure initial conditions do not contain the constant functions, since with periodic boundary conditions the constant functions introduce a null space */ 356 ierr = VecSum(u,&sum);CHKERRQ(ierr); 357 ierr = VecShift(u,-sum/lenglob);CHKERRQ(ierr); 358 PetscFunctionReturn(0); 359 } 360 361 362 /* 363 TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function. 364 365 InitialConditions() computes the initial conditions for the begining of the Tao iterations 366 367 Input Parameter: 368 u - uninitialized solution vector (global) 369 appctx - user-defined application context 370 371 Output Parameter: 372 u - vector with solution at initial time (global) 373 */ 374 PetscErrorCode TrueSolution(Vec u,AppCtx *appctx) 375 { 376 PetscScalar *s; 377 const PetscScalar *xg; 378 PetscErrorCode ierr; 379 PetscInt i,j,lenglob; 380 PetscReal sum; 381 382 PetscFunctionBegin; 383 ierr = DMDAVecGetArray(appctx->da,u,&s);CHKERRQ(ierr); 384 ierr = DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg);CHKERRQ(ierr); 385 lenglob = appctx->param.E*(appctx->param.N-1); 386 for (i=0; i<lenglob; i++) { 387 s[i]= 0; 388 for (j=0; j<appctx->ncoeff; j++) { 389 s[i] += appctx->solutioncoefficients[j]*PetscSinScalar(2*(j+1)*PETSC_PI*xg[i]); 390 } 391 } 392 ierr = DMDAVecRestoreArray(appctx->da,u,&s);CHKERRQ(ierr); 393 ierr = DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg);CHKERRQ(ierr); 394 /* make sure initial conditions do not contain the constant functions, since with periodic boundary conditions the constant functions introduce a null space */ 395 ierr = VecSum(u,&sum);CHKERRQ(ierr); 396 ierr = VecShift(u,-sum/lenglob);CHKERRQ(ierr); 397 PetscFunctionReturn(0); 398 } 399 /* --------------------------------------------------------------------- */ 400 /* 401 Sets the desired profile for the final end time 402 403 Input Parameters: 404 t - final time 405 obj - vector storing the desired profile 406 appctx - user-defined application context 407 408 */ 409 PetscErrorCode ComputeReference(TS ts,PetscReal t,Vec obj,AppCtx *appctx) 410 { 411 PetscScalar *s,tc; 412 const PetscScalar *xg; 413 PetscErrorCode ierr; 414 PetscInt i, j,lenglob; 415 416 PetscFunctionBegin; 417 ierr = DMDAVecGetArray(appctx->da,obj,&s);CHKERRQ(ierr); 418 ierr = DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg);CHKERRQ(ierr); 419 lenglob = appctx->param.E*(appctx->param.N-1); 420 for (i=0; i<lenglob; i++) { 421 s[i]= 0; 422 for (j=0; j<appctx->ncoeff; j++) { 423 tc = -appctx->param.mu*(j+1)*(j+1)*4.0*PETSC_PI*PETSC_PI*t; 424 s[i] += appctx->solutioncoefficients[j]*PetscSinScalar(2*(j+1)*PETSC_PI*(xg[i] + appctx->param.a*t))*PetscExpReal(tc); 425 } 426 } 427 ierr = DMDAVecRestoreArray(appctx->da,obj,&s);CHKERRQ(ierr); 428 ierr = DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg);CHKERRQ(ierr); 429 PetscFunctionReturn(0); 430 } 431 432 PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx) 433 { 434 PetscErrorCode ierr; 435 AppCtx *appctx = (AppCtx*)ctx; 436 437 PetscFunctionBegin; 438 ierr = MatMult(appctx->SEMop.keptstiff,globalin,globalout);CHKERRQ(ierr); 439 PetscFunctionReturn(0); 440 } 441 442 PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec globalin,Mat A, Mat B,void *ctx) 443 { 444 PetscErrorCode ierr; 445 AppCtx *appctx = (AppCtx*)ctx; 446 447 PetscFunctionBegin; 448 ierr = MatCopy(appctx->SEMop.keptstiff,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 449 PetscFunctionReturn(0); 450 } 451 452 /* --------------------------------------------------------------------- */ 453 454 /* 455 RHSLaplacian - matrix for diffusion 456 457 Input Parameters: 458 ts - the TS context 459 t - current time (ignored) 460 X - current solution (ignored) 461 dummy - optional user-defined context, as set by TSetRHSJacobian() 462 463 Output Parameters: 464 AA - Jacobian matrix 465 BB - optionally different matrix from which the preconditioner is built 466 str - flag indicating matrix structure 467 468 Scales by the inverse of the mass matrix (perhaps that should be pulled out) 469 470 */ 471 PetscErrorCode RHSLaplacian(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx) 472 { 473 PetscReal **temp; 474 PetscReal vv; 475 AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 476 PetscErrorCode ierr; 477 PetscInt i,xs,xn,l,j; 478 PetscInt *rowsDM; 479 480 PetscFunctionBegin; 481 /* 482 Creates the element stiffness matrix for the given gll 483 */ 484 ierr = PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr); 485 486 /* scale by the size of the element */ 487 for (i=0; i<appctx->param.N; i++) { 488 vv=-appctx->param.mu*2.0/appctx->param.Le; 489 for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv; 490 } 491 492 ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 493 ierr = DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr); 494 495 if (appctx->param.N-1 < 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Polynomial order must be at least 2"); 496 xs = xs/(appctx->param.N-1); 497 xn = xn/(appctx->param.N-1); 498 499 ierr = PetscMalloc1(appctx->param.N,&rowsDM);CHKERRQ(ierr); 500 /* 501 loop over local elements 502 */ 503 for (j=xs; j<xs+xn; j++) { 504 for (l=0; l<appctx->param.N; l++) { 505 rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l; 506 } 507 ierr = MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES);CHKERRQ(ierr); 508 } 509 ierr = PetscFree(rowsDM);CHKERRQ(ierr); 510 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 511 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 512 ierr = VecReciprocal(appctx->SEMop.mass);CHKERRQ(ierr); 513 ierr = MatDiagonalScale(A,appctx->SEMop.mass,0);CHKERRQ(ierr); 514 ierr = VecReciprocal(appctx->SEMop.mass);CHKERRQ(ierr); 515 516 ierr = PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr); 517 PetscFunctionReturn(0); 518 } 519 520 /* 521 Almost identical to Laplacian 522 523 Note that the element matrix is NOT scaled by the size of element like the Laplacian term. 524 */ 525 PetscErrorCode RHSAdvection(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx) 526 { 527 PetscReal **temp; 528 PetscReal vv; 529 AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 530 PetscErrorCode ierr; 531 PetscInt i,xs,xn,l,j; 532 PetscInt *rowsDM; 533 534 PetscFunctionBegin; 535 /* 536 Creates the element stiffness matrix for the given gll 537 */ 538 ierr = PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr); 539 540 /* scale by the size of the element */ 541 for (i=0; i<appctx->param.N; i++) { 542 vv = -appctx->param.a; 543 for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv; 544 } 545 546 ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 547 ierr = DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr); 548 549 if (appctx->param.N-1 < 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Polynomial order must be at least 2"); 550 xs = xs/(appctx->param.N-1); 551 xn = xn/(appctx->param.N-1); 552 553 ierr = PetscMalloc1(appctx->param.N,&rowsDM);CHKERRQ(ierr); 554 /* 555 loop over local elements 556 */ 557 for (j=xs; j<xs+xn; j++) { 558 for (l=0; l<appctx->param.N; l++) { 559 rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l; 560 } 561 ierr = MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES);CHKERRQ(ierr); 562 } 563 ierr = PetscFree(rowsDM);CHKERRQ(ierr); 564 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 565 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 566 ierr = VecReciprocal(appctx->SEMop.mass);CHKERRQ(ierr); 567 ierr = MatDiagonalScale(A,appctx->SEMop.mass,0);CHKERRQ(ierr); 568 ierr = VecReciprocal(appctx->SEMop.mass);CHKERRQ(ierr); 569 570 ierr = PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr); 571 PetscFunctionReturn(0); 572 } 573 574 /* ------------------------------------------------------------------ */ 575 /* 576 FormFunctionGradient - Evaluates the function and corresponding gradient. 577 578 Input Parameters: 579 tao - the Tao context 580 ic - the input vector 581 ctx - optional user-defined context, as set when calling TaoSetObjectiveAndGradientRoutine() 582 583 Output Parameters: 584 f - the newly evaluated function 585 G - the newly evaluated gradient 586 587 Notes: 588 589 The forward equation is 590 M u_t = F(U) 591 which is converted to 592 u_t = M^{-1} F(u) 593 in the user code since TS has no direct way of providing a mass matrix. The Jacobian of this is 594 M^{-1} J 595 where J is the Jacobian of F. Now the adjoint equation is 596 M v_t = J^T v 597 but TSAdjoint does not solve this since it can only solve the transposed system for the 598 Jacobian the user provided. Hence TSAdjoint solves 599 w_t = J^T M^{-1} w (where w = M v) 600 since there is no way to indicate the mass matrix as a separate entitity to TS. Thus one 601 must be careful in initializing the "adjoint equation" and using the result. This is 602 why 603 G = -2 M(u(T) - u_d) 604 below (instead of -2(u(T) - u_d) 605 606 607 */ 608 PetscErrorCode FormFunctionGradient(Tao tao,Vec ic,PetscReal *f,Vec G,void *ctx) 609 { 610 AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 611 PetscErrorCode ierr; 612 Vec temp; 613 614 PetscFunctionBegin; 615 ierr = TSSetTime(appctx->ts,0.0);CHKERRQ(ierr); 616 ierr = TSSetStepNumber(appctx->ts,0);CHKERRQ(ierr); 617 ierr = TSResetTrajectory(appctx->ts);CHKERRQ(ierr); 618 ierr = TSSetTimeStep(appctx->ts,appctx->initial_dt);CHKERRQ(ierr); 619 ierr = VecCopy(ic,appctx->dat.curr_sol);CHKERRQ(ierr); 620 621 ierr = TSSolve(appctx->ts,appctx->dat.curr_sol);CHKERRQ(ierr); 622 ierr = VecCopy(appctx->dat.curr_sol,appctx->dat.joe);CHKERRQ(ierr); 623 624 /* Compute the difference between the current ODE solution and target ODE solution */ 625 ierr = VecWAXPY(G,-1.0,appctx->dat.curr_sol,appctx->dat.reference);CHKERRQ(ierr); 626 627 /* Compute the objective/cost function */ 628 ierr = VecDuplicate(G,&temp);CHKERRQ(ierr); 629 ierr = VecPointwiseMult(temp,G,G);CHKERRQ(ierr); 630 ierr = VecDot(temp,appctx->SEMop.mass,f);CHKERRQ(ierr); 631 ierr = VecDestroy(&temp);CHKERRQ(ierr); 632 633 /* Compute initial conditions for the adjoint integration. See Notes above */ 634 ierr = VecScale(G, -2.0);CHKERRQ(ierr); 635 ierr = VecPointwiseMult(G,G,appctx->SEMop.mass);CHKERRQ(ierr); 636 ierr = TSSetCostGradients(appctx->ts,1,&G,NULL);CHKERRQ(ierr); 637 638 ierr = TSAdjointSolve(appctx->ts);CHKERRQ(ierr); 639 /* ierr = VecPointwiseDivide(G,G,appctx->SEMop.mass);CHKERRQ(ierr);*/ 640 PetscFunctionReturn(0); 641 } 642 643 PetscErrorCode MonitorError(Tao tao,void *ctx) 644 { 645 AppCtx *appctx = (AppCtx*)ctx; 646 Vec temp,grad; 647 PetscReal nrm; 648 PetscErrorCode ierr; 649 PetscInt its; 650 PetscReal fct,gnorm; 651 652 PetscFunctionBegin; 653 ierr = VecDuplicate(appctx->dat.ic,&temp);CHKERRQ(ierr); 654 ierr = VecWAXPY(temp,-1.0,appctx->dat.ic,appctx->dat.true_solution);CHKERRQ(ierr); 655 ierr = VecPointwiseMult(temp,temp,temp);CHKERRQ(ierr); 656 ierr = VecDot(temp,appctx->SEMop.mass,&nrm);CHKERRQ(ierr); 657 nrm = PetscSqrtReal(nrm); 658 ierr = TaoGetGradientVector(tao,&grad);CHKERRQ(ierr); 659 ierr = VecPointwiseMult(temp,temp,temp);CHKERRQ(ierr); 660 ierr = VecDot(temp,appctx->SEMop.mass,&gnorm);CHKERRQ(ierr); 661 gnorm = PetscSqrtReal(gnorm); 662 ierr = VecDestroy(&temp);CHKERRQ(ierr); 663 ierr = TaoGetIterationNumber(tao,&its);CHKERRQ(ierr); 664 ierr = TaoGetObjective(tao,&fct);CHKERRQ(ierr); 665 if (!its) { 666 ierr = PetscPrintf(PETSC_COMM_WORLD,"%% Iteration Error Objective Gradient-norm\n");CHKERRQ(ierr); 667 ierr = PetscPrintf(PETSC_COMM_WORLD,"history = [\n");CHKERRQ(ierr); 668 } 669 ierr = PetscPrintf(PETSC_COMM_WORLD,"%3D %g %g %g\n",its,(double)nrm,(double)fct,(double)gnorm);CHKERRQ(ierr); 670 PetscFunctionReturn(0); 671 } 672 673 PetscErrorCode MonitorDestroy(void **ctx) 674 { 675 PetscErrorCode ierr; 676 677 PetscFunctionBegin; 678 ierr = PetscPrintf(PETSC_COMM_WORLD,"];\n");CHKERRQ(ierr); 679 PetscFunctionReturn(0); 680 } 681 682 683 /*TEST 684 685 build: 686 requires: !complex 687 688 test: 689 requires: !single 690 args: -ts_adapt_dt_max 3.e-3 -E 10 -N 8 -ncoeff 5 -tao_bqnls_mat_lmvm_scale_type none 691 692 test: 693 suffix: cn 694 requires: !single 695 args: -ts_type cn -ts_dt .003 -pc_type lu -E 10 -N 8 -ncoeff 5 -tao_bqnls_mat_lmvm_scale_type none 696 697 test: 698 suffix: 2 699 requires: !single 700 args: -ts_adapt_dt_max 3.e-3 -E 10 -N 8 -ncoeff 5 -a .1 -tao_bqnls_mat_lmvm_scale_type none 701 702 703 TEST*/ 704