ProblemVec.cc 29.6 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
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
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
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
#include "ProblemVec.h"
#include "RecoveryEstimator.h"
#include "Serializer.h"
#include "AbstractFunction.h"
#include "Operator.h"
#include "SystemVector.h"
#include "DOFMatrix.h"
#include "FiniteElemSpace.h"
#include "Estimator.h"
#include "Marker.h"
#include "AdaptInfo.h"
#include "FileWriter.h"
#include "CoarseningManager.h"
#include "RefinementManager.h"
#include "Mesh.h"
#include "OEMSolver.h"
#include "Preconditioner.h"
#include "MatVecMultiplier.h"
#include "DirichletBC.h"
#include "RobinBC.h"
#include "PeriodicBC.h"
#include "Lagrange.h"

namespace AMDiS {

  ProblemVec *ProblemVec::traversePtr_ = NULL;

  void ProblemVec::initialize(Flag initFlag,
			      ProblemVec *adoptProblem,
			      Flag adoptFlag)
  {
    FUNCNAME("ProblemVec::initialize()");
    
    // === create meshes ===
    if (meshes_.size() != 0) { 
      WARNING("meshes already created\n");
    } else {
      if (initFlag.isSet(CREATE_MESH) || 
	  ((!adoptFlag.isSet(INIT_MESH))&&
	   (initFlag.isSet(INIT_SYSTEM) || initFlag.isSet(INIT_FE_SPACE)))) {
	createMesh();
      } 
      if (adoptProblem && 
	  (adoptFlag.isSet(INIT_MESH) || 
	   adoptFlag.isSet(INIT_SYSTEM) ||
	   adoptFlag.isSet(INIT_FE_SPACE))) {
	meshes_ = adoptProblem->getMeshes();
	componentMeshes_ = adoptProblem->componentMeshes_;
	refinementManager_ = adoptProblem->refinementManager_;
	coarseningManager_ = adoptProblem->coarseningManager_;
      }
    }

    if (meshes_.size() == 0) 
      WARNING("no mesh created\n");

    // === create fespace ===
    if (feSpaces_.size() != 0) {
      WARNING("feSpaces already created\n");
    } else {
      if (initFlag.isSet(INIT_FE_SPACE) || 
	  (initFlag.isSet(INIT_SYSTEM)&&!adoptFlag.isSet(INIT_FE_SPACE))) {
	createFESpace();
      } 
      if (adoptProblem &&
	  (adoptFlag.isSet(INIT_FE_SPACE) || adoptFlag.isSet(INIT_SYSTEM))) {
	feSpaces_ = adoptProblem->getFESpaces();
	componentSpaces_ = adoptProblem->componentSpaces_;
      }
    }

    if (feSpaces_.size() == 0) 
      WARNING("no feSpace created\n");

    // === create system ===
    if (initFlag.isSet(INIT_SYSTEM)) {
      createMatricesAndVectors();
    } 
    if (adoptProblem && adoptFlag.isSet(INIT_SYSTEM)) {
      solution_ = adoptProblem->getSolution();
      rhs_ = adoptProblem->getRHS();
      systemMatrix_ = adoptProblem->getSystemMatrix();
    }

    // === create solver ===
    if (solver_) {
      WARNING("solver already created\n");
    } else {
      if (initFlag.isSet(INIT_SOLVER)) {
	createSolver();
      } 
      if (adoptProblem && adoptFlag.isSet(INIT_SOLVER)) {
	TEST_EXIT(!solver_)("solver already created\n");
	solver_ = adoptProblem->getSolver();
      }
    }

    if (!solver_) 
      WARNING("no solver created\n");

    // === create estimator ===
    if (initFlag.isSet(INIT_ESTIMATOR)) {
      createEstimator();
    } 
    if (adoptProblem && adoptFlag.isSet(INIT_ESTIMATOR)) {
      estimator_ = adoptProblem->getEstimator();
    } 

    // === create marker ===
    if (initFlag.isSet(INIT_MARKER)) {
      createMarker();
    } 
    if (adoptProblem && adoptFlag.isSet(INIT_MARKER)) {
      marker_ = adoptProblem->getMarker();
    } 


    // === create file writer ===
    if (initFlag.isSet(INIT_FILEWRITER)) {
      createFileWriter();
    }

    
    // === read serialization and init mesh ===
    
    // There are two possiblities where the user can define a serialization
    // to be read from disk. Either by providing the parameter -rs when executing
    // the program or in the init file. The -rs parameter is always checked first,
    // because it can be added automatically when  rescheduling the program
    // before timeout of the runqueue.

    int readSerialization = 0;
    ::std::string serializationFilename = "";
    GET_PARAMETER(0, "argv->rs", &serializationFilename);

    // If the parameter -rs is set, we do nothing here, because the problem will be
    // deserialized in the constructor of a following AdaptInstationary initialization.
    if (!serializationFilename.compare("")) {
      int readSerializationWithAdaptInfo = 0;

      GET_PARAMETER(0, name_ + "->input->read serialization", "%d", 
		    &readSerialization);
      GET_PARAMETER(0, name_ + "->input->serialization with adaptinfo", "%d",
		    &readSerializationWithAdaptInfo);

      // The serialization file is only read, if the adaptInfo part should not be used.
      // If the adaptInfo part should be also read, the serialization file will be read
      // in the constructor of the AdaptInstationary problem, because we do not have here
      // the adaptInfo object.
      if (readSerialization && !readSerializationWithAdaptInfo) {
	GET_PARAMETER(0, name_ + "->input->serialization filename", 
		      &serializationFilename);
	TEST_EXIT(serializationFilename != "")("no serialization file\n");

	MSG("Deserialization from file: %s\n", serializationFilename.c_str());
	::std::ifstream in(serializationFilename.c_str());
	deserialize(in);
	in.close();
      } else {
	// Initialize the meshes if there is no serialization file.
	for (int i = 0; i < static_cast<int>(meshes_.size()); i++) {
	  if (initFlag.isSet(INIT_MESH) && 
	      meshes_[i] && 
	      !(meshes_[i]->isInitialized())) {
	    meshes_[i]->initialize();
	    
	    // do global refinements
	    int globalRefinements = 0;
	    GET_PARAMETER(0, meshes_[0]->getName() + "->global refinements", "%d", 
			  &globalRefinements);
	    for (int gr = 0; gr < static_cast<int>(meshes_.size()); gr++) {
	      refinementManager_->globalRefine(meshes_[gr], globalRefinements);
	    }
	  }
	}	
      }
    }

    doOtherStuff();
  }

