DOFVector.hh 38.7 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
#include <list>

#include "FixVec.h"
#include "Boundary.h"
#include "DOFAdmin.h"
#include "ElInfo.h"
#include "Error.h"
#include "FiniteElemSpace.h"
#include "Global.h"
#include "Mesh.h"
#include "Quadrature.h"
#include "AbstractFunction.h"
#include "BoundaryManager.h"
#include "ElementVector.h"
#include "Assembler.h"

namespace AMDiS {

  template<typename T>
  void DOFVectorBase<T>::addElementVector(T factor, 
					  const ElementVector &elVec, 
					  const BoundaryType *bound,
					  bool add)
  {
25
    FUNCNAME("DOFVector::addElementVector()");
26
27
28

    int n_row = elVec.getSize();

29
    for (DegreeOfFreedom i = 0; i < n_row; i++) {
30
31
32
33
      BoundaryCondition *condition = 
	bound ? this->getBoundaryManager()->getBoundaryCondition(bound[i]) : NULL;

      if(!(condition && condition->isDirichlet())) {
34
	DegreeOfFreedom irow = elVec.dofIndices[i];
35
36
37
38
39
40
41
42
43
	(*this)[irow] = (add ? (*this)[irow] : 0.0);
	(*this)[irow] += factor * elVec[i];
      }
    }
  }

  template<typename T>
  DOFVector<T>::DOFVector(const FiniteElemSpace* f,::std::string n)
    : DOFVectorBase<T>(f, n),
44
45
      refineInter(false), 
      coarsenOperation(NO_OPERATION)
46
47
48
49
50
51
52
53
54
  {
    init(f, n);
  } 

  template<typename T>
  void DOFVector<T>::init(const FiniteElemSpace* f,::std::string n)
  {
    name = n;
    feSpace = f;
55
    if (feSpace && feSpace->getAdmin()) {
56
57
58
59
60
61
62
63
      (feSpace->getAdmin())->addDOFIndexed(this);
    }  
    this->boundaryManager = NEW BoundaryManager;
  }

  template<typename T>
  DOFVector<T>::~DOFVector()
  {
64
    if (feSpace && feSpace->getAdmin()) {
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
      (feSpace->getAdmin())->removeDOFIndexed(this);
    }

    if (this->boundaryManager) {
      DELETE this->boundaryManager;
    }
  }

  template<typename T>
  DOFVector<T> * DOFVector<T>::traverseVector = NULL;

  template<typename T>
  FastQuadrature *DOFVector<T>::quad_fast = NULL;

  template<typename T>
  double DOFVector<T>::norm = 0.0;

  template<typename T>
  int DOFVector<T>::dim = 0;

  template<typename T>
  double DOFVector<T>::nrm2() const
  {
    FUNCNAME("DOFVector<T>::nrm2()");

90
91
92
    TEST_EXIT_DBG(feSpace)("feSpace is NULL\n");
    TEST_EXIT_DBG(feSpace->getAdmin())("admin is NULL\n");
    TEST_EXIT_DBG(static_cast<int>(vec.size()) >= feSpace->getAdmin()->getUsedSize())
93
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(),
94
95
96
       feSpace->getAdmin()->getUsedSize());
    
    double nrm = 0.0;
97
98
99
100
101
102
103
104
105
106
107
108
    Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(this)), USED_DOFS);
    for(vecIterator.reset(); !vecIterator.end(); ++vecIterator)
      nrm += (*vecIterator) * (*vecIterator);

    return(sqrt(nrm));
  }

  template<typename T>
  double DOFVector<T>::squareNrm2() const
  {
    FUNCNAME("DOFVector<T>::nrm2()");

109
110
111
    TEST_EXIT_DBG(feSpace)("feSpace is NULL\n");
    TEST_EXIT_DBG(feSpace->getAdmin())("admin is NULL\n");
    TEST_EXIT_DBG(static_cast<int>(vec.size()) >= feSpace->getAdmin()->getUsedSize())
112
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(),
113
114
115
       feSpace->getAdmin()->getUsedSize());
    
    double nrm = 0.0;
116
117
118
119
120
121
122
123
124
125
    Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(this)), USED_DOFS);
    for(vecIterator.reset(); !vecIterator.end(); ++vecIterator)
      nrm += (*vecIterator) * (*vecIterator);

    return nrm;
  }

  template<typename T>
  T DOFVector<T>::asum() const
  {
126
    FUNCNAME("DOFVector<T>::asum()");
127

128
129
130
131
132
    TEST_EXIT_DBG(feSpace)("feSpace is NULL\n");
    TEST_EXIT_DBG(feSpace->getAdmin())("admin is NULL\n");
    TEST_EXIT_DBG(static_cast<int>(vec.size()) >= feSpace->getAdmin()->getUsedSize())
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(),
       feSpace->getAdmin()->getUsedSize());
133

134
    double nrm = 0.0;
135
136
137
138
139
140
141
142
143
144
    Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(this)), USED_DOFS);
    for(vecIterator.reset(); !vecIterator.end(); ++vecIterator)
      nrm += abs(*vecIterator);

    return(nrm);
  }

  template<typename T>
  T DOFVector<T>::sum() const
  {
145
    FUNCNAME("DOFVector<T>::sum()");
146

147
148
149
150
151
    TEST_EXIT_DBG(feSpace)("feSpace is NULL\n");
    TEST_EXIT_DBG(feSpace->getAdmin())("admin is NULL\n");
    TEST_EXIT_DBG(static_cast<int>(vec.size()) >= feSpace->getAdmin()->getUsedSize())
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(),
       feSpace->getAdmin()->getUsedSize());
152

153
    double nrm = 0.0;    
