/**
     * Simple Parabolic PDE u' = del squared u
     *
     * With u = 0 on the boundaries of the unit cube. Subject to the initial
     * condition u(0,x,y,z)=sin( PI x)sin( PI y)sin( PI z).
     */
    void TestSimpleLinearParabolicSolver3DZeroDirich()
    {
        // read mesh on [0,1]x[0,1]x[0,1]
        TrianglesMeshReader<3,3> mesh_reader("mesh/test/data/cube_136_elements");
        TetrahedralMesh<3,3> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        // Instantiate PDE object
        HeatEquation<3> pde;

        // Boundary conditions - zero dirichlet everywhere on boundary
        BoundaryConditionsContainer<3,3,1> bcc;
        bcc.DefineZeroDirichletOnMeshBoundary(&mesh);

        // Solver
        SimpleLinearParabolicSolver<3,3> solver(&mesh,&pde,&bcc);

        /*
         * Choose initial condition sin(x*pi)*sin(y*pi)*sin(z*pi) as
         * this is an eigenfunction of the heat equation.
         */
        std::vector<double> init_cond(mesh.GetNumNodes());
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];
            double y = mesh.GetNode(i)->GetPoint()[1];
            double z = mesh.GetNode(i)->GetPoint()[2];
            init_cond[i] = sin(x*M_PI)*sin(y*M_PI)*sin(z*M_PI);
        }
        Vec initial_condition = PetscTools::CreateVec(init_cond);

        double t_end = 0.1;
        solver.SetTimes(0, t_end);
        solver.SetTimeStep(0.001);

        solver.SetInitialCondition(initial_condition);

        Vec result = solver.Solve();
        ReplicatableVector result_repl(result);

        // Check solution is u = e^{-3*t*pi*pi} sin(x*pi)*sin(y*pi)*sin(z*pi), t=0.1
        for (unsigned i=0; i<result_repl.GetSize(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];
            double y = mesh.GetNode(i)->GetPoint()[1];
            double z = mesh.GetNode(i)->GetPoint()[2];
            double u = exp(-3*t_end*M_PI*M_PI)*sin(x*M_PI)*sin(y*M_PI)*sin(z*M_PI);
            TS_ASSERT_DELTA(result_repl[i], u, 0.1);
        }

        PetscTools::Destroy(initial_condition);
        PetscTools::Destroy(result);
    }
    void TestValidate()
    {
        // Load a 2D square mesh with 1 central non-boundary node
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/square_4_elements");
        TetrahedralMesh<2,2> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        BoundaryConditionsContainer<2,2,1> bcc;

        // No BCs yet, so shouldn't validate
        TS_ASSERT(!bcc.Validate(&mesh));

        // Add some BCs
        ConstBoundaryCondition<2> *bc = new ConstBoundaryCondition<2>(0.0);
        bcc.AddDirichletBoundaryCondition(mesh.GetNode(0), bc);
        bcc.AddDirichletBoundaryCondition(mesh.GetNode(1), bc);
        bcc.AddDirichletBoundaryCondition(mesh.GetNode(3), bc);
        TetrahedralMesh<2,2>::BoundaryElementIterator iter
        = mesh.GetBoundaryElementIteratorEnd();
        iter--;
        bcc.AddNeumannBoundaryCondition(*iter, bc); // 2 to 3
        iter--;
        bcc.AddNeumannBoundaryCondition(*iter, bc); // 1 to 2

        TS_ASSERT(bcc.Validate(&mesh));
    }
    void TestDefineZeroDirichletOnMeshBoundary()
    {
        // Load a 2D square mesh with 1 central non-boundary node
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/square_4_elements");
        TetrahedralMesh<2,2> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        BoundaryConditionsContainer<2,2,1> bcc;

        bcc.DefineZeroDirichletOnMeshBoundary(&mesh);

        // Check boundary nodes have the right condition
        for (int i=0; i<4; i++)
        {
            double value = bcc.GetDirichletBCValue(mesh.GetNode(i));
            TS_ASSERT_DELTA(value, 0.0, 1e-12);
        }

        // Check non-boundary node has no condition
        TS_ASSERT(!bcc.HasDirichletBoundaryCondition(mesh.GetNode(4)));
    }
    void TestAnyNonZeroNeumannConditionsAndApplyNeumannToMeshBoundary()
    {
        // Load a 2D square mesh with 1 central non-boundary node
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/square_4_elements");
        TetrahedralMesh<2,2> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        BoundaryConditionsContainer<2,2,1> bcc;
        BoundaryConditionsContainer<2,2,2> bcc_2unknowns;

        TS_ASSERT_EQUALS(bcc.AnyNonZeroNeumannConditions(), false);

        bcc.DefineZeroNeumannOnMeshBoundary(&mesh);

        TetrahedralMesh<2,2>::BoundaryElementIterator iter;
        iter = mesh.GetBoundaryElementIteratorBegin();
        while (iter != mesh.GetBoundaryElementIteratorEnd())
        {
            TS_ASSERT(bcc.HasNeumannBoundaryCondition(*iter));
            double value = bcc.GetNeumannBCValue(*iter, (*iter)->GetNode(0)->GetPoint());
            TS_ASSERT_DELTA(value, 0.0, 1e-8);

            iter++;
        }
        TS_ASSERT_EQUALS(bcc.AnyNonZeroNeumannConditions(), false);

        iter = mesh.GetBoundaryElementIteratorBegin();

        ConstBoundaryCondition<2>* p_boundary_condition2 = new ConstBoundaryCondition<2>(-1);

        bcc_2unknowns.AddNeumannBoundaryCondition(*iter, p_boundary_condition2);
        TS_ASSERT_EQUALS(bcc_2unknowns.AnyNonZeroNeumannConditions(), true);
    }
    void TestAddNeumannBoundaryConditions()
    {
          // Load a 2D square mesh with 1 central non-boundary node
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/square_4_elements");
        TetrahedralMesh<2,2> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        BoundaryConditionsContainer<2,2,2> bcc;

        // No BCs yet, so shouldn't validate
        TS_ASSERT(!bcc.Validate(&mesh));

        // Add some BCs
        ConstBoundaryCondition<2> *bc1 = new ConstBoundaryCondition<2>(2.0);
        ConstBoundaryCondition<2> *bc2 = new ConstBoundaryCondition<2>(-3.0);

        TetrahedralMesh<2,2>::BoundaryElementIterator iter = mesh.GetBoundaryElementIteratorEnd();
        iter--;
        bcc.AddNeumannBoundaryCondition(*iter, bc1, 0);
        bcc.AddNeumannBoundaryCondition(*iter, bc2, 1);
        iter--;
        bcc.AddNeumannBoundaryCondition(*iter, bc1, 0);
        iter--;
        bcc.AddNeumannBoundaryCondition(*iter, bc2, 1);

        iter = mesh.GetBoundaryElementIteratorEnd();
        iter--;
        TS_ASSERT_DELTA(bcc.GetNeumannBCValue(*iter, ChastePoint<2>(), 0), 2.0, 1e-9);
        TS_ASSERT_DELTA(bcc.GetNeumannBCValue(*iter, ChastePoint<2>(), 1), -3.0, 1e-9);

        iter--;
        TS_ASSERT_DELTA(bcc.GetNeumannBCValue(*iter, ChastePoint<2>(), 0), 2.0, 1e-9);
        TS_ASSERT_DELTA(bcc.GetNeumannBCValue(*iter, ChastePoint<2>(), 1), 0.0, 1e-9);

        iter--;
        TS_ASSERT_DELTA(bcc.GetNeumannBCValue(*iter, ChastePoint<2>(), 0), 0.0, 1e-9);
        TS_ASSERT_DELTA(bcc.GetNeumannBCValue(*iter, ChastePoint<2>(), 1), -3.0, 1e-9);
    }
    /*
     * Define a particular test.
     */
    void TestSchnackenbergSystemOnButterflyMesh() throw (Exception)
    {
        /* As usual, we first create a mesh. Here we are using a 2d mesh of a butterfly-shaped domain. */
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/butterfly");
        TetrahedralMesh<2,2> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        /* We scale the mesh to an appropriate size. */
        mesh.Scale(0.2, 0.2);

        /* Next, we instantiate the PDE system to be solved. We pass the parameter values into the
         * constructor.  (The order is D,,1,,  D,,2,,  k,,1,,  k,,-1,,  k,,2,,  k,,3,,) */
        SchnackenbergCoupledPdeSystem<2> pde(1e-4, 1e-2, 0.1, 0.2, 0.3, 0.1);

        /*
         * Then we have to define the boundary conditions. As we are in 2d, {{{SPACE_DIM}}}=2 and
         * {{{ELEMENT_DIM}}}=2. We also have two unknowns u and v,
         * so in this case {{{PROBLEM_DIM}}}=2. The value of each boundary condition is
         * given by the spatially uniform steady state solution of the Schnackenberg system,
         * given by u = (k,,1,, + k,,2,,)/k,,-1,,, v = k,,2,,k,,-1,,^2^/k,,3,,(k,,1,, + k,,2,,)^2^.
         */
        BoundaryConditionsContainer<2,2,2> bcc;
        ConstBoundaryCondition<2>* p_bc_for_u = new ConstBoundaryCondition<2>(2.0);
        ConstBoundaryCondition<2>* p_bc_for_v = new ConstBoundaryCondition<2>(0.75);
        for (TetrahedralMesh<2,2>::BoundaryNodeIterator node_iter = mesh.GetBoundaryNodeIteratorBegin();
             node_iter != mesh.GetBoundaryNodeIteratorEnd();
             ++node_iter)
        {
            bcc.AddDirichletBoundaryCondition(*node_iter, p_bc_for_u, 0);
            bcc.AddDirichletBoundaryCondition(*node_iter, p_bc_for_v, 1);
        }

        /* This is the solver for solving coupled systems of linear parabolic PDEs and ODEs,
         * which takes in the mesh, the PDE system, the boundary conditions and optionally
         * a vector of ODE systems (one for each node in the mesh). Since in this example
         * we are solving a system of coupled PDEs only, we do not supply this last argument. */
        LinearParabolicPdeSystemWithCoupledOdeSystemSolver<2,2,2> solver(&mesh, &pde, &bcc);

        /* Then we set the end time and time step and the output directory to which results will be written. */
        double t_end = 10;
        solver.SetTimes(0, t_end);
        solver.SetTimeStep(1e-1);
        solver.SetSamplingTimeStep(1);
        solver.SetOutputDirectory("TestSchnackenbergSystemOnButterflyMesh");

        /* We create a vector of initial conditions for u and v that are random perturbations
         * of the spatially uniform steady state and pass this to the solver. */
        std::vector<double> init_conds(2*mesh.GetNumNodes());
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            init_conds[2*i] = fabs(2.0 + RandomNumberGenerator::Instance()->ranf());
            init_conds[2*i + 1] = fabs(0.75 + RandomNumberGenerator::Instance()->ranf());
        }
        Vec initial_condition = PetscTools::CreateVec(init_conds);
        solver.SetInitialCondition(initial_condition);

        /* We now solve the PDE system and write results to VTK files, for
         * visualization using Paraview.  Results will be written to CHASTE_TEST_OUTPUT/TestSchnackenbergSystemOnButterflyMesh
         * as a results.pvd file and several results_[time].vtu files.
         * You should see something like [[Image(u.png, 350px)]] for u and [[Image(v.png, 350px)]] for v.
         */
        solver.SolveAndWriteResultsToFile();

        /*
         * All PETSc {{{Vec}}}s should be destroyed when they are no longer needed.
         */
        PetscTools::Destroy(initial_condition);
    }
    /* Define a particular test. Note the {{{throw(Exception)}}} at the end of the
     * declaration. This causes {{{Exception}}} messages to be printed out if an
     * {{{Exception}}} is thrown, rather than just getting the message "terminate
     * called after throwing an instance of 'Exception' " */
    void TestSolvingNonlinearEllipticPde() throw(Exception)
    {
        /* As usual, first create a mesh. */
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/square_128_elements");
        TetrahedralMesh<2,2> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        /* Next, instantiate the PDE to be solved. */
        MyNonlinearPde pde;

        /*
         * Then we have to define the boundary conditions. First, the Dirichlet boundary
         * condition, u=0 on x=0, using the boundary node iterator.
         */
        BoundaryConditionsContainer<2,2,1> bcc;
        ConstBoundaryCondition<2>* p_zero_bc = new ConstBoundaryCondition<2>(0.0);
        for (TetrahedralMesh<2,2>::BoundaryNodeIterator node_iter = mesh.GetBoundaryNodeIteratorBegin();
             node_iter != mesh.GetBoundaryNodeIteratorEnd();
             node_iter++)
        {
            if (fabs((*node_iter)->GetPoint()[1]) < 1e-12)
            {
                bcc.AddDirichletBoundaryCondition(*node_iter, p_zero_bc);
            }
        }

        /* And then the Neumman conditions. Neumann boundary condition are defined on
         * surface elements, and for this problem, the Neumman boundary value depends
         * on the position in space, so we make use of the {{{FunctionalBoundaryCondition}}}
         * object, which contains a pointer to a function, and just returns the value
         * of that function for the required point when the {{{GetValue}}} method is called.
         */
        FunctionalBoundaryCondition<2>* p_functional_bc = new FunctionalBoundaryCondition<2>(&MyNeummanFunction);
        /* Loop over surface elements. */
        for (TetrahedralMesh<2,2>::BoundaryElementIterator elt_iter = mesh.GetBoundaryElementIteratorBegin();
             elt_iter != mesh.GetBoundaryElementIteratorEnd();
             elt_iter++)
        {
            /* Get the y value of any node (here, the zero-th). */
            double y = (*elt_iter)->GetNodeLocation(0,1);
            /* If y=1... */
            if (fabs(y-1.0) < 1e-12)
            {
                /* ... then associate the functional boundary condition, (Dgradu).n = y,
                 *  with the surface element... */
                bcc.AddNeumannBoundaryCondition(*elt_iter, p_functional_bc);
            }
            else
            {
                /* ...else associate the zero boundary condition (i.e. zero flux) with this
                 * element. */
                bcc.AddNeumannBoundaryCondition(*elt_iter, p_zero_bc);
            }
        }
        /* Note that in the above loop, the zero Neumman boundary condition was applied
         * to all surface elements for which y!=1, which included the Dirichlet surface
         * y=0. This is OK, as Dirichlet boundary conditions are applied to the finite
         * element matrix after Neumman boundary conditions, where the appropriate rows
         * in the matrix are overwritten.
         *
         * This is the solver for solving nonlinear problems, which, as usual,
         * takes in the mesh, the PDE, and the boundary conditions. */
        SimpleNonlinearEllipticSolver<2,2> solver(&mesh, &pde, &bcc);

        /* The solver also needs to be given an initial guess, which will be
         * a PETSc vector. We can make use of a helper method to create it.
         */
        Vec initial_guess = PetscTools::CreateAndSetVec(mesh.GetNumNodes(), 0.25);

        /* '''Optional:''' To use Chaste's Newton solver to solve nonlinear vector equations that are
         * assembled, rather than the default PETSc nonlinear solvers, we can
         * do the following: */
        SimpleNewtonNonlinearSolver newton_solver;
        solver.SetNonlinearSolver(&newton_solver);
        /* '''Optional:''' We can also manually set tolerances, and whether to print statistics, with
         * this nonlinear vector equation solver */
        newton_solver.SetTolerance(1e-10);
        newton_solver.SetWriteStats();

        /* Now call {{{Solve}}}, passing in the initial guess */
        Vec answer = solver.Solve(initial_guess);

        /* Note that we could have got the solver to not use an analytical Jacobian
         * and use a numerically-calculated Jacobian instead, by passing in false as a second
         * parameter:
         */
        //Vec answer = solver.Solve(initial_guess, false);

        /* Once solved, we can check the obtained solution against the analytical
         * solution. */
        ReplicatableVector answer_repl(answer);
        for (unsigned i=0; i<answer_repl.GetSize(); i++)
        {
            double y = mesh.GetNode(i)->GetPoint()[1];
            double exact_u = sqrt(y*(4-y));
            TS_ASSERT_DELTA(answer_repl[i], exact_u, 0.15);
        }

        /* Finally, we have to remember to destroy the PETSc {{{Vec}}}s. */
        PetscTools::Destroy(initial_guess);
        PetscTools::Destroy(answer);
    }
    /**
     * This test provides an example of how to solve a coupled PDE system
     * where there is no coupled ODE system, and can be used as a template
     * for solving standard reaction-diffusion problems arising in the
     * study of pattern formation on fixed domains.
     */
    void TestSchnackenbergCoupledPdeSystemIn1dWithNonZeroDirichlet()
    {
        // Create mesh of domain [0,1]
        TrianglesMeshReader<1,1> mesh_reader("mesh/test/data/1D_0_to_1_1000_elements");
        TetrahedralMesh<1,1> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        // Create PDE system object
        SchnackenbergCoupledPdeSystem<1> pde(1e-4, 1e-2, 0.1, 0.2, 0.3, 0.1);

        // Create non-zero Dirichlet boundary conditions for each state variable
        BoundaryConditionsContainer<1,1,2> bcc;
        ConstBoundaryCondition<1>* p_bc_for_u = new ConstBoundaryCondition<1>(2.0);
        ConstBoundaryCondition<1>* p_bc_for_v = new ConstBoundaryCondition<1>(0.75);
        bcc.AddDirichletBoundaryCondition(mesh.GetNode(0), p_bc_for_u, 0);
        bcc.AddDirichletBoundaryCondition(mesh.GetNode(0), p_bc_for_v, 1);
        bcc.AddDirichletBoundaryCondition(mesh.GetNode(mesh.GetNumNodes()-1), p_bc_for_u, 0);
        bcc.AddDirichletBoundaryCondition(mesh.GetNode(mesh.GetNumNodes()-1), p_bc_for_v, 1);

        // Create PDE system solver
        LinearParabolicPdeSystemWithCoupledOdeSystemSolver<1,1,2> solver(&mesh, &pde, &bcc);

        // Set end time and time step
        double t_end = 10;
        solver.SetTimes(0, t_end);
        solver.SetTimeStep(1e-1);

        // Create initial conditions that are random perturbations of the uniform steady state
        std::vector<double> init_conds(2*mesh.GetNumNodes());
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            init_conds[2*i] = fabs(2.0 + RandomNumberGenerator::Instance()->ranf());
            init_conds[2*i + 1] = fabs(0.75 + RandomNumberGenerator::Instance()->ranf());
        }
        Vec initial_condition = PetscTools::CreateVec(init_conds);
        solver.SetInitialCondition(initial_condition);

        // Solve PDE system and store result
        Vec solution = solver.Solve();
        ReplicatableVector solution_repl(solution);

        // Write results for visualization in gnuplot
        OutputFileHandler handler("TestSchnackenbergCoupledPdeSystemIn1dWithNonZeroDirichlet", false);
        out_stream results_file = handler.OpenOutputFile("schnackenberg.dat");
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            double x = mesh.GetNode(i)->rGetLocation()[0];
            double u = solution_repl[2*i];
            double v = solution_repl[2*i + 1];
            (*results_file) << x << "\t" << u << "\t" << v << "\n" << std::flush;
        }
        results_file->close();

        std::string results_filename = handler.GetOutputDirectoryFullPath() + "schnackenberg.dat";
        NumericFileComparison comp_results(results_filename, "pde/test/data/schnackenberg.dat");
        TS_ASSERT(comp_results.CompareFiles(1e-3));

        // Tidy up
        PetscTools::Destroy(initial_condition);
        PetscTools::Destroy(solution);
    }
    void TestHeatEquationWithSourceWithCoupledOdeSystemIn1dWithZeroNeumann()
    {
        // Create mesh of domain [0,1]
        TrianglesMeshReader<1,1> mesh_reader("mesh/test/data/1D_0_to_1_100_elements");
        TetrahedralMesh<1,1> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        // Create PDE system object
        HeatEquationWithSourceForCoupledOdeSystem<1> pde;

        // Define zero Neumann boundary conditions
        BoundaryConditionsContainer<1,1,1> bcc;
        ConstBoundaryCondition<1>* p_boundary_condition = new ConstBoundaryCondition<1>(0.0);
        TetrahedralMesh<1,1>::BoundaryElementIterator iter = mesh.GetBoundaryElementIteratorBegin();
        bcc.AddNeumannBoundaryCondition(*iter, p_boundary_condition);
        iter = mesh.GetBoundaryElementIteratorEnd();
        iter--;
        bcc.AddNeumannBoundaryCondition(*iter, p_boundary_condition);

        // Create the correct number of ODE systems
        double a = 5.0;
        std::vector<AbstractOdeSystemForCoupledPdeSystem*> ode_systems;
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            ode_systems.push_back(new OdeSystemForCoupledHeatEquationWithSource(a));
        }

        // Create PDE system solver
        LinearParabolicPdeSystemWithCoupledOdeSystemSolver<1,1,1> solver(&mesh, &pde, &bcc, ode_systems);

        // Test setting end time and timestep
        TS_ASSERT_THROWS_THIS(solver.SetTimes(1.0, 0.0), "Start time has to be less than end time");
        TS_ASSERT_THROWS_THIS(solver.SetTimeStep(0.0), "Time step has to be greater than zero");

        // Set end time and timestep
        double t_end = 0.1;
        solver.SetTimes(0, t_end);
        solver.SetTimeStep(0.001);

        // Set initial condition u(x,0) = 1 + cos(pi*x)
        std::vector<double> init_cond(mesh.GetNumNodes());
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];
            init_cond[i] = 1 + cos(M_PI*x);
        }
        Vec initial_condition = PetscTools::CreateVec(init_cond);
        solver.SetInitialCondition(initial_condition);

        // Solve PDE system and store result
        Vec result = solver.Solve();
        ReplicatableVector result_repl(result);

        /*
         * Test that solution is given by
         *
         * u(x,t) = 1 + (1 - exp(-a*t))/a + exp(-pi*pi*t)*cos(pi*x),
         * v(x,t) = exp(-a*t),
         *
         * with t = t_end.
         */
        for (unsigned i=0; i<result_repl.GetSize(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];

            double u = 1 + (1 - exp(-a*t_end))/a + exp(-M_PI*M_PI*t_end)*cos(M_PI*x);
            TS_ASSERT_DELTA(result_repl[i], u, 0.1);

            double u_from_v = solver.GetOdeSystemAtNode(i)->rGetPdeSolution()[0];
            TS_ASSERT_DELTA(result_repl[i], u_from_v, 0.1);

            double v = exp(-a*t_end);
            TS_ASSERT_DELTA(ode_systems[i]->rGetStateVariables()[0], v, 0.1);
        }

        // Test the method GetOdeSystemAtNode()
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            TS_ASSERT(solver.GetOdeSystemAtNode(i) != NULL);
            TS_ASSERT_DELTA(static_cast<OdeSystemForCoupledHeatEquationWithSource*>(solver.GetOdeSystemAtNode(i))->GetA(), 5.0, 1e-6);
        }

        // Tidy up
        PetscTools::Destroy(initial_condition);
        PetscTools::Destroy(result);
    }
    void TestHeatEquationWithCoupledOdeSystemIn2dWithZeroDirichlet()
    {
        // Create mesh of the domain [0,1]x[0,1]
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/square_4096_elements");
        TetrahedralMesh<2,2> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        // Create PDE system object
        HeatEquationForCoupledOdeSystem<2> pde;

        // Define zero Dirichlet boundary conditions on entire boundary
        BoundaryConditionsContainer<2,2,1> bcc;
        bcc.DefineZeroDirichletOnMeshBoundary(&mesh);

        // Create the correct number of ODE systems
        double a = 5.0;
        std::vector<AbstractOdeSystemForCoupledPdeSystem*> ode_systems;
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            ode_systems.push_back(new OdeSystemForCoupledHeatEquation(a));
        }

        // Create PDE system solver
        LinearParabolicPdeSystemWithCoupledOdeSystemSolver<2,2,1> solver(&mesh, &pde, &bcc, ode_systems);

        // Set end time and timestep
        double t_end = 0.01;
        solver.SetTimes(0, t_end);
        solver.SetTimeStep(0.001);

        /*
         * Set initial condition
         *
         * u(x,y,0) = sin(pi*x)*sin(pi*y),
         *
         * which is an eigenfunction of the heat equation.
         */
        std::vector<double> init_cond(mesh.GetNumNodes());
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];
            double y = mesh.GetNode(i)->GetPoint()[1];
            init_cond[i] = sin(M_PI*x)*sin(M_PI*y);
        }
        Vec initial_condition = PetscTools::CreateVec(init_cond);
        solver.SetInitialCondition(initial_condition);

        // Solve PDE system and store result
        Vec result = solver.Solve();
        ReplicatableVector result_repl(result);

        /*
         * Test that solution is given by
         *
         * u(x,y,t) = e^{-2*pi*pi*t} sin(pi*x)*sin(pi*y),
         * v(x,y,t) = 1 + (1 - e^{-2*pi*pi*t})*sin(pi*x)*sin(pi*y)*a/(2*pi*pi),
         *
         * with t = t_end.
         */
        for (unsigned i=0; i<result_repl.GetSize(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];
            double y = mesh.GetNode(i)->GetPoint()[1];

            double u = exp(-2*M_PI*M_PI*t_end)*sin(M_PI*x)*sin(M_PI*y);
            double v = 1.0 + (a/(2*M_PI*M_PI))*(1 - exp(-2*M_PI*M_PI*t_end))*sin(M_PI*x)*sin(M_PI*y);

            TS_ASSERT_DELTA(result_repl[i], u, 0.01);
            TS_ASSERT_DELTA(ode_systems[i]->rGetStateVariables()[0], v, 0.01);
        }

        // Tidy up
        PetscTools::Destroy(initial_condition);
        PetscTools::Destroy(result);
    }
    void TestSolveAndWriteResultsToFileMethod()
    {
        // Create mesh of the domain [0,1]x[0,1]
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/square_128_elements");
        TetrahedralMesh<2,2> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        // Create PDE system object
        HeatEquationForCoupledOdeSystem<2> pde;

        // Define zero Dirichlet boundary conditions on entire boundary
        BoundaryConditionsContainer<2,2,1> bcc;
        bcc.DefineZeroDirichletOnMeshBoundary(&mesh);

        // Create the correct number of ODE systems
        double a = 5.0;
        std::vector<AbstractOdeSystemForCoupledPdeSystem*> ode_systems;
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            ode_systems.push_back(new OdeSystemForCoupledHeatEquation(a));
        }

        // Create PDE system solver
        LinearParabolicPdeSystemWithCoupledOdeSystemSolver<2,2,1> solver(&mesh, &pde, &bcc, ode_systems);

        // Set end time and timestep (end time is not a multiple of timestep, for coverage)
        double t_end = 0.105;

        /*
         * Set initial condition
         *
         * u(x,y,0) = sin(pi*x)*sin(pi*y),
         *
         * which is an eigenfunction of the heat equation.
         */
        std::vector<double> init_cond(mesh.GetNumNodes());
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];
            double y = mesh.GetNode(i)->GetPoint()[1];
            init_cond[i] = sin(M_PI*x)*sin(M_PI*y);
        }
        Vec initial_condition = PetscTools::CreateVec(init_cond);

        // Need an output folder
        TS_ASSERT_THROWS_THIS(solver.SolveAndWriteResultsToFile(),
                "SetOutputDirectory() must be called prior to SolveAndWriteResultsToFile()");
        solver.SetOutputDirectory("TestHeatEquationForCoupledOdeSystemIn2dWithZeroDirichletWithOutput");

        // Need a time interval
        TS_ASSERT_THROWS_THIS(solver.SolveAndWriteResultsToFile(),
                "SetTimes() must be called prior to SolveAndWriteResultsToFile()");
        solver.SetTimes(0, t_end);

        // Need a timestep
        TS_ASSERT_THROWS_THIS(solver.SolveAndWriteResultsToFile(),
                "SetTimeStep() must be called prior to SolveAndWriteResultsToFile()");
        solver.SetTimeStep(0.01);

        // Need sampling interval
        TS_ASSERT_THROWS_THIS(solver.SolveAndWriteResultsToFile(),
                "SetSamplingTimeStep() must be called prior to SolveAndWriteResultsToFile()");
        solver.SetSamplingTimeStep(0.1);

        // Need initial condition
        TS_ASSERT_THROWS_THIS(solver.SolveAndWriteResultsToFile(),
                "SetInitialCondition() must be called prior to SolveAndWriteResultsToFile()");
        solver.SetInitialCondition(initial_condition);

        solver.SolveAndWriteResultsToFile();
