~ubuntu-branches/ubuntu/trusty/rheolef/trusty

1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
///
/// This file is part of Rheolef.
///
/// Copyright (C) 2000-2009 Pierre Saramito <Pierre.Saramito@imag.fr>
///
/// Rheolef is free software; you can redistribute it and/or modify
/// it under the terms of the GNU General Public License as published by
/// the Free Software Foundation; either version 2 of the License, or
/// (at your option) any later version.
///
/// Rheolef is distributed in the hope that it will be useful,
/// but WITHOUT ANY WARRANTY; without even the implied warranty of
/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
/// GNU General Public License for more details.
///
/// You should have received a copy of the GNU General Public License
/// along with Rheolef; if not, write to the Free Software
/// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
///
/// =========================================================================
#include "rheolef/quadrature.h"
#include "rheolef/gauss_jacobi.h"
using namespace rheolef;
using namespace std;

/* ------------------------------------------
 * Gauss quadrature formulae on the triangle
 * order r : exact for all polynoms of P_r
 *
 * References for optimal formulae :
 *   1)  G. Dhatt and G. Touzot,
 *       Une presentation de la methode des elements finis,
 *       Maloine editeur, Paris
 *       Deuxieme edition,
 *       1984,
 *       page 297
 *   2)  M. Crouzeix et A. L. Mignot
 *       Exercices d'analyse des equations differentielles,
 *       Masson,
 *       1986,
 *       p. 61
 * PROBLEM: direct optimized formulae have hard-coded
 *   coefficients in double precision.
 * TODO: compute numeric value at the fly why full 
 *   precision.
 * ------------------------------------------ 
 */

