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

/** \file BasisFunction.h */

#ifndef AMDIS_BASISFUNCTION_H
#define AMDIS_BASISFUNCTION_H

// ============================================================================
// ===== includes =============================================================
// ============================================================================

#include <string>
#include "Global.h"
#include "Boundary.h"
Thomas Witkowski's avatar
Thomas Witkowski committed
32
#include "MatrixVector.h"
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55

namespace AMDiS {

  // ============================================================================
  // ===== forward declarations =================================================
  // ============================================================================

  class DOFAdmin;
  class Element;
  class ElInfo;
  class RCNeighbourList;
  template<typename T> class WorldVector;
  template<typename T> class WorldMatrix;
  class Quadrature;

  template <typename ReturnType, typename ArgumentType> class AbstractFunction;
  template <typename T> class DOFVector;
  template <typename T> class DOFIndexed;
  template <typename T> class DimVec;
  template <typename T> class DimMat;
  template <typename T, GeoIndex d> class FixVec;
  template <typename T> class VectorOfFixVecs;

56
57
58
59
60
61
62
63
64
65
66

  // ============================================================================
  // ===== function interfaces===================================================
  // ============================================================================

  /** \brief
   * Function interface for evaluating basis functions.
   */
  class BasFctType
  {
  public:
67
    BasFctType() {}
68

69
    virtual ~BasFctType() {}
70
71
72
73
74
75
76
77
78
79
80

    virtual double operator()(const DimVec<double>&) const = 0;
  };


  /** \brief
   * Function interface for evaluating gradients of basis functions.
   */   
  class GrdBasFctType
  {
  public:
81
    GrdBasFctType() {}
82

83
    virtual ~GrdBasFctType() {}
84
85
86
87
88
89
90
91
92
93
94
95

    virtual void operator()(const DimVec<double>&,
			    DimVec<double>&) const = 0;
  };

  
  /** \brief
   * Function interface for evaluating second derivative of basis functions.
   */
  class D2BasFctType
  {
  public:
96
    D2BasFctType() {}
97

98
    virtual ~D2BasFctType() {}
99
100
101
102
103
104

    virtual void operator()(const DimVec<double>&,
			    DimMat<double>&) const = 0;
  };
			    

105
106
107
108
109
110
111
112
  // ============================================================================
  // ===== typedefs =============================================================
  // ============================================================================

  typedef BasFctType *BFptr;
  typedef GrdBasFctType *GBFptr;
  typedef D2BasFctType *DBFptr;

113

114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
  // ============================================================================
  // ===== class BasisFunction ==================================================
  // ============================================================================

  /** \ingroup FEMSpace
   * \brief
   * Base class for finite element basis functions. In order to build up a
   * finite element space, we have to specify a set of local basis functions.
   * Together with the correspondig DOF administration and the underlying mesh,
   * the finite element space is given. 
   * This class holds the local basis functions and their derivatives of the
   * reference element. They are evaluated at barycentric coordinates, so they
   * can be used on every element of the mesh.  
   */
  class BasisFunction
  {  
  protected:
    /** \brief
     * Creates a BasisFunction object of given dim and degree 
     */
134
    BasisFunction(const std::string& name, int dim, int degree);
135
136
137
138
139
140
141
142
143
144
145

    /** \brief
     * destructor
     */
    virtual ~BasisFunction();

  public:
    /** \brief
     * compares two BasisFunction objects.
     */
    virtual bool operator==(const BasisFunction& a) const {
146
      return a.getName() == name;
147
    }
148
149
150
151
152