  void ProblemVec::createMesh() 
  {
    FUNCNAME("ProblemVec::createMesh()");

    componentMeshes_.resize(numComponents_);
    ::std::map<int, Mesh*> meshForRefinementSet;
    char number[3];

    ::std::string meshName("");
    GET_PARAMETER(0, name_ + "->mesh", &meshName);
    TEST_EXIT(meshName != "")("no mesh name spezified\n");
    int dim = 0;
    GET_PARAMETER(0, name_ + "->dim", "%d", &dim);
    TEST_EXIT(dim)("no problem dimension spezified!\n");

197
    for (int i = 0; i < numComponents_; i++) {
198
      sprintf(number, "%d", i);
199
200
201
      int refSet = -1;
      GET_PARAMETER(0, name_ + "->refinement set[" + number + "]", "%d", &refSet);
      if (refSet < 0) {
202
203
204
	WARNING("refinement set for component %d not set\n", i);
	refSet = 0;
      }
205
      if (meshForRefinementSet[refSet] == NULL) {
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
	Mesh *newMesh = NEW Mesh(meshName, dim);
	meshForRefinementSet[refSet] = newMesh;
	meshes_.push_back(newMesh);
      }
      componentMeshes_[i] = meshForRefinementSet[refSet];
    }
    switch(dim) {
    case 1:
      coarseningManager_ = NEW CoarseningManager1d();
      refinementManager_ = NEW RefinementManager1d();
      break;
    case 2:
      coarseningManager_ = NEW CoarseningManager2d();
      refinementManager_ = NEW RefinementManager2d();
      break;
    case 3:
      coarseningManager_ = NEW CoarseningManager3d();
      refinementManager_ = NEW RefinementManager3d();
      break;
    default:
      ERROR_EXIT("invalid dim!\n");
    }
  }

  void ProblemVec::createFESpace()
  {
    FUNCNAME("ProblemVec::createFESpace()");

    int degree = 1;
    char number[3];

    ::std::map< ::std::pair<Mesh*, int>, FiniteElemSpace*> feSpaceMap;
238
    int dim = -1;
239
240
241
242
243
    GET_PARAMETER(0, name_ + "->dim", "%d", &dim);
    TEST_EXIT(dim != -1)("no problem dimension spezified!\n");

    componentSpaces_.resize(numComponents_, NULL);

244
    for (int i = 0; i < numComponents_; i++) {
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
      sprintf(number, "%d", i);
      GET_PARAMETER(0, name_ + "->polynomial degree[" + number + "]","%d", &degree);

      TEST_EXIT(componentSpaces_[i] == NULL)("feSpace already created\n");

      if (feSpaceMap[::std::pair<Mesh*, int>(componentMeshes_[i], degree)] == NULL) {
	FiniteElemSpace *newFESpace = 
	  FiniteElemSpace::provideFESpace(NULL,
					  Lagrange::getLagrange(dim, degree),
					  componentMeshes_[i],
					  name_ + "->feSpace");
	feSpaceMap[::std::pair<Mesh*, int>(componentMeshes_[i], degree)] = newFESpace;
	feSpaces_.push_back(newFESpace);
      }
      componentSpaces_[i] = 
	feSpaceMap[::std::pair<Mesh*, int>(componentMeshes_[i], degree)];
    }

    // create dof admin for vertex dofs if neccessary
264
    for (int i = 0; i < static_cast<int>(meshes_.size()); i++) {
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
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
431
432
433
434
435
436
437
438
      if (meshes_[i]->getNumberOfDOFs(VERTEX) == 0) {
	DimVec<int> ln_dof(meshes_[i]->getDim(), DEFAULT_VALUE, 0);
	ln_dof[VERTEX]= 1;
	meshes_[i]->createDOFAdmin("vertex dofs", ln_dof);      
      }
    }
  }