//#ifdef CHASTE_VTK
///\todo #1967 Check that the file was output and has expected content
//#endif // CHASTE_VTK

        // Tidy up
        PetscTools::Destroy(initial_condition);
    }
    void TestPerdiodicBoundaryConditions()
    {
        const int SIZE = 5;

        DistributedVectorFactory factory(SIZE);
        Vec template_vec = factory.CreateVec(2);
        LinearSystem linear_system(template_vec, 2*SIZE);
        PetscTools::Destroy(template_vec);

        for (int i = 0; i < 2*SIZE; i++)
        {
            for (int j = 0; j < 2*SIZE; j++)
            {
                // LHS matrix is all 2s
                linear_system.SetMatrixElement(i,j,2);
            }
            // RHS vector is all 3s
            linear_system.SetRhsVectorElement(i,3);
        }

        linear_system.AssembleIntermediateLinearSystem();

        Node<3>* nodes[SIZE];
        BoundaryConditionsContainer<3,3,2> bcc;

        for (unsigned i=0; i<SIZE-1; i++)
        {
            nodes[i] = new Node<3>(i,true);
        }

        bcc.AddPeriodicBoundaryCondition(nodes[0], nodes[1]);
        bcc.AddPeriodicBoundaryCondition(nodes[2], nodes[3]);

        bcc.ApplyPeriodicBcsToLinearProblem(linear_system, true, true);

        linear_system.AssembleFinalLinearSystem();

        ReplicatableVector rhs_repl(linear_system.rGetRhsVector());
        TS_ASSERT_DELTA(rhs_repl[0], 0.0, 1e-12); // node 0, variable 0
        TS_ASSERT_DELTA(rhs_repl[1], 0.0, 1e-12); // node 0, variable 1
        TS_ASSERT_DELTA(rhs_repl[2], 3.0, 1e-12);
        TS_ASSERT_DELTA(rhs_repl[3], 3.0, 1e-12);
        TS_ASSERT_DELTA(rhs_repl[4], 0.0, 1e-12); // node 2, variable 0
        TS_ASSERT_DELTA(rhs_repl[5], 0.0, 1e-12); // node 2, variable 1
        TS_ASSERT_DELTA(rhs_repl[6], 3.0, 1e-12);
        TS_ASSERT_DELTA(rhs_repl[7], 3.0, 1e-12);
        TS_ASSERT_DELTA(rhs_repl[8], 3.0, 1e-12);
        TS_ASSERT_DELTA(rhs_repl[9], 3.0, 1e-12);


        //
        //  Matrix should have
        //   row 0 altered to be [1, 0 -1, 0, 0, ..., 0]
        //   row 1 altered to be [0, 1, 0 -1, 0, ..., 0]
        //   row 4 altered to be [0, 0, 0, 0, 1, 0, -1, 0, 0, 0]
        //   row 5 altered to be [0, 0, 0, 0, 0, 1, 0, -1, 0, 0]
        //   All other rows just [2, 2, ..., 2]

        Mat& r_mat = linear_system.rGetLhsMatrix();

        PetscInt lo, hi;
        PetscMatTools::GetOwnershipRange(r_mat, lo, hi);
        for(int i=lo; i<hi; i++)
        {
            if(i==0 || i==1 || i==4 || i==5)
            {
                unsigned col_one = i;
                unsigned col_minus_one = i+2;

                for(unsigned j=0; j<2*SIZE; j++)
                {
                    double val = 0.0;
                    if(j==col_one)
                    {
                        val = 1.0;
                    }
                    if(j==col_minus_one)
                    {
                        val = -1.0;
                    }
                    TS_ASSERT_DELTA( PetscMatTools::GetElement(r_mat, i, j), val, 1e-12);
                }
            }
            else
            {
                for(unsigned j=0; j<2*SIZE; j++)
                {
                    TS_ASSERT_DELTA( PetscMatTools::GetElement(r_mat, i, j), 2.0, 1e-12);
                }
            }
        }

        for (unsigned i=0; i<SIZE-1; i++)
        {
            delete nodes[i];
        }
    }
    void TestHeatEquationWithCoupledOdeSystemIn1dWithMixed()
    {
        // Create mesh of domain [0,1]
        TrianglesMeshReader<1,1> mesh_reader("mesh/test/data/1D_0_to_1_100_elements");
        TetrahedralMesh<1,1> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        // Create PDE system object
        HeatEquationForCoupledOdeSystem<1> pde;

        // Define non-zero Neumann boundary condition at x=0
        BoundaryConditionsContainer<1,1,1> bcc;
        ConstBoundaryCondition<1>* p_boundary_condition = new ConstBoundaryCondition<1>(1.0);
        TetrahedralMesh<1,1>::BoundaryElementIterator iter = mesh.GetBoundaryElementIteratorBegin();
        bcc.AddNeumannBoundaryCondition(*iter, p_boundary_condition);

        // Define zero Dirichlet boundary condition at x=1
        ConstBoundaryCondition<1>* p_boundary_condition2 = new ConstBoundaryCondition<1>(0.0);
        TetrahedralMesh<1,1>::BoundaryNodeIterator node_iter = mesh.GetBoundaryNodeIteratorEnd();
        --node_iter;
        bcc.AddDirichletBoundaryCondition(*node_iter, p_boundary_condition2);

        // Create the correct number of ODE systems
        double a = 5.0;
        std::vector<AbstractOdeSystemForCoupledPdeSystem*> ode_systems;
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            ode_systems.push_back(new OdeSystemForCoupledHeatEquation(a));
        }

        // Create PDE system solver
        LinearParabolicPdeSystemWithCoupledOdeSystemSolver<1,1,1> solver(&mesh, &pde, &bcc, ode_systems);

        // Set end time and timestep
        double t_end = 0.1;
        solver.SetTimes(0, t_end);
        solver.SetTimeStep(0.001);

        // Set initial condition u(x,0) = 1 - x
        std::vector<double> init_cond(mesh.GetNumNodes());
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];
            init_cond[i] = 1 - x;
        }
        Vec initial_condition = PetscTools::CreateVec(init_cond);
        solver.SetInitialCondition(initial_condition);

        // Solve PDE system and store result
        Vec result = solver.Solve();
        ReplicatableVector result_repl(result);

        /*
         * Test that solution is given by
         *
         * u(x,t) = 1 - x,
         * v(x,t) = 1 + a*(1-x)*t,
         *
         * with t = t_end.
         */
        for (unsigned i=0; i<result_repl.GetSize(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];
            double u = 1 - x;
            TS_ASSERT_DELTA(result_repl[i], u, 0.1);

            double v = 1 + a*(1-x)*t_end;
            TS_ASSERT_DELTA(ode_systems[i]->rGetStateVariables()[0], v, 0.1);
        }

        // Tidy up
        PetscTools::Destroy(initial_condition);
        PetscTools::Destroy(result);
    }
    /**
     * Simple Parabolic PDE u' = del squared u + 1
     *
     * With u = -(1/6)(x^2+y^2+z^2) on the boundaries of the unit cube.
     *
     * Subject to the initial condition
     * u(0,x,y,z)=sin( PI x)sin( PI y)sin( PI z) - (1/6)(x^2+y^2+z^2)
     *
     */
    void TestSimpleLinearParabolicSolver3DZeroDirichWithSourceTerm()
    {
        // Create mesh from mesh reader
        TrianglesMeshReader<3,3> mesh_reader("mesh/test/data/cube_136_elements");
        TetrahedralMesh<3,3> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        // Instantiate PDE object
        HeatEquationWithSourceTerm<3> pde;

        // Boundary conditions
        BoundaryConditionsContainer<3,3,1> bcc;
        TetrahedralMesh<3,3>::BoundaryNodeIterator iter = mesh.GetBoundaryNodeIteratorBegin();

        while (iter < mesh.GetBoundaryNodeIteratorEnd())
        {
            double x = (*iter)->GetPoint()[0];
            double y = (*iter)->GetPoint()[1];
            double z = (*iter)->GetPoint()[2];
            ConstBoundaryCondition<3>* p_dirichlet_boundary_condition =
                new ConstBoundaryCondition<3>(-1.0/6*(x*x+y*y+z*z));
            bcc.AddDirichletBoundaryCondition(*iter, p_dirichlet_boundary_condition);
            iter++;
        }

        // Solver
        SimpleLinearParabolicSolver<3,3> solver(&mesh,&pde,&bcc);

        // Initial condition u(0,x) = sin(x*pi)*sin(y*pi)*sin(z*pi)-1/6*(x^2+y^2+z^2)
        Vec initial_condition = PetscTools::CreateVec(mesh.GetNumNodes());

        double* p_initial_condition;
        VecGetArray(initial_condition, &p_initial_condition);

        int lo, hi;
        VecGetOwnershipRange(initial_condition, &lo, &hi);
        for (int global_index = lo; global_index < hi; global_index++)
        {
            int local_index = global_index - lo;
            double x = mesh.GetNode(global_index)->GetPoint()[0];
            double y = mesh.GetNode(global_index)->GetPoint()[1];
            double z = mesh.GetNode(global_index)->GetPoint()[2];
            p_initial_condition[local_index] = sin(x*M_PI)*sin(y*M_PI)*sin(z*M_PI)-1.0/6*(x*x+y*y+z*z);
        }
        VecRestoreArray(initial_condition, &p_initial_condition);

        double t_end = 0.1;
        solver.SetTimes(0, t_end);
        solver.SetTimeStep(0.01);

        solver.SetInitialCondition(initial_condition);
        Vec result = solver.Solve();

        // Check result
        double* p_result;
        VecGetArray(result, &p_result);

        // Solution should be u = e^{-t*2*pi*pi} sin(x*pi) sin(y*pi) sin(z*pi) - 1/6(x^2+y^2+z^2), t=0.1
        for (int global_index = lo; global_index < hi; global_index++)
        {
            int local_index = global_index - lo;
            double x = mesh.GetNode(global_index)->GetPoint()[0];
            double y = mesh.GetNode(global_index)->GetPoint()[1];
            double z = mesh.GetNode(global_index)->GetPoint()[2];
            double u = exp(-t_end*3*M_PI*M_PI)*sin(x*M_PI)*sin(y*M_PI)*sin(z*M_PI)-1.0/6*(x*x+y*y+z*z);
            TS_ASSERT_DELTA(p_result[local_index], u, 0.1);
        }
        VecRestoreArray(result, &p_result);
        PetscTools::Destroy(initial_condition);
        PetscTools::Destroy(result);
    }
    /**
     * Simple Parabolic PDE u' = del squared u
     *
     * With u = x on 5 boundaries of the unit cube, and
     * u_n = 1 on the x face of the cube.
     *
     * Subject to the initial condition
     * u(0,x,y,z)=sin( PI x)sin( PI y)sin( PI z) + x
     */
    void TestSimpleLinearParabolicSolver3DNeumannOnCoarseMesh()
    {
        // Create mesh from mesh reader
        TrianglesMeshReader<3,3> mesh_reader("mesh/test/data/cube_136_elements");

        TetrahedralMesh<3,3> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        // Instantiate PDE object
        HeatEquation<3> pde;

        // Boundary conditions
        BoundaryConditionsContainer<3,3,1> bcc;
        TetrahedralMesh<3,3>::BoundaryNodeIterator iter = mesh.GetBoundaryNodeIteratorBegin();

        while (iter != mesh.GetBoundaryNodeIteratorEnd())
        {
            double x = (*iter)->GetPoint()[0];
            double y = (*iter)->GetPoint()[1];
            double z = (*iter)->GetPoint()[2];


            if ((fabs(y) < 0.01) || (fabs(y - 1.0) < 0.01) ||
                (fabs(x) < 0.01) ||
                (fabs(z) < 0.01) || (fabs(z - 1.0) < 0.01) )
            {
                ConstBoundaryCondition<3>* p_dirichlet_boundary_condition =
                    new ConstBoundaryCondition<3>(x);
                bcc.AddDirichletBoundaryCondition(*iter, p_dirichlet_boundary_condition);
            }

            iter++;
        }

        TetrahedralMesh<3,3>::BoundaryElementIterator surf_iter = mesh.GetBoundaryElementIteratorBegin();
        ConstBoundaryCondition<3>* p_neumann_boundary_condition =
            new ConstBoundaryCondition<3>(1.0);

        while (surf_iter != mesh.GetBoundaryElementIteratorEnd())
        {
            int node = (*surf_iter)->GetNodeGlobalIndex(0);
            double x = mesh.GetNode(node)->GetPoint()[0];

            if (fabs(x - 1.0) < 0.01)
            {
                bcc.AddNeumannBoundaryCondition(*surf_iter, p_neumann_boundary_condition);
            }

            surf_iter++;
        }

        // Solver
        SimpleLinearParabolicSolver<3,3> solver(&mesh,&pde,&bcc);

        // Initial condition u(0,x,y) = sin(0.5*PI*x)*sin(PI*y)+x
        Vec initial_condition = PetscTools::CreateVec(mesh.GetNumNodes());

        double* p_initial_condition;
        VecGetArray(initial_condition, &p_initial_condition);

        int lo, hi;
        VecGetOwnershipRange(initial_condition, &lo, &hi);
        for (int global_index = lo; global_index < hi; global_index++)
        {
            int local_index = global_index - lo;
            double x = mesh.GetNode(global_index)->GetPoint()[0];
            double y = mesh.GetNode(global_index)->GetPoint()[1];
            double z = mesh.GetNode(global_index)->GetPoint()[2];

            p_initial_condition[local_index] = sin(0.5*M_PI*x)*sin(M_PI*y)*sin(M_PI*z)+x;
        }
        VecRestoreArray(initial_condition, &p_initial_condition);

        solver.SetTimes(0, 0.1);
        solver.SetTimeStep(0.01);

        solver.SetInitialCondition(initial_condition);
        Vec result = solver.Solve();

        // Check result
        double* p_result;
        VecGetArray(result, &p_result);

        // Solution should be u = e^{-5/2*PI*PI*t} sin(0.5*PI*x)*sin(PI*y)*sin(PI*z)+x, t=0.1
        for (int global_index = lo; global_index < hi; global_index++)
        {
            int local_index = global_index - lo;
            double x = mesh.GetNode(global_index)->GetPoint()[0];
            double y = mesh.GetNode(global_index)->GetPoint()[1];
            double z = mesh.GetNode(global_index)->GetPoint()[2];

            double u = exp((-5/2)*M_PI*M_PI*0.1) * sin(0.5*M_PI*x) * sin(M_PI*y)* sin(M_PI*z) + x;
            TS_ASSERT_DELTA(p_result[local_index], u, u*0.15);
        }
        VecRestoreArray(result, &p_result);
        PetscTools::Destroy(initial_condition);
        PetscTools::Destroy(result);
    }
    // test 2D problem - takes a long time to run.
    // solution is incorrect to specified tolerance.
    void xTestSimpleLinearParabolicSolver2DNeumannWithSmallTimeStepAndFineMesh()
    {
        // Create mesh from mesh reader
        FemlabMeshReader<2,2> mesh_reader("mesh/test/data/",
                                          "femlab_fine_square_nodes.dat",
                                          "femlab_fine_square_elements.dat",
                                          "femlab_fine_square_edges.dat");

        TetrahedralMesh<2,2> mesh;
        mesh.ConstructFromMeshReader(mesh_reader);

        // Instantiate PDE object
        HeatEquation<2> pde;

        // Boundary conditions - zero dirichlet on boundary;
        BoundaryConditionsContainer<2,2,1> bcc;
        TetrahedralMesh<2,2>::BoundaryNodeIterator iter = mesh.GetBoundaryNodeIteratorBegin();

        while (iter != mesh.GetBoundaryNodeIteratorEnd())
        {
            double x = (*iter)->GetPoint()[0];
            double y = (*iter)->GetPoint()[1];

            ConstBoundaryCondition<2>* p_dirichlet_boundary_condition =
                new ConstBoundaryCondition<2>(x);

            if (fabs(y) < 0.01)
            {
                bcc.AddDirichletBoundaryCondition(*iter, p_dirichlet_boundary_condition);
            }

            if (fabs(y - 1.0) < 0.01)
            {
                bcc.AddDirichletBoundaryCondition(*iter, p_dirichlet_boundary_condition);
            }

            if (fabs(x) < 0.01)
            {
                bcc.AddDirichletBoundaryCondition(*iter, p_dirichlet_boundary_condition);
            }

            iter++;
        }

        TetrahedralMesh<2,2>::BoundaryElementIterator surf_iter = mesh.GetBoundaryElementIteratorBegin();
        ConstBoundaryCondition<2>* p_neumann_boundary_condition =
            new ConstBoundaryCondition<2>(1.0);

        while (surf_iter != mesh.GetBoundaryElementIteratorEnd())
        {
            int node = (*surf_iter)->GetNodeGlobalIndex(0);
            double x = mesh.GetNode(node)->GetPoint()[0];

            if (fabs(x - 1.0) < 0.01)
            {
                bcc.AddNeumannBoundaryCondition(*surf_iter, p_neumann_boundary_condition);
            }

            surf_iter++;
        }

        // Solver
        SimpleLinearParabolicSolver<2,2> solver(&mesh,&pde,&bcc);

        // Initial condition u(0,x,y) = sin(0.5*M_PI*x)*sin(M_PI*y)+x
        std::vector<double> init_cond(mesh.GetNumNodes());
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];
            double y = mesh.GetNode(i)->GetPoint()[1];
            init_cond[i] = sin(0.5*M_PI*x)*sin(M_PI*y)+x;
        }
        Vec initial_condition = PetscTools::CreateVec(init_cond);

        double t_end = 0.1;
        solver.SetTimes(0, t_end);
        solver.SetTimeStep(0.001);

        solver.SetInitialCondition(initial_condition);

        Vec result = solver.Solve();
        ReplicatableVector result_repl(result);

        // Check solution is u = e^{-5/4*M_PI*M_PI*t} sin(0.5*M_PI*x)*sin(M_PI*y)+x, t=0.1
        for (unsigned i=0; i<result_repl.GetSize(); i++)
        {
            double x = mesh.GetNode(i)->GetPoint()[0];
            double y = mesh.GetNode(i)->GetPoint()[1];
            double u = exp((-5/4)*M_PI*M_PI*t_end) * sin(0.5*M_PI*x) * sin(M_PI*y) + x;
            TS_ASSERT_DELTA(result_repl[i], u, 0.001);
        }

        PetscTools::Destroy(result);
        PetscTools::Destroy(initial_condition);
    }
    void TestSolvingEllipticPde() throw(Exception)
    {
        /* First we declare a mesh reader which reads mesh data files of the 'Triangle'
         * format. The path given is relative to the main Chaste directory. As we are in 2d,
         * the reader will look for three datafiles, [name].nodes, [name].ele and [name].edge.
         * Note that the first template argument here is the spatial dimension of the
         * elements in the mesh ({{{ELEMENT_DIM}}}), and the second is the dimension of the nodes,
         * i.e. the dimension of the space the mesh lives in ({{{SPACE_DIM}}}). Usually
         * {{{ELEMENT_DIM}}} and {{{SPACE_DIM}}} will be equal. */
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/square_128_elements");
        /* Now declare a tetrahedral mesh with the same dimensions... */
        TetrahedralMesh<2,2> mesh;
        /* ... and construct the mesh using the mesh reader. */
        mesh.ConstructFromMeshReader(mesh_reader);

        /* Next we instantiate an instance of our PDE we wish to solve. */
        MyPde pde;

        /* A set of boundary conditions are stored in a {{{BoundaryConditionsContainer}}}. The
         * three template arguments are ELEMENT_DIM, SPACE_DIM and PROBLEM_DIM, the latter being
         * the number of unknowns we are solving for. We have one unknown (ie u is a scalar, not
         * a vector), so in this case {{{PROBLEM_DIM}}}=1. */
        BoundaryConditionsContainer<2,2,1> bcc;

        /* Defining the boundary conditions is the only particularly fiddly part of solving PDEs,
         * unless they are very simple, such as u=0 on the boundary, which could be done
         * as follows: */
        //bcc.DefineZeroDirichletOnMeshBoundary(&mesh);

        /* We want to specify u=0 on x=0 and y=0.  To do this, we first create the boundary condition
         * object saying what the value of the condition is at any particular point in space.  Here
         * we use the class `ConstBoundaryCondition`, a subclass of `AbstractBoundaryCondition` that
         * yields the same constant value (0.0 here) everywhere it is used.
         *
         * Note that the object is allocated with `new`.  The `BoundaryConditionsContainer` object deals
         * with deleting its associated boundary condition objects.  Note too that we could allocate a
         * separate condition object for each boundary node, but using the same object where possible is
         * more memory efficient.
         */
        ConstBoundaryCondition<2>* p_zero_boundary_condition = new ConstBoundaryCondition<2>(0.0);
        /* We then get a boundary node iterator from the mesh... */
        TetrahedralMesh<2,2>::BoundaryNodeIterator iter = mesh.GetBoundaryNodeIteratorBegin();
        /* ...and loop over the boundary nodes, getting the x and y values. */
        while (iter < mesh.GetBoundaryNodeIteratorEnd())
        {
            double x = (*iter)->GetPoint()[0];
            double y = (*iter)->GetPoint()[1];
            /* If x=0 or y=0... */
            if ((x==0) || (y==0))
            {
                /* ...associate the zero boundary condition created above with this boundary node
                 * ({{{*iter}}} being a pointer to a {{{Node<2>}}}).
                 */
                bcc.AddDirichletBoundaryCondition(*iter, p_zero_boundary_condition);
            }
            iter++;
        }

        /* Now we create Neumann boundary conditions for the ''surface elements'' on x=1 and y=1. Note that
         * Dirichlet boundary conditions are defined on nodes, whereas Neumann boundary conditions are
         * defined on surface elements. Note also that the natural boundary condition statement for this
         * PDE is (D grad u).n = g(x) (where n is the outward-facing surface normal), and g(x) is a prescribed
         * function, ''not'' something like du/dn=g(x). Hence the boundary condition we are specifying is
         * (D grad u).n = 0.
         *
         * '''Important note for 1D:''' This means that if we were solving 2u,,xx,,=f(x) in 1D, and
         * wanted to specify du/dx=1 on the LHS boundary, the Neumann boundary value we have to specify is
         * -2, as n=-1 (outward facing normal) so (D gradu).n = -2 when du/dx=1.
         *
         * To define Neumann bcs, we reuse the zero boundary condition object defined above, but apply it
         * at surface elements.  We loop over these using another iterator provided by the mesh class.
         */
        TetrahedralMesh<2,2>::BoundaryElementIterator surf_iter
            = mesh.GetBoundaryElementIteratorBegin();
        while (surf_iter != mesh.GetBoundaryElementIteratorEnd())
        {
            /* Get the x and y values of any node (here, the 0th) in the element. */
            unsigned node_index = (*surf_iter)->GetNodeGlobalIndex(0);
            double x = mesh.GetNode(node_index)->GetPoint()[0];
            double y = mesh.GetNode(node_index)->GetPoint()[1];

            /* If x=1 or y=1... */
            if ( (fabs(x-1.0) < 1e-6) || (fabs(y-1.0) < 1e-6) )
            {
                /* ...associate the boundary condition with the surface element. */
                bcc.AddNeumannBoundaryCondition(*surf_iter, p_zero_boundary_condition);
            }

            /* Finally increment the iterator. */
            surf_iter++;
        }

        /* Next we define the solver of the PDE.
         * To solve an {{{AbstractLinearEllipticPde}}} (which is the type of PDE {{{MyPde}}} is),
         * we use a {{{SimpleLinearEllipticSolver}}}. The solver, again templated over
         * {{{ELEMENT_DIM}}} and {{{SPACE_DIM}}}, needs to be given (pointers to) the mesh,
         * pde and boundary conditions.
         */
        SimpleLinearEllipticSolver<2,2> solver(&mesh, &pde, &bcc);

        /* To solve, just call {{{Solve()}}}. A PETSc vector is returned. */
        Vec result = solver.Solve();

        /* It is a pain to access the individual components of a PETSc vector, even when running only on
         * one process. A helper class called {{{ReplicatableVector}}} has been created. Create
         * an instance of one of these, using the PETSc {{{Vec}}} as the data. The ''i''th
         * component of {{{result}}} can now be obtained by simply doing {{{result_repl[i]}}}.
         */
        ReplicatableVector result_repl(result);

        /* Let us write out the solution to a file. To do this, create an
         * {{{OutputFileHandler}}}, passing in the directory we want files written to.
         * This is relative to the directory defined by the CHASTE_TEST_OUTPUT environment
         * variable - usually `/tmp/$USER/testoutput`. Note by default the output directory
         * passed in is emptied by this command. To avoid this, {{{false}}} can be passed in as a second
         * parameter.
         */
        OutputFileHandler output_file_handler("TestSolvingLinearPdeTutorial");

        /* Create an {{{out_stream}}}, which is a stream to a particular file. An {{{out_stream}}}
         * is a smart pointer to a `std::ofstream`. */
        out_stream p_file = output_file_handler.OpenOutputFile("linear_solution.txt");

        /* Loop over the entries of the solution. */
        for (unsigned i=0; i<result_repl.GetSize(); i++)
        {
            /* Get the x and y-values of the node corresponding to this entry. The method
             * {{{GetNode}}} on the mesh class returns a pointer to a {{{Node}}}. */
            double x = mesh.GetNode(i)->rGetLocation()[0];
            double y = mesh.GetNode(i)->rGetLocation()[1];

            /* Get the computed solution at this node from the {{{ReplicatableVector}}}. */
            double u = result_repl[i];

            /* Finally, write x, y and u to the output file. The solution could then be
             * visualised in (eg) matlab, using the commands:
             * {{{sol=load('linear_solution.txt'); plot3(sol(:,1),sol(:,2),sol(:,3),'.');}}}*/
            (*p_file) << x << " " << y << " " << u << "\n";
        }

        /* All PETSc {{{Vec}}}s should be destroyed when they are no longer needed, or you will have a memory leak. */
        PetscTools::Destroy(result);
    }
    /*
     * == Test 2: Solving a linear parabolic PDE ==
     *
     * Now we solve a parabolic PDE. We choose a simple problem so that the code changes
     * needed from the elliptic case are clearer. We will solve
     * du/dt = div(grad u) + u, in 3d, with boundary conditions u=1 on the boundary, and initial
     * conditions u=1.
     *
     */
    void TestSolvingParabolicPde() throw(Exception)
    {
        /* Create a 10 by 10 by 10 mesh in 3D, this time using the {{{ConstructRegularSlabMesh}}} method
         * on the mesh. The first parameter is the cartesian space-step and the other three parameters are the width, height and depth of the mesh.*/
        TetrahedralMesh<3,3> mesh;
        mesh.ConstructRegularSlabMesh(0.1, 1.0, 1.0, 1.0);

        /* Our PDE object should be a class that is derived from the {{{AbstractLinearParabolicPde}}}.
         * We could write it ourselves as in the previous test, but since the PDE we want to solve is
         * so simple, it has already been defined (look it up! - it is located in pde/test/pdes).
         */
        HeatEquationWithSourceTerm<3> pde;

        /* Create a new boundary conditions container and specify u=1.0 on the boundary. */
        BoundaryConditionsContainer<3,3,1> bcc;
        bcc.DefineConstantDirichletOnMeshBoundary(&mesh, 1.0);

        /* Create an instance of the solver, passing in the mesh, pde and boundary conditions.
         */
        SimpleLinearParabolicSolver<3,3> solver(&mesh,&pde,&bcc);

        /* For parabolic problems, initial conditions are also needed. The solver will expect
         * a PETSc vector, where the i-th entry is the initial solution at node i, to be passed
         * in. To create this PETSc {{{Vec}}}, we will use a helper function in the {{{PetscTools}}}
         * class to create a {{{Vec}}} of size num_nodes, with each entry set to 1.0. Then we
         * set the initial condition on the solver. */
        Vec initial_condition = PetscTools::CreateAndSetVec(mesh.GetNumNodes(), 1.0);
        solver.SetInitialCondition(initial_condition);

        /* Next define the start time, end time, and timestep, and set them. */
        double t_start = 0;
        double t_end = 1;
        double dt = 0.01;
        solver.SetTimes(t_start, t_end);
        solver.SetTimeStep(dt);


        /* HOW_TO_TAG PDE
         * Output results to file for time-dependent PDE solvers
         */


        /* When we call Solve() below we will just get the solution at the final time. If we want
         * to have intermediate solutions written to file, we do the following. We start by
         * specifying an output directory and filename prefix for our results file:
         */
        solver.SetOutputDirectoryAndPrefix("ParabolicSolverTutorial","results");
        /* When an output directory has been specified, the solver writes output in HDF5 format. To
         * convert this to another output format, we call the relevant method. Here, we convert
         * the output to plain text files (VTK or cmgui formats are also possible). We also say how
         * often to write the data, telling the solver to output results to file every tenth timestep.
         * The solver will create one file for each variable (in this case there is only one variable)
         * and for each time, so for example, the file
         * `results_Variable_0_10` is the results for u, over all nodes, at the 11th printed time.
         * Have a look in the output directory after running the test. (For comments on visualising the data in
         * matlab or octave, see the end of the tutorial UserTutorials/WritingPdeSolvers.)
         */
        solver.SetOutputToTxt(true);
        solver.SetPrintingTimestepMultiple(10);

        /* Now we can solve the problem. The {{{Vec}}} that is returned can be passed into a
         * {{{ReplicatableVector}}} as before.
         */
        Vec solution = solver.Solve();
        ReplicatableVector solution_repl(solution);

        /* Let's also solve the equivalent static PDE, i.e. set du/dt=0, so 0=div(gradu) + u. This
         * is easy, as the PDE class has already been defined. */
        SimplePoissonEquation<3,3> static_pde;
        SimpleLinearEllipticSolver<3,3> static_solver(&mesh, &static_pde, &bcc);
        Vec static_solution = static_solver.Solve();
        ReplicatableVector static_solution_repl(static_solution);

        /* We can now compare the solution of the parabolic PDE at t=1 with the static solution,
         * to see if the static equilibrium solution was reached in the former. (Ideally we should
         * compute some relative error, but we just compute an absolute error for simplicity.) */
        for (unsigned i=0; i<static_solution_repl.GetSize(); i++)
        {
            TS_ASSERT_DELTA( solution_repl[i], static_solution_repl[i], 1e-3);
        }

        /* All PETSc vectors should be destroyed when they are no longer needed. */
        PetscTools::Destroy(initial_condition);
        PetscTools::Destroy(solution);
        PetscTools::Destroy(static_solution);
    }