Hello -

Does anyone have any experience with getting SLEPc to work efficiently in 
parallel? For some reason, SLEPc in parallel takes roughly 140 times more 
time to find one eigenvalue and eigenvector than the same solution in 
serial (i.e. the same matrices). For example, I have a problem with 4,905 
degrees of freedom. It takes the serial solve (assembly and solution of one 
eigenpair) 0.16667 seconds, while the parallel solution (2 processors) 
takes > 24 seconds. For larger problems I haven't even decided to wait on 
the solution of the parallel problem to see how long it would actually take 
(the longest I waited was 45 minutes for a problem that took 1.2 minutes in 
serial). 

For the parallel solution I am using the AMG preconditioner but maybe this 
is a poor choice for this problem? I'm really not sure. The eigenvalue 
problem I am solving is:
K x = w M x , where K is the elastic stiffness matrix, M is the mass 
matrix, w is an eigen-frequency of the system, and x is the eigenvector. 
Typically, the shift - invert method works very well (as it certainly does 
in serial) but I can't quite figure out why SLEPc is so helplessly slow in 
parallel. If anyone has any experience with this and would be willing to 
give me some advice, I would really appreciate it!

Here is the parallel version of the "solve" function.

unsigned int EigenvalueProblem<dim>::solve ()
  {
    TimerOutput::Scope t(computing_timer, "solve");
    pcout << "   Number of eigenvectors requested: "
              << eigenvectors.size() << std::endl;

    PETScWrappers::PreconditionBoomerAMG::AdditionalData data;
    data.symmetric_operator = true;
    PETScWrappers::PreconditionBoomerAMG preconditioner(mpi_communicator, 
data);
    SolverControl linear_solver_control(dof_handler.n_dofs(),1.0e-12,false,
false);
    PETScWrappers::SolverCG linear_solver(linear_solver_control,
mpi_communicator);
    linear_solver.initialize(preconditioner);

    SolverControl solver_control (2000, 1e-9, false, false);
    SLEPcWrappers::SolverLanczos eigensolver(solver_control, 
mpi_communicator);

    double shift_freq = parameters.get_double ("Shift frequency");
    double eigen_shift = std::pow( 2.0 * PI * shift_freq, 2.0);
    SLEPcWrappers::TransformationShiftInvert::AdditionalData additional_data
(eigen_shift);
    SLEPcWrappers::TransformationShiftInvert shift (mpi_communicator, 
additional_data);
    shift.set_solver(linear_solver);
    eigensolver.set_transformation (shift);
    eigensolver.set_which_eigenpairs (EPS_SMALLEST_REAL);
    eigensolver.set_problem_type (EPS_GHEP);

    eigensolver.solve (stiffness_matrix,mass_matrix,eigenvalues,eigenvectors
,
                       eigenvectors.size());

    for (unsigned int i=0; i<eigenvectors.size(); ++i)
        eigenvectors[i] /= eigenvectors[i].linfty_norm ();

    // Finally return the number of iterations it took to converge:
    return solver_control.last_step ();
  }

 And here is the serial version of the solve function:
unsigned int EigenvalueProblem<dim>::solve ()
  {
      
    std::cout << "   Number of eigenvectors requested: "
              << eigenvectors.size() << std::endl;

    SolverControl solver_control (2000, 1e-9);
    SLEPcWrappers::SolverKrylovSchur eigensolver(solver_control);
    eigensolver.set_which_eigenpairs (EPS_SMALLEST_REAL);
    eigensolver.set_problem_type (EPS_GHEP);
    
    double eigen_shift = std::pow( 2.0 * PI * parameters.get_double ("Shift 
frequency"), 2.0);
    SLEPcWrappers::TransformationShiftInvert::AdditionalData additional_data
(eigen_shift);
    SLEPcWrappers::TransformationShiftInvert shift (MPI_COMM_WORLD, 
additional_data);
    eigensolver.set_transformation (shift);
      
    eigensolver.solve (stiffness_matrix,mass_matrix,eigenvalues,eigenvectors
,
                       eigenvectors.size());

    for (unsigned int i=0; i<eigenvectors.size(); ++i)
      eigenvectors[i] /= eigenvectors[i].linfty_norm ();

    // Finally return the number of iterations it took to converge:
    return solver_control.last_step ();
  }

Thank you again for taking a look at my problem,
Jonathan

-- 
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 [email protected].
For more options, visit https://groups.google.com/d/optout.

Reply via email to