  void ProblemVec::createMatricesAndVectors()
  {
    FUNCNAME("ProblemVec::createMatricesAndVectors()");

    int i;

    // === create vectors and system matrix ===

    systemMatrix_ = NEW Matrix<DOFMatrix*>(numComponents_, numComponents_);
    systemMatrix_->set(NULL);
    rhs_ = NEW SystemVector("rhs", componentSpaces_, numComponents_);
    solution_ = NEW SystemVector("solution", componentSpaces_, numComponents_);

    char number[10];
    ::std::string numberedName;
    for (i = 0; i < numComponents_; i++) {
      (*systemMatrix_)[i][i] = NEW DOFMatrix(componentSpaces_[i], 
					     componentSpaces_[i], "A_ii");
      (*systemMatrix_)[i][i]->setCoupleMatrix(false);
      sprintf(number, "[%d]", i);
      numberedName = "rhs" + ::std::string(number);
      rhs_->setDOFVector(i, NEW DOFVector<double>(componentSpaces_[i], numberedName));
      numberedName = name_ + ::std::string(number);
      solution_->setDOFVector(i, NEW DOFVector<double>(componentSpaces_[i], 
						       numberedName));
      solution_->getDOFVector(i)->refineInterpol(true);
      solution_->getDOFVector(i)->setCoarsenOperation(COARSE_INTERPOL);
      solution_->getDOFVector(i)->set(0.0);
    }

    // === create matVec ===
    matVec_ = NEW StandardMatVec<Matrix<DOFMatrix*>, SystemVector>(systemMatrix_);
  }

  void ProblemVec::createSolver()
  {
    FUNCNAME("ProblemVec::createSolver()");

    // === create solver ===
    ::std::string solverType("no");
    GET_PARAMETER(0, name_ + "->solver", &solverType);
    OEMSolverCreator<SystemVector> *solverCreator = 
      dynamic_cast<OEMSolverCreator<SystemVector>*>(
						    CreatorMap<OEMSolver<SystemVector> >
						    ::getCreator(solverType)
						    );
    TEST_EXIT(solverCreator)("no solver type\n");
    solverCreator->setName(name_ + "->solver");
    solver_ = solverCreator->create();
    solver_->initParameters();

    // === create preconditioners ===
    ::std::string preconType("no");

    PreconditionerScal *scalPrecon;
    PreconditionerVec *vecPrecon = NEW PreconditionerVec(numComponents_);

    GET_PARAMETER(0, name_ + "->solver->left precon", &preconType);
    CreatorInterface<PreconditionerScal> *preconCreator =
      CreatorMap<PreconditionerScal>::getCreator(preconType);

    int i, j;

    if (!preconCreator->isNullCreator()) {
      dynamic_cast<PreconditionerScalCreator*>(preconCreator)->
	setName(name_ + "->solver->left precon");

      for(i = 0; i < numComponents_; i++) {
	dynamic_cast<PreconditionerScalCreator*>(preconCreator)->
	  setSizeAndRow(numComponents_, i);
    
	scalPrecon = preconCreator->create();
	for(j = 0; j < numComponents_; j++) {
	  scalPrecon->setMatrix(&(*systemMatrix_)[i][j], j);
	}
	vecPrecon->setScalarPrecon(i, scalPrecon);
      }
      leftPrecon_ = vecPrecon;
    }


    vecPrecon = NEW PreconditionerVec(numComponents_);

    GET_PARAMETER(0, name_ + "->solver->right precon", &preconType);
    preconCreator = 
      CreatorMap<PreconditionerScal>::getCreator(preconType);

    if(!preconCreator->isNullCreator()) {
      dynamic_cast<PreconditionerScalCreator*>(preconCreator)->
	setName(name_ + "->solver->left precon");


      for(i = 0; i < numComponents_; i++) {
	dynamic_cast<PreconditionerScalCreator*>(preconCreator)->
	  setSizeAndRow(numComponents_, i);
    
	scalPrecon = preconCreator->create();
	for(j = 0; j < numComponents_; j++) {
	  scalPrecon->setMatrix(&(*systemMatrix_)[i][j], j);
	}
	vecPrecon->setScalarPrecon(i, scalPrecon);
      }
      rightPrecon_ = vecPrecon;
    }


    // === create vector creator ===
    solver_->setVectorCreator(NEW SystemVector::Creator("temp",
							componentSpaces_, 
							numComponents_));
  }

