#include #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 "DirichletBC.h" #include "RobinBC.h" #include "PeriodicBC.h" #include "Lagrange.h" #include "Flag.h" #include "TraverseParallel.h" #include "VtkWriter.h" #include "ValueReader.h" namespace AMDiS { 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 the adopt problem has fewer components than this problem, but only one // mesh for all component, than scal up the componentMeshes array. if (adoptProblem->getNumComponents() < nComponents) { TEST_EXIT(meshes.size() == 1)("Daran muss ich noch arbeiten!\n"); componentMeshes.resize(nComponents); for (int i = adoptProblem->getNumComponents(); i < nComponents; i++) { componentMeshes[i] = componentMeshes[0]; } } } } 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; traverseInfo = adoptProblem->traverseInfo; // If the adopt problem has fewer components than this problem, but only one // fe space for all component, than scal up the componentSpaces array. if (adoptProblem->getNumComponents() < nComponents) { TEST_EXIT(feSpaces.size() == 1)("Daran muss ich noch arbeiten!\n"); componentSpaces.resize(nComponents); for (int i = adoptProblem->getNumComponents(); i < nComponents; i++) { componentSpaces[i] = componentSpaces[0]; } } } } 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(meshes.size()); i++) if (initFlag.isSet(INIT_MESH) && meshes[i] && !(meshes[i]->isInitialized())) meshes[i]->initialize(); // === read value file and use it for the mesh values === std::string valueFilename(""); GET_PARAMETER(0, meshes[0]->getName() + "->value file name", &valueFilename); if (valueFilename.length()) { ValueReader::readValue(valueFilename, meshes[0], solution->getDOFVector(0), meshes[0]->getMacroFileInfo()); meshes[0]->clearMacroFileInfo(); } // === do global refinements === for (int i = 0; i < static_cast(meshes.size()); i++) if (initFlag.isSet(INIT_MESH) && meshes[i]) refinementManager->globalRefine(meshes[i], globalRefinements); } } doOtherStuff(); } void ProblemVec::createMesh() { FUNCNAME("ProblemVec::createMesh()"); componentMeshes.resize(nComponents); std::map 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); nMeshes++; } 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()"); std::map< std::pair, 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); traverseInfo.resize(nComponents); for (int i = 0; i < nComponents; i++) { char number[3]; sprintf(number, "%d", i); int degree = 1; GET_PARAMETER(0, name + "->polynomial degree[" + number + "]","%d", °ree); TEST_EXIT(componentSpaces[i] == NULL)("feSpace already created\n"); if (feSpaceMap[std::pair(componentMeshes[i], degree)] == NULL) { stringstream s; s << name << "->feSpace[" << i << "]"; FiniteElemSpace *newFESpace = FiniteElemSpace::provideFESpace(NULL, Lagrange::getLagrange(dim, degree), componentMeshes[i], s.str()); feSpaceMap[std::pair(componentMeshes[i], degree)] = newFESpace; feSpaces.push_back(newFESpace); } componentSpaces[i] = feSpaceMap[std::pair(componentMeshes[i], degree)]; } for (int i = 0; i < nComponents; i++) { for (int j = 0; j < nComponents; j++) { traverseInfo.getMatrix(i, j).setFESpace(componentSpaces[i], componentSpaces[j]); } traverseInfo.getVector(i).setFESpace(componentSpaces[i]); } // create dof admin for vertex dofs if neccessary for (int i = 0; i < static_cast(meshes.size()); i++) { if (meshes[i]->getNumberOfDOFs(VERTEX) == 0) { DimVec 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(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(componentSpaces[i], numberedName)); numberedName = name + std::string(number); solution->setDOFVector(i, NEW DOFVector(componentSpaces[i], numberedName)); solution->getDOFVector(i)->setCoarsenOperation(COARSE_INTERPOL); solution->getDOFVector(i)->set(0.0); } } void ProblemVec::createSolver() { FUNCNAME("ProblemVec::createSolver()"); // === create solver === std::string solverType("0"); GET_PARAMETER(0, name + "->solver", &solverType); OEMSolverCreator *solverCreator = dynamic_cast(CreatorMap::getCreator(solverType)); TEST_EXIT(solverCreator)("no solver type\n"); solverCreator->setName(name + "->solver"); solver = solverCreator->create(); solver->initParameters(); } void ProblemVec::createEstimator() { FUNCNAME("ProblemVec::createEstimator()"); // create and set leaf data prototype for (int i = 0; i < static_cast(meshes.size()); i++) { meshes[i]->setElementDataPrototype (NEW LeafDataEstimatableVec(NEW LeafDataCoarsenableVec)); } char number[3]; std::string estName; for (int 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("0"); GET_PARAMETER(0, estName, &estimatorType); EstimatorCreator *estimatorCreator = dynamic_cast( CreatorMap::getCreator(estimatorType)); if (estimatorCreator) { estimatorCreator->setName(estName); estimatorCreator->setRow(i); if (estimatorType == "recovery") { dynamic_cast(estimatorCreator)-> setSolution(solution->getDOFVector(i)); } estimator[i] = estimatorCreator->create(); } if (estimator[i]) { for (int 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 nMarkersCreated = 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]) { nMarkersCreated++; // If there is more than one marker, and all components are defined // on the same mesh, the maximum marking has to be enabled. if ((nMarkersCreated > 1) && (nMeshes == 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* > 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(this)); } } void ProblemVec::doOtherStuff() { } void ProblemVec::solve(AdaptInfo *adaptInfo, bool fixedMatrix) { 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->solveSystem(solverMatrix, *solution, *rhs); #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()"); VtkWriter::writeFile(solution->getDOFVector(0), "test.vtu"); 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(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(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, bool asmMatrix, bool asmVector) { FUNCNAME("ProblemVec::buildAfterCoarsen()"); clock_t first = clock(); #ifdef _OPENMP double wtime = omp_get_wtime(); #endif for (int i = 0; i < static_cast(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]) { DOFMatrix* dof_matrix= (*systemMatrix)[i][j]; DOFMatrix::base_matrix_type& base_matrix= dof_matrix->getBaseMatrix(); int nnz_per_row= 0; if (num_rows(base_matrix) != 0) nnz_per_row= int(double(base_matrix.nnz()) / num_rows(base_matrix) * 1.2); if (nnz_per_row < 5) nnz_per_row= 5; // Correct dimensionality of matrix base_matrix.change_dim(componentSpaces[i]->getAdmin()->getUsedSize(), componentSpaces[j]->getAdmin()->getUsedSize()); // Reuse old sparsity information (if available) or default // dof_matrix->startInsertion(nnz_per_row); } } } traverseInfo.updateStatus(); 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) 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 (traverseInfo.getStatus(i, j) == SingleComponentInfo::EQ_SPACES_NO_AUX || traverseInfo.getStatus(i, j) == SingleComponentInfo::EQ_SPACES_WITH_AUX) { // Row fe space and col fe space are both equal if (traverseInfo.getStatus(i) == SingleComponentInfo::EQ_SPACES_NO_AUX || traverseInfo.getStatus(i) == SingleComponentInfo::EQ_SPACES_WITH_AUX) { // The simplest case: either the right hand side has no operaters, no aux // fe spaces, or all aux fe spaces are equal to the row and col fe space. assembleOnOneMesh(componentSpaces[i], assembleFlag, assembleMatrix ? matrix : NULL, ((i == j) && asmVector) ? rhs->getDOFVector(i) : NULL); } else if (traverseInfo.getStatus(i) == SingleComponentInfo::EQ_SPACES_WITH_DIF_AUX) { // Row fe space and col fe space are both equal, but right hand side has at // least one another aux fe space. assembleOnOneMesh(componentSpaces[i], assembleFlag, assembleMatrix ? matrix : NULL, ((i == j) && asmVector) ? rhs->getDOFVector(i) : NULL); assembleOnDifMeshes2(componentSpaces[i], traverseInfo.getAuxFESpace(i, j), assembleFlag, NULL, ((i == j) && asmVector) ? rhs->getDOFVector(i) : NULL); } else { ERROR_EXIT("Possible? If yes, not yet implemented!\n"); } } else if (traverseInfo.getStatus(i, j) == SingleComponentInfo::EQ_SPACES_WITH_DIF_AUX) { assembleOnOneMesh(componentSpaces[i], assembleFlag, assembleMatrix ? matrix : NULL, ((i == j) && asmVector) ? rhs->getDOFVector(i) : NULL); assembleOnDifMeshes2(componentSpaces[i], traverseInfo.getAuxFESpace(i, j), assembleFlag, assembleMatrix ? matrix : NULL, ((i == j) && asmVector) ? rhs->getDOFVector(i) : NULL); } else if (traverseInfo.getStatus(i, j) == SingleComponentInfo::DIF_SPACES_NO_AUX || traverseInfo.getStatus(i, j) == SingleComponentInfo::DIF_SPACES_WITH_AUX) { assembleOnDifMeshes(componentSpaces[i], componentSpaces[j], assembleFlag, assembleMatrix ? matrix : NULL, ((i == j) && asmVector) ? rhs->getDOFVector(i) : NULL); } else { ERROR_EXIT("Not yet implemented!\n"); } // TODO: ExitMatrix should be called after finishInsertion! if (assembleMatrix && matrix->getBoundaryManager()) matrix->getBoundaryManager()->exitMatrix(matrix); assembledMatrix[i][j] = true; } // And now assemble boundary conditions on the vectors assembleBoundaryConditions(rhs->getDOFVector(i), solution->getDOFVector(i), componentMeshes[i], assembleFlag); } int nnz = 0; // Finish insertion for (int i = 0; i < nComponents; i++) for (int j = 0; j < nComponents; j++) if ((*systemMatrix)[i][j]) { (*systemMatrix)[i][j]->finishInsertion(); nnz += (*systemMatrix)[i][j]->getBaseMatrix().nnz(); } // clock_t first1 = clock(); solverMatrix.setMatrix(*systemMatrix); // clock_t first2 = clock(); createPrecon(); // clock_t first3 = clock(); // std::cout << "T1 = " << TIME_USED(first1, first2) << std::endl; // std::cout << "T2 = " << TIME_USED(first2, first3) << std::endl; INFO(info, 8)("fillin of assembled matrix: %d\n", nnz); #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::createPrecon() { std::string preconType("no"); GET_PARAMETER(0, name + "->solver->left precon", &preconType); CreatorInterface *preconCreator = CreatorMap::getCreator(preconType); solver->setLeftPrecon( preconCreator->create(solverMatrix.getMatrix()) ); preconType= "no"; GET_PARAMETER(0, name + "->solver->right precon", &preconType); preconCreator = CreatorMap::getCreator(preconType); solver->setRightPrecon( preconCreator->create(solverMatrix.getMatrix()) ); } 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(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::interpolInitialSolution(std::vector >*> *fct) { FUNCNAME("ProblemVec::interpolInitialSolution()"); solution->interpol(fct); } void ProblemVec::addMatrixOperator(Operator *op, int i, int j, double *factor, double *estFactor, bool fake) { 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); if (!fake) { traverseInfo.getMatrix(i, j).setAuxFESpaces(op->getAuxFESpaces()); for (int k = 0; k < static_cast(op->getAuxFESpaces().size()); k++) { if ((op->getAuxFESpaces())[k]->getMesh() != componentSpaces[i]->getMesh() || (op->getAuxFESpaces())[k]->getMesh() != componentSpaces[j]->getMesh()) { op->setNeedDualTraverse(true); break; } } } } void ProblemVec::addVectorOperator(Operator *op, int i, double *factor, double *estFactor, bool fake) { FUNCNAME("ProblemVec::addVectorOperator()"); rhs->getDOFVector(i)->addOperator(op, factor, estFactor); if (!fake) { traverseInfo.getVector(i).setAuxFESpaces(op->getAuxFESpaces()); for (int j = 0; j < static_cast(op->getAuxFESpaces().size()); j++) { if ((op->getAuxFESpaces())[j]->getMesh() != componentSpaces[i]->getMesh()) { op->setNeedDualTraverse(true); break; } } } } void ProblemVec::addDirichletBC(BoundaryType type, int system, AbstractFunction >* 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 > *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 > *n, AbstractFunction > *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 *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 *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(); tmpMatrix->getBaseMatrix().change_dim(matrix->getRowFESpace()->getAdmin()->getUsedSize(), matrix->getColFESpace()->getAdmin()->getUsedSize()); tmpMatrix->startInsertion(); } if (vector) { tmpVector = NEW DOFVector(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); // After creating privat copies of the DOFMatrix and the DOFVector, all threads // have to wait at this barrier. Especially for small problems this is required, // because otherwise one thread may be finished with assembling, before another // has made his private copy. #ifdef _OPENMP #pragma omp barrier #endif 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, NULL); elInfo = stack.traverseNext(elInfo); } tmpMatrix->finishInsertion(); // 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 matrix->getBaseMatrix() += tmpMatrix->getBaseMatrix(); } #ifdef _OPENMP #pragma omp barrier #endif #ifdef _OPENMP #pragma omp master #endif { if (matrix) matrix->startInsertion(); } if (matrix) { // Remove rows corresponding to DOFs on a Dirichlet boundary. #ifdef _OPENMP #pragma omp critical #endif 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 *vector) { const BasisFunction *basisFcts = rowFeSpace->getBasisFcts(); BoundaryType *bound = NULL; if (useGetBound) bound = GET_MEMORY(BoundaryType, basisFcts->getNumber()); if (matrix) matrix->startInsertion(); DualTraverse dualTraverse; ElInfo *rowElInfo, *colElInfo; ElInfo *largeElInfo, *smallElInfo; dualTraverse.setFillSubElemMat(true, basisFcts); bool cont = dualTraverse.traverseFirst(rowFeSpace->getMesh(), colFeSpace->getMesh(), -1, -1, assembleFlag, assembleFlag, &rowElInfo, &colElInfo, &smallElInfo, &largeElInfo); while (cont) { if (useGetBound) basisFcts->getBound(rowElInfo, bound); if (matrix) { matrix->assemble(1.0, rowElInfo, colElInfo, smallElInfo, largeElInfo, 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::assembleOnDifMeshes2(const FiniteElemSpace *mainFeSpace, const FiniteElemSpace *auxFeSpace, Flag assembleFlag, DOFMatrix *matrix, DOFVector *vector) { Mesh *mainMesh = mainFeSpace->getMesh(); Mesh *auxMesh = auxFeSpace->getMesh(); const BasisFunction *basisFcts = mainFeSpace->getBasisFcts(); BoundaryType *bound = NULL; if (useGetBound) bound = GET_MEMORY(BoundaryType, basisFcts->getNumber()); if (matrix) matrix->startInsertion(); DualTraverse dualTraverse; ElInfo *mainElInfo, *auxElInfo; ElInfo *largeElInfo, *smallElInfo; dualTraverse.setFillSubElemMat(true, basisFcts); bool cont = dualTraverse.traverseFirst(mainMesh, auxMesh, -1, -1, assembleFlag, assembleFlag, &mainElInfo, &auxElInfo, &smallElInfo, &largeElInfo); while (cont) { if (useGetBound) basisFcts->getBound(mainElInfo, bound); if (matrix) matrix->assemble2(1.0, mainElInfo, auxElInfo, smallElInfo, largeElInfo, bound); if (vector) vector->assemble2(1.0, mainElInfo, auxElInfo, smallElInfo, largeElInfo, bound); cont = dualTraverse.traverseNext(&mainElInfo, &auxElInfo, &smallElInfo, &largeElInfo); } if (useGetBound) FREE_MEMORY(bound, BoundaryType, basisFcts->getNumber()); } void ProblemVec::assembleBoundaryConditions(DOFVector *rhs, DOFVector *solution, Mesh *mesh, Flag assembleFlag) { /* ================ Initialization of vectors ==================== */ if (rhs->getBoundaryManager()) rhs->getBoundaryManager()->initVector(rhs); if (solution->getBoundaryManager()) solution->getBoundaryManager()->initVector(solution); #ifdef _OPENMP TraverseParallelStack stack; #else TraverseStack stack; #endif /* ================= Parallel Boundary Assemblage ================= */ #ifdef _OPENMP #pragma omp parallel #endif { // Each thread assembles on its own dof-vectors. DOFVector *tmpRhsVec = NEW DOFVector(rhs->getFESpace(), "tmpRhs"); DOFVector *tmpSolVec = NEW DOFVector(solution->getFESpace(), "tmpSol"); tmpRhsVec->set(0.0); tmpSolVec->set(0.0); // (Parallel) traverse of mesh. ElInfo *elInfo = stack.traverseFirst(mesh, -1, assembleFlag); while (elInfo) { if (rhs->getBoundaryManager()) rhs->getBoundaryManager()-> fillBoundaryConditions(elInfo, tmpRhsVec); if (solution->getBoundaryManager()) solution->getBoundaryManager()->fillBoundaryConditions(elInfo, tmpSolVec); elInfo = stack.traverseNext(elInfo); } // After (parallel) mesh traverse, the result is applied to the final // vectors. This section is not allowed to be executed by more than one // thread at the same time. #ifdef _OPENMP #pragma omp critical #endif { DOFVector::Iterator rhsIt(rhs, USED_DOFS); DOFVector::Iterator solIt(solution, USED_DOFS); DOFVector::Iterator tmpRhsIt(tmpRhsVec, USED_DOFS); DOFVector::Iterator tmpSolIt(tmpSolVec, USED_DOFS); for (rhsIt.reset(), solIt.reset(), tmpRhsIt.reset(), tmpSolIt.reset(); !rhsIt.end(); ++rhsIt, ++solIt, ++tmpRhsIt, ++tmpSolIt) { if (*tmpRhsIt != 0.0) *rhsIt = *tmpRhsIt; if (*tmpSolIt != 0.0) *solIt = *tmpSolIt; } } // pragma omp critical DELETE tmpRhsVec; DELETE tmpSolVec; } // pragma omp parallel /* ======================= Finalize vectors ================== */ if (rhs->getBoundaryManager()) rhs->getBoundaryManager()->exitVector(rhs); if (solution->getBoundaryManager()) solution->getBoundaryManager()->exitVector(solution); } void ProblemVec::writeResidualMesh(int comp, AdaptInfo *adaptInfo, const std::string name) { FUNCNAME("ProblemVec::writeResidualMesh()"); std::map vec; TraverseStack stack; ElInfo *elInfo = stack.traverseFirst(this->getMesh(comp), -1, Mesh::CALL_LEAF_EL | Mesh::FILL_COORDS); while (elInfo) { vec[elInfo->getElement()->getIndex()] = elInfo->getElement()->getEstimation(comp); elInfo = stack.traverseNext(elInfo); } ElementFileWriter fw(name, this->getFESpace(comp), 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(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(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 *tmp = NEW DOFVector(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 *sol = NEW DOFVector(componentSpaces[i], "tmp"); sol->interpol(exactSolutionFcts[i]); DOFVector::Iterator it1(sol, USED_DOFS); DOFVector::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 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; } } }