    /** \brief
     * Returns !(*this == b)
     */
    inline bool operator!=(const BasisFunction& b) const {
153
      return !(operator == (b));
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

    /** \brief
     * Used by \ref getDOFIndices and \ref getVec
     */
    virtual int* orderOfPositionIndices(const Element* el, GeoIndex position, 
					int positionIndex) const = 0;

    /** \brief
     * Pointer to a function which connects the set of local basis functions
     * with its global DOFs.
     * getDOFIndices(el, admin, dof) returns a pointer to a const vector of 
     * length \ref nBasFcts where the i-th entry is the index of the DOF 
     * associated to the i-th basis function; arguments are the actual element 
     * el and the DOF admin admin of the corresponding finite element space 
     * (these indices depend on all defined DOF admins and thus on the 
     * corresponding admin); if the last argument dof is NULL, getDOFndices 
     * has to provide memory for storing this vector, which is overwritten on the
     * next call of getDOFIndices; if dof is not NULL, dof is a pointer to a 
     * vector which has to be filled;   
     */
    virtual const DegreeOfFreedom* getDOFIndices(const Element*,
						 const DOFAdmin&, 
						 DegreeOfFreedom*) const = 0;

    /** \brief
Thomas Witkowski's avatar
Thomas Witkowski committed
180
181
182
183
184
185
186
187
188
     * The second argument 'bound' has to be a pointer to a vector which has 
     * to be filled. Its length is \ref nBasFcts (the number of basis functions
     * in the used finite element space). After calling this function, the i-th 
     * entry of the array is the boundary type of the i-th basis function of this
     * element.
     * 
     * This function needs boundary information within the ElInfo object; thus, 
     * all routines using this function on the elements need the FILL_BOUND 
     * flag during mesh traversal;
189
     */
Thomas Witkowski's avatar
Thomas Witkowski committed
190
    virtual void getBound(const ElInfo*, BoundaryType *) const {};
191
192
193
194
195
196

    /** \brief
     * Returns \ref degree of BasisFunction
     */
    inline const int getDegree() const { 
      return degree; 
197
    }
198
199
200
201
202
203

    /** \brief
     * Returns \ref dim of BasisFunction
     */
    inline const int getDim() const { 
      return dim; 
204
    }
205
206
207
208
209
210

    /** \brief
     * Returns \ref nBasFcts which is the number of local basis functions
     */
    inline const int getNumber() const { 
      return nBasFcts; 
211
    }
212
213
214
215

    /** \brief
     * Returns \ref name of BasisFunction
     */
216
    inline const std::string& getName() const { 
217
      return name; 
218
    }
219
220
221
222
223
224
225
226
227
228
229

    /** \brief
     * Returns \ref nDOF[i]
     */
    int getNumberOfDOFs(int i) const;

    /** \brief
     * Returns \ref nDOF
     */
    inline DimVec<int>* getNumberOfDOFs() const { 
      return nDOF; 
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

    /** \brief
     * Initialisation of the \ref nDOF vector. Must be implemented by sub classes
     */
    virtual void setNDOF() = 0;

    /** \brief
     * Returns the barycentric coordinates of the i-th basis function.
     */
    virtual DimVec<double> *getCoords(int i) const = 0;

    /** \brief
     * Returns a pointer to a const vector with interpolation coefficients of the
     * function f; if indices is a pointer to NULL, the coefficient for all 
     * basis functions are calculated and the i-th entry in the vector is the 
     * coefficient of the i-th basis function; if indices is non NULL, only the 
     * coefficients for a subset of the local basis functions has to be 
     * calculated; n is the number of those basis functions, indices[0], . . . 
     * , indices[n-1] are the local indices of the basis functions where the
     * coefficients have to be calculated, and the i-th entry in the return 
     * vector is then the coefficient of the indices[i]-th basis function; coeff 
     * may be a pointer to a vector which has to be filled 
     * (compare the dof argument of \ref getDOFIndices());
     * such a function usually needs vertex coordinate information; thus, all 
     * routines using this function on the elements need the FILL COORDS flag 
     * during mesh traversal.
     * Must be implemented by sub classes.
     */
    virtual const double* interpol(const ElInfo *el_info, 
				   int n, const int *indices, 
				   AbstractFunction<double, WorldVector<double> > *f,
				   double *coeff) = 0;


    /** \brief
     * WorldVector<double> valued interpol function.
     */
    virtual const WorldVector<double>* 
    interpol(const ElInfo *el_info, int no, 
	     const int *b_no,
	     AbstractFunction<WorldVector<double>,WorldVector<double> > *f, 
	     WorldVector<double> *vec) = 0;

    /** \brief
     * Returns the i-th local basis function
     */
    inline BasFctType *getPhi(int i) const { 
      return (*phi)[i]; 
279
    }
280
281
282
283
284
285

    /** \brief
     * Returns the gradient of the i-th local basis function
     */
    inline GrdBasFctType *getGrdPhi(int i) const { 
      return  (*grdPhi)[i]; 
286
    }
287
288
289
290
291
292

    /** \brief
     * Returns the second derivative of the i-th local basis function
     */
    inline D2BasFctType *getD2Phi(int i) const { 
      return (*d2Phi)[i]; 
293
    }
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316

    /** \brief
     * Approximates the L2 scalar products of a given function with all basis 
     * functions by numerical quadrature and adds the corresponding values to a 
     * DOF vector;
     * f is a pointer for the evaluation of the given function in world 
     * coordinates x and returns the value of that function at x; if f is a NULL
     *  pointer, nothing is done;
     * fh is the DOF vector where at the i-th entry the approximation of the L2 
     * scalar product of the given function with the i-th global basis function 
     * of fh->feSpace is stored;
     * quad is the quadrature for the approximation of the integral on each leaf 
     * element of fh->feSpace->mesh; if quad is a NULL pointer, a default 
     * quadrature which is exact of degree 2*fh->feSpace->basFcts->degree-2 is 
     * used.
     * The integrals are approximated by looping over all leaf elements, 
     * computing the approximations to the element contributions and adding 
     * these values to the vector fh by add element vec().
     * The vector fh is not initialized with 0.0; only the new contributions are 
     * added
     */
    virtual void l2ScpFctBas(Quadrature*,
			     AbstractFunction<double, WorldVector<double> >* /*f*/,
317
318
			     DOFVector<double>* /*fh*/)
    {}
319
320
321
322
323
324

    /** \brief
     * WorldVector<double> valued l2ScpFctBas function
     */
    virtual void l2ScpFctBas(Quadrature* ,
			     AbstractFunction<WorldVector<double>, WorldVector<double> >* /*f*/,
325
326
			     DOFVector<WorldVector<double> >* /*fh*/) 
    {}
327
328
329
330
331


    /** \brief
     * Interpolates a DOFIndexed<double> after refinement
     */
332
333
    virtual void  refineInter(DOFIndexed<double> *, RCNeighbourList*, int)
    {}
334
335
336
337

    /** \brief
     * Interpolates a DOFIndexed<double> after coarsening
     */
338
339
    virtual void  coarseInter(DOFIndexed<double> *, RCNeighbourList*, int)
    {}
340
341
342
343

    /** \brief
     * Restricts a DOFIndexed<double> after coarsening
     */
344
345
    virtual void  coarseRestr(DOFIndexed<double> *, RCNeighbourList*, int)
    {}
346
347
348
349

    /** \brief
     * Interpolates a DOFVector<WorldVector<double> > after refinement
     */
350
351
    virtual void  refineInter(DOFVector<WorldVector<double> >*, RCNeighbourList*, int)
    {}
352
353
354
355

    /** \brief
     * Interpolates a DOFVector<WorldVector<double> > after coarsening
     */
356
357
    virtual void  coarseInter(DOFVector<WorldVector<double> >*, RCNeighbourList*, int)
    {}
358
359
360
361

    /** \brief
     * Restricts a DOFVector<WorldVector<double> > after coarsening
     */
362
363
    virtual void  coarseRestr(DOFVector<WorldVector<double> >*, RCNeighbourList*, int)
    {}
364
365
366
367
368
369
370
371
372

    /** \brief
     * Returns local dof indices of the element for the given fe space.
     */
    virtual const DegreeOfFreedom *getLocalIndices(const Element*,
						   const DOFAdmin*,
						   DegreeOfFreedom*) const
    {
      return NULL;
373
    }
374

Thomas Witkowski's avatar
Thomas Witkowski committed
375
376
377
378
379
380
    /** \brief
     * Returns local dof indices of the element for the given fe space.
     */
    virtual void getLocalIndicesVec(const Element*,
				    const DOFAdmin*,
				    Vector<DegreeOfFreedom>*) const
381
    {}
Thomas Witkowski's avatar
Thomas Witkowski committed
382

383
384
385
386
387
388
389

    /** \brief
     * Evaluates elements value at barycentric coordinates lambda with local 
     * coefficient vector uh.
     */
    double evalUh(const DimVec<double>& lambda, const double* uh) const;
  
390
391
392
393
394
395
396
397
398
    /** \brief
     * Evaluates elements value at barycentric coordinates lambda with local 
     * coefficient vector uh. If val is not NULL the result will be stored 
     * there, otherwise a pointer to a static local variable is returned which 
     * will be overwritten after the next call.
     */
    const WorldVector<double>& evalUh(const DimVec<double>& lambda, 
				      const WorldVector<double>* uh, WorldVector<double>* val) const;
    
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
    /** \brief
     * Evaluates the gradient at barycentric coordinates lambda. Lambda is the
     * Jacobian of the barycentric coordinates. uh is the local coefficient
     * vector. If val is not NULL the result will be stored 
     * there, otherwise a pointer to a static local variable is returned which 
     * will be overwritten after the next call.
     */
    const WorldVector<double>& evalGrdUh(const DimVec<double>& lambda,
					 const DimVec<WorldVector<double> >& Lambda,
					 const double* uh,  WorldVector<double>* val) const;

    /** \brief
     * Evaluates the second derivative at barycentric coordinates lambda. 
     * Lambda is the Jacobian of the barycentric coordinates. uh is the local 
     * coefficient vector. If val is not NULL the result will be stored 
     * there, otherwise a pointer to a static local variable is returned which 
     * will be overwritten after the next call.
     */
    const WorldMatrix<double>& evalD2Uh(const DimVec<double>& lambda,
					const DimVec<WorldVector<double> >& Lambda,
					const double* uh, WorldMatrix<double>* val) const;

  protected:
    /** \brief
     * Textual description
     */
425
    std::string name;     
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441

    /** \brief
     * Number of basisfunctions on one Element
     */                    
    int nBasFcts;

    /** \brief
     * Maximal degree of the basis functions
     */                    
    int degree;

    /** \brief
     * Dimension of the basis functions
     */                    
    int dim;

Thomas Witkowski's avatar
Thomas Witkowski committed
442
443
444
445
446
    /** \brief
     * Dimension of the world.
     */
    int dow;

447
448
449
450
451
452
453
454
    /** \brief
     * Number of DOFs at the different positions
     */                    
    DimVec<int> *nDOF;

    /** \brief
     * Vector of the local functions
     */
455
    std::vector<BasFctType*> *phi;
456
457
458
459

    /** \brief
     * Vector of gradients
     */
460
    std::vector<GrdBasFctType*> *grdPhi;
461
462
463
464

    /** \brief
     * Vector of second derivatives
     */
465
    std::vector<D2BasFctType*> *d2Phi;
Thomas Witkowski's avatar
Thomas Witkowski committed
466
467
468
469
470
471


    /** \brief
     * Is used by function evalGrdUh. To make it thread safe, we need a
     * temporary DimVec for every thread.
     */
472
    std::vector<DimVec<double>* > grdTmpVec1;
Thomas Witkowski's avatar
Thomas Witkowski committed
473
474
475
476
477

    /** \brief
     * Is used by function evalGrdUh. To make it thread safe, we need a
     * temporary DimVec for every thread.
     */
478
    std::vector<DimVec<double>* > grdTmpVec2;
479
480
481
482
483
  };

}

#endif  // AMDIS_BASISFUNCTION_H