Quadrature.h 17.6 KB
Newer Older
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
// ============================================================================
// ==                                                                        ==
// == AMDiS - Adaptive multidimensional simulations                          ==
// ==                                                                        ==
// ============================================================================
// ==                                                                        ==
// ==  crystal growth group                                                  ==
// ==                                                                        ==
// ==  Stiftung caesar                                                       ==
// ==  Ludwig-Erhard-Allee 2                                                 ==
// ==  53175 Bonn                                                            ==
// ==  germany                                                               ==
// ==                                                                        ==
// ============================================================================
// ==                                                                        ==
// ==  http://www.caesar.de/cg/AMDiS                                         ==
// ==                                                                        ==
// ============================================================================

/** \file Quadrature.h */

#ifndef AMDIS_QUADRATURE_H
#define AMDIS_QUADRATURE_H

#include "BasisFunction.h"
#include "Flag.h"
#include "MemoryManager.h"
#include "FixVec.h"
#include <list>

namespace AMDiS {

  template<typename T> class WorldVector;
  template<typename T> class DimVec;
  template<typename T> class VectorOfFixVecs;
  template<typename T> class MatrixOfFixVecs;

  // ============================================================================
  // ===== class Quadrature =====================================================
  // ============================================================================

  /** 
   * \ingroup Assembler
   *
   * \brief
   * For the assemblage of the system matrix and right hand side vector of the 
   * linear system, we have to compute integrals, for example: 
   * \f[ \int_{\Omega} f(x)\varphi_i(x) dx \f]
   * For general data A, b, c, and f, these integrals can not be calculated 
   * exactly. Quadrature formulas have to be used in order to calculate the 
   * integrals approximately. Numerical integration in finite element methods is
   * done by looping over all grid elements and using a quadrature formula on 
   * each element.
   */
  class Quadrature
  {
  public:
    MEMORY_MANAGED(Quadrature);

  protected:
    /** \brief
     * Avoids call of default constructor
     */
    Quadrature();

    /** \brief
     * Constructs a Quadrature with name name_ of degree degree_ for dim dim_.
     * The Quadrature has n_points_ quadrature points with barycentric 
     * coordinates lambda_ and weights w_. The constructor is protected because
     * access to a Quadrature should be done via \ref provideQuadrature.
     */
    Quadrature(const char* name_,
73
74
75
	       int degree_,
	       int dim_,
	       int n_points_,
76
	       VectorOfFixVecs<DimVec<double> > *lambda_,
77
	       double* w_)
78
79
80
81
82
83
      : name(name_),
	degree(degree_),
	dim(dim_),
	n_points(n_points_),
	lambda(lambda_),
	w(w_) 
84
    {}
85
86

  public:
87
    /// Copy constructor
88
89
    Quadrature(const Quadrature&);

90
    /// Destructor
91
92
    ~Quadrature();

93
    /// Returns a Quadrature for dimension dim exact for degree degree.
Thomas Witkowski's avatar
Thomas Witkowski committed
94
    static Quadrature *provideQuadrature(int dim, int degree);
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110

    /** \brief
     * Approximates an integral by the numerical quadrature described by quad;
     * f is a pointer to an AbstractFunction to be integrated, evaluated in 
     * barycentric coordinates; the return value is 
     * \f[ \sum_{k = 0}^{n_points-1} w[k] * (*f)(lambda[k]) \f]
     * For the approximation of \f$ \int_S f\f$ we have to multiply this value 
     * with d!|S| for a simplex S; for a parametric simplex f should be a pointer
     * to a function which calculates 
     * \f$ f(\lambda)|det DF_S(\hat{x}(\lambda))| \f$.
     */
    double integrateStdSimplex(AbstractFunction<double, DimVec<double> > *f);
  
    /** \brief
     * Returns \ref name
     */
111
    inline const std::string& getName() { 
Thomas Witkowski's avatar
Thomas Witkowski committed
112
      return name; 
113
    }
114

115
    /// Returns \ref n_points
Thomas Witkowski's avatar
Thomas Witkowski committed
116
117
    inline int getNumPoints() const {
      return n_points;
118
    }
119
120
121
122

