File: cs_quadrature.c

package info (click to toggle)
code-saturne 4.3.3%2Brepack-1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 77,992 kB
  • sloc: ansic: 281,257; f90: 122,305; python: 56,490; makefile: 3,915; xml: 3,285; cpp: 3,183; sh: 1,139; lex: 176; yacc: 101; sed: 16
file content (400 lines) | stat: -rw-r--r-- 13,813 bytes parent folder | download
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
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
/*============================================================================
 * Routines to handle quadrature rules
 *============================================================================*/

/*
  This file is part of Code_Saturne, a general-purpose CFD tool.

  Copyright (C) 1998-2016 EDF S.A.

  This program 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.

  This program 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
  this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
  Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/

/*----------------------------------------------------------------------------*/

#include "cs_defs.h"

/*----------------------------------------------------------------------------
 * Standard C library headers
 *----------------------------------------------------------------------------*/

#include <math.h>
#include <assert.h>

/*----------------------------------------------------------------------------
 *  Local headers
 *----------------------------------------------------------------------------*/

#include "cs_cdo.h"

/*----------------------------------------------------------------------------
 * Header for the current file
 *----------------------------------------------------------------------------*/

#include "cs_quadrature.h"

/*----------------------------------------------------------------------------*/

BEGIN_C_DECLS

/*=============================================================================
 * Local Macro definitions
 *============================================================================*/

static const double  _quad_25ov48 = 25./48.;
static const double  _quad_9ov16 = 9./16.;
static const double  _quad_over18 = 1./18.;
static const double  _quad_over6 = 1./6.;
static const double  _quad_over3 = 1./3.;
static const double  _quad_9ov40 = 9./40.;
static const double  _quad_31ov240 = 31./240.;

/* Constant quadrature weights */
static double  _edge_quad2c1;
static double  _edge_quad2c2;
static double  _edge_quad3c1;
static double  _edge_quad3c2;
static double  _tria_quad7c1;
static double  _tria_quad7c2;
static double  _tria_quad7c3;
static double  _tria_quad7c4;
static double  _tetr_quad4c1;
static double  _tetr_quad4c2;

static const char
cs_quadrature_type_name[CS_QUADRATURE_N_TYPES][CS_BASE_STRING_LEN] =
  { N_("none"),
    N_("barycentric"),
    N_("higher (single weight)"),
    N_("highest") };

/*============================================================================
 * Private function prototypes
 *============================================================================*/

/*============================================================================
 * Public function prototypes
 *============================================================================*/

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Compute constant weights for all quadratures used
 */
/*----------------------------------------------------------------------------*/