154
155
156
157
158
159
160
161
162
163
    Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(this)), USED_DOFS);
    for(vecIterator.reset(); !vecIterator.end(); ++vecIterator)
      nrm += *vecIterator;

    return(nrm);
  }

  template<typename T>
  void DOFVector<T>::set(T alpha)
  {
164
    FUNCNAME("DOFVector<T>::set()");
165

166
167
168
    TEST_EXIT_DBG(feSpace)("feSpace is NULL\n");
    TEST_EXIT_DBG(feSpace->getAdmin())("admin is NULL\n");
    TEST_EXIT_DBG(static_cast<int>(vec.size()) >= feSpace->getAdmin()->getUsedSize())
169
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(),
170
171
       feSpace->getAdmin()->getUsedSize());
    
172
173
174
175
176
177
178
179
180
181
    Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(this), USED_DOFS);
    for(vecIterator.reset(); !vecIterator.end(); ++vecIterator) {
      *vecIterator = alpha ; 
    };
  }


  template<typename T>
  void DOFVector<T>::copy(const DOFVector<T>& x)
  {
182
183
184
185
186
    FUNCNAME("DOFVector<T>::copy()");
    
    TEST_EXIT_DBG(feSpace && x.feSpace)
      ("feSpace is NULL: %8X, %8X\n", feSpace, x.feSpace);
    TEST_EXIT_DBG(feSpace->getAdmin() && (feSpace->getAdmin() == x.feSpace->getAdmin()))
187
188
      ("no admin or different admins: %8X, %8X\n",
       feSpace->getAdmin(), x.feSpace->getAdmin());
189
    TEST_EXIT_DBG(static_cast<int>(vec.size()) >= feSpace->getAdmin()->getUsedSize())
190
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), 
191
192
       feSpace->getAdmin()->getUsedSize());
    TEST_EXIT_DBG(static_cast<int>(x.vec.size()) >= feSpace->getAdmin()->getUsedSize())
193
      ("x.size = %d too small: admin->sizeUsed = %d\n", x.vec.size(), 
194
195
       feSpace->getAdmin()->getUsedSize());
    
196
197
    Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(this), USED_DOFS);
    Iterator xIterator(dynamic_cast<DOFVector<T>*>(const_cast<DOFVector<T>*>(&x)), USED_DOFS);
198
199
200
201
202
203
    for (vecIterator.reset(), xIterator.reset();
	 !vecIterator.end();
	 ++vecIterator, ++xIterator) {
      
      *vecIterator = *xIterator; 
    };
204
205
206
207
208
209
  }


  template<typename T>
  T DOFVector<T>::min() const
  {
210
211
212
213
214
    FUNCNAME("DOFVector<T>::min()");
    
    TEST_EXIT_DBG(feSpace && feSpace->getAdmin())
      ("pointer is NULL: %8X, %8X\n", this, feSpace->getAdmin());
    TEST_EXIT_DBG((static_cast<int>(vec.size())) >= feSpace->getAdmin()->getUsedSize())
215
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), 
216
217
218
       feSpace->getAdmin()->getUsedSize());

    T m;    
219
    Iterator vecIterator(const_cast<DOFIndexed<T>*>(dynamic_cast<const DOFIndexed<T>*>(this)), USED_DOFS);
220
221
222
223
    for (vecIterator.reset(), m = *vecIterator; !vecIterator.end(); ++vecIterator) {
      m = ::std::min(m, *vecIterator);
    }

224
225
226
227
228
229
    return m;
  }

  template<typename T>
  T DOFVector<T>::max() const 
  {
230
231
232
233
234
    FUNCNAME("DOFVector<T>::max()");
    
    TEST_EXIT_DBG(feSpace && feSpace->getAdmin())
      ("pointer is NULL: %8X, %8X\n", this, feSpace->getAdmin());
    TEST_EXIT_DBG((static_cast<int>(vec.size())) >= feSpace->getAdmin()->getUsedSize())
235
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), 
236
       feSpace->getAdmin()->getUsedSize());
237

238
    T m;    
239
    Iterator vecIterator(const_cast<DOFIndexed<T>*>(dynamic_cast<const DOFIndexed<T>*>(this)), USED_DOFS);
240
241
242
243
    for (vecIterator.reset(), m = *vecIterator; !vecIterator.end(); ++vecIterator) {
      m = ::std::max(m, *vecIterator);
    }

