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

/** \file DOFVector.h */

#ifndef AMDIS_DOFVECTOR_H 
#define AMDIS_DOFVECTOR_H 
 
// ===========================================================================
// ===== includes ============================================================
// ===========================================================================

29
30
31
32
#include <vector> 
#include <memory> 
#include <list> 

33
34
35
36
37
38
39
40
41
42
43
44
#include "FixVec.h"
#include "Global.h" 
#include "Flag.h" 
#include "RCNeighbourList.h" 
#include "DOFIterator.h"
#include "DOFIndexed.h"
#include "DOFContainer.h"
#include "MemoryManager.h"
#include "Boundary.h"
#include "CreatorInterface.h"
#include "Serializable.h"
#include "DOFMatrix.h" 
Thomas Witkowski's avatar
Thomas Witkowski committed
45
#include "BasisFunction.h"
46
47
48
49

#ifdef HAVE_PARALLEL_AMDIS
#include "petscao.h"
#endif
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
 
namespace AMDiS {

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

  class Mesh; 
  class FiniteElemSpace; 
  class ElInfo; 
  class DOFAdmin; 
  class BasisFunction; 
  class FillInfo; 
  class Quadrature; 
  class FastQuadrature;
  class DOFMatrix; 
  class MultiGridSortSolver; 
  class Operator;
  class ElementVector;
  class BoundaryManager;

71
  template<typename ResultType, typename ArgumentType> class AbstractFunction;
72
73
74
75
76
77


  template<typename T> 
  class DOFVectorBase : public DOFIndexed<T>
  {
  public:
Thomas Witkowski's avatar
Thomas Witkowski committed
78

79
80
81
    DOFVectorBase() 
      : feSpace(NULL),
	elementVector(NULL),
Thomas Witkowski's avatar
Thomas Witkowski committed
82
83
        boundaryManager(NULL),
        nBasFcts(0)
84
    {}
Thomas Witkowski's avatar
Thomas Witkowski committed
85
    
86
    DOFVectorBase(const FiniteElemSpace *f, std::string n);
87

Thomas Witkowski's avatar
Thomas Witkowski committed
88
    virtual ~DOFVectorBase();
89

Thomas Witkowski's avatar
Thomas Witkowski committed
90
91
92
93
94
    /** \brief 
     * For the given element, this function returns an array of all DOFs of this
     * DOFVector that are defined on this element.
     */
    virtual const T *getLocalVector(const Element *el, T* localVec) const;
95

Thomas Witkowski's avatar
Thomas Witkowski committed
96
97
98
99
    /** \brief
     * Evaluates the DOF vector at a set of quadrature points defined on the 
     * given element.
     */
100
101
    const T *getVecAtQPs(const ElInfo *elInfo, 
			 const Quadrature *quad,
102
			 const FastQuadrature *quadFast,
103
			 T *vecAtQPs) const;
104

Thomas Witkowski's avatar
Thomas Witkowski committed
105
106
107
108
109
110
    const T *getVecAtQPs(const ElInfo *smallElInfo, 
			 const ElInfo *largeElInfo,
			 const Quadrature *quad,
			 const FastQuadrature *quadFast,
			 T *vecAtQPs) const;

111
112
    const WorldVector<T> *getGrdAtQPs(const ElInfo *elInfo, 
				      const Quadrature *quad,
113
				      const FastQuadrature *quadFast,
114
				      WorldVector<T> *grdAtQPs) const;
115

Thomas Witkowski's avatar
Thomas Witkowski committed
116
117
118
119
120
121
    const WorldVector<T> *getGrdAtQPs(const ElInfo *smallElInfo, 
				      const ElInfo *largeElInfo,
				      const Quadrature *quad,
				      const FastQuadrature *quadFast,
				      WorldVector<T> *grdAtQPs) const;

122
123
    const WorldMatrix<T> *getD2AtQPs(const ElInfo *elInfo, 
				     const Quadrature *quad,
124
				     const FastQuadrature *quadFast,
125
				     WorldMatrix<T> *d2AtQPs) const;
126

Thomas Witkowski's avatar
Thomas Witkowski committed
127
    /// Returns the FE space the DOF vector is defined on.
128
    inline const FiniteElemSpace* getFESpace() const {
129
      return feSpace;
130
    }
131

Thomas Witkowski's avatar
Thomas Witkowski committed
132
133
134
135
136
137
138
139
140
141
142
143
144
    /** \brief
     * Assembles the element vector for the given ellement and adds the
     * element matrix to the current DOF vector.
     */ 
    void assemble(T factor, ElInfo *elInfo,			    
		  const BoundaryType *bound, 
		  Operator *op = NULL);

