DOFVector.h 22.4 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
// ============================================================================
// ==                                                                        ==
// == 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 ============================================================
// ===========================================================================

#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
41
#include "BasisFunction.h"
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
#include <vector> 
#include <memory> 
#include <list> 
 
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;

66
  template<typename ResultType, typename ArgumentType> class AbstractFunction;
67
68
69
70
71
72


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

74
75
76
    DOFVectorBase() 
      : feSpace(NULL),
	elementVector(NULL),
Thomas Witkowski's avatar
Thomas Witkowski committed
77
78
        boundaryManager(NULL),
        nBasFcts(0)
79
    {}
Thomas Witkowski's avatar
Thomas Witkowski committed
80
    
81
    DOFVectorBase(const FiniteElemSpace *f, std::string n);
82

Thomas Witkowski's avatar
Thomas Witkowski committed
83
    virtual ~DOFVectorBase();
84

85
86
    virtual const T *getLocalVector(const Element *el, 
				    T* localVec) const;
87

88
89
    const T *getVecAtQPs(const ElInfo *elInfo, 
			 const Quadrature *quad,
90
			 const FastQuadrature *quadFast,
91
			 T *vecAtQPs) const;
92

93
94
    const WorldVector<T> *getGrdAtQPs(const ElInfo *elInfo, 
				      const Quadrature *quad,
95
				      const FastQuadrature *quadFast,
96
				      WorldVector<T> *grdAtQPs) const;
97

98
99
    const WorldMatrix<T> *getD2AtQPs(const ElInfo *elInfo, 
				     const Quadrature *quad,
100
				     const FastQuadrature *quadFast,
101
				     WorldMatrix<T> *d2AtQPs) const;
102

103
    inline const FiniteElemSpace* getFESpace() const {
104
      return feSpace;
105
    }
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122

    ElementVector *assemble(T factor, ElInfo *elInfo,
			    const BoundaryType *bound, 
			    Operator *op = NULL);

    void addElementVector(T sign,
			  const ElementVector &elVec, 
			  const BoundaryType *bound,
			  bool add = true); 
 
    inline void addOperator(Operator* op, 
			    double *factor = NULL,
			    double *estFactor = NULL) 
    {
      operators.push_back(op);
      operatorFactor.push_back(factor);
      operatorEstFactor.push_back(estFactor);
123
    }
124

125
    inline std::vector<double*>::iterator getOperatorFactorBegin() {
126
      return operatorFactor.begin();
127
    }
128

129
    inline std::vector<double*>::iterator getOperatorFactorEnd() {
130
      return operatorFactor.end();
131
    }
132

133
    inline std::vector<double*>::iterator getOperatorEstFactorBegin() {
134
      return operatorEstFactor.begin();
135
    }
136

137
    inline std::vector<double*>::iterator getOperatorEstFactorEnd() {
138
      return operatorEstFactor.end();
139
    }
140
141


142
    inline std::vector<Operator*>::iterator getOperatorsBegin() {
143
      return operators.begin();
144
    }
145

146
    inline std::vector<Operator*>::iterator getOperatorsEnd() {
147
      return operators.end();
148
    }
149
150
151
152
153
154
155
156
157
158
159

    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);

160
    inline std::vector<Operator*>& getOperators() { 
Thomas Witkowski's avatar
Thomas Witkowski committed
161
      return operators; 
162
    }
Thomas Witkowski's avatar
Thomas Witkowski committed
163

164
    inline std::vector<double*>& getOperatorFactor() { 
Thomas Witkowski's avatar
Thomas Witkowski committed
165
      return operatorFactor; 
166
    }
167

168
    inline std::vector<double*>& getOperatorEstFactor() { 
169
      return operatorEstFactor; 
170
    }
171
172
173
174

    /** \brief
     * Returns \ref name
     */
175
    inline const std::string& getName() const { 
Thomas Witkowski's avatar
Thomas Witkowski committed
176
      return name; 
177
    } 
178

Thomas Witkowski's avatar
Thomas Witkowski committed
179
180
    inline BoundaryManager* getBoundaryManager() const { 
      return boundaryManager; 
181
    }
182
183
184

    inline void setBoundaryManager(BoundaryManager *bm) {
      boundaryManager = bm;
185
    }
186
187

  protected:
Thomas Witkowski's avatar
Thomas Witkowski committed
188
189
190
    /** \brief
     *
     */
191
192
    const FiniteElemSpace *feSpace;