244
245
246
247
248
249
250
251
    return m;
  }

  template<typename T>
  void gemv(MatrixTranspose transpose, T alpha,
	    const DOFMatrix& a, const DOFVector<T>& x,
	    T beta, DOFVector<T>& y)
  {
252
253
254
255
    FUNCNAME("DOFVector<T>::gemv()");

    int j, jcol, ysize;
    T sum, ax;
256
    const DOFMatrix::MatrixRow *row;
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
    
    TEST_EXIT_DBG(a.getRowFESpace() && 
		  a.getColFESpace() && 
		  x.getFESpace() && 
		  y.getFESpace())
      ("feSpace is NULL: %8X, %8X, %8X, %8X\n", 
       a.getRowFESpace(), a.getColFESpace(), x.getFESpace(), y.getFESpace());
    TEST_EXIT_DBG(a.getRowFESpace()->getAdmin() && 
		  a.getColFESpace()->getAdmin() && 
		  (((transpose == NoTranspose) && 
		    (a.getColFESpace()->getAdmin() == x.getFESpace()->getAdmin()) && 
		    (a.getRowFESpace()->getAdmin() == y.getFESpace()->getAdmin()))||
		   ((transpose == Transpose) && 
		    (a.getRowFESpace()->getAdmin() == x.getFESpace()->getAdmin()) && 
		    (a.getColFESpace()->getAdmin() == y.getFESpace()->getAdmin()))))
272
273
274
      ("no admin or different admins: %8X, %8X, %8X, %8X\n",
       a.getRowFESpace()->getAdmin(), a.getColFESpace()->getAdmin(), 
       x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin());
275
    
276
    if (transpose == NoTranspose) {
277
      TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize())
278
	("x.size = %d too small: admin->sizeUsed = %d\n",
279
280
	 x.getSize(), a.getColFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>(y.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
281
	("y.size = %d too small: admin->sizeUsed = %d\n",
282
283
	 y.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
284
	("a.size = %d too small: admin->sizeUsed = %d\n",
285
	a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
286
287
    }
    else  if (transpose == Transpose) {
288
      TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
289
	("x.size = %d too small: admin->sizeUsed = %d\n",
290
291
	 x.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>(y.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize())
292
	("y.size = %d too small: admin->sizeUsed = %d\n",
293
294
	 y.getSize(), a.getColFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
295
	("a.size = %d too small: admin->sizeUsed = %d\n",
296
	a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
    }

    ysize = y.getSize();

    typename DOFVector<T>::Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(&y), FREE_DOFS);
    for(vecIterator.reset();
	!vecIterator.end(); ++vecIterator) { 
      *vecIterator = 0;
    };

    if (transpose == NoTranspose) {
      typename DOFVector<T>::Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(&y), USED_DOFS); 
      DOFMatrix::Iterator rowIterator(const_cast<DOFMatrix*>(&a), USED_DOFS);
      for(vecIterator.reset(), rowIterator.reset();
	  !rowIterator.end();
	  ++rowIterator, ++vecIterator) { 
	sum = 0;
	row = &(a[rowIterator.getDOFIndex()]);
	for(::std::vector<MatEntry>::iterator colIterator = rowIterator->begin();
	    colIterator != rowIterator->end();
	    colIterator++) {
	  jcol = colIterator->col;
	  if (jcol >= 0) { // entry used? 
	    sum += (static_cast<T>(colIterator->entry)) * x[jcol];
	  } else {
	    if (jcol == DOFMatrix::NO_MORE_ENTRIES)
	      break;
	  }
	}
	*vecIterator *= beta;
	*vecIterator += alpha * sum;
      };
    } else if (transpose == Transpose) {
      typename DOFVector<T>::Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(&y), USED_DOFS);
      for(vecIterator.reset();
	  !vecIterator.end(); 
	  ++vecIterator) {
	*vecIterator  *= beta ; 
      };
    
      typename DOFVector<T>::Iterator xIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(&x)), USED_DOFS);
      DOFMatrix::Iterator rowIterator(const_cast<DOFMatrix*>(&a), USED_DOFS);
      for(xIterator.reset(), rowIterator.reset();
	  !rowIterator.end();
	  ++rowIterator, ++xIterator) { 
	ax = alpha * (*xIterator);
	row = &(a[rowIterator.getDOFIndex()]);
	for(::std::vector<MatEntry>::iterator colIterator = rowIterator->begin();
	    colIterator != rowIterator->end();
	    colIterator++) {
	  jcol = colIterator->col;
	  if (jcol >= 0) // entry used?
	    y[jcol] += ax * (static_cast<T>(colIterator->entry));
	  else 
	    if (jcol == DOFMatrix::NO_MORE_ENTRIES) break;
	}
      }
    }
    else {
      ERROR_EXIT("transpose=%d\n", transpose);
    }
  }

  template<typename T>
  void DOFVector<T>::print() const
  {
    FUNCNAME("DOFVector<T>::print");
    int  i, j;
    const DOFAdmin *admin = NULL;
    const char     *format;

    if (feSpace) admin = feSpace->getAdmin();

    MSG("Vec `%s':\n", name.c_str());
    j = 0;
    if (admin) {
      if (admin->getUsedSize() > 100)
	format = "%s(%3d,%10.5le)";
      else if (admin->getUsedSize() > 10)
	format = "%s(%2d,%10.5le)";
      else
	format = "%s(%1d,%10.5le)";

      Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(this)), USED_DOFS);
      for(vecIterator.reset();
	  !vecIterator.end(); ++vecIterator) {
	if ((j % 3) == 0) {
	  if (j) Msg::print("\n");
	  MSG(format, "", vecIterator.getDOFIndex(), *vecIterator);
	}
	else 
	  Msg::print(format, " ", vecIterator.getDOFIndex(), *vecIterator);
	j++;
      };
      Msg::print("\n");
    }
    else
      {
	MSG("no DOFAdmin, print whole vector.\n");
    
	for (i = 0; i < static_cast<int>( vec.size()); i++) {
	  if ((j % 3) == 0)
	    {
	      if (j) Msg::print("\n");
	      MSG("(%d,%10.5e)",i,vec[i]);
	    }
	  else 
	    Msg::print(" (%d,%10.5e)",i,vec[i]);
	  j++;
	}
	Msg::print("\n");
      }
    return;
  }


  template<typename T>
  T DOFVectorBase<T>::evalUh(const DimVec<double>& lambda, 
			     DegreeOfFreedom* dof_indices)
  {
    int i;
    BasisFunction* phi = const_cast<BasisFunction*>(this->getFESpace()->getBasisFcts());
    int numberOfBasFcts = phi->getNumber();
    T val = 0.0;

    for (i = 0; i < numberOfBasFcts; i++)
      val += (*this)[dof_indices[i]]*(*phi->getPhi(i))(lambda);

    return val;
  }

  template<typename T>
  void DOFVector<T>::interpol(AbstractFunction<T, WorldVector<double> > *fct)
  {
431
    FUNCNAME("DOFVector<T>::interpol()");
432

433
434
435
    TEST_EXIT_DBG(fct)("no function to interpolate\n");

    interFct = fct;
436

437
438
439
440
441
    if (!this->getFESpace()) {
      MSG("no dof admin in vec %s, skipping interpolation\n",
	  this->getName().c_str());
      return;
    }
442

443
444
445
446
447
    if (!(this->getFESpace()->getAdmin())) {
      MSG("no dof admin in feSpace %s, skipping interpolation\n", 
	  this->getFESpace()->getName().c_str());
      return;
    }
448
  
449
450
451
452
453
    if (!(this->getFESpace()->getBasisFcts())) {
      MSG("no basis functions in admin of vec %s, skipping interpolation\n",
	  this->getName().c_str());
      return;
    }
454

455
456
457
458
459
    if (!(fct)) {
      MSG("function that should be interpolated only pointer to NULL, ");
      Msg::print("skipping interpolation\n");
      return;
    }
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
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
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
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
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752

    traverseVector = this;
    this->getFESpace()->getMesh()->traverse(-1, 
					    Mesh::CALL_LEAF_EL |
					    Mesh::FILL_COORDS, 
					    interpolFct);

    return;
  }

  template<typename T>
  int DOFVector<T>::interpolFct(ElInfo* elinfo)
  {
    int i;
    const BasisFunction *basFct = traverseVector->getFESpace()->getBasisFcts();
    const DOFAdmin* admin = traverseVector->getFESpace()->getAdmin();
    const DegreeOfFreedom   *dof = 
      basFct->getLocalIndices(const_cast<Element *>(elinfo->getElement()), 
			      admin, NULL);
    const T *inter_val = 
      const_cast<BasisFunction*>(basFct)->interpol(elinfo, 
						   0, 
						   NULL,
						   traverseVector->interFct, 
						   NULL);

    int number = basFct->getNumber();
    for (i = 0; i < number; i++)
      (*traverseVector)[dof[i]] = inter_val[i];

    return 0;
  }

  // das hat Andreas eingefuegt: Integral...

  template<typename T>
  double DOFVector<T>::Int(Quadrature* q) const
  {
    FUNCNAME("DOFVector::Int");
  
    Mesh* mesh = feSpace->getMesh();

    int deg;
    dim = mesh->getDim();

    if (!q)
      {
	deg = 2*feSpace->getBasisFcts()->getDegree();
	q = Quadrature::provideQuadrature(dim, deg);
      }

    quad_fast = FastQuadrature::provideFastQuadrature(feSpace->getBasisFcts(), *q, INIT_PHI);

    norm = 0.0;
    traverseVector = const_cast<DOFVector<T>*>(this);

    mesh->traverse(-1, 
		   Mesh::CALL_LEAF_EL|
		   Mesh::FILL_COORDS |
		   Mesh::FILL_DET, 
		   Int_fct);

    return norm;  
  }

  template<typename T>
  int DOFVector<T>::Int_fct(ElInfo *elinfo)
  {
    double det, normT;
    const T *uh_vec;
    int iq;

    det = elinfo->getDet();
    //uh_loc = traverseVector->getLocalVector(elinfo->getElement(), NULL);
    //uh_vec = quad_fast->uhAtQp(uh_loc, NULL);

    uh_vec = traverseVector->getVecAtQPs(elinfo,
					 NULL,
					 quad_fast,
					 NULL);

    int numPoints = quad_fast->getNumPoints();
    for (normT = iq = 0; iq < numPoints; iq++)
      {
	normT += quad_fast->getWeight(iq)*(uh_vec[iq]);
      }
    norm += det*normT;

    return 0;
  }

  // ... und die L1-Norm ...

  template<typename T>
  double DOFVector<T>::L1Norm(Quadrature* q) const
  {
    FUNCNAME("DOFVector::L1Norm");
  
    Mesh* mesh = feSpace->getMesh();

    int deg;
    dim = mesh->getDim();

    if (!q)
      {
	deg = 2*feSpace->getBasisFcts()->getDegree();
	q = Quadrature::provideQuadrature(dim, deg);
      }

    quad_fast = FastQuadrature::provideFastQuadrature(feSpace->getBasisFcts(),*q,INIT_PHI);

    norm = 0.0;
    traverseVector = const_cast<DOFVector<T>*>(this);

    mesh->traverse(-1, 
		   Mesh::CALL_LEAF_EL|
		   Mesh::FILL_COORDS |
		   Mesh::FILL_DET, 
		   L1Norm_fct);

    return norm;  
  }

  template<typename T>
  int DOFVector<T>::L1Norm_fct(ElInfo *elinfo)
  {
    double det, normT;
    const T *uh_loc, *uh_vec;
    int iq;

    det = elinfo->getDet();
    //   uh_loc = 
    //     traverseVector->getLocalVector(elinfo->getElement(), 
    // 				   NULL);
    //   uh_vec = quad_fast->uhAtQp(uh_loc, NULL);

    uh_vec = traverseVector->getVecAtQPs(elinfo,
					 NULL,
					 quad_fast,
					 NULL);

    int numPoints = quad_fast->getNumPoints();
    for (normT = iq = 0; iq < numPoints; iq++)
      {
	normT += quad_fast->getWeight(iq)*abs(uh_vec[iq]);
      }
    norm += det*normT;

    return 0;
  }


  // bis hierhin gehen Andreas Ergaenzungen...

  template<typename T>
  double DOFVector<T>::L2NormSquare(Quadrature* q) const
  {
    FUNCNAME("DOFVector::L2NormSquare");
  
    Mesh* mesh = feSpace->getMesh();

    int deg;
    dim = mesh->getDim();

    if (!q)
      {
	deg = 2*feSpace->getBasisFcts()->getDegree();
	q = Quadrature::provideQuadrature(dim, deg);
      }

    quad_fast = FastQuadrature::provideFastQuadrature(feSpace->getBasisFcts(), 
						      *q, 
						      INIT_PHI);

    norm = 0.0;
    traverseVector = const_cast<DOFVector<T>*>(this);

    mesh->traverse(-1, Mesh::CALL_LEAF_EL|Mesh::FILL_COORDS|Mesh::FILL_DET, L2Norm_fct);

    return norm;  
  }

  template<typename T>
  int DOFVector<T>::L2Norm_fct(ElInfo *elinfo)
  {
    double det, normT;
    const T *uh_vec;
    int iq;

    det = elinfo->getDet();
    //   uh_loc = 
    //     traverseVector->getLocalVector(elinfo->getElement(), 
    // 				   NULL);
    //   uh_vec = quad_fast->uhAtQp(uh_loc, NULL);

    uh_vec = traverseVector->getVecAtQPs(elinfo,
					 NULL,
					 quad_fast,
					 NULL);

    int numPoints = quad_fast->getNumPoints();
    for (normT = iq = 0; iq < numPoints; iq++)
      {
	normT += quad_fast->getWeight(iq)*sqr(uh_vec[iq]);
      }
    norm += det*normT;

    return 0;
  }

  template<typename T>
  int DOFVector<T>::H1Norm_fct(ElInfo *elinfo)
  {
    double norm2, normT;
    //const T *uh_loc;
    const WorldVector<T> *grduh_vec;
    int iq, j;

    double det = elinfo->getDet();
    //const DimVec<WorldVector<double> > &Lambda = elinfo->getGrdLambda();
  
    //det = elinfo->calcGrdLambda(Lambda);
    //uh_loc = 
    //  traverseVector->getLocalVector(elinfo->getElement(), NULL);
    //grduh_vec = 
    //  quad_fast->grdUhAtQp(const_cast<const DimVec<WorldVector<T> >&>(Lambda), uh_loc, NULL);

    grduh_vec = traverseVector->getGrdAtQPs(elinfo, NULL, quad_fast, NULL);

    int dimOfWorld = Global::getGeo(WORLD);

    int numPoints = quad_fast->getNumPoints();
    for (normT = iq = 0; iq < numPoints; iq++)
      {
	for (norm2 = j = 0; j < dimOfWorld; j++)
	  norm2 += sqr(grduh_vec[iq][j]);

	normT += quad_fast->getWeight(iq)*norm2;
      }
    norm += det*normT;

    return 0;
  }


  template<typename T>
  double DOFVector<T>::H1NormSquare(Quadrature *q) const
  {
    FUNCNAME("DOFVector::H1NormSquare");
    int deg;
  
    Mesh *mesh = feSpace->getMesh();
    dim = mesh->getDim();

    if (!q)
      {
	deg = 2*feSpace->getBasisFcts()->getDegree()-2;
	q = Quadrature::provideQuadrature(dim, deg);
      }
    quad_fast = 
      FastQuadrature::provideFastQuadrature(feSpace->getBasisFcts(), 
					    *q, 
					    INIT_GRD_PHI);

    norm = 0.0;
    traverseVector = const_cast<DOFVector<T>*>(this);

    mesh->traverse(-1, Mesh::CALL_LEAF_EL | Mesh::FILL_COORDS | 
		   Mesh::FILL_DET | Mesh::FILL_GRD_LAMBDA,
		   H1Norm_fct);

    return norm;
  }

  template<typename T>
  void DOFVector<T>::compressDOFIndexed(int first, int last, 
					::std::vector<DegreeOfFreedom> &newDOF)
  {
    int i, j;
    for(i = first; i <= last; i++) {
      if((j = newDOF[i]) >= 0) {
	vec[j] = vec[i];
      }
    }
  }

  template<typename T>
  ElementVector* DOFVectorBase<T>::assemble(T factor, ElInfo *elInfo,
					    const BoundaryType *bound, 
					    Operator *op)
  {
    FUNCNAME("DOFVector::assemble");

753
754
    if (!(op || this->operators.size())) 
      return NULL;
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
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

    Operator *operat = op ? op : this->operators[0];

    this->elementVector = 
      operat->getAssembler()->initElementVector(this->elementVector, elInfo);

    if(op) {
      op->getElementVector(elInfo, this->elementVector);
    } else {
      ::std::vector<Operator*>::iterator it;
      ::std::vector<double*>::iterator factorIt;
      for(it = this->operators.begin(), factorIt = this->operatorFactor.begin();	
	  it != this->operators.end(); 
	  ++it, ++factorIt) 
	{
	  (*it)->getElementVector(elInfo, 
				  this->elementVector, 
				  *factorIt ? **factorIt : 1.0);
	}
    }

    addElementVector(factor,
		     *this->elementVector, 
		     bound);

    return this->elementVector;
  }

  template<typename T>
  Flag DOFVectorBase<T>::getAssembleFlag()
  {
    Flag fillFlag(0);
    ::std::vector<Operator*>::iterator op;
    for(op = this->operators.begin(); op != this->operators.end(); ++op) {
      fillFlag |= (*op)->getFillFlag();
    }
    return fillFlag;
  }

  template<typename T>
  DOFVector<T>& DOFVector<T>::operator=(const DOFVector<T>& rhs )
  {
    feSpace=rhs.feSpace;
    vec=rhs.vec;
    this->elementVector=NULL;
    interFct=rhs.interFct;
    refineInter=rhs.refineInter;
    coarsenOperation=rhs.coarsenOperation;
    this->operators=rhs.operators;
    this->operatorFactor=rhs.operatorFactor;
    if (rhs.boundaryManager) {
      if(this->boundaryManager) delete this->boundaryManager; 
      this->boundaryManager = new BoundaryManager(*rhs.boundaryManager);
      //    boundaryManager->setDOFVector(this);
    }
    else 
      this->boundaryManager=NULL;
    return *this;
  }

  template<typename T>
  const DOFVector<T>& operator*=(DOFVector<T>& x, T scal)
  {
818
    FUNCNAME("DOFVector<T>::operator*=(DOFVector<T>& x, T scal)");
819

820
821
    TEST_EXIT_DBG(x.getFESpace() && x.getFESpace()->getAdmin())
      ("pointer is NULL: %8X, %8X\n", x.getFESpace(), x.getFESpace()->getAdmin());
822
823

    typename DOFVector<T>::Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(&x), USED_DOFS);