  void ProblemVec::createEstimator()
  {
    FUNCNAME("ProblemVec::createEstimator()");

    int i, j;

    // create and set leaf data prototype
    for(i = 0; i < static_cast<int>(meshes_.size()); i++) {
      meshes_[i]->setElementDataPrototype
	(NEW LeafDataEstimatableVec(NEW LeafDataCoarsenableVec));
    }  

    char number[3];
    ::std::string estName;

    for(i = 0; i < numComponents_; i++) {
      TEST_EXIT(estimator_[i] == NULL)("estimator already created\n");
      sprintf(number, "%d", i);
      estName = name_ + "->estimator[" + ::std::string(number) + "]";

      // === create estimator ===
      ::std::string estimatorType("no");
      GET_PARAMETER(0, estName, &estimatorType);
      EstimatorCreator *estimatorCreator = 
	dynamic_cast<EstimatorCreator*>(
					CreatorMap<Estimator>::getCreator(estimatorType));
      if(estimatorCreator) {
	estimatorCreator->setName(estName);
	estimatorCreator->setRow(i);
	if(estimatorType == "recovery") {
	  dynamic_cast<RecoveryEstimator::Creator*>(estimatorCreator)->
	    setSolution(solution_->getDOFVector(i));
	}
	estimator_[i] = estimatorCreator->create();
      }


      if(estimator_[i]) {
	for(j=0; j < numComponents_; j++) {
	  estimator_[i]->addSystem((*systemMatrix_)[i][j], 
				   solution_->getDOFVector(j), 
				   rhs_->getDOFVector(j));
	}
      }
    }
  }

  void ProblemVec::createMarker()
  {
    FUNCNAME("ProblemVec::createMarker()");

    ::std::string numberedName;
    char number[10];
    int numMarkersCreated = 0;
439
440

    for (int i = 0; i < numComponents_; i++) {
441
442
      sprintf(number, "[%d]", i);
      numberedName = name_ + "->marker" + ::std::string(number);
443
444
445
446
      marker_[i] = Marker::createMarker(numberedName, i);
      if (marker_[i]) {
	numMarkersCreated++;
	if (numMarkersCreated > 1)
447
448
449
450
451
452
453
454
455
456
457
458
459
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
	  marker_[i]->setMaximumMarking(true);
      }
    }
  }

  void ProblemVec::createFileWriter()
  {
    FUNCNAME("ProblemVec::createFileWriter()");
  

    // Create one filewriter for all components of the problem
    ::std::string numberedName  = name_ + "->output";
    ::std::string filename = "";
    GET_PARAMETER(0, numberedName + "->filename", &filename);

    if (filename != "") {
      ::std::vector< DOFVector<double>* > solutionList(numComponents_);

      for (int i = 0; i < numComponents_; i++) {
	TEST_EXIT(componentMeshes_[0] == componentMeshes_[i])
	  ("All Meshes have to be equal to write a vector file.\n");

	solutionList[i] = solution_->getDOFVector(i);
      }

      fileWriters_.push_back(NEW FileWriter(numberedName,
					    componentMeshes_[0],
					    solutionList));
    }


    // Create own filewriters for each components of the problem
    char number[10];
    for (int i = 0; i < numComponents_; i++) {
      sprintf(number, "[%d]", i);
      numberedName  = name_ + "->output" + ::std::string(number);
      filename = "";
      GET_PARAMETER(0, numberedName + "->filename", &filename);

      if (filename != "") {
	fileWriters_.push_back(NEW FileWriter(numberedName, 
					      componentMeshes_[i], 
					      solution_->getDOFVector(i)));
      }
    }


    // Check for serializer
    int writeSerialization = 0;
    GET_PARAMETER(0, name_ + "->write serialization", "%d", &writeSerialization);
    if (writeSerialization) {
      MSG("Use are using the obsolete parameter: %s->write serialization\n", name_.c_str());
      MSG("Please use instead the following parameter: %s->output->write serialization\n", name_.c_str());
      ERROR_EXIT("Usage of an obsolete parameter (see message above)!\n");
    }

    GET_PARAMETER(0, name_ + "->output->write serialization", "%d", &writeSerialization);
    if (writeSerialization) {
      fileWriters_.push_back(NEW Serializer<ProblemVec>(this));
    }
  }

  void ProblemVec::doOtherStuff()
  {
  }

  void ProblemVec::solve(AdaptInfo *adaptInfo) 
  {
    FUNCNAME("Problem::solve()");

    if (!solver_) {
      WARNING("no solver\n");
      return;
    }

#ifdef _OPENMP
    double wtime = omp_get_wtime();
#endif

    clock_t first = clock();
    int iter = solver_->solve(matVec_, solution_, rhs_, leftPrecon_, rightPrecon_);   
528
    
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
#ifdef _OPENMP
    INFO(info_, 8)("solution of discrete system needed %.5f seconds system time / %.5f seconds wallclock time\n",
		   TIME_USED(first, clock()),
		   omp_get_wtime() - wtime);
#else
    INFO(info_, 8)("solution of discrete system needed %.5f seconds\n",
		   TIME_USED(first, clock()));
#endif


    adaptInfo->setSolverIterations(iter);
    adaptInfo->setMaxSolverIterations(solver_->getMaxIterations());
    adaptInfo->setSolverTolerance(solver_->getTolerance());
    adaptInfo->setSolverResidual(solver_->getResidual());
  }

