Hello,
I am developing code in the Step-40 like setting - Cahn-
Hilliard Navier Stokes system. I am witnessing very poor
performance of startup for large problems.
For the system with
1458 dof for Cahn Hilliard
47000 dofs for Navier Stoke
I got the following results (4 computer cores)
Setup CH, wall time: 0.1601s.
Setup NS, wall time: 29.9953s.
Setup, wall time: 30.1752s.
The times deteriorate if I solve larger problems.
I guess, that there is something wrong with my setup function
for the Navier-Stokes part of the system.
Could You help me please to find out, where is the performance flaw?
I use namespace ::LinearAlgebraTrilinos for backend of computations.
LA::MPI::Vector locally_relevant_solution_nse;
LA::MPI::Vector old_solution_nse;
LA::MPI::Vector solution_nse;
IndexSet locally_owned_dofs_nse;
IndexSet locally_relevant_dofs_nse;
{
computing_timer.enter_subsection("Setup NS");
dof_handler_nse.distribute_dofs(fe_nse);
DoFRenumbering::Cuthill_McKee(dof_handler_nse);
locally_owned_dofs_nse = dof_handler_nse.locally_owned_dofs();
DoFTools::extract_locally_relevant_dofs(dof_handler_nse,
locally_relevant_dofs_nse);
locally_relevant_solution_nse.reinit(locally_owned_dofs_nse,
locally_relevant_dofs_nse, mpi_communicator);
old_solution_nse.reinit(locally_owned_dofs_nse,
locally_relevant_dofs_nse, mpi_communicator);
old_solution_nse = 0;
solution_nse.reinit(locally_owned_dofs_nse,
locally_relevant_dofs_nse,
mpi_communicator);
system_rhs_nse.reinit(locally_owned_dofs_nse, mpi_communicator);
system_rhs_nse = 0;
constraint_matrix_nse.clear();
constraint_matrix_nse.reinit(locally_relevant_dofs_nse);
DoFTools::make_hanging_node_constraints(dof_handler_nse,
constraint_matrix_nse);
std::vector<bool> component_mask(dim + 1, false);
for (int i = 0; i < dim; ++i)
component_mask[i] = true; // velocities
VectorTools::interpolate_boundary_values(dof_handler_nse, 0,
ConstantFunction<dim>(0., dim + 1), constraint_matrix_nse,
component_mask);
VectorTools::interpolate_boundary_values(dof_handler_nse, 1,
InflowVelocityBoundaryValues<dim>(), constraint_matrix_nse,
component_mask);
constraint_matrix_nse.close();
LA::Vector<double> vec_old_solution(dof_handler_nse.n_dofs());
VectorTools::interpolate(dof_handler_nse, ZeroFunction<3>(dim + 1),
vec_old_solution);
old_solution_nse = vec_old_solution;
DynamicSparsityPattern csp(locally_relevant_dofs_nse);
DoFTools::make_sparsity_pattern(dof_handler_nse, csp,
constraint_matrix_nse, false);
SparsityTools::distribute_sparsity_pattern(csp,
dof_handler_nse.n_locally_owned_dofs_per_processor(),
mpi_communicator, locally_relevant_dofs_nse);
system_matrix_nse.reinit(locally_owned_dofs_nse,
locally_owned_dofs_nse,
csp, mpi_communicator);
computing_timer.leave_subsection();
}
Thank You
Marek C
--
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.