824
    for (vecIterator.reset(); !vecIterator.end(); ++vecIterator) {
825
826
      (*vecIterator) *= scal; 
    };
827

828
829
830
831
832
833
834
    return x;
  }


  template<typename T>
  const DOFVector<T>& operator+=(DOFVector<T>& x, const DOFVector<T>& y)
  {
835
836
837
    FUNCNAME("DOFVector<T>::operator+=(DOFVector<T>& x, const DOFVector<T>& y)");
    
    TEST_EXIT_DBG(x.getFESpace() && y.getFESpace())
838
      ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace());
839
840
    TEST_EXIT_DBG(x.getFESpace()->getAdmin() &&
	      (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin()))
841
842
      ("no admin or different admins: %8X, %8X\n",
       x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin());
843
844
    TEST_EXIT_DBG(x.getSize() == y.getSize())("different sizes\n");
    
845
846
    typename DOFVector<T>::Iterator xIterator(dynamic_cast<DOFIndexed<T>*>(&x), USED_DOFS);
    typename DOFVector<T>::Iterator yIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(&y)), USED_DOFS);
847
848
849
850
851
852
    for (xIterator.reset(), yIterator.reset();
	 !xIterator.end();
	 ++xIterator, ++yIterator) {     
      *xIterator += *yIterator; 
    }