  void ProblemVec::estimate(AdaptInfo *adaptInfo) 
  {
    FUNCNAME("ProblemVec::estimate()");

    clock_t first = clock();

551
552
553
554
#ifdef _OPENMP
    double wtime = omp_get_wtime();
#endif

555
556
557
558
559
560
561
562
563
564
565
566
567
568
    for (int i = 0; i < numComponents_; i++) {
      Estimator *scalEstimator = estimator_[i];

      if (scalEstimator) {
	scalEstimator->estimate(adaptInfo->getTimestep());
	adaptInfo->setEstSum(scalEstimator->getErrorSum(), i);
	adaptInfo->setEstMax(scalEstimator->getErrorMax(), i);
	adaptInfo->setTimeEstSum(scalEstimator->getTimeEst(), i);
	adaptInfo->setTimeEstMax(scalEstimator->getTimeEstMax(), i);
      } else {
	WARNING("no estimator for component %d\n" , i);
      }
    }

569
570
571
572
573
574
575
576
577
578
#ifdef _OPENMP
    INFO(info_, 8)("estimation of the error needed %.5f seconds system time / %.5f seconds wallclock time\n",
		   TIME_USED(first, clock()),
		   omp_get_wtime() - wtime);
#else
    INFO(info_, 8)("estimation of the error needed %.5f seconds\n",
		   TIME_USED(first, clock()));

#endif

579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
  }

  Flag ProblemVec::markElements(AdaptInfo *adaptInfo) 
  {
    FUNCNAME("ProblemVec::markElements()");

    // to enforce albert-like behavior: refinement even if space tolerance
    // here is reached already because of time adaption
    allowFirstRefinement();

    Flag markFlag = 0;
    for (int i = 0; i < numComponents_; i++) {
      if (marker_[i]) {
	markFlag |= marker_[i]->markMesh(adaptInfo, componentMeshes_[i]);
      } else {
	WARNING("no marker for component %d\n", i);
      }
    }
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
    return markFlag;
  }

  Flag ProblemVec::refineMesh(AdaptInfo *adaptInfo) 
  {
    FUNCNAME("ProblemVec::refineMesh()");

    int numMeshes = static_cast<int>(meshes_.size());
    Flag refineFlag = 0;
    for (int i = 0; i < numMeshes; i++) {
      refineFlag |= refinementManager_->refineMesh(meshes_[i]);
    }
    return refineFlag;
  }

  Flag ProblemVec::coarsenMesh(AdaptInfo *adaptInfo) 
  {
    FUNCNAME("ProblemVec::coarsenMesh()");

    int i, numMeshes = static_cast<int>(meshes_.size());
    Flag coarsenFlag = 0;
    for(i = 0; i < numMeshes; i++) {
      if(adaptInfo->isCoarseningAllowed(i)) {
	coarsenFlag |= coarseningManager_->coarsenMesh(meshes_[i]);

	WARNING("coarsening for component %d no allowed\n", i);
      }
    }
    return coarsenFlag;
  }

  Flag ProblemVec::oneIteration(AdaptInfo *adaptInfo, Flag toDo)
  {
    FUNCNAME("ProblemVec::oneIteration()");

    if (allowFirstRef_) {
      for (int i = 0; i < numComponents_; i++) {
	adaptInfo->allowRefinement(true, i);
      }
      allowFirstRef_ = false;
    } else {
      for (int i = 0; i < numComponents_; i++) {
	if (adaptInfo->spaceToleranceReached(i)) {
	  adaptInfo->allowRefinement(false, i);
	} else {
	  adaptInfo->allowRefinement(true, i);	
	}
      }
    }

    return StandardProblemIteration::oneIteration(adaptInfo, toDo);
  }

