DOFVector.hh 37.9 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
#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"
16
#include "OpenMP.h"
17
18
19
20

namespace AMDiS {

  template<typename T>
21
  DOFVectorBase<T>::DOFVectorBase(const FiniteElemSpace *f, std::string n)
Thomas Witkowski's avatar
Thomas Witkowski committed
22
23
24
25
    : feSpace(f),
      name(n),
      elementVector(NULL),
      boundaryManager(NULL)
26
  {
Thomas Witkowski's avatar
Thomas Witkowski committed
27
    nBasFcts = feSpace->getBasisFcts()->getNumber();
28
    int dim = feSpace->getMesh()->getDim();
29

Thomas Witkowski's avatar
Thomas Witkowski committed
30
    localIndices.resize(omp_get_max_threads());
31
32
    grdPhis.resize(omp_get_max_threads());
    D2Phis.resize(omp_get_max_threads());
33

Thomas Witkowski's avatar
Thomas Witkowski committed
34
35
    for (int i = 0; i < omp_get_max_threads(); i++) {
      localIndices[i] = GET_MEMORY(DegreeOfFreedom, this->nBasFcts);      
36
37
      grdPhis[i] = NEW DimVec<double>(dim, DEFAULT_VALUE, 0.0);
      D2Phis[i] = NEW DimMat<double>(dim, NO_INIT);
38
39
    }
  }
Thomas Witkowski's avatar
Thomas Witkowski committed
40
  
41
  template<typename T>
Thomas Witkowski's avatar
Thomas Witkowski committed
42
43
44
45
  DOFVectorBase<T>::~DOFVectorBase()
  {
    for (int i = 0; i < static_cast<int>(localIndices.size()); i++) {
      FREE_MEMORY(localIndices[i], DegreeOfFreedom, this->nBasFcts);
46
47
      DELETE grdPhis[i];
      DELETE D2Phis[i];
Thomas Witkowski's avatar
Thomas Witkowski committed
48
49
50
51
    }
  }
  
  template<typename T>
52
  DOFVector<T>::DOFVector(const FiniteElemSpace* f, std::string n)
53
    : DOFVectorBase<T>(f, n),
54
55
      refineInter(false), 
      coarsenOperation(NO_OPERATION)
56
57
58
59
60
  {
    init(f, n);
  } 

  template<typename T>
61
  void DOFVector<T>::init(const FiniteElemSpace* f, std::string n)
62
63
64
  {
    name = n;
    feSpace = f;
65
    if (feSpace && feSpace->getAdmin()) {
66
      (feSpace->getAdmin())->addDOFIndexed(this);
Thomas Witkowski's avatar
Thomas Witkowski committed
67
68
    }
      
69
70
71
72
73
74
    this->boundaryManager = NEW BoundaryManager;
  }

  template<typename T>
  DOFVector<T>::~DOFVector()
  {
75
    if (feSpace && feSpace->getAdmin()) {
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
      (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;

Thomas Witkowski's avatar
Thomas Witkowski committed
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
  template<typename T>
  void DOFVectorBase<T>::addElementVector(T factor, 
					  const ElementVector &elVec, 
					  const BoundaryType *bound,
					  bool add)
  {
    FUNCNAME("DOFVector::addElementVector()");

    int n_row = elVec.getSize();

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

      if(!(condition && condition->isDirichlet())) {
	DegreeOfFreedom irow = elVec.dofIndices[i];
	(*this)[irow] = (add ? (*this)[irow] : 0.0);
	(*this)[irow] += factor * elVec[i];
      }
    }
  }

118
119
120
121
122
  template<typename T>
  double DOFVector<T>::nrm2() const
  {
    FUNCNAME("DOFVector<T>::nrm2()");

123
124
125
    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())
126
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(),
127
128
129
       feSpace->getAdmin()->getUsedSize());
    
    double nrm = 0.0;
130
131
132
133
134
135
136
137
138
139
140
141
    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()");

142
143
144
    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())
145
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(),
146
147
148
       feSpace->getAdmin()->getUsedSize());
    
    double nrm = 0.0;