853
854
855
856
857
858
    return x;
  }

  template<typename T>
  const DOFVector<T>& operator-=(DOFVector<T>& x, const DOFVector<T>& y)
  {
859
    FUNCNAME("DOFVector<T>::operator-=(DOFVector<T>& x, const DOFVector<T>& y)");
860

861
    TEST_EXIT_DBG(x.getFESpace() && y.getFESpace())
862
      ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace());
863
864
    TEST_EXIT_DBG(x.getFESpace()->getAdmin() &&
	      (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin()))
865
866
      ("no admin or different admins: %8X, %8X\n",
       x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin());
867
868
    TEST_EXIT_DBG(x.getSize() == y.getSize())("different sizes\n");
    
869
870
    typename DOFVector<T>::Iterator xIterator(dynamic_cast<DOFIndexed<T>*>(&x), USED_DOFS);
    typename DOFVector<T>::Iterator yIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(&y)), USED_DOFS);
871
872
873
874
875
    for (xIterator.reset(), yIterator.reset();
	 !xIterator.end();
	 ++xIterator, ++yIterator) {
      *xIterator -= *yIterator; 
    }
876
877
878
879
880
881
    return x;
  }

  template<typename T>
  const DOFVector<T>& operator*=(DOFVector<T>& x, const DOFVector<T>& y)
  {
882
883
884
    FUNCNAME("DOFVector<T>::operator*=(DOFVector<T>& x, const DOFVector<T>& y)");
    
    TEST_EXIT_DBG(x.getFESpace() && y.getFESpace())
885
      ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace());