    /** \brief
     * Returns \ref w[p]
     */
Thomas Witkowski's avatar
Thomas Witkowski committed
123
124
    inline double getWeight(int p) const {
      return w[p];
125
    }
126
127
128
129

    /** \brief
     * Returns \ref w.
     */
Thomas Witkowski's avatar
Thomas Witkowski committed
130
131
    inline double* getWeight() const { 
      return w; 
132
    }
133
134
135
136

    /** \brief
     * Returns \ref dim
     */
Thomas Witkowski's avatar
Thomas Witkowski committed
137
138
    inline int getDim() const { 
      return dim; 
139
    }
140
141
142
143

    /** \brief
     * Returns \ref degree
     */
Thomas Witkowski's avatar
Thomas Witkowski committed
144
145
    inline int getDegree() const { 
      return degree; 
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

    /** \brief
     * Returns a pointer to a vector storing the values of a doubled valued 
     * function at all quadrature points; f is that AbstractFunction
     * , evaluated in barycentric coordinates; if vec is not NULL, the values are
     * stored in this vector, otherwise the values are stored in some static 
     * local vector, which is overwritten on the next call
     */
    const double *fAtQp(const AbstractFunction<double, DimVec<double> >& f,
			double *vec) const ;

    /** \brief
     * Returns a pointer to a vector storing the gradient (with respect to world
     * coordinates) of a double valued function at all quadrature points;
     * grdF is a pointer to a AbstractFunction, evaluated in barycentric
     * coordinates and returning a pointer to a WorldVector storing the gradient;
     * if vec is not NULL, the values are stored in this vector, otherwise the 
     * values are stored in some local static vector, which is overwritten on the
     * next call
     */
    const WorldVector<double> *grdFAtQp(const AbstractFunction<WorldVector<double>, 
					DimVec<double> >& grdF, 
					WorldVector<double>* vec) const;
  


    /** \brief
     * Returns \ref lambda[a][b] which is the b-th coordinate entry of the a-th
     * quadrature point
     */
Thomas Witkowski's avatar
Thomas Witkowski committed
177
    inline double getLambda(int a, int b) const {
178
      return (lambda ? (*lambda)[a][b] : 0.0);
179
    }
180
181
182
183
184
185
186

    /** \brief
     * Returns \ref lambda[a] which is a DimVec<double> containing the 
     * coordiantes of the a-th quadrature point
     */
    inline const DimVec<double>& getLambda(int a) const {
      return (*lambda)[a]; 
187
    }
188
189
190
191

    /** \brief
     * Returns \ref lambda which is a VectorOfFixvecs<DimVec<double> >.
     */
Thomas Witkowski's avatar
Thomas Witkowski committed
192
193
    VectorOfFixVecs<DimVec<double> > *getLambda() const { 
      return lambda; 
194
    }
195
196
197
198
199
200
201
202
203
204
205
206


  public:
    /** \brief
     * Maximal number of quadrature points for the different dimensions
     */
    static const int maxNQuadPoints[4];

  protected:
    /** \brief
     * Name of this Quadrature
     */
207
    std::string name;
208
209
210
211
212
213
214
215
216
217
218

    /** \brief
     * Quadrature is exact of this degree
     */
    int degree;

    /** \brief
     * Quadrature for dimension dim
     */
    int dim;

219
    /// Number of quadrature points
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
    int n_points;

    /** \brief
     * Vector of quadrature points given in barycentric coordinates
     */
    VectorOfFixVecs<DimVec<double> > *lambda; 

    /** \brief
     * Vector of quadrature weights
     */
    double *w;

  protected:
    /** \brief
     * Initialisation of all static Quadrature objects which will be returned
     * by \ref provideQuadrature()
     */
    static void initStaticQuadratures();

    /** \name static quadratures, used weights, and barycentric coords
     * \{
     */
    static Quadrature **quad_nd[4];
    static Quadrature *quad_0d[1];
    static Quadrature *quad_1d[20];
    static Quadrature *quad_2d[18];
    static Quadrature *quad_3d[8];

    static VectorOfFixVecs<DimVec<double> > *x_0d;
    static double *w_0d;

    static VectorOfFixVecs<DimVec<double> > *x0_1d;
    static VectorOfFixVecs<DimVec<double> > *x1_1d;
    static VectorOfFixVecs<DimVec<double> > *x2_1d;
    static VectorOfFixVecs<DimVec<double> > *x3_1d;
    static VectorOfFixVecs<DimVec<double> > *x4_1d;
    static VectorOfFixVecs<DimVec<double> > *x5_1d;
    static VectorOfFixVecs<DimVec<double> > *x6_1d;
    static VectorOfFixVecs<DimVec<double> > *x7_1d;
    static VectorOfFixVecs<DimVec<double> > *x8_1d;
    static VectorOfFixVecs<DimVec<double> > *x9_1d;
    static double *w0_1d;
    static double *w1_1d;
    static double *w2_1d;
    static double *w3_1d;
    static double *w4_1d;
    static double *w5_1d;
    static double *w6_1d;
    static double *w7_1d;
    static double *w8_1d;
    static double *w9_1d;

    static VectorOfFixVecs<DimVec<double> > *x1_2d;
    static VectorOfFixVecs<DimVec<double> > *x2_2d;
    static VectorOfFixVecs<DimVec<double> > *x3_2d;
    static VectorOfFixVecs<DimVec<double> > *x4_2d;
    static VectorOfFixVecs<DimVec<double> > *x5_2d;
    static VectorOfFixVecs<DimVec<double> > *x7_2d;
    static VectorOfFixVecs<DimVec<double> > *x8_2d;
    static VectorOfFixVecs<DimVec<double> > *x9_2d;
    static VectorOfFixVecs<DimVec<double> > *x10_2d;
    static VectorOfFixVecs<DimVec<double> > *x11_2d;
    static VectorOfFixVecs<DimVec<double> > *x12_2d;
    static VectorOfFixVecs<DimVec<double> > *x17_2d;
    static double *w1_2d;
    static double *w2_2d;
    static double *w3_2d;
    static double *w4_2d;
    static double *w5_2d;
    static double *w7_2d;
    static double *w8_2d;
    static double *w9_2d;
    static double *w10_2d;
    static double *w11_2d;
    static double *w12_2d;
    static double *w17_2d;

    static VectorOfFixVecs<DimVec<double> > *x1_3d;
    static VectorOfFixVecs<DimVec<double> > *x2_3d;
    static VectorOfFixVecs<DimVec<double> > *x3_3d;
    static VectorOfFixVecs<DimVec<double> > *x4_3d;
    static VectorOfFixVecs<DimVec<double> > *x5_3d;
    static VectorOfFixVecs<DimVec<double> > *x7_3d;
    static double *w1_3d;
    static double *w2_3d;
    static double *w3_3d;
    static double *w4_3d;
    static double *w5_3d;
    static double *w7_3d;

    /** \} */
  };