149
150
151
152
153
154
155
156
157
158
    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
  {
159
    FUNCNAME("DOFVector<T>::asum()");
160

161
162
163
164
165
    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());
166

167
    double nrm = 0.0;
168
169
170
171
172
173
174
175
176
177
    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
  {
178
    FUNCNAME("DOFVector<T>::sum()");
179

180
181
182
183
184
    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());
185

186
    double nrm = 0.0;    
187
188
189
190
191
192
193
194
195
196
    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)
  {
197
    FUNCNAME("DOFVector<T>::set()");
198

199
200
201
    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())
202
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(),
203
204
       feSpace->getAdmin()->getUsedSize());
    
205
206
207
208
209
210
211
212
213
214
    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)
  {
215
216
217
218
219
    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()))
220
221
      ("no admin or different admins: %8X, %8X\n",
       feSpace->getAdmin(), x.feSpace->getAdmin());
222
    TEST_EXIT_DBG(static_cast<int>(vec.size()) >= feSpace->getAdmin()->getUsedSize())
223
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), 
224
225
       feSpace->getAdmin()->getUsedSize());
    TEST_EXIT_DBG(static_cast<int>(x.vec.size()) >= feSpace->getAdmin()->getUsedSize())
226
      ("x.size = %d too small: admin->sizeUsed = %d\n", x.vec.size(), 
227
228
       feSpace->getAdmin()->getUsedSize());
    
229
230
    Iterator vecIterator(dynamic_cast<DOFIndexed<T>*>(this), USED_DOFS);
    Iterator xIterator(dynamic_cast<DOFVector<T>*>(const_cast<DOFVector<T>*>(&x)), USED_DOFS);
231
232
233
234
235
236
    for (vecIterator.reset(), xIterator.reset();
	 !vecIterator.end();
	 ++vecIterator, ++xIterator) {
      
      *vecIterator = *xIterator; 
    };
237
238
239
240
241
242
  }


  template<typename T>
  T DOFVector<T>::min() const
  {
243
244
245
246
247
    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())
248
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), 
249
250
251
       feSpace->getAdmin()->getUsedSize());

    T m;    
252
    Iterator vecIterator(const_cast<DOFIndexed<T>*>(dynamic_cast<const DOFIndexed<T>*>(this)), USED_DOFS);
253
    for (vecIterator.reset(), m = *vecIterator; !vecIterator.end(); ++vecIterator) {
254
      m = std::min(m, *vecIterator);
255
256
    }

257
258
259
260
261
262
    return m;
  }

  template<typename T>
  T DOFVector<T>::max() const 
  {
263
264
265
266
267
    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())
268
      ("size = %d too small: admin->sizeUsed = %d\n", vec.size(), 
269
       feSpace->getAdmin()->getUsedSize());
270

271
    T m;    
272
    Iterator vecIterator(const_cast<DOFIndexed<T>*>(dynamic_cast<const DOFIndexed<T>*>(this)), USED_DOFS);
273
    for (vecIterator.reset(), m = *vecIterator; !vecIterator.end(); ++vecIterator) {
274
      m = std::max(m, *vecIterator);
275
276
    }

277
278
279
280
281
282
283
284
    return m;
  }

  template<typename T>
  void gemv(MatrixTranspose transpose, T alpha,
	    const DOFMatrix& a, const DOFVector<T>& x,
	    T beta, DOFVector<T>& y)
  {
285
286
287
288
    FUNCNAME("DOFVector<T>::gemv()");

    int j, jcol, ysize;
    T sum, ax;
289
    const DOFMatrix::MatrixRow *row;
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
    
    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()))))
305
306
307
      ("no admin or different admins: %8X, %8X, %8X, %8X\n",
       a.getRowFESpace()->getAdmin(), a.getColFESpace()->getAdmin(), 
       x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin());