  void ProblemVec::buildAfterCoarsen(AdaptInfo *adaptInfo, Flag flag) 
  {
    FUNCNAME("ProblemVec::buildAfterCoarsen()");

    clock_t first = clock();

657
658
659
660
#ifdef _OPENMP
    double wtime = omp_get_wtime();
#endif

Thomas Witkowski's avatar
Thomas Witkowski committed
661
    for (int i = 0; i < static_cast<int>(meshes_.size()); i++) {
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
      meshes_[i]->dofCompress();
    }

    Flag assembleFlag = 
      flag | 
      (*systemMatrix_)[0][0]->getAssembleFlag() | 
      rhs_->getDOFVector(0)->getAssembleFlag()   |
      Mesh::CALL_LEAF_EL                        | 
      Mesh::FILL_COORDS                         |
      Mesh::FILL_DET                            |
      Mesh::FILL_GRD_LAMBDA |
      Mesh::FILL_NEIGH;

    if (useGetBound_) {
      assembleFlag |= Mesh::FILL_BOUND;
    }

Thomas Witkowski's avatar
Thomas Witkowski committed
679

680
681
682
683
684
685
686
    for (int i = 0; i < numComponents_; i++) {
      MSG("%d DOFs for %s\n", 
	  componentSpaces_[i]->getAdmin()->getUsedSize(), 
	  componentSpaces_[i]->getName().c_str());

      rhs_->getDOFVector(i)->set(0.0);
      for (int j = 0; j < numComponents_; j++) {
687
688
689
690
691
692
693
	if ((*systemMatrix_)[i][j]) {
	  // The matrix should not be deleted, if it was assembled before
	  // and it is marked to be assembled only once.
	  if (!(assembleMatrixOnlyOnce_[i][j] && assembledMatrix_[i][j])) {
	    (*systemMatrix_)[i][j]->clear();
	  }
	}
694
695
696
      }
    }

697
698
699
700
701
    int i;
#ifdef _OPENMP
#pragma omp parallel for 
#endif
    for (i = 0; i < numComponents_; i++) {
702
703
      const BasisFunction *basisFcts = componentSpaces_[i]->getBasisFcts();

704
      for (int j = 0; j < numComponents_; j++) {
705
706
707
708
	// Only if this variable is true, the current matrix will be assembled.	
	bool assembleMatrix = true;
	// The DOFMatrix which should be assembled (or not, if assembleMatrix
	// will be set to false).
709
710
	DOFMatrix *matrix = (*systemMatrix_)[i][j];

711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
	// If the matrix was assembled before and it is marked to be assembled
	// only once, it will not be assembled.
	if (assembleMatrixOnlyOnce_[i][j] && assembledMatrix_[i][j]) {
	  assembleMatrix = false;
	}
	// If there is no DOFMatrix (e.g. if it is completly 0), do not assemble.
	if (!matrix) {
	  assembleMatrix = false;
	}

	// If the matrix should not be assembled, the rhs vector has to be considered.
	// This will be only done, if i == j. So, if both is not true, we can jump
	// to the next matrix.
	if (!assembleMatrix && i != j) {
	  continue;
	}

728
729
730
	if (assembleMatrix && matrix->getBoundaryManager())
	  matrix->getBoundaryManager()->initMatrix(matrix);

Thomas Witkowski's avatar
Thomas Witkowski committed
731
	
732
733
734
	if (componentMeshes_[i] != componentMeshes_[j]) {
	  ERROR_EXIT("not yet\n");
	} else {
735
736
	  BoundaryType *bound = NULL;
	  if (useGetBound_) {
737
	    bound = GET_MEMORY(BoundaryType, basisFcts->getNumber());
738
	  }
Thomas Witkowski's avatar
Thomas Witkowski committed
739
	  
Thomas Witkowski's avatar
Thomas Witkowski committed
740
741
742
	  TraverseStack stack;
	  ElInfo *elInfo = stack.traverseFirst(componentMeshes_[i], -1, assembleFlag);
	  
743
	  while (elInfo) {
744
	    if (useGetBound_) {
745
	      basisFcts->getBound(elInfo, bound);
746
747
	    }
	          
748
	    if (assembleMatrix) {
749
	      matrix->assemble(1.0, elInfo, bound);
750
	      
751
752
	      if (matrix->getBoundaryManager()) {
		matrix->
Thomas Witkowski's avatar
Thomas Witkowski committed
753
		  getBoundaryManager()->
754
755
		  fillBoundaryConditions(elInfo, matrix);		
	      }		      
756
	    }
757
	    
758
759
760
	    if (i == j) {
	      rhs_->getDOFVector(i)->assemble(1.0, elInfo, bound);
	    }
Thomas Witkowski's avatar
Thomas Witkowski committed
761
	    
762
763
	    elInfo = stack.traverseNext(elInfo);
	  }
Thomas Witkowski's avatar
Thomas Witkowski committed
764
	  
765

766
	  if (assembleMatrix && matrix->getBoundaryManager())
Thomas Witkowski's avatar
Thomas Witkowski committed
767
	  	    matrix->getBoundaryManager()->exitMatrix(matrix);	  
768
769

	  if (useGetBound_) {
770
	    FREE_MEMORY(bound, BoundaryType, basisFcts->getNumber());
Thomas Witkowski's avatar
Thomas Witkowski committed
771
	  }	  
772
	}
Thomas Witkowski's avatar
Thomas Witkowski committed
773
	
774
	assembledMatrix_[i][j] = true;
775
776
777
778
      }

      // fill boundary conditions
      if (rhs_->getDOFVector(i)->getBoundaryManager())
779
	rhs_->getDOFVector(i)->getBoundaryManager()->initVector(rhs_->getDOFVector(i));     
780
      
781
      if (solution_->getDOFVector(i)->getBoundaryManager())
782
      	solution_->getDOFVector(i)->getBoundaryManager()->initVector(solution_->getDOFVector(i));
783

Thomas Witkowski's avatar
Thomas Witkowski committed
784
785
      TraverseStack stack;
      ElInfo *elInfo = stack.traverseFirst(componentMeshes_[i], -1, assembleFlag);
786
      while (elInfo) {
787
	if (rhs_->getDOFVector(i)->getBoundaryManager())
788
789
	  rhs_->getDOFVector(i)->getBoundaryManager()->
	    fillBoundaryConditions(elInfo, rhs_->getDOFVector(i));
790
791

	if (solution_->getDOFVector(i)->getBoundaryManager())
792
793
794
795
	  solution_->getDOFVector(i)->getBoundaryManager()->
	    fillBoundaryConditions(elInfo, solution_->getDOFVector(i));
	elInfo = stack.traverseNext(elInfo);
      }
796
      
797
798
799
      if (rhs_->getDOFVector(i)->getBoundaryManager())
	rhs_->getDOFVector(i)->getBoundaryManager()->exitVector(rhs_->getDOFVector(i));
      if (solution_->getDOFVector(i)->getBoundaryManager())
Thomas Witkowski's avatar
Thomas Witkowski committed
800
      solution_->getDOFVector(i)->getBoundaryManager()->exitVector(solution_->getDOFVector(i));    
801
    }
802
803
804
805
806
807

#ifdef _OPENMP
    INFO(info_, 8)("buildAfterCoarsen needed %.5f seconds system time / %.5f seconds wallclock time\n",
		   TIME_USED(first, clock()),
		   omp_get_wtime() - wtime);
#else
808
    INFO(info_, 8)("buildAfterCoarsen needed %.5f seconds\n",
Thomas Witkowski's avatar
Thomas Witkowski committed
809
		   TIME_USED(first, clock()));
810
#endif
811
812
813
814
815
816
  }

