1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] ="Solves a simple data assimilation problem with one dimensional advection diffusion equation using TSAdjoint\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown 6c4762a1bSJed Brown Not yet tested in parallel 7c4762a1bSJed Brown 8c4762a1bSJed Brown */ 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* ------------------------------------------------------------------------ 11c4762a1bSJed Brown 12c4762a1bSJed Brown This program uses the one-dimensional advection-diffusion equation), 13c4762a1bSJed Brown u_t = mu*u_xx - a u_x, 14c4762a1bSJed Brown on the domain 0 <= x <= 1, with periodic boundary conditions 15c4762a1bSJed Brown 16c4762a1bSJed Brown to demonstrate solving a data assimilation problem of finding the initial conditions 17c4762a1bSJed Brown to produce a given solution at a fixed time. 18c4762a1bSJed Brown 19c4762a1bSJed Brown The operators are discretized with the spectral element method 20c4762a1bSJed Brown 21c4762a1bSJed Brown ------------------------------------------------------------------------- */ 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* 24c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this file 25c4762a1bSJed Brown automatically includes: 26c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 27c4762a1bSJed Brown petscmat.h - matrices 28c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 29c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 30c4762a1bSJed Brown petscksp.h - linear solvers petscsnes.h - nonlinear solvers 31c4762a1bSJed Brown */ 32c4762a1bSJed Brown 33c4762a1bSJed Brown #include <petsctao.h> 34c4762a1bSJed Brown #include <petscts.h> 35c4762a1bSJed Brown #include <petscdt.h> 36c4762a1bSJed Brown #include <petscdraw.h> 37c4762a1bSJed Brown #include <petscdmda.h> 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* 40c4762a1bSJed Brown User-defined application context - contains data needed by the 41c4762a1bSJed Brown application-provided call-back routines. 42c4762a1bSJed Brown */ 43c4762a1bSJed Brown 44c4762a1bSJed Brown typedef struct { 45c4762a1bSJed Brown PetscInt n; /* number of nodes */ 46c4762a1bSJed Brown PetscReal *nodes; /* GLL nodes */ 47c4762a1bSJed Brown PetscReal *weights; /* GLL weights */ 48c4762a1bSJed Brown } PetscGLL; 49c4762a1bSJed Brown 50c4762a1bSJed Brown typedef struct { 51c4762a1bSJed Brown PetscInt N; /* grid points per elements*/ 52c4762a1bSJed Brown PetscInt E; /* number of elements */ 53c4762a1bSJed Brown PetscReal tol_L2,tol_max; /* error norms */ 54c4762a1bSJed Brown PetscInt steps; /* number of timesteps */ 55c4762a1bSJed Brown PetscReal Tend; /* endtime */ 56c4762a1bSJed Brown PetscReal mu; /* viscosity */ 57c4762a1bSJed Brown PetscReal a; /* advection speed */ 58c4762a1bSJed Brown PetscReal L; /* total length of domain */ 59c4762a1bSJed Brown PetscReal Le; 60c4762a1bSJed Brown PetscReal Tadj; 61c4762a1bSJed Brown } PetscParam; 62c4762a1bSJed Brown 63c4762a1bSJed Brown typedef struct { 64c4762a1bSJed Brown Vec reference; /* desired end state */ 65c4762a1bSJed Brown Vec grid; /* total grid */ 66c4762a1bSJed Brown Vec grad; 67c4762a1bSJed Brown Vec ic; 68c4762a1bSJed Brown Vec curr_sol; 69c4762a1bSJed Brown Vec joe; 70c4762a1bSJed Brown Vec true_solution; /* actual initial conditions for the final solution */ 71c4762a1bSJed Brown } PetscData; 72c4762a1bSJed Brown 73c4762a1bSJed Brown typedef struct { 74c4762a1bSJed Brown Vec grid; /* total grid */ 75c4762a1bSJed Brown Vec mass; /* mass matrix for total integration */ 76c4762a1bSJed Brown Mat stiff; /* stifness matrix */ 77c4762a1bSJed Brown Mat advec; 78c4762a1bSJed Brown Mat keptstiff; 79c4762a1bSJed Brown PetscGLL gll; 80c4762a1bSJed Brown } PetscSEMOperators; 81c4762a1bSJed Brown 82c4762a1bSJed Brown typedef struct { 83c4762a1bSJed Brown DM da; /* distributed array data structure */ 84c4762a1bSJed Brown PetscSEMOperators SEMop; 85c4762a1bSJed Brown PetscParam param; 86c4762a1bSJed Brown PetscData dat; 87c4762a1bSJed Brown TS ts; 88c4762a1bSJed Brown PetscReal initial_dt; 89c4762a1bSJed Brown PetscReal *solutioncoefficients; 90c4762a1bSJed Brown PetscInt ncoeff; 91c4762a1bSJed Brown } AppCtx; 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* 94c4762a1bSJed Brown User-defined routines 95c4762a1bSJed Brown */ 96c4762a1bSJed Brown extern PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*); 97c4762a1bSJed Brown extern PetscErrorCode RHSLaplacian(TS,PetscReal,Vec,Mat,Mat,void*); 98c4762a1bSJed Brown extern PetscErrorCode RHSAdvection(TS,PetscReal,Vec,Mat,Mat,void*); 99c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*); 100c4762a1bSJed Brown extern PetscErrorCode ComputeReference(TS,PetscReal,Vec,AppCtx*); 101c4762a1bSJed Brown extern PetscErrorCode MonitorError(Tao,void*); 102c4762a1bSJed Brown extern PetscErrorCode MonitorDestroy(void**); 103c4762a1bSJed Brown extern PetscErrorCode ComputeSolutionCoefficients(AppCtx*); 104c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 105c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 106c4762a1bSJed Brown 107c4762a1bSJed Brown int main(int argc,char **argv) 108c4762a1bSJed Brown { 109c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 110c4762a1bSJed Brown Tao tao; 111c4762a1bSJed Brown Vec u; /* approximate solution vector */ 112c4762a1bSJed Brown PetscInt i, xs, xm, ind, j, lenglob; 113c4762a1bSJed Brown PetscReal x, *wrk_ptr1, *wrk_ptr2; 114c4762a1bSJed Brown MatNullSpace nsp; 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117c4762a1bSJed Brown Initialize program and set problem parameters 118c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 119c4762a1bSJed Brown PetscFunctionBegin; 120c4762a1bSJed Brown 121*327415f7SBarry Smith PetscFunctionBeginUser; 1229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 123c4762a1bSJed Brown 124c4762a1bSJed Brown /*initialize parameters */ 125c4762a1bSJed Brown appctx.param.N = 10; /* order of the spectral element */ 126c4762a1bSJed Brown appctx.param.E = 8; /* number of elements */ 127c4762a1bSJed Brown appctx.param.L = 1.0; /* length of the domain */ 128c4762a1bSJed Brown appctx.param.mu = 0.00001; /* diffusion coefficient */ 129c4762a1bSJed Brown appctx.param.a = 0.0; /* advection speed */ 130c4762a1bSJed Brown appctx.initial_dt = 1e-4; 131c4762a1bSJed Brown appctx.param.steps = PETSC_MAX_INT; 132c4762a1bSJed Brown appctx.param.Tend = 0.01; 133c4762a1bSJed Brown appctx.ncoeff = 2; 134c4762a1bSJed Brown 1359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&appctx.param.N,NULL)); 1369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-E",&appctx.param.E,NULL)); 1379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-ncoeff",&appctx.ncoeff,NULL)); 1389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-Tend",&appctx.param.Tend,NULL)); 1399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&appctx.param.mu,NULL)); 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-a",&appctx.param.a,NULL)); 141c4762a1bSJed Brown appctx.param.Le = appctx.param.L/appctx.param.E; 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 144c4762a1bSJed Brown Create GLL data structures 145c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1469566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(appctx.param.N,&appctx.SEMop.gll.nodes,appctx.param.N,&appctx.SEMop.gll.weights)); 1479566063dSJacob Faibussowitsch PetscCall(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N,PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights)); 148c4762a1bSJed Brown appctx.SEMop.gll.n = appctx.param.N; 149c4762a1bSJed Brown lenglob = appctx.param.E*(appctx.param.N-1); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* 152c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 153c4762a1bSJed Brown and to set up the ghost point communication pattern. There are E*(Nl-1)+1 154c4762a1bSJed Brown total grid values spread equally among all the processors, except first and last 155c4762a1bSJed Brown */ 156c4762a1bSJed Brown 1579566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,lenglob,1,1,NULL,&appctx.da)); 1589566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(appctx.da)); 1599566063dSJacob Faibussowitsch PetscCall(DMSetUp(appctx.da)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* 162c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the 163c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that 164c4762a1bSJed Brown have the same types. 165c4762a1bSJed Brown */ 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(appctx.da,&u)); 1689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&appctx.dat.ic)); 1699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&appctx.dat.true_solution)); 1709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&appctx.dat.reference)); 1719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&appctx.SEMop.grid)); 1729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&appctx.SEMop.mass)); 1739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&appctx.dat.curr_sol)); 1749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&appctx.dat.joe)); 175c4762a1bSJed Brown 1769566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx.da,&xs,NULL,NULL,&xm,NULL,NULL)); 1779566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1)); 1789566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 181c4762a1bSJed Brown 182c4762a1bSJed Brown xs=xs/(appctx.param.N-1); 183c4762a1bSJed Brown xm=xm/(appctx.param.N-1); 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* 186c4762a1bSJed Brown Build total grid and mass over entire mesh (multi-elemental) 187c4762a1bSJed Brown */ 188c4762a1bSJed Brown 189c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 190c4762a1bSJed Brown for (j=0; j<appctx.param.N-1; j++) { 191c4762a1bSJed Brown x = (appctx.param.Le/2.0)*(appctx.SEMop.gll.nodes[j]+1.0)+appctx.param.Le*i; 192c4762a1bSJed Brown ind=i*(appctx.param.N-1)+j; 193c4762a1bSJed Brown wrk_ptr1[ind]=x; 194c4762a1bSJed Brown wrk_ptr2[ind]=.5*appctx.param.Le*appctx.SEMop.gll.weights[j]; 195c4762a1bSJed Brown if (j==0) wrk_ptr2[ind]+=.5*appctx.param.Le*appctx.SEMop.gll.weights[j]; 196c4762a1bSJed Brown } 197c4762a1bSJed Brown } 1989566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1)); 1999566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2)); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 202c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 203c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2049566063dSJacob Faibussowitsch PetscCall(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE)); 2059566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(appctx.da,&appctx.SEMop.stiff)); 2069566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(appctx.da,&appctx.SEMop.advec)); 207c4762a1bSJed Brown 208c4762a1bSJed Brown /* 209c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 210c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 211c4762a1bSJed Brown as a time-dependent matrix. 212c4762a1bSJed Brown */ 2139566063dSJacob Faibussowitsch PetscCall(RHSLaplacian(appctx.ts,0.0,u,appctx.SEMop.stiff,appctx.SEMop.stiff,&appctx)); 2149566063dSJacob Faibussowitsch PetscCall(RHSAdvection(appctx.ts,0.0,u,appctx.SEMop.advec,appctx.SEMop.advec,&appctx)); 2159566063dSJacob Faibussowitsch PetscCall(MatAXPY(appctx.SEMop.stiff,-1.0,appctx.SEMop.advec,DIFFERENT_NONZERO_PATTERN)); 2169566063dSJacob Faibussowitsch PetscCall(MatDuplicate(appctx.SEMop.stiff,MAT_COPY_VALUES,&appctx.SEMop.keptstiff)); 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */ 2199566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nsp)); 2209566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.stiff,nsp)); 2219566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(nsp,appctx.SEMop.stiff,NULL)); 2229566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* Create the TS solver that solves the ODE and its adjoint; set its options */ 2259566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&appctx.ts)); 2269566063dSJacob Faibussowitsch PetscCall(TSSetSolutionFunction(appctx.ts,(PetscErrorCode (*)(TS,PetscReal,Vec, void *))ComputeReference,&appctx)); 2279566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(appctx.ts,TS_LINEAR)); 2289566063dSJacob Faibussowitsch PetscCall(TSSetType(appctx.ts,TSRK)); 2299566063dSJacob Faibussowitsch PetscCall(TSSetDM(appctx.ts,appctx.da)); 2309566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx.ts,0.0)); 2319566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(appctx.ts,appctx.initial_dt)); 2329566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(appctx.ts,appctx.param.steps)); 2339566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(appctx.ts,appctx.param.Tend)); 2349566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(appctx.ts,TS_EXACTFINALTIME_MATCHSTEP)); 2359566063dSJacob Faibussowitsch PetscCall(TSSetTolerances(appctx.ts,1e-7,NULL,1e-7,NULL)); 2369566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx.ts)); 237c4762a1bSJed Brown /* Need to save initial timestep user may have set with -ts_dt so it can be reset for each new TSSolve() */ 2389566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(appctx.ts,&appctx.initial_dt)); 2399566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(appctx.ts,NULL,TSComputeRHSFunctionLinear,&appctx)); 2409566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,TSComputeRHSJacobianConstant,&appctx)); 2419566063dSJacob Faibussowitsch /* PetscCall(TSSetRHSFunction(appctx.ts,NULL,RHSFunction,&appctx)); 2429566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,RHSJacobian,&appctx)); */ 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* Set random initial conditions as initial guess, compute analytic reference solution and analytic (true) initial conditions */ 2459566063dSJacob Faibussowitsch PetscCall(ComputeSolutionCoefficients(&appctx)); 2469566063dSJacob Faibussowitsch PetscCall(InitialConditions(appctx.dat.ic,&appctx)); 2479566063dSJacob Faibussowitsch PetscCall(ComputeReference(appctx.ts,appctx.param.Tend,appctx.dat.reference,&appctx)); 2489566063dSJacob Faibussowitsch PetscCall(ComputeReference(appctx.ts,0.0,appctx.dat.true_solution,&appctx)); 249c4762a1bSJed Brown 250f32d6360SSatish Balay /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */ 2519566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(appctx.ts)); 2529566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx.ts)); 253f32d6360SSatish Balay 254c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 2559566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao)); 2569566063dSJacob Faibussowitsch PetscCall(TaoSetMonitor(tao,MonitorError,&appctx,MonitorDestroy)); 2579566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOBQNLS)); 2589566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,appctx.dat.ic)); 259c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 2609566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&appctx)); 261c4762a1bSJed Brown /* Check for any TAO command line options */ 2629566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(tao,1e-8,PETSC_DEFAULT,PETSC_DEFAULT)); 2639566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 2649566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 265c4762a1bSJed Brown 2669566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 2679566063dSJacob Faibussowitsch PetscCall(PetscFree(appctx.solutioncoefficients)); 2689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.advec)); 2699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.stiff)); 2709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.keptstiff)); 2719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 2729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.ic)); 2739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.joe)); 2749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.true_solution)); 2759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.reference)); 2769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.SEMop.grid)); 2779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.SEMop.mass)); 2789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.curr_sol)); 2799566063dSJacob Faibussowitsch PetscCall(PetscFree2(appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights)); 2809566063dSJacob Faibussowitsch PetscCall(DMDestroy(&appctx.da)); 2819566063dSJacob Faibussowitsch PetscCall(TSDestroy(&appctx.ts)); 282c4762a1bSJed Brown 283c4762a1bSJed Brown /* 284c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 285c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 286c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 287c4762a1bSJed Brown options are chosen (e.g., -log_summary). 288c4762a1bSJed Brown */ 2899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 290b122ec5aSJacob Faibussowitsch return 0; 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293c4762a1bSJed Brown /* 294c4762a1bSJed Brown Computes the coefficients for the analytic solution to the PDE 295c4762a1bSJed Brown */ 296c4762a1bSJed Brown PetscErrorCode ComputeSolutionCoefficients(AppCtx *appctx) 297c4762a1bSJed Brown { 298c4762a1bSJed Brown PetscRandom rand; 299c4762a1bSJed Brown PetscInt i; 300c4762a1bSJed Brown 301c4762a1bSJed Brown PetscFunctionBegin; 3029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->ncoeff,&appctx->solutioncoefficients)); 3039566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 3049566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand,.9,1.0)); 305c4762a1bSJed Brown for (i=0; i<appctx->ncoeff; i++) { 3069566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&appctx->solutioncoefficients[i])); 307c4762a1bSJed Brown } 3089566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 309c4762a1bSJed Brown PetscFunctionReturn(0); 310c4762a1bSJed Brown } 311c4762a1bSJed Brown 312c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 313c4762a1bSJed Brown /* 314c4762a1bSJed Brown InitialConditions - Computes the (random) initial conditions for the Tao optimization solve (these are also initial conditions for the first TSSolve() 315c4762a1bSJed Brown 316c4762a1bSJed Brown Input Parameter: 317c4762a1bSJed Brown u - uninitialized solution vector (global) 318c4762a1bSJed Brown appctx - user-defined application context 319c4762a1bSJed Brown 320c4762a1bSJed Brown Output Parameter: 321c4762a1bSJed Brown u - vector with solution at initial time (global) 322c4762a1bSJed Brown */ 323c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 324c4762a1bSJed Brown { 325c4762a1bSJed Brown PetscScalar *s; 326c4762a1bSJed Brown const PetscScalar *xg; 327c4762a1bSJed Brown PetscInt i,j,lenglob; 328c4762a1bSJed Brown PetscReal sum,val; 329c4762a1bSJed Brown PetscRandom rand; 330c4762a1bSJed Brown 331c4762a1bSJed Brown PetscFunctionBegin; 3329566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 3339566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand,.9,1.0)); 3349566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da,u,&s)); 3359566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 336c4762a1bSJed Brown lenglob = appctx->param.E*(appctx->param.N-1); 337c4762a1bSJed Brown for (i=0; i<lenglob; i++) { 338c4762a1bSJed Brown s[i]= 0; 339c4762a1bSJed Brown for (j=0; j<appctx->ncoeff; j++) { 3409566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&val)); 341c4762a1bSJed Brown s[i] += val*PetscSinScalar(2*(j+1)*PETSC_PI*xg[i]); 342c4762a1bSJed Brown } 343c4762a1bSJed Brown } 3449566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 3459566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da,u,&s)); 3469566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 347c4762a1bSJed Brown /* make sure initial conditions do not contain the constant functions, since with periodic boundary conditions the constant functions introduce a null space */ 3489566063dSJacob Faibussowitsch PetscCall(VecSum(u,&sum)); 3499566063dSJacob Faibussowitsch PetscCall(VecShift(u,-sum/lenglob)); 350c4762a1bSJed Brown PetscFunctionReturn(0); 351c4762a1bSJed Brown } 352c4762a1bSJed Brown 353c4762a1bSJed Brown /* 354c4762a1bSJed Brown TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function. 355c4762a1bSJed Brown 356a5b23f4aSJose E. Roman InitialConditions() computes the initial conditions for the beginning of the Tao iterations 357c4762a1bSJed Brown 358c4762a1bSJed Brown Input Parameter: 359c4762a1bSJed Brown u - uninitialized solution vector (global) 360c4762a1bSJed Brown appctx - user-defined application context 361c4762a1bSJed Brown 362c4762a1bSJed Brown Output Parameter: 363c4762a1bSJed Brown u - vector with solution at initial time (global) 364c4762a1bSJed Brown */ 365c4762a1bSJed Brown PetscErrorCode TrueSolution(Vec u,AppCtx *appctx) 366c4762a1bSJed Brown { 367c4762a1bSJed Brown PetscScalar *s; 368c4762a1bSJed Brown const PetscScalar *xg; 369c4762a1bSJed Brown PetscInt i,j,lenglob; 370c4762a1bSJed Brown PetscReal sum; 371c4762a1bSJed Brown 372c4762a1bSJed Brown PetscFunctionBegin; 3739566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da,u,&s)); 3749566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 375c4762a1bSJed Brown lenglob = appctx->param.E*(appctx->param.N-1); 376c4762a1bSJed Brown for (i=0; i<lenglob; i++) { 377c4762a1bSJed Brown s[i]= 0; 378c4762a1bSJed Brown for (j=0; j<appctx->ncoeff; j++) { 379c4762a1bSJed Brown s[i] += appctx->solutioncoefficients[j]*PetscSinScalar(2*(j+1)*PETSC_PI*xg[i]); 380c4762a1bSJed Brown } 381c4762a1bSJed Brown } 3829566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da,u,&s)); 3839566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 384c4762a1bSJed Brown /* make sure initial conditions do not contain the constant functions, since with periodic boundary conditions the constant functions introduce a null space */ 3859566063dSJacob Faibussowitsch PetscCall(VecSum(u,&sum)); 3869566063dSJacob Faibussowitsch PetscCall(VecShift(u,-sum/lenglob)); 387c4762a1bSJed Brown PetscFunctionReturn(0); 388c4762a1bSJed Brown } 389c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 390c4762a1bSJed Brown /* 391c4762a1bSJed Brown Sets the desired profile for the final end time 392c4762a1bSJed Brown 393c4762a1bSJed Brown Input Parameters: 394c4762a1bSJed Brown t - final time 395c4762a1bSJed Brown obj - vector storing the desired profile 396c4762a1bSJed Brown appctx - user-defined application context 397c4762a1bSJed Brown 398c4762a1bSJed Brown */ 399c4762a1bSJed Brown PetscErrorCode ComputeReference(TS ts,PetscReal t,Vec obj,AppCtx *appctx) 400c4762a1bSJed Brown { 401c4762a1bSJed Brown PetscScalar *s,tc; 402c4762a1bSJed Brown const PetscScalar *xg; 403c4762a1bSJed Brown PetscInt i, j,lenglob; 404c4762a1bSJed Brown 405c4762a1bSJed Brown PetscFunctionBegin; 4069566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da,obj,&s)); 4079566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 408c4762a1bSJed Brown lenglob = appctx->param.E*(appctx->param.N-1); 409c4762a1bSJed Brown for (i=0; i<lenglob; i++) { 410c4762a1bSJed Brown s[i]= 0; 411c4762a1bSJed Brown for (j=0; j<appctx->ncoeff; j++) { 412c4762a1bSJed Brown tc = -appctx->param.mu*(j+1)*(j+1)*4.0*PETSC_PI*PETSC_PI*t; 413c4762a1bSJed Brown s[i] += appctx->solutioncoefficients[j]*PetscSinScalar(2*(j+1)*PETSC_PI*(xg[i] + appctx->param.a*t))*PetscExpReal(tc); 414c4762a1bSJed Brown } 415c4762a1bSJed Brown } 4169566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da,obj,&s)); 4179566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 418c4762a1bSJed Brown PetscFunctionReturn(0); 419c4762a1bSJed Brown } 420c4762a1bSJed Brown 421c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx) 422c4762a1bSJed Brown { 423c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 424c4762a1bSJed Brown 425c4762a1bSJed Brown PetscFunctionBegin; 4269566063dSJacob Faibussowitsch PetscCall(MatMult(appctx->SEMop.keptstiff,globalin,globalout)); 427c4762a1bSJed Brown PetscFunctionReturn(0); 428c4762a1bSJed Brown } 429c4762a1bSJed Brown 430c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec globalin,Mat A, Mat B,void *ctx) 431c4762a1bSJed Brown { 432c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 433c4762a1bSJed Brown 434c4762a1bSJed Brown PetscFunctionBegin; 4359566063dSJacob Faibussowitsch PetscCall(MatCopy(appctx->SEMop.keptstiff,A,DIFFERENT_NONZERO_PATTERN)); 436c4762a1bSJed Brown PetscFunctionReturn(0); 437c4762a1bSJed Brown } 438c4762a1bSJed Brown 439c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 440c4762a1bSJed Brown 441c4762a1bSJed Brown /* 442c4762a1bSJed Brown RHSLaplacian - matrix for diffusion 443c4762a1bSJed Brown 444c4762a1bSJed Brown Input Parameters: 445c4762a1bSJed Brown ts - the TS context 446c4762a1bSJed Brown t - current time (ignored) 447c4762a1bSJed Brown X - current solution (ignored) 448c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 449c4762a1bSJed Brown 450c4762a1bSJed Brown Output Parameters: 451c4762a1bSJed Brown AA - Jacobian matrix 452c4762a1bSJed Brown BB - optionally different matrix from which the preconditioner is built 453c4762a1bSJed Brown str - flag indicating matrix structure 454c4762a1bSJed Brown 455c4762a1bSJed Brown Scales by the inverse of the mass matrix (perhaps that should be pulled out) 456c4762a1bSJed Brown 457c4762a1bSJed Brown */ 458c4762a1bSJed Brown PetscErrorCode RHSLaplacian(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx) 459c4762a1bSJed Brown { 460c4762a1bSJed Brown PetscReal **temp; 461c4762a1bSJed Brown PetscReal vv; 462c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 463c4762a1bSJed Brown PetscInt i,xs,xn,l,j; 464c4762a1bSJed Brown PetscInt *rowsDM; 465c4762a1bSJed Brown 466c4762a1bSJed Brown PetscFunctionBegin; 467c4762a1bSJed Brown /* 468c4762a1bSJed Brown Creates the element stiffness matrix for the given gll 469c4762a1bSJed Brown */ 4709566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 471c4762a1bSJed Brown 472c4762a1bSJed Brown /* scale by the size of the element */ 473c4762a1bSJed Brown for (i=0; i<appctx->param.N; i++) { 474c4762a1bSJed Brown vv=-appctx->param.mu*2.0/appctx->param.Le; 475c4762a1bSJed Brown for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv; 476c4762a1bSJed Brown } 477c4762a1bSJed Brown 4789566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 4799566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL)); 480c4762a1bSJed Brown 4813c859ba3SBarry Smith PetscCheck(appctx->param.N-1 >= 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Polynomial order must be at least 2"); 482c4762a1bSJed Brown xs = xs/(appctx->param.N-1); 483c4762a1bSJed Brown xn = xn/(appctx->param.N-1); 484c4762a1bSJed Brown 4859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->param.N,&rowsDM)); 486c4762a1bSJed Brown /* 487c4762a1bSJed Brown loop over local elements 488c4762a1bSJed Brown */ 489c4762a1bSJed Brown for (j=xs; j<xs+xn; j++) { 490c4762a1bSJed Brown for (l=0; l<appctx->param.N; l++) { 491c4762a1bSJed Brown rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l; 492c4762a1bSJed Brown } 4939566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES)); 494c4762a1bSJed Brown } 4959566063dSJacob Faibussowitsch PetscCall(PetscFree(rowsDM)); 4969566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 4979566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 4989566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 4999566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,appctx->SEMop.mass,0)); 5009566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 501c4762a1bSJed Brown 5029566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 503c4762a1bSJed Brown PetscFunctionReturn(0); 504c4762a1bSJed Brown } 505c4762a1bSJed Brown 506c4762a1bSJed Brown /* 507c4762a1bSJed Brown Almost identical to Laplacian 508c4762a1bSJed Brown 509c4762a1bSJed Brown Note that the element matrix is NOT scaled by the size of element like the Laplacian term. 510c4762a1bSJed Brown */ 511c4762a1bSJed Brown PetscErrorCode RHSAdvection(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx) 512c4762a1bSJed Brown { 513c4762a1bSJed Brown PetscReal **temp; 514c4762a1bSJed Brown PetscReal vv; 515c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 516c4762a1bSJed Brown PetscInt i,xs,xn,l,j; 517c4762a1bSJed Brown PetscInt *rowsDM; 518c4762a1bSJed Brown 519c4762a1bSJed Brown PetscFunctionBegin; 520c4762a1bSJed Brown /* 521c4762a1bSJed Brown Creates the element stiffness matrix for the given gll 522c4762a1bSJed Brown */ 5239566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 524c4762a1bSJed Brown 525c4762a1bSJed Brown /* scale by the size of the element */ 526c4762a1bSJed Brown for (i=0; i<appctx->param.N; i++) { 527c4762a1bSJed Brown vv = -appctx->param.a; 528c4762a1bSJed Brown for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv; 529c4762a1bSJed Brown } 530c4762a1bSJed Brown 5319566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 5329566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL)); 533c4762a1bSJed Brown 5343c859ba3SBarry Smith PetscCheck(appctx->param.N-1 >= 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Polynomial order must be at least 2"); 535c4762a1bSJed Brown xs = xs/(appctx->param.N-1); 536c4762a1bSJed Brown xn = xn/(appctx->param.N-1); 537c4762a1bSJed Brown 5389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->param.N,&rowsDM)); 539c4762a1bSJed Brown /* 540c4762a1bSJed Brown loop over local elements 541c4762a1bSJed Brown */ 542c4762a1bSJed Brown for (j=xs; j<xs+xn; j++) { 543c4762a1bSJed Brown for (l=0; l<appctx->param.N; l++) { 544c4762a1bSJed Brown rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l; 545c4762a1bSJed Brown } 5469566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES)); 547c4762a1bSJed Brown } 5489566063dSJacob Faibussowitsch PetscCall(PetscFree(rowsDM)); 5499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 5509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 5519566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 5529566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,appctx->SEMop.mass,0)); 5539566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 554c4762a1bSJed Brown 5559566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 556c4762a1bSJed Brown PetscFunctionReturn(0); 557c4762a1bSJed Brown } 558c4762a1bSJed Brown 559c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 560c4762a1bSJed Brown /* 561c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 562c4762a1bSJed Brown 563c4762a1bSJed Brown Input Parameters: 564c4762a1bSJed Brown tao - the Tao context 565c4762a1bSJed Brown ic - the input vector 566a82e8c82SStefano Zampini ctx - optional user-defined context, as set when calling TaoSetObjectiveAndGradient() 567c4762a1bSJed Brown 568c4762a1bSJed Brown Output Parameters: 569c4762a1bSJed Brown f - the newly evaluated function 570c4762a1bSJed Brown G - the newly evaluated gradient 571c4762a1bSJed Brown 572c4762a1bSJed Brown Notes: 573c4762a1bSJed Brown 574c4762a1bSJed Brown The forward equation is 575c4762a1bSJed Brown M u_t = F(U) 576c4762a1bSJed Brown which is converted to 577c4762a1bSJed Brown u_t = M^{-1} F(u) 578c4762a1bSJed Brown in the user code since TS has no direct way of providing a mass matrix. The Jacobian of this is 579c4762a1bSJed Brown M^{-1} J 580c4762a1bSJed Brown where J is the Jacobian of F. Now the adjoint equation is 581c4762a1bSJed Brown M v_t = J^T v 582c4762a1bSJed Brown but TSAdjoint does not solve this since it can only solve the transposed system for the 583c4762a1bSJed Brown Jacobian the user provided. Hence TSAdjoint solves 584c4762a1bSJed Brown w_t = J^T M^{-1} w (where w = M v) 585a5b23f4aSJose E. Roman since there is no way to indicate the mass matrix as a separate entity to TS. Thus one 586c4762a1bSJed Brown must be careful in initializing the "adjoint equation" and using the result. This is 587c4762a1bSJed Brown why 588c4762a1bSJed Brown G = -2 M(u(T) - u_d) 589c4762a1bSJed Brown below (instead of -2(u(T) - u_d) 590c4762a1bSJed Brown 591c4762a1bSJed Brown */ 592c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec ic,PetscReal *f,Vec G,void *ctx) 593c4762a1bSJed Brown { 594c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 595c4762a1bSJed Brown Vec temp; 596c4762a1bSJed Brown 597c4762a1bSJed Brown PetscFunctionBegin; 5989566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx->ts,0.0)); 5999566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(appctx->ts,0)); 6009566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(appctx->ts,appctx->initial_dt)); 6019566063dSJacob Faibussowitsch PetscCall(VecCopy(ic,appctx->dat.curr_sol)); 602c4762a1bSJed Brown 6039566063dSJacob Faibussowitsch PetscCall(TSSolve(appctx->ts,appctx->dat.curr_sol)); 6049566063dSJacob Faibussowitsch PetscCall(VecCopy(appctx->dat.curr_sol,appctx->dat.joe)); 605c4762a1bSJed Brown 606c4762a1bSJed Brown /* Compute the difference between the current ODE solution and target ODE solution */ 6079566063dSJacob Faibussowitsch PetscCall(VecWAXPY(G,-1.0,appctx->dat.curr_sol,appctx->dat.reference)); 608c4762a1bSJed Brown 609c4762a1bSJed Brown /* Compute the objective/cost function */ 6109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(G,&temp)); 6119566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(temp,G,G)); 6129566063dSJacob Faibussowitsch PetscCall(VecDot(temp,appctx->SEMop.mass,f)); 6139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&temp)); 614c4762a1bSJed Brown 615c4762a1bSJed Brown /* Compute initial conditions for the adjoint integration. See Notes above */ 6169566063dSJacob Faibussowitsch PetscCall(VecScale(G, -2.0)); 6179566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(G,G,appctx->SEMop.mass)); 6189566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(appctx->ts,1,&G,NULL)); 619c4762a1bSJed Brown 6209566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(appctx->ts)); 6219566063dSJacob Faibussowitsch /* PetscCall(VecPointwiseDivide(G,G,appctx->SEMop.mass));*/ 622c4762a1bSJed Brown PetscFunctionReturn(0); 623c4762a1bSJed Brown } 624c4762a1bSJed Brown 625c4762a1bSJed Brown PetscErrorCode MonitorError(Tao tao,void *ctx) 626c4762a1bSJed Brown { 627c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 628c4762a1bSJed Brown Vec temp,grad; 629c4762a1bSJed Brown PetscReal nrm; 630c4762a1bSJed Brown PetscInt its; 631c4762a1bSJed Brown PetscReal fct,gnorm; 632c4762a1bSJed Brown 633c4762a1bSJed Brown PetscFunctionBegin; 6349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx->dat.ic,&temp)); 6359566063dSJacob Faibussowitsch PetscCall(VecWAXPY(temp,-1.0,appctx->dat.ic,appctx->dat.true_solution)); 6369566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(temp,temp,temp)); 6379566063dSJacob Faibussowitsch PetscCall(VecDot(temp,appctx->SEMop.mass,&nrm)); 638c4762a1bSJed Brown nrm = PetscSqrtReal(nrm); 6399566063dSJacob Faibussowitsch PetscCall(TaoGetGradient(tao,&grad,NULL,NULL)); 6409566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(temp,temp,temp)); 6419566063dSJacob Faibussowitsch PetscCall(VecDot(temp,appctx->SEMop.mass,&gnorm)); 642c4762a1bSJed Brown gnorm = PetscSqrtReal(gnorm); 6439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&temp)); 6449566063dSJacob Faibussowitsch PetscCall(TaoGetIterationNumber(tao,&its)); 6459566063dSJacob Faibussowitsch PetscCall(TaoGetSolutionStatus(tao,NULL,&fct,NULL,NULL,NULL,NULL)); 646c4762a1bSJed Brown if (!its) { 6479566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%% Iteration Error Objective Gradient-norm\n")); 6489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"history = [\n")); 649c4762a1bSJed Brown } 65063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%3" PetscInt_FMT " %g %g %g\n",its,(double)nrm,(double)fct,(double)gnorm)); 651c4762a1bSJed Brown PetscFunctionReturn(0); 652c4762a1bSJed Brown } 653c4762a1bSJed Brown 654c4762a1bSJed Brown PetscErrorCode MonitorDestroy(void **ctx) 655c4762a1bSJed Brown { 656c4762a1bSJed Brown PetscFunctionBegin; 6579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"];\n")); 658c4762a1bSJed Brown PetscFunctionReturn(0); 659c4762a1bSJed Brown } 660c4762a1bSJed Brown 661c4762a1bSJed Brown /*TEST 662c4762a1bSJed Brown 663c4762a1bSJed Brown build: 664c4762a1bSJed Brown requires: !complex 665c4762a1bSJed Brown 666c4762a1bSJed Brown test: 667c4762a1bSJed Brown requires: !single 668c4762a1bSJed Brown args: -ts_adapt_dt_max 3.e-3 -E 10 -N 8 -ncoeff 5 -tao_bqnls_mat_lmvm_scale_type none 669c4762a1bSJed Brown 670c4762a1bSJed Brown test: 671c4762a1bSJed Brown suffix: cn 672c4762a1bSJed Brown requires: !single 673c4762a1bSJed Brown args: -ts_type cn -ts_dt .003 -pc_type lu -E 10 -N 8 -ncoeff 5 -tao_bqnls_mat_lmvm_scale_type none 674c4762a1bSJed Brown 675c4762a1bSJed Brown test: 676c4762a1bSJed Brown suffix: 2 677c4762a1bSJed Brown requires: !single 678c4762a1bSJed Brown args: -ts_adapt_dt_max 3.e-3 -E 10 -N 8 -ncoeff 5 -a .1 -tao_bqnls_mat_lmvm_scale_type none 679c4762a1bSJed Brown 680c4762a1bSJed Brown TEST*/ 681