    void assemble2(T factor, 
		   ElInfo *mainElInfo, ElInfo *auxElInfo,
		   ElInfo *smallElInfo, ElInfo *largeElInfo,
		   const BoundaryType *bound, 
		   Operator *op = NULL);
145
146
147
148
149

    void addElementVector(T sign,
			  const ElementVector &elVec, 
			  const BoundaryType *bound,
			  bool add = true); 
Thomas Witkowski's avatar
Thomas Witkowski committed
150
151
152
153
154
155
156

    /* \brief
     * That function must be called after the matrix assembling has been finished. 
     * This makes it possible to start some cleanup or matrix data compressing 
     * procedures.
     */
    void finishAssembling();
157
158
159
160
161
162
163
164
 
    inline void addOperator(Operator* op, 
			    double *factor = NULL,
			    double *estFactor = NULL) 
    {
      operators.push_back(op);
      operatorFactor.push_back(factor);
      operatorEstFactor.push_back(estFactor);
165
    }
166

167
    inline std::vector<double*>::iterator getOperatorFactorBegin() {
168
      return operatorFactor.begin();
169
    }
170

171
    inline std::vector<double*>::iterator getOperatorFactorEnd() {
172
      return operatorFactor.end();
173
    }
174

175
    inline std::vector<double*>::iterator getOperatorEstFactorBegin() {
176
      return operatorEstFactor.begin();
177
    }
178

179
    inline std::vector<double*>::iterator getOperatorEstFactorEnd() {
180
      return operatorEstFactor.end();
181
    }
182
183


184
    inline std::vector<Operator*>::iterator getOperatorsBegin() {
185
      return operators.begin();
186
    }
187

188
    inline std::vector<Operator*>::iterator getOperatorsEnd() {
189
      return operators.end();
190
    }
191
192
193
194
195
196
197
198
199
200
201

    Flag getAssembleFlag();

    /** \brief
     * Evaluates \f[ u_h(x(\lambda)) = \sum_{i=0}^{m-1} vec[ind[i]] * 
     * \varphi^i(\lambda) \f] where \f$ \varphi^i \f$ is the i-th basis function,
     * \f$ x(\lambda) \f$ are the world coordinates of lambda and
     * \f$ m \f$ is the number of basis functions
     */
    T evalUh(const DimVec<double>& lambda, DegreeOfFreedom* ind);

202
    inline std::vector<Operator*>& getOperators() { 
Thomas Witkowski's avatar
Thomas Witkowski committed
203
      return operators; 
204
    }
Thomas Witkowski's avatar
Thomas Witkowski committed
205

206
    inline std::vector<double*>& getOperatorFactor() { 
Thomas Witkowski's avatar
Thomas Witkowski committed
207
      return operatorFactor; 
208
    }
209

210
    inline std::vector<double*>& getOperatorEstFactor() { 
211
      return operatorEstFactor; 
212
    }
213

Thomas Witkowski's avatar
Thomas Witkowski committed
214
    /// Returns \ref name
215
    inline const std::string& getName() const { 
Thomas Witkowski's avatar
Thomas Witkowski committed
216
      return name; 
217
    } 
218

Thomas Witkowski's avatar
Thomas Witkowski committed
219
220
    inline BoundaryManager* getBoundaryManager() const { 
      return boundaryManager; 
221
    }
222
223
224