  // ============================================================================
  // ===== class FastQuadrature =================================================
  // ============================================================================

  /** \brief
   * Pre-compute the values of all basis functions at all quadrature nodes;  
   */
  const Flag INIT_PHI=1;

  /** \brief
   * Pre-compute the gradients (with respect to the barycentric coordinates) of
   * all basis functions at all quadrature nodes
   */
  const Flag INIT_GRD_PHI=2; 

  /** \brief
   * pre-compute all 2nd derivatives (with respect to the barycentric 
   * coordinates) of all basis functions at all quadrature nodes;
   */
  const Flag INIT_D2_PHI=4;

  // ============================================================================

  /** 
   * \ingroup Integration
   *
   *\brief
   * Often numerical integration involves basis functions, such as the assembling
   * of the system matrix and right hand side, or the integration of finite 
   * element functions. Since numerical quadrature involves only the values at 
   * the quadrature points and the values of basis functions and its derivatives
   * are the same at these points for all elements of the grid, such routines can
   * be much more efficient, if they can use pre-computed values of the basis
   * functions at the quadrature points. In this case the basis functions do not 
   * have to be evaluated for each quadrature point on every element newly.
   * Information that should be pre-computed can be specified by the following 
   * symbolic constants:
   * \ref INIT_PHI, \ref INIT_GRD_PHI, \ref INIT_D2_PHI
   */
  class FastQuadrature
  {
  public:
    MEMORY_MANAGED(FastQuadrature);

