Hello,
I have the feeling, that my setup is the root cause of poor performance
of both solve and assemble (3s for assemble in release mode for 2000 dofs
is according to my opinion slow [17s in debug mode]).


Could You look at following stripped snippets of setup, solve and assemble, and
give me advice concerning effective desing?

template<int dim>
void CahnHilliard<dim>::setup_systems() {

        dof_handler_phase.distribute_dofs(fe_phase_pot);

        locally_owned_dofs_phase = dof_handler_phase.locally_owned_dofs();
        DoFTools::extract_locally_relevant_dofs(dof_handler_phase,
                locally_relevant_dofs_phase);

        locally_relevant_solution_phase.reinit(locally_owned_dofs_phase,
                locally_relevant_dofs_phase, mpi_communicator);

        old_solution_phase.reinit(locally_owned_dofs_phase,
                locally_relevant_dofs_phase, mpi_communicator);
        old_solution_phase = 0;

        current_solution_phase.reinit(locally_owned_dofs_phase,
                locally_relevant_dofs_phase, mpi_communicator);
        current_solution_phase = 0;

        solution_update.reinit(locally_owned_dofs_phase,
                locally_relevant_dofs_phase, mpi_communicator);
        solution_update = 0;

        solution_phase.reinit(locally_owned_dofs_phase,
                locally_relevant_dofs_phase, mpi_communicator);

        system_rhs_phase.reinit(locally_owned_dofs_phase, mpi_communicator);
        system_rhs_phase = 0;

        residuum_phase.reinit(locally_owned_dofs_phase, mpi_communicator);
        residuum_phase = 0;

       constraint_matrix_phase.clear();
        constraint_matrix_phase.reinit(locally_relevant_dofs_phase);
        DoFTools::make_hanging_node_constraints(dof_handler_phase,
                constraint_matrix_phase);
        constraint_matrix_phase.close();

        constraint_matrix_phase_constant.clear();
        constraint_matrix_phase_constant.reinit(locally_relevant_dofs_phase);
        DoFTools::make_hanging_node_constraints(dof_handler_phase,
                constraint_matrix_phase_constant);
        constraint_matrix_phase_constant.close();

        LA::Vector<double> vec_old_solution(dof_handler_phase.n_dofs());
        VectorTools::interpolate(dof_handler_phase, InitialValuePhase<dim>(),
                vec_old_solution);

        old_solution_phase = vec_old_solution;

        DynamicSparsityPattern csp(dof_handler_phase.n_dofs(),
                dof_handler_phase.n_dofs());

        DoFTools::make_sparsity_pattern(dof_handler_phase, csp,
                constraint_matrix_phase, false);

        SparsityTools::distribute_sparsity_pattern(csp,
                dof_handler_phase.n_locally_owned_dofs_per_processor(),
                mpi_communicator, locally_relevant_dofs_phase);

        system_matrix_phase_constant.reinit(locally_owned_dofs_phase,
                locally_owned_dofs_phase, csp, mpi_communicator);

        jacobian_phase.reinit(locally_owned_dofs_phase,
                locally_owned_dofs_phase, csp, mpi_communicator);

        completely_distributed_solution_phase.reinit(locally_owned_dofs_phase,
                mpi_communicator);

    }



my solve method


template<int dim>
void CahnHilliard<dim>::solve_direct_phase() {

    SolverControl cn(100, 1.e-10, false, false);
    TrilinosWrappers::SolverDirect::AdditionalData add_data(
            "Amesos_Superludist");
    TrilinosWrappers::SolverDirect solver(cn, add_data);
    solver.solve(system_matrix_phase_constant,
            completely_distributed_solution_phase, system_rhs_phase);

    constraint_matrix_phase.distribute(completely_distributed_solution_phase);

    solution_phase = completely_distributed_solution_phase;
}





My back end data structures use TrilinosWrappers:

    LA::MPI::Vector solution_phase;
    LA::MPI::Vector solution_update;
    LA::MPI::Vector current_solution_phase;


LA::MPI::SparseMatrix  system_matrix_phase_constant,
    LA::MPI::Vector  system_rhs_phase;
LA::MPI::Vector locally_relevant_solution_phase;
    LA::MPI::Vector old_solution_phase;
    LA::MPI::Vector completely_distributed_solution_phase;



my assembly function looks like



template<int dim>
void CahnHilliard<dim>::assemble_system_phase() {

....initial setup for assemble

    for (; cell != endc; ++cell, ++cell_nse) {
        if (cell->is_locally_owned()) {
            fe_v.reinit(cell);
            fe_v_nse.reinit(cell_nse);
            cell->get_dof_indices(dof_indices);
            cell_nse->get_dof_indices(dof_indices_nse);
            assemble_cell_term_phase(fe_v, fe_v_nse, dof_indices,
                    dof_indices_nse);
        }
    }

    system_matrix_phase_constant.compress(VectorOperation::add);
    system_rhs_phase.compress(VectorOperation::add);

}



Maybe I am missing some point in the effective design of the distributed
computation like the one in Step-40.

Any advice will be appreciated

Thank You

Marek





2016-07-22 8:42 GMT+02:00 Vinetou Incucuna <[email protected]>:
> Hello,
> once more again
>
>>
>> 71874 dofs for Cahn-Hilliard part of system, cca 2000-3000 dofs per computer 
>> core
>
> I have solved the  Navier-Stokes-Cahn-Hilliard system on 32 cores as 
> decoupled.
> Algorithm:
> 1)assemble, solver Cahn-Hilliard    ->71874 dofs
> 2)assemble,solve Navier-Stokes    -> 2773956 dofs
>
> Therefore ad1)  is split into small subsystems (approximately 2000
> unknowns for each core) and therefore the Amesos_Superludis should
> have worked nicely?
> Hasn't it? I admit  that concerning  the ad)2 is your remark sound.
>
> Marek
>
>
> 2016-07-21 22:12 GMT+02:00 Bruno Turcksin <[email protected]>:
>>
>> Marek,
>>
>> 2016-07-21 15:59 GMT-04:00 Vinetou Incucuna <[email protected]>:
>> > I suppose, that the next
>> > step will be an   solver for general systems,
>> > like
>> > TrilinosWrappers::SolverBicgstab  or TrilinosWrappers::SolverGMRES
>> > ?
>> Yes that's right. However, the performance of Krylov solvers depends
>> on good preconditioning. So you will need to find a preconditioner for
>> a your problem.
>>
>> Best,
>>
>> Bruno
>>
>> --
>> 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.

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