    inline void setBoundaryManager(BoundaryManager *bm) {
      boundaryManager = bm;
225
    }
226
227

  protected:
Thomas Witkowski's avatar
Thomas Witkowski committed
228
    ///
229
230
    const FiniteElemSpace *feSpace;

Thomas Witkowski's avatar
Thomas Witkowski committed
231
    ///
232
    std::string name;
233

Thomas Witkowski's avatar
Thomas Witkowski committed
234
    ///
235
236
    ElementVector *elementVector;

Thomas Witkowski's avatar
Thomas Witkowski committed
237
    ///
238
    std::vector<Operator*> operators;
239

Thomas Witkowski's avatar
Thomas Witkowski committed
240
    ///
241
    std::vector<double*> operatorFactor;
242

Thomas Witkowski's avatar
Thomas Witkowski committed
243
    ///
244
    std::vector<double*> operatorEstFactor;
245

Thomas Witkowski's avatar
Thomas Witkowski committed
246
    ///
247
248
    BoundaryManager *boundaryManager;

Thomas Witkowski's avatar
Thomas Witkowski committed
249
    /// Number of basis functions of the used finite element space.
Thomas Witkowski's avatar
Thomas Witkowski committed
250
    int nBasFcts;
251

Thomas Witkowski's avatar
Thomas Witkowski committed
252
    /// Are used to store temporary local dofs of an element. Thread safe.
253
254
    std::vector<DegreeOfFreedom*> localIndices;

Thomas Witkowski's avatar
Thomas Witkowski committed
255
    /// Are used to store temporary local values of an element. Thread safe.
Thomas Witkowski's avatar
Thomas Witkowski committed
256
257
    std::vector<T*> localVectors;

Thomas Witkowski's avatar
Thomas Witkowski committed
258
    /// Temporarly used in \ref getGrdAtQPs. Thread safe.
Thomas Witkowski's avatar
Thomas Witkowski committed
259
260
    std::vector<DimVec<double>*> grdTmp;

Thomas Witkowski's avatar
Thomas Witkowski committed
261
    /// Temporarly used in \ref getGrdAtQPs. Thread safe.
262
263
    std::vector<DimVec<double>*> grdPhis;

Thomas Witkowski's avatar
Thomas Witkowski committed
264
    /// Temporarly used in \ref getD2AtQPs. Thread safe.
265
    std::vector<DimMat<double>*> D2Phis;
Thomas Witkowski's avatar
Thomas Witkowski committed
266

Thomas Witkowski's avatar
Thomas Witkowski committed
267
    /// Dimension of the mesh this DOFVectorBase belongs to
Thomas Witkowski's avatar
Thomas Witkowski committed
268
    int dim;
Thomas Witkowski's avatar
Thomas Witkowski committed
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


  // ===========================================================================
  // ===== defs ================================================================
  // ===========================================================================

  /** \brief
   * Specifies which operation should be done after coarsening
   */
  typedef enum{
    NO_OPERATION = 0,   
    COARSE_RESTRICT = 1,
    COARSE_INTERPOL = 2 
  } CoarsenOperation;
 
  // ===========================================================================
  // ===== class DOFVector =====================================================
  // ===========================================================================

  /** \ingroup DOFAdministration 
   * \brief
   * The DOFs described above are just integers that can be used as indices into 
   * vectors and matrices. During refinement and coarsening of the mesh, the 
   * number of used DOFs, the meaning of one integer index, and even the total 
   * range of DOFs change. To be able to handle these changes automatically for 
   * all vectors, which are indexed by the DOFs, special data structures are 
   * used which contain such vector data. Lists of these structures are kept in 
   * DOFAdmin, so that all vectors in the lists can be resized together with the
   * range of DOFs. During refinement and coarsening of elements, values can be
   * interpolated automatically to new DOFs, and restricted from old DOFs.
   */
  template<typename T> 
  class DOFVector : public DOFVectorBase<T>, public Serializable
  {  
  public:
    MEMORY_MANAGED(DOFVector<T>);