Thomas Witkowski's avatar
Thomas Witkowski committed
193
194
195
    /** \brief
     *
     */
196
    std::string name;
197

Thomas Witkowski's avatar
Thomas Witkowski committed
198
199
200
    /** \brief
     *
     */
201
202
    ElementVector *elementVector;

Thomas Witkowski's avatar
Thomas Witkowski committed
203
204
205
    /** \brief
     *
     */
206
    std::vector<Operator*> operators;
207

Thomas Witkowski's avatar
Thomas Witkowski committed
208
209
210
    /** \brief
     *
     */
211
    std::vector<double*> operatorFactor;
212

Thomas Witkowski's avatar
Thomas Witkowski committed
213
214
215
    /** \brief
     *
     */
216
    std::vector<double*> operatorEstFactor;
217

Thomas Witkowski's avatar
Thomas Witkowski committed
218
219
220
    /** \brief
     *
     */
221
222
    BoundaryManager *boundaryManager;

Thomas Witkowski's avatar
Thomas Witkowski committed
223
224
225
226
    /** \brief
     * Number of basis functions of the used finite element space.
     */
    int nBasFcts;
227

Thomas Witkowski's avatar
Thomas Witkowski committed
228
229
230
231
    /** \brief
     * Are used to store temporary local dofs of an element. 
     * Thread safe.
     */
232
233
    std::vector<DegreeOfFreedom*> localIndices;

Thomas Witkowski's avatar
Thomas Witkowski committed
234
235
236
237
238
239
240
241
242
243
244
    /** \brief
     * Are used to store temporary local values of an element. 
     * Thread safe.
     */
    std::vector<T*> localVectors;

    /** \brief
     * Temporarly used in \ref getGrdAtQPs. Thread safe.
     */
    std::vector<DimVec<double>*> grdTmp;

245
246
247
248
249
250
251
252
253
    /** \brief
     * Temporarly used in \ref getGrdAtQPs. Thread safe.
     */
    std::vector<DimVec<double>*> grdPhis;

    /** \brief
     * Temporarly used in \ref getD2AtQPs. Thread safe.
     */
    std::vector<DimMat<double>*> D2Phis;
Thomas Witkowski's avatar
Thomas Witkowski committed
254
255
256
257
258

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


  // ===========================================================================
  // ===== 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)
305
      {}
306
307
308

      Iterator(DOFAdmin *admin, DOFIndexed<T> *c, DOFIteratorType type)
	: DOFIterator<T>(admin, c, type)
309
      {}
310
311
312
313
314
315
    };

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

316
317
318
      Creator(FiniteElemSpace *feSpace_) 
        : feSpace(feSpace_) 
      {}
319
320
321

      DOFVector<T> *create() { 
	return NEW DOFVector<T>(feSpace, ""); 
322
      }
323
324
325

      void free(DOFVector<T> *vec) { 
	DELETE vec; 
326
      }
327
328
329
330
331
332
333
334
335
336
337
338
339

    private:
      FiniteElemSpace *feSpace;
    };

  public:
    /** \brief
     * Empty constructor. No initialization!
     */
    DOFVector() 
      : DOFVectorBase<T>(), 
	feSpace(NULL),
	coarsenOperation(NO_OPERATION)
340
    {}
341
342
343
344

    /** \brief
     * Constructs a DOFVector with name n belonging to FiniteElemSpace f
     */
345
    DOFVector(const FiniteElemSpace* f, std::string n); 
346
347
348
349

    /** \brief
     * Initialization.
     */
350
    void init(const FiniteElemSpace* f, std::string n);
351
352
353
354
355

    /** \brief
     * Copy Constructor
     */
    DOFVector(const DOFVector& rhs) {
Thomas Witkowski's avatar
Thomas Witkowski committed
356
357
358
      *this = rhs;   
      name = rhs.name + "copy";
      if (feSpace && feSpace->getAdmin()) {
359
360
	(dynamic_cast<DOFAdmin*>(feSpace->getAdmin()))->addDOFIndexed(this);
      }
361
    }
362
363
364
365
366
367
368
369
370

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

    /** \brief
     * Returns iterator to the begin of \ref vec
     */
371
    typename std::vector<T>::iterator begin() { 
372
      return vec.begin(); 
373
    }
374
375
376
377

    /** \brief
     * Returns iterator to the end of \ref vec
     */
378
    typename std::vector<T>::iterator end() { 
379
      return vec.end(); 
380
    }
