Reference documentation for deal.II version 8.5.1
mg_tools.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/logstream.h>
17 #include <deal.II/base/thread_management.h>
18 #include <deal.II/lac/sparsity_pattern.h>
19 #include <deal.II/lac/block_sparsity_pattern.h>
20 #include <deal.II/lac/dynamic_sparsity_pattern.h>
21 #include <deal.II/lac/sparsity_pattern.h>
22 #include <deal.II/lac/block_vector.h>
23 #include <deal.II/lac/sparse_matrix.h>
24 #include <deal.II/lac/block_sparse_matrix.h>
25 #include <deal.II/lac/block_vector.h>
26 #include <deal.II/grid/tria.h>
27 #include <deal.II/grid/tria_iterator.h>
28 #include <deal.II/dofs/dof_accessor.h>
29 #include <deal.II/multigrid/mg_tools.h>
30 #include <deal.II/multigrid/mg_base.h>
31 #include <deal.II/base/mg_level_object.h>
32 #include <deal.II/dofs/dof_handler.h>
33 #include <deal.II/dofs/dof_tools.h>
34 #include <deal.II/fe/fe.h>
35 #include <deal.II/fe/component_mask.h>
36 #include <deal.II/numerics/matrix_tools.h>
37 
38 #include <vector>
39 #include <algorithm>
40 #include <numeric>
41 
42 DEAL_II_NAMESPACE_OPEN
43 
44 
45 namespace MGTools
46 {
47 
48 
49  // specializations for 1D
50  template <>
51  void
53  const DoFHandler<1,1> &,
54  const unsigned int,
55  std::vector<unsigned int> &,
56  const DoFTools::Coupling)
57  {
58  Assert(false, ExcNotImplemented());
59  }
60 
61 
62 
63  template <>
64  void
66  const DoFHandler<1,1> &,
67  const unsigned int,
68  std::vector<unsigned int> &,
71  {
72  Assert(false, ExcNotImplemented());
73  }
74 
75 
76 
77  template <>
78  void
80  const DoFHandler<1,2> &,
81  const unsigned int,
82  std::vector<unsigned int> &,
83  const DoFTools::Coupling)
84  {
85  Assert(false, ExcNotImplemented());
86  }
87 
88 
89  template <>
90  void
92  const DoFHandler<1,2> &,
93  const unsigned int,
94  std::vector<unsigned int> &,
97  {
98  Assert(false, ExcNotImplemented());
99  }
100 
101 
102 
103 // Template for 2D and 3D. For 1D see specialization above
104  template <int dim, int spacedim>
105  void
107  const DoFHandler<dim,spacedim> &dofs,
108  const unsigned int level,
109  std::vector<unsigned int> &row_lengths,
110  const DoFTools::Coupling flux_coupling)
111  {
112  Assert (row_lengths.size() == dofs.n_dofs(),
113  ExcDimensionMismatch(row_lengths.size(), dofs.n_dofs()));
114 
115  // Function starts here by
116  // resetting the counters.
117  std::fill(row_lengths.begin(), row_lengths.end(), 0);
118  // We need the user flags, so we
119  // save them for later restoration
120  std::vector<bool> old_flags;
121  // We need a non-constant
122  // triangulation for the user
123  // flags. Since we restore them in
124  // the end, this cast is safe.
125  Triangulation<dim,spacedim> &user_flags_triangulation =
126  const_cast<Triangulation<dim,spacedim>&> (dofs.get_triangulation());
127  user_flags_triangulation.save_user_flags(old_flags);
128  user_flags_triangulation.clear_user_flags();
129 
130  const typename DoFHandler<dim,spacedim>::cell_iterator end = dofs.end(level);
132  std::vector<types::global_dof_index> cell_indices;
133  std::vector<types::global_dof_index> neighbor_indices;
134 
135  // We loop over cells and go from
136  // cells to lower dimensional
137  // objects. This is the only way to
138  // cope with the fact, that an
139  // unknown number of cells may
140  // share an object of dimension
141  // smaller than dim-1.
142  for (cell = dofs.begin(level); cell != end; ++cell)
143  {
144  const FiniteElement<dim> &fe = cell->get_fe();
145  cell_indices.resize(fe.dofs_per_cell);
146  cell->get_mg_dof_indices(cell_indices);
147  unsigned int i = 0;
148  // First, dofs on
149  // vertices. We assume that
150  // each vertex dof couples
151  // with all dofs on
152  // adjacent grid cells.
153 
154  // Adding all dofs of the cells
155  // will add dofs of the faces
156  // of the cell adjacent to the
157  // vertex twice. Therefore, we
158  // subtract these here and add
159  // them in a loop over the
160  // faces below.
161 
162  // in 1d, faces and vertices
163  // are identical. Nevertheless,
164  // this will only work if
165  // dofs_per_face is zero and
166  // dofs_per_vertex is
167  // arbitrary, not the other way
168  // round.
169 //TODO: This assumes that even in hp context, the dofs per face coincide!
170  unsigned int increment = fe.dofs_per_cell - dim * fe.dofs_per_face;
171  while (i < fe.first_line_index)
172  row_lengths[cell_indices[i++]] += increment;
173  // From now on, if an object is
174  // a cell, its dofs only couple
175  // inside the cell. Since the
176  // faces are handled below, we
177  // have to subtract ALL faces
178  // in this case.
179 
180  // In all other cases we
181  // subtract adjacent faces to be
182  // added in the loop below.
183  increment = (dim>1)
184  ? fe.dofs_per_cell - (dim-1) * fe.dofs_per_face
186  while (i < fe.first_quad_index)
187  row_lengths[cell_indices[i++]] += increment;
188 
189  // Now quads in 2D and 3D
190  increment = (dim>2)
191  ? fe.dofs_per_cell - (dim-2) * fe.dofs_per_face
193  while (i < fe.first_hex_index)
194  row_lengths[cell_indices[i++]] += increment;
195  // Finally, cells in 3D
197  while (i < fe.dofs_per_cell)
198  row_lengths[cell_indices[i++]] += increment;
199 
200  // At this point, we have
201  // counted all dofs
202  // contributiong from cells
203  // coupled topologically to the
204  // adjacent cells, but we
205  // subtracted some faces.
206 
207  // Now, let's go by the faces
208  // and add the missing
209  // contribution as well as the
210  // flux contributions.
211  for (unsigned int iface=0; iface<GeometryInfo<dim>::faces_per_cell; ++iface)
212  {
213  bool level_boundary = cell->at_boundary(iface);
214  typename DoFHandler<dim,spacedim>::cell_iterator neighbor;
215  if (!level_boundary)
216  {
217  neighbor = cell->neighbor(iface);
218  if (static_cast<unsigned int>(neighbor->level()) != level)
219  level_boundary = true;
220  }
221 
222  if (level_boundary)
223  {
224  for (unsigned int local_dof=0; local_dof<fe.dofs_per_cell; ++local_dof)
225  row_lengths[cell_indices[local_dof]] += fe.dofs_per_face;
226  continue;
227  }
228 
229  const FiniteElement<dim> &nfe = neighbor->get_fe();
230  typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(iface);
231 
232  // Flux couplings are
233  // computed from both sides
234  // for simplicity.
235 
236  // The dofs on the common face
237  // will be handled below,
238  // therefore, we subtract them
239  // here.
240  if (flux_coupling != DoFTools::none)
241  {
242  const unsigned int dof_increment = nfe.dofs_per_cell - nfe.dofs_per_face;
243  for (unsigned int local_dof=0; local_dof<fe.dofs_per_cell; ++local_dof)
244  row_lengths[cell_indices[local_dof]] += dof_increment;
245  }
246 
247  // Do this only once per
248  // face.
249  if (face->user_flag_set())
250  continue;
251  face->set_user_flag();
252  // At this point, we assume
253  // that each cell added its
254  // dofs minus the face to
255  // the couplings of the
256  // face dofs. Since we
257  // subtracted two faces, we
258  // have to re-add one.
259 
260  // If one side of the face
261  // is refined, all the fine
262  // face dofs couple with
263  // the coarse one.
264  neighbor_indices.resize(nfe.dofs_per_cell);
265  neighbor->get_mg_dof_indices(neighbor_indices);
266  for (unsigned int local_dof=0; local_dof<fe.dofs_per_cell; ++local_dof)
267  row_lengths[cell_indices[local_dof]] += nfe.dofs_per_face;
268  for (unsigned int local_dof=0; local_dof<nfe.dofs_per_cell; ++local_dof)
269  row_lengths[neighbor_indices[local_dof]] += fe.dofs_per_face;
270  }
271  }
272  user_flags_triangulation.load_user_flags(old_flags);
273  }
274 
275 
276 // This is the template for 2D and 3D. See version for 1D above
277  template <int dim, int spacedim>
278  void
280  const DoFHandler<dim,spacedim> &dofs,
281  const unsigned int level,
282  std::vector<unsigned int> &row_lengths,
283  const Table<2,DoFTools::Coupling> &couplings,
284  const Table<2,DoFTools::Coupling> &flux_couplings)
285  {
286  Assert (row_lengths.size() == dofs.n_dofs(),
287  ExcDimensionMismatch(row_lengths.size(), dofs.n_dofs()));
288 
289  // Function starts here by
290  // resetting the counters.
291  std::fill(row_lengths.begin(), row_lengths.end(), 0);
292  // We need the user flags, so we
293  // save them for later restoration
294  std::vector<bool> old_flags;
295  // We need a non-constant
296  // triangulation for the user
297  // flags. Since we restore them in
298  // the end, this cast is safe.
299  Triangulation<dim,spacedim> &user_flags_triangulation =
300  const_cast<Triangulation<dim,spacedim>&> (dofs.get_triangulation());
301  user_flags_triangulation.save_user_flags(old_flags);
302  user_flags_triangulation.clear_user_flags();
303 
304  const typename DoFHandler<dim,spacedim>::cell_iterator end = dofs.end(level);
306  std::vector<types::global_dof_index> cell_indices;
307  std::vector<types::global_dof_index> neighbor_indices;
308 
309  // We have to translate the
310  // couplings from components to
311  // blocks, so this works for
312  // nonprimitive elements as well.
313  std::vector<Table<2, DoFTools::Coupling> > couple_cell;
314  std::vector<Table<2, DoFTools::Coupling> > couple_face;
315  DoFTools::convert_couplings_to_blocks(dofs, couplings, couple_cell);
316  DoFTools::convert_couplings_to_blocks(dofs, flux_couplings, couple_face);
317 
318  // We loop over cells and go from
319  // cells to lower dimensional
320  // objects. This is the only way to
321  // cope withthe fact, that an
322  // unknown number of cells may
323  // share an object of dimension
324  // smaller than dim-1.
325  for (cell = dofs.begin_active(); cell != end; ++cell)
326  {
327  const FiniteElement<dim> &fe = cell->get_fe();
328  const unsigned int fe_index = cell->active_fe_index();
329 
330  Assert (couplings.n_rows()==fe.n_components(),
331  ExcDimensionMismatch(couplings.n_rows(), fe.n_components()));
332  Assert (couplings.n_cols()==fe.n_components(),
333  ExcDimensionMismatch(couplings.n_cols(), fe.n_components()));
334  Assert (flux_couplings.n_rows()==fe.n_components(),
335  ExcDimensionMismatch(flux_couplings.n_rows(), fe.n_components()));
336  Assert (flux_couplings.n_cols()==fe.n_components(),
337  ExcDimensionMismatch(flux_couplings.n_cols(), fe.n_components()));
338 
339  cell_indices.resize(fe.dofs_per_cell);
340  cell->get_mg_dof_indices(cell_indices);
341  unsigned int i = 0;
342  // First, dofs on
343  // vertices. We assume that
344  // each vertex dof couples
345  // with all dofs on
346  // adjacent grid cells.
347 
348  // Adding all dofs of the cells
349  // will add dofs of the faces
350  // of the cell adjacent to the
351  // vertex twice. Therefore, we
352  // subtract these here and add
353  // them in a loop over the
354  // faces below.
355 
356  // in 1d, faces and vertices
357  // are identical. Nevertheless,
358  // this will only work if
359  // dofs_per_face is zero and
360  // dofs_per_vertex is
361  // arbitrary, not the other way
362  // round.
363  unsigned int increment;
364  while (i < fe.first_line_index)
365  {
366  for (unsigned int base=0; base<fe.n_base_elements(); ++base)
367  for (unsigned int mult=0; mult<fe.element_multiplicity(base); ++mult)
368  if (couple_cell[fe_index](fe.system_to_block_index(i).first,
369  fe.first_block_of_base(base) + mult) != DoFTools::none)
370  {
371  increment = fe.base_element(base).dofs_per_cell
372  - dim * fe.base_element(base).dofs_per_face;
373  row_lengths[cell_indices[i]] += increment;
374  }
375  ++i;
376  }
377  // From now on, if an object is
378  // a cell, its dofs only couple
379  // inside the cell. Since the
380  // faces are handled below, we
381  // have to subtract ALL faces
382  // in this case.
383 
384  // In all other cases we
385  // subtract adjacent faces to be
386  // added in the loop below.
387  while (i < fe.first_quad_index)
388  {
389  for (unsigned int base=0; base<fe.n_base_elements(); ++base)
390  for (unsigned int mult=0; mult<fe.element_multiplicity(base); ++mult)
391  if (couple_cell[fe_index](fe.system_to_block_index(i).first,
392  fe.first_block_of_base(base) + mult) != DoFTools::none)
393  {
394  increment = fe.base_element(base).dofs_per_cell
395  - ((dim>1)
396  ? (dim-1)
398  * fe.base_element(base).dofs_per_face;
399  row_lengths[cell_indices[i]] += increment;
400  }
401  ++i;
402  }
403 
404  // Now quads in 2D and 3D
405  while (i < fe.first_hex_index)
406  {
407  for (unsigned int base=0; base<fe.n_base_elements(); ++base)
408  for (unsigned int mult=0; mult<fe.element_multiplicity(base); ++mult)
409  if (couple_cell[fe_index](fe.system_to_block_index(i).first,
410  fe.first_block_of_base(base) + mult) != DoFTools::none)
411  {
412  increment = fe.base_element(base).dofs_per_cell
413  - ((dim>2)
414  ? (dim-2)
416  * fe.base_element(base).dofs_per_face;
417  row_lengths[cell_indices[i]] += increment;
418  }
419  ++i;
420  }
421 
422  // Finally, cells in 3D
423  while (i < fe.dofs_per_cell)
424  {
425  for (unsigned int base=0; base<fe.n_base_elements(); ++base)
426  for (unsigned int mult=0; mult<fe.element_multiplicity(base); ++mult)
427  if (couple_cell[fe_index](fe.system_to_block_index(i).first,
428  fe.first_block_of_base(base) + mult) != DoFTools::none)
429  {
430  increment = fe.base_element(base).dofs_per_cell
432  * fe.base_element(base).dofs_per_face;
433  row_lengths[cell_indices[i]] += increment;
434  }
435  ++i;
436  }
437 
438  // At this point, we have
439  // counted all dofs
440  // contributiong from cells
441  // coupled topologically to the
442  // adjacent cells, but we
443  // subtracted some faces.
444 
445  // Now, let's go by the faces
446  // and add the missing
447  // contribution as well as the
448  // flux contributions.
449  for (unsigned int iface=0; iface<GeometryInfo<dim>::faces_per_cell; ++iface)
450  {
451  bool level_boundary = cell->at_boundary(iface);
452  typename DoFHandler<dim,spacedim>::cell_iterator neighbor;
453  if (!level_boundary)
454  {
455  neighbor = cell->neighbor(iface);
456  if (static_cast<unsigned int>(neighbor->level()) != level)
457  level_boundary = true;
458  }
459 
460  if (level_boundary)
461  {
462  for (unsigned int local_dof=0; local_dof<fe.dofs_per_cell; ++local_dof)
463  row_lengths[cell_indices[local_dof]] += fe.dofs_per_face;
464  continue;
465  }
466 
467  const FiniteElement<dim> &nfe = neighbor->get_fe();
468  typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(iface);
469 
470  // Flux couplings are
471  // computed from both sides
472  // for simplicity.
473 
474  // The dofs on the common face
475  // will be handled below,
476  // therefore, we subtract them
477  // here.
478  for (unsigned int base=0; base<nfe.n_base_elements(); ++base)
479  for (unsigned int mult=0; mult<nfe.element_multiplicity(base); ++mult)
480  for (unsigned int local_dof=0; local_dof<fe.dofs_per_cell; ++local_dof)
481  if (couple_face[fe_index](fe.system_to_block_index(local_dof).first,
482  nfe.first_block_of_base(base) + mult) != DoFTools::none)
483  {
484  const unsigned int dof_increment = nfe.base_element(base).dofs_per_cell
485  - nfe.base_element(base).dofs_per_face;
486  row_lengths[cell_indices[local_dof]] += dof_increment;
487  }
488 
489  // Do this only once per
490  // face and not on the
491  // hanging faces.
492  if (face->user_flag_set())
493  continue;
494  face->set_user_flag();
495  // At this point, we assume
496  // that each cell added its
497  // dofs minus the face to
498  // the couplings of the
499  // face dofs. Since we
500  // subtracted two faces, we
501  // have to re-add one.
502 
503  // If one side of the face
504  // is refined, all the fine
505  // face dofs couple with
506  // the coarse one.
507 
508  // Wolfgang, do they couple
509  // with each other by
510  // constraints?
511 
512  // This will not work with
513  // different couplings on
514  // different cells.
515  neighbor_indices.resize(nfe.dofs_per_cell);
516  neighbor->get_mg_dof_indices(neighbor_indices);
517  for (unsigned int base=0; base<nfe.n_base_elements(); ++base)
518  for (unsigned int mult=0; mult<nfe.element_multiplicity(base); ++mult)
519  for (unsigned int local_dof=0; local_dof<fe.dofs_per_cell; ++local_dof)
520  if (couple_cell[fe_index](fe.system_to_component_index(local_dof).first,
521  nfe.first_block_of_base(base) + mult) != DoFTools::none)
522  row_lengths[cell_indices[local_dof]]
523  += nfe.base_element(base).dofs_per_face;
524  for (unsigned int base=0; base<fe.n_base_elements(); ++base)
525  for (unsigned int mult=0; mult<fe.element_multiplicity(base); ++mult)
526  for (unsigned int local_dof=0; local_dof<nfe.dofs_per_cell; ++local_dof)
527  if (couple_cell[fe_index](nfe.system_to_component_index(local_dof).first,
528  fe.first_block_of_base(base) + mult) != DoFTools::none)
529  row_lengths[neighbor_indices[local_dof]]
530  += fe.base_element(base).dofs_per_face;
531  }
532  }
533  user_flags_triangulation.load_user_flags(old_flags);
534  }
535 
536 
537 
538  template <typename DoFHandlerType, typename SparsityPatternType>
539  void make_sparsity_pattern (const DoFHandlerType &dof,
540  SparsityPatternType &sparsity,
541  const unsigned int level)
542  {
543  const types::global_dof_index n_dofs = dof.n_dofs(level);
544  (void)n_dofs;
545 
546  Assert (sparsity.n_rows() == n_dofs,
547  ExcDimensionMismatch (sparsity.n_rows(), n_dofs));
548  Assert (sparsity.n_cols() == n_dofs,
549  ExcDimensionMismatch (sparsity.n_cols(), n_dofs));
550 
551  const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
552  std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
553  typename DoFHandlerType::cell_iterator cell = dof.begin(level),
554  endc = dof.end(level);
555  for (; cell!=endc; ++cell)
556  if (dof.get_triangulation().locally_owned_subdomain()==numbers::invalid_subdomain_id
557  || cell->level_subdomain_id()==dof.get_triangulation().locally_owned_subdomain())
558  {
559  cell->get_mg_dof_indices (dofs_on_this_cell);
560  // make sparsity pattern for this cell
561  for (unsigned int i=0; i<dofs_per_cell; ++i)
562  for (unsigned int j=0; j<dofs_per_cell; ++j)
563  sparsity.add (dofs_on_this_cell[i],
564  dofs_on_this_cell[j]);
565  }
566  }
567 
568 
569 
570  template <int dim, typename SparsityPatternType, int spacedim>
571  void
573  SparsityPatternType &sparsity,
574  const unsigned int level)
575  {
576  const types::global_dof_index n_dofs = dof.n_dofs(level);
577  (void)n_dofs;
578 
579  Assert (sparsity.n_rows() == n_dofs,
580  ExcDimensionMismatch (sparsity.n_rows(), n_dofs));
581  Assert (sparsity.n_cols() == n_dofs,
582  ExcDimensionMismatch (sparsity.n_cols(), n_dofs));
583 
584  const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
585  std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
586  std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
587  typename DoFHandler<dim,spacedim>::cell_iterator cell = dof.begin(level),
588  endc = dof.end(level);
589  for (; cell!=endc; ++cell)
590  {
591  if (!cell->is_locally_owned_on_level()) continue;
592 
593  cell->get_mg_dof_indices (dofs_on_this_cell);
594  // make sparsity pattern for this cell
595  for (unsigned int i=0; i<dofs_per_cell; ++i)
596  for (unsigned int j=0; j<dofs_per_cell; ++j)
597  sparsity.add (dofs_on_this_cell[i],
598  dofs_on_this_cell[j]);
599 
600  // Loop over all interior neighbors
601  for (unsigned int face = 0;
602  face < GeometryInfo<dim>::faces_per_cell;
603  ++face)
604  {
605  bool use_face = false;
606  if ((! cell->at_boundary(face)) &&
607  (static_cast<unsigned int>(cell->neighbor_level(face)) == level))
608  use_face = true;
609  else if (cell->has_periodic_neighbor(face) &&
610  (static_cast<unsigned int>(cell->periodic_neighbor_level(face)) == level))
611  use_face = true;
612 
613  if (use_face)
614  {
616  neighbor = cell->neighbor_or_periodic_neighbor(face);
617  neighbor->get_mg_dof_indices (dofs_on_other_cell);
618  // only add one direction The other is taken care of by
619  // neighbor (except when the neighbor is not owned by the same
620  // processor)
621  for (unsigned int i=0; i<dofs_per_cell; ++i)
622  {
623  for (unsigned int j=0; j<dofs_per_cell; ++j)
624  {
625  sparsity.add (dofs_on_this_cell[i],
626  dofs_on_other_cell[j]);
627  }
628  }
629  if (neighbor->is_locally_owned_on_level() == false)
630  for (unsigned int i=0; i<dofs_per_cell; ++i)
631  for (unsigned int j=0; j<dofs_per_cell; ++j)
632  {
633  sparsity.add (dofs_on_other_cell[i],
634  dofs_on_other_cell[j]);
635  sparsity.add (dofs_on_other_cell[i],
636  dofs_on_this_cell[j]);
637  }
638  }
639  }
640  }
641  }
642 
643 
644 
645  template <int dim, typename SparsityPatternType, int spacedim>
646  void
648  SparsityPatternType &sparsity,
649  const unsigned int level)
650  {
651  Assert ((level>=1) && (level<dof.get_triangulation().n_global_levels()),
652  ExcIndexRange(level, 1, dof.get_triangulation().n_global_levels()));
653 
654  const types::global_dof_index fine_dofs = dof.n_dofs(level);
655  const types::global_dof_index coarse_dofs = dof.n_dofs(level-1);
656  (void)fine_dofs;
657  (void)coarse_dofs;
658 
659  // Matrix maps from fine level to coarse level
660 
661  Assert (sparsity.n_rows() == coarse_dofs,
662  ExcDimensionMismatch (sparsity.n_rows(), coarse_dofs));
663  Assert (sparsity.n_cols() == fine_dofs,
664  ExcDimensionMismatch (sparsity.n_cols(), fine_dofs));
665 
666  const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
667  std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
668  std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
669  typename DoFHandler<dim,spacedim>::cell_iterator cell = dof.begin(level),
670  endc = dof.end(level);
671  for (; cell!=endc; ++cell)
672  {
673  if (!cell->is_locally_owned_on_level()) continue;
674 
675  cell->get_mg_dof_indices (dofs_on_this_cell);
676  // Loop over all interior neighbors
677  for (unsigned int face = 0;
678  face < GeometryInfo<dim>::faces_per_cell;
679  ++face)
680  {
681  // Neighbor is coarser
682  bool use_face = false;
683  if ((! cell->at_boundary(face)) &&
684  (static_cast<unsigned int>(cell->neighbor_level(face)) != level))
685  use_face = true;
686  else if (cell->has_periodic_neighbor(face) &&
687  (static_cast<unsigned int>(cell->periodic_neighbor_level(face)) != level))
688  use_face = true;
689 
690  if (use_face)
691  {
693  neighbor = cell->neighbor_or_periodic_neighbor(face);
694  neighbor->get_mg_dof_indices (dofs_on_other_cell);
695 
696  for (unsigned int i=0; i<dofs_per_cell; ++i)
697  {
698  for (unsigned int j=0; j<dofs_per_cell; ++j)
699  {
700  sparsity.add (dofs_on_other_cell[i],
701  dofs_on_this_cell[j]);
702  sparsity.add (dofs_on_other_cell[j],
703  dofs_on_this_cell[i]);
704  }
705  }
706  }
707  }
708  }
709  }
710 
711 
712 
713  template <int dim, typename SparsityPatternType, int spacedim>
714  void
716  SparsityPatternType &sparsity,
717  const unsigned int level,
718  const Table<2,DoFTools::Coupling> &int_mask,
719  const Table<2,DoFTools::Coupling> &flux_mask)
720  {
721  const FiniteElement<dim> &fe = dof.get_fe();
722  const types::global_dof_index n_dofs = dof.n_dofs(level);
723  const unsigned int n_comp = fe.n_components();
724  (void)n_dofs;
725  (void)n_comp;
726 
727  Assert (sparsity.n_rows() == n_dofs,
728  ExcDimensionMismatch (sparsity.n_rows(), n_dofs));
729  Assert (sparsity.n_cols() == n_dofs,
730  ExcDimensionMismatch (sparsity.n_cols(), n_dofs));
731  Assert (int_mask.n_rows() == n_comp,
732  ExcDimensionMismatch (int_mask.n_rows(), n_comp));
733  Assert (int_mask.n_cols() == n_comp,
734  ExcDimensionMismatch (int_mask.n_cols(), n_comp));
735  Assert (flux_mask.n_rows() == n_comp,
736  ExcDimensionMismatch (flux_mask.n_rows(), n_comp));
737  Assert (flux_mask.n_cols() == n_comp,
738  ExcDimensionMismatch (flux_mask.n_cols(), n_comp));
739 
740  const unsigned int total_dofs = fe.dofs_per_cell;
741  std::vector<types::global_dof_index> dofs_on_this_cell(total_dofs);
742  std::vector<types::global_dof_index> dofs_on_other_cell(total_dofs);
743  Table<2,bool> support_on_face(total_dofs, GeometryInfo<dim>::faces_per_cell);
744 
745  typename DoFHandler<dim,spacedim>::cell_iterator cell = dof.begin(level),
746  endc = dof.end(level);
747 
749  int_dof_mask = DoFTools::dof_couplings_from_component_couplings(fe, int_mask),
750  flux_dof_mask = DoFTools::dof_couplings_from_component_couplings(fe, flux_mask);
751 
752  for (unsigned int i=0; i<total_dofs; ++i)
753  for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
754  support_on_face(i,f) = fe.has_support_on_face(i,f);
755 
756  // Clear user flags because we will
757  // need them. But first we save
758  // them and make sure that we
759  // restore them later such that at
760  // the end of this function the
761  // Triangulation will be in the
762  // same state as it was at the
763  // beginning of this function.
764  std::vector<bool> user_flags;
765  dof.get_triangulation().save_user_flags(user_flags);
766  const_cast<Triangulation<dim,spacedim> &>(dof.get_triangulation()).clear_user_flags ();
767 
768  for (; cell!=endc; ++cell)
769  {
770  if (!cell->is_locally_owned_on_level()) continue;
771 
772  cell->get_mg_dof_indices (dofs_on_this_cell);
773  // make sparsity pattern for this cell
774  for (unsigned int i=0; i<total_dofs; ++i)
775  for (unsigned int j=0; j<total_dofs; ++j)
776  if (int_dof_mask[i][j] != DoFTools::none)
777  sparsity.add (dofs_on_this_cell[i],
778  dofs_on_this_cell[j]);
779 
780  // Loop over all interior neighbors
781  for (unsigned int face = 0;
782  face < GeometryInfo<dim>::faces_per_cell;
783  ++face)
784  {
785  typename DoFHandler<dim,spacedim>::face_iterator cell_face = cell->face(face);
786  if (cell_face->user_flag_set ())
787  continue;
788 
789  if (cell->at_boundary (face) && !cell->has_periodic_neighbor(face))
790  {
791  for (unsigned int i=0; i<total_dofs; ++i)
792  {
793  const bool i_non_zero_i = support_on_face (i, face);
794  for (unsigned int j=0; j<total_dofs; ++j)
795  {
796  const bool j_non_zero_i = support_on_face (j, face);
797 
798  if (flux_dof_mask(i,j) == DoFTools::always)
799  sparsity.add (dofs_on_this_cell[i],
800  dofs_on_this_cell[j]);
801  if (flux_dof_mask(i,j) == DoFTools::nonzero
802  && i_non_zero_i && j_non_zero_i)
803  sparsity.add (dofs_on_this_cell[i],
804  dofs_on_this_cell[j]);
805  }
806  }
807  }
808  else
809  {
811  neighbor = cell->neighbor_or_periodic_neighbor(face);
812 
813  if (neighbor->level() < cell->level())
814  continue;
815 
816  unsigned int neighbor_face = cell->has_periodic_neighbor(face)?
817  cell->periodic_neighbor_of_periodic_neighbor(face):
818  cell->neighbor_of_neighbor(face);
819 
820  neighbor->get_mg_dof_indices (dofs_on_other_cell);
821  for (unsigned int i=0; i<total_dofs; ++i)
822  {
823  const bool i_non_zero_i = support_on_face (i, face);
824  const bool i_non_zero_e = support_on_face (i, neighbor_face);
825  for (unsigned int j=0; j<total_dofs; ++j)
826  {
827  const bool j_non_zero_i = support_on_face (j, face);
828  const bool j_non_zero_e = support_on_face (j, neighbor_face);
829  if (flux_dof_mask(i,j) == DoFTools::always)
830  {
831  sparsity.add (dofs_on_this_cell[i],
832  dofs_on_other_cell[j]);
833  sparsity.add (dofs_on_other_cell[i],
834  dofs_on_this_cell[j]);
835  sparsity.add (dofs_on_this_cell[i],
836  dofs_on_this_cell[j]);
837  sparsity.add (dofs_on_other_cell[i],
838  dofs_on_other_cell[j]);
839  }
840  if (flux_dof_mask(i,j) == DoFTools::nonzero)
841  {
842  if (i_non_zero_i && j_non_zero_e)
843  sparsity.add (dofs_on_this_cell[i],
844  dofs_on_other_cell[j]);
845  if (i_non_zero_e && j_non_zero_i)
846  sparsity.add (dofs_on_other_cell[i],
847  dofs_on_this_cell[j]);
848  if (i_non_zero_i && j_non_zero_i)
849  sparsity.add (dofs_on_this_cell[i],
850  dofs_on_this_cell[j]);
851  if (i_non_zero_e && j_non_zero_e)
852  sparsity.add (dofs_on_other_cell[i],
853  dofs_on_other_cell[j]);
854  }
855 
856  if (flux_dof_mask(j,i) == DoFTools::always)
857  {
858  sparsity.add (dofs_on_this_cell[j],
859  dofs_on_other_cell[i]);
860  sparsity.add (dofs_on_other_cell[j],
861  dofs_on_this_cell[i]);
862  sparsity.add (dofs_on_this_cell[j],
863  dofs_on_this_cell[i]);
864  sparsity.add (dofs_on_other_cell[j],
865  dofs_on_other_cell[i]);
866  }
867  if (flux_dof_mask(j,i) == DoFTools::nonzero)
868  {
869  if (j_non_zero_i && i_non_zero_e)
870  sparsity.add (dofs_on_this_cell[j],
871  dofs_on_other_cell[i]);
872  if (j_non_zero_e && i_non_zero_i)
873  sparsity.add (dofs_on_other_cell[j],
874  dofs_on_this_cell[i]);
875  if (j_non_zero_i && i_non_zero_i)
876  sparsity.add (dofs_on_this_cell[j],
877  dofs_on_this_cell[i]);
878  if (j_non_zero_e && i_non_zero_e)
879  sparsity.add (dofs_on_other_cell[j],
880  dofs_on_other_cell[i]);
881  }
882  }
883  }
884  neighbor->face(neighbor_face)->set_user_flag ();
885  }
886  }
887  }
888 
889  // finally restore the user flags
890  const_cast<Triangulation<dim,spacedim> &>(dof.get_triangulation()).load_user_flags(user_flags);
891  }
892 
893 
894 
895  template <int dim, typename SparsityPatternType, int spacedim>
896  void
898  SparsityPatternType &sparsity,
899  const unsigned int level,
900  const Table<2,DoFTools::Coupling> &flux_mask)
901  {
902  const FiniteElement<dim> &fe = dof.get_fe();
903  const unsigned int n_comp = fe.n_components();
904  (void)n_comp;
905 
906  Assert ((level>=1) && (level<dof.get_triangulation().n_global_levels()),
907  ExcIndexRange(level, 1, dof.get_triangulation().n_global_levels()));
908 
909  const types::global_dof_index fine_dofs = dof.n_dofs(level);
910  const types::global_dof_index coarse_dofs = dof.n_dofs(level-1);
911  (void)fine_dofs;
912  (void)coarse_dofs;
913 
914  // Matrix maps from fine level to coarse level
915 
916  Assert (sparsity.n_rows() == coarse_dofs,
917  ExcDimensionMismatch (sparsity.n_rows(), coarse_dofs));
918  Assert (sparsity.n_cols() == fine_dofs,
919  ExcDimensionMismatch (sparsity.n_cols(), fine_dofs));
920  Assert (flux_mask.n_rows() == n_comp,
921  ExcDimensionMismatch (flux_mask.n_rows(), n_comp));
922  Assert (flux_mask.n_cols() == n_comp,
923  ExcDimensionMismatch (flux_mask.n_cols(), n_comp));
924 
925  const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
926  std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
927  std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
928  Table<2,bool> support_on_face(dofs_per_cell, GeometryInfo<dim>::faces_per_cell);
929 
930  typename DoFHandler<dim,spacedim>::cell_iterator cell = dof.begin(level),
931  endc = dof.end(level);
932 
933  const Table<2,DoFTools::Coupling> flux_dof_mask
935 
936  for (unsigned int i=0; i<dofs_per_cell; ++i)
937  for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
938  support_on_face(i,f) = fe.has_support_on_face(i,f);
939 
940  for (; cell!=endc; ++cell)
941  {
942  if (!cell->is_locally_owned_on_level()) continue;
943 
944  cell->get_mg_dof_indices (dofs_on_this_cell);
945  // Loop over all interior neighbors
946  for (unsigned int face = 0;
947  face < GeometryInfo<dim>::faces_per_cell;
948  ++face)
949  {
950  // Neighbor is coarser
951  bool use_face = false;
952  if ((! cell->at_boundary(face)) &&
953  (static_cast<unsigned int>(cell->neighbor_level(face)) != level))
954  use_face = true;
955  else if (cell->has_periodic_neighbor(face) &&
956  (static_cast<unsigned int>(cell->periodic_neighbor_level(face)) != level))
957  use_face = true;
958 
959  if (use_face)
960  {
962  neighbor = cell->neighbor_or_periodic_neighbor(face);
963  neighbor->get_mg_dof_indices (dofs_on_other_cell);
964 
965  for (unsigned int i=0; i<dofs_per_cell; ++i)
966  {
967  for (unsigned int j=0; j<dofs_per_cell; ++j)
968  {
969  if (flux_dof_mask(i,j) != DoFTools::none)
970  {
971  sparsity.add (dofs_on_other_cell[i],
972  dofs_on_this_cell[j]);
973  sparsity.add (dofs_on_other_cell[j],
974  dofs_on_this_cell[i]);
975  }
976  }
977  }
978  }
979  }
980  }
981  }
982 
983 
984 
985  template <int dim, int spacedim>
986  void
988  std::vector<std::vector<types::global_dof_index> > &result,
989  bool only_once,
990  std::vector<unsigned int> target_component)
991  {
992  const FiniteElement<dim> &fe = dof_handler.get_fe();
993  const unsigned int n_components = fe.n_components();
994  const unsigned int nlevels = dof_handler.get_triangulation().n_global_levels();
995 
996  Assert (result.size() == nlevels,
997  ExcDimensionMismatch(result.size(), nlevels));
998 
999  if (target_component.size() == 0)
1000  {
1001  target_component.resize(n_components);
1002  for (unsigned int i=0; i<n_components; ++i)
1003  target_component[i] = i;
1004  }
1005 
1006  Assert(target_component.size() == n_components,
1007  ExcDimensionMismatch(target_component.size(), n_components));
1008 
1009  for (unsigned int l=0; l<nlevels; ++l)
1010  {
1011  result[l].resize (n_components);
1012  std::fill (result[l].begin(),result[l].end(), 0U);
1013 
1014  // special case for only one
1015  // component. treat this first
1016  // since it does not require any
1017  // computations
1018  if (n_components == 1)
1019  {
1020  result[l][0] = dof_handler.n_dofs(l);
1021  }
1022  else
1023  {
1024  // otherwise determine the number
1025  // of dofs in each component
1026  // separately. do so in parallel
1027  std::vector<std::vector<bool> >
1028  dofs_in_component (n_components,
1029  std::vector<bool>(dof_handler.n_dofs(l),
1030  false));
1031  std::vector<ComponentMask> component_select (n_components);
1032  Threads::TaskGroup<> tasks;
1033  for (unsigned int i=0; i<n_components; ++i)
1034  {
1035  void (*fun_ptr) (const unsigned int level,
1036  const DoFHandler<dim,spacedim> &,
1037  const ComponentMask &,
1038  std::vector<bool> &)
1040 
1041  std::vector<bool> tmp(n_components, false);
1042  tmp[i] = true;
1043  component_select[i] = ComponentMask(tmp);
1044 
1045  tasks += Threads::new_task (fun_ptr,
1046  l, dof_handler,
1047  component_select[i],
1048  dofs_in_component[i]);
1049  }
1050  tasks.join_all();
1051 
1052  // next count what we got
1053  unsigned int component = 0;
1054  for (unsigned int b=0; b<fe.n_base_elements(); ++b)
1055  {
1056  const FiniteElement<dim> &base = fe.base_element(b);
1057  // Dimension of base element
1058  unsigned int d = base.n_components();
1059 
1060  for (unsigned int m=0; m<fe.element_multiplicity(b); ++m)
1061  {
1062  for (unsigned int dd=0; dd<d; ++dd)
1063  {
1064  if (base.is_primitive() || (!only_once || dd==0))
1065  result[l][target_component[component]]
1066  += std::count(dofs_in_component[component].begin(),
1067  dofs_in_component[component].end(),
1068  true);
1069  ++component;
1070  }
1071  }
1072  }
1073  // finally sanity check
1074  Assert (!dof_handler.get_fe().is_primitive()
1075  ||
1076  std::accumulate (result[l].begin(),
1077  result[l].end(), 0U)
1078  ==
1079  dof_handler.n_dofs(l),
1080  ExcInternalError());
1081  }
1082  }
1083  }
1084 
1085 
1086 
1087  template <typename DoFHandlerType>
1088  void
1090  (const DoFHandlerType &dof_handler,
1091  std::vector<std::vector<types::global_dof_index> > &dofs_per_block,
1092  std::vector<unsigned int> target_block)
1093  {
1095  const unsigned int n_blocks = fe.n_blocks();
1096  const unsigned int n_levels = dof_handler.get_triangulation().n_global_levels();
1097 
1098  AssertDimension (dofs_per_block.size(), n_levels);
1099 
1100  for (unsigned int l=0; l<n_levels; ++l)
1101  std::fill (dofs_per_block[l].begin(), dofs_per_block[l].end(), 0U);
1102  // If the empty vector was given as
1103  // default argument, set up this
1104  // vector as identity.
1105  if (target_block.size()==0)
1106  {
1107  target_block.resize(n_blocks);
1108  for (unsigned int i=0; i<n_blocks; ++i)
1109  target_block[i] = i;
1110  }
1111  Assert(target_block.size()==n_blocks,
1112  ExcDimensionMismatch(target_block.size(),n_blocks));
1113 
1114  const unsigned int max_block
1115  = *std::max_element (target_block.begin(),
1116  target_block.end());
1117  const unsigned int n_target_blocks = max_block + 1;
1118  (void)n_target_blocks;
1119 
1120  for (unsigned int l=0; l<n_levels; ++l)
1121  AssertDimension (dofs_per_block[l].size(), n_target_blocks);
1122 
1123  // special case for only one
1124  // block. treat this first
1125  // since it does not require any
1126  // computations
1127  if (n_blocks == 1)
1128  {
1129  for (unsigned int l=0; l<n_levels; ++l)
1130  dofs_per_block[l][0] = dof_handler.n_dofs(l);
1131  return;
1132  }
1133  // otherwise determine the number
1134  // of dofs in each block
1135  // separately. do so in parallel
1136  for (unsigned int l=0; l<n_levels; ++l)
1137  {
1138  std::vector<std::vector<bool> >
1139  dofs_in_block (n_blocks, std::vector<bool>(dof_handler.n_dofs(l), false));
1140  std::vector<BlockMask> block_select (n_blocks);
1141  Threads::TaskGroup<> tasks;
1142  for (unsigned int i=0; i<n_blocks; ++i)
1143  {
1144  void (*fun_ptr) (const unsigned int level,
1145  const DoFHandlerType &,
1146  const BlockMask &,
1147  std::vector<bool> &)
1148  = &DoFTools::extract_level_dofs<DoFHandlerType>;
1149 
1150  std::vector<bool> tmp(n_blocks, false);
1151  tmp[i] = true;
1152  block_select[i] = tmp;
1153 
1154  tasks += Threads::new_task (fun_ptr,
1155  l, dof_handler, block_select[i],
1156  dofs_in_block[i]);
1157  }
1158  tasks.join_all ();
1159 
1160  // next count what we got
1161  for (unsigned int block=0; block<fe.n_blocks(); ++block)
1162  dofs_per_block[l][target_block[block]]
1163  += std::count(dofs_in_block[block].begin(),
1164  dofs_in_block[block].end(),
1165  true);
1166  }
1167  }
1168 
1169 
1170 
1171  template <int dim, int spacedim>
1172  void
1174  const DoFHandler<dim,spacedim> &dof,
1175  const typename FunctionMap<dim>::type &function_map,
1176  std::vector<std::set<types::global_dof_index> > &boundary_indices,
1177  const ComponentMask &component_mask)
1178  {
1179  Assert (boundary_indices.size() == dof.get_triangulation().n_global_levels(),
1180  ExcDimensionMismatch (boundary_indices.size(),
1181  dof.get_triangulation().n_global_levels()));
1182 
1183  std::set<types::boundary_id> boundary_ids;
1184  typename std::map<types::boundary_id, const Function<dim>* >::const_iterator it
1185  = function_map.begin();
1186  for (; it!=function_map.end(); ++it)
1187  boundary_ids.insert(it->first);
1188 
1189  std::vector<IndexSet> boundary_indexset;
1190  make_boundary_list(dof, boundary_ids, boundary_indexset, component_mask);
1191  for (unsigned int i=0; i<dof.get_triangulation().n_global_levels(); ++i)
1192  boundary_indices[i].insert(boundary_indexset[i].begin(),
1193  boundary_indexset[i].end());
1194  }
1195 
1196 
1197  template <int dim, int spacedim>
1198  void
1200  const typename FunctionMap<dim>::type &function_map,
1201  std::vector<IndexSet> &boundary_indices,
1202  const ComponentMask &component_mask)
1203  {
1204  Assert (boundary_indices.size() == dof.get_triangulation().n_global_levels(),
1205  ExcDimensionMismatch (boundary_indices.size(),
1206  dof.get_triangulation().n_global_levels()));
1207 
1208  std::set<types::boundary_id> boundary_ids;
1209  typename std::map<types::boundary_id, const Function<dim>* >::const_iterator it = function_map.begin();
1210  for (; it!=function_map.end(); ++it)
1211  boundary_ids.insert(it->first);
1212 
1213  make_boundary_list (dof, boundary_ids, boundary_indices, component_mask);
1214  }
1215 
1216 
1217 
1218  template <int dim, int spacedim>
1219  void
1221  const std::set<types::boundary_id> &boundary_ids,
1222  std::vector<IndexSet> &boundary_indices,
1223  const ComponentMask &component_mask)
1224  {
1225  boundary_indices.resize(dof.get_triangulation().n_global_levels());
1226 
1227  // if for whatever reason we were passed an empty set, return immediately
1228  if (boundary_ids.size() == 0)
1229  return;
1230 
1231  for (unsigned int i=0; i<dof.get_triangulation().n_global_levels(); ++i)
1232  if (boundary_indices[i].size() == 0)
1233  boundary_indices[i] = IndexSet (dof.n_dofs(i));
1234 
1235  const unsigned int n_components = DoFTools::n_components(dof);
1236  const bool fe_is_system = (n_components != 1);
1237 
1238  std::vector<types::global_dof_index> local_dofs;
1239  local_dofs.reserve (DoFTools::max_dofs_per_face(dof));
1240  std::fill (local_dofs.begin (),
1241  local_dofs.end (),
1243 
1244  // First, deal with the simpler case when we have to identify all boundary dofs
1245  if (component_mask.n_selected_components(n_components) == n_components)
1246  {
1248  cell = dof.begin(),
1249  endc = dof.end();
1250  for (; cell!=endc; ++cell)
1251  {
1252  if (dof.get_triangulation().locally_owned_subdomain()!=numbers::invalid_subdomain_id
1253  && cell->level_subdomain_id()==numbers::artificial_subdomain_id)
1254  continue;
1255  const FiniteElement<dim> &fe = cell->get_fe();
1256  const unsigned int level = cell->level();
1257  local_dofs.resize(fe.dofs_per_face);
1258 
1259  for (unsigned int face_no = 0; face_no < GeometryInfo<dim>::faces_per_cell;
1260  ++face_no)
1261  if (cell->at_boundary(face_no) == true)
1262  {
1263  const typename DoFHandler<dim,spacedim>::face_iterator
1264  face = cell->face(face_no);
1265  const types::boundary_id bi = face->boundary_id();
1266  // Face is listed in boundary map
1267  if (boundary_ids.find(bi) != boundary_ids.end())
1268  {
1269  face->get_mg_dof_indices(level, local_dofs);
1270  boundary_indices[level].add_indices(local_dofs.begin(), local_dofs.end());
1271  }
1272  }
1273  }
1274  }
1275  else
1276  {
1277  Assert (component_mask.n_selected_components(n_components) > 0,
1278  ExcMessage("It's probably worthwhile to select at least one component."));
1279 
1281  cell = dof.begin(),
1282  endc = dof.end();
1283  for (; cell!=endc; ++cell)
1284  if (dof.get_triangulation().locally_owned_subdomain()==numbers::invalid_subdomain_id
1285  || cell->level_subdomain_id()!=numbers::artificial_subdomain_id)
1286  for (unsigned int face_no = 0; face_no < GeometryInfo<dim>::faces_per_cell;
1287  ++face_no)
1288  {
1289  if (cell->at_boundary(face_no) == false)
1290  continue;
1291 
1292  const FiniteElement<dim> &fe = cell->get_fe();
1293  const unsigned int level = cell->level();
1294 
1295  typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(face_no);
1296  const types::boundary_id boundary_component = face->boundary_id();
1297  if (boundary_ids.find(boundary_component) != boundary_ids.end())
1298  // we want to constrain this boundary
1299  {
1300 
1301  for (unsigned int i=0; i<cell->get_fe().dofs_per_cell; ++i)
1302  {
1303  const ComponentMask &nonzero_component_array
1304  = cell->get_fe().get_nonzero_components (i);
1305  // if we want to constrain one of the nonzero components,
1306  // we have to constrain all of them
1307 
1308  bool selected = false;
1309  for (unsigned int c=0; c<n_components; ++c)
1310  if (nonzero_component_array[c] == true
1311  && component_mask[c]== true)
1312  {
1313  selected = true;
1314  break;
1315  }
1316  if (selected)
1317  for (unsigned int c=0; c<n_components; ++c)
1318  Assert (nonzero_component_array[c] == false || component_mask[c] == true,
1319  ExcMessage ("You are using a non-primitive FiniteElement "
1320  "and try to constrain just some of its components!"));
1321  }
1322 
1323  // get indices, physical location and boundary values of dofs on this face
1324  local_dofs.resize (fe.dofs_per_face);
1325  face->get_mg_dof_indices (level, local_dofs);
1326  if (fe_is_system)
1327  {
1328  for (unsigned int i=0; i<local_dofs.size(); ++i)
1329  {
1330  unsigned int component = numbers::invalid_unsigned_int;
1331  if (fe.is_primitive())
1332  component = fe.face_system_to_component_index(i).first;
1333  else
1334  {
1335  // Just pick the first of the components
1336  // We already know that either all or none
1337  // of the components are selected
1338  const ComponentMask &nonzero_component_array
1339  = cell->get_fe().get_nonzero_components (i);
1340  for (unsigned int c=0; c<n_components; ++c)
1341  if (nonzero_component_array[c] == true)
1342  {
1343  component = c;
1344  break;
1345  }
1346  }
1348  if (component_mask[component] == true)
1349  boundary_indices[level].add_index(local_dofs[i]);
1350  }
1351  }
1352  else
1353  boundary_indices[level].add_indices(local_dofs.begin(), local_dofs.end());
1354  }
1355  }
1356  }
1357  }
1358 
1359 
1360  template <int dim, int spacedim>
1361  void
1362  extract_non_interface_dofs (const DoFHandler<dim,spacedim> &mg_dof_handler,
1363  std::vector<std::set<types::global_dof_index> > &non_interface_dofs)
1364  {
1365  Assert (non_interface_dofs.size() == mg_dof_handler.get_triangulation().n_global_levels(),
1366  ExcDimensionMismatch (non_interface_dofs.size(),
1367  mg_dof_handler.get_triangulation().n_global_levels()));
1368 
1369  const FiniteElement<dim,spacedim> &fe = mg_dof_handler.get_fe();
1370 
1371  const unsigned int dofs_per_cell = fe.dofs_per_cell;
1372  const unsigned int dofs_per_face = fe.dofs_per_face;
1373 
1374  std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
1375  std::vector<bool> cell_dofs(dofs_per_cell, false);
1376  std::vector<bool> cell_dofs_interface(dofs_per_cell, false);
1377 
1378  typename DoFHandler<dim>::cell_iterator cell = mg_dof_handler.begin(),
1379  endc = mg_dof_handler.end();
1380 
1381 
1382  for (; cell!=endc; ++cell)
1383  {
1384  if (mg_dof_handler.get_triangulation().locally_owned_subdomain()!=numbers::invalid_subdomain_id
1385  && cell->level_subdomain_id()!=mg_dof_handler.get_triangulation().locally_owned_subdomain())
1386  continue;
1387 
1388  std::fill (cell_dofs.begin(), cell_dofs.end(), false);
1389  std::fill (cell_dofs_interface.begin(), cell_dofs_interface.end(), false);
1390 
1391  for (unsigned int face_nr=0; face_nr<GeometryInfo<dim>::faces_per_cell; ++face_nr)
1392  {
1393  const typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(face_nr);
1394  if (!face->at_boundary())
1395  {
1396  //interior face
1397  const typename DoFHandler<dim>::cell_iterator
1398  neighbor = cell->neighbor(face_nr);
1399 
1400  if ((neighbor->level() < cell->level()))
1401  {
1402  for (unsigned int j=0; j<dofs_per_face; ++j)
1403  cell_dofs_interface[fe.face_to_cell_index(j,face_nr)] = true;
1404  }
1405  else
1406  {
1407  for (unsigned int j=0; j<dofs_per_face; ++j)
1408  cell_dofs[fe.face_to_cell_index(j,face_nr)] = true;
1409  }
1410  }
1411  else
1412  {
1413  //boundary face
1414  for (unsigned int j=0; j<dofs_per_face; ++j)
1415  cell_dofs[fe.face_to_cell_index(j,face_nr)] = true;
1416  }
1417  }
1418 
1419  const unsigned int level = cell->level();
1420  cell->get_mg_dof_indices (local_dof_indices);
1421 
1422  for (unsigned int i=0; i<dofs_per_cell; ++i)
1423  if (cell_dofs[i] && !cell_dofs_interface[i])
1424  non_interface_dofs[level].insert(local_dof_indices[i]);
1425  }
1426  }
1427 
1428 
1429  template <int dim, int spacedim>
1430  void
1432  std::vector<IndexSet> &interface_dofs)
1433  {
1434  Assert (interface_dofs.size() == mg_dof_handler.get_triangulation().n_global_levels(),
1435  ExcDimensionMismatch (interface_dofs.size(),
1436  mg_dof_handler.get_triangulation().n_global_levels()));
1437 
1438  std::vector<std::vector<types::global_dof_index> >
1439  tmp_interface_dofs(interface_dofs.size());
1440 
1441  const FiniteElement<dim,spacedim> &fe = mg_dof_handler.get_fe();
1442 
1443  const unsigned int dofs_per_cell = fe.dofs_per_cell;
1444  const unsigned int dofs_per_face = fe.dofs_per_face;
1445 
1446  std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
1447 
1448  std::vector<bool> cell_dofs(dofs_per_cell, false);
1449 
1450  typename DoFHandler<dim>::cell_iterator cell = mg_dof_handler.begin(),
1451  endc = mg_dof_handler.end();
1452 
1453  for (; cell!=endc; ++cell)
1454  {
1455  // Do not look at artificial level cells (in a serial computation we
1456  // need to ignore the level_subdomain_id() because it is never set).
1457  if (mg_dof_handler.get_triangulation().locally_owned_subdomain()!=numbers::invalid_subdomain_id
1458  && cell->level_subdomain_id()==numbers::artificial_subdomain_id)
1459  continue;
1460 
1461  bool has_coarser_neighbor = false;
1462 
1463  std::fill (cell_dofs.begin(), cell_dofs.end(), false);
1464 
1465  for (unsigned int face_nr=0; face_nr<GeometryInfo<dim>::faces_per_cell; ++face_nr)
1466  {
1467  const typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(face_nr);
1468  if (!face->at_boundary())
1469  {
1470  //interior face
1471  const typename DoFHandler<dim>::cell_iterator
1472  neighbor = cell->neighbor(face_nr);
1473 
1474  // only process cell pairs if one or both of them are owned by me (ignore if running in serial)
1475  if (mg_dof_handler.get_triangulation().locally_owned_subdomain()!=numbers::invalid_subdomain_id
1476  &&
1477  neighbor->level_subdomain_id()==numbers::artificial_subdomain_id)
1478  continue;
1479 
1480  // Do refinement face from the coarse side
1481  if (neighbor->level() < cell->level())
1482  {
1483  for (unsigned int j=0; j<dofs_per_face; ++j)
1484  cell_dofs[fe.face_to_cell_index(j,face_nr)] = true;
1485 
1486  has_coarser_neighbor = true;
1487  }
1488  }
1489  }
1490 
1491  if (has_coarser_neighbor == false)
1492  continue;
1493 
1494  const unsigned int level = cell->level();
1495  cell->get_mg_dof_indices (local_dof_indices);
1496 
1497  for (unsigned int i=0; i<dofs_per_cell; ++i)
1498  {
1499  if (cell_dofs[i])
1500  tmp_interface_dofs[level].push_back(local_dof_indices[i]);
1501  }
1502  }
1503 
1504  for (unsigned int l=0; l<mg_dof_handler.get_triangulation().n_global_levels(); ++l)
1505  {
1506  interface_dofs[l].clear();
1507  std::sort(tmp_interface_dofs[l].begin(), tmp_interface_dofs[l].end());
1508  interface_dofs[l].add_indices(tmp_interface_dofs[l].begin(),
1509  std::unique(tmp_interface_dofs[l].begin(),
1510  tmp_interface_dofs[l].end()));
1511  interface_dofs[l].compress();
1512  }
1513 
1514  }
1515 }
1516 
1517 
1518 // explicit instantiations
1519 #include "mg_tools.inst"
1520 
1521 DEAL_II_NAMESPACE_CLOSE
const unsigned int first_hex_index
Definition: fe_base.h:271
static const unsigned int invalid_unsigned_int
Definition: types.h:170
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
const types::subdomain_id invalid_subdomain_id
Definition: types.h:245
void clear_user_flags()
Definition: tria.cc:9889
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:728
void load_user_flags(std::istream &in)
Definition: tria.cc:9946
cell_iterator end() const
Definition: dof_handler.cc:756
void make_flux_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity, const unsigned int level)
Definition: mg_tools.cc:572
void make_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity, const unsigned int level)
Definition: mg_tools.cc:539
ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:228
bool is_primitive() const
Definition: fe.h:3034
const FiniteElement< dim, spacedim > & get_fe() const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
Table< 2, Coupling > dof_couplings_from_component_couplings(const FiniteElement< dim, spacedim > &fe, const Table< 2, Coupling > &component_couplings)
unsigned int n_blocks() const
void make_boundary_list(const DoFHandler< dim, spacedim > &mg_dof, const typename FunctionMap< dim >::type &function_map, std::vector< std::set< types::global_dof_index > > &boundary_indices, const ComponentMask &component_mask=ComponentMask())
Definition: mg_tools.cc:1173
void count_dofs_per_component(const DoFHandler< dim, spacedim > &mg_dof, std::vector< std::vector< types::global_dof_index > > &result, const bool only_once=false, std::vector< unsigned int > target_component=std::vector< unsigned int >())
Definition: mg_tools.cc:987
static ::ExceptionBase & ExcMessage(std::string arg1)
ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:256
const unsigned int first_quad_index
Definition: fe_base.h:266
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:313
unsigned int element_multiplicity(const unsigned int index) const
Definition: fe.h:2863
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index) const
Definition: fe.h:2895
types::global_dof_index n_dofs() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe.cc:1251
std::map< types::boundary_id, const Function< dim, Number > * > type
Definition: function_map.h:81
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
Definition: fe.cc:528
const unsigned int dofs_per_cell
Definition: fe_base.h:295
types::global_dof_index first_block_of_base(const unsigned int b) const
Definition: fe.h:2949
void convert_couplings_to_blocks(const hp::DoFHandler< dim, spacedim > &dof_handler, const Table< 2, Coupling > &table_by_component, std::vector< Table< 2, Coupling > > &tables_by_block)
Definition: dof_tools.cc:2110
const types::subdomain_id artificial_subdomain_id
Definition: types.h:261
void compute_row_length_vector(const DoFHandler< dim, spacedim > &dofs, const unsigned int level, std::vector< unsigned int > &row_lengths, const DoFTools::Coupling flux_couplings=DoFTools::none)
Definition: mg_tools.cc:106
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
Definition: fe.h:2839
unsigned int n_components() const
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
Definition: fe.h:2982
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
const unsigned int dofs_per_face
Definition: fe_base.h:288
const unsigned int first_line_index
Definition: fe_base.h:261
void count_dofs_per_block(const DoFHandlerType &dof_handler, std::vector< std::vector< types::global_dof_index > > &dofs_per_block, std::vector< unsigned int > target_block=std::vector< unsigned int >())
Definition: mg_tools.cc:1090
void extract_inner_interface_dofs(const DoFHandler< dim, spacedim > &mg_dof_handler, std::vector< IndexSet > &interface_dofs)
Definition: mg_tools.cc:1431
void save_user_flags(std::ostream &out) const
Definition: tria.cc:9899
static ::ExceptionBase & ExcNotImplemented()
const Triangulation< dim, spacedim > & get_triangulation() const
Definition: table.h:33
unsigned char boundary_id
Definition: types.h:110
void extract_level_dofs(const unsigned int level, const DoFHandlerType &dof, const ComponentMask &component_mask, std::vector< bool > &selected_dofs)
Definition: dof_tools.cc:489
void make_flux_sparsity_pattern_edge(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity, const unsigned int level)
Definition: mg_tools.cc:647
unsigned int n_base_elements() const
Definition: fe.h:2853
Task< RT > new_task(const std_cxx11::function< RT()> &function)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe.cc:1079
static ::ExceptionBase & ExcInternalError()