Hello Lailai

Here is an example from my code where I use MUMPS via Trilinos

      static TrilinosWrappers::SolverDirect::AdditionalData data
(false, "Amesos_Mumps");      static SolverControl solver_control (1,
0);      // solve for x      {         TrilinosWrappers::MPI::Vector
tmp (locally_owned_dofs, mpi_communicator);
TrilinosWrappers::SolverDirect solver (solver_control, data);
solver.solve (system_matrix, tmp, rhs);         constraints.distribute
(tmp);         x = tmp;      }


If you want to do the LU decomposition only once and reuse it, e.g., in
case of mass matrix, you can do something like this

      static TrilinosWrappers::SolverDirect::AdditionalData data
(false, "Amesos_Mumps");      static SolverControl solver_control (1,
0);

      // put following in your class, not sure that you really need a shared_ptr

std_cxx11::shared_ptr<TrilinosWrappers::SolverDirect> mumps_solver;
// If it is first time, compute LU decomposition if(!mumps_solver) {
mumps_solver = std_cxx11::shared_ptr<TrilinosWrappers::SolverDirect> (new
TrilinosWrappers::SolverDirect(solver_control, data)); mumps_solver->
initialize (mass_matrix); } // solve for ax { TrilinosWrappers::MPI::Vector
tmp (locally_owned_dofs, mpi_communicator); mumps_solver->solve (tmp, rhs);
constraints.distribute (tmp); ax = tmp; }


Best
praveen

On Sat, Feb 11, 2017 at 9:19 AM, Lailai Zhu <lailaizh...@gmail.com> wrote:

> hi, dear all,
>
> I am looking for an dealII example that uses trilino's  amesos direct
> solver with MPI,
> i guess more specifically, i mean how to call amesos_mumps solver based on
> MPI. I have installed
> them correctly with dealII. I think i need to use the subroutine
>
> void solve
> <https://dealii.org/8.4.1/doxygen/deal.II/classTrilinosWrappers_1_1SolverDirect.html#a6ce3a3768442425f3b32997b58600eb7>
> (const SparseMatrix
> <https://dealii.org/8.4.1/doxygen/deal.II/classTrilinosWrappers_1_1SparseMatrix.html>
> &A,::parallel::distributed::Vector
> <https://dealii.org/8.4.1/doxygen/deal.II/classparallel_1_1distributed_1_1Vector.html><
> double > &x, const ::parallel::distributed::Vector
> <https://dealii.org/8.4.1/doxygen/deal.II/classparallel_1_1distributed_1_1Vector.html><
> double > &b)
>
> May I ask is there any direct example or something close? Thanks in
> advance,
>
> best,
> lailai
>
> --
> 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.
>

-- 
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.

Reply via email to