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 Concepts: TS^time-dependent linear problems 11c4762a1bSJed Brown Concepts: TS^heat equation 12c4762a1bSJed Brown Concepts: TS^diffusion equation 13c4762a1bSJed Brown Concepts: adjoints 14c4762a1bSJed Brown Processors: n 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* ------------------------------------------------------------------------ 18c4762a1bSJed Brown 19c4762a1bSJed Brown This program uses the one-dimensional advection-diffusion equation), 20c4762a1bSJed Brown u_t = mu*u_xx - a u_x, 21c4762a1bSJed Brown on the domain 0 <= x <= 1, with periodic boundary conditions 22c4762a1bSJed Brown 23c4762a1bSJed Brown to demonstrate solving a data assimilation problem of finding the initial conditions 24c4762a1bSJed Brown to produce a given solution at a fixed time. 25c4762a1bSJed Brown 26c4762a1bSJed Brown The operators are discretized with the spectral element method 27c4762a1bSJed Brown 28c4762a1bSJed Brown ------------------------------------------------------------------------- */ 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* 31c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this file 32c4762a1bSJed Brown automatically includes: 33c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 34c4762a1bSJed Brown petscmat.h - matrices 35c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 36c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 37c4762a1bSJed Brown petscksp.h - linear solvers petscsnes.h - nonlinear solvers 38c4762a1bSJed Brown */ 39c4762a1bSJed Brown 40c4762a1bSJed Brown #include <petsctao.h> 41c4762a1bSJed Brown #include <petscts.h> 42c4762a1bSJed Brown #include <petscdt.h> 43c4762a1bSJed Brown #include <petscdraw.h> 44c4762a1bSJed Brown #include <petscdmda.h> 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* 47c4762a1bSJed Brown User-defined application context - contains data needed by the 48c4762a1bSJed Brown application-provided call-back routines. 49c4762a1bSJed Brown */ 50c4762a1bSJed Brown 51c4762a1bSJed Brown typedef struct { 52c4762a1bSJed Brown PetscInt n; /* number of nodes */ 53c4762a1bSJed Brown PetscReal *nodes; /* GLL nodes */ 54c4762a1bSJed Brown PetscReal *weights; /* GLL weights */ 55c4762a1bSJed Brown } PetscGLL; 56c4762a1bSJed Brown 57c4762a1bSJed Brown typedef struct { 58c4762a1bSJed Brown PetscInt N; /* grid points per elements*/ 59c4762a1bSJed Brown PetscInt E; /* number of elements */ 60c4762a1bSJed Brown PetscReal tol_L2,tol_max; /* error norms */ 61c4762a1bSJed Brown PetscInt steps; /* number of timesteps */ 62c4762a1bSJed Brown PetscReal Tend; /* endtime */ 63c4762a1bSJed Brown PetscReal mu; /* viscosity */ 64c4762a1bSJed Brown PetscReal a; /* advection speed */ 65c4762a1bSJed Brown PetscReal L; /* total length of domain */ 66c4762a1bSJed Brown PetscReal Le; 67c4762a1bSJed Brown PetscReal Tadj; 68c4762a1bSJed Brown } PetscParam; 69c4762a1bSJed Brown 70c4762a1bSJed Brown typedef struct { 71c4762a1bSJed Brown Vec reference; /* desired end state */ 72c4762a1bSJed Brown Vec grid; /* total grid */ 73c4762a1bSJed Brown Vec grad; 74c4762a1bSJed Brown Vec ic; 75c4762a1bSJed Brown Vec curr_sol; 76c4762a1bSJed Brown Vec joe; 77c4762a1bSJed Brown Vec true_solution; /* actual initial conditions for the final solution */ 78c4762a1bSJed Brown } PetscData; 79c4762a1bSJed Brown 80c4762a1bSJed Brown typedef struct { 81c4762a1bSJed Brown Vec grid; /* total grid */ 82c4762a1bSJed Brown Vec mass; /* mass matrix for total integration */ 83c4762a1bSJed Brown Mat stiff; /* stifness matrix */ 84c4762a1bSJed Brown Mat advec; 85c4762a1bSJed Brown Mat keptstiff; 86c4762a1bSJed Brown PetscGLL gll; 87c4762a1bSJed Brown } PetscSEMOperators; 88c4762a1bSJed Brown 89c4762a1bSJed Brown typedef struct { 90c4762a1bSJed Brown DM da; /* distributed array data structure */ 91c4762a1bSJed Brown PetscSEMOperators SEMop; 92c4762a1bSJed Brown PetscParam param; 93c4762a1bSJed Brown PetscData dat; 94c4762a1bSJed Brown TS ts; 95c4762a1bSJed Brown PetscReal initial_dt; 96c4762a1bSJed Brown PetscReal *solutioncoefficients; 97c4762a1bSJed Brown PetscInt ncoeff; 98c4762a1bSJed Brown } AppCtx; 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* 101c4762a1bSJed Brown User-defined routines 102c4762a1bSJed Brown */ 103c4762a1bSJed Brown extern PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*); 104c4762a1bSJed Brown extern PetscErrorCode RHSLaplacian(TS,PetscReal,Vec,Mat,Mat,void*); 105c4762a1bSJed Brown extern PetscErrorCode RHSAdvection(TS,PetscReal,Vec,Mat,Mat,void*); 106c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*); 107c4762a1bSJed Brown extern PetscErrorCode ComputeReference(TS,PetscReal,Vec,AppCtx*); 108c4762a1bSJed Brown extern PetscErrorCode MonitorError(Tao,void*); 109c4762a1bSJed Brown extern PetscErrorCode MonitorDestroy(void**); 110c4762a1bSJed Brown extern PetscErrorCode ComputeSolutionCoefficients(AppCtx*); 111c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 112c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 113c4762a1bSJed Brown 114c4762a1bSJed Brown int main(int argc,char **argv) 115c4762a1bSJed Brown { 116c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 117c4762a1bSJed Brown Tao tao; 118c4762a1bSJed Brown Vec u; /* approximate solution vector */ 119c4762a1bSJed Brown PetscInt i, xs, xm, ind, j, lenglob; 120c4762a1bSJed Brown PetscReal x, *wrk_ptr1, *wrk_ptr2; 121c4762a1bSJed Brown MatNullSpace nsp; 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124c4762a1bSJed Brown Initialize program and set problem parameters 125c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 126c4762a1bSJed Brown PetscFunctionBegin; 127c4762a1bSJed Brown 128*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /*initialize parameters */ 131c4762a1bSJed Brown appctx.param.N = 10; /* order of the spectral element */ 132c4762a1bSJed Brown appctx.param.E = 8; /* number of elements */ 133c4762a1bSJed Brown appctx.param.L = 1.0; /* length of the domain */ 134c4762a1bSJed Brown appctx.param.mu = 0.00001; /* diffusion coefficient */ 135c4762a1bSJed Brown appctx.param.a = 0.0; /* advection speed */ 136c4762a1bSJed Brown appctx.initial_dt = 1e-4; 137c4762a1bSJed Brown appctx.param.steps = PETSC_MAX_INT; 138c4762a1bSJed Brown appctx.param.Tend = 0.01; 139c4762a1bSJed Brown appctx.ncoeff = 2; 140c4762a1bSJed Brown 1415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&appctx.param.N,NULL)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-E",&appctx.param.E,NULL)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ncoeff",&appctx.ncoeff,NULL)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-Tend",&appctx.param.Tend,NULL)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&appctx.param.mu,NULL)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-a",&appctx.param.a,NULL)); 147c4762a1bSJed Brown appctx.param.Le = appctx.param.L/appctx.param.E; 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 150c4762a1bSJed Brown Create GLL data structures 151c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(appctx.param.N,&appctx.SEMop.gll.nodes,appctx.param.N,&appctx.SEMop.gll.weights)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N,PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights)); 154c4762a1bSJed Brown appctx.SEMop.gll.n = appctx.param.N; 155c4762a1bSJed Brown lenglob = appctx.param.E*(appctx.param.N-1); 156c4762a1bSJed Brown 157c4762a1bSJed Brown /* 158c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 159c4762a1bSJed Brown and to set up the ghost point communication pattern. There are E*(Nl-1)+1 160c4762a1bSJed Brown total grid values spread equally among all the processors, except first and last 161c4762a1bSJed Brown */ 162c4762a1bSJed Brown 1635f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,lenglob,1,1,NULL,&appctx.da)); 1645f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(appctx.da)); 1655f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(appctx.da)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* 168c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the 169c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that 170c4762a1bSJed Brown have the same types. 171c4762a1bSJed Brown */ 172c4762a1bSJed Brown 1735f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(appctx.da,&u)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.dat.ic)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.dat.true_solution)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.dat.reference)); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.SEMop.grid)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.SEMop.mass)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.dat.curr_sol)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.dat.joe)); 181c4762a1bSJed Brown 1825f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(appctx.da,&xs,NULL,NULL,&xm,NULL,NULL)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 187c4762a1bSJed Brown 188c4762a1bSJed Brown xs=xs/(appctx.param.N-1); 189c4762a1bSJed Brown xm=xm/(appctx.param.N-1); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* 192c4762a1bSJed Brown Build total grid and mass over entire mesh (multi-elemental) 193c4762a1bSJed Brown */ 194c4762a1bSJed Brown 195c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 196c4762a1bSJed Brown for (j=0; j<appctx.param.N-1; j++) { 197c4762a1bSJed Brown x = (appctx.param.Le/2.0)*(appctx.SEMop.gll.nodes[j]+1.0)+appctx.param.Le*i; 198c4762a1bSJed Brown ind=i*(appctx.param.N-1)+j; 199c4762a1bSJed Brown wrk_ptr1[ind]=x; 200c4762a1bSJed Brown wrk_ptr2[ind]=.5*appctx.param.Le*appctx.SEMop.gll.weights[j]; 201c4762a1bSJed Brown if (j==0) wrk_ptr2[ind]+=.5*appctx.param.Le*appctx.SEMop.gll.weights[j]; 202c4762a1bSJed Brown } 203c4762a1bSJed Brown } 2045f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2)); 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 208c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 209c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2105f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE)); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(appctx.da,&appctx.SEMop.stiff)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(appctx.da,&appctx.SEMop.advec)); 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* 215c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 216c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 217c4762a1bSJed Brown as a time-dependent matrix. 218c4762a1bSJed Brown */ 2195f80ce2aSJacob Faibussowitsch CHKERRQ(RHSLaplacian(appctx.ts,0.0,u,appctx.SEMop.stiff,appctx.SEMop.stiff,&appctx)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(RHSAdvection(appctx.ts,0.0,u,appctx.SEMop.advec,appctx.SEMop.advec,&appctx)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(appctx.SEMop.stiff,-1.0,appctx.SEMop.advec,DIFFERENT_NONZERO_PATTERN)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(appctx.SEMop.stiff,MAT_COPY_VALUES,&appctx.SEMop.keptstiff)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */ 2255f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nsp)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetNullSpace(appctx.SEMop.stiff,nsp)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceTest(nsp,appctx.SEMop.stiff,NULL)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceDestroy(&nsp)); 229c4762a1bSJed Brown 230c4762a1bSJed Brown /* Create the TS solver that solves the ODE and its adjoint; set its options */ 2315f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&appctx.ts)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolutionFunction(appctx.ts,(PetscErrorCode (*)(TS,PetscReal,Vec, void *))ComputeReference,&appctx)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(appctx.ts,TS_LINEAR)); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(appctx.ts,TSRK)); 2355f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(appctx.ts,appctx.da)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTime(appctx.ts,0.0)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(appctx.ts,appctx.initial_dt)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(appctx.ts,appctx.param.steps)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(appctx.ts,appctx.param.Tend)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(appctx.ts,TS_EXACTFINALTIME_MATCHSTEP)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTolerances(appctx.ts,1e-7,NULL,1e-7,NULL)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(appctx.ts)); 243c4762a1bSJed Brown /* Need to save initial timestep user may have set with -ts_dt so it can be reset for each new TSSolve() */ 2445f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTimeStep(appctx.ts,&appctx.initial_dt)); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(appctx.ts,NULL,TSComputeRHSFunctionLinear,&appctx)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,TSComputeRHSJacobianConstant,&appctx)); 2475f80ce2aSJacob Faibussowitsch /* CHKERRQ(TSSetRHSFunction(appctx.ts,NULL,RHSFunction,&appctx)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,RHSJacobian,&appctx)); */ 249c4762a1bSJed Brown 250c4762a1bSJed Brown /* Set random initial conditions as initial guess, compute analytic reference solution and analytic (true) initial conditions */ 2515f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeSolutionCoefficients(&appctx)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(InitialConditions(appctx.dat.ic,&appctx)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeReference(appctx.ts,appctx.param.Tend,appctx.dat.reference,&appctx)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeReference(appctx.ts,0.0,appctx.dat.true_solution,&appctx)); 255c4762a1bSJed Brown 256f32d6360SSatish Balay /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */ 2575f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSaveTrajectory(appctx.ts)); 2585f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(appctx.ts)); 259f32d6360SSatish Balay 260c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 2615f80ce2aSJacob Faibussowitsch CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetMonitor(tao,MonitorError,&appctx,MonitorDestroy)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(tao,TAOBQNLS)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetSolution(tao,appctx.dat.ic)); 265c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 2665f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&appctx)); 267c4762a1bSJed Brown /* Check for any TAO command line options */ 2685f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetTolerances(tao,1e-8,PETSC_DEFAULT,PETSC_DEFAULT)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetFromOptions(tao)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSolve(tao)); 271c4762a1bSJed Brown 2725f80ce2aSJacob Faibussowitsch CHKERRQ(TaoDestroy(&tao)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(appctx.solutioncoefficients)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&appctx.SEMop.advec)); 2755f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&appctx.SEMop.stiff)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&appctx.SEMop.keptstiff)); 2775f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 2785f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.dat.ic)); 2795f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.dat.joe)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.dat.true_solution)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.dat.reference)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.SEMop.grid)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.SEMop.mass)); 2845f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.dat.curr_sol)); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&appctx.da)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&appctx.ts)); 288c4762a1bSJed Brown 289c4762a1bSJed Brown /* 290c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 291c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 292c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 293c4762a1bSJed Brown options are chosen (e.g., -log_summary). 294c4762a1bSJed Brown */ 295*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 296*b122ec5aSJacob Faibussowitsch return 0; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299c4762a1bSJed Brown /* 300c4762a1bSJed Brown Computes the coefficients for the analytic solution to the PDE 301c4762a1bSJed Brown */ 302c4762a1bSJed Brown PetscErrorCode ComputeSolutionCoefficients(AppCtx *appctx) 303c4762a1bSJed Brown { 304c4762a1bSJed Brown PetscRandom rand; 305c4762a1bSJed Brown PetscInt i; 306c4762a1bSJed Brown 307c4762a1bSJed Brown PetscFunctionBegin; 3085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(appctx->ncoeff,&appctx->solutioncoefficients)); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 3105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand,.9,1.0)); 311c4762a1bSJed Brown for (i=0; i<appctx->ncoeff; i++) { 3125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&appctx->solutioncoefficients[i])); 313c4762a1bSJed Brown } 3145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 315c4762a1bSJed Brown PetscFunctionReturn(0); 316c4762a1bSJed Brown } 317c4762a1bSJed Brown 318c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 319c4762a1bSJed Brown /* 320c4762a1bSJed Brown InitialConditions - Computes the (random) initial conditions for the Tao optimization solve (these are also initial conditions for the first TSSolve() 321c4762a1bSJed Brown 322c4762a1bSJed Brown Input Parameter: 323c4762a1bSJed Brown u - uninitialized solution vector (global) 324c4762a1bSJed Brown appctx - user-defined application context 325c4762a1bSJed Brown 326c4762a1bSJed Brown Output Parameter: 327c4762a1bSJed Brown u - vector with solution at initial time (global) 328c4762a1bSJed Brown */ 329c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 330c4762a1bSJed Brown { 331c4762a1bSJed Brown PetscScalar *s; 332c4762a1bSJed Brown const PetscScalar *xg; 333c4762a1bSJed Brown PetscInt i,j,lenglob; 334c4762a1bSJed Brown PetscReal sum,val; 335c4762a1bSJed Brown PetscRandom rand; 336c4762a1bSJed Brown 337c4762a1bSJed Brown PetscFunctionBegin; 3385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 3395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand,.9,1.0)); 3405f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(appctx->da,u,&s)); 3415f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 342c4762a1bSJed Brown lenglob = appctx->param.E*(appctx->param.N-1); 343c4762a1bSJed Brown for (i=0; i<lenglob; i++) { 344c4762a1bSJed Brown s[i]= 0; 345c4762a1bSJed Brown for (j=0; j<appctx->ncoeff; j++) { 3465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&val)); 347c4762a1bSJed Brown s[i] += val*PetscSinScalar(2*(j+1)*PETSC_PI*xg[i]); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown } 3505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 3515f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(appctx->da,u,&s)); 3525f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 353c4762a1bSJed Brown /* make sure initial conditions do not contain the constant functions, since with periodic boundary conditions the constant functions introduce a null space */ 3545f80ce2aSJacob Faibussowitsch CHKERRQ(VecSum(u,&sum)); 3555f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(u,-sum/lenglob)); 356c4762a1bSJed Brown PetscFunctionReturn(0); 357c4762a1bSJed Brown } 358c4762a1bSJed Brown 359c4762a1bSJed Brown /* 360c4762a1bSJed Brown TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function. 361c4762a1bSJed Brown 362a5b23f4aSJose E. Roman InitialConditions() computes the initial conditions for the beginning of the Tao iterations 363c4762a1bSJed Brown 364c4762a1bSJed Brown Input Parameter: 365c4762a1bSJed Brown u - uninitialized solution vector (global) 366c4762a1bSJed Brown appctx - user-defined application context 367c4762a1bSJed Brown 368c4762a1bSJed Brown Output Parameter: 369c4762a1bSJed Brown u - vector with solution at initial time (global) 370c4762a1bSJed Brown */ 371c4762a1bSJed Brown PetscErrorCode TrueSolution(Vec u,AppCtx *appctx) 372c4762a1bSJed Brown { 373c4762a1bSJed Brown PetscScalar *s; 374c4762a1bSJed Brown const PetscScalar *xg; 375c4762a1bSJed Brown PetscInt i,j,lenglob; 376c4762a1bSJed Brown PetscReal sum; 377c4762a1bSJed Brown 378c4762a1bSJed Brown PetscFunctionBegin; 3795f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(appctx->da,u,&s)); 3805f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 381c4762a1bSJed Brown lenglob = appctx->param.E*(appctx->param.N-1); 382c4762a1bSJed Brown for (i=0; i<lenglob; i++) { 383c4762a1bSJed Brown s[i]= 0; 384c4762a1bSJed Brown for (j=0; j<appctx->ncoeff; j++) { 385c4762a1bSJed Brown s[i] += appctx->solutioncoefficients[j]*PetscSinScalar(2*(j+1)*PETSC_PI*xg[i]); 386c4762a1bSJed Brown } 387c4762a1bSJed Brown } 3885f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(appctx->da,u,&s)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 390c4762a1bSJed Brown /* make sure initial conditions do not contain the constant functions, since with periodic boundary conditions the constant functions introduce a null space */ 3915f80ce2aSJacob Faibussowitsch CHKERRQ(VecSum(u,&sum)); 3925f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(u,-sum/lenglob)); 393c4762a1bSJed Brown PetscFunctionReturn(0); 394c4762a1bSJed Brown } 395c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 396c4762a1bSJed Brown /* 397c4762a1bSJed Brown Sets the desired profile for the final end time 398c4762a1bSJed Brown 399c4762a1bSJed Brown Input Parameters: 400c4762a1bSJed Brown t - final time 401c4762a1bSJed Brown obj - vector storing the desired profile 402c4762a1bSJed Brown appctx - user-defined application context 403c4762a1bSJed Brown 404c4762a1bSJed Brown */ 405c4762a1bSJed Brown PetscErrorCode ComputeReference(TS ts,PetscReal t,Vec obj,AppCtx *appctx) 406c4762a1bSJed Brown { 407c4762a1bSJed Brown PetscScalar *s,tc; 408c4762a1bSJed Brown const PetscScalar *xg; 409c4762a1bSJed Brown PetscInt i, j,lenglob; 410c4762a1bSJed Brown 411c4762a1bSJed Brown PetscFunctionBegin; 4125f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(appctx->da,obj,&s)); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 414c4762a1bSJed Brown lenglob = appctx->param.E*(appctx->param.N-1); 415c4762a1bSJed Brown for (i=0; i<lenglob; i++) { 416c4762a1bSJed Brown s[i]= 0; 417c4762a1bSJed Brown for (j=0; j<appctx->ncoeff; j++) { 418c4762a1bSJed Brown tc = -appctx->param.mu*(j+1)*(j+1)*4.0*PETSC_PI*PETSC_PI*t; 419c4762a1bSJed Brown s[i] += appctx->solutioncoefficients[j]*PetscSinScalar(2*(j+1)*PETSC_PI*(xg[i] + appctx->param.a*t))*PetscExpReal(tc); 420c4762a1bSJed Brown } 421c4762a1bSJed Brown } 4225f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(appctx->da,obj,&s)); 4235f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 424c4762a1bSJed Brown PetscFunctionReturn(0); 425c4762a1bSJed Brown } 426c4762a1bSJed Brown 427c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx) 428c4762a1bSJed Brown { 429c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 430c4762a1bSJed Brown 431c4762a1bSJed Brown PetscFunctionBegin; 4325f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(appctx->SEMop.keptstiff,globalin,globalout)); 433c4762a1bSJed Brown PetscFunctionReturn(0); 434c4762a1bSJed Brown } 435c4762a1bSJed Brown 436c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec globalin,Mat A, Mat B,void *ctx) 437c4762a1bSJed Brown { 438c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 439c4762a1bSJed Brown 440c4762a1bSJed Brown PetscFunctionBegin; 4415f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(appctx->SEMop.keptstiff,A,DIFFERENT_NONZERO_PATTERN)); 442c4762a1bSJed Brown PetscFunctionReturn(0); 443c4762a1bSJed Brown } 444c4762a1bSJed Brown 445c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 446c4762a1bSJed Brown 447c4762a1bSJed Brown /* 448c4762a1bSJed Brown RHSLaplacian - matrix for diffusion 449c4762a1bSJed Brown 450c4762a1bSJed Brown Input Parameters: 451c4762a1bSJed Brown ts - the TS context 452c4762a1bSJed Brown t - current time (ignored) 453c4762a1bSJed Brown X - current solution (ignored) 454c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 455c4762a1bSJed Brown 456c4762a1bSJed Brown Output Parameters: 457c4762a1bSJed Brown AA - Jacobian matrix 458c4762a1bSJed Brown BB - optionally different matrix from which the preconditioner is built 459c4762a1bSJed Brown str - flag indicating matrix structure 460c4762a1bSJed Brown 461c4762a1bSJed Brown Scales by the inverse of the mass matrix (perhaps that should be pulled out) 462c4762a1bSJed Brown 463c4762a1bSJed Brown */ 464c4762a1bSJed Brown PetscErrorCode RHSLaplacian(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx) 465c4762a1bSJed Brown { 466c4762a1bSJed Brown PetscReal **temp; 467c4762a1bSJed Brown PetscReal vv; 468c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 469c4762a1bSJed Brown PetscInt i,xs,xn,l,j; 470c4762a1bSJed Brown PetscInt *rowsDM; 471c4762a1bSJed Brown 472c4762a1bSJed Brown PetscFunctionBegin; 473c4762a1bSJed Brown /* 474c4762a1bSJed Brown Creates the element stiffness matrix for the given gll 475c4762a1bSJed Brown */ 4765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 477c4762a1bSJed Brown 478c4762a1bSJed Brown /* scale by the size of the element */ 479c4762a1bSJed Brown for (i=0; i<appctx->param.N; i++) { 480c4762a1bSJed Brown vv=-appctx->param.mu*2.0/appctx->param.Le; 481c4762a1bSJed Brown for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv; 482c4762a1bSJed Brown } 483c4762a1bSJed Brown 4845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 4855f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL)); 486c4762a1bSJed Brown 4873c859ba3SBarry Smith PetscCheck(appctx->param.N-1 >= 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Polynomial order must be at least 2"); 488c4762a1bSJed Brown xs = xs/(appctx->param.N-1); 489c4762a1bSJed Brown xn = xn/(appctx->param.N-1); 490c4762a1bSJed Brown 4915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(appctx->param.N,&rowsDM)); 492c4762a1bSJed Brown /* 493c4762a1bSJed Brown loop over local elements 494c4762a1bSJed Brown */ 495c4762a1bSJed Brown for (j=xs; j<xs+xn; j++) { 496c4762a1bSJed Brown for (l=0; l<appctx->param.N; l++) { 497c4762a1bSJed Brown rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l; 498c4762a1bSJed Brown } 4995f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES)); 500c4762a1bSJed Brown } 5015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rowsDM)); 5025f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 5035f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 5045f80ce2aSJacob Faibussowitsch CHKERRQ(VecReciprocal(appctx->SEMop.mass)); 5055f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,appctx->SEMop.mass,0)); 5065f80ce2aSJacob Faibussowitsch CHKERRQ(VecReciprocal(appctx->SEMop.mass)); 507c4762a1bSJed Brown 5085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 509c4762a1bSJed Brown PetscFunctionReturn(0); 510c4762a1bSJed Brown } 511c4762a1bSJed Brown 512c4762a1bSJed Brown /* 513c4762a1bSJed Brown Almost identical to Laplacian 514c4762a1bSJed Brown 515c4762a1bSJed Brown Note that the element matrix is NOT scaled by the size of element like the Laplacian term. 516c4762a1bSJed Brown */ 517c4762a1bSJed Brown PetscErrorCode RHSAdvection(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx) 518c4762a1bSJed Brown { 519c4762a1bSJed Brown PetscReal **temp; 520c4762a1bSJed Brown PetscReal vv; 521c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 522c4762a1bSJed Brown PetscInt i,xs,xn,l,j; 523c4762a1bSJed Brown PetscInt *rowsDM; 524c4762a1bSJed Brown 525c4762a1bSJed Brown PetscFunctionBegin; 526c4762a1bSJed Brown /* 527c4762a1bSJed Brown Creates the element stiffness matrix for the given gll 528c4762a1bSJed Brown */ 5295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 530c4762a1bSJed Brown 531c4762a1bSJed Brown /* scale by the size of the element */ 532c4762a1bSJed Brown for (i=0; i<appctx->param.N; i++) { 533c4762a1bSJed Brown vv = -appctx->param.a; 534c4762a1bSJed Brown for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv; 535c4762a1bSJed Brown } 536c4762a1bSJed Brown 5375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 5385f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL)); 539c4762a1bSJed Brown 5403c859ba3SBarry Smith PetscCheck(appctx->param.N-1 >= 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Polynomial order must be at least 2"); 541c4762a1bSJed Brown xs = xs/(appctx->param.N-1); 542c4762a1bSJed Brown xn = xn/(appctx->param.N-1); 543c4762a1bSJed Brown 5445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(appctx->param.N,&rowsDM)); 545c4762a1bSJed Brown /* 546c4762a1bSJed Brown loop over local elements 547c4762a1bSJed Brown */ 548c4762a1bSJed Brown for (j=xs; j<xs+xn; j++) { 549c4762a1bSJed Brown for (l=0; l<appctx->param.N; l++) { 550c4762a1bSJed Brown rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l; 551c4762a1bSJed Brown } 5525f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES)); 553c4762a1bSJed Brown } 5545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rowsDM)); 5555f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 5565f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 5575f80ce2aSJacob Faibussowitsch CHKERRQ(VecReciprocal(appctx->SEMop.mass)); 5585f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,appctx->SEMop.mass,0)); 5595f80ce2aSJacob Faibussowitsch CHKERRQ(VecReciprocal(appctx->SEMop.mass)); 560c4762a1bSJed Brown 5615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 562c4762a1bSJed Brown PetscFunctionReturn(0); 563c4762a1bSJed Brown } 564c4762a1bSJed Brown 565c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 566c4762a1bSJed Brown /* 567c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 568c4762a1bSJed Brown 569c4762a1bSJed Brown Input Parameters: 570c4762a1bSJed Brown tao - the Tao context 571c4762a1bSJed Brown ic - the input vector 572a82e8c82SStefano Zampini ctx - optional user-defined context, as set when calling TaoSetObjectiveAndGradient() 573c4762a1bSJed Brown 574c4762a1bSJed Brown Output Parameters: 575c4762a1bSJed Brown f - the newly evaluated function 576c4762a1bSJed Brown G - the newly evaluated gradient 577c4762a1bSJed Brown 578c4762a1bSJed Brown Notes: 579c4762a1bSJed Brown 580c4762a1bSJed Brown The forward equation is 581c4762a1bSJed Brown M u_t = F(U) 582c4762a1bSJed Brown which is converted to 583c4762a1bSJed Brown u_t = M^{-1} F(u) 584c4762a1bSJed Brown in the user code since TS has no direct way of providing a mass matrix. The Jacobian of this is 585c4762a1bSJed Brown M^{-1} J 586c4762a1bSJed Brown where J is the Jacobian of F. Now the adjoint equation is 587c4762a1bSJed Brown M v_t = J^T v 588c4762a1bSJed Brown but TSAdjoint does not solve this since it can only solve the transposed system for the 589c4762a1bSJed Brown Jacobian the user provided. Hence TSAdjoint solves 590c4762a1bSJed Brown w_t = J^T M^{-1} w (where w = M v) 591a5b23f4aSJose E. Roman since there is no way to indicate the mass matrix as a separate entity to TS. Thus one 592c4762a1bSJed Brown must be careful in initializing the "adjoint equation" and using the result. This is 593c4762a1bSJed Brown why 594c4762a1bSJed Brown G = -2 M(u(T) - u_d) 595c4762a1bSJed Brown below (instead of -2(u(T) - u_d) 596c4762a1bSJed Brown 597c4762a1bSJed Brown */ 598c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec ic,PetscReal *f,Vec G,void *ctx) 599c4762a1bSJed Brown { 600c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 601c4762a1bSJed Brown Vec temp; 602c4762a1bSJed Brown 603c4762a1bSJed Brown PetscFunctionBegin; 6045f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTime(appctx->ts,0.0)); 6055f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetStepNumber(appctx->ts,0)); 6065f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(appctx->ts,appctx->initial_dt)); 6075f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(ic,appctx->dat.curr_sol)); 608c4762a1bSJed Brown 6095f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(appctx->ts,appctx->dat.curr_sol)); 6105f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(appctx->dat.curr_sol,appctx->dat.joe)); 611c4762a1bSJed Brown 612c4762a1bSJed Brown /* Compute the difference between the current ODE solution and target ODE solution */ 6135f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(G,-1.0,appctx->dat.curr_sol,appctx->dat.reference)); 614c4762a1bSJed Brown 615c4762a1bSJed Brown /* Compute the objective/cost function */ 6165f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(G,&temp)); 6175f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(temp,G,G)); 6185f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(temp,appctx->SEMop.mass,f)); 6195f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&temp)); 620c4762a1bSJed Brown 621c4762a1bSJed Brown /* Compute initial conditions for the adjoint integration. See Notes above */ 6225f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(G, -2.0)); 6235f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(G,G,appctx->SEMop.mass)); 6245f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetCostGradients(appctx->ts,1,&G,NULL)); 625c4762a1bSJed Brown 6265f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdjointSolve(appctx->ts)); 6275f80ce2aSJacob Faibussowitsch /* CHKERRQ(VecPointwiseDivide(G,G,appctx->SEMop.mass));*/ 628c4762a1bSJed Brown PetscFunctionReturn(0); 629c4762a1bSJed Brown } 630c4762a1bSJed Brown 631c4762a1bSJed Brown PetscErrorCode MonitorError(Tao tao,void *ctx) 632c4762a1bSJed Brown { 633c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 634c4762a1bSJed Brown Vec temp,grad; 635c4762a1bSJed Brown PetscReal nrm; 636c4762a1bSJed Brown PetscInt its; 637c4762a1bSJed Brown PetscReal fct,gnorm; 638c4762a1bSJed Brown 639c4762a1bSJed Brown PetscFunctionBegin; 6405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(appctx->dat.ic,&temp)); 6415f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(temp,-1.0,appctx->dat.ic,appctx->dat.true_solution)); 6425f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(temp,temp,temp)); 6435f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(temp,appctx->SEMop.mass,&nrm)); 644c4762a1bSJed Brown nrm = PetscSqrtReal(nrm); 6455f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetGradient(tao,&grad,NULL,NULL)); 6465f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(temp,temp,temp)); 6475f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(temp,appctx->SEMop.mass,&gnorm)); 648c4762a1bSJed Brown gnorm = PetscSqrtReal(gnorm); 6495f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&temp)); 6505f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetIterationNumber(tao,&its)); 6515f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetSolutionStatus(tao,NULL,&fct,NULL,NULL,NULL,NULL)); 652c4762a1bSJed Brown if (!its) { 6535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%% Iteration Error Objective Gradient-norm\n")); 6545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"history = [\n")); 655c4762a1bSJed Brown } 6565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%3D %g %g %g\n",its,(double)nrm,(double)fct,(double)gnorm)); 657c4762a1bSJed Brown PetscFunctionReturn(0); 658c4762a1bSJed Brown } 659c4762a1bSJed Brown 660c4762a1bSJed Brown PetscErrorCode MonitorDestroy(void **ctx) 661c4762a1bSJed Brown { 662c4762a1bSJed Brown PetscFunctionBegin; 6635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"];\n")); 664c4762a1bSJed Brown PetscFunctionReturn(0); 665c4762a1bSJed Brown } 666c4762a1bSJed Brown 667c4762a1bSJed Brown /*TEST 668c4762a1bSJed Brown 669c4762a1bSJed Brown build: 670c4762a1bSJed Brown requires: !complex 671c4762a1bSJed Brown 672c4762a1bSJed Brown test: 673c4762a1bSJed Brown requires: !single 674c4762a1bSJed Brown args: -ts_adapt_dt_max 3.e-3 -E 10 -N 8 -ncoeff 5 -tao_bqnls_mat_lmvm_scale_type none 675c4762a1bSJed Brown 676c4762a1bSJed Brown test: 677c4762a1bSJed Brown suffix: cn 678c4762a1bSJed Brown requires: !single 679c4762a1bSJed Brown args: -ts_type cn -ts_dt .003 -pc_type lu -E 10 -N 8 -ncoeff 5 -tao_bqnls_mat_lmvm_scale_type none 680c4762a1bSJed Brown 681c4762a1bSJed Brown test: 682c4762a1bSJed Brown suffix: 2 683c4762a1bSJed Brown requires: !single 684c4762a1bSJed Brown args: -ts_adapt_dt_max 3.e-3 -E 10 -N 8 -ncoeff 5 -a .1 -tao_bqnls_mat_lmvm_scale_type none 685c4762a1bSJed Brown 686c4762a1bSJed Brown TEST*/ 687