    /** \ingroup DOFAdministration
     * \brief
     * Enables the access of DOFVector<T>::Iterator. Alias for DOFIterator<T>
     */
    class Iterator : public DOFIterator<T> {
    public:
      Iterator(DOFIndexed<T> *c, DOFIteratorType type)
	: DOFIterator<T>(c, type)
315
      {}
316
317
318

      Iterator(DOFAdmin *admin, DOFIndexed<T> *c, DOFIteratorType type)
	: DOFIterator<T>(admin, c, type)
319
      {}
320
321
322
323
324
325
    };

    class Creator : public CreatorInterface<DOFVector<T> > {
    public:
      MEMORY_MANAGED(Creator);

326
327
328
      Creator(FiniteElemSpace *feSpace_) 
        : feSpace(feSpace_) 
      {}
329
330
331

      DOFVector<T> *create() { 
	return NEW DOFVector<T>(feSpace, ""); 
332
      }
333
334
335

      void free(DOFVector<T> *vec) { 
	DELETE vec; 
336
      }
337
338
339
340
341
342
343
344
345
346
347
348
349

    private:
      FiniteElemSpace *feSpace;
    };

  public:
    /** \brief
     * Empty constructor. No initialization!
     */
    DOFVector() 
      : DOFVectorBase<T>(), 
	feSpace(NULL),
	coarsenOperation(NO_OPERATION)
350
    {}
351
352
353
354

    /** \brief
     * Constructs a DOFVector with name n belonging to FiniteElemSpace f
     */
355
    DOFVector(const FiniteElemSpace* f, std::string n); 
356
357
358
359

    /** \brief
     * Initialization.
     */
360
    void init(const FiniteElemSpace* f, std::string n);
361
362
363
364
365

    /** \brief
     * Copy Constructor
     */
    DOFVector(const DOFVector& rhs) {
Thomas Witkowski's avatar
Thomas Witkowski committed
366
367
368
      *this = rhs;   
      name = rhs.name + "copy";
      if (feSpace && feSpace->getAdmin()) {
369
370
	(dynamic_cast<DOFAdmin*>(feSpace->getAdmin()))->addDOFIndexed(this);
      }
371
    }
372
373
374
375
376
377
378
379
380

    /** \brief
     * Destructor
     */
    virtual ~DOFVector();

    /** \brief
     * Returns iterator to the begin of \ref vec
     */
381
    typename std::vector<T>::iterator begin() { 
382
      return vec.begin(); 
383
    }
384
385
386
387

    /** \brief
     * Returns iterator to the end of \ref vec
     */
388
    typename std::vector<T>::iterator end() { 
389
      return vec.end(); 
390
    }
391
392
393
394
395
396
  
    /** \brief
     * Used by DOFAdmin to compress this DOFVector. Implementation of
     * DOFIndexedBase::compress()
     */
    virtual void compressDOFIndexed(int first, int last,
397
				    std::vector<DegreeOfFreedom> &newDof);
398
399
400
401
402
403

    /** \brief
     * Sets \ref coarsenOperation to op
     */
    inline void setCoarsenOperation(CoarsenOperation op) { 
      coarsenOperation = op; 
404
    }
405
406
407
408

    /** \brief
     * Returns \ref coarsenOperation
     */
409
410
    inline CoarsenOperation getCoarsenOperation() { 
      return coarsenOperation; 
411
    }
412
413
414
415

    /** \brief
     * Restriction after coarsening. Implemented for DOFVector<double>
     */
416
    inline void coarseRestrict(RCNeighbourList&, int) {}
417
418
419
420

