// // Software License for AMDiS // // Copyright (c) 2010 Dresden University of Technology // All rights reserved. // Authors: Simon Vey, Thomas Witkowski et al. // // This file is part of AMDiS // // See also license.opensource.txt in the distribution. #include "AMDiS.h" #include "parallel/PetscSolverGlobalMatrix.h" #include "parallel/StdMpi.h" #include "parallel/MpiHelper.h" namespace AMDiS { void PetscSolverGlobalMatrix::fillPetscMatrix(Matrix *mat) { FUNCNAME("PetscSolverGlobalMatrix::fillPetscMatrix()"); TEST_EXIT_DBG(meshDistributor)("No mesh distributor object defined!\n"); TEST_EXIT_DBG(mat)("No DOF matrix defined!\n"); double wtime = MPI::Wtime(); vector feSpaces = getFeSpaces(mat); int nRankRows = meshDistributor->getNumberRankDofs(feSpaces); int nOverallRows = meshDistributor->getNumberOverallDofs(feSpaces); // === Create PETSc vector (solution and a temporary vector). === VecCreateMPI(PETSC_COMM_WORLD, nRankRows, nOverallRows, &petscSolVec); VecCreateMPI(PETSC_COMM_WORLD, nRankRows, nOverallRows, &petscTmpVec); int testddd = 1; Parameters::get("block size", testddd); if (testddd > 1) { VecSetBlockSize(petscSolVec, testddd); VecSetBlockSize(petscTmpVec, testddd); } int recvAllValues = 0; int sendValue = static_cast(meshDistributor->getLastMeshChangeIndex() != lastMeshNnz); meshDistributor->getMpiComm().Allreduce(&sendValue, &recvAllValues, 1, MPI_INT, MPI_SUM); recvAllValues = 1; if (!d_nnz || recvAllValues != 0 || alwaysCreateNnzStructure) { // Global DOF mapping must only be recomputed, if the NNZ structure has // changed (thus, also the mesh has changed). createGlobalDofMapping(feSpaces); if (d_nnz) { delete [] d_nnz; d_nnz = NULL; delete [] o_nnz; o_nnz = NULL; } createPetscNnzStructure(mat); lastMeshNnz = meshDistributor->getLastMeshChangeIndex(); } // === Create PETSc matrix with the computed nnz data structure. === MatCreateMPIAIJ(PETSC_COMM_WORLD, nRankRows, nRankRows, nOverallRows, nOverallRows, 0, d_nnz, 0, o_nnz, &petscMatrix); if (testddd > 1) { MatSetBlockSize(petscMatrix, testddd); MSG("MAT SET BLOCK SIZE: %d\n", testddd); } #if (DEBUG != 0) MSG("Fill petsc matrix 1 needed %.5f seconds\n", MPI::Wtime() - wtime); #endif #if (DEBUG != 0) int a, b; MatGetOwnershipRange(petscMatrix, &a, &b); TEST_EXIT(a == meshDistributor->getStartDofs(feSpaces)) ("Wrong matrix ownership range!\n"); TEST_EXIT(b == meshDistributor->getStartDofs(feSpaces) + nRankRows) ("Wrong matrix ownership range!\n"); #endif // === Transfer values from DOF matrices to the PETSc matrix. === int nComponents = mat->getNumRows(); for (int i = 0; i < nComponents; i++) for (int j = 0; j < nComponents; j++) if ((*mat)[i][j]) setDofMatrix((*mat)[i][j], i, j); #if (DEBUG != 0) MSG("Fill petsc matrix 2 needed %.5f seconds\n", MPI::Wtime() - wtime); #endif MatAssemblyBegin(petscMatrix, MAT_FINAL_ASSEMBLY); MatAssemblyEnd(petscMatrix, MAT_FINAL_ASSEMBLY); // === Init PETSc solver. === KSPCreate(PETSC_COMM_WORLD, &solver); KSPGetPC(solver, &pc); KSPSetOperators(solver, petscMatrix, petscMatrix, SAME_NONZERO_PATTERN); KSPSetTolerances(solver, 0.0, 1e-8, PETSC_DEFAULT, PETSC_DEFAULT); KSPSetType(solver, KSPBCGS); KSPSetOptionsPrefix(solver, kspPrefix.c_str()); KSPSetFromOptions(solver); PCSetFromOptions(pc); // Do not delete the solution vector, use it for the initial guess. if (!zeroStartVector) KSPSetInitialGuessNonzero(solver, PETSC_TRUE); MSG("Fill petsc matrix needed %.5f seconds\n", MPI::Wtime() - wtime); } void PetscSolverGlobalMatrix::fillPetscRhs(SystemVector *vec) { FUNCNAME("PetscSolverGlobalMatrix::fillPetscRhs()"); TEST_EXIT_DBG(vec)("No DOF vector defined!\n"); TEST_EXIT_DBG(meshDistributor)("No mesh distributor defined!\n"); vector feSpaces = getFeSpaces(vec); int nRankRows = meshDistributor->getNumberRankDofs(feSpaces); int nOverallRows = meshDistributor->getNumberOverallDofs(feSpaces); VecCreateMPI(PETSC_COMM_WORLD, nRankRows, nOverallRows, &petscRhsVec); int testddd = 1; Parameters::get("block size", testddd); if (testddd > 1) VecSetBlockSize(petscRhsVec, testddd); // === Transfer values from DOF vector to the PETSc vector. === for (int i = 0; i < vec->getSize(); i++) setDofVector(petscRhsVec, vec->getDOFVector(i), i); VecAssemblyBegin(petscRhsVec); VecAssemblyEnd(petscRhsVec); if (removeRhsNullSpace) { MSG("Remove constant null space from the RHS!\n"); MatNullSpace sp; MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, PETSC_NULL, &sp); MatNullSpaceRemove(sp, petscRhsVec, PETSC_NULL); MatNullSpaceDestroy(&sp); } } void PetscSolverGlobalMatrix::solvePetscMatrix(SystemVector &vec, AdaptInfo *adaptInfo) { FUNCNAME("PetscSolverGlobalMatrix::solvePetscMatrix()"); int nComponents = vec.getSize(); // === Set old solution to be initiual guess for PETSc solver. === if (!zeroStartVector) { VecSet(petscSolVec, 0.0); for (int i = 0; i < nComponents; i++) setDofVector(petscSolVec, vec.getDOFVector(i), i, true); VecAssemblyBegin(petscSolVec); VecAssemblyEnd(petscSolVec); } // PETSc. PetscErrorCode solverError = KSPSolve(solver, petscRhsVec, petscSolVec); if (solverError != 0) { AMDiS::finalize(); exit(-1); } // === Transfere values from PETSc's solution vectors to the DOF vectors. === PetscScalar *vecPointer; VecGetArray(petscSolVec, &vecPointer); int c = 0; for (int i = 0; i < nComponents; i++) { DOFVector &dv = *(vec.getDOFVector(i)); const FiniteElemSpace *feSpace = dv.getFeSpace(); int nRankDofs = meshDistributor->getNumberRankDofs(feSpace); for (int j = 0; j < nRankDofs; j++) dv[meshDistributor->mapLocalToDof(feSpace, j)] = vecPointer[c++]; } VecRestoreArray(petscSolVec, &vecPointer); // === Synchronize DOFs at common DOFs, i.e., DOFs that correspond to === // === more than one partition. === meshDistributor->synchVector(vec); // Print iteration counter and residual norm of the solution. printSolutionInfo(adaptInfo); // === Destroy PETSc's variables. === VecDestroy(&petscRhsVec); } void PetscSolverGlobalMatrix::destroyMatrixData() { FUNCNAME("PetscSolverGlobalMatrix::destroyMatrixData()"); MatDestroy(&petscMatrix); KSPDestroy(&solver); VecDestroy(&petscSolVec); VecDestroy(&petscTmpVec); } void PetscSolverGlobalMatrix::setDofMatrix(DOFMatrix* mat, int nRowMat, int nColMat) { FUNCNAME("PetscSolverGlobalMatrix::setDofMatrix()"); TEST_EXIT(mat)("No DOFMatrix!\n"); using mtl::tag::row; using mtl::tag::nz; using mtl::begin; using mtl::end; namespace traits = mtl::traits; typedef DOFMatrix::base_matrix_type Matrix; traits::col::type col(mat->getBaseMatrix()); traits::const_value::type value(mat->getBaseMatrix()); typedef traits::range_generator::type cursor_type; typedef traits::range_generator::type icursor_type; vector cols; vector values; cols.reserve(300); values.reserve(300); vector globalCols; // Get periodic mapping object PeriodicMap &perMap = meshDistributor->getPeriodicMap(); // === Traverse all rows of the dof matrix and insert row wise the values === // === to the PETSc matrix. === for (cursor_type cursor = begin(mat->getBaseMatrix()), cend = end(mat->getBaseMatrix()); cursor != cend; ++cursor) { const FiniteElemSpace *rowFe = mat->getRowFeSpace(); const FiniteElemSpace *colFe = mat->getColFeSpace(); // Global index of the current row DOF. int globalRowDof = meshDistributor->mapDofToGlobal(rowFe, *cursor); // Test if the current row DOF is a periodic DOF. bool periodicRow = perMap.isPeriodic(rowFe, globalRowDof); if (!periodicRow) { // === Row DOF index is not periodic. === // Get PETSc's mat row index. int rowIndex = dofToMatIndex.get(nRowMat, globalRowDof); cols.clear(); values.clear(); for (icursor_type icursor = begin(cursor), icend = end(cursor); icursor != icend; ++icursor) { // Global index of the current column index. int globalColDof = meshDistributor->mapDofToGlobal(colFe, col(*icursor)); // Test if the current col dof is a periodic dof. bool periodicCol = perMap.isPeriodic(colFe, globalColDof); // Get PETSc's mat col index. int colIndex = dofToMatIndex.get(nColMat, globalColDof); // Ignore all zero entries, expect it is a diagonal entry. if (value(*icursor) == 0.0 && rowIndex != colIndex) continue; if (!periodicCol) { // Calculate the exact position of the column index in the PETSc matrix. cols.push_back(colIndex); values.push_back(value(*icursor)); } else { // === Row index is not periodic, but column index is. === // Create set of all periodic associations of the column index. std::set perAsc; std::set& perColAsc = perMap.getAssociations(colFe, globalColDof); for (std::set::iterator it = perColAsc.begin(); it != perColAsc.end(); ++it) if (*it >= -3) perAsc.insert(*it); // Scale value to the number of periodic associations of the column index. double scaledValue = value(*icursor) * pow(0.5, static_cast(perAsc.size())); // === Create set of all matrix column indices due to the periodic === // === associations of the column DOF index. === vector newCols; // First, add the original matrix index. newCols.push_back(globalColDof); // And add all periodic matrix indices. for (std::set::iterator it = perAsc.begin(); it != perAsc.end(); ++it) { int nCols = static_cast(newCols.size()); for (int i = 0; i < nCols; i++) { TEST_EXIT_DBG(perMap.isPeriodic(colFe, *it, newCols[i])) ("Wrong periodic DOF associations at boundary %d with DOF %d!\n", *it, newCols[i]); newCols.push_back(perMap.map(colFe, *it, newCols[i])); } } for (unsigned int i = 0; i < newCols.size(); i++) { cols.push_back(dofToMatIndex.get(nColMat, newCols[i])); values.push_back(scaledValue); } } } MatSetValues(petscMatrix, 1, &rowIndex, cols.size(), &(cols[0]), &(values[0]), ADD_VALUES); } else { // === Row DOF index is periodic. === // Because this row is periodic, we will have to add the entries of this // matrix row to multiple rows. The following maps store to each row an // array of column indices and values of the entries that must be added to // the PETSc matrix. map > colsMap; map > valsMap; // Traverse all column entries. for (icursor_type icursor = begin(cursor), icend = end(cursor); icursor != icend; ++icursor) { // Global index of the current column index. int globalColDof = meshDistributor->mapDofToGlobal(colFe, col(*icursor)); // Ignore all zero entries, expect it is a diagonal entry. if (value(*icursor) == 0.0 && globalRowDof != globalColDof) continue; // === Add all periodic associations of both, the row and the column === // === indices to the set perAsc. === std::set perAsc; if (perMap.isPeriodic(colFe, globalColDof)) { std::set& perColAsc = perMap.getAssociations(colFe, globalColDof); for (std::set::iterator it = perColAsc.begin(); it != perColAsc.end(); ++it) if (*it >= -3) perAsc.insert(*it); } std::set& perRowAsc = perMap.getAssociations(rowFe, globalRowDof); for (std::set::iterator it = perRowAsc.begin(); it != perRowAsc.end(); ++it) if (*it >= -3) perAsc.insert(*it); // Scale the value with respect to the number of periodic associations. double scaledValue = value(*icursor) * pow(0.5, static_cast(perAsc.size())); // === Create all matrix entries with respect to the periodic === // === associations of the row and column indices. === vector > entry; // First, add the original entry. entry.push_back(make_pair(globalRowDof, globalColDof)); // Then, traverse the periodic associations of the row and column // indices and create the corresponding entries. for (std::set::iterator it = perAsc.begin(); it != perAsc.end(); ++it) { int nEntry = static_cast(entry.size()); for (int i = 0; i < nEntry; i++) { int perRowDof = 0; if (perMap.isPeriodic(rowFe, *it, entry[i].first)) perRowDof = perMap.map(rowFe, *it, entry[i].first); else perRowDof = entry[i].first; int perColDof; if (perMap.isPeriodic(colFe, *it, entry[i].second)) perColDof = perMap.map(colFe, *it, entry[i].second); else perColDof = entry[i].second; entry.push_back(make_pair(perRowDof, perColDof)); } } // === Translate the matrix entries to PETSc's matrix. for (unsigned int i = 0; i < entry.size(); i++) { int rowIdx = dofToMatIndex.get(nRowMat, entry[i].first); int colIdx = dofToMatIndex.get(nColMat, entry[i].second); colsMap[rowIdx].push_back(colIdx); valsMap[rowIdx].push_back(scaledValue); } } // === Finally, add all periodic rows to the PETSc matrix. === for (map >::iterator rowIt = colsMap.begin(); rowIt != colsMap.end(); ++rowIt) { TEST_EXIT_DBG(rowIt->second.size() == valsMap[rowIt->first].size()) ("Should not happen!\n"); int rowIndex = rowIt->first; MatSetValues(petscMatrix, 1, &rowIndex, rowIt->second.size(), &(rowIt->second[0]), &(valsMap[rowIt->first][0]), ADD_VALUES); } } } } void PetscSolverGlobalMatrix::setDofVector(Vec& petscVec, DOFVector* vec, int nRowVec, bool rankOnly) { FUNCNAME("PetscSolverGlobalMatrix::setDofVector()"); const FiniteElemSpace *feSpace = vec->getFeSpace(); PeriodicMap &perMap = meshDistributor->getPeriodicMap(); // Traverse all used DOFs in the dof vector. DOFVector::Iterator dofIt(vec, USED_DOFS); for (dofIt.reset(); !dofIt.end(); ++dofIt) { if (rankOnly && !meshDistributor->getIsRankDof(feSpace, dofIt.getDOFIndex())) continue; // Calculate global row index of the DOF. DegreeOfFreedom globalRowDof = meshDistributor->mapDofToGlobal(feSpace, dofIt.getDOFIndex()); // Get PETSc's mat index of the row DOF. int index = dofToMatIndex.get(nRowVec, globalRowDof); if (perMap.isPeriodic(feSpace, globalRowDof)) { std::set& perAsc = perMap.getAssociations(feSpace, globalRowDof); double value = *dofIt / (perAsc.size() + 1.0); VecSetValues(petscVec, 1, &index, &value, ADD_VALUES); for (std::set::iterator perIt = perAsc.begin(); perIt != perAsc.end(); ++perIt) { int mappedDof = perMap.map(feSpace, *perIt, globalRowDof); int mappedIndex = dofToMatIndex.get(nRowVec, mappedDof); VecSetValues(petscVec, 1, &mappedIndex, &value, ADD_VALUES); } } else { // The DOF index is not periodic. double value = *dofIt; VecSetValues(petscVec, 1, &index, &value, ADD_VALUES); } } } void PetscSolverGlobalMatrix::createPetscNnzStructure(Matrix *mat) { FUNCNAME("PetscSolverGlobalMatrix::createPetscNnzStructure()"); TEST_EXIT_DBG(!d_nnz)("There is something wrong!\n"); TEST_EXIT_DBG(!o_nnz)("There is something wrong!\n"); vector feSpaces = getFeSpaces(mat); int nRankRows = meshDistributor->getNumberRankDofs(feSpaces); int rankStartIndex = meshDistributor->getStartDofs(feSpaces); d_nnz = new int[nRankRows]; o_nnz = new int[nRankRows]; for (int i = 0; i < nRankRows; i++) { d_nnz[i] = 0; o_nnz[i] = 0; } using mtl::tag::row; using mtl::tag::nz; using mtl::begin; using mtl::end; namespace traits = mtl::traits; typedef DOFMatrix::base_matrix_type Matrix; typedef vector > MatrixNnzEntry; typedef map RankToDofContainer; // Stores to each rank a list of nnz entries (i.e. pairs of row and column // index) that this rank will send to. These nnz entries will be assembled // on this rank, but because the row DOFs are not DOFs of this rank they // will be send to the owner of the row DOFs. map sendMatrixEntry; // Maps to each DOF that must be send to another rank the rank number of the // receiving rank. map, int> sendDofToRank; // First, create for all ranks, to which we send data to, MatrixNnzEntry // object with 0 entries. for (unsigned int i = 0; i < feSpaces.size(); i++) { for (DofComm::Iterator it(meshDistributor->getRecvDofs(), feSpaces[i]); !it.end(); it.nextRank()) { sendMatrixEntry[it.getRank()].resize(0); for (; !it.endDofIter(); it.nextDof()) sendDofToRank[make_pair(it.getDofIndex(), i)] = it.getRank(); } } // Create list of ranks from which we receive data from. std::set recvFromRank; for (unsigned int i = 0; i < feSpaces.size(); i++) for (DofComm::Iterator it(meshDistributor->getSendDofs(), feSpaces[i]); !it.end(); it.nextRank()) recvFromRank.insert(it.getRank()); // === Traverse matrices to create nnz data. === int nComponents = mat->getNumRows(); for (int i = 0; i < nComponents; i++) { for (int j = 0; j < nComponents; j++) { if (!(*mat)[i][j]) continue; TEST_EXIT_DBG((*mat)[i][j]->getRowFeSpace() == feSpaces[i]) ("Should not happen!\n"); TEST_EXIT_DBG((*mat)[i][j]->getColFeSpace() == feSpaces[j]) ("Should not happen!\n"); Matrix bmat = (*mat)[i][j]->getBaseMatrix(); traits::col::type col(bmat); traits::const_value::type value(bmat); typedef traits::range_generator::type cursor_type; typedef traits::range_generator::type icursor_type; for (cursor_type cursor = begin(bmat), cend = end(bmat); cursor != cend; ++cursor) { int globalRowDof = meshDistributor->mapDofToGlobal(feSpaces[i], *cursor); // The corresponding global matrix row index of the current row DOF. int petscRowIdx = dofToMatIndex.get(i, globalRowDof); if (meshDistributor->getIsRankDof(feSpaces[i], *cursor)) { // === The current row DOF is a rank DOF, so create the === // === corresponding nnz values directly on rank's nnz data. === // This is the local row index of the local PETSc matrix. int localPetscRowIdx = petscRowIdx - rankStartIndex; TEST_EXIT_DBG(localPetscRowIdx >= 0 && localPetscRowIdx < nRankRows) ("Should not happen! \n Debug info: localRowIdx = %d globalRowIndx = %d petscRowIdx = %d localPetscRowIdx = %d rStart = %d nCompontens = %d nRankRows = %d\n", *cursor, meshDistributor->mapDofToGlobal(feSpaces[i], *cursor), petscRowIdx, localPetscRowIdx, rankStartIndex, nComponents, nRankRows); // Traverse all non zero entries in this row. for (icursor_type icursor = begin(cursor), icend = end(cursor); icursor != icend; ++icursor) { int globalColDof = meshDistributor->mapDofToGlobal(feSpaces[j], col(*icursor)); int petscColIdx = dofToMatIndex.get(j, globalColDof); if (value(*icursor) != 0.0 || petscRowIdx == petscColIdx) { // The row DOF is a rank DOF, if also the column is a rank DOF, // increment the d_nnz values for this row, otherwise the // o_nnz value. if (petscColIdx >= rankStartIndex && petscColIdx < rankStartIndex + nRankRows) d_nnz[localPetscRowIdx]++; else o_nnz[localPetscRowIdx]++; } } } else { // === The current row DOF is not a rank DOF, i.e., its values === // === are also created on this rank, but afterthere they will === // === be send to another rank. So we need to send also the === // === corresponding nnz structure of this row to the corres- === // === ponding rank. === // Send all non zero entries to the member of the row DOF. int sendToRank = sendDofToRank[make_pair(*cursor, i)]; for (icursor_type icursor = begin(cursor), icend = end(cursor); icursor != icend; ++icursor) { if (value(*icursor) != 0.0) { int globalColDof = meshDistributor->mapDofToGlobal(feSpaces[j], col(*icursor)); int petscColIdx = dofToMatIndex.get(j, globalColDof); sendMatrixEntry[sendToRank]. push_back(make_pair(petscRowIdx, petscColIdx)); } } } // if (isRankDof[*cursor]) ... else ... } // for each row in mat[i][j] } } // === Send and recv the nnz row structure to/from other ranks. === StdMpi stdMpi(meshDistributor->getMpiComm(), true); stdMpi.send(sendMatrixEntry); for (std::set::iterator it = recvFromRank.begin(); it != recvFromRank.end(); ++it) stdMpi.recv(*it); stdMpi.startCommunication(); // === Evaluate the nnz structure this rank got from other ranks and add === // === it to the PETSc nnz data structure. === for (map::iterator it = stdMpi.getRecvData().begin(); it != stdMpi.getRecvData().end(); ++it) { if (it->second.size() > 0) { for (unsigned int i = 0; i < it->second.size(); i++) { int r = it->second[i].first; int c = it->second[i].second; int localRowIdx = r - rankStartIndex; TEST_EXIT_DBG(localRowIdx >= 0 && localRowIdx < nRankRows) ("Got row index %d/%d (nRankRows = %d) from rank %d. Should not happen!\n", r, localRowIdx, nRankRows, it->first); if (c < rankStartIndex || c >= rankStartIndex + nRankRows) o_nnz[localRowIdx]++; else d_nnz[localRowIdx]++; } } } // The above algorithm for calculating the number of nnz per row over- // approximates the value, i.e., the number is always equal or larger to // the real number of nnz values in the global parallel matrix. For small // matrices, the problem may arise, that the result is larger than the // number of elements in a row. This is fixed in the following. if (nRankRows < 100) for (int i = 0; i < nRankRows; i++) d_nnz[i] = std::min(d_nnz[i], nRankRows); } void PetscSolverGlobalMatrix::createGlobalDofMapping(vector &feSpaces) { FUNCNAME("PetscSolverGlobalMatrix::createGlobalDofMapping()"); int offset = meshDistributor->getStartDofs(feSpaces); Mesh *mesh = meshDistributor->getMesh(); dofToMatIndex.clear(); for (unsigned int i = 0; i < feSpaces.size(); i++) { // === Create indices for all DOFs in rank' domain. === std::set rankDofSet; mesh->getAllDofs(feSpaces[i], rankDofSet); for (std::set::iterator it = rankDofSet.begin(); it != rankDofSet.end(); ++it) if (meshDistributor->getIsRankDof(feSpaces[i], **it)) { int globalIndex = meshDistributor->mapDofToGlobal(feSpaces[i], **it); int globalMatIndex = globalIndex - meshDistributor->getStartDofs(feSpaces[i]) + offset; dofToMatIndex.add(i, globalIndex, globalMatIndex); } // === Communicate interior boundary DOFs between domains. === StdMpi > stdMpi(meshDistributor->getMpiComm()); for (DofComm::Iterator it(meshDistributor->getSendDofs(), feSpaces[i]); !it.end(); it.nextRank()) { vector sendGlobalDofs; for (; !it.endDofIter(); it.nextDof()) { int globalIndex = meshDistributor->mapDofToGlobal(feSpaces[i], it.getDofIndex()); int globalMatIndex = dofToMatIndex.get(i, globalIndex); sendGlobalDofs.push_back(globalMatIndex); } stdMpi.send(it.getRank(), sendGlobalDofs); } for (DofComm::Iterator it(meshDistributor->getRecvDofs(), feSpaces[i]); !it.end(); it.nextRank()) stdMpi.recv(it.getRank()); stdMpi.startCommunication(); for (DofComm::Iterator it(meshDistributor->getRecvDofs(), feSpaces[i]); !it.end(); it.nextRank()) for (; !it.endDofIter(); it.nextDof()) { int globalIndex = meshDistributor->mapDofToGlobal(feSpaces[i], it.getDofIndex()); int globalMatIndex = stdMpi.getRecvData(it.getRank())[it.getDofCounter()]; dofToMatIndex.add(i, globalIndex, globalMatIndex); } // === Communicate DOFs on periodic boundaries. === stdMpi.clear(); for (DofComm::Iterator it(meshDistributor->getPeriodicDofs(), feSpaces[i]); !it.end(); it.nextRank()) { vector sendGlobalDofs; for (; !it.endDofIter(); it.nextDof()) { int ind0 = meshDistributor->mapDofToGlobal(feSpaces[i], it.getDofIndex()); int ind1 = dofToMatIndex.get(i, ind0); sendGlobalDofs.push_back(ind0); sendGlobalDofs.push_back(ind1); } stdMpi.send(it.getRank(), sendGlobalDofs); stdMpi.recv(it.getRank()); } stdMpi.startCommunication(); for (DofComm::Iterator it(meshDistributor->getPeriodicDofs(), feSpaces[i]); !it.end(); it.nextRank()) for (; !it.endDofIter(); it.nextDof()) { int ind0 = stdMpi.getRecvData(it.getRank())[it.getDofCounter() * 2]; int ind1 = stdMpi.getRecvData(it.getRank())[it.getDofCounter() * 2 + 1]; dofToMatIndex.add(i, ind0, ind1); } // === Update offset. === offset += meshDistributor->getNumberRankDofs(feSpaces[i]); } } }