Flow123d  jenkins-Flow123d-windows32-release-multijob-51
quadrature_lib.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Definitions of particular quadrature rules on simplices.
27  * @author Jan Stebel
28  */
29 
30 #include "system/global_defs.h"
31 #include "system/system.hh"
33 
34 
35 #define FLOAT double
36 #define SHORT int
37 #define _F(n) n
38 #include "quad.c"
39 #undef FLOAT
40 #undef SHORT
41 #undef _F
42 
43 
44 using namespace arma;
45 
46 
47 template<unsigned int dim>
48 QGauss<dim>::QGauss(const unsigned int order)
49 {
50  typedef QUAD* pQUAD;
51 
52  static const pQUAD
53  q1d[] = { QUAD_1D_P1,
61  },
62  q2d[] = { QUAD_2D_P1,
69  },
70  q3d[] = { QUAD_3D_P1,
75  };
76  static const double unit_cell_volume[] = { 1, 1, 0.5, 1./6 };
77  const pQUAD *q;
78  unsigned int nquads;
79  vec::fixed<dim> p;
80 
81  switch (dim)
82  {
83  case 0:
84  this->quadrature_points.push_back(p);
85  this->weights.push_back(1);
86  return;
87  case 1:
88  q = q1d;
89  nquads = sizeof(q1d) / sizeof(pQUAD);
90  break;
91  case 2:
92  q = q2d;
93  nquads = sizeof(q2d) / sizeof(pQUAD);
94  break;
95  case 3:
96  q = q3d;
97  nquads = sizeof(q3d) / sizeof(pQUAD);
98  break;
99  }
100 
101  ASSERT(order < nquads, "Quadrature of given order is not implemented.");
102 
103  for (int i=0; i<q[order]->npoints; i++)
104  {
105  for (unsigned int j=0; j<dim; j++)
106  p(j) = q[order]->points[i*(dim+1)+j];
107 
108  this->quadrature_points.push_back(p);
109  // The weights must be adjusted according to the volume of the unit cell:
110  // 1D cell: volume 1
111  // 2D cell: volume 1/2
112  // 3D cell: volume 1/6
113  this->weights.push_back(q[order]->weights[i]*unit_cell_volume[dim]);
114  }
115 }
116 
117 
118 template class QGauss<0>;
119 template class QGauss<1>;
120 template class QGauss<2>;
121 template class QGauss<3>;
#define QUAD_2D_P20
Definition: quad.h:115
#define QUAD_1D_P14
Definition: quad.h:63
#define QUAD_3D_P5
Definition: quad.h:129
#define QUAD_1D_P10
Definition: quad.h:57
#define QUAD_3D_P7
Definition: quad.h:133
#define QUAD_3D_P6
Definition: quad.h:131
#define QUAD_1D_P15
Definition: quad.h:64
#define QUAD_3D_P8
Definition: quad.h:135
#define QUAD_2D_P9
Definition: quad.h:93
#define QUAD_1D_P11
Definition: quad.h:58
#define QUAD_1D_P2
Definition: quad.h:45
#define QUAD_2D_P21
Definition: quad.h:117
#define QUAD_2D_P17
Definition: quad.h:109
#define QUAD_2D_P10
Definition: quad.h:95
Definition: quad.h:29
#define QUAD_1D_P18
Definition: quad.h:69
#define QUAD_2D_P19
Definition: quad.h:113
#define QUAD_2D_P2
Definition: quad.h:79
#define QUAD_1D_P5
Definition: quad.h:49
#define QUAD_1D_P21
Definition: quad.h:73
#define QUAD_3D_P13
Definition: quad.h:145
Symmetric Gauss-Legendre quadrature formulae on simplices.
#define QUAD_2D_P14
Definition: quad.h:103
#define QUAD_2D_P6
Definition: quad.h:87
#define QUAD_2D_P12
Definition: quad.h:99
Global macros to enhance readability and debugging, general constants.
#define QUAD_3D_P2
Definition: quad.h:123
#define ASSERT(...)
Definition: global_defs.h:121
#define QUAD_2D_P8
Definition: quad.h:91
#define QUAD_2D_P16
Definition: quad.h:107
#define QUAD_2D_P1
Definition: quad.h:77
#define QUAD_1D_P13
Definition: quad.h:61
#define QUAD_1D_P9
Definition: quad.h:55
#define QUAD_2D_P5
Definition: quad.h:85
#define QUAD_2D_P18
Definition: quad.h:111
#define QUAD_2D_P4
Definition: quad.h:83
#define QUAD_1D_P16
Definition: quad.h:66
#define QUAD_1D_P8
Definition: quad.h:54
#define QUAD_3D_P12
Definition: quad.h:143
#define QUAD_1D_P17
Definition: quad.h:67
#define QUAD_2D_P13
Definition: quad.h:101
#define QUAD_3D_P14
Definition: quad.h:147
#define QUAD_3D_P11
Definition: quad.h:141
#define QUAD_1D_P3
Definition: quad.h:46
#define QUAD_3D_P10
Definition: quad.h:139
Definitions of particular quadrature rules on simplices.
#define QUAD_1D_P19
Definition: quad.h:70
#define QUAD_3D_P1
Definition: quad.h:121
#define QUAD_2D_P7
Definition: quad.h:89
#define QUAD_1D_P4
Definition: quad.h:48
#define QUAD_1D_P1
Definition: quad.h:43
#define QUAD_1D_P12
Definition: quad.h:60
#define QUAD_2D_P11
Definition: quad.h:97
#define QUAD_2D_P15
Definition: quad.h:105
#define QUAD_1D_P20
Definition: quad.h:72
#define QUAD_2D_P3
Definition: quad.h:81
#define QUAD_1D_P7
Definition: quad.h:52
QGauss(const unsigned int order)
Create a formula of given order.
#define QUAD_3D_P4
Definition: quad.h:127
#define QUAD_3D_P9
Definition: quad.h:137
#define QUAD_3D_P3
Definition: quad.h:125
#define QUAD_1D_P6
Definition: quad.h:51