381
382
383
384
385
386
  
    /** \brief
     * Used by DOFAdmin to compress this DOFVector. Implementation of
     * DOFIndexedBase::compress()
     */
    virtual void compressDOFIndexed(int first, int last,
387
				    std::vector<DegreeOfFreedom> &newDof);
388
389
390
391
392
393

    /** \brief
     * Sets \ref coarsenOperation to op
     */
    inline void setCoarsenOperation(CoarsenOperation op) { 
      coarsenOperation = op; 
394
    }
395
396
397
398

    /** \brief
     * Returns \ref coarsenOperation
     */
399
400
    inline CoarsenOperation getCoarsenOperation() { 
      return coarsenOperation; 
401
    }
402
403
404
405

    /** \brief
     * Restriction after coarsening. Implemented for DOFVector<double>
     */
406
    inline void coarseRestrict(RCNeighbourList&, int) {}
407
408
409
410

    /** \brief
     * Interpolation after refinement.
     */
411
    inline void refineInterpol(RCNeighbourList&, int) {}
412
413
414
415

    /** \brief
     * Returns \ref vec
     */
416
    std::vector<T>& getVector() { 
Thomas Witkowski's avatar
Thomas Witkowski committed
417
      return vec;
418
    }
419
420
421
422

    /** \brief
     * Returns size of \ref vec
     */
423
424
    inline int getSize() const { 
      return vec.size();
425
    } 
426
427
428
429

    /** \brief
     * Returns used size of the vector.
     */
430
431
    inline int getUsedSize() const { 
      return feSpace->getAdmin()->getUsedSize(); 
432
    }
433
434
435
436
437

    /** \brief
     * Resizes \ref vec to n
     */
    inline void resize(int n) { 
438
439
      FUNCNAME("DOFVector<T>::resize()");
      TEST_EXIT_DBG((n >= 0))("Can't resize DOFVector to negative size\n"); 
440
      vec.resize(n);
441
    } 
442
443
444
445
446

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

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

462
463
    inline int getSize() { 
      return vec.size(); 
464
    }
465
466
467
468
469
470

    /** \brief
     * Returns \ref vec[i]
     */
    inline T& operator[](DegreeOfFreedom i) {
      FUNCNAME("DOFVector<T>::operator[]");
471
472
      TEST_EXIT_DBG((i >= 0) && (i < static_cast<int>(vec.size())))
	("Illegal vector index %d.\n",i);
473
      return vec[i];
474
    }
Thomas Witkowski's avatar
Thomas Witkowski committed
475
 
476
477
478
479
480
481
482
483
484
    /** \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
485
 
486
487
488
489
490
    /** \brief
     * Calculates L2 norm of this DOFVector
     */
    inline double L2Norm(Quadrature* q = NULL) const {
      return sqrt(L2NormSquare());
491
    }
492
493
494
495
496
497
498
499
500
501
502

    /** \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());
503
    }; 
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524

    /** \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();
525
    }
526
527
528
529
530
531
532
533
534
535
536

    /** \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();
537
    } 
538
539
540
541
542

    /** \brief
     * Calculates doublewell of this DOFVector
     */
    double DoubleWell(Quadrature* q = NULL) const;
Thomas Witkowski's avatar
Thomas Witkowski committed
543
 
544
545
546
547
    /** \brief
     * Calculates the sum of this DOFVector
     */
    T sum() const; 
Thomas Witkowski's avatar
Thomas Witkowski committed
548
 
549
550
551
552
553
554
555
556
557
558
559
    /** \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;
560
    } 
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579

    /** \brief
     * Assignment operator between two vectors
     */
    DOFVector<T>& operator=(const DOFVector<T>& ); 

    /** \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
580
581
582
583
584
585
    T max() const;

    /** \brief
     * Returns absolute maximum of DOFVector
     */
    T absMax() const;
586
587
588
589
590

    /** \brief
     * Used by interpol while mesh traversal
     */
    static int interpolFct(ElInfo* elinfo);
591
 
592
593
594
595
596
    /** \brief
     * Prints \ref vec to stdout
     */
    void print() const; 

597
598
599
    /** \brief
     *
     */
600
601
    int calcMemoryUsage() const;

602
603
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);

    // ===== Serializable implementation =====

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

Thomas Witkowski's avatar
Thomas Witkowski committed
633
  protected: 