  protected:
    /** \brief
     * Constructs a FastQuadrature for the given Quadrature, BasisFunction, and
     * flag.
     */
    FastQuadrature(BasisFunction* basFcts, 
		   Quadrature* quad, 
		   Flag flag)
      : init_flag(flag), 
	phi(NULL), 
	grdPhi(NULL), 
	D2Phi(NULL), 
	quadrature(quad), 
	basisFunctions(basFcts) 
371
    {}
372
373
374
375
376
377
378
379
380

    /** \brief
     * Copy constructor
     */
    FastQuadrature(const FastQuadrature&);

    /** \brief
     * Extended copy constructor
     */
381
    FastQuadrature(const FastQuadrature&, const Flag);
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402

    /** \brief
     * Destructor
     */
    ~FastQuadrature();

  public:
    /** \brief
     * Returns a FastQuadrature for the given BasisFunction, Quadrature, and
     * flags
     */
    static FastQuadrature* provideFastQuadrature(const BasisFunction*,
						 const Quadrature&,
						 Flag);

    /** \brief
     * inits FastQuadrature like speciefied in flag
     */
    void init(Flag init_flag);

    inline bool initialized(Flag flag) {
403
      if (flag == INIT_PHI) {
404
405
406
	return (phi != NULL);
      }

407
      if (flag == INIT_GRD_PHI) {
408
409
410
	return (grdPhi != NULL);
      }

411
      if (flag == INIT_D2_PHI) {
412
413
414
415
416
	return (D2Phi != NULL);
      }

      ERROR_EXIT("invalid flag\n");
      return false;
417
    }
418
419
420
421

    /** \brief
     * Returns \ref quadrature
     */
422
423
    inline const Quadrature* getQuadrature() const { 
      return quadrature; 
424
    }
425
426
427
428

    /** \brief
     * Returns \ref max_points
     */
429
430
    inline int getMaxQuadPoints() { 
      return max_points; 
431
    }
432
433
434
435

    /** \brief
     * Returns (*\ref D2Phi)[q][i][j][m]
     */
436
    const double getSecDer(int q, int i, int j, int m) const;
437
438
439
440
441
442
443
444
445

    /** \brief
     * Returns (*\ref D2Phi)[q]
     */
    const VectorOfFixVecs<DimMat<double> > *getSecDer(int q) const;

    /** \brief
     * Returns (*\ref grdPhi)[q][i][j]
     */
446
447
    inline const double getGradient(int q, int i ,int j) const {
      return (grdPhi) ? (*grdPhi)[q][i][j] : 0.0;
448
    }
449
450
451
452

    /** \brief
     * Returns (*\ref grdPhi)[q]
     */
453
454
    inline VectorOfFixVecs<DimVec<double> >* getGradient(int q) const {
      return (grdPhi) ? &((*grdPhi)[q]) : NULL;
455
    }
456
457
458
459

    /** \brief
     * Returns \ref phi[q][i]
     */
460
461
    inline const double getPhi(int q, int i) const {
      return (phi) ? phi[q][i] : 0.0;
462
    }
463
464
465
466

    /** \brief
     * Returns \ref phi[q]
     */
467
468
    inline const double *getPhi(int q) const {
      return (phi) ? phi[q] : NULL;
469
    }
470
471
472
473

    /** \brief
     * Returns \ref quadrature ->integrateStdSimplex(f)
     */
474
    inline double integrateStdSimplex(AbstractFunction<double, DimVec<double> > *f) {
475
      return quadrature->integrateStdSimplex(f); 
476
    }
477
478
479
480
  
    /** \brief
     * Returns \ref quadrature ->getNumPoints()
     */
481
482
    inline int getNumPoints() const { 
      return quadrature->getNumPoints();
483
    }
484
485
486
487