    /** \brief
     * Interpolation after refinement.
     */
421
    inline void refineInterpol(RCNeighbourList&, int) {}
422
423
424
425

    /** \brief
     * Returns \ref vec
     */
426
    std::vector<T>& getVector() { 
Thomas Witkowski's avatar
Thomas Witkowski committed
427
      return vec;
428
    }
429
430
431
432

    /** \brief
     * Returns size of \ref vec
     */
433
434
    inline int getSize() const { 
      return vec.size();
435
    } 
436
437
438
439

    /** \brief
     * Returns used size of the vector.
     */
440
441
    inline int getUsedSize() const { 
      return feSpace->getAdmin()->getUsedSize(); 
442
    }
443
444
445
446
447

    /** \brief
     * Resizes \ref vec to n
     */
    inline void resize(int n) { 
448
449
      FUNCNAME("DOFVector<T>::resize()");
      TEST_EXIT_DBG((n >= 0))("Can't resize DOFVector to negative size\n"); 
450
      vec.resize(n);
451
    } 
452
453
454
455
456

    /** \brief
     * Resizes \ref vec to n and inits new values with init
     */
    inline void resize(int n, T init) { 
457
458
      FUNCNAME("DOFVector<T>::resize()");
      TEST_EXIT_DBG((n >= 0))("Can't resize DOFVector to negative size\n"); 
459
      vec.resize(n, init);
460
    } 
461
462
463
464
465
466

    /** \brief
     * Returns \ref vec[i]
     */
    inline const T& operator[](DegreeOfFreedom i) const {
      FUNCNAME("DOFVector<T>::operator[]");
467
468
      TEST_EXIT_DBG((i>= 0) && (i < static_cast<int>(vec.size())))
	("Illegal vector index %d.\n",i);
469
      return vec[i];
470
    } 
471

472
473
    inline int getSize() { 
      return vec.size(); 
474
    }
475
476
477
478
479
480

    /** \brief
     * Returns \ref vec[i]
     */
    inline T& operator[](DegreeOfFreedom i) {
      FUNCNAME("DOFVector<T>::operator[]");
481
482
      TEST_EXIT_DBG((i >= 0) && (i < static_cast<int>(vec.size())))
	("Illegal vector index %d.\n",i);
483
      return vec[i];
484
    }
Thomas Witkowski's avatar
Thomas Witkowski committed
485
 
486
487
488
489
490
491
492
493
494
    /** \brief
     * Calculates Integral of this DOFVector
     */
    double Int(Quadrature* q = NULL) const;

    /** \brief
     * Calculates L1 norm of this DOFVector
     */
    double L1Norm(Quadrature* q = NULL) const;
Thomas Witkowski's avatar
Thomas Witkowski committed
495
 
496
497
498
499
500
    /** \brief
     * Calculates L2 norm of this DOFVector
     */
    inline double L2Norm(Quadrature* q = NULL) const {
      return sqrt(L2NormSquare());
501
    }
502
503
504
505
506
507
508
509
510
511
512

    /** \brief
     * Calculates square of L2 norm of this DOFVector
     */
    double L2NormSquare(Quadrature* q = NULL) const;

    /** \brief
     * Calculates H1 norm of this DOFVector
     */
    inline double H1Norm(Quadrature* q = NULL) const {
      return sqrt(H1NormSquare());
513
    }
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534

    /** \brief
     * Calculates square of H1 norm of this DOFVector
     */
    double H1NormSquare(Quadrature* q = NULL) const;  

    /** \brief
     * Calculates euclidian norm of this DOFVector
     */
    double nrm2() const; 

    /** \brief
     * Returns square of the euclidian norm.
     */
    double squareNrm2() const;

    /** \brief
     * Calculates l2 norm of this DOFVector
     */
    inline double l2norm() const { 
      return nrm2();
535
    }
536
537
538
539
540
541
542
543
544
545
546

