Hello, 

The VectorTools::interpolate_boundary_values function works for 
interpolating fe_Q<dim> functions but not fe_DGQ<dim> . 
Was this meant to be the case? 
I attached a minimal code that illustrates this.
   

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/0a8d1da6-e564-45de-bc3d-1f5a9d5dfad9n%40googlegroups.com.
#include <deal.II/grid/tria.h>
#include <deal.II/grid/grid_generator.h>

#include <deal.II/fe/fe_q.h>
#include <deal.II/fe/fe_dgq.h>

#include <deal.II/dofs/dof_handler.h>
#include <deal.II/dofs/dof_tools.h>

#include <deal.II/numerics/data_out.h>
#include <deal.II/base/function.h>
#include <deal.II/numerics/vector_tools.h>

#include <fstream>
using namespace dealii;
int main()

{
    // create an feq obkect and interpolate a value of 1 at boundary 0
    {
        Triangulation<2> tria;
        GridGenerator::hyper_cube(tria, 0, 1, true);
        tria.refine_global(4);
        FE_Q<2> fe(1);
        DoFHandler<2> dof(tria);

        Vector<double> solution;

        dof.distribute_dofs(fe);
        solution.reinit(dof.n_dofs());

        std::map<types::global_dof_index, double> boundary_values;
        VectorTools::interpolate_boundary_values(dof,
                                                 0,
                                                 Functions::ConstantFunction<2>(1),
                                                 boundary_values);
        for (auto &boundary_value : boundary_values)
            solution(boundary_value.first) = boundary_value.second;

        DataOut<2> data_out;

        data_out.attach_dof_handler(dof);
        data_out.add_data_vector(solution, "solution");
        data_out.build_patches();

        const std::string filename = "solution-feq.vtu";
        std::ofstream output(filename);
        data_out.write_vtu(output);
    }

    // Same as above but
    // we create an feGQq obkect and interpolate
    {
        Triangulation<2> tria;
        GridGenerator::hyper_cube(tria, 0, 1, true);
        tria.refine_global(4);
        FE_DGQ<2> fe(1);                 // only difference between this and the above code 
        DoFHandler<2> dof(tria);

        Vector<double> solution;

        dof.distribute_dofs(fe);
        solution.reinit(dof.n_dofs());

        std::map<types::global_dof_index, double> boundary_values;
        VectorTools::interpolate_boundary_values(dof,
                                                 0,
                                                 Functions::ConstantFunction<2>(1),
                                                 boundary_values);
        for (auto &boundary_value : boundary_values)
            solution(boundary_value.first) = boundary_value.second;

        DataOut<2> data_out;

        data_out.attach_dof_handler(dof);
        data_out.add_data_vector(solution, "solution");
        data_out.build_patches();

        const std::string filename = "solution-fedgq.vtu";
        std::ofstream output(filename);
        data_out.write_vtu(output);
    }
}

Reply via email to