Reference documentation for deal.II version 9.0.0
quadrature.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 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/geometry_info.h>
17 #include <deal.II/base/quadrature.h>
18 #include <deal.II/base/qprojector.h>
19 #include <deal.II/base/memory_consumption.h>
20 #include <deal.II/base/utilities.h>
21 #include <deal.II/base/std_cxx14/memory.h>
22 
23 #include <cmath>
24 #include <cstdlib>
25 #include <iterator>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 
30 template <>
31 Quadrature<0>::Quadrature (const unsigned int n_q)
32  :
33  quadrature_points (n_q),
34  weights (n_q, 0),
35  is_tensor_product_flag (false)
36 {}
37 
38 
39 
40 template <int dim>
41 Quadrature<dim>::Quadrature (const unsigned int n_q)
42  :
43  quadrature_points (n_q, Point<dim>()),
44  weights (n_q, 0),
45  is_tensor_product_flag (dim==1)
46 {}
47 
48 
49 
50 template <int dim>
51 void
52 Quadrature<dim>::initialize (const std::vector<Point<dim> > &p,
53  const std::vector<double> &w)
54 {
55  AssertDimension (w.size(), p.size());
56  quadrature_points = p;
57  weights = w;
58 }
59 
60 
61 
62 template <int dim>
63 Quadrature<dim>::Quadrature (const std::vector<Point<dim> > &points,
64  const std::vector<double> &weights)
65  :
66  quadrature_points(points),
67  weights(weights),
68  is_tensor_product_flag (dim==1)
69 {
70  Assert (weights.size() == points.size(),
71  ExcDimensionMismatch(weights.size(), points.size()));
72 }
73 
74 
75 
76 template <int dim>
77 Quadrature<dim>::Quadrature (const std::vector<Point<dim> > &points)
78  :
79  quadrature_points(points),
80  weights(points.size(), std::numeric_limits<double>::infinity()),
81  is_tensor_product_flag (dim==1)
82 {
83  Assert(weights.size() == points.size(),
84  ExcDimensionMismatch(weights.size(), points.size()));
85 }
86 
87 
88 
89 template <int dim>
91  :
92  quadrature_points(std::vector<Point<dim> > (1, point)),
93  weights(std::vector<double> (1, 1.)),
94  is_tensor_product_flag (true),
95  tensor_basis (new std::array<Quadrature<1>, dim>())
96 {
97  for (unsigned int i=0; i<dim; ++i)
98  {
99  const std::vector<Point<1> > quad_vec_1d (1, Point<1>(point[i]));
100  (*tensor_basis)[i] = Quadrature<1>(quad_vec_1d, weights);
101  }
102 }
103 
104 
105 
106 template <>
107 Quadrature<1>::Quadrature (const Point<1> &point)
108  :
109  quadrature_points(std::vector<Point<1> > (1, point)),
110  weights(std::vector<double> (1, 1.)),
111  is_tensor_product_flag (true)
112 {
113 }
114 
115 
116 
117 template <>
118 Quadrature<0>::Quadrature (const SubQuadrature &,
119  const Quadrature<1> &)
120 {
121  Assert(false, ExcImpossibleInDim(0));
122 }
123 
124 
125 
126 template <int dim>
128  const Quadrature<1> &q2)
129  :
130  quadrature_points (q1.size() * q2.size()),
131  weights (q1.size() * q2.size()),
132  is_tensor_product_flag (q1.is_tensor_product())
133 {
134  unsigned int present_index = 0;
135  for (unsigned int i2=0; i2<q2.size(); ++i2)
136  for (unsigned int i1=0; i1<q1.size(); ++i1)
137  {
138  // compose coordinates of
139  // new quadrature point by tensor
140  // product in the last component
141  for (unsigned int d=0; d<dim-1; ++d)
142  quadrature_points[present_index](d)
143  = q1.point(i1)(d);
144  quadrature_points[present_index](dim-1)
145  = q2.point(i2)(0);
146 
147  weights[present_index] = q1.weight(i1) * q2.weight(i2);
148 
149  ++present_index;
150  };
151 
152 #ifdef DEBUG
153  if (size() > 0)
154  {
155  double sum = 0;
156  for (unsigned int i=0; i<size(); ++i)
157  sum += weights[i];
158  // we cannot guarantee the sum of weights
159  // to be exactly one, but it should be
160  // near that.
161  Assert ((sum>0.999999) && (sum<1.000001), ExcInternalError());
162  }
163 #endif
164 
166  {
167  tensor_basis = std_cxx14::make_unique<std::array<Quadrature<1>, dim>> ();
168  for (unsigned int i=0; i<dim-1; ++i)
169  (*tensor_basis)[i] = q1.get_tensor_basis()[i];
170  (*tensor_basis)[dim-1] = q2;
171  }
172 }
173 
174 
175 
176 template <>
177 Quadrature<1>::Quadrature (const SubQuadrature &,
178  const Quadrature<1> &q2)
179  :
180  quadrature_points (q2.size()),
181  weights (q2.size()),
182  is_tensor_product_flag (true)
183 {
184  unsigned int present_index = 0;
185  for (unsigned int i2=0; i2<q2.size(); ++i2)
186  {
187  // compose coordinates of
188  // new quadrature point by tensor
189  // product in the last component
190  quadrature_points[present_index](0)
191  = q2.point(i2)(0);
192 
193  weights[present_index] = q2.weight(i2);
194 
195  ++present_index;
196  }
197 
198 #ifdef DEBUG
199  if (size() > 0)
200  {
201  double sum = 0;
202  for (unsigned int i=0; i<size(); ++i)
203  sum += weights[i];
204  // we cannot guarantee the sum of weights
205  // to be exactly one, but it should be
206  // near that.
207  Assert ((sum>0.999999) && (sum<1.000001), ExcInternalError());
208  }
209 #endif
210 }
211 
212 
213 
214 template <>
216  :
217  Subscriptor(),
218 // quadrature_points(1),
219  weights(1,1.),
220  is_tensor_product_flag (false)
221 {}
222 
223 
224 template <>
226  :
227  Subscriptor()
228 {
229  // this function should never be
230  // called -- this should be the
231  // copy constructor in 1d...
232  Assert (false, ExcImpossibleInDim(1));
233 }
234 
235 
236 
237 template <int dim>
239 :
240 Subscriptor(),
241  quadrature_points (Utilities::fixed_power<dim>(q.size())),
242  weights (Utilities::fixed_power<dim>(q.size())),
243  is_tensor_product_flag(true)
244 {
245  Assert (dim <= 3, ExcNotImplemented());
246 
247  const unsigned int n0 = q.size();
248  const unsigned int n1 = (dim>1) ? n0 : 1;
249  const unsigned int n2 = (dim>2) ? n0 : 1;
250 
251  unsigned int k=0;
252  for (unsigned int i2=0; i2<n2; ++i2)
253  for (unsigned int i1=0; i1<n1; ++i1)
254  for (unsigned int i0=0; i0<n0; ++i0)
255  {
256  quadrature_points[k](0) = q.point(i0)(0);
257  if (dim>1)
258  quadrature_points[k](1) = q.point(i1)(0);
259  if (dim>2)
260  quadrature_points[k](2) = q.point(i2)(0);
261  weights[k] = q.weight(i0);
262  if (dim>1)
263  weights[k] *= q.weight(i1);
264  if (dim>2)
265  weights[k] *= q.weight(i2);
266  ++k;
267  }
268 
269  tensor_basis = std_cxx14::make_unique<std::array<Quadrature<1>, dim>> ();
270  for (unsigned int i=0; i<dim; ++i)
271  (*tensor_basis)[i] = q;
272 }
273 
274 
275 
276 template <int dim>
278  :
279  Subscriptor(),
280  quadrature_points (q.quadrature_points),
281  weights (q.weights),
282  is_tensor_product_flag (q.is_tensor_product_flag)
283 {
284  if (dim>1 && is_tensor_product_flag)
285  tensor_basis = std_cxx14::make_unique<std::array<Quadrature<1>, dim>> (*q.tensor_basis);
286 }
287 
288 
289 
290 template <int dim>
293 {
294  weights = q.weights;
295  quadrature_points = q.quadrature_points;
296  is_tensor_product_flag = q.is_tensor_product_flag;
297  if (dim >1 && is_tensor_product_flag)
298  {
299  if (tensor_basis == nullptr)
300  tensor_basis = std_cxx14::make_unique<std::array<Quadrature<1>, dim>> (*q.tensor_basis);
301  else
302  *tensor_basis = *q.tensor_basis;
303  }
304  return *this;
305 }
306 
307 
308 
309 template <int dim>
310 bool
312 {
313  return ((quadrature_points == q.quadrature_points)
314  &&
315  (weights == q.weights));
316 }
317 
318 
319 
320 template <int dim>
321 std::size_t
323 {
324  return (MemoryConsumption::memory_consumption (quadrature_points) +
326 }
327 
328 
329 
330 template <int dim>
331 typename std::conditional<dim==1, std::array<Quadrature<1>, dim>,const std::array<Quadrature<1>,dim>&>::type
333 {
334  Assert (this->is_tensor_product_flag == true,
335  ExcMessage("This function only makes sense if "
336  "this object represents a tensor product!"));
337  Assert (tensor_basis != nullptr, ExcInternalError());
338 
339  return *tensor_basis;
340 }
341 
342 
343 
344 template <>
345 std::array<Quadrature<1>, 1>
347 {
348  Assert (this->is_tensor_product_flag == true,
349  ExcMessage("This function only makes sense if "
350  "this object represents a tensor product!"));
351 
352  return std::array<Quadrature<1>, 1> {{*this}};
353 }
354 
355 
356 
357 //---------------------------------------------------------------------------
358 template <int dim>
360  : Quadrature<dim>(qx.size())
361 {
362  Assert (dim==1, ExcImpossibleInDim(dim));
363  unsigned int k=0;
364  for (unsigned int k1=0; k1<qx.size(); ++k1)
365  {
366  this->quadrature_points[k](0) = qx.point(k1)(0);
367  this->weights[k++] = qx.weight(k1);
368  }
369  Assert (k==this->size(), ExcInternalError());
370  this->is_tensor_product_flag = true;
371 }
372 
373 
374 
375 template <int dim>
377  const Quadrature<1> &qy)
378  : Quadrature<dim>(qx.size()*qy.size())
379 {
380  Assert (dim==2, ExcImpossibleInDim(dim));
381 }
382 
383 
384 
385 template <>
387  const Quadrature<1> &qy)
388  : Quadrature<2>(qx.size()*qy.size())
389 {
390  unsigned int k=0;
391  for (unsigned int k2=0; k2<qy.size(); ++k2)
392  for (unsigned int k1=0; k1<qx.size(); ++k1)
393  {
394  this->quadrature_points[k](0) = qx.point(k1)(0);
395  this->quadrature_points[k](1) = qy.point(k2)(0);
396  this->weights[k++] = qx.weight(k1) * qy.weight(k2);
397  }
398  Assert (k==this->size(), ExcInternalError());
399  this->is_tensor_product_flag = true;
400  const std::array<Quadrature<1>, 2> q_array {{qx, qy}};
401  this->tensor_basis = std_cxx14::make_unique<std::array<Quadrature<1>, 2>>(q_array);
402 }
403 
404 
405 
406 template <int dim>
408  const Quadrature<1> &qy,
409  const Quadrature<1> &qz)
410  : Quadrature<dim>(qx.size()*qy.size()*qz.size())
411 {
412  Assert (dim==3, ExcImpossibleInDim(dim));
413 }
414 
415 
416 
417 template <>
419  const Quadrature<1> &qy,
420  const Quadrature<1> &qz)
421  : Quadrature<3>(qx.size()*qy.size()*qz.size())
422 {
423  unsigned int k=0;
424  for (unsigned int k3=0; k3<qz.size(); ++k3)
425  for (unsigned int k2=0; k2<qy.size(); ++k2)
426  for (unsigned int k1=0; k1<qx.size(); ++k1)
427  {
428  this->quadrature_points[k](0) = qx.point(k1)(0);
429  this->quadrature_points[k](1) = qy.point(k2)(0);
430  this->quadrature_points[k](2) = qz.point(k3)(0);
431  this->weights[k++] = qx.weight(k1) * qy.weight(k2) * qz.weight(k3);
432  }
433  Assert (k==this->size(), ExcInternalError());
434  this->is_tensor_product_flag = true;
435  const std::array<Quadrature<1>, 3> q_array {{qx, qy, qz}};
436  this->tensor_basis = std_cxx14::make_unique<std::array<Quadrature<1>, 3>>(q_array);
437 }
438 
439 
440 
441 //---------------------------------------------------------------------------
442 
443 
444 
445 template <int dim>
448 {
449  std::vector<Point<2> > q_points (q.size());
450  std::vector<double> weights (q.size());
451  for (unsigned int i=0; i<q.size(); ++i)
452  {
453  q_points[i][0] = q.point(i)[1];
454  q_points[i][1] = q.point(i)[0];
455 
456  weights[i] = q.weight(i);
457  }
458 
459  return Quadrature<2> (q_points, weights);
460 }
461 
462 
463 template <int dim>
466  const unsigned int n_times)
467 {
468  std::vector<Point<2> > q_points (q.size());
469  std::vector<double> weights (q.size());
470  for (unsigned int i=0; i<q.size(); ++i)
471  {
472  switch (n_times%4)
473  {
474  case 0:
475  // 0 degree
476  q_points[i][0] = q.point(i)[0];
477  q_points[i][1] = q.point(i)[1];
478  break;
479  case 1:
480  // 90 degree counterclockwise
481  q_points[i][0] = 1.0 - q.point(i)[1];
482  q_points[i][1] = q.point(i)[0];
483  break;
484  case 2:
485  // 180 degree counterclockwise
486  q_points[i][0] = 1.0 - q.point(i)[0];
487  q_points[i][1] = 1.0 - q.point(i)[1];
488  break;
489  case 3:
490  // 270 degree counterclockwise
491  q_points[i][0] = q.point(i)[1];
492  q_points[i][1] = 1.0 - q.point(i)[0];
493  break;
494  }
495 
496  weights[i] = q.weight(i);
497  }
498 
499  return Quadrature<2> (q_points, weights);
500 }
501 
502 
503 template <>
504 void
506  const unsigned int face_no,
507  std::vector<Point<1> > &q_points)
508 {
509  const unsigned int dim=1;
511  AssertDimension (q_points.size(), 1);
512 
513  q_points[0] = Point<dim>((double) face_no);
514 }
515 
516 
517 
518 template <>
519 void
521  const unsigned int face_no,
522  std::vector<Point<2> > &q_points)
523 {
524  const unsigned int dim=2;
526  Assert (q_points.size() == quadrature.size(),
527  ExcDimensionMismatch (q_points.size(), quadrature.size()));
528 
529  for (unsigned int p=0; p<quadrature.size(); ++p)
530  switch (face_no)
531  {
532  case 0:
533  q_points[p] = Point<dim>(0,quadrature.point(p)(0));
534  break;
535  case 1:
536  q_points[p] = Point<dim>(1,quadrature.point(p)(0));
537  break;
538  case 2:
539  q_points[p] = Point<dim>(quadrature.point(p)(0),0);
540  break;
541  case 3:
542  q_points[p] = Point<dim>(quadrature.point(p)(0),1);
543  break;
544  default:
545  Assert (false, ExcInternalError());
546  };
547 }
548 
549 
550 
551 template <>
552 void
554  const unsigned int face_no,
555  std::vector<Point<3> > &q_points)
556 {
557  const unsigned int dim=3;
559  Assert (q_points.size() == quadrature.size(),
560  ExcDimensionMismatch (q_points.size(), quadrature.size()));
561 
562  for (unsigned int p=0; p<quadrature.size(); ++p)
563  switch (face_no)
564  {
565  case 0:
566  q_points[p] = Point<dim>(0,
567  quadrature.point(p)(0),
568  quadrature.point(p)(1));
569  break;
570  case 1:
571  q_points[p] = Point<dim>(1,
572  quadrature.point(p)(0),
573  quadrature.point(p)(1));
574  break;
575  case 2:
576  q_points[p] = Point<dim>(quadrature.point(p)(1),
577  0,
578  quadrature.point(p)(0));
579  break;
580  case 3:
581  q_points[p] = Point<dim>(quadrature.point(p)(1),
582  1,
583  quadrature.point(p)(0));
584  break;
585  case 4:
586  q_points[p] = Point<dim>(quadrature.point(p)(0),
587  quadrature.point(p)(1),
588  0);
589  break;
590  case 5:
591  q_points[p] = Point<dim>(quadrature.point(p)(0),
592  quadrature.point(p)(1),
593  1);
594  break;
595 
596  default:
597  Assert (false, ExcInternalError());
598  };
599 }
600 
601 
602 
603 template <>
604 void
606  const unsigned int face_no,
607  const unsigned int,
608  std::vector<Point<1> > &q_points,
609  const RefinementCase<0> &)
610 {
611  const unsigned int dim=1;
613  AssertDimension (q_points.size(), 1);
614 
615  q_points[0] = Point<dim>((double) face_no);
616 }
617 
618 
619 
620 template <>
621 void
623  const unsigned int face_no,
624  const unsigned int subface_no,
625  std::vector<Point<2> > &q_points,
626  const RefinementCase<1> &)
627 {
628  const unsigned int dim=2;
631 
632  Assert (q_points.size() == quadrature.size(),
633  ExcDimensionMismatch (q_points.size(), quadrature.size()));
634 
635  for (unsigned int p=0; p<quadrature.size(); ++p)
636  switch (face_no)
637  {
638  case 0:
639  switch (subface_no)
640  {
641  case 0:
642  q_points[p] = Point<dim>(0,quadrature.point(p)(0)/2);
643  break;
644  case 1:
645  q_points[p] = Point<dim>(0,quadrature.point(p)(0)/2+0.5);
646  break;
647  default:
648  Assert (false, ExcInternalError());
649  };
650  break;
651  case 1:
652  switch (subface_no)
653  {
654  case 0:
655  q_points[p] = Point<dim>(1,quadrature.point(p)(0)/2);
656  break;
657  case 1:
658  q_points[p] = Point<dim>(1,quadrature.point(p)(0)/2+0.5);
659  break;
660  default:
661  Assert (false, ExcInternalError());
662  };
663  break;
664  case 2:
665  switch (subface_no)
666  {
667  case 0:
668  q_points[p]
669  = Point<dim>(quadrature.point(p)(0)/2,0);
670  break;
671  case 1:
672  q_points[p]
673  = Point<dim>(quadrature.point(p)(0)/2+0.5,0);
674  break;
675  default:
676  Assert (false, ExcInternalError());
677  };
678  break;
679  case 3:
680  switch (subface_no)
681  {
682  case 0:
683  q_points[p] = Point<dim>(quadrature.point(p)(0)/2,1);
684  break;
685  case 1:
686  q_points[p] = Point<dim>(quadrature.point(p)(0)/2+0.5,1);
687  break;
688  default:
689  Assert (false, ExcInternalError());
690  };
691  break;
692 
693  default:
694  Assert (false, ExcInternalError());
695  };
696 }
697 
698 
699 
700 template <>
701 void
703  const unsigned int face_no,
704  const unsigned int subface_no,
705  std::vector<Point<3> > &q_points,
706  const RefinementCase<2> &ref_case)
707 {
708  const unsigned int dim=3;
711  Assert (q_points.size() == quadrature.size(),
712  ExcDimensionMismatch (q_points.size(), quadrature.size()));
713 
714  // one coordinate is at a const value. for
715  // faces 0, 2 and 4 this value is 0.0, for
716  // faces 1, 3 and 5 it is 1.0
717  double const_value=face_no%2;
718  // local 2d coordinates are xi and eta,
719  // global 3d coordinates are x, y and
720  // z. those have to be mapped. the following
721  // indices tell, which global coordinate
722  // (0->x, 1->y, 2->z) corresponds to which
723  // local one
724  unsigned int xi_index = numbers::invalid_unsigned_int,
725  eta_index = numbers::invalid_unsigned_int,
726  const_index = face_no/2;
727  // the xi and eta values have to be scaled
728  // (by factor 0.5 or factor 1.0) depending on
729  // the refinement case and translated (by 0.0
730  // or 0.5) depending on the refinement case
731  // and subface_no.
732  double xi_scale=1.0,
733  eta_scale=1.0,
734  xi_translation=0.0,
735  eta_translation=0.0;
736  // set the index mapping between local and
737  // global coordinates
738  switch (face_no/2)
739  {
740  case 0:
741  xi_index=1;
742  eta_index=2;
743  break;
744  case 1:
745  xi_index=2;
746  eta_index=0;
747  break;
748  case 2:
749  xi_index=0;
750  eta_index=1;
751  break;
752  }
753  // set the scale and translation parameter
754  // for individual subfaces
755  switch ((unsigned char)ref_case)
756  {
757  case RefinementCase<dim-1>::cut_x:
758  xi_scale=0.5;
759  xi_translation = subface_no%2 * 0.5;
760  break;
761  case RefinementCase<dim-1>::cut_y:
762  eta_scale=0.5;
763  eta_translation = subface_no%2 * 0.5;
764  break;
765  case RefinementCase<dim-1>::cut_xy:
766  xi_scale= 0.5;
767  eta_scale=0.5;
768  xi_translation = int(subface_no%2) * 0.5;
769  eta_translation = int(subface_no/2) * 0.5;
770  break;
771  default:
772  Assert(false,ExcInternalError());
773  break;
774  }
775  // finally, compute the scaled, translated,
776  // projected quadrature points
777  for (unsigned int p=0; p<quadrature.size(); ++p)
778  {
779  q_points[p][xi_index] = xi_scale * quadrature.point(p)(0) + xi_translation;
780  q_points[p][eta_index] = eta_scale * quadrature.point(p)(1) + eta_translation;
781  q_points[p][const_index] = const_value;
782  }
783 }
784 
785 
786 template <>
789 {
790  const unsigned int dim = 1;
791 
792  const unsigned int n_points = 1,
794 
795  // first fix quadrature points
796  std::vector<Point<dim> > q_points;
797  q_points.reserve(n_points * n_faces);
798  std::vector <Point<dim> > help(n_points);
799 
800 
801  // project to each face and append
802  // results
803  for (unsigned int face=0; face<n_faces; ++face)
804  {
805  project_to_face(quadrature, face, help);
806  std::copy (help.begin(), help.end(),
807  std::back_inserter (q_points));
808  }
809 
810  // next copy over weights
811  std::vector<double> weights;
812  weights.reserve (n_points * n_faces);
813  for (unsigned int face=0; face<n_faces; ++face)
814  std::copy (quadrature.get_weights().begin(),
815  quadrature.get_weights().end(),
816  std::back_inserter (weights));
817 
818  Assert (q_points.size() == n_points * n_faces,
819  ExcInternalError());
820  Assert (weights.size() == n_points * n_faces,
821  ExcInternalError());
822 
823  return Quadrature<dim>(q_points, weights);
824 }
825 
826 
827 
828 template <>
830 QProjector<2>::project_to_all_faces (const SubQuadrature &quadrature)
831 {
832  const unsigned int dim = 2;
833 
834  const unsigned int n_points = quadrature.size(),
836 
837  // first fix quadrature points
838  std::vector<Point<dim> > q_points;
839  q_points.reserve(n_points * n_faces);
840  std::vector <Point<dim> > help(n_points);
841 
842  // project to each face and append
843  // results
844  for (unsigned int face=0; face<n_faces; ++face)
845  {
846  project_to_face(quadrature, face, help);
847  std::copy (help.begin(), help.end(),
848  std::back_inserter (q_points));
849  }
850 
851  // next copy over weights
852  std::vector<double> weights;
853  weights.reserve (n_points * n_faces);
854  for (unsigned int face=0; face<n_faces; ++face)
855  std::copy (quadrature.get_weights().begin(),
856  quadrature.get_weights().end(),
857  std::back_inserter (weights));
858 
859  Assert (q_points.size() == n_points * n_faces,
860  ExcInternalError());
861  Assert (weights.size() == n_points * n_faces,
862  ExcInternalError());
863 
864  return Quadrature<dim>(q_points, weights);
865 }
866 
867 
868 
869 template <>
871 QProjector<3>::project_to_all_faces (const SubQuadrature &quadrature)
872 {
873  const unsigned int dim = 3;
874 
875  SubQuadrature q_reflected=reflect (quadrature);
876  SubQuadrature q[8]=
877  {
878  quadrature,
879  rotate (quadrature,1),
880  rotate (quadrature,2),
881  rotate (quadrature,3),
882  q_reflected,
883  rotate (q_reflected,3),
884  rotate (q_reflected,2),
885  rotate (q_reflected,1)
886  };
887 
888 
889 
890  const unsigned int n_points = quadrature.size(),
892 
893  // first fix quadrature points
894  std::vector<Point<dim> > q_points;
895  q_points.reserve(n_points * n_faces * 8);
896  std::vector <Point<dim> > help(n_points);
897 
898  std::vector<double> weights;
899  weights.reserve (n_points * n_faces * 8);
900 
901  // do the following for all possible
902  // mutations of a face (mutation==0
903  // corresponds to a face with standard
904  // orientation, no flip and no rotation)
905  for (unsigned int mutation=0; mutation<8; ++mutation)
906  {
907  // project to each face and append
908  // results
909  for (unsigned int face=0; face<n_faces; ++face)
910  {
911  project_to_face(q[mutation], face, help);
912  std::copy (help.begin(), help.end(),
913  std::back_inserter (q_points));
914  }
915 
916  // next copy over weights
917  for (unsigned int face=0; face<n_faces; ++face)
918  std::copy (q[mutation].get_weights().begin(),
919  q[mutation].get_weights().end(),
920  std::back_inserter (weights));
921  }
922 
923 
924  Assert (q_points.size() == n_points * n_faces * 8,
925  ExcInternalError());
926  Assert (weights.size() == n_points * n_faces * 8,
927  ExcInternalError());
928 
929  return Quadrature<dim>(q_points, weights);
930 }
931 
932 
933 
934 template <>
937 {
938  const unsigned int dim = 1;
939 
940  const unsigned int n_points = 1,
942  subfaces_per_face = GeometryInfo<dim>::max_children_per_face;
943 
944  // first fix quadrature points
945  std::vector<Point<dim> > q_points;
946  q_points.reserve (n_points * n_faces * subfaces_per_face);
947  std::vector <Point<dim> > help(n_points);
948 
949  // project to each face and copy
950  // results
951  for (unsigned int face=0; face<n_faces; ++face)
952  for (unsigned int subface=0; subface<subfaces_per_face; ++subface)
953  {
954  project_to_subface(quadrature, face, subface, help);
955  std::copy (help.begin(), help.end(),
956  std::back_inserter (q_points));
957  };
958 
959  // next copy over weights
960  std::vector<double> weights;
961  weights.reserve (n_points * n_faces * subfaces_per_face);
962  for (unsigned int face=0; face<n_faces; ++face)
963  for (unsigned int subface=0; subface<subfaces_per_face; ++subface)
964  std::copy (quadrature.get_weights().begin(),
965  quadrature.get_weights().end(),
966  std::back_inserter (weights));
967 
968  Assert (q_points.size() == n_points * n_faces * subfaces_per_face,
969  ExcInternalError());
970  Assert (weights.size() == n_points * n_faces * subfaces_per_face,
971  ExcInternalError());
972 
973  return Quadrature<dim>(q_points, weights);
974 }
975 
976 
977 
978 template <>
980 QProjector<2>::project_to_all_subfaces (const SubQuadrature &quadrature)
981 {
982  const unsigned int dim = 2;
983 
984  const unsigned int n_points = quadrature.size(),
986  subfaces_per_face = GeometryInfo<dim>::max_children_per_face;
987 
988  // first fix quadrature points
989  std::vector<Point<dim> > q_points;
990  q_points.reserve (n_points * n_faces * subfaces_per_face);
991  std::vector <Point<dim> > help(n_points);
992 
993  // project to each face and copy
994  // results
995  for (unsigned int face=0; face<n_faces; ++face)
996  for (unsigned int subface=0; subface<subfaces_per_face; ++subface)
997  {
998  project_to_subface(quadrature, face, subface, help);
999  std::copy (help.begin(), help.end(),
1000  std::back_inserter (q_points));
1001  };
1002 
1003  // next copy over weights
1004  std::vector<double> weights;
1005  weights.reserve (n_points * n_faces * subfaces_per_face);
1006  for (unsigned int face=0; face<n_faces; ++face)
1007  for (unsigned int subface=0; subface<subfaces_per_face; ++subface)
1008  std::copy (quadrature.get_weights().begin(),
1009  quadrature.get_weights().end(),
1010  std::back_inserter (weights));
1011 
1012  Assert (q_points.size() == n_points * n_faces * subfaces_per_face,
1013  ExcInternalError());
1014  Assert (weights.size() == n_points * n_faces * subfaces_per_face,
1015  ExcInternalError());
1016 
1017  return Quadrature<dim>(q_points, weights);
1018 }
1019 
1020 
1021 
1022 template <>
1024 QProjector<3>::project_to_all_subfaces (const SubQuadrature &quadrature)
1025 {
1026  const unsigned int dim = 3;
1027  SubQuadrature q_reflected=reflect (quadrature);
1028  SubQuadrature q[8]=
1029  {
1030  quadrature,
1031  rotate (quadrature,1),
1032  rotate (quadrature,2),
1033  rotate (quadrature,3),
1034  q_reflected,
1035  rotate (q_reflected,3),
1036  rotate (q_reflected,2),
1037  rotate (q_reflected,1)
1038  };
1039 
1040  const unsigned int n_points = quadrature.size(),
1042  total_subfaces_per_face = 2 + 2 + 4;
1043 
1044  // first fix quadrature points
1045  std::vector<Point<dim> > q_points;
1046  q_points.reserve (n_points * n_faces * total_subfaces_per_face * 8);
1047  std::vector <Point<dim> > help(n_points);
1048 
1049  std::vector<double> weights;
1050  weights.reserve (n_points * n_faces * total_subfaces_per_face * 8);
1051 
1052  // do the following for all possible
1053  // mutations of a face (mutation==0
1054  // corresponds to a face with standard
1055  // orientation, no flip and no rotation)
1056  for (unsigned int mutation=0; mutation<8; ++mutation)
1057  {
1058  // project to each face and copy
1059  // results
1060  for (unsigned int face=0; face<n_faces; ++face)
1061  for (unsigned int ref_case=RefinementCase<dim-1>::cut_xy;
1062  ref_case>=RefinementCase<dim-1>::cut_x;
1063  --ref_case)
1064  for (unsigned int subface=0; subface<GeometryInfo<dim-1>::n_children(RefinementCase<dim-1>(ref_case)); ++subface)
1065  {
1066  project_to_subface(q[mutation], face, subface, help,
1067  RefinementCase<dim-1>(ref_case));
1068  std::copy (help.begin(), help.end(),
1069  std::back_inserter (q_points));
1070  }
1071 
1072  // next copy over weights
1073  for (unsigned int face=0; face<n_faces; ++face)
1074  for (unsigned int ref_case=RefinementCase<dim-1>::cut_xy;
1075  ref_case>=RefinementCase<dim-1>::cut_x;
1076  --ref_case)
1077  for (unsigned int subface=0; subface<GeometryInfo<dim-1>::n_children(RefinementCase<dim-1>(ref_case)); ++subface)
1078  std::copy (q[mutation].get_weights().begin(),
1079  q[mutation].get_weights().end(),
1080  std::back_inserter (weights));
1081  }
1082 
1083  Assert (q_points.size() == n_points * n_faces * total_subfaces_per_face * 8,
1084  ExcInternalError());
1085  Assert (weights.size() == n_points * n_faces * total_subfaces_per_face * 8,
1086  ExcInternalError());
1087 
1088  return Quadrature<dim>(q_points, weights);
1089 }
1090 
1091 
1092 
1093 // This function is not used in the library
1094 template <int dim>
1097  const unsigned int child_no)
1098 {
1101 
1102  const unsigned int n_q_points = quadrature.size();
1103 
1104  std::vector<Point<dim> > q_points(n_q_points);
1105  for (unsigned int i=0; i<n_q_points; ++i)
1107  quadrature.point(i), child_no);
1108 
1109  // for the weights, things are
1110  // equally simple: copy them and
1111  // scale them
1112  std::vector<double> weights = quadrature.get_weights ();
1113  for (unsigned int i=0; i<n_q_points; ++i)
1114  weights[i] *= (1./GeometryInfo<dim>::max_children_per_cell);
1115 
1116  return Quadrature<dim> (q_points, weights);
1117 }
1118 
1119 
1120 template <int dim>
1123 {
1124  const unsigned int n_points = quadrature.size(),
1126 
1127  std::vector<Point<dim> > q_points(n_points * n_children);
1128  std::vector<double> weights(n_points * n_children);
1129 
1130  // project to each child and copy
1131  // results
1132  for (unsigned int child=0; child<n_children; ++child)
1133  {
1134  Quadrature<dim> help = project_to_child(quadrature, child);
1135  for (unsigned int i=0; i<n_points; ++i)
1136  {
1137  q_points[child*n_points+i] = help.point(i);
1138  weights[child*n_points+i] = help.weight(i);
1139  }
1140  }
1141  return Quadrature<dim>(q_points, weights);
1142 }
1143 
1144 
1145 
1146 
1147 template <int dim>
1150  const Quadrature<1> &quadrature,
1151  const Point<dim> &p1,
1152  const Point<dim> &p2)
1153 {
1154  const unsigned int n = quadrature.size();
1155  std::vector<Point<dim> > points(n);
1156  std::vector<double> weights(n);
1157  const double length = p1.distance(p2);
1158 
1159  for (unsigned int k=0; k<n; ++k)
1160  {
1161  const double alpha = quadrature.point(k)(0);
1162  points[k] = alpha * p2;
1163  points[k] += (1.-alpha) * p1;
1164  weights[k] = length * quadrature.weight(k);
1165  }
1166  return Quadrature<dim> (points, weights);
1167 }
1168 
1169 
1170 
1171 template <int dim>
1174 face (const unsigned int face_no,
1175  const bool face_orientation,
1176  const bool face_flip,
1177  const bool face_rotation,
1178  const unsigned int n_quadrature_points)
1179 {
1181  ExcInternalError());
1182 
1183  switch (dim)
1184  {
1185  case 1:
1186  case 2:
1187  return face_no * n_quadrature_points;
1188 
1189 
1190  case 3:
1191  {
1192  // in 3d, we have to account for faces that
1193  // have non-standard face orientation, flip
1194  // and rotation. thus, we have to store
1195  // _eight_ data sets per face or subface
1196 
1197  // set up a table with the according offsets
1198  // for non-standard orientation, first index:
1199  // face_orientation (standard true=1), second
1200  // index: face_flip (standard false=0), third
1201  // index: face_rotation (standard false=0)
1202  //
1203  // note, that normally we should use the
1204  // obvious offsets 0,1,2,3,4,5,6,7. However,
1205  // prior to the changes enabling flipped and
1206  // rotated faces, in many places of the
1207  // library the convention was used, that the
1208  // first dataset with offset 0 corresponds to
1209  // a face in standard orientation. therefore
1210  // we use the offsets 4,5,6,7,0,1,2,3 here to
1211  // stick to that (implicit) convention
1212  static const unsigned int offset[2][2][2]=
1213  {
1214  { {4*GeometryInfo<dim>::faces_per_cell, 5*GeometryInfo<dim>::faces_per_cell}, // face_orientation=false; face_flip=false; face_rotation=false and true
1216  }, // face_orientation=false; face_flip=true; face_rotation=false and true
1217  { {0*GeometryInfo<dim>::faces_per_cell, 1*GeometryInfo<dim>::faces_per_cell}, // face_orientation=true; face_flip=false; face_rotation=false and true
1219  }
1220  }; // face_orientation=true; face_flip=true; face_rotation=false and true
1221 
1222  return ((face_no
1223  + offset[face_orientation][face_flip][face_rotation])
1224  * n_quadrature_points);
1225  }
1226 
1227  default:
1228  Assert (false, ExcInternalError());
1229  }
1231 }
1232 
1233 
1234 
1235 template <>
1238 subface (const unsigned int face_no,
1239  const unsigned int subface_no,
1240  const bool,
1241  const bool,
1242  const bool,
1243  const unsigned int n_quadrature_points,
1245 {
1247  ExcInternalError());
1249  ExcInternalError());
1250 
1251  return ((face_no * GeometryInfo<1>::max_children_per_face +
1252  subface_no)
1253  * n_quadrature_points);
1254 }
1255 
1256 
1257 
1258 template <>
1261 subface (const unsigned int face_no,
1262  const unsigned int subface_no,
1263  const bool,
1264  const bool,
1265  const bool,
1266  const unsigned int n_quadrature_points,
1268 {
1270  ExcInternalError());
1272  ExcInternalError());
1273 
1274  return ((face_no * GeometryInfo<2>::max_children_per_face +
1275  subface_no)
1276  * n_quadrature_points);
1277 }
1278 
1279 
1280 template <>
1283 subface (const unsigned int face_no,
1284  const unsigned int subface_no,
1285  const bool face_orientation,
1286  const bool face_flip,
1287  const bool face_rotation,
1288  const unsigned int n_quadrature_points,
1289  const internal::SubfaceCase<3> ref_case)
1290 {
1291  const unsigned int dim = 3;
1292 
1294  ExcInternalError());
1296  ExcInternalError());
1297 
1298  // As the quadrature points created by
1299  // QProjector are on subfaces in their
1300  // "standard location" we have to use a
1301  // permutation of the equivalent subface
1302  // number in order to respect face
1303  // orientation, flip and rotation. The
1304  // information we need here is exactly the
1305  // same as the
1306  // GeometryInfo<3>::child_cell_on_face info
1307  // for the bottom face (face 4) of a hex, as
1308  // on this the RefineCase of the cell matches
1309  // that of the face and the subfaces are
1310  // numbered in the same way as the child
1311  // cells.
1312 
1313  // in 3d, we have to account for faces that
1314  // have non-standard face orientation, flip
1315  // and rotation. thus, we have to store
1316  // _eight_ data sets per face or subface
1317  // already for the isotropic
1318  // case. Additionally, we have three
1319  // different refinement cases, resulting in
1320  // <tt>4 + 2 + 2 = 8</tt> different subfaces
1321  // for each face.
1322  const unsigned int total_subfaces_per_face=8;
1323 
1324  // set up a table with the according offsets
1325  // for non-standard orientation, first index:
1326  // face_orientation (standard true=1), second
1327  // index: face_flip (standard false=0), third
1328  // index: face_rotation (standard false=0)
1329  //
1330  // note, that normally we should use the
1331  // obvious offsets 0,1,2,3,4,5,6,7. However,
1332  // prior to the changes enabling flipped and
1333  // rotated faces, in many places of the
1334  // library the convention was used, that the
1335  // first dataset with offset 0 corresponds to
1336  // a face in standard orientation. therefore
1337  // we use the offsets 4,5,6,7,0,1,2,3 here to
1338  // stick to that (implicit) convention
1339  static const unsigned int orientation_offset[2][2][2]=
1340  {
1341  {
1342  // face_orientation=false; face_flip=false; face_rotation=false and true
1343  {
1344  4*GeometryInfo<dim>::faces_per_cell*total_subfaces_per_face,
1345  5*GeometryInfo<dim>::faces_per_cell *total_subfaces_per_face
1346  },
1347  // face_orientation=false; face_flip=true; face_rotation=false and true
1348  {
1349  6*GeometryInfo<dim>::faces_per_cell*total_subfaces_per_face,
1350  7*GeometryInfo<dim>::faces_per_cell *total_subfaces_per_face
1351  }
1352  },
1353  {
1354  // face_orientation=true; face_flip=false; face_rotation=false and true
1355  {
1356  0*GeometryInfo<dim>::faces_per_cell*total_subfaces_per_face,
1357  1*GeometryInfo<dim>::faces_per_cell *total_subfaces_per_face
1358  },
1359  // face_orientation=true; face_flip=true; face_rotation=false and true
1360  {
1361  2*GeometryInfo<dim>::faces_per_cell*total_subfaces_per_face,
1362  3*GeometryInfo<dim>::faces_per_cell *total_subfaces_per_face
1363  }
1364  }
1365  };
1366 
1367  // set up a table with the offsets for a
1368  // given refinement case respecting the
1369  // corresponding number of subfaces. the
1370  // index corresponds to (RefineCase::Type - 1)
1371 
1372  // note, that normally we should use the
1373  // obvious offsets 0,2,6. However, prior to
1374  // the implementation of anisotropic
1375  // refinement, in many places of the library
1376  // the convention was used, that the first
1377  // dataset with offset 0 corresponds to a
1378  // standard (isotropic) face
1379  // refinement. therefore we use the offsets
1380  // 6,4,0 here to stick to that (implicit)
1381  // convention
1382  static const unsigned int ref_case_offset[3]=
1383  {
1384  6, //cut_x
1385  4, //cut_y
1386  0 //cut_xy
1387  };
1388 
1389 
1390  // for each subface of a given FaceRefineCase
1391  // there is a corresponding equivalent
1392  // subface number of one of the "standard"
1393  // RefineCases (cut_x, cut_y, cut_xy). Map
1394  // the given values to those equivalent
1395  // ones.
1396 
1397  // first, define an invalid number
1398  static const unsigned int e = numbers::invalid_unsigned_int;
1399 
1400  static const RefinementCase<dim-1>
1402  =
1403  {
1404  // case_none. there should be only
1405  // invalid values here. However, as
1406  // this function is also called (in
1407  // tests) for cells which have no
1408  // refined faces, use isotropic
1409  // refinement instead
1410  {
1411  RefinementCase<dim-1>::cut_xy,
1412  RefinementCase<dim-1>::cut_xy,
1413  RefinementCase<dim-1>::cut_xy,
1414  RefinementCase<dim-1>::cut_xy
1415  },
1416  // case_x
1417  {
1418  RefinementCase<dim-1>::cut_x,
1419  RefinementCase<dim-1>::cut_x,
1420  RefinementCase<dim-1>::no_refinement,
1421  RefinementCase<dim-1>::no_refinement
1422  },
1423  // case_x1y
1424  {
1425  RefinementCase<dim-1>::cut_xy,
1426  RefinementCase<dim-1>::cut_xy,
1427  RefinementCase<dim-1>::cut_x,
1428  RefinementCase<dim-1>::no_refinement
1429  },
1430  // case_x2y
1431  {
1432  RefinementCase<dim-1>::cut_x,
1433  RefinementCase<dim-1>::cut_xy,
1434  RefinementCase<dim-1>::cut_xy,
1435  RefinementCase<dim-1>::no_refinement
1436  },
1437  // case_x1y2y
1438  {
1439  RefinementCase<dim-1>::cut_xy,
1440  RefinementCase<dim-1>::cut_xy,
1441  RefinementCase<dim-1>::cut_xy,
1442  RefinementCase<dim-1>::cut_xy
1443  },
1444  // case_y
1445  {
1446  RefinementCase<dim-1>::cut_y,
1447  RefinementCase<dim-1>::cut_y,
1448  RefinementCase<dim-1>::no_refinement,
1449  RefinementCase<dim-1>::no_refinement
1450  },
1451  // case_y1x
1452  {
1453  RefinementCase<dim-1>::cut_xy,
1454  RefinementCase<dim-1>::cut_xy,
1455  RefinementCase<dim-1>::cut_y,
1456  RefinementCase<dim-1>::no_refinement
1457  },
1458  // case_y2x
1459  {
1460  RefinementCase<dim-1>::cut_y,
1461  RefinementCase<dim-1>::cut_xy,
1462  RefinementCase<dim-1>::cut_xy,
1463  RefinementCase<dim-1>::no_refinement
1464  },
1465  // case_y1x2x
1466  {
1467  RefinementCase<dim-1>::cut_xy,
1468  RefinementCase<dim-1>::cut_xy,
1469  RefinementCase<dim-1>::cut_xy,
1470  RefinementCase<dim-1>::cut_xy
1471  },
1472  // case_xy (case_isotropic)
1473  {
1474  RefinementCase<dim-1>::cut_xy,
1475  RefinementCase<dim-1>::cut_xy,
1476  RefinementCase<dim-1>::cut_xy,
1477  RefinementCase<dim-1>::cut_xy
1478  }
1479  };
1480 
1481  static const unsigned int
1483  =
1484  {
1485  // case_none, see above
1486  {0,1,2,3},
1487  // case_x
1488  {0,1,e,e},
1489  // case_x1y
1490  {0,2,1,e},
1491  // case_x2y
1492  {0,1,3,e},
1493  // case_x1y2y
1494  {0,2,1,3},
1495  // case_y
1496  {0,1,e,e},
1497  // case_y1x
1498  {0,1,1,e},
1499  // case_y2x
1500  {0,2,3,e},
1501  // case_y1x2x
1502  {0,1,2,3},
1503  // case_xy (case_isotropic)
1504  {0,1,2,3}
1505  };
1506 
1507  // If face-orientation or face_rotation are
1508  // non-standard, cut_x and cut_y have to be
1509  // exchanged.
1510  static const RefinementCase<dim-1> ref_case_permutation[4]
1511  = {RefinementCase<dim-1>::no_refinement,
1512  RefinementCase<dim-1>::cut_y,
1513  RefinementCase<dim-1>::cut_x,
1514  RefinementCase<dim-1>::cut_xy
1515  };
1516 
1517  // set a corresponding (equivalent)
1518  // RefineCase and subface number
1519  const RefinementCase<dim-1> equ_ref_case=equivalent_refine_case[ref_case][subface_no];
1520  const unsigned int equ_subface_no=equivalent_subface_number[ref_case][subface_no];
1521  // make sure, that we got a valid subface and RefineCase
1523  Assert(equ_subface_no!=e, ExcInternalError());
1524  // now, finally respect non-standard faces
1525  const RefinementCase<dim-1>
1526  final_ref_case = (face_orientation==face_rotation
1527  ?
1528  ref_case_permutation[equ_ref_case]
1529  :
1530  equ_ref_case);
1531 
1532  // what we have now is the number of
1533  // the subface in the natural
1534  // orientation of the *face*. what we
1535  // need to know is the number of the
1536  // subface concerning the standard face
1537  // orientation as seen from the *cell*.
1538 
1539  // this mapping is not trivial, but we
1540  // have done exactly this stuff in the
1541  // child_cell_on_face function. in
1542  // order to reduce the amount of code
1543  // as well as to make maintaining the
1544  // functionality easier we want to
1545  // reuse that information. So we note
1546  // that on the bottom face (face 4) of
1547  // a hex cell the local x and y
1548  // coordinates of the face and the cell
1549  // coincide, thus also the refinement
1550  // case of the face corresponds to the
1551  // refinement case of the cell
1552  // (ignoring cell refinement along the
1553  // z direction). Using this knowledge
1554  // we can (ab)use the
1555  // child_cell_on_face function to do
1556  // exactly the transformation we are in
1557  // need of now
1558  const unsigned int
1559  final_subface_no = GeometryInfo<dim>::child_cell_on_face(RefinementCase<dim>(final_ref_case),
1560  4,
1561  equ_subface_no,
1562  face_orientation,
1563  face_flip,
1564  face_rotation,
1565  equ_ref_case);
1566 
1567  return (((face_no * total_subfaces_per_face
1568  + ref_case_offset[final_ref_case-1]
1569  + final_subface_no)
1570  + orientation_offset[face_orientation][face_flip][face_rotation]
1571  )
1572  * n_quadrature_points);
1573 }
1574 
1575 
1576 
1577 template <int dim>
1580  const unsigned int face_no)
1581 {
1582  std::vector<Point<dim> > points(quadrature.size());
1583  project_to_face(quadrature, face_no, points);
1584  return Quadrature<dim>(points, quadrature.get_weights());
1585 }
1586 
1587 
1588 template <int dim>
1591  const unsigned int face_no,
1592  const unsigned int subface_no,
1593  const RefinementCase<dim-1> &ref_case)
1594 {
1595  std::vector<Point<dim> > points(quadrature.size());
1596  project_to_subface(quadrature, face_no, subface_no, points, ref_case);
1597  return Quadrature<dim>(points, quadrature.get_weights());
1598 }
1599 
1600 
1601 // ------------------------------------------------------------ //
1602 
1603 namespace internal
1604 {
1605  namespace QIteratedImplementation
1606  {
1607  namespace
1608  {
1609  bool uses_both_endpoints(const Quadrature<1> &base_quadrature)
1610  {
1611  bool at_left = false,
1612  at_right = false;
1613  for (unsigned int i = 0; i<base_quadrature.size(); ++i)
1614  {
1615  if (base_quadrature.point(i) == Point<1>(0.0))
1616  at_left = true;
1617  if (base_quadrature.point(i) == Point<1>(1.0))
1618  at_right = true;
1619  }
1620 
1621  return (at_left && at_right);
1622  }
1623  }
1624  }
1625 }
1626 
1627 // template <>
1628 // void
1629 // QIterated<1>::fill(Quadrature<1>& dst,
1630 // const Quadrature<1> &base_quadrature,
1631 // const unsigned int n_copies)
1632 // {
1633 // Assert (n_copies > 0, ExcZero());
1634 // Assert (base_quadrature.size() > 0, ExcZero());
1635 
1636 // const unsigned int np =
1637 // uses_both_endpoints(base_quadrature)
1638 // ? (base_quadrature.size()-1) * n_copies + 1
1639 // : base_quadrature.size() * n_copies;
1640 
1641 // dst.quadrature_points.resize(np);
1642 // dst.weights.resize(np);
1643 
1644 // if (!uses_both_endpoints(base_quadrature))
1645 // // we don't have to skip some
1646 // // points in order to get a
1647 // // reasonable quadrature formula
1648 // {
1649 // unsigned int next_point = 0;
1650 // for (unsigned int copy=0; copy<n_copies; ++copy)
1651 // for (unsigned int q_point=0; q_point<base_quadrature.size(); ++q_point)
1652 // {
1653 // dst.quadrature_points[next_point](0)
1654 // = (copy + base_quadrature.point(q_point)(0)) / n_copies;
1655 // dst.weights[next_point]
1656 // = base_quadrature.weight(q_point) / n_copies;
1657 // ++next_point;
1658 // };
1659 // }
1660 // else
1661 // // skip doubly available points
1662 // {
1663 // unsigned int next_point = 0;
1664 
1665 // // first find out the weights of
1666 // // the left and the right boundary
1667 // // points. note that these usually
1668 // // are but need not necessarily be
1669 // // the same
1670 // double double_point_weight = 0;
1671 // unsigned int n_end_points = 0;
1672 // for (unsigned int i=0; i<base_quadrature.size(); ++i)
1673 // // add up the weight if this
1674 // // is an endpoint
1675 // if ((base_quadrature.point(i)(0) == 0.) ||
1676 // (base_quadrature.point(i)(0) == 1.))
1677 // {
1678 // double_point_weight += base_quadrature.weight(i);
1679 // ++n_end_points;
1680 // };
1681 // // scale the weight correctly
1682 // double_point_weight /= n_copies;
1683 
1684 // // make sure the base quadrature formula
1685 // // has only one quadrature point
1686 // // per end point
1687 // Assert (n_end_points == 2, ExcInvalidQuadratureFormula());
1688 
1689 
1690 // for (unsigned int copy=0; copy<n_copies; ++copy)
1691 // for (unsigned int q_point=0; q_point<base_quadrature.size(); ++q_point)
1692 // {
1693 // // skip the left point of
1694 // // this copy since we
1695 // // have already entered
1696 // // it the last time
1697 // if ((copy > 0) &&
1698 // (base_quadrature.point(q_point)(0) == 0.))
1699 // continue;
1700 
1701 // dst.quadrature_points[next_point](0)
1702 // = (copy+base_quadrature.point(q_point)(0)) / n_copies;
1703 
1704 // // if this is the
1705 // // rightmost point of one
1706 // // of the non-last
1707 // // copies: give it the
1708 // // double weight
1709 // if ((copy != n_copies-1) &&
1710 // (base_quadrature.point(q_point)(0) == 1.))
1711 // dst.weights[next_point] = double_point_weight;
1712 // else
1713 // dst.weights[next_point] = base_quadrature.weight(q_point) /
1714 // n_copies;
1715 
1716 // ++next_point;
1717 // };
1718 // };
1719 
1720 // #if DEBUG
1721 // double sum_of_weights = 0;
1722 // for (unsigned int i=0; i<dst.size(); ++i)
1723 // sum_of_weights += dst.weight(i);
1724 // Assert (std::fabs(sum_of_weights-1) < 1e-15,
1725 // ExcInternalError());
1726 // #endif
1727 
1728 // }
1729 
1730 
1731 template <>
1733  const unsigned int )
1734  :
1735  Quadrature<0>()
1736 {
1737  Assert (false, ExcNotImplemented());
1738 }
1739 
1740 
1741 
1742 template <>
1743 QIterated<1>::QIterated (const Quadrature<1> &base_quadrature,
1744  const unsigned int n_copies)
1745  :
1746  Quadrature<1> (internal::QIteratedImplementation::uses_both_endpoints(base_quadrature) ?
1747  (base_quadrature.size()-1) * n_copies + 1 :
1748  base_quadrature.size() * n_copies)
1749 {
1750 // fill(*this, base_quadrature, n_copies);
1751  Assert (base_quadrature.size() > 0, ExcNotInitialized());
1752  Assert (n_copies > 0, ExcZero());
1753 
1754  if (!internal::QIteratedImplementation::uses_both_endpoints(base_quadrature))
1755  // we don't have to skip some
1756  // points in order to get a
1757  // reasonable quadrature formula
1758  {
1759  unsigned int next_point = 0;
1760  for (unsigned int copy=0; copy<n_copies; ++copy)
1761  for (unsigned int q_point=0; q_point<base_quadrature.size(); ++q_point)
1762  {
1763  this->quadrature_points[next_point]
1764  = Point<1>(base_quadrature.point(q_point)(0) / n_copies
1765  +
1766  (1.0*copy)/n_copies);
1767  this->weights[next_point]
1768  = base_quadrature.weight(q_point) / n_copies;
1769 
1770  ++next_point;
1771  };
1772  }
1773  else
1774  // skip doubly available points
1775  {
1776  unsigned int next_point = 0;
1777 
1778  // first find out the weights of
1779  // the left and the right boundary
1780  // points. note that these usually
1781  // are but need not necessarily be
1782  // the same
1783  double double_point_weight = 0;
1784  unsigned int n_end_points = 0;
1785  for (unsigned int i=0; i<base_quadrature.size(); ++i)
1786  // add up the weight if this
1787  // is an endpoint
1788  if ((base_quadrature.point(i) == Point<1>(0.0)) ||
1789  (base_quadrature.point(i) == Point<1>(1.0)))
1790  {
1791  double_point_weight += base_quadrature.weight(i);
1792  ++n_end_points;
1793  };
1794  // scale the weight correctly
1795  double_point_weight /= n_copies;
1796 
1797  // make sure the base quadrature formula
1798  // has only one quadrature point
1799  // per end point
1800  Assert (n_end_points == 2, ExcInvalidQuadratureFormula());
1801 
1802 
1803  for (unsigned int copy=0; copy<n_copies; ++copy)
1804  for (unsigned int q_point=0; q_point<base_quadrature.size(); ++q_point)
1805  {
1806  // skip the left point of
1807  // this copy since we
1808  // have already entered
1809  // it the last time
1810  if ((copy > 0) &&
1811  (base_quadrature.point(q_point) == Point<1>(0.0)))
1812  continue;
1813 
1814  this->quadrature_points[next_point]
1815  = Point<1>(base_quadrature.point(q_point)(0) / n_copies
1816  +
1817  (1.0*copy)/n_copies);
1818 
1819  // if this is the
1820  // rightmost point of one
1821  // of the non-last
1822  // copies: give it the
1823  // double weight
1824  if ((copy != n_copies-1) &&
1825  (base_quadrature.point(q_point) == Point<1>(1.0)))
1826  this->weights[next_point] = double_point_weight;
1827  else
1828  this->weights[next_point] = base_quadrature.weight(q_point) /
1829  n_copies;
1830 
1831  ++next_point;
1832  };
1833  };
1834 
1835 #if DEBUG
1836  double sum_of_weights = 0;
1837  for (unsigned int i=0; i<this->size(); ++i)
1838  sum_of_weights += this->weight(i);
1839  Assert (std::fabs(sum_of_weights-1) < 1e-13,
1840  ExcInternalError());
1841 #endif
1842 }
1843 
1844 
1845 // template <int dim>
1846 // void
1847 // QIterated<dim>::fill(Quadrature<dim>&, const Quadrature<1>&, unsigned int)
1848 // {
1849 // Assert(false, ExcNotImplemented());
1850 // }
1851 
1852 
1853 // construct higher dimensional quadrature formula by tensor product
1854 // of lower dimensional iterated quadrature formulae
1855 template <int dim>
1857  const unsigned int N)
1858  :
1859  Quadrature<dim> (QIterated<dim-1>(base_quadrature, N),
1860  QIterated<1>(base_quadrature, N))
1861 {}
1862 
1863 
1864 
1865 // explicit instantiations; note: we need them all for all dimensions
1866 template class Quadrature<0>;
1867 template class Quadrature<1>;
1868 template class Quadrature<2>;
1869 template class Quadrature<3>;
1870 template class QAnisotropic<1>;
1871 template class QAnisotropic<2>;
1872 template class QAnisotropic<3>;
1873 template class QIterated<1>;
1874 template class QIterated<2>;
1875 template class QIterated<3>;
1876 template class QProjector<1>;
1877 template class QProjector<2>;
1878 template class QProjector<3>;
1879 
1880 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
Definition: types.h:173
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
Definition: quadrature.cc:1096
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
std::vector< double > weights
Definition: quadrature.h:258
const std::vector< double > & get_weights() const
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Quadrature(const unsigned int n_quadrature_points=0)
Definition: quadrature.cc:41
static Quadrature< 2 > rotate(const Quadrature< 2 > &q, const unsigned int n_times)
Definition: quadrature.cc:465
#define AssertIndexRange(index, range)
Definition: exceptions.h:1284
const Point< dim > & point(const unsigned int i) const
STL namespace.
std::conditional< dim==1, std::array< Quadrature< 1 >, dim >, const std::array< Quadrature< 1 >, dim > & >::type get_tensor_basis() const
Definition: quadrature.cc:332
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
Definition: point.h:104
Quadrature & operator=(const Quadrature< dim > &)
Definition: quadrature.cc:292
static Quadrature< dim > project_to_all_subfaces(const SubQuadrature &quadrature)
static Quadrature< dim > project_to_line(const Quadrature< 1 > &quadrature, const Point< dim > &p1, const Point< dim > &p2)
Definition: quadrature.cc:1149
QAnisotropic(const Quadrature< 1 > &qx)
Definition: quadrature.cc:359
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
void initialize(const std::vector< Point< dim > > &points, const std::vector< double > &weights)
Definition: quadrature.cc:52
T sum(const T &t, const MPI_Comm &mpi_communicator)
static unsigned int child_cell_on_face(const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim-1 > &face_refinement_case=RefinementCase< dim-1 >::isotropic_refinement)
#define Assert(cond, exc)
Definition: exceptions.h:1142
static void project_to_subface(const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim > > &q_points, const RefinementCase< dim-1 > &ref_case=RefinementCase< dim-1 >::isotropic_refinement)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::size_t memory_consumption() const
Definition: quadrature.cc:322
void rotate(const double angle, Triangulation< 2 > &triangulation)
Definition: grid_tools.cc:661
std::unique_ptr< std::array< Quadrature< 1 >, dim > > tensor_basis
Definition: quadrature.h:273
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
bool is_tensor_product_flag
Definition: quadrature.h:267
std::vector< Point< dim > > quadrature_points
Definition: quadrature.h:252
unsigned int size() const
static Quadrature< 2 > reflect(const Quadrature< 2 > &q)
Definition: quadrature.cc:447
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
Definition: cuda.h:31
bool operator==(const Quadrature< dim > &p) const
Definition: quadrature.cc:311
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
QIterated(const Quadrature< 1 > &base_quadrature, const unsigned int n_copies)
Definition: quadrature.cc:1856
static DataSetDescriptor subface(const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points, const internal::SubfaceCase< dim > ref_case=internal::SubfaceCase< dim >::case_isotropic)
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcZero()
static Quadrature< dim > project_to_all_children(const Quadrature< dim > &quadrature)
Definition: quadrature.cc:1122
double weight(const unsigned int i) const
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcInternalError()
static DataSetDescriptor face(const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
Definition: quadrature.cc:1174