void
quadrature_on_geo::init_triangle (quadrature_option_type opt)
{
  // -------------------------------------------------------------------------
  // special case : superconvergent patch recovery nodes & weights
  // -------------------------------------------------------------------------
  quadrature_option_type::family_type f = opt.get_family();
  if (f == quadrature_option_type::superconvergent) {
    switch (opt.get_order()) {
      case 0:
      case 1:
        f =  quadrature_option_type::gauss;	
        break;
      case 2:
        f =  quadrature_option_type::middle_edge;	
        break;
      default:
	error_macro ("unsupported superconvergent("<<opt.get_order()<<")");
    }
  }
  // -------------------------------------------------------------------------
  // special case : Middle Edge formulae
  // e.g. when using special option in riesz_representer
  // -------------------------------------------------------------------------
  if (f == quadrature_option_type::middle_edge) {
    switch (opt.get_order()) {
     case 0 :
     case 1 :
     case 2 :
	    // middle edge formulae:
	    wx(x(Float(0.5), Float(0.0)), 1./6);
	    wx(x(Float(0.0), Float(0.5)), 1./6);
	    wx(x(Float(0.5), Float(0.5)), 1./6);
	    break;
     default:
	    error_macro ("unsupported Middle-Edge("<<opt.get_order()<<")");
    }
    return;
  }
  // -------------------------------------------------------------------------
  // special case : Gauss-Lobatto points
  // e.g. when using special option in riesz() function
  // -------------------------------------------------------------------------
  if (f == quadrature_option_type::gauss_lobatto) {
    switch (opt.get_order()) {
     case 0 :
     case 1 :
	    // trapezes:
	    wx(x(Float(0), Float(0)), 1./6);
	    wx(x(Float(1), Float(0)), 1./6);
	    wx(x(Float(0), Float(1)), 1./6);
	    break;
     case 2 :
	    // simpson:
	    wx(x(Float(0), Float(0)), 3/Float(120));
	    wx(x(Float(1), Float(0)), 3/Float(120));
	    wx(x(Float(0), Float(1)), 3/Float(120));
	    wx(x(Float(0.5), Float(0.0)), 8/Float(120));
	    wx(x(Float(0.0), Float(0.5)), 8/Float(120));
	    wx(x(Float(0.5), Float(0.5)), 8/Float(120));
	    wx(x(1/Float(3), 1/Float(3)), 27/Float(120));
	    break;
     default:
	    error_macro ("unsupported Gauss-Lobatto("<<opt.get_order()<<")");
    }
    return;
  }
  // -------------------------------------------------------------------------
  // default & general case : Gauss points
  // -------------------------------------------------------------------------
  assert_macro (f == quadrature_option_type::gauss,
        "unsupported quadrature family \"" << opt.get_family_name() << "\"");

  switch (opt.get_order()) {
   case 0 :
   case 1 :
	    // better than the general case: 
	    //  r=0 => 2 nodes
	    //  r=1 => 4 nodes
	    //  here : 1 node
	    wx(x(1/Float(3), 1/Float(3)), 1./2);
	    break;
   case 2 :
	    // better than the general case: 
	    //  r=2 => 6 nodes
	    //  here : 3 node
	    wx(x(1/Float(6), 1/Float(6)), 1/Float(6));
	    wx(x(4/Float(6), 1/Float(6)), 1/Float(6));
	    wx(x(1/Float(6), 4/Float(6)), 1/Float(6));
	    break;
#if !(defined(_RHEOLEF_HAVE_CLN) || defined(_RHEOLEF_HAVE_DOUBLEDOUBLE) || defined(_RHEOLEF_HAVE_LONG_DOUBLE))
	    // numerical values are inlined in double precision : 15 digits !
   case 3 :
   case 4 : {
	    // better than the general case: 
	    //  r=3 => 9 nodes
	    //  r=4 => 12 nodes
	    //  here : 6 node
            Float a  = 0.445948490915965;
            Float b  = 0.091576213509771; 
	    Float wa = 0.111690794839005;
	    Float wb = 0.054975871827661;
	    wx(x(a, a),     wa);
  	    wx(x(1-2*a, a), wa);
  	    wx(x(a, 1-2*a), wa);
  	    wx(x(b, b),     wb);
  	    wx(x(1-2*b, b), wb);
  	    wx(x(b, 1-2*b), wb);
	    break;
	  }
   case 5 :
   case 6 : {
	    // better than the general case: 
	    //  r=5 => 16 nodes
	    //  r=6 => 20 nodes
	    //  here : 12 node
	    Float a  = 0.063089014491502;
	    Float b  = 0.249286745170910;
	    Float c  = 0.310352451033785;
	    Float d  = 0.053145049844816;
	    Float wa = 0.025422453185103;
	    Float wb = 0.058393137863189;
	    Float wc = 0.041425537809187;
	    wx(x(a, a),         wa);
	    wx(x(1-2*a, a),     wa);
	    wx(x(a, 1-2*a),     wa);
	    wx(x(b, b),         wb);
	    wx(x(1-2*b, b),     wb);
	    wx(x(b, 1-2*b),     wb);
	    wx(x(c, d),         wc);
	    wx(x(d, c),         wc);
	    wx(x(1-(c+d), c),   wc);
	    wx(x(1-(c+d), d),   wc);
	    wx(x(c, 1-(c+d)),   wc);
	    wx(x(d, 1-(c+d)),   wc);
	    break;
	  }
#endif // BIGFLOAT
   default: {
    // Gauss-Legendre quadrature formulae 
    //  where Legendre = Jacobi(alpha=0,beta=0) polynoms
    size_t r = opt.get_order();
    size_type n0 = n_node_gauss(r);
    size_type n1 = n_node_gauss(r+1);
    vector<Float> zeta0(n0), omega0(n0);
    vector<Float> zeta1(n1), omega1(n1);
    gauss_jacobi (n0, 0, 0, zeta0.begin(), omega0.begin());
    gauss_jacobi (n1, 0, 0, zeta1.begin(), omega1.begin());

    for (size_type i = 0; i < n0; i++) {
      for (size_type j = 0; j < n1; j++) {
        // we transform the square into the triangle 
	Float eta_0 = (1+zeta0[i])*(1-zeta1[j])/4;
	Float eta_1 =              (1+zeta1[j])/2;
	Float J     =              (1-zeta1[j])/8;
        wx (x(eta_0,eta_1), J*omega0[i]*omega1[j]);
      }
    }
   }
  }
}