308
    
309
    if (transpose == NoTranspose) {
310
      TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize())
311
	("x.size = %d too small: admin->sizeUsed = %d\n",
312
313
	 x.getSize(), a.getColFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>(y.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
314
	("y.size = %d too small: admin->sizeUsed = %d\n",
315
316
	 y.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
317
	("a.size = %d too small: admin->sizeUsed = %d\n",
318
	a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
319
320
    }
    else  if (transpose == Transpose) {
321
      TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
322
	("x.size = %d too small: admin->sizeUsed = %d\n",
323
324
	 x.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>(y.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize())
325
	("y.size = %d too small: admin->sizeUsed = %d\n",
326
327
	 y.getSize(), a.getColFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
328
	("a.size = %d too small: admin->sizeUsed = %d\n",
329
	a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
    }

    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()]);
348
	for(std::vector<MatEntry>::iterator colIterator = rowIterator->begin();
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
	    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()]);
377
	for(std::vector<MatEntry>::iterator colIterator = rowIterator->begin();
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
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
	    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)
  {
464
    FUNCNAME("DOFVector<T>::interpol()");
465

466
467
468
    TEST_EXIT_DBG(fct)("no function to interpolate\n");

    interFct = fct;
469

470
471
472
473
474
    if (!this->getFESpace()) {
      MSG("no dof admin in vec %s, skipping interpolation\n",
	  this->getName().c_str());
      return;
    }
475

476
477
478
479
480
    if (!(this->getFESpace()->getAdmin())) {
      MSG("no dof admin in feSpace %s, skipping interpolation\n", 
	  this->getFESpace()->getName().c_str());
      return;
    }
481
  
482
483
484
485
486
    if (!(this->getFESpace()->getBasisFcts())) {
      MSG("no basis functions in admin of vec %s, skipping interpolation\n",
	  this->getName().c_str());
      return;
    }
487

488
489
490
491
492
    if (!(fct)) {
      MSG("function that should be interpolated only pointer to NULL, ");
      Msg::print("skipping interpolation\n");
      return;
    }
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

    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();
566
    uh_vec = traverseVector->getVecAtQPs(elinfo, NULL, quad_fast, NULL);
567
568

    int numPoints = quad_fast->getNumPoints();
569
570
571
    for (normT = iq = 0; iq < numPoints; iq++) {
      normT += quad_fast->getWeight(iq)*(uh_vec[iq]);
    }
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
    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();
617
    uh_vec = traverseVector->getVecAtQPs(elinfo, NULL, quad_fast, NULL);
618
619

    int numPoints = quad_fast->getNumPoints();
620
621
622
    for (normT = iq = 0; iq < numPoints; iq++) {
      normT += quad_fast->getWeight(iq)*abs(uh_vec[iq]);
    }
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
    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();
667
    uh_vec = traverseVector->getVecAtQPs(elinfo, NULL, quad_fast, NULL);
668
669

    int numPoints = quad_fast->getNumPoints();
670
671
672
    for (normT = iq = 0; iq < numPoints; iq++) {
      normT += quad_fast->getWeight(iq)*sqr(uh_vec[iq]);
    }
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
    norm += det*normT;

    return 0;
  }

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

    double det = elinfo->getDet();

    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, 
736
					std::vector<DegreeOfFreedom> &newDOF)
737
738
739
740
741
742
743
744
745
746
747
748
749
750
  {
    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)
  {
751
    FUNCNAME("DOFVector::assemble()");
752

753
754
    if (!(op || this->operators.size())) 
      return NULL;
755
756
757
758

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

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

761
    if (op) {
762
763
      op->getElementVector(elInfo, this->elementVector);
    } else {
764
765
      std::vector<Operator*>::iterator it;
      std::vector<double*>::iterator factorIt;
766
767
768
769
770
771
772
      for (it = this->operators.begin(), factorIt = this->operatorFactor.begin();	
	   it != this->operators.end(); 
	   ++it, ++factorIt) {
	(*it)->getElementVector(elInfo, 
				this->elementVector, 
				*factorIt ? **factorIt : 1.0);
      }
773
774
775
    }

    addElementVector(factor,
776
777
    		     *this->elementVector, 
    		     bound);
778
779
780
781
782
783
784
785

    return this->elementVector;
  }

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

  template<typename T>
794
795
796
797
798
799
800
801
802
803
  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;
804
    if (rhs.boundaryManager) {
805
      if(this->boundaryManager) delete this->boundaryManager; 
806
807
808
      this->boundaryManager = new BoundaryManager(*rhs.boundaryManager);
      //    boundaryManager->setDOFVector(this);
    }
809
810
    else 
      this->boundaryManager=NULL;
811
812
813
814
815
816
    return *this;
  }

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

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

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

827
828
829
830
831
832
833
    return x;
  }


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