634
635
636
637
638
639
640
641
642
    /** \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
643
 
644
645
646
647
648
649
650
651
652
653
    /** \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);

654
655
656
657
658
    /** \brief
     * Used by DoubleWell while mesh traversal
     */
    static int DoubleWell_fct(ElInfo* elinfo);

659
660
661
662
  protected: 
    /** \brief
     * Name of this DOFVector
     */
663
    std::string name; 
664
665
666
667
668
669
670
671
672

    /** \brief
     * FiniteElemSpace of the vector
     */
    const FiniteElemSpace *feSpace; 

    /** \brief
     * Data container
     */
673
    std::vector<T> vec; 
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
 
    /** \brief
     * Specifies what operation should be performed after coarsening
     */
    CoarsenOperation coarsenOperation;

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

    /** \brief
     * Used for mesh traversal
     */
    static DOFVector<T> *traverseVector;

  protected:
    /** \brief
     * Used while calculating vector norms
     */
    static FastQuadrature *quad_fast;

    /** \brief
     * Stores the last calculated vector norm
     */
    static double norm;

    /** \brief
     * Dimension of the mesh this DOFVector belongs to
     */
    static int dim;
  }; 

Thomas Witkowski's avatar
Thomas Witkowski committed
707
708
  // ===============================================================================

709
710
711
712
713
714
  template<>
  void DOFVector<double>::refineInterpol(RCNeighbourList&, int); 

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

Thomas Witkowski's avatar
Thomas Witkowski committed
715
716
  inline double min(const DOFVector<double>& v) {
    return v.min();
717
  } 
Thomas Witkowski's avatar
Thomas Witkowski committed
718
719
720

  inline double max(const DOFVector<double>& v) {
    return v.max();
721
  }
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740

  // ===========================================================================
  // ===== 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
     */
741
    DOFVectorDOF(const FiniteElemSpace* feSpace_, std::string name_)
742
743
744
      : DOFVector<DegreeOfFreedom>(feSpace_, name_)
    {
      feSpace->getAdmin()->addDOFContainer(this);
745
    }
746
747
748
749
750
751
  
    /** \brief
     * Deregisters itself at DOFAdmin.
     */
    ~DOFVectorDOF() {
      feSpace->getAdmin()->removeDOFContainer(this);
752
    }
753
754
755
756
757
758
759

    /** \brief
     * Implements DOFContainer::operator[]() by calling 
     * DOFVEctor<DegreeOfFreedom>::operator[]()
     */
    DegreeOfFreedom& operator[](DegreeOfFreedom i) {
      return DOFVector<DegreeOfFreedom>::operator[](i);
760
    }
761
762
763
764
765
766

    /** \brief
     * Implements DOFIndexedBase::getSize()
     */
    int getSize() const {
      return DOFVector<DegreeOfFreedom>::getSize();
767
    }
768
769
770
771
772
773
774
775
776
777
778
779
780
781

    /** \brief
     * Implements DOFIndexedBase::resize()
     */
    void resize(int size) {
      DOFVector<DegreeOfFreedom>::resize(size);
    }

    void freeDOFContent(DegreeOfFreedom dof);

  protected:
    DOFVectorDOF();
  };

Thomas Witkowski's avatar
Thomas Witkowski committed
782
783
784
785
786


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


787
788
789
  template<typename T>
  double norm(DOFVector<T> *vec) {
    return vec->nrm2();
790
  }
791
792
793
794

  template<typename T>
  double L2Norm(DOFVector<T> *vec) {
    return vec->L2Norm();
795
  }
796
797
798
799

  template<typename T>
  double H1Norm(DOFVector<T> *vec) {
    return vec->H1Norm();
800
  }
801
802
803
804

  template<typename T>
  void print(DOFVector<T> *vec) {
    vec->print();
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
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858

  // 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>
859
860
861
  inline void scal(T a, DOFVector<T>& y) {
    y *= a;
  }
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877

  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);

878
879
880
881
882
883
884
885
886
  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);


887
  template<typename T>
888
  inline void set(DOFVector<T>& vec, T d) {
889
    vec.set(d);
890
  }
891
892

  template<typename T>
893
  inline void setValue(DOFVector<T>& vec, T d) {
894
    vec.set(d);
895
  }
896
897
898
899

  template<typename T>
  inline int size(DOFVector<T> *vec) {
    return vec->getUsedSize();
900
  }
901
902
903
904
905
906
907
908

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

#include "DOFVector.hh"

#endif  // !_DOFVECTOR_H_