886
887
    TEST_EXIT_DBG(x.getFESpace()->getAdmin() &&
	      (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin()))
888
889
      ("no admin or different admins: %8X, %8X\n",
       x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin());
890
891
    TEST_EXIT_DBG(x.getSize() == y.getSize())("different sizes\n");
    
892
893
    typename DOFVector<T>::Iterator xIterator(dynamic_cast<DOFIndexed<T>*>(&x), USED_DOFS);
    typename DOFVector<T>::Iterator yIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(&y)), USED_DOFS);
894
895
896
897
898
    for (xIterator.reset(), yIterator.reset();
	 !xIterator.end();
	 ++xIterator, ++yIterator) {
      *xIterator *= *yIterator; 
    }
899
900
901
902
903
904
905
906
907

    return x;
  }

  template<typename T>
  T operator*(DOFVector<T>& x, DOFVector<T>& y)
  {
    FUNCNAME("DOFVector<T>::operator*(DOFVector<T>& x, DOFVector<T>& y)");

908
    TEST_EXIT_DBG(x.getFESpace() && y.getFESpace())
909
      ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace());
910
911
    TEST_EXIT_DBG(x.getFESpace()->getAdmin() && 
	      (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin()))
912
913
      ("no admin or different admins: %8X, %8X\n",
       x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin());
914
    TEST_EXIT_DBG(x.getSize() == y.getSize())("different sizes\n");
915

916
    T dot = 0;    
917
918
919
920
    typename DOFVector<T>::Iterator xIterator(dynamic_cast<DOFIndexed<T>*>(&x), USED_DOFS);
    typename DOFVector<T>::Iterator yIterator(dynamic_cast<DOFIndexed<T>*>(&y), USED_DOFS);
    for (xIterator.reset(), yIterator.reset();
	 !xIterator.end();
921
922
923
	 ++xIterator, ++yIterator) {
      dot += (*xIterator) * (*yIterator);
    };