  void ProblemVec::writeFiles(AdaptInfo *adaptInfo, bool force) 
  {
    FUNCNAME("ProblemVec::writeFiles()");

817
818
819
820
821
822
823
824
825
826
827
828
829
    clock_t first = clock();

#ifdef _OPENMP
    double wtime = omp_get_wtime();
#endif

    int i;
    int size = static_cast<int>(fileWriters_.size());
#ifdef _OPENMP
#pragma omp parallel for schedule(static, 1)
#endif
    for (i = 0; i < size; i++) {
      fileWriters_[i]->writeFiles(adaptInfo, force);
830
    }
831
832
833
834
835
836
837
838
839
    
#ifdef _OPENMP
    INFO(info_, 8)("writeFiles needed %.5f seconds system time / %.5f seconds wallclock time\n",
		   TIME_USED(first, clock()),
		   omp_get_wtime() - wtime);
#else
    INFO(info_, 8)("writeFiles needed %.5f seconds\n",
		   TIME_USED(first, clock()));
#endif
840
841
  }

Thomas Witkowski's avatar
Thomas Witkowski committed
842
843
844
845
846
847
848
  void ProblemVec::writeDelayedFiles()
  {
    for (int i = 0; i < static_cast<int>(fileWriters_.size()); i++) {
      fileWriters_[i]->writeDelayedFiles();
    }
  }

Thomas Witkowski's avatar
Thomas Witkowski committed
849
850
851
852
853
854
855
856
857
858
  bool ProblemVec::existsDelayedCalculation()
  {
    for (int i = 0; i < static_cast<int>(fileWriters_.size()); i++) {
      if (fileWriters_[i]->isWritingDelayed())
	return true;   
    }

    return false;
  }

859
860
861
862
863
864
865
866
867
868
869
870
871
872
  void ProblemVec::interpolInitialSolution(::std::vector<AbstractFunction<double, WorldVector<double> >*> *fct) 
  {
    FUNCNAME("ProblemVec::interpolInitialSolution()");

    solution_->interpol(fct);
  }

  void ProblemVec::addMatrixOperator(Operator *op, 
				     int i, int j,
				     double *factor,
				     double *estFactor)
  {
    FUNCNAME("ProblemVec::addMatrixOperator()");

873
    if (!(*systemMatrix_)[i][j]) {
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
      TEST_EXIT(i != j)("should have been created already\n");
      (*systemMatrix_)[i][j] = NEW DOFMatrix(componentSpaces_[i],
					     componentSpaces_[j],
					     "");
      (*systemMatrix_)[i][j]->setCoupleMatrix(true);

      (*systemMatrix_)[i][j]->getBoundaryManager()->
	setBoundaryConditionMap((*systemMatrix_)[i][i]->getBoundaryManager()->
				getBoundaryConditionMap());
    }    
    (*systemMatrix_)[i][j]->addOperator(op, factor, estFactor);
  }

