#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 "DualTraverse.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"
#include "Flag.h"
#include "TraverseParallel.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 {
	int globalRefinements = 0;
	GET_PARAMETER(0, meshes_[0]->getName() + "->global refinements", "%d", 
		      &globalRefinements);

	// 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();	    
	    refinementManager_->globalRefine(meshes_[i], globalRefinements);
	  }
	}	
      }
    }

    doOtherStuff();
  }

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

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

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

    for (int i = 0; i < nComponents; i++) {
      sprintf(number, "%d", i);
      int refSet = -1;
      GET_PARAMETER(0, name_ + "->refinement set[" + number + "]", "%d", &refSet);
      if (refSet < 0) {
	refSet = 0;
      }
      if (meshForRefinementSet[refSet] == NULL) {
	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;
    int dim = -1;
    GET_PARAMETER(0, name_ + "->dim", "%d", &dim);
    TEST_EXIT(dim != -1)("no problem dimension specified!\n");

    componentSpaces.resize(nComponents, NULL);

    for (int i = 0; i < nComponents; i++) {
      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
    for (int i = 0; i < static_cast<int>(meshes_.size()); i++) {
      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()");

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

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

    char number[10];
    std::string numberedName;
    for (int i = 0; i < nComponents; 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(nComponents);

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


    vecPrecon = NEW PreconditionerVec(nComponents);

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


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

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

    for (int i = 0; i < nComponents; i++) {
      sprintf(number, "[%d]", i);
      numberedName = name_ + "->marker" + std::string(number);
      marker[i] = Marker::createMarker(numberedName, i);
      if (marker[i]) {
	numMarkersCreated++;
	if (numMarkersCreated > 1)
	  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(nComponents);

      for (int i = 0; i < nComponents; 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 < nComponents; 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_);   
    
#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();

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

    if (computeExactError) {
      computeError(adaptInfo);
    } else {
      for (int i = 0; i < nComponents; 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);
	}
      }
    }

#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

  }

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

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

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

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

    int nMeshes = static_cast<int>(meshes_.size());
    Flag coarsenFlag = 0;
    for (int i = 0; i < nMeshes; 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 < nComponents; i++) {
	adaptInfo->allowRefinement(true, i);
      }
      allowFirstRef_ = false;
    } else {
      for (int i = 0; i < nComponents; 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();

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

    for (int i = 0; i < static_cast<int>(meshes_.size()); i++) {
      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;
    }


    for (int i = 0; i < nComponents; 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 < nComponents; j++) {
	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();
	  }
	}
      }
    }

    for (int i = 0; i < nComponents; i++) {
      for (int j = 0; j < nComponents; j++) {
	// 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).
	DOFMatrix *matrix = (*systemMatrix_)[i][j];

	// 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;
	}

	if (assembleMatrix && matrix->getBoundaryManager())
	  matrix->getBoundaryManager()->initMatrix(matrix);

	if (componentSpaces[i] == componentSpaces[j]) {
	  //	  std::cout << "--- assembleOnOneMesh ---\n";
	  assembleOnOneMesh(componentSpaces[i],
			    assembleFlag,
			    assembleMatrix ? matrix : NULL,
			    (i == j) ? rhs_->getDOFVector(i) : NULL);
	  //	  std::cout << "--- Finished ---\n";
	  //	  std::cout << "--- Dim: " << matrix->getUsedSize() << "x" << matrix->getNumCols() << "---\n";
	} else {
	  //	  std::cout << "--- assembleOnDifMesh ---\n";
	  assembleOnDifMeshes(componentSpaces[i], componentSpaces[j],
			      assembleFlag,
			      assembleMatrix ? matrix : NULL,
			      (i == j) ? rhs_->getDOFVector(i) : NULL);	  
	  //	  std::cout << "--- Finished ---\n";
	  //	  std::cout << "--- Dim: " << matrix->getUsedSize() << "x" << matrix->getNumCols() << "---\n";
	}

	if (assembleMatrix && matrix->getBoundaryManager())
	  matrix->getBoundaryManager()->exitMatrix(matrix);	  
	
	assembledMatrix_[i][j] = true;
      }

      // fill boundary conditions
      if (rhs_->getDOFVector(i)->getBoundaryManager())
	rhs_->getDOFVector(i)->getBoundaryManager()->initVector(rhs_->getDOFVector(i));     
      
      if (solution_->getDOFVector(i)->getBoundaryManager())
      	solution_->getDOFVector(i)->getBoundaryManager()->initVector(solution_->getDOFVector(i));

      TraverseStack stack;
      ElInfo *elInfo = stack.traverseFirst(componentMeshes[i], -1, assembleFlag);
      while (elInfo) {
	if (rhs_->getDOFVector(i)->getBoundaryManager())
	  rhs_->getDOFVector(i)->getBoundaryManager()->
	    fillBoundaryConditions(elInfo, rhs_->getDOFVector(i));

	if (solution_->getDOFVector(i)->getBoundaryManager())
	  solution_->getDOFVector(i)->getBoundaryManager()->
	    fillBoundaryConditions(elInfo, solution_->getDOFVector(i));

	elInfo = stack.traverseNext(elInfo);
      }

      if (rhs_->getDOFVector(i)->getBoundaryManager())
	rhs_->getDOFVector(i)->getBoundaryManager()->exitVector(rhs_->getDOFVector(i));
      if (solution_->getDOFVector(i)->getBoundaryManager())
	solution_->getDOFVector(i)->getBoundaryManager()->exitVector(solution_->getDOFVector(i));    
    }

#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
    INFO(info_, 8)("buildAfterCoarsen needed %.5f seconds\n",
		   TIME_USED(first, clock()));
#endif
  }

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

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

  void ProblemVec::writeDelayedFiles()
  {
    for (int i = 0; i < static_cast<int>(fileWriters_.size()); i++) {
      fileWriters_[i]->writeDelayedFiles();
    }
  }

  bool ProblemVec::existsDelayedCalculation()
  {
    for (int i = 0; i < static_cast<int>(fileWriters_.size()); i++) {
      if (fileWriters_[i]->isWritingDelayed())
	return true;   
    }

    return false;
  }

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

    if (!(*systemMatrix_)[i][j]) {
      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 < nComponents; i++) {
      if (systemMatrix_ && (*systemMatrix_)[system][i]) {
	(*systemMatrix_)[system][i]->getBoundaryManager()->addBoundaryCondition(dirichlet);
      }
    }

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

    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]);
    if (rhs_)
      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]);
    if (rhs_)
      rhs_->getDOFVector(row)->getBoundaryManager()->addBoundaryCondition(robin);

    if (systemMatrix_ && (*systemMatrix_)[row][col]) {
      (*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);

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

    if (rhs_) 
      rhs_->getDOFVector(row)->getBoundaryManager()->
	addBoundaryCondition(periodic);
  }

  void ProblemVec::assembleOnOneMesh(FiniteElemSpace *feSpace, Flag assembleFlag,
				     DOFMatrix *matrix, DOFVector<double> *vector)
  {
    Mesh *mesh = feSpace->getMesh();
    const BasisFunction *basisFcts = feSpace->getBasisFcts();

#ifdef _OPENMP
    TraverseParallelStack stack;
#else
    TraverseStack stack;
#endif   

#ifdef _OPENMP
#pragma omp parallel
#endif
    {
      BoundaryType *bound = useGetBound_ ? GET_MEMORY(BoundaryType, basisFcts->getNumber()) : NULL;

      // Create for every thread its private matrix and vector, on that
      // the thread will assemble its part of the mesh.
      DOFMatrix *tmpMatrix = NULL;
      DOFVector<double> *tmpVector = NULL; 

      if (matrix) {
	tmpMatrix = NEW DOFMatrix(matrix->getRowFESpace(),
				  matrix->getColFESpace(),
				  "tmp");

	// Copy the global matrix to the private matrix, because we need the
	// operators defined on the global matrix in the private one. Only the
	// values have to be set to zero.
	*tmpMatrix = *matrix;
	tmpMatrix->clear();
      }

      if (vector) {
	tmpVector = NEW DOFVector<double>(vector->getFESpace(), "tmp");

	// Copy the global vector to the private vector, because we need the
	// operatirs defined on the global vector in the private one. But set
	// the values to zero of the private vector after copying.
	*tmpVector = *vector;
	tmpVector->set(0.0);
      }

      // Because we are using the parallel traverse stack, each thread will
      // traverse only a part of the mesh.
      ElInfo *elInfo = stack.traverseFirst(mesh, -1, assembleFlag);
      while (elInfo) {
	if (useGetBound_) {
	  basisFcts->getBound(elInfo, bound);
	}
	
	if (matrix) {
	  tmpMatrix->assemble(1.0, elInfo, bound);
	  
	  // Take the matrix boundary manager from the public matrix,
	  // but assemble the boundary conditions on the thread private matrix.
	  if (matrix->getBoundaryManager()) {
	    matrix->getBoundaryManager()->
	      fillBoundaryConditions(elInfo, tmpMatrix);
	  }		      
	}
	
	if (vector) {
	  tmpVector->assemble(1.0, elInfo, bound);
	}
	
	elInfo = stack.traverseNext(elInfo);
      }

      // After mesh traverse, all thread have to added their private matrices and
      // vectors to the global public matrix and public vector. Therefore, this is 
      // a critical section, which is allowed to be executed by on thread only at 
      // the same time.

      if (matrix) {
#ifdef _OPENMP
#pragma omp critical
#endif
	{
	  addDOFMatrix(matrix, tmpMatrix);

	  // Remove rows corresponding to DOFs on a Dirichlet boundary.
	  matrix->removeRowsWithDBC(tmpMatrix->getApplyDBCs());
	}

	DELETE tmpMatrix;
      }

      if (vector) {
#ifdef _OPENMP
#pragma omp critical
#endif
	*vector += *tmpVector;

	DELETE tmpVector;
      }

      if (useGetBound_) {
	FREE_MEMORY(bound, BoundaryType, basisFcts->getNumber());
      }	      

    } // pragma omp parallel

  }

  void ProblemVec::assembleOnDifMeshes(FiniteElemSpace *rowFeSpace, FiniteElemSpace *colFeSpace,
				       Flag assembleFlag,
				       DOFMatrix *matrix, DOFVector<double> *vector)
  {
    Mesh *rowMesh = rowFeSpace->getMesh();
    Mesh *colMesh = colFeSpace->getMesh();

    const BasisFunction *basisFcts = rowFeSpace->getBasisFcts();
    BoundaryType *bound = NULL;
    if (useGetBound_) {
      bound = GET_MEMORY(BoundaryType, basisFcts->getNumber());
    }

    DualTraverse dualTraverse;
    ElInfo *rowElInfo, *colElInfo;
    ElInfo *largeElInfo, *smallElInfo;

    dualTraverse.setFillSubElemInfo(true);
    bool cont = dualTraverse.traverseFirst(rowMesh, colMesh, -1, -1,
					   assembleFlag, assembleFlag,
					   &rowElInfo, &colElInfo,
					   &smallElInfo, &largeElInfo);
    while (cont) {
      Element *rowElem = rowElInfo->getElement();
      Element *colElem = colElInfo->getElement();

      if (rowElInfo->getLevel() != colElInfo->getLevel()) {
	if (smallElInfo == rowElInfo) {
	  std::cout << "Row = small\n";
	  ERROR_EXIT("NOCH EIN PAAR GEDANKEN MACHEN!\n");
	} else {
	  std::cout << "Row = large\n";
	}

	if (largeElInfo == colElInfo) 
	  std::cout << "Col = large\n";
	else
	  std::cout << "Col = small\n";

	if (useGetBound_) {
	  basisFcts->getBound(rowElInfo, bound);
	}

	if (matrix) {
	  matrix->assemble(1.0, rowElInfo, colElInfo, bound);
	}
      } else {
	if (useGetBound_) {
	  basisFcts->getBound(rowElInfo, bound);
	}
      
	if (matrix) {
	  matrix->assemble(1.0, rowElInfo, bound);
	  
	  if (matrix->getBoundaryManager()) {
	    matrix->getBoundaryManager()->
	      fillBoundaryConditions(rowElInfo, matrix);
	  }		      
	}
      }
      
      if (vector) {
	vector->assemble(1.0, rowElInfo, bound);
      }

      cont = dualTraverse.traverseNext(&rowElInfo, &colElInfo,
				       &smallElInfo, &largeElInfo);
    }

    if (useGetBound_) {
      FREE_MEMORY(bound, BoundaryType, basisFcts->getNumber());
    }
  }

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

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

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

    for (int i = 0; i < nComponents; i++) {		
      TEST_EXIT(exactSolutionFcts[i])("No solution function given!\n");

      // Compute the difference between exact and computed solution
      DOFVector<double> *tmp = NEW DOFVector<double>(componentSpaces[i], "tmp");
      tmp->interpol(exactSolutionFcts[i]);
      double solMax = tmp->absMax();
      *tmp -= *(solution_->getDOFVector(i));
      
      MSG("L2    error = %.8e\n", tmp->L2Norm());
      MSG("L-inf error = %.8e\n", tmp->absMax() / solMax);
      
      adaptInfo->setEstSum(tmp->absMax() / solMax, i);
      adaptInfo->setEstMax(tmp->absMax() / solMax, i);
      
      // To set element estimates, compute a vector with the difference
      // between exact and computed solution for each DOF.
      DOFVector<double> *sol = NEW DOFVector<double>(componentSpaces[i], "tmp");
      sol->interpol(exactSolutionFcts[i]);
      DOFVector<double>::Iterator it1(sol, USED_DOFS);
      DOFVector<double>::Iterator it2(tmp, USED_DOFS);
      for (it1.reset(), it2.reset(); !it1.end(); ++it1, ++it2) {
	if ((abs(*it1) <= DBL_TOL) || (abs(*it2) <= DBL_TOL)) {
	  *it2 = 0.0;
	} else {
	  *it2 = abs(*it2 / *it1);
	}
      }

      // Compute estimate for every mesh element
      Vector<DegreeOfFreedom> locInd(componentSpaces[i]->getBasisFcts()->getNumber());
      TraverseStack stack;
      ElInfo *elInfo = stack.traverseFirst(componentMeshes[i], -1, Mesh::CALL_LEAF_EL);
      while (elInfo) {
	componentSpaces[i]->getBasisFcts()->getLocalIndicesVec(elInfo->getElement(),
							       componentSpaces[i]->getAdmin(),
							       &locInd);
	double estimate = 0.0;
	for (int j = 0; j < componentSpaces[i]->getBasisFcts()->getNumber(); j++) {
	  estimate += (*tmp)[locInd[j]];
	}
	elInfo->getElement()->setEstimation(estimate, i);
	elInfo->getElement()->setMark(0);
								
	elInfo = stack.traverseNext(elInfo);
      }  
      
      DELETE tmp;	
      DELETE sol;
    }						           
  }
}