924
925
926
927
928
929
930
931
932
933
934
935
936

    return(dot);
  }

  template<typename T>
  void mv(MatrixTranspose transpose,
	  const DOFMatrix &a, 
	  const DOFVector<T>&x,
	  DOFVector<T> &result,
	  bool add)
  {
    FUNCNAME("DOFVector<T>::mv()");

937
    TEST_EXIT_DBG(a.getRowFESpace() && a.getColFESpace() && x.getFESpace() && result.getFESpace())
938
939
      ("getFESpace() is NULL: %8X, %8X, %8X, %8X\n", 
       a.getRowFESpace(), a.getColFESpace(), x.getFESpace(), result.getFESpace());
940
941
942
943
944
945
946
947
    
    TEST_EXIT_DBG((a.getRowFESpace()->getAdmin() && a.getColFESpace()->getAdmin()) && 
	      (((transpose == NoTranspose) && 
		(a.getColFESpace()->getAdmin() == x.getFESpace()->getAdmin()) && 
		(a.getRowFESpace()->getAdmin() == result.getFESpace()->getAdmin())) ||
	       ((transpose == Transpose) && 
		(a.getRowFESpace()->getAdmin() == x.getFESpace()->getAdmin()) && 
		(a.getColFESpace()->getAdmin() == result.getFESpace()->getAdmin()))))
948
949
950
      ("no admin or different admins: %8X, %8X, %8X, %8X\n",
       a.getRowFESpace()->getAdmin(), a.getColFESpace()->getAdmin(), 
       x.getFESpace()->getAdmin(), result.getFESpace()->getAdmin());
951
    
952
    if (transpose == NoTranspose) {
953
      TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize())
954
	("x.size = %d too small: admin->sizeUsed = %d\n",
955
956
	 x.getSize(), a.getColFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>( result.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
957
	("size = %d too small: admin->sizeUsed = %d\n",
958
959
	 result.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
960
	("a.size = %d too small: admin->sizeUsed = %d\n",
961
962
	 a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
      
963

Thomas Witkowski's avatar
Thomas Witkowski committed
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
      typename DOFVector<T>::Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(&result), USED_DOFS); 
      DOFMatrix::Iterator rowIterator(const_cast<DOFMatrix*>(&a), USED_DOFS);
      for(vecIterator.reset(), rowIterator.reset();
      	  !rowIterator.end();
      	  ++rowIterator, ++vecIterator) { 
      	
	double sum = 0;
      	if (!add) 
	  *vecIterator = 0.0;
      	for(::std::vector<MatEntry>::iterator colIterator = rowIterator->begin();
      	    colIterator != rowIterator->end();
      	    colIterator++) {
      	  
	  int jcol = colIterator->col;
      	  if (jcol >= 0) { // entry used? 
      	    sum += (static_cast<double>(colIterator->entry)) * x[jcol];
      	  } else {
      	    if (jcol == DOFMatrix::NO_MORE_ENTRIES)
      	      break;
      	  }
      	}
      	*vecIterator += sum;
986
987
988
      }

    } else if (transpose == Transpose) {
989
      TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
990
	("x.size = %d too small: admin->sizeUsed = %d\n",
991
992
	 x.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>( result.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize())
993
	("size = %d too small: admin->sizeUsed = %d\n",
994
995
	 result.getSize(), a.getColFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
996
	("a.size = %d too small: admin->sizeUsed = %d\n",
997
998
	 a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
      
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
      if(!add) {
	typename DOFVector<T>::Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(&result), USED_DOFS);
	for(vecIterator.reset();
	    !vecIterator.end(); ++vecIterator) {
	  *vecIterator = 0; 
	};
      }

      typename DOFVector<T>::Iterator xIterator(const_cast<DOFIndexed<T>*>(dynamic_cast<const DOFIndexed<T>*>(&x)), USED_DOFS);
      DOFMatrix::Iterator rowIterator(const_cast<DOFMatrix*>(&a), USED_DOFS);
      for(xIterator.reset(), rowIterator.reset();
	  !rowIterator.end();
	  ++rowIterator, ++xIterator) { 
	T ax = (*xIterator);
	for(::std::vector<MatEntry>::iterator colIterator = rowIterator->begin();
	    colIterator != rowIterator->end();
	    colIterator++) {
	  int jcol = colIterator->col;
	  if (jcol >= 0) // entry used?
	    result[jcol] += (static_cast<double>(colIterator->entry)) * ax;
	  else 
	    if (jcol == DOFMatrix::NO_MORE_ENTRIES) break;
	}
      }
    } else {
      ERROR_EXIT("transpose = %d\n", transpose);
    }
  }

  template<typename T>
  void axpy(double alpha, 
	    const DOFVector<T>& x, 
	    DOFVector<T>& y)
  {
    FUNCNAME("DOFVector<T>::axpy()");

1035
    TEST_EXIT_DBG(x.getFESpace() && y.getFESpace())
1036
1037
      ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace());

1038
1039
    TEST_EXIT_DBG(x.getFESpace()->getAdmin() && 
	      (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin()))
1040
1041
      ("no admin or different admins: %8X, %8X\n",
       x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin());
1042
    TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize())
1043
      ("size = %d too small: admin->size = %d\n", x.getSize(), 
1044
1045
       x.getFESpace()->getAdmin()->getUsedSize());
    TEST_EXIT_DBG(static_cast<int>(y.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize())
1046
      ("y.size = %d too small: admin->size = %d\n", y.getSize(), 
1047
1048
       x.getFESpace()->getAdmin()->getUsedSize());
    
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
      // This is the old implementation of the mv-multiplication. It have been changed
      // because of the OpenMP-parallelization:
//     typename DOFVector<T>::Iterator xIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(&x)), USED_DOFS);
//     typename DOFVector<T>::Iterator yIterator(dynamic_cast<DOFIndexed<T>*>(&y), USED_DOFS);
//     for(xIterator.reset(), yIterator.reset();
// 	!xIterator.end();
// 	++xIterator, ++yIterator) 
//       {  
// 	*yIterator += alpha * (*xIterator); 
//       };

      int i;
      int maxI = y.getSize();

#ifdef _OPENMP
#pragma omp parallel for schedule(dynamic, 25000) default(shared) private(i)
#endif
      for (i = 0; i < maxI; i++) {
1067
	if (!x.getFESpace()->getAdmin()->isDOFFree(i)) {
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
	  y[i] += alpha * x[i];
	}
      }
  }

  template<typename T>
  const DOFVector<T>& operator*(const DOFVector<T>& v, double d)
  {
    static DOFVector<T> result;
    return mult(d, v, result);
  }

  template<typename T>
  const DOFVector<T>& operator*(double d, const DOFVector<T>& v)
  {
    static DOFVector<T> result;
    return mult(d, v, result);
  }

  template<typename T>
  const DOFVector<T>& operator+(const DOFVector<T> &v1 , const DOFVector<T> &v2)
  {
    static DOFVector<T> result;
    return add(v1, v2, result);
  }


  template<typename T>
  void xpay(double alpha, 
	    const DOFVector<T>& x, 
	    DOFVector<T>& y)
  {
    FUNCNAME("DOFVector<T>::xpay()");

1102
    TEST_EXIT_DBG(x.getFESpace() && y.getFESpace())
1103
1104
      ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace());

1105
1106
    TEST_EXIT_DBG(x.getFESpace()->getAdmin() && 
	      (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin()))
1107
1108
      ("no admin or different admins: %8X, %8X\n",
       x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin());
1109
    TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize())
1110
      ("size = %d too small: admin->size = %d\n", x.getSize(), 
1111
1112
       x.getFESpace()->getAdmin()->getUsedSize());
    TEST_EXIT_DBG(static_cast<int>(y.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize())
1113
      ("y.size = %d too small: admin->size = %d\n", y.getSize(), 
1114
1115
       x.getFESpace()->getAdmin()->getUsedSize());
    
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
    // This is the old implementation of the mv-multiplication. It have been changed
    // because of the OpenMP-parallelization:
    //     typename DOFVector<T>::Iterator xIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(&x)), USED_DOFS);
    //     typename DOFVector<T>::Iterator yIterator(dynamic_cast<DOFIndexed<T>*>(&y), USED_DOFS);
    //     for(xIterator.reset(), yIterator.reset();
    // 	!xIterator.end();
    // 	++xIterator, ++yIterator) 
    //       {  
    // 	*yIterator = alpha *(*yIterator)+ (*xIterator); 
    //       };

      int i;
      int maxI = y.getSize();

1130
1131
      const DOFAdmin* admin = x.getFESpace()->getAdmin();

1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
#ifdef _OPENMP
#pragma omp parallel for schedule(dynamic, 25000) default(shared) private(i)
#endif
      for (i = 0; i < maxI; i++) {
	if (!admin->isDOFFree(i)) {
	  y[i] = alpha * y[i] + x[i];
	}
      }
  }

  template<typename T>
  inline const DOFVector<T>& mult(double scal,
				  const DOFVector<T>& v,
				  DOFVector<T>& result)
  {
    typename DOFVector<T>::Iterator vIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(&v)), USED_DOFS);
    typename DOFVector<T>::Iterator rIterator(dynamic_cast<DOFIndexed<T>*>(&result), USED_DOFS);
    for(vIterator.reset(), rIterator.reset();
	!vIterator.end();
	++vIterator, ++rIterator) 
      {  
	*rIterator = scal * (*vIterator);
      };

    return result;
  }

  template<typename T>
  inline const DOFVector<T>& add(const DOFVector<T>& v,
				 double scal,
				 DOFVector<T>& result)
  {
    typename DOFVector<T>::Iterator vIterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(&v)), USED_DOFS);
    typename DOFVector<T>::Iterator rIterator(dynamic_cast<DOFIndexed<T>*>(&result), USED_DOFS);
    for(vIterator.reset(), rIterator.reset();
	!vIterator.end();
	++vIterator, ++rIterator) 
      {  
	*rIterator = (*vIterator) + scal;
      };
    return result;
  }

  template<typename T>
  inline const DOFVector<T>& add(const DOFVector<T>& v1,
				 const DOFVector<T>& v2,
				 DOFVector<T>& result)
  {
    typename DOFVector<T>::Iterator v1Iterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(&v1)), USED_DOFS);
    typename DOFVector<T>::Iterator v2Iterator(dynamic_cast<DOFIndexed<T>*>(const_cast<DOFVector<T>*>(&v2)), USED_DOFS);
    typename DOFVector<T>::Iterator rIterator(dynamic_cast<DOFIndexed<T>*>(&result), USED_DOFS);
    for(v1Iterator.reset(), v2Iterator.reset(), rIterator.reset();
	!v1Iterator.end();
	++v1Iterator, ++v2Iterator, ++rIterator) 
      {  
	*rIterator = (*v1Iterator) + (*v2Iterator);
      };
    return result;

  }

  template<typename T>
  const T *DOFVectorBase<T>::getLocalVector(const Element *el, T *d) const
  {
    static T* localVec = NULL;
    const DOFAdmin* admin = feSpace->getAdmin();

    int i;
    int nBasFcts = feSpace->getBasisFcts()->getNumber();

    T *result;
  
    if(d) {
      result = d;
    } else {
      if(localVec) delete [] localVec;
      localVec = new T[nBasFcts]; 
      result = localVec;
    }

    const DegreeOfFreedom *localIndices = 
      feSpace->getBasisFcts()->getLocalIndices(el, admin, NULL);

    for(i = 0; i < nBasFcts; i++) {
      result[i] = (*this)[localIndices[i]];
    }

    return result;
  }

  template<typename T>
  const T *DOFVectorBase<T>::getVecAtQPs(const ElInfo         *elInfo, 
					 const Quadrature     *quad,
					 const FastQuadrature *quadFast,
					 T                    *vecAtQPs) const
  {
    FUNCNAME("DOFVector<T>::getVecAtQPs()");
  
1230
    TEST_EXIT_DBG(quad || quadFast)("neither quad nor quadFast defined\n");
1231

1232
1233
1234
    if (quad && quadFast) {
      TEST_EXIT_DBG(quad == quadFast->getQuadrature())
      	("quad != quadFast->quadrature\n");
1235
1236
    }

1237
    TEST_EXIT_DBG(!quadFast || quadFast->getBasisFunctions() == feSpace->getBasisFcts())
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
      ("invalid basis functions");

    Element *el = elInfo->getElement();

    const Quadrature *quadrature = quadFast ? quadFast->getQuadrature() : quad;

    const BasisFunction *basFcts = feSpace->getBasisFcts();

    int numPoints = quadrature->getNumPoints();
    int nBasFcts  = basFcts->getNumber();
    int i, j;

    static T *localvec = NULL;

    T *result;

    if (vecAtQPs) {
      result = vecAtQPs;
    } else {
      if(localvec) delete [] localvec;
      localvec = new T[numPoints];
      for(i = 0; i < numPoints; i++) {
	localvec[i] = 0.0;
      }
      result = localvec;
    }
  
    const T *localVec = getLocalVector(el, NULL);

    for(i = 0; i < numPoints; i++) {
      for(result[i] = j = 0; j < nBasFcts; j++) {
	result[i] +=  localVec[j] * 
	  (quadFast ? 
	   (quadFast->getPhi(i, j)) : 
	   ((*(basFcts->getPhi(j)))(quad->getLambda(i))));
      }
    }

    return const_cast<const T*>(result);
  }

}