  void ProblemVec::addVectorOperator(Operator *op, int i,
				     double *factor,
				     double *estFactor)
  {
    FUNCNAME("ProblemVec::addVectorOperator()");

    rhs_->getDOFVector(i)->addOperator(op, factor, estFactor);
  }

  void ProblemVec::addDirichletBC(BoundaryType type, int system,
				  AbstractFunction<double, WorldVector<double> >* b)
  {
    FUNCNAME("ProblemVec::addDirichletBC()");

    DirichletBC *dirichlet = new DirichletBC(type, 
					     b, 
					     componentSpaces_[system]);
    for (int i = 0; i < numComponents_; i++) {
      if (systemMatrix_ && (*systemMatrix_)[system][i]) {
	(*systemMatrix_)[system][i]->getBoundaryManager()->addBoundaryCondition(dirichlet);
      }
    }
909

910
911
    if (rhs_)
      rhs_->getDOFVector(system)->getBoundaryManager()->addBoundaryCondition(dirichlet);
912

913
914
915
916
917
918
919
920
921
922
923
924
925
    if (solution_)
      solution_->getDOFVector(system)->getBoundaryManager()->addBoundaryCondition(dirichlet);
  }

  void ProblemVec::addNeumannBC(BoundaryType type, int row, int col, 
				AbstractFunction<double, WorldVector<double> > *n)
  {
    FUNCNAME("ProblemVec::addNeumannBC()");

    NeumannBC *neumann = 
      new NeumannBC(type, n, 
		    componentSpaces_[row], 
		    componentSpaces_[col]);
926
    if (rhs_)
927
928
929
930
931
932
933
934
935
936
937
938
939
      rhs_->getDOFVector(row)->getBoundaryManager()->addBoundaryCondition(neumann);
  }

  void ProblemVec::addRobinBC(BoundaryType type, int row, int col, 
			      AbstractFunction<double, WorldVector<double> > *n,
			      AbstractFunction<double, WorldVector<double> > *r)
  {
    FUNCNAME("ProblemVec::addRobinBC()");

    RobinBC *robin = 
      new RobinBC(type, n, r, 
		  componentSpaces_[row], 
		  componentSpaces_[col]);
940
    if (rhs_)
941
      rhs_->getDOFVector(row)->getBoundaryManager()->addBoundaryCondition(robin);
942
943

    if (systemMatrix_ && (*systemMatrix_)[row][col]) {
944
945
946
947
948
949
950
951
952
953
954
955
      (*systemMatrix_)[row][col]->getBoundaryManager()->addBoundaryCondition(robin);
    }
  }

  void ProblemVec::addPeriodicBC(BoundaryType type, int row, int col) 
  {
    FUNCNAME("ProblemVec::addPeriodicBC()");

    FiniteElemSpace *feSpace = componentSpaces_[row];

    PeriodicBC *periodic = new PeriodicBC(type, feSpace);

956
    if (systemMatrix_ && (*systemMatrix_)[row][col]) 
957
      (*systemMatrix_)[row][col]->getBoundaryManager()->addBoundaryCondition(periodic);
958
959

    if (rhs_) 
960
961
962
963
      rhs_->getDOFVector(row)->getBoundaryManager()->
	addBoundaryCondition(periodic);
  }

964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
  void ProblemVec::writeResidualMesh(AdaptInfo *adaptInfo, const ::std::string name)
  {
    FUNCNAME("ProblemVec::writeResidualMesh()");

    Mesh *mesh = this->getMesh(0);
    FiniteElemSpace *fe = this->getFESpace(0);
    
    std::map<int, double> vec;
    
    TraverseStack stack;
    ElInfo *elInfo = stack.traverseFirst(mesh,
					 -1, 
					 Mesh::CALL_LEAF_EL | 
					 Mesh::FILL_COORDS);
    
    while (elInfo) {		  
      Element *el = elInfo->getElement();
      double lError = el->getEstimation(0);
      
      vec[elInfo->getElement()->getIndex()] = lError;
      elInfo = stack.traverseNext(elInfo);
    }
    
    ElementFileWriter fw(name, mesh, fe, vec);
    fw.writeFiles(adaptInfo, true);    
  }

991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
  void ProblemVec::serialize(::std::ostream &out) 
  {
    FUNCNAME("ProblemVec::serialize()");

    SerializerUtil::serializeBool(out, &allowFirstRef_);
    
    for (int i = 0; i < static_cast<int>(meshes_.size()); i++) {
      meshes_[i]->serialize(out);
    }

    solution_->serialize(out);
  }

  void ProblemVec::deserialize(::std::istream &in) 
  {
    FUNCNAME("ProblemVec::deserialize()");

    SerializerUtil::deserializeBool(in, &allowFirstRef_);

    for (int i = 0; i < static_cast<int>(meshes_.size()); i++) {
      meshes_[i]->deserialize(in);
    }

    solution_->deserialize(in);
  }
}