Hello all, I'm trying to add a single DoF to my Trilinos system (not associated with the FE system, similar to the question asked here <https://groups.google.com/d/msg/dealii/FpDMRyhsmqI/XqiHH2qZpWEJ>). I followed the advice there, and set it up as a block system, where block (0,0) is the usual FE matrix, and block (1,1) is a 1x1 "matrix". Setting this up to work in parallel with Trilinos was tricky, but I got it working by using the usual partitioning IndexSet and an IndexSet with 1 entry owned by processor 0, which I use to reinit the new matrix blocks (the 1x1 block, and the 1xn and nx1 off-diagonal blocks).
The code I've written works, however, as I've scaled up the system it's become a huge bottleneck. Definitely this is because I'm using the deal.ii BlockSparsityPattern to reinitialize the matrix blocks in a not efficient way, but I haven't been able to find a better way to do it. I've attached the relevant piece of code below, if anyone could suggest a more functional way to do this, I'd appreciate it. // Initialize the block sparse matrix, where block (0,0) is the original matrix TrilinosWrappers::BlockSparseMatrix added_Jacobian_matrix; added_Jacobian_matrix.reinit(2,2); added_Jacobian_matrix.block(0,0).reinit(Jacobian_matrix); // Initialize the block sparsity pattern BlockSparsityPattern added_sp; added_sp.reinit(2,2); added_sp.block(0,0).reinit(Jacobian_matrix.m(), Jacobian_matrix.n(), 0); added_sp.block(1, 1).reinit(1,1,0); // This sets the off-diagonal blocks to the correct size FullMatrix<double> col(Jacobian_matrix.n(), 1), row(1, Jacobian_matrix.m()); for (unsigned int j=0; j<Jacobian_matrix.n(); ++j) { col(j,0) = 1; row(0,j) = 1; } added_sp.block(0, 1).copy_from(col); added_sp.block(1, 0).copy_from(row); added_sp.collect_sizes(); added_sp.compress(); // Set up a 1-entry IndexSet // IndexSet partitioning is what's used by the original matrix IndexSet single; single.clear(); single.set_size(1); if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0) single.add_index(0); // Reinitialize the matrix blocks using the deal.ii SparsityPattern added_Jacobian_matrix.block(1, 0).reinit( single, partitioning, added_sp.block(1, 0), MPI_COMM_WORLD, false); added_Jacobian_matrix.block(0, 1).reinit( partitioning, single, added_sp.block(0, 1), MPI_COMM_WORLD, false); added_Jacobian_matrix.block(1, 1).reinit( single, single, added_sp.block(1, 1), MPI_COMM_WORLD, false); -- The deal.II project is located at http://www.dealii.org/ For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en --- You received this message because you are subscribed to the Google Groups "deal.II User Group" group. To unsubscribe from this group and stop receiving emails from it, send an email to dealii+unsubscr...@googlegroups.com. For more options, visit https://groups.google.com/d/optout.