Reference documentation for deal.II version 9.1.0-pre
step-45.h
1  0,
109  /*b_id2*/ 1,
110  /*direction*/ 0,
111  matched_pairs);
112 @endcode
113 would yield periodicity constraints such that @f$u(0,y)=u(1,y)@f$ for all
114 @f$y\in[0,1]@f$.
115 
116 If we instead consider the parallelogram given by the convex hull of
117 @f$(0,0)@f$, @f$(1,1)@f$, @f$(1,2)@f$, @f$(0,1)@f$ we can achieve the constraints
118 @f$u(0,y)=u(1,y+1)@f$ by specifying an @p offset:
119 @code
121  /*b_id1*/ 0,
122  /*b_id2*/ 1,
123  /*direction*/ 0,
124  matched_pairs,
125  Tensor<1, 2>(0.,1.));
126 @endcode
127 or
128 @code
130  /*b_id1*/ 0,
131  /*b_id2*/ 1,
132  /*arbitrary direction*/ 0,
133  matched_pairs,
134  Tensor<1, 2>(1.,1.));
135 @endcode
136 Here, again, the assignment of boundary indicators 0 and 1 stems from
137 what GridGenerator::parallelogram() documents.
138 
139 The resulting @p matched_pairs can be used in
140 DoFTools::make_periodicity_constraints for populating a AffineConstraints
141 with periodicity constraints:
142 @code
143 DoFTools::make_periodicity_constraints(matched_pairs, constraints);
144 @endcode
145 
146 Apart from this high level interface there are also variants of
147 DoFTools::make_periodicity_constraints available that combine those two
148 steps (see the variants of DofTools::make_periodicity_constraints).
149 
150 There is also a low level interface to
151 DoFTools::make_periodicity_constraints if more flexibility is needed. The
152 low level variant allows to directly specify two faces that shall be
153 constrained:
154 @code
155 using namespace DoFTools;
156 make_periodicity_constraints(face_1,
157  face_2,
158  affine_constraints,
159  component_mask = <default value>;
160  face_orientation = <default value>,
161  face_flip = <default value>,
162  face_rotation = <default value>,
163  matrix = <default value>);
164 @endcode
165 Here, we need to specify the orientation of the two faces using
166 @p face_orientation, @p face_flip and @p face_orientation. For a closer description
167 have a look at the documentation of DoFTools::make_periodicity_constraints.
168 The remaining parameters are the same as for the high level interface apart
169 from the self-explaining @p component_mask and @p affine_constraints.
170 
171 
172 <a name="problem"></a>
173 <a name="Apracticalexample"></a><h1>A practical example</h1>
174 
175 
176 In the following, we show how to use the above functions in a more involved
177 example. The task is to enforce rotated periodicity constraints for the
178 velocity component of a Stokes flow.
179 
180 On a quarter-circle defined by @f$\Omega=\{{\bf x}\in(0,1)^2:\|{\bf x}\|\in (0.5,1)\}@f$ we are
181 going to solve the Stokes problem
182 @f{eqnarray*}
183  -\Delta \; \textbf{u} + \nabla p &=& (\exp(-100*\|{\bf x}-(.75,0.1)^T\|^2),0)^T, \\
184  -\textrm{div}\; \textbf{u}&=&0,\\
185  \textbf{u}|_{\Gamma_1}&=&{\bf 0},
186 @f}
187 where the boundary @f$\Gamma_1@f$ is defined as @f$\Gamma_1:=\{x\in \partial\Omega: \|x\|\in\{0.5,1\}\}@f$.
188 For the remaining parts of the boundary we are going to use periodic boundary conditions, i.e.
189 @f{align*}
190  u_x(0,\nu)&=-u_y(\nu,0)&\nu&\in[0,1]\\
191  u_y(0,\nu)&=u_x(\nu,0)&\nu&\in[0,1].
192 @f}
193 
194 The mesh will be generated by GridGenerator::quarter_hyper_shell(),
195 which also documents how it assigns boundary indicators to its various
196 boundaries if its `colorize` argument is set to `true`.
197  * <a name="CommProg"></a>
198  * <h1> The commented program</h1>
199  *
200  * This example program is a slight modification of @ref step_22 "step-22" running in parallel
201  * using Trilinos to demonstrate the usage of periodic boundary conditions in
202  * deal.II. We thus omit to discuss the majority of the source code and only
203  * comment on the parts that deal with periodicity constraints. For the rest
204  * have a look at @ref step_22 "step-22" and the full source code at the bottom.
205  *
206 
207  *
208  * In order to implement periodic boundary conditions only two functions
209  * have to be modified:
210  * - <code>StokesProblem<dim>::setup_dofs()</code>:
211  * To populate a AffineConstraints object with periodicity constraints
212  * - <code>StokesProblem<dim>::run()</code>:
213  * To supply a distributed triangulation with periodicity information.
214  *
215  *
216  * The rest of the program is identical to @ref step_22 "step-22", so let us skip this part
217  * and only show these two functions in the following. (The full program can be
218  * found in the "Plain program" section below, though.)
219  *
220 
221  *
222  *
223 
224  *
225  *
226  *
227  *
228  * <a name="Settingupperiodicityconstraintsondistributedtriangulations"></a>
229  * <h3>Setting up periodicity constraints on distributed triangulations</h3>
230  *
231  * @code
232  * template <int dim>
233  * void StokesProblem<dim>::create_mesh()
234  * {
235  * Point<dim> center;
236  * const double inner_radius = .5;
237  * const double outer_radius = 1.;
238  *
240  * triangulation, center, inner_radius, outer_radius, 0, true);
241  *
242  * @endcode
243  *
244  * Before we can prescribe periodicity constraints, we need to ensure that
245  * cells on opposite sides of the domain but connected by periodic faces are
246  * part of the ghost layer if one of them is stored on the local processor.
247  * At this point we need to think about how we want to prescribe
248  * periodicity. The vertices @f$\text{vertices}_2@f$ of a face on the left
249  * boundary should be matched to the vertices @f$\text{vertices}_1@f$ of a face
250  * on the lower boundary given by @f$\text{vertices}_2=R\cdot
251  * \text{vertices}_1+b@f$ where the rotation matrix @f$R@f$ and the offset @f$b@f$ are
252  * given by
253  * @f{align*}
254  * R=\begin{pmatrix}
255  * 0&1\\-1&0
256  * \end{pmatrix},
257  * \quad
258  * b=\begin{pmatrix}0&0\end{pmatrix}.
259  * @f}
260  * The data structure we are saving the resulting information into is here
261  * based on the Triangulation.
262  *
263  * @code
264  * std::vector<GridTools::PeriodicFacePair<
266  * periodicity_vector;
267  *
268  * FullMatrix<double> rotation_matrix(dim);
269  * rotation_matrix[0][1] = 1.;
270  * rotation_matrix[1][0] = -1.;
271  *
272  * GridTools::collect_periodic_faces(triangulation,
273  * 2,
274  * 3,
275  * 1,
276  * periodicity_vector,
277  * Tensor<1, dim>(),
278  * rotation_matrix);
279  *
280  * @endcode
281  *
282  * Now telling the triangulation about the desired periodicity is
283  * particularly easy by just calling
285  *
286  * @code
287  * triangulation.add_periodicity(periodicity_vector);
288  *
289  * triangulation.refine_global(4 - dim);
290  * }
291  *
292  *
293  * @endcode
294  *
295  *
296  * <a name="Settingupperiodicityconstraintsondistributedtriangulations"></a>
297  * <h3>Setting up periodicity constraints on distributed triangulations</h3>
298  *
299  * @code
300  * template <int dim>
301  * void StokesProblem<dim>::setup_dofs()
302  * {
303  * dof_handler.distribute_dofs(fe);
304  *
305  * std::vector<unsigned int> block_component(dim + 1, 0);
306  * block_component[dim] = 1;
307  * DoFRenumbering::component_wise(dof_handler, block_component);
308  *
309  * std::vector<types::global_dof_index> dofs_per_block(2);
310  * DoFTools::count_dofs_per_block(dof_handler,
311  * dofs_per_block,
312  * block_component);
313  * const unsigned int n_u = dofs_per_block[0], n_p = dofs_per_block[1];
314  *
315  * {
316  * owned_partitioning.clear();
317  * IndexSet locally_owned_dofs = dof_handler.locally_owned_dofs();
318  * owned_partitioning.push_back(locally_owned_dofs.get_view(0, n_u));
319  * owned_partitioning.push_back(locally_owned_dofs.get_view(n_u, n_u + n_p));
320  *
321  * relevant_partitioning.clear();
322  * IndexSet locally_relevant_dofs;
324  * locally_relevant_dofs);
325  * relevant_partitioning.push_back(locally_relevant_dofs.get_view(0, n_u));
326  * relevant_partitioning.push_back(
327  * locally_relevant_dofs.get_view(n_u, n_u + n_p));
328  *
329  * constraints.clear();
330  * constraints.reinit(locally_relevant_dofs);
331  *
332  * FEValuesExtractors::Vector velocities(0);
333  *
334  * DoFTools::make_hanging_node_constraints(dof_handler, constraints);
336  * dof_handler,
337  * 0,
338  * BoundaryValues<dim>(),
339  * constraints,
340  * fe.component_mask(velocities));
342  * dof_handler,
343  * 1,
344  * BoundaryValues<dim>(),
345  * constraints,
346  * fe.component_mask(velocities));
347  *
348  * @endcode
349  *
350  * After we provided the mesh with the necessary information for the
351  * periodicity constraints, we are now able to actual create them. For
352  * describing the matching we are using the same approach as before, i.e.,
353  * the @f$\text{vertices}_2@f$ of a face on the left boundary should be
354  * matched to the vertices
355  * @f$\text{vertices}_1@f$ of a face on the lower boundary given by
356  * @f$\text{vertices}_2=R\cdot \text{vertices}_1+b@f$ where the rotation
357  * matrix @f$R@f$ and the offset @f$b@f$ are given by
358  * @f{align*}
359  * R=\begin{pmatrix}
360  * 0&1\\-1&0
361  * \end{pmatrix},
362  * \quad
363  * b=\begin{pmatrix}0&0\end{pmatrix}.
364  * @f}
365  * These two objects not only describe how faces should be matched but
366  * also in which sense the solution should be transformed from
367  * @f$\text{face}_2@f$ to
368  * @f$\text{face}_1@f$.
369  *
370  * @code
371  * FullMatrix<double> rotation_matrix(dim);
372  * rotation_matrix[0][1] = 1.;
373  * rotation_matrix[1][0] = -1.;
374  *
375  * Tensor<1, dim> offset;
376  *
377  * @endcode
378  *
379  * For setting up the constraints, we first store the periodicity
380  * information in an auxiliary object of type
381  * <code>std::vector@<GridTools::PeriodicFacePair<typename
382  * DoFHandler@<dim@>::cell_iterator@> </code>. The periodic boundaries
383  * have the boundary indicators 2 (x=0) and 3 (y=0). All the other
384  * parameters we have set up before. In this case the direction does not
385  * matter. Due to @f$\text{vertices}_2=R\cdot \text{vertices}_1+b@f$ this is
386  * exactly what we want.
387  *
388  * @code
389  * std::vector<
391  * periodicity_vector;
392  *
393  * const unsigned int direction = 1;
394  *
395  * GridTools::collect_periodic_faces(dof_handler,
396  * 2,
397  * 3,
398  * direction,
399  * periodicity_vector,
400  * offset,
401  * rotation_matrix);
402  *
403  * @endcode
404  *
405  * Next, we need to provide information on which vector valued components
406  * of the solution should be rotated. Since we choose here to just
407  * constraint the velocity and this starts at the first component of the
408  * solution vector, we simply insert a 0:
409  *
410  * @code
411  * std::vector<unsigned int> first_vector_components;
412  * first_vector_components.push_back(0);
413  *
414  * @endcode
415  *
416  * After setting up all the information in periodicity_vector all we have
417  * to do is to tell make_periodicity_constraints to create the desired
418  * constraints.
419  *
420  * @code
421  * DoFTools::make_periodicity_constraints<DoFHandler<dim>>(
422  * periodicity_vector,
423  * constraints,
424  * fe.component_mask(velocities),
425  * first_vector_components);
426  *
428  * dof_handler,
429  * 0,
430  * BoundaryValues<dim>(),
431  * constraints,
432  * fe.component_mask(velocities));
434  * dof_handler,
435  * 1,
436  * BoundaryValues<dim>(),
437  * constraints,
438  * fe.component_mask(velocities));
439  * }
440  *
441  * constraints.close();
442  *
443  * {
444  * TrilinosWrappers::BlockSparsityPattern bsp(owned_partitioning,
445  * owned_partitioning,
446  * relevant_partitioning,
447  * mpi_communicator);
448  *
449  * Table<2, DoFTools::Coupling> coupling(dim + 1, dim + 1);
450  * for (unsigned int c = 0; c < dim + 1; ++c)
451  * for (unsigned int d = 0; d < dim + 1; ++d)
452  * if (!((c == dim) && (d == dim)))
453  * coupling[c][d] = DoFTools::always;
454  * else
455  * coupling[c][d] = DoFTools::none;
456  *
457  * DoFTools::make_sparsity_pattern(dof_handler,
458  * coupling,
459  * bsp,
460  * constraints,
461  * false,
463  * mpi_communicator));
464  *
465  * bsp.compress();
466  *
467  * system_matrix.reinit(bsp);
468  * }
469  *
470  * {
471  * TrilinosWrappers::BlockSparsityPattern preconditioner_bsp(
472  * owned_partitioning,
473  * owned_partitioning,
474  * relevant_partitioning,
475  * mpi_communicator);
476  *
477  * Table<2, DoFTools::Coupling> preconditioner_coupling(dim + 1, dim + 1);
478  * for (unsigned int c = 0; c < dim + 1; ++c)
479  * for (unsigned int d = 0; d < dim + 1; ++d)
480  * if ((c == dim) && (d == dim))
481  * preconditioner_coupling[c][d] = DoFTools::always;
482  * else
483  * preconditioner_coupling[c][d] = DoFTools::none;
484  *
485  * DoFTools::make_sparsity_pattern(dof_handler,
486  * preconditioner_coupling,
487  * preconditioner_bsp,
488  * constraints,
489  * false,
491  * mpi_communicator));
492  *
493  * preconditioner_bsp.compress();
494  *
495  * preconditioner_matrix.reinit(preconditioner_bsp);
496  * }
497  *
498  * system_rhs.reinit(owned_partitioning, mpi_communicator);
499  * solution.reinit(owned_partitioning,
500  * relevant_partitioning,
501  * mpi_communicator);
502  * }
503  *
504  * @endcode
505  *
506  * The rest of the program is then again identical to @ref step_22 "step-22". We will omit
507  * it here now, but as before, you can find these parts in the "Plain program"
508  * section below.
509  *
510 
511  *
512 <a name="Results"></a><h1>Results</h1>
513 
514 
515 The created output is not very surprising. We simply see that the solution is
516 periodic with respect to the left and lower boundary:
517 
518 <img src="https://www.dealii.org/images/steps/developer/step-45.periodic.png" alt="">
519 
520 Without the periodicity constraints we would have ended up with the following solution:
521 
522 <img src="https://www.dealii.org/images/steps/developer/step-45.non_periodic.png" alt="">
523 <a name="PlainProg"></a>
524 <h1> The plain program</h1>
525 @include "step-45.cc"
526  */
void interpolate_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > * > &function_map, std::map< types::global_dof_index, number > &boundary_values, const ComponentMask &component_mask=ComponentMask())
Contents is actually a matrix.
Point< spacedim > point(const gp_Pnt &p, const double &tolerance=1e-10)
Definition: utilities.cc:180
void make_periodicity_constraints(const FaceIterator &face_1, const typename identity< FaceIterator >::type &face_2, AffineConstraints< double > &constraints, const ComponentMask &component_mask=ComponentMask(), const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const FullMatrix< double > &matrix=FullMatrix< double >(), const std::vector< unsigned int > &first_vector_components=std::vector< unsigned int >())
void quarter_hyper_shell(Triangulation< dim > &tria, const Point< dim > &center, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, const bool colorize=false)
void component_wise(DoFHandlerType &dof_handler, const std::vector< unsigned int > &target_component=std::vector< unsigned int >())
void parallelogram(Triangulation< dim > &tria, const Point< dim >(&corners)[dim], const bool colorize=false)
VectorizedArray< Number > exp(const ::VectorizedArray< Number > &x)
void make_hanging_node_constraints(const DoFHandlerType &dof_handler, AffineConstraints< number > &constraints)
void clear()
Definition: index_set.h:1587
void extract_locally_relevant_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:1143
IndexSet get_view(const size_type begin, const size_type end) const
Definition: index_set.cc:236
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
void count_dofs_per_block(const DoFHandlerType &dof, std::vector< types::global_dof_index > &dofs_per_block, const std::vector< unsigned int > &target_block=std::vector< unsigned int >())
Definition: dof_tools.cc:1969
Definition: mpi.h:55
void collect_periodic_faces(const MeshType &mesh, const types::boundary_id b_id1, const types::boundary_id b_id2, const int direction, std::vector< PeriodicFacePair< typename MeshType::cell_iterator >> &matched_pairs, const Tensor< 1, MeshType::space_dimension > &offset=::Tensor< 1, MeshType::space_dimension >(), const FullMatrix< double > &matrix=FullMatrix< double >())
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:80
Definition: table.h:37
void make_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
typename::Triangulation< dim, spacedim >::cell_iterator cell_iterator
Definition: tria.h:269
virtual void add_periodicity(const std::vector<::GridTools::PeriodicFacePair< cell_iterator >> &) override
Definition: tria.cc:4563