void VtkNonlinearElasticitySolutionWriter<DIM>::Write()
{
    if (mpSolver->mOutputDirectory=="")
    {
        EXCEPTION("No output directory was given to the mechanics solver");
    }

#ifdef CHASTE_VTK
    VtkMeshWriter<DIM, DIM> mesh_writer(mpSolver->mOutputDirectory + "/vtk", "solution", true);

    // write the displacement
    std::vector<c_vector<double,DIM> > displacement(mpSolver->mrQuadMesh.GetNumNodes());
    std::vector<c_vector<double,DIM> >& r_spatial_solution = mpSolver->rGetSpatialSolution();
    for(unsigned i=0; i<mpSolver->mrQuadMesh.GetNumNodes(); i++)
    {
        for(unsigned j=0; j<DIM; j++)
        {
            displacement[i](j) = r_spatial_solution[i](j)- mpSolver->mrQuadMesh.GetNode(i)->rGetLocation()[j];
        }
    }
    mesh_writer.AddPointData("Displacement", displacement);

    // write pressures
    if (mpSolver->mCompressibilityType==INCOMPRESSIBLE)
    {
        mesh_writer.AddPointData("Pressure", mpSolver->rGetPressures());
    }

    // write the element attribute as cell data.
    std::vector<double> element_attribute;
    for(typename QuadraticMesh<DIM>::ElementIterator iter = mpSolver->mrQuadMesh.GetElementIteratorBegin();
        iter != mpSolver->mrQuadMesh.GetElementIteratorEnd();
        ++iter)
    {
        element_attribute.push_back(iter->GetAttribute());
    }
    mesh_writer.AddCellData("Attribute", element_attribute);

    // write strains if requested
    if (mWriteElementWiseStrains)
    {
        mTensorData.clear();
        mTensorData.resize(mpSolver->mrQuadMesh.GetNumElements());

        std::string name;
        switch(mElementWiseStrainType)
        {
            case DEFORMATION_GRADIENT_F:
            {
                name = "deformation_gradient_F";
                break;
            }
            case DEFORMATION_TENSOR_C:
            {
                name = "deformation_tensor_C";
                break;
            }
            case LAGRANGE_STRAIN_E:
            {
                name = "Lagrange_strain_E";
                break;
            }
            default:
            {
                NEVER_REACHED;
            }
        }

        for (typename AbstractTetrahedralMesh<DIM,DIM>::ElementIterator iter = mpSolver->mrQuadMesh.GetElementIteratorBegin();
             iter != mpSolver->mrQuadMesh.GetElementIteratorEnd();
             ++iter)
        {
            mpSolver->GetElementCentroidStrain(mElementWiseStrainType, *iter, mTensorData[iter->GetIndex()]);
        }

        mesh_writer.AddTensorCellData(name, mTensorData);
    }
//// Future..
//        if (mWriteNodeWiseStresses)
//        {
//            std::vector<c_matrix<double,DIM,DIM> > tensor_data;
//            // use recoverer
//            mesh_writer.AddTensorCellData("Stress_NAME_ME", tensor_data);
//        }

    // final write
    mesh_writer.WriteFilesUsingMesh(mpSolver->mrQuadMesh);
#endif // CHASTE_VTK
}
double ElectrodesStimulusFactory<DIM>::ComputeElectrodeTotalFlux(AbstractChasteRegion<DIM>* pRegion, double stimulusMagnitude)
{
    GaussianQuadratureRule<DIM>* pQuadRule = new GaussianQuadratureRule<DIM>(2);

    //the basis functions
    c_vector<double, DIM+1> phi;

    double total_electrode_flux = 0.0;
    double ret;

    for (typename AbstractTetrahedralMesh<DIM,DIM>::NodeIterator node_iter=this->mpMesh->GetNodeIteratorBegin();
         node_iter != this->mpMesh->GetNodeIteratorEnd();
         ++node_iter)
    {
        if ( pRegion->DoesContain( (*node_iter).GetPoint() ) )
        {
            unsigned node_index = node_iter->GetIndex();
            assert(node_index < this->mpMesh->GetNumNodes());
            double contribution_of_this_node = 0.0;
            //loop over the elements where this node is contained
            for(std::set<unsigned>::iterator iter = this->mpMesh->GetNode(node_index)->rGetContainingElementIndices().begin();
                iter != this->mpMesh->GetNode(node_index)->rGetContainingElementIndices().end();
                ++iter)
            {
                Element<DIM,DIM>* p_element = this->mpMesh->GetElement(*iter);

                /*Determine jacobian for this element*/
                c_matrix<double, DIM, DIM> jacobian;
                c_matrix<double, DIM, DIM> inverse_jacobian;//unused here, but needed for the function below
                double jacobian_determinant;
                this->mpMesh->GetInverseJacobianForElement(p_element->GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);

                double contribution_of_this_element = 0.0;//...to this node
                 // loop over Gauss points
                for (unsigned quad_index=0; quad_index < pQuadRule->GetNumQuadPoints(); quad_index++)
                {
                    const ChastePoint<DIM>& quad_point = pQuadRule->rGetQuadPoint(quad_index);
                    BasisFunction::ComputeBasisFunctions(quad_point, phi);

                    double interpolated_stimulus = 0.0;
                    //loop over nodes in this element
                    for (unsigned node_index_in_element = 0; node_index_in_element < p_element->GetNumNodes(); node_index_in_element++)
                    {
                        //const Node<DIM>* p_node = p_element->GetNode(node_index_in_element);
                        assert(p_element->GetNumNodes() == DIM+1);
                        interpolated_stimulus += stimulusMagnitude*phi(node_index_in_element);
                        contribution_of_this_element += interpolated_stimulus * phi(node_index_in_element) * jacobian_determinant * pQuadRule->GetWeight(quad_index);
                    }

                }/*end of loop over gauss points*/
                contribution_of_this_node += contribution_of_this_element;

            }/*end of loop over elements where the node is contained*/
            total_electrode_flux += contribution_of_this_node;
        }/* end of if that checks if node is in the electrode*/
    }/* end of loop over nodes in the mesh*/
#ifndef NDEBUG
    int mpi_ret = MPI_Allreduce(&total_electrode_flux, &ret, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
    assert(mpi_ret == MPI_SUCCESS);
#else
    MPI_Allreduce(&total_electrode_flux, &ret, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
#endif

    //clear up memory
    delete pQuadRule;

    assert(ret < DBL_MAX);
    return ret;
}