    /** \brief
     * Calculates the absolute sum of this DOFVector
     */
    T asum() const; 

    /** \brief
     * Calculates the l1 norm of this DOFVector
     */
    inline double l1norm() const { 
      return asum();
547
    } 
548
549
550
551
552

    /** \brief
     * Calculates doublewell of this DOFVector
     */
    double DoubleWell(Quadrature* q = NULL) const;
Thomas Witkowski's avatar
Thomas Witkowski committed
553
 
554
555
556
557
    /** \brief
     * Calculates the sum of this DOFVector
     */
    T sum() const; 
Thomas Witkowski's avatar
Thomas Witkowski committed
558
 
559
560
561
562
563
564
565
566
567
568
569
    /** \brief
     * Sets \ref vec[i] = val, i=0 , ... , size
     */
    void set(T val); 

    /** \brief
     * Assignment operator for setting complete vector to a certain value d
     */
    inline DOFVector<T>& operator=(T d) {
      set(d); 
      return *this;
570
    } 
571
572
573
574

    /** \brief
     * Assignment operator between two vectors
     */
575
    DOFVector<T>& operator=(const DOFVector<T>&); 
576
577
578
579
580
581
582
583
584
585
586
587
588
589

    /** \brief
     * vec[i] = v.vec[i]
     */
    void copy(const DOFVector<T>& v); 

    /** \brief
     * Returns minimum of DOFVector
     */
    T min() const; 

    /** \brief
     * Returns maximum of DOFVector
     */
Thomas Witkowski's avatar
Thomas Witkowski committed
590
591
    T max() const;

592
    /// Returns absolute maximum of DOFVector
Thomas Witkowski's avatar
Thomas Witkowski committed
593
    T absMax() const;
594

595
    /// Used by interpol while mesh traversal
596
    static int interpolFct(ElInfo* elinfo);
597
 
598
    /// Prints \ref vec to stdout
599
600
    void print() const; 

601
    ///
602
603
    int calcMemoryUsage() const;

604
605
606
607
608
609
610
611
    /** \brief
     * Computes the coefficients of the interpolant of the function fct and
     * stores these in the DOFVector
     */
    void interpol(AbstractFunction<T, WorldVector<double> > *fct);

    void interpol(DOFVector<T> *v, double factor);

612
    void serialize(std::ostream &out) {
613
614
615
      unsigned int size = vec.size();
      out.write(reinterpret_cast<const char*>(&size), sizeof(unsigned int));
      out.write(reinterpret_cast<const char*>(&(vec[0])), size * sizeof(T));
616
    }
617

618
    void deserialize(std::istream &in) {
619
620
621
622
      unsigned int size;
      in.read(reinterpret_cast<char*>(&size), sizeof(unsigned int));
      vec.resize(size);
      in.read(reinterpret_cast<char*>(&(vec[0])), size * sizeof(T));
623
    }
624
625
626
627
628
629
630

    DOFVector<WorldVector<T> > *getGradient(DOFVector<WorldVector<T> >*) const;

    WorldVector<DOFVector<T>*> *getGradient(WorldVector<DOFVector<T>*> *grad) const;

    DOFVector<WorldVector<T> >*

631
    getRecoveryGradient(DOFVector<WorldVector<T> >*) const;
632

633
634
635
636
637
638
639
#ifdef HAVE_PARALLEL_AMDIS
    /// Sets the petsc application ordering object to map dof indices.
    void useApplicationOrdering(AO *ao) {
      applicationOrdering = ao;
    }
#endif

Thomas Witkowski's avatar
Thomas Witkowski committed
640
  protected: 
641
642
643
644
645
646
647
648
649
    /** \brief
     * Used by Int while mesh traversal
     */
    static int Int_fct(ElInfo* elinfo);