852
853
854
855
856
857
    return x;
  }

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

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

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

    return x;
  }

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

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

915
    T dot = 0;    
916
917
918
919
    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();
920
921
922
	 ++xIterator, ++yIterator) {
      dot += (*xIterator) * (*yIterator);
    };
923
924
925
926
927
928
929
930
931
932
933
934
935

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

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

Thomas Witkowski's avatar
Thomas Witkowski committed
963
964
965
966
967
968
969
970
971
      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;
972
      	for(std::vector<MatEntry>::iterator colIterator = rowIterator->begin();
Thomas Witkowski's avatar
Thomas Witkowski committed
973
974
975
976
977
978
979
980
981
982
983
984
      	    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;
985
986
987
      }

    } else if (transpose == Transpose) {
988
      TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
989
	("x.size = %d too small: admin->sizeUsed = %d\n",
990
991
	 x.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>( result.getSize()) >= a.getColFESpace()->getAdmin()->getUsedSize())
992
	("size = %d too small: admin->sizeUsed = %d\n",
993
994
	 result.getSize(), a.getColFESpace()->getAdmin()->getUsedSize());
      TEST_EXIT_DBG(static_cast<int>( a.getSize()) >= a.getRowFESpace()->getAdmin()->getUsedSize())
995
	("a.size = %d too small: admin->sizeUsed = %d\n",
996
997
	 a.getSize(), a.getRowFESpace()->getAdmin()->getUsedSize());
      
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
      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);
1012
	for (std::vector<MatEntry>::iterator colIterator = rowIterator->begin();
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
	    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()");

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

1037
1038
    TEST_EXIT_DBG(x.getFESpace()->getAdmin() && 
	      (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin()))
1039
1040
      ("no admin or different admins: %8X, %8X\n",
       x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin());
1041
    TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize())
1042
      ("size = %d too small: admin->size = %d\n", x.getSize(), 
1043
1044
       x.getFESpace()->getAdmin()->getUsedSize());
    TEST_EXIT_DBG(static_cast<int>(y.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize())
1045
      ("y.size = %d too small: admin->size = %d\n", y.getSize(), 
1046
1047
       x.getFESpace()->getAdmin()->getUsedSize());
    
1048
1049
1050
1051
1052
1053
1054
    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); 
    };
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
  }

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

1086
    TEST_EXIT_DBG(x.getFESpace() && y.getFESpace())
1087
1088
      ("feSpace is NULL: %8X, %8X\n", x.getFESpace(), y.getFESpace());

1089
1090
    TEST_EXIT_DBG(x.getFESpace()->getAdmin() && 
	      (x.getFESpace()->getAdmin() == y.getFESpace()->getAdmin()))
1091
1092
      ("no admin or different admins: %8X, %8X\n",
       x.getFESpace()->getAdmin(), y.getFESpace()->getAdmin());
