Dear all,

I am having some trouble with two commands.   The first one I am trying to use 
is based off an e-mail from Sept 4th and Sept 6, 2009 between Martin and 
Wolfgang.  Part of it was 

------------------------

  std::map<unsigned int, double> boundary_values;
  for (cell = dof_handler.begin_active(); ...)
    for (f=0; f<GeometryInfo<dim>::faces_per_cell;...)
      if (cell->face(f)->at_boundary())
        {
          cell->face(f)->get_dof_indices(face_dof_indices);
          for (unsigned int i=0; ...)
            if (fe.face_system_to_component_index(i) belongs to one of
                the vector components for which you use the RT element)
              boundary_values[face_dof_indices[i]] = 0;
        }

Does that make sense?

W.

----------------------------

I am trying to use  "cell->face(f)->get_dof_indices(face_dof_indices); " and 
have not had any success.  I was hoping to get these global indices for 
evaluating my right hand side in the Navier-Stokes code I am trying to write.  
The error that I am getting is 

" 'struct dealii::TriaAccessor < 1, 2, 2 >' has no member named 
'get_dof_indices' "

Since this was not working, I was looking at a second command.  

DoFTools::extract_boundary_dofs (dof_handler, component_select, selected_dofs, 
boundary_indicators);   

However, I have been having trouble with the boundary_indicators argument.  I 
do not get any errors when I pass in the 'boundary_indicators' argument with a 
set<unsigned char> , but I also do not get any information back.  The code I am 
running to test out both of these commands is below (minus the includes).  The 
two lines using these commands are commented out.

Thanks for the help,
Dan Brauss


 
using namespace dealii;

int main ()
{
    const unsigned int degree = 1;
    const unsigned int dim = 2;

    Triangulation<dim>      triangulation;

    SparsityPattern       sparsity_pattern;

    FESystem<dim>           fe (FE_Q<dim>(degree+1), dim,
                                FE_Q<dim>(degree), 1);
    DoFHandler<dim>         dof_handler(triangulation);


    std::vector<unsigned int> subdivisions (dim,1);
    subdivisions[0] = 4;

    const Point<dim> bottom_left = (dim == 2 ?
                                     Point<dim> (-1,-1) : 
                                     Point<dim> (-1,0,-1));
    const Point<dim> top_right   = (dim == 2 ?
                                     Point<dim> (1,1) : 
                                     Point<dim> (1,1,0));
    GridGenerator::subdivided_hyper_rectangle (triangulation,
                                               subdivisions,
                                               bottom_left,
                                               top_right);

    triangulation.refine_global (0);

    for ( Triangulation<dim>::active_cell_iterator
          cell = triangulation.begin_active();
        cell != triangulation.end(); ++cell)
      for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
        {
          if (cell->face(f)->center()[0] == -1)  // x-component is -1
              cell->face(f)->set_all_boundary_indicators(1);

          if (cell->face(f)->center()[1] == 1)   // y-component is 1
              cell->face(f)->set_all_boundary_indicators(1);

          if (cell->face(f)->center()[1] == -1)   // y-component is -1
              cell->face(f)->set_all_boundary_indicators(1);

          if (cell->face(f)->center()[0] == 1)   // x-component is 1
              cell->face(f)->set_all_boundary_indicators(0);
 
        }  

    dof_handler.distribute_dofs(fe);


    std::vector<unsigned int> face_dof_indices (fe.dofs_per_face,0);

    for ( Triangulation<dim>::active_cell_iterator
          cell = triangulation.begin_active();
        cell != triangulation.end(); ++cell)
      for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
        {
          if (cell->face(f)->center()[0] == -1)  //x -component is -1
            {
/*              cell->face(f)->get_dof_indices(face_dof_indices);               
    */
              for (unsigned int k=0; k<fe.dofs_per_face; ++k)
                {
                  std::cout << "local dof on face = "
                            << k 
                            << "and global is "
                            << face_dof_indices[k]
                            << std::endl;
                }
            }
        }

    std::vector<bool> selected_dofs1 (dof_handler.n_dofs());
    std::vector<bool> selected_dofs2 (dof_handler.n_dofs());
    std::vector<bool> component_select (dim+1, false);
    component_select[dim] = true;

    std::set<unsigned char> boundary_indicators;
    unsigned char b = '1';    
    boundary_indicators.insert(b);

    DoFTools::extract_dofs (dof_handler,
                                     component_select,
                                     selected_dofs1);
                                     

    DoFTools::extract_boundary_dofs (dof_handler,
                                     component_select,
                                     selected_dofs2);

/*    DoFTools::extract_boundary_dofs (dof_handler,
                                     component_select,
                                     selected_dofs2,
                                     boundary_indicators);                      
         */

    for (unsigned int i=0; i<dof_handler.n_dofs(); ++i)
      {
        std::cout << "degree of freedom = " 
                  << i 
                  << " "
                  << "pressure basis fcn = "
                  << selected_dofs1[i]
                  << " "
                  << "boundary dof = "
                  << selected_dofs2[i]
                  << "(true=1 / false=0)"
                  << std::endl;
      }

 


_______________________________________________
dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii

Reply via email to