    /** \brief
     * Used by L1Norm while mesh traversal
     */
    static int L1Norm_fct(ElInfo* elinfo);
Thomas Witkowski's avatar
Thomas Witkowski committed
650
 
651
652
653
654
655
656
657
658
659
660
    /** \brief
     * Used by L2Norm while mesh traversal
     */
    static int L2Norm_fct(ElInfo* elinfo);

    /** \brief
     * Used by H1Norm while mesh traversal
     */
    static int H1Norm_fct(ElInfo* elinfo);

661
662
663
664
665
    /** \brief
     * Used by DoubleWell while mesh traversal
     */
    static int DoubleWell_fct(ElInfo* elinfo);

666
  protected: 
667
    /// Name of this DOFVector
668
    std::string name; 
669

670
    /// FiniteElemSpace of the vector
671
672
    const FiniteElemSpace *feSpace; 

673
    /// Data container
674
    std::vector<T> vec; 
675
 
676
    /// Specifies what operation should be performed after coarsening
677
678
    CoarsenOperation coarsenOperation;

679
    /// Used by \ref interpol
680
681
    AbstractFunction<T, WorldVector<double> > *interFct;

682
    /// Used for mesh traversal
683
684
    static DOFVector<T> *traverseVector;

685
686
687
688
689
#ifdef HAVE_PARALLEL_AMDIS
    /// Petsc application ordering to map dof indices.
    AO *applicationOrdering;
#endif

690
  protected:
691
    /// Used while calculating vector norms
692
693
    static FastQuadrature *quad_fast;

694
    /// Stores the last calculated vector norm
695
696
697
    static double norm;
  }; 

Thomas Witkowski's avatar
Thomas Witkowski committed
698
699
  // ===============================================================================

700
701
702
703
704
705
  template<>
  void DOFVector<double>::refineInterpol(RCNeighbourList&, int); 

  template<>
  void DOFVector<double>::coarseRestrict(RCNeighbourList&, int); 

Thomas Witkowski's avatar
Thomas Witkowski committed
706
707
  inline double min(const DOFVector<double>& v) {
    return v.min();
708
  } 
Thomas Witkowski's avatar
Thomas Witkowski committed
709
710
711

  inline double max(const DOFVector<double>& v) {
    return v.max();
712
  }
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731

  // ===========================================================================
  // ===== class DOFVectorDOF ==================================================
  // ===========================================================================

  /** \ingroup DOFAdministration
   * \brief
   * A DOFVector that stores DOF indices.
   */
  class DOFVectorDOF : public DOFVector<DegreeOfFreedom>,
		       public DOFContainer
  {
  public:  
    MEMORY_MANAGED(DOFVectorDOF);

    /** \brief
     * Calls constructor of DOFVector<DegreeOfFreedom> and registers itself
     * as DOFContainer at DOFAdmin
     */
732
    DOFVectorDOF(const FiniteElemSpace* feSpace_, std::string name_)
733
734
735
      : DOFVector<DegreeOfFreedom>(feSpace_, name_)
    {
      feSpace->getAdmin()->addDOFContainer(this);
736
    }
737
  
738
    /// Deregisters itself at DOFAdmin.
739
740
    ~DOFVectorDOF() {
      feSpace->getAdmin()->removeDOFContainer(this);
741
    }
742
743
744
745
746
747
748

    /** \brief
     * Implements DOFContainer::operator[]() by calling 
     * DOFVEctor<DegreeOfFreedom>::operator[]()
     */
    DegreeOfFreedom& operator[](DegreeOfFreedom i) {
      return DOFVector<DegreeOfFreedom>::operator[](i);
749
    }
750

751
    /// Implements DOFIndexedBase::getSize()
752
753
    int getSize() const {
      return DOFVector<DegreeOfFreedom>::getSize();
754
    }
755

756
    /// Implements DOFIndexedBase::resize()
757
758
759
760
761
762
763
764
765
766
    void resize(int size) {
      DOFVector<DegreeOfFreedom>::resize(size);
    }

    void freeDOFContent(DegreeOfFreedom dof);

  protected:
    DOFVectorDOF();
  };

Thomas Witkowski's avatar
Thomas Witkowski committed
767
768
769
770
771


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


772
773
774
  template<typename T>
  double norm(DOFVector<T> *vec) {
    return vec->nrm2();
775
  }
776
777
778
779