    /** \brief
     * Returns \ref quadrature ->getWeight(p)
     */
488
489
    inline double getWeight(int p) const {
      return quadrature->getWeight(p);
490
    }
491
492
493
494

    /** \brief
     * Returns \ref quadrature ->getDim()
     */
495
496
    inline int getDim() const { 
      return quadrature->getDim(); 
497
    }
498
499
500
501

    /** \brief
     * Returns \ref quadrature ->getDegree()
     */
502
503
    inline int getDegree() const { 
      return quadrature->getDegree();
504
    }
505
506
507
508
509
510
511
512
513
514

    /** \brief
     * Returns \ref quadrature ->grdFAtQp(f, vec)
     */
    inline const WorldVector<double> 
    *grdFAtQp(const AbstractFunction<WorldVector<double>, 
	      DimVec<double> >& f,
	      WorldVector<double>* vec) const 
    { 
      return quadrature->grdFAtQp(f, vec);
515
    }
516
517
518
519
520
  
    /** \brief
     * Returns \ref quadrature ->fAtQp(f, vec)
     */
    inline const double *fAtQp(const AbstractFunction<double, 
521
			       DimVec<double> >& f, double *vec) const 
522
523
    {
      return quadrature->fAtQp(f, vec);    
524
    }
525
526
527
528
529
530

    /** \brief
     * Returns \ref quadrature ->getLambda(a,b)
     */
    inline double getLambda(int a,int b) const {
      return quadrature->getLambda(a,b);
531
    }
532
533
534
535
536
537

    /** \brief
     * Returns \ref quadrature ->getLambda(a)
     */
    inline const DimVec<double>& getLambda(int a) const { 
      return quadrature->getLambda(a); 
538
    }
539
540
541
542

    /** \brief
     * Returns \ref basisFunctions
     */
543
544
    inline BasisFunction* getBasisFunctions() const { 
      return basisFunctions; 
545
    }
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582

  protected:
    /** \brief
     * Specifies which information should be pre-computed. Can be \ref INIT_PHI, 
     * \ref INIT_GRD_PHI, or \ref INIT_D2_PHI
     */
    Flag init_flag;

    /** \brief
     * Matrix storing function values if the flag \ref INIT_PHI is set;
     * phi[i][j] stores the value \ref basisFunctions->phi[j]
     * (quadrature->lambda[i]), 0 <= j < basisFunctions->getNumber()  and 
     * 0 <= i < n_points
     */
    double **phi;

    /** \brief
     * Matrix storing all gradients (with respect to the barycentric coordinates)
     * if the flag \ref INIT_GRD_PHI is set; grdPhi[i][j][k] stores the value 
     * basisFunctions->grdPhi[j](quadrature->lambda[i])[k] 
     * for 0 <= j < basisFunctions->getNumber(),
     * 0 <= i < . . . , n_points, and 0 <= k < DIM
     */
    MatrixOfFixVecs<DimVec<double> > *grdPhi;

    /** \brief
     * Matrix storing all second derivatives (with respect to the barycentric
     * coordinates) if the flag \ref INIT_D2_PHI is set; D2Phi[i][j][k][l] stores
     * the value basisFunctions->D2Phi[j](quadrature->lambda[i])[k][l] 
     * for 0 <= j < basisFunctions->getNumber(), 
     * 0 <= i < n_points, and 0 <= k,l < DIM
     */
    MatrixOfFixVecs<DimMat<double> > *D2Phi;

    /** \brief
     * List of all used FastQuadratures
     */
583
    static std::list<FastQuadrature*> fastQuadList;
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606

    /** \brief
     * Maximal number of quadrature points for all yet initialised FastQuadrature
     * objects. This value may change after a new initialisation of a
     * FastQuadrature
     */
    static int max_points;

    /** \brief
     * This FastQuadrature stores values for Quadrature quadrature
     */
    Quadrature* quadrature;

    /** \brief
     * Values stored for basis functions basisFunctions
     */
    BasisFunction* basisFunctions;

  };

}

#endif  // AMDIS_QUADRATURE_H