void
cs_quadrature_setup(void)
{
  /* Quadrature on edge with 2 points */
  _edge_quad2c1 = 0.5*(1 + _quad_over3*sqrt(3.));
  _edge_quad2c2 = 0.5*(1 - _quad_over3*sqrt(3.));

  /* Quadrature on edge with 3 points */
  _edge_quad3c1 = 0.5*(1 + sqrt(0.6));
  _edge_quad3c2 = 0.5*(1 - sqrt(0.6));

  /* Quadrature on triangles with 7 points */
  _tria_quad7c1 = (6. - sqrt(15.))/21.;
  _tria_quad7c2 = (6. + sqrt(15.))/21.;
  _tria_quad7c3 = _quad_31ov240 - sqrt(15.)/1200.;
  _tria_quad7c4 = _quad_31ov240 + sqrt(15.)/1200.;

  /* Quadrature on tetrahedron with 4 points */
  _tetr_quad4c1 = 0.05*(5. - sqrt(5));
  _tetr_quad4c2 = 1. -3.*_tetr_quad4c1;
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Compute quadrature points for an edge from v1 -> v2 (2 points)
 *          Exact for polynomial function up to order 3
 *
 * \param[in]      v1       first vertex
 * \param[in]      v2       second vertex
 * \param[in]      len      length of edge [v1, v2]
 * \param[in, out] gpts     gauss points
 * \param[in, out] w        weight (same weight for the two points)
 */
/*----------------------------------------------------------------------------*/

void
cs_quadrature_edge_2pts(const cs_real_3_t  v1,
                        const cs_real_3_t  v2,
                        double             len,
                        cs_real_3_t        gpts[],
                        double            *w)
{
  int  k;

  /* Compute quadrature points */
  for (k = 0; k < 3; k++) {
    gpts[0][k] = v1[k]*_edge_quad2c2 + v2[k]*_edge_quad2c1;
    gpts[1][k] = v1[k]*_edge_quad2c1 + v2[k]*_edge_quad2c2;
  }

  /* Compute weights */
  *w= 0.5*len;
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Compute quadrature points for an edge from v1 -> v2 (3 points)
 *          Exact for polynomial function up to order 5
 *
 * \param[in]      v1       first vertex
 * \param[in]      v2       second vertex
 * \param[in]      len      length of edge [v1, v2]
 * \param[in ,out] gpts     gauss points
 * \param[in, out] w        weights
 */
/*----------------------------------------------------------------------------*/

void
cs_quadrature_edge_3pts(const cs_real_3_t  v1,
                        const cs_real_3_t  v2,
                        double             len,
                        cs_real_3_t        gpts[],
                        double             w[])
{
  int  k;

  const double  b = len * _quad_over18;

  /* Compute quadrature points */
  for (k = 0; k < 3; k++) {
    gpts[0][k] = 0.5*( v1[k] + v2[k] );
    gpts[1][k] = v1[k]*_edge_quad3c1 + v2[k]*_edge_quad3c2;
    gpts[2][k] = v1[k]*_edge_quad3c2 + v2[k]*_edge_quad3c1;
  }

  /* Compute weights */
  w[0]= 8*b, w[1] = w[2] = 5*b;
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Compute quadrature points for a triangle (3 points)
 *          Exact for polynomial function up to order 2
 *
 * \param[in]      v1       first vertex
 * \param[in]      v2       second vertex
 * \param[in]      v3       third vertex
 * \param[in]      area     area of triangle {v1, v2, v3}
 * \param[in, out] gpts     gauss points
 * \param[in, out] w        weight (same weight for the three points)
 */
/*----------------------------------------------------------------------------*/

void
cs_quadrature_tria_3pts(const cs_real_3_t   v1,
                        const cs_real_3_t   v2,
                        const cs_real_3_t   v3,
                        double              area,
                        cs_real_3_t         gpts[],
                        double             *w)
{
  int  k;

  /* Compute quadrature points */
  for (k = 0; k < 3; k++) {
    gpts[0][k] = 0.5*( v1[k] + v2[k] );
    gpts[1][k] = 0.5*( v1[k] + v3[k] );
    gpts[2][k] = 0.5*( v2[k] + v3[k] );
  }

  /* Compute weight */
  *w = _quad_over3 * area;
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Compute quadrature points for a triangle (4 points)
 *          Exact for polynomial function up to order 3
 *
 * \param[in]      v1       first vertex
 * \param[in]      v2       second vertex
 * \param[in]      v3       third vertex
 * \param[in]      area     area of triangle {v1, v2, v3}
 * \param[in, out] gpts     gauss points
 * \param[in, out] w        weights
 */
/*----------------------------------------------------------------------------*/

void
cs_quadrature_tria_4pts(const cs_real_3_t   v1,
                        const cs_real_3_t   v2,
                        const cs_real_3_t   v3,
                        double              area,
                        cs_real_3_t         gpts[],
                        double              w[])
{
  int  k;

  /* Compute quadrature points */
  for (k = 0; k < 3; k++) {
    gpts[0][k] = _quad_over3*( v1[k] + v2[k] + v3[k] );
    gpts[1][k] = 0.2*( v1[k] + v2[k] ) + 0.6*v3[k];
    gpts[2][k] = 0.2*( v1[k] + v3[k] ) + 0.6*v2[k];
    gpts[3][k] = 0.2*( v2[k] + v3[k] ) + 0.6*v1[k];
  }

  /* Compute weight */
  w[0] = -_quad_9ov16 * area;
  w[1] = w[2] = w[3] = _quad_25ov48 * area;
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief   Compute quadrature points for a triangle (7 points)
 *          Exact for polynomial function up to order 5
 *
 * \param[in]      v1       first vertex
 * \param[in]      v2       second vertex
 * \param[in]      v3       third vertex
 * \param[in]      area     area of triangle {v1, v2, v3}
 * \param[in, out] gpts     gauss points
 * \param[in, out] w        weights
 */
/*----------------------------------------------------------------------------*/

void
cs_quadrature_tria_7pts(const cs_real_3_t   v1,
                        const cs_real_3_t   v2,
                        const cs_real_3_t   v3,
                        double              area,
                        cs_real_3_t         gpts[],
                        double              w[])
{
  int  k;

  /* Compute quadrature points */
  for (k = 0; k < 3; k++) {
    gpts[0][k] = _quad_over3*( v1[k] + v2[k] + v3[k] );
    gpts[1][k] = _tria_quad7c1* (v1[k] + v2[k]) + (1-2*_tria_quad7c1)*v3[k];
    gpts[2][k] = _tria_quad7c1* (v3[k] + v1[k]) + (1-2*_tria_quad7c1)*v2[k];
    gpts[3][k] = _tria_quad7c1* (v2[k] + v3[k]) + (1-2*_tria_quad7c1)*v1[k];
    gpts[4][k] = _tria_quad7c2* (v1[k] + v2[k]) + (1-2*_tria_quad7c2)*v3[k];
    gpts[5][k] = _tria_quad7c2* (v3[k] + v1[k]) + (1-2*_tria_quad7c2)*v2[k];
    gpts[6][k] = _tria_quad7c2* (v2[k] + v3[k]) + (1-2*_tria_quad7c2)*v1[k];
  }

  /* Compute weights */
  w[0] = _quad_9ov40 * area;
  w[1] = w[2] = w[3] = _tria_quad7c3 * area;
  w[4] = w[5] = w[6] = _tria_quad7c4 * area;

}

/*----------------------------------------------------------------------------*/
/*!
 * \brief  Compute the quadrature in a tetrehedra. Exact for 2nd order
 *         polynomials (order 3).
 *
 * \param[in]       xv       first vertex
 * \param[in]       xe       second vertex
 * \param[in]       xf       third vertex
 * \param[in]       xc       fourth vertex
 * \param[in]       vol      volume of tetrahedron {xv, xe, xf, xc}
 * \param[in, out]  gpts     4 Gauss points (size = 3*4)
 * \param[in, out]  w        weight (same value for all points)
 */
/*----------------------------------------------------------------------------*/

void
cs_quadrature_tet_4pts(const cs_real_3_t  xv,
                       const cs_real_3_t  xe,
                       const cs_real_3_t  xf,
                       const cs_real_3_t  xc,
                       double             vol,
                       cs_real_3_t        gpts[],
                       double            *w)
{
  int  k;

  /* Compute Gauss points */
  for (k = 0; k < 3; k++) {
    gpts[0][k] = _tetr_quad4c1*(xv[k] + xe[k] + xf[k]) + _tetr_quad4c2*xc[k];
    gpts[1][k] = _tetr_quad4c1*(xe[k] + xf[k] + xc[k]) + _tetr_quad4c2*xv[k];
    gpts[2][k] = _tetr_quad4c1*(xf[k] + xc[k] + xv[k]) + _tetr_quad4c2*xe[k];
    gpts[3][k] = _tetr_quad4c1*(xc[k] + xv[k] + xe[k]) + _tetr_quad4c2*xf[k];
  }

  /* Compute weights (multiplicity 4) */
  *w = 0.25*vol;
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief  Compute the quadrature in a tetrehedra. Exact for 3rd order
 *         polynomials (order 4).
 *
 * \param[in]       xv       first vertex
 * \param[in]       xe       second vertex
 * \param[in]       xf       third vertex
 * \param[in]       xc       fourth vertex
 * \param[in]       vol      volume of tetrahedron {xv, xe, xf, xc}
 * \param[in, out]  gpts     5 Gauss points (size = 3*5)
 * \param[in, out]  weights  5 weigths related to each Gauss point
 */
/*----------------------------------------------------------------------------*/

void
cs_quadrature_tet_5pts(const cs_real_3_t  xv,
                       const cs_real_3_t  xe,
                       const cs_real_3_t  xf,
                       const cs_real_3_t  xc,
                       double             vol,
                       cs_real_3_t        gpts[],
                       double             weights[])
{
  int  k;

  const double  wv1 = -0.8*vol;  /* multiplicity 1 */
  const double  wv2 = 0.45*vol;  /* multiplicity 4 */

  /* compute Gauss points */
  for (k = 0; k < 3; k++) {
    gpts[0][k] = _quad_over6*(xv[k] + xe[k] + xf[k]) + 0.5*xc[k];
    gpts[1][k] = _quad_over6*(xe[k] + xf[k] + xc[k]) + 0.5*xv[k];
    gpts[2][k] = _quad_over6*(xf[k] + xc[k] + xv[k]) + 0.5*xe[k];
    gpts[3][k] = _quad_over6*(xc[k] + xv[k] + xe[k]) + 0.5*xf[k];
    gpts[4][k] = 0.25*(xv[k] + xe[k] + xf[k] + xc[k]);
  }

  /* Compute weights */
  for (k = 0; k < 4; k++)
    weights[k] = wv2;
  weights[4] = wv1;
}

/*----------------------------------------------------------------------------*/
/*!
 * \brief  Return th name associated to a type of quadrature
 *
 * \param[in]     type     cs_quadra_type_t
 *
 * \return the name associated to a given type of quadrature
 */
/*----------------------------------------------------------------------------*/

const char *
cs_quadrature_get_type_name(const cs_quadra_type_t  type)
{
  return cs_quadrature_type_name[type];
}

/*----------------------------------------------------------------------------*/

END_C_DECLS