  template<typename T>
  double L2Norm(DOFVector<T> *vec) {
    return vec->L2Norm();
780
  }
781
782
783
784

  template<typename T>
  double H1Norm(DOFVector<T> *vec) {
    return vec->H1Norm();
785
  }
786
787
788
789

  template<typename T>
  void print(DOFVector<T> *vec) {
    vec->print();
790
  }
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843

  // point wise multiplication
  template<typename T>
  const DOFVector<T>& operator*=(DOFVector<T>& x, const DOFVector<T>& y);

  // multiplication with scalar

  template<typename T>
  const DOFVector<T>& operator*=(DOFVector<T>& x, T scal);

  // scalar product

  template<typename T>
  T operator*(DOFVector<T>& x, DOFVector<T>& y);

  // addition

  template<typename T>
  const DOFVector<T>& operator+=(DOFVector<T>& x, const DOFVector<T>& y);

  // subtraction

  template<typename T>
  const DOFVector<T>& operator-=(DOFVector<T>& x, const DOFVector<T>& y);

  template<typename T>
  const DOFVector<T>& operator*(const DOFVector<T>& v, double d);

  template<typename T>
  const DOFVector<T>& operator*(double d, const DOFVector<T>& v);

  template<typename T>
  const DOFVector<T>& operator+(const DOFVector<T>&v1 , const DOFVector<T>& v2);



  // y = a*x + y

  template<typename T>
  void axpy(double a,const DOFVector<T>& x, DOFVector<T>& y);

  // matrix vector product
  template<typename T>
  void mv(MatrixTranspose transpose, 
	  const DOFMatrix &a, 
	  const DOFVector<T> &x,
	  DOFVector<T> &result,
	  bool add = false); 

  template<typename T>
  void xpay(double a,const DOFVector<T>& x,DOFVector<T>& y);

  template<typename T>
844
845
846
  inline void scal(T a, DOFVector<T>& y) {
    y *= a;
  }
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862

  template<typename T>
  inline const DOFVector<T>& mult(double scal,
				  const DOFVector<T>& v,
				  DOFVector<T>& result);

  template<typename T>
  inline const DOFVector<T>& add(const DOFVector<T>& v,
				 double scal,
				 DOFVector<T>& result);

  template<typename T>
  inline const DOFVector<T>& add(const DOFVector<T>& v1,
				 const DOFVector<T>& v2,
				 DOFVector<T>& result);

863
864
865
866
867
868
869
870
871
  template<typename T>
  inline const DOFVector<T>& mod(const DOFVector<T>& v,
				 DOFVector<T>& result);

  template<typename T>
  inline const DOFVector<T>& Tanh(const DOFVector<T>& v,
				  DOFVector<T>& result);


872
  template<typename T>
873
  inline void set(DOFVector<T>& vec, T d) {
874
    vec.set(d);
875
  }
876
877

  template<typename T>
878
  inline void setValue(DOFVector<T>& vec, T d) {
879
    vec.set(d);
880
  }
881
882
883
884

  template<typename T>
  inline int size(DOFVector<T> *vec) {
    return vec->getUsedSize();
885
  }
886
887
888
889
890
891
892
893

  WorldVector<DOFVector<double>*> *transform(DOFVector<WorldVector<double> > *vec,
					     WorldVector<DOFVector<double>*> *result);
}

#include "DOFVector.hh"

#endif  // !_DOFVECTOR_H_