1093
    TEST_EXIT_DBG(static_cast<int>(x.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize())
1094
      ("size = %d too small: admin->size = %d\n", x.getSize(), 
1095
1096
       x.getFESpace()->getAdmin()->getUsedSize());
    TEST_EXIT_DBG(static_cast<int>(y.getSize()) >= x.getFESpace()->getAdmin()->getUsedSize())
1097
      ("y.size = %d too small: admin->size = %d\n", y.getSize(), 
1098
1099
       x.getFESpace()->getAdmin()->getUsedSize());
    
1100
1101
1102
1103
1104
1105
1106
    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); 
    };
1107
1108
1109
1110
1111
1112
1113
1114
1115
  }

  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);
1116
1117
1118
1119
1120
    for (vIterator.reset(), rIterator.reset();
	 !vIterator.end();
	 ++vIterator, ++rIterator) {
      *rIterator = scal * (*vIterator);
    };
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
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

    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;
Thomas Witkowski's avatar
Thomas Witkowski committed
1163
    T *result;   
1164
    if (d) {
1165
1166
      result = d;
    } else {
1167
1168
      if (localVec) 
	delete [] localVec;
1169
      localVec = new T[nBasFcts]; 
1170
      result = localVec;     
1171
    }
1172
    
Thomas Witkowski's avatar
Thomas Witkowski committed
1173
1174
1175
    DegreeOfFreedom *myLocalIndices = localIndices[omp_get_thread_num()];
    feSpace->getBasisFcts()->getLocalIndices(el, feSpace->getAdmin(), 
					     myLocalIndices);
1176

1177
    for (int i = 0; i < nBasFcts; i++) {
Thomas Witkowski's avatar
Thomas Witkowski committed
1178
      result[i] = (*this)[myLocalIndices[i]];
1179
    }
Thomas Witkowski's avatar
Thomas Witkowski committed
1180
       
1181
1182
1183
1184
    return result;
  }

  template<typename T>
1185
1186
  const T *DOFVectorBase<T>::getVecAtQPs(const ElInfo *elInfo, 
					 const Quadrature *quad,
1187
					 const FastQuadrature *quadFast,
1188
					 T *vecAtQPs) const
1189
1190
1191
  {
    FUNCNAME("DOFVector<T>::getVecAtQPs()");
  
1192
    TEST_EXIT_DBG(quad || quadFast)("neither quad nor quadFast defined\n");
1193

1194
1195
1196
    if (quad && quadFast) {
      TEST_EXIT_DBG(quad == quadFast->getQuadrature())
      	("quad != quadFast->quadrature\n");
1197
1198
    }

1199
    TEST_EXIT_DBG(!quadFast || quadFast->getBasisFunctions() == feSpace->getBasisFcts())
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
      ("invalid basis functions");

    Element *el = elInfo->getElement();

    const Quadrature *quadrature = quadFast ? quadFast->getQuadrature() : quad;
    const BasisFunction *basFcts = feSpace->getBasisFcts();

    int numPoints = quadrature->getNumPoints();

    static T *localvec = NULL;

    T *result;

    if (vecAtQPs) {
      result = vecAtQPs;
    } else {
1216
1217
      if (localvec) 
	delete [] localvec;
1218
      localvec = new T[numPoints];
1219
      for (int i = 0; i < numPoints; i++) {
1220
1221
1222
1223
	localvec[i] = 0.0;
      }
      result = localvec;
    }
1224
1225
1226
1227
1228
1229
1230
1231
      
    T *localVec = new T[nBasFcts];
    getLocalVector(el, localVec);
    
    for (int i = 0; i < numPoints; i++) {
      result[i] = 0.0;
      for (int j = 0; j < nBasFcts; j++) {
	result[i] += localVec[j] * 
1232
1233
1234
1235
1236
	  (quadFast ? 
	   (quadFast->getPhi(i, j)) : 
	   ((*(basFcts->getPhi(j)))(quad->getLambda(i))));
      }
    }
1237
1238
1239
    
    delete [] localVec;
    
1240
1241
1242
1243
    return const_cast<const T*>(result);
  }

}