Dear Bastian and Daniel, > That's what I was assuming. So yesterday I did it with this code : > // This block add to the periodicity constraint the little > compression we want > { > IndexSet selected_dofs_left; > std::set< types::boundary_id > boundary_ids_left= > std::set<types::boundary_id>(); > boundary_ids_left.insert(1); > DoFTools::extract_boundary_dofs(dof_handler_ref, > fe.component_mask(x_displacement), selected_dofs_left, > boundary_ids_left); > unsigned int nb_dofs_face_left = selected_dofs_left.n_elements(); > IndexSet::ElementIterator dofs_left = selected_dofs_left.begin(); > for(unsigned int i = 0; i < nb_dofs_face_left; i++) > { > //constraints.add_line(*dofs_left ); > constraints.set_inhomogeneity(*dofs_left, apply_dirichlet_bc ? > displacement_side_1 : 0.0); > dofs_left++; > } > } > And it works well. Thanks a lot for your help.I think this topic can be > set to *Answered*! >
I am going to do the same as Bastian, namely periodic boundary condition for displacement such that u( 0, y ) = u( 1, y ) + *lambda ,* * so I used the recommended above code by Bastian for applying inhomogenity to predefined periodicity, in order to have reletive dispacement between right and left faces,* *but my results shows that both left and right faces are moving equally to the same direction, not relative to each other.* *In addition, I activated the uncommented line, " * //constraints.add_line(*dofs_left );" otherwise I got this error: "call add_line() before calling set_inhomogeneity()" This is my periodic constraints( left side boundry_id = 0 and right side boundry_id = 1) : { std::vector<GridTools::PeriodicFacePair<typename DoFHandler<dim>::cell_iterator> > periodicity_vector_x; GridTools::collect_periodic_faces(dof_handler, /*b_id1*/ 0, /*b_id2*/ 1, /*direction*/ 0, periodicity_vector_x); DoFTools::make_periodicity_constraints<DoFHandler<dim> > (periodicity_vector_x, constraints, fe.component_mask(x_displacement)); } { IndexSet selected_dofs_right; std::set< types::boundary_id > boundary_ids_right= std::set<types::boundary_id>(); boundary_ids_right.insert(1); DoFTools::extract_boundary_dofs(dof_handler, fe.component_mask(x_displacement), selected_dofs_right, boundary_ids_right); unsigned int nb_dofs_face_right = selected_dofs_right.n_elements(); IndexSet::ElementIterator dofs_right = selected_dofs_right.begin(); for(unsigned int i = 0; i < nb_dofs_face_right; i++) { // constraints.add_line(*dofs_right ); constraints.set_inhomogeneity(*dofs_right, (apply_dirichlet_bc ? 0.01e-9 : 0.0)); dofs_right++; } } Any hint would be appreciated. Thanks, Hamed -- 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.