And for the record, attached is a corrected version of the MWE illustrating the correct and intended behaviour of constraint creation for pure p-refinement.
Number of cells: 2 Number of degrees of freedom: 11 Number of constraints: 1 On Saturday, August 27, 2016 at 9:41:41 AM UTC+2, Jean-Paul Pelteret wrote: > > Ha ha. Dammit, what a silly mistake :-) Thanks Wolfgang! > > On Friday, August 26, 2016 at 11:08:17 PM UTC+2, Wolfgang Bangerth wrote: >> >> On 08/25/2016 02:34 AM, Jean-Paul Pelteret wrote: >> > Hi Deepak, >> > >> > Ok, so I've created a minimal working example to demonstrate what >> you're >> > seeing here. I'm not sufficiently familiar with the p-refinement part >> of the >> > hp functionality/implementation to tell whether this is truly the >> expected >> > behaviour or not. I see that the examples and description >> > < >> https://www.dealii.org/8.4.1/doxygen/deal.II/group__constraints.html#ga3eaa31a679484e80c193e74e8a967dc8> >> >> >> > in the documentation mix h- and p-refinement, so are _perhaps_ not >> > sufficiently clear to describe the expected outcome of this >> configuration. Or >> > maybe its a bug arising from the adjacent cells being on the same >> refinement >> > level (i.e. there are no actual hanging _nodes_). Perhaps someone else >> could >> > clarify this? >> >> It's not a bug but a feature :-) In your code, you call >> >> >> // Build constraints >> ConstraintMatrix constraints; >> DoFTools::make_hanging_node_constraints (dof_handler, >> constraints); >> constraints.clear(); // ************* >> >> I suspect what you really wanted to call in the last line is >> 'constraints.close()' -- because if you do so, the result is pretty much >> what >> I think you probably wanted ;-) >> >> Cheers >> W. >> >> -- >> ------------------------------------------------------------------------ >> Wolfgang Bangerth email: bange...@colostate.edu >> www: http://www.math.colostate.edu/~bangerth/ >> >> -- 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. For more options, visit https://groups.google.com/d/optout.
#include <deal.II/dofs/dof_handler.h> #include <deal.II/dofs/dof_tools.h> #include <deal.II/fe/fe_q.h> #include <deal.II/fe/fe_values.h> #include <deal.II/grid/grid_generator.h> #include <deal.II/grid/tria.h> #include <deal.II/grid/tria_accessor.h> #include <deal.II/grid/tria_iterator.h> #include <deal.II/hp/dof_handler.h> #include <deal.II/hp/fe_collection.h> #include <deal.II/hp/fe_values.h> #include <deal.II/lac/constraint_matrix.h> #include <deal.II/numerics/data_out.h> #include <iostream> #include <fstream> using namespace dealii; int main() { const int dim = 2; Triangulation<dim> tria; const hp::FECollection<dim> fe_collection (FE_Q<dim>(1), FE_Q<dim>(2)); hp::DoFHandler<dim> dof_handler(tria); // Two cell domain std::vector<unsigned int> divisions (dim,1); divisions[0] = 2; Point<dim> pt_1, pt_2; for (unsigned int d=0; d<dim; ++d) pt_2[d] = 1.0; GridGenerator::subdivided_hyper_rectangle(tria, divisions, pt_1, pt_2); // Set active FE index (p-refinement level; default index is 0) dof_handler.begin_active()->set_active_fe_index(1); dof_handler.distribute_dofs(fe_collection); // Build constraints ConstraintMatrix constraints; DoFTools::make_hanging_node_constraints (dof_handler, constraints); constraints.close(); // Print to screen std::cout << "Number of cells: " << tria.n_active_cells() << "\n" << "Number of degrees of freedom: " << dof_handler.n_dofs() << "\n" << "Number of constraints: " << constraints.n_constraints() << "\n"; // Output p-refinement level Vector<double> p_refinement (tria.n_active_cells()); unsigned int cell_no = 0; for (typename hp::DoFHandler<dim>::active_cell_iterator cell = dof_handler.begin_active(); cell != dof_handler.end(); ++cell, ++cell_no) { p_refinement[cell_no] = cell->active_fe_index(); } DataOut<dim,hp::DoFHandler<dim> > data_out; data_out.attach_dof_handler (dof_handler); data_out.add_data_vector (p_refinement, "p_refinement_level"); data_out.build_patches (); std::ostringstream filename; filename << "solution.vtk"; std::ofstream output (filename.str().c_str()); data_out.write_vtk (output); return 0; }