Reference documentation for deal.II version 8.4.2
mg_tools.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2015 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  if ( (! cell->at_boundary(face)) &&
606  (static_cast<unsigned int>(cell->neighbor_level(face)) == level) )
607  {
609  neighbor = cell->neighbor(face);
610  neighbor->get_mg_dof_indices (dofs_on_other_cell);
611  // only add one direction The other is taken care of by
612  // neighbor (except when the neighbor is not owned by the same
613  // processor)
614  for (unsigned int i=0; i<dofs_per_cell; ++i)
615  {
616  for (unsigned int j=0; j<dofs_per_cell; ++j)
617  {
618  sparsity.add (dofs_on_this_cell[i],
619  dofs_on_other_cell[j]);
620  }
621  }
622  if (neighbor->is_locally_owned_on_level() == false)
623  for (unsigned int i=0; i<dofs_per_cell; ++i)
624  for (unsigned int j=0; j<dofs_per_cell; ++j)
625  {
626  sparsity.add (dofs_on_other_cell[i],
627  dofs_on_other_cell[j]);
628  sparsity.add (dofs_on_other_cell[i],
629  dofs_on_this_cell[j]);
630  }
631  }
632  }
633  }
634  }
635 
636 
637 
638  template <int dim, typename SparsityPatternType, int spacedim>
639  void
641  SparsityPatternType &sparsity,
642  const unsigned int level)
643  {
644  Assert ((level>=1) && (level<dof.get_triangulation().n_global_levels()),
645  ExcIndexRange(level, 1, dof.get_triangulation().n_global_levels()));
646 
647  const types::global_dof_index fine_dofs = dof.n_dofs(level);
648  const types::global_dof_index coarse_dofs = dof.n_dofs(level-1);
649  (void)fine_dofs;
650  (void)coarse_dofs;
651 
652  // Matrix maps from fine level to coarse level
653 
654  Assert (sparsity.n_rows() == coarse_dofs,
655  ExcDimensionMismatch (sparsity.n_rows(), coarse_dofs));
656  Assert (sparsity.n_cols() == fine_dofs,
657  ExcDimensionMismatch (sparsity.n_cols(), fine_dofs));
658 
659  const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
660  std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
661  std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
662  typename DoFHandler<dim,spacedim>::cell_iterator cell = dof.begin(level),
663  endc = dof.end(level);
664  for (; cell!=endc; ++cell)
665  {
666  if (!cell->is_locally_owned_on_level()) continue;
667 
668  cell->get_mg_dof_indices (dofs_on_this_cell);
669  // Loop over all interior neighbors
670  for (unsigned int face = 0;
671  face < GeometryInfo<dim>::faces_per_cell;
672  ++face)
673  {
674  // Neighbor is coarser
675 
676  if ( (! cell->at_boundary(face)) &&
677  (static_cast<unsigned int>(cell->neighbor_level(face)) != level) )
678  {
680  neighbor = cell->neighbor(face);
681  neighbor->get_mg_dof_indices (dofs_on_other_cell);
682 
683  for (unsigned int i=0; i<dofs_per_cell; ++i)
684  {
685  for (unsigned int j=0; j<dofs_per_cell; ++j)
686  {
687  sparsity.add (dofs_on_other_cell[i],
688  dofs_on_this_cell[j]);
689  sparsity.add (dofs_on_other_cell[j],
690  dofs_on_this_cell[i]);
691  }
692  }
693  }
694  }
695  }
696  }
697 
698 
699 
700  template <int dim, typename SparsityPatternType, int spacedim>
701  void
703  SparsityPatternType &sparsity,
704  const unsigned int level,
705  const Table<2,DoFTools::Coupling> &int_mask,
706  const Table<2,DoFTools::Coupling> &flux_mask)
707  {
708  const FiniteElement<dim> &fe = dof.get_fe();
709  const types::global_dof_index n_dofs = dof.n_dofs(level);
710  const unsigned int n_comp = fe.n_components();
711  (void)n_dofs;
712  (void)n_comp;
713 
714  Assert (sparsity.n_rows() == n_dofs,
715  ExcDimensionMismatch (sparsity.n_rows(), n_dofs));
716  Assert (sparsity.n_cols() == n_dofs,
717  ExcDimensionMismatch (sparsity.n_cols(), n_dofs));
718  Assert (int_mask.n_rows() == n_comp,
719  ExcDimensionMismatch (int_mask.n_rows(), n_comp));
720  Assert (int_mask.n_cols() == n_comp,
721  ExcDimensionMismatch (int_mask.n_cols(), n_comp));
722  Assert (flux_mask.n_rows() == n_comp,
723  ExcDimensionMismatch (flux_mask.n_rows(), n_comp));
724  Assert (flux_mask.n_cols() == n_comp,
725  ExcDimensionMismatch (flux_mask.n_cols(), n_comp));
726 
727  const unsigned int total_dofs = fe.dofs_per_cell;
728  std::vector<types::global_dof_index> dofs_on_this_cell(total_dofs);
729  std::vector<types::global_dof_index> dofs_on_other_cell(total_dofs);
730  Table<2,bool> support_on_face(total_dofs, GeometryInfo<dim>::faces_per_cell);
731 
732  typename DoFHandler<dim,spacedim>::cell_iterator cell = dof.begin(level),
733  endc = dof.end(level);
734 
736  int_dof_mask = DoFTools::dof_couplings_from_component_couplings(fe, int_mask),
737  flux_dof_mask = DoFTools::dof_couplings_from_component_couplings(fe, flux_mask);
738 
739  for (unsigned int i=0; i<total_dofs; ++i)
740  for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
741  support_on_face(i,f) = fe.has_support_on_face(i,f);
742 
743  // Clear user flags because we will
744  // need them. But first we save
745  // them and make sure that we
746  // restore them later such that at
747  // the end of this function the
748  // Triangulation will be in the
749  // same state as it was at the
750  // beginning of this function.
751  std::vector<bool> user_flags;
752  dof.get_triangulation().save_user_flags(user_flags);
753  const_cast<Triangulation<dim,spacedim> &>(dof.get_triangulation()).clear_user_flags ();
754 
755  for (; cell!=endc; ++cell)
756  {
757  if (!cell->is_locally_owned_on_level()) continue;
758 
759  cell->get_mg_dof_indices (dofs_on_this_cell);
760  // make sparsity pattern for this cell
761  for (unsigned int i=0; i<total_dofs; ++i)
762  for (unsigned int j=0; j<total_dofs; ++j)
763  if (int_dof_mask[i][j] != DoFTools::none)
764  sparsity.add (dofs_on_this_cell[i],
765  dofs_on_this_cell[j]);
766 
767  // Loop over all interior neighbors
768  for (unsigned int face = 0;
769  face < GeometryInfo<dim>::faces_per_cell;
770  ++face)
771  {
772  typename DoFHandler<dim,spacedim>::face_iterator cell_face = cell->face(face);
773  if (cell_face->user_flag_set ())
774  continue;
775 
776  if (cell->at_boundary (face) )
777  {
778  for (unsigned int i=0; i<total_dofs; ++i)
779  {
780  const bool i_non_zero_i = support_on_face (i, face);
781  for (unsigned int j=0; j<total_dofs; ++j)
782  {
783  const bool j_non_zero_i = support_on_face (j, face);
784 
785  if (flux_dof_mask(i,j) == DoFTools::always)
786  sparsity.add (dofs_on_this_cell[i],
787  dofs_on_this_cell[j]);
788  if (flux_dof_mask(i,j) == DoFTools::nonzero
789  && i_non_zero_i && j_non_zero_i)
790  sparsity.add (dofs_on_this_cell[i],
791  dofs_on_this_cell[j]);
792  }
793  }
794  }
795  else
796  {
798  neighbor = cell->neighbor(face);
799 
800  if (neighbor->level() < cell->level())
801  continue;
802 
803  unsigned int neighbor_face = cell->neighbor_of_neighbor(face);
804 
805  neighbor->get_mg_dof_indices (dofs_on_other_cell);
806  for (unsigned int i=0; i<total_dofs; ++i)
807  {
808  const bool i_non_zero_i = support_on_face (i, face);
809  const bool i_non_zero_e = support_on_face (i, neighbor_face);
810  for (unsigned int j=0; j<total_dofs; ++j)
811  {
812  const bool j_non_zero_i = support_on_face (j, face);
813  const bool j_non_zero_e = support_on_face (j, neighbor_face);
814  if (flux_dof_mask(i,j) == DoFTools::always)
815  {
816  sparsity.add (dofs_on_this_cell[i],
817  dofs_on_other_cell[j]);
818  sparsity.add (dofs_on_other_cell[i],
819  dofs_on_this_cell[j]);
820  sparsity.add (dofs_on_this_cell[i],
821  dofs_on_this_cell[j]);
822  sparsity.add (dofs_on_other_cell[i],
823  dofs_on_other_cell[j]);
824  }
825  if (flux_dof_mask(i,j) == DoFTools::nonzero)
826  {
827  if (i_non_zero_i && j_non_zero_e)
828  sparsity.add (dofs_on_this_cell[i],
829  dofs_on_other_cell[j]);
830  if (i_non_zero_e && j_non_zero_i)
831  sparsity.add (dofs_on_other_cell[i],
832  dofs_on_this_cell[j]);
833  if (i_non_zero_i && j_non_zero_i)
834  sparsity.add (dofs_on_this_cell[i],
835  dofs_on_this_cell[j]);
836  if (i_non_zero_e && j_non_zero_e)
837  sparsity.add (dofs_on_other_cell[i],
838  dofs_on_other_cell[j]);
839  }
840 
841  if (flux_dof_mask(j,i) == DoFTools::always)
842  {
843  sparsity.add (dofs_on_this_cell[j],
844  dofs_on_other_cell[i]);
845  sparsity.add (dofs_on_other_cell[j],
846  dofs_on_this_cell[i]);
847  sparsity.add (dofs_on_this_cell[j],
848  dofs_on_this_cell[i]);
849  sparsity.add (dofs_on_other_cell[j],
850  dofs_on_other_cell[i]);
851  }
852  if (flux_dof_mask(j,i) == DoFTools::nonzero)
853  {
854  if (j_non_zero_i && i_non_zero_e)
855  sparsity.add (dofs_on_this_cell[j],
856  dofs_on_other_cell[i]);
857  if (j_non_zero_e && i_non_zero_i)
858  sparsity.add (dofs_on_other_cell[j],
859  dofs_on_this_cell[i]);
860  if (j_non_zero_i && i_non_zero_i)
861  sparsity.add (dofs_on_this_cell[j],
862  dofs_on_this_cell[i]);
863  if (j_non_zero_e && i_non_zero_e)
864  sparsity.add (dofs_on_other_cell[j],
865  dofs_on_other_cell[i]);
866  }
867  }
868  }
869  neighbor->face(neighbor_face)->set_user_flag ();
870  }
871  }
872  }
873 
874  // finally restore the user flags
875  const_cast<Triangulation<dim,spacedim> &>(dof.get_triangulation()).load_user_flags(user_flags);
876  }
877 
878 
879 
880  template <int dim, typename SparsityPatternType, int spacedim>
881  void
883  SparsityPatternType &sparsity,
884  const unsigned int level,
885  const Table<2,DoFTools::Coupling> &flux_mask)
886  {
887  const FiniteElement<dim> &fe = dof.get_fe();
888  const unsigned int n_comp = fe.n_components();
889  (void)n_comp;
890 
891  Assert ((level>=1) && (level<dof.get_triangulation().n_global_levels()),
892  ExcIndexRange(level, 1, dof.get_triangulation().n_global_levels()));
893 
894  const types::global_dof_index fine_dofs = dof.n_dofs(level);
895  const types::global_dof_index coarse_dofs = dof.n_dofs(level-1);
896  (void)fine_dofs;
897  (void)coarse_dofs;
898 
899  // Matrix maps from fine level to coarse level
900 
901  Assert (sparsity.n_rows() == coarse_dofs,
902  ExcDimensionMismatch (sparsity.n_rows(), coarse_dofs));
903  Assert (sparsity.n_cols() == fine_dofs,
904  ExcDimensionMismatch (sparsity.n_cols(), fine_dofs));
905  Assert (flux_mask.n_rows() == n_comp,
906  ExcDimensionMismatch (flux_mask.n_rows(), n_comp));
907  Assert (flux_mask.n_cols() == n_comp,
908  ExcDimensionMismatch (flux_mask.n_cols(), n_comp));
909 
910  const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell;
911  std::vector<types::global_dof_index> dofs_on_this_cell(dofs_per_cell);
912  std::vector<types::global_dof_index> dofs_on_other_cell(dofs_per_cell);
913  Table<2,bool> support_on_face(dofs_per_cell, GeometryInfo<dim>::faces_per_cell);
914 
915  typename DoFHandler<dim,spacedim>::cell_iterator cell = dof.begin(level),
916  endc = dof.end(level);
917 
918  const Table<2,DoFTools::Coupling> flux_dof_mask
920 
921  for (unsigned int i=0; i<dofs_per_cell; ++i)
922  for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
923  support_on_face(i,f) = fe.has_support_on_face(i,f);
924 
925  for (; cell!=endc; ++cell)
926  {
927  if (!cell->is_locally_owned_on_level()) continue;
928 
929  cell->get_mg_dof_indices (dofs_on_this_cell);
930  // Loop over all interior neighbors
931  for (unsigned int face = 0;
932  face < GeometryInfo<dim>::faces_per_cell;
933  ++face)
934  {
935  // Neighbor is coarser
936 
937  if ( (! cell->at_boundary(face)) &&
938  (static_cast<unsigned int>(cell->neighbor_level(face)) != level) )
939  {
941  neighbor = cell->neighbor(face);
942  neighbor->get_mg_dof_indices (dofs_on_other_cell);
943 
944  for (unsigned int i=0; i<dofs_per_cell; ++i)
945  {
946  for (unsigned int j=0; j<dofs_per_cell; ++j)
947  {
948  if (flux_dof_mask(i,j) != DoFTools::none)
949  {
950  sparsity.add (dofs_on_other_cell[i],
951  dofs_on_this_cell[j]);
952  sparsity.add (dofs_on_other_cell[j],
953  dofs_on_this_cell[i]);
954  }
955  }
956  }
957  }
958  }
959  }
960  }
961 
962 
963 
964  template <int dim, int spacedim>
965  void
967  std::vector<std::vector<types::global_dof_index> > &result,
968  bool only_once,
969  std::vector<unsigned int> target_component)
970  {
971  const FiniteElement<dim> &fe = dof_handler.get_fe();
972  const unsigned int n_components = fe.n_components();
973  const unsigned int nlevels = dof_handler.get_triangulation().n_global_levels();
974 
975  Assert (result.size() == nlevels,
976  ExcDimensionMismatch(result.size(), nlevels));
977 
978  if (target_component.size() == 0)
979  {
980  target_component.resize(n_components);
981  for (unsigned int i=0; i<n_components; ++i)
982  target_component[i] = i;
983  }
984 
985  Assert(target_component.size() == n_components,
986  ExcDimensionMismatch(target_component.size(), n_components));
987 
988  for (unsigned int l=0; l<nlevels; ++l)
989  {
990  result[l].resize (n_components);
991  std::fill (result[l].begin(),result[l].end(), 0U);
992 
993  // special case for only one
994  // component. treat this first
995  // since it does not require any
996  // computations
997  if (n_components == 1)
998  {
999  result[l][0] = dof_handler.n_dofs(l);
1000  }
1001  else
1002  {
1003  // otherwise determine the number
1004  // of dofs in each component
1005  // separately. do so in parallel
1006  std::vector<std::vector<bool> >
1007  dofs_in_component (n_components,
1008  std::vector<bool>(dof_handler.n_dofs(l),
1009  false));
1010  std::vector<ComponentMask> component_select (n_components);
1011  Threads::TaskGroup<> tasks;
1012  for (unsigned int i=0; i<n_components; ++i)
1013  {
1014  void (*fun_ptr) (const unsigned int level,
1015  const DoFHandler<dim,spacedim> &,
1016  const ComponentMask &,
1017  std::vector<bool> &)
1018  = &DoFTools::extract_level_dofs<DoFHandler<dim,spacedim> >;
1019 
1020  std::vector<bool> tmp(n_components, false);
1021  tmp[i] = true;
1022  component_select[i] = ComponentMask(tmp);
1023 
1024  tasks += Threads::new_task (fun_ptr,
1025  l, dof_handler,
1026  component_select[i],
1027  dofs_in_component[i]);
1028  }
1029  tasks.join_all();
1030 
1031  // next count what we got
1032  unsigned int component = 0;
1033  for (unsigned int b=0; b<fe.n_base_elements(); ++b)
1034  {
1035  const FiniteElement<dim> &base = fe.base_element(b);
1036  // Dimension of base element
1037  unsigned int d = base.n_components();
1038 
1039  for (unsigned int m=0; m<fe.element_multiplicity(b); ++m)
1040  {
1041  for (unsigned int dd=0; dd<d; ++dd)
1042  {
1043  if (base.is_primitive() || (!only_once || dd==0))
1044  result[l][target_component[component]]
1045  += std::count(dofs_in_component[component].begin(),
1046  dofs_in_component[component].end(),
1047  true);
1048  ++component;
1049  }
1050  }
1051  }
1052  // finally sanity check
1053  Assert (!dof_handler.get_fe().is_primitive()
1054  ||
1055  std::accumulate (result[l].begin(),
1056  result[l].end(), 0U)
1057  ==
1058  dof_handler.n_dofs(l),
1059  ExcInternalError());
1060  }
1061  }
1062  }
1063 
1064 
1065 
1066  template <typename DoFHandlerType>
1067  void
1069  (const DoFHandlerType &dof_handler,
1070  std::vector<std::vector<types::global_dof_index> > &dofs_per_block,
1071  std::vector<unsigned int> target_block)
1072  {
1074  const unsigned int n_blocks = fe.n_blocks();
1075  const unsigned int n_levels = dof_handler.get_triangulation().n_global_levels();
1076 
1077  AssertDimension (dofs_per_block.size(), n_levels);
1078 
1079  for (unsigned int l=0; l<n_levels; ++l)
1080  std::fill (dofs_per_block[l].begin(), dofs_per_block[l].end(), 0U);
1081  // If the empty vector was given as
1082  // default argument, set up this
1083  // vector as identity.
1084  if (target_block.size()==0)
1085  {
1086  target_block.resize(n_blocks);
1087  for (unsigned int i=0; i<n_blocks; ++i)
1088  target_block[i] = i;
1089  }
1090  Assert(target_block.size()==n_blocks,
1091  ExcDimensionMismatch(target_block.size(),n_blocks));
1092 
1093  const unsigned int max_block
1094  = *std::max_element (target_block.begin(),
1095  target_block.end());
1096  const unsigned int n_target_blocks = max_block + 1;
1097  (void)n_target_blocks;
1098 
1099  for (unsigned int l=0; l<n_levels; ++l)
1100  AssertDimension (dofs_per_block[l].size(), n_target_blocks);
1101 
1102  // special case for only one
1103  // block. treat this first
1104  // since it does not require any
1105  // computations
1106  if (n_blocks == 1)
1107  {
1108  for (unsigned int l=0; l<n_levels; ++l)
1109  dofs_per_block[l][0] = dof_handler.n_dofs(l);
1110  return;
1111  }
1112  // otherwise determine the number
1113  // of dofs in each block
1114  // separately. do so in parallel
1115  for (unsigned int l=0; l<n_levels; ++l)
1116  {
1117  std::vector<std::vector<bool> >
1118  dofs_in_block (n_blocks, std::vector<bool>(dof_handler.n_dofs(l), false));
1119  std::vector<BlockMask> block_select (n_blocks);
1120  Threads::TaskGroup<> tasks;
1121  for (unsigned int i=0; i<n_blocks; ++i)
1122  {
1123  void (*fun_ptr) (const unsigned int level,
1124  const DoFHandlerType &,
1125  const BlockMask &,
1126  std::vector<bool> &)
1127  = &DoFTools::extract_level_dofs<DoFHandlerType>;
1128 
1129  std::vector<bool> tmp(n_blocks, false);
1130  tmp[i] = true;
1131  block_select[i] = tmp;
1132 
1133  tasks += Threads::new_task (fun_ptr,
1134  l, dof_handler, block_select[i],
1135  dofs_in_block[i]);
1136  }
1137  tasks.join_all ();
1138 
1139  // next count what we got
1140  for (unsigned int block=0; block<fe.n_blocks(); ++block)
1141  dofs_per_block[l][target_block[block]]
1142  += std::count(dofs_in_block[block].begin(),
1143  dofs_in_block[block].end(),
1144  true);
1145  }
1146  }
1147 
1148 
1149  template <>
1150  void
1152  const DoFHandler<1,1> &,
1153  const FunctionMap<1>::type &,
1154  std::vector<std::set<types::global_dof_index> > &,
1155  const ComponentMask &)
1156  {
1157  Assert(false, ExcNotImplemented());
1158  }
1159 
1160 
1161 
1162  template <>
1163  void
1165  const DoFHandler<1,2> &,
1166  const FunctionMap<1>::type &,
1167  std::vector<std::set<types::global_dof_index> > &,
1168  const ComponentMask &)
1169  {
1170  Assert(false, ExcNotImplemented());
1171  }
1172 
1173 
1174 
1175  template <int dim, int spacedim>
1176  void
1178  const DoFHandler<dim,spacedim> &dof,
1179  const typename FunctionMap<dim>::type &function_map,
1180  std::vector<std::set<types::global_dof_index> > &boundary_indices,
1181  const ComponentMask &component_mask)
1182  {
1183  // if for whatever reason we were
1184  // passed an empty map, return
1185  // immediately
1186  if (function_map.size() == 0)
1187  return;
1188 
1189  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
1190 
1191  (void)n_levels;
1192 
1193 
1194  const unsigned int n_components = DoFTools::n_components(dof);
1195  const bool fe_is_system = (n_components != 1);
1196 
1197  AssertDimension (boundary_indices.size(), n_levels);
1198 
1199  std::vector<types::global_dof_index> local_dofs;
1200  local_dofs.reserve (DoFTools::max_dofs_per_face(dof));
1201  std::fill (local_dofs.begin (),
1202  local_dofs.end (),
1204 
1205  // First, deal with the simpler
1206  // case when we have to identify
1207  // all boundary dofs
1208  if (component_mask.n_selected_components(n_components) == n_components)
1209  {
1211  cell = dof.begin(),
1212  endc = dof.end();
1213  for (; cell!=endc; ++cell)
1214  {
1215  if (dof.get_triangulation().locally_owned_subdomain()!=numbers::invalid_subdomain_id
1216  && cell->level_subdomain_id()==numbers::artificial_subdomain_id)
1217  continue;
1218  const FiniteElement<dim> &fe = cell->get_fe();
1219  const unsigned int level = cell->level();
1220  local_dofs.resize(fe.dofs_per_face);
1221 
1222  for (unsigned int face_no = 0; face_no < GeometryInfo<dim>::faces_per_cell;
1223  ++face_no)
1224  if (cell->at_boundary(face_no) == true)
1225  {
1226  const typename DoFHandler<dim,spacedim>::face_iterator
1227  face = cell->face(face_no);
1228  const types::boundary_id bi = face->boundary_id();
1229  // Face is listed in
1230  // boundary map
1231  if (function_map.find(bi) != function_map.end())
1232  {
1233  face->get_mg_dof_indices(level, local_dofs);
1234  for (unsigned int i=0; i<fe.dofs_per_face; ++i)
1235  boundary_indices[level].insert(local_dofs[i]);
1236  }
1237  }
1238  }
1239  }
1240  else
1241  {
1242  Assert (component_mask.n_selected_components(n_components) > 0,
1243  ExcMessage("It's probably worthwhile to select at least one component."));
1244 
1246  cell = dof.begin(),
1247  endc = dof.end();
1248  for (; cell!=endc; ++cell)
1249  if (dof.get_triangulation().locally_owned_subdomain()==numbers::invalid_subdomain_id
1250  || cell->level_subdomain_id()!=numbers::artificial_subdomain_id)
1251  for (unsigned int face_no = 0; face_no < GeometryInfo<dim>::faces_per_cell;
1252  ++face_no)
1253  {
1254  if (cell->at_boundary(face_no) == false)
1255  continue;
1256 
1257  const FiniteElement<dim> &fe = cell->get_fe();
1258  const unsigned int level = cell->level();
1259 
1260  // we can presently deal only with
1261  // primitive elements for boundary
1262  // values. this does not preclude
1263  // us using non-primitive elements
1264  // in components that we aren't
1265  // interested in, however. make
1266  // sure that all shape functions
1267  // that are non-zero for the
1268  // components we are interested in,
1269  // are in fact primitive
1270  for (unsigned int i=0; i<cell->get_fe().dofs_per_cell; ++i)
1271  {
1272  const ComponentMask &nonzero_component_array
1273  = cell->get_fe().get_nonzero_components (i);
1274  for (unsigned int c=0; c<n_components; ++c)
1275  if ((nonzero_component_array[c] == true)
1276  &&
1277  (component_mask[c] == true))
1278  Assert (cell->get_fe().is_primitive (i),
1279  ExcMessage ("This function can only deal with requested boundary "
1280  "values that correspond to primitive (scalar) base "
1281  "elements"));
1282  }
1283 
1284  typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(face_no);
1285  const types::boundary_id boundary_component = face->boundary_id();
1286  if (function_map.find(boundary_component) != function_map.end())
1287  // face is of the right component
1288  {
1289  // get indices, physical location and
1290  // boundary values of dofs on this
1291  // face
1292  local_dofs.resize (fe.dofs_per_face);
1293  face->get_mg_dof_indices (level, local_dofs);
1294  if (fe_is_system)
1295  {
1296  // enter those dofs
1297  // into the list that
1298  // match the
1299  // component
1300  // signature. avoid
1301  // the usual
1302  // complication that
1303  // we can't just use
1304  // *_system_to_component_index
1305  // for non-primitive
1306  // FEs
1307  for (unsigned int i=0; i<local_dofs.size(); ++i)
1308  {
1309  unsigned int component;
1310  if (fe.is_primitive())
1311  component = fe.face_system_to_component_index(i).first;
1312  else
1313  {
1314  // non-primitive
1315  // case. make
1316  // sure that
1317  // this
1318  // particular
1319  // shape
1320  // function
1321  // _is_
1322  // primitive,
1323  // and get at
1324  // it's
1325  // component. use
1326  // usual
1327  // trick to
1328  // transfer
1329  // face dof
1330  // index to
1331  // cell dof
1332  // index
1333  const unsigned int cell_i
1334  = (dim == 1 ?
1335  i
1336  :
1337  (dim == 2 ?
1338  (i<2*fe.dofs_per_vertex ? i : i+2*fe.dofs_per_vertex)
1339  :
1340  (dim == 3 ?
1341  (i<4*fe.dofs_per_vertex ?
1342  i
1343  :
1344  (i<4*fe.dofs_per_vertex+4*fe.dofs_per_line ?
1345  i+4*fe.dofs_per_vertex
1346  :
1347  i+4*fe.dofs_per_vertex+8*fe.dofs_per_line))
1348  :
1350  Assert (cell_i < fe.dofs_per_cell, ExcInternalError());
1351 
1352  // make sure
1353  // that if
1354  // this is
1355  // not a
1356  // primitive
1357  // shape function,
1358  // then all
1359  // the
1360  // corresponding
1361  // components
1362  // in the
1363  // mask are
1364  // not set
1365 // if (!fe.is_primitive(cell_i))
1366 // for (unsigned int c=0; c<n_components; ++c)
1367 // if (fe.get_nonzero_components(cell_i)[c])
1368 // Assert (component_mask[c] == false,
1369 // ExcFENotPrimitive());
1370 
1371 // let's pick the first of possibly more than one non-zero
1372 // components. if shape function is non-primitive, then we will ignore
1373 // the result in the following anyway, otherwise there's only one
1374 // non-zero component which we will use
1375  component = fe.get_nonzero_components(cell_i).first_selected_component();
1376  }
1377 
1378  if (component_mask[component] == true)
1379  boundary_indices[level].insert(local_dofs[i]);
1380  }
1381  }
1382  else
1383  for (unsigned int i=0; i<local_dofs.size(); ++i)
1384  boundary_indices[level].insert(local_dofs[i]);
1385  }
1386  }
1387  }
1388  }
1389 
1390 
1391  template <int dim, int spacedim>
1392  void
1394  const typename FunctionMap<dim>::type &function_map,
1395  std::vector<IndexSet> &boundary_indices,
1396  const ComponentMask &component_mask)
1397  {
1398  Assert (boundary_indices.size() == dof.get_triangulation().n_global_levels(),
1399  ExcDimensionMismatch (boundary_indices.size(),
1400  dof.get_triangulation().n_global_levels()));
1401 
1402  std::vector<std::set<types::global_dof_index> >
1403  my_boundary_indices (dof.get_triangulation().n_global_levels());
1404  make_boundary_list (dof, function_map, my_boundary_indices, component_mask);
1405  for (unsigned int i=0; i<dof.get_triangulation().n_global_levels(); ++i)
1406  {
1407  boundary_indices[i] = IndexSet (dof.n_dofs(i));
1408  boundary_indices[i].add_indices (my_boundary_indices[i].begin(),
1409  my_boundary_indices[i].end());
1410  }
1411  }
1412 
1413 
1414  template <int dim, int spacedim>
1415  void
1416  extract_non_interface_dofs (const DoFHandler<dim,spacedim> &mg_dof_handler,
1417  std::vector<std::set<types::global_dof_index> > &non_interface_dofs)
1418  {
1419  Assert (non_interface_dofs.size() == mg_dof_handler.get_triangulation().n_global_levels(),
1420  ExcDimensionMismatch (non_interface_dofs.size(),
1421  mg_dof_handler.get_triangulation().n_global_levels()));
1422 
1423  const FiniteElement<dim,spacedim> &fe = mg_dof_handler.get_fe();
1424 
1425  const unsigned int dofs_per_cell = fe.dofs_per_cell;
1426  const unsigned int dofs_per_face = fe.dofs_per_face;
1427 
1428  std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
1429  std::vector<bool> cell_dofs(dofs_per_cell, false);
1430  std::vector<bool> cell_dofs_interface(dofs_per_cell, false);
1431 
1432  typename DoFHandler<dim>::cell_iterator cell = mg_dof_handler.begin(),
1433  endc = mg_dof_handler.end();
1434 
1435 
1436  for (; cell!=endc; ++cell)
1437  {
1438  if (mg_dof_handler.get_triangulation().locally_owned_subdomain()!=numbers::invalid_subdomain_id
1439  && cell->level_subdomain_id()!=mg_dof_handler.get_triangulation().locally_owned_subdomain())
1440  continue;
1441 
1442  std::fill (cell_dofs.begin(), cell_dofs.end(), false);
1443  std::fill (cell_dofs_interface.begin(), cell_dofs_interface.end(), false);
1444 
1445  for (unsigned int face_nr=0; face_nr<GeometryInfo<dim>::faces_per_cell; ++face_nr)
1446  {
1447  const typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(face_nr);
1448  if (!face->at_boundary())
1449  {
1450  //interior face
1451  const typename DoFHandler<dim>::cell_iterator
1452  neighbor = cell->neighbor(face_nr);
1453 
1454  if ((neighbor->level() < cell->level()))
1455  {
1456  for (unsigned int j=0; j<dofs_per_face; ++j)
1457  cell_dofs_interface[fe.face_to_cell_index(j,face_nr)] = true;
1458  }
1459  else
1460  {
1461  for (unsigned int j=0; j<dofs_per_face; ++j)
1462  cell_dofs[fe.face_to_cell_index(j,face_nr)] = true;
1463  }
1464  }
1465  else
1466  {
1467  //boundary face
1468  for (unsigned int j=0; j<dofs_per_face; ++j)
1469  cell_dofs[fe.face_to_cell_index(j,face_nr)] = true;
1470  }
1471  }
1472 
1473  const unsigned int level = cell->level();
1474  cell->get_mg_dof_indices (local_dof_indices);
1475 
1476  for (unsigned int i=0; i<dofs_per_cell; ++i)
1477  if (cell_dofs[i] && !cell_dofs_interface[i])
1478  non_interface_dofs[level].insert(local_dof_indices[i]);
1479  }
1480  }
1481 
1482 
1483  template <int dim, int spacedim>
1484  void
1486  std::vector<IndexSet> &interface_dofs)
1487  {
1488  Assert (interface_dofs.size() == mg_dof_handler.get_triangulation().n_global_levels(),
1489  ExcDimensionMismatch (interface_dofs.size(),
1490  mg_dof_handler.get_triangulation().n_global_levels()));
1491 
1492  std::vector<std::vector<types::global_dof_index> >
1493  tmp_interface_dofs(interface_dofs.size());
1494 
1495  const FiniteElement<dim,spacedim> &fe = mg_dof_handler.get_fe();
1496 
1497  const unsigned int dofs_per_cell = fe.dofs_per_cell;
1498  const unsigned int dofs_per_face = fe.dofs_per_face;
1499 
1500  std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
1501 
1502  std::vector<bool> cell_dofs(dofs_per_cell, false);
1503 
1504  typename DoFHandler<dim>::cell_iterator cell = mg_dof_handler.begin(),
1505  endc = mg_dof_handler.end();
1506 
1507  for (; cell!=endc; ++cell)
1508  {
1509  // Do not look at artificial level cells (in a serial computation we
1510  // need to ignore the level_subdomain_id() because it is never set).
1511  if (mg_dof_handler.get_triangulation().locally_owned_subdomain()!=numbers::invalid_subdomain_id
1512  && cell->level_subdomain_id()==numbers::artificial_subdomain_id)
1513  continue;
1514 
1515  bool has_coarser_neighbor = false;
1516 
1517  std::fill (cell_dofs.begin(), cell_dofs.end(), false);
1518 
1519  for (unsigned int face_nr=0; face_nr<GeometryInfo<dim>::faces_per_cell; ++face_nr)
1520  {
1521  const typename DoFHandler<dim,spacedim>::face_iterator face = cell->face(face_nr);
1522  if (!face->at_boundary())
1523  {
1524  //interior face
1525  const typename DoFHandler<dim>::cell_iterator
1526  neighbor = cell->neighbor(face_nr);
1527 
1528  // only process cell pairs if one or both of them are owned by me (ignore if running in serial)
1529  if (mg_dof_handler.get_triangulation().locally_owned_subdomain()!=numbers::invalid_subdomain_id
1530  &&
1531  neighbor->level_subdomain_id()==numbers::artificial_subdomain_id)
1532  continue;
1533 
1534  // Do refinement face from the coarse side
1535  if (neighbor->level() < cell->level())
1536  {
1537  for (unsigned int j=0; j<dofs_per_face; ++j)
1538  cell_dofs[fe.face_to_cell_index(j,face_nr)] = true;
1539 
1540  has_coarser_neighbor = true;
1541  }
1542  }
1543  }
1544 
1545  if (has_coarser_neighbor == false)
1546  continue;
1547 
1548  const unsigned int level = cell->level();
1549  cell->get_mg_dof_indices (local_dof_indices);
1550 
1551  for (unsigned int i=0; i<dofs_per_cell; ++i)
1552  {
1553  if (cell_dofs[i])
1554  tmp_interface_dofs[level].push_back(local_dof_indices[i]);
1555  }
1556  }
1557 
1558  for (unsigned int l=0; l<mg_dof_handler.get_triangulation().n_global_levels(); ++l)
1559  {
1560  interface_dofs[l].clear();
1561  std::sort(tmp_interface_dofs[l].begin(), tmp_interface_dofs[l].end());
1562  interface_dofs[l].add_indices(tmp_interface_dofs[l].begin(),
1563  std::unique(tmp_interface_dofs[l].begin(),
1564  tmp_interface_dofs[l].end()));
1565  interface_dofs[l].compress();
1566  }
1567 
1568  }
1569 }
1570 
1571 
1572 // explicit instantiations
1573 #include "mg_tools.inst"
1574 
1575 DEAL_II_NAMESPACE_CLOSE
const unsigned int first_hex_index
Definition: fe_base.h:259
static const unsigned int invalid_unsigned_int
Definition: types.h:164
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
const types::subdomain_id invalid_subdomain_id
Definition: types.h:239
void clear_user_flags()
Definition: tria.cc:9618
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:833
void load_user_flags(std::istream &in)
Definition: tria.cc:9675
::ExceptionBase & ExcMessage(std::string arg1)
cell_iterator end() const
Definition: dof_handler.cc:861
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:221
const FiniteElement< dim, spacedim > & get_fe() const
const unsigned int dofs_per_line
Definition: fe_base.h:232
Table< 2, Coupling > dof_couplings_from_component_couplings(const FiniteElement< dim, spacedim > &fe, const Table< 2, Coupling > &component_couplings)
unsigned int n_blocks() const
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: dof_handler.cc:845
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:1177
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:966
ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:249
const unsigned int first_quad_index
Definition: fe_base.h:254
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:294
unsigned int element_multiplicity(const unsigned int index) const
Definition: fe.h:2766
bool is_primitive(const unsigned int i) const
Definition: fe.h:2937
const ComponentMask & get_nonzero_components(const unsigned int i) const
Definition: fe.h:2915
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index) const
Definition: fe.h:2798
types::global_dof_index n_dofs() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe.cc:1236
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:529
const unsigned int dofs_per_cell
Definition: fe_base.h:283
types::global_dof_index first_block_of_base(const unsigned int b) const
Definition: fe.h:2852
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:2007
const types::subdomain_id artificial_subdomain_id
Definition: types.h:255
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:2742
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:2885
unsigned int first_selected_component(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
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:276
const unsigned int first_line_index
Definition: fe_base.h:249
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:1069
void extract_inner_interface_dofs(const DoFHandler< dim, spacedim > &mg_dof_handler, std::vector< IndexSet > &interface_dofs)
Definition: mg_tools.cc:1485
void save_user_flags(std::ostream &out) const
Definition: tria.cc:9628
const Triangulation< dim, spacedim > & get_triangulation() const
const unsigned int dofs_per_vertex
Definition: fe_base.h:226
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:453
void make_flux_sparsity_pattern_edge(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternType &sparsity, const unsigned int level)
Definition: mg_tools.cc:640
unsigned int n_base_elements() const
Definition: fe.h:2756
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:1080