I think the reason is that step-17 uses an inefficient constructor for
the system matrix:
    system_matrix.reinit (mpi_communicator,
                          dof_handler.n_dofs(),
                          dof_handler.n_dofs(),
                          n_local_dofs,
                          n_local_dofs,
                          dof_handler.max_couplings_between_dofs());

This will lead to very slow setup and large memory requirements. You
can try something closer to what is used in step-40 (using
DynamicSparsityPattern).


On Fri, Jul 29, 2016 at 3:41 PM, Pete Griffin <[email protected]> wrote:
> I did another plot (see attached) of Memory usage vs. #DOF with step-18.
> This is, as the other two, a 3d elasticity problem. The results in terms of
> memory usage were in line with step-8 and contrasted step-17. I will try to
> understand step-18 well enough to transfer the MPI/PETSc stuff to step-17
> with the goal of getting memory usage of step-17 down to expected. Whether
> the changes it will work on multiprocessor system I will have no way of
> knowing since I run only on a single processor.
>
> Pete Griffin
>
>
> On Monday, July 25, 2016 at 1:59:53 PM UTC-4, Pete Griffin wrote:
>>
>> I found that running step-8 and step-17 on a single processor Intel® Core™
>> i7-3630QM CPU @ 2.40GHz × 8 used substantially more Peak resident memory (>
>> 5x) than I thought it would. This surprised me since I thought from reading
>> step-17 that the memory increase was on the order of the solution vector
>> which should have been << 2x greater. I verified some of the larger memory
>> usage numbers using top.
>>
>> Is my assumption correct that 5x Peak resident memory is more than it
>> should be?
>>
>> The results of other simulations with a beam with body-force load and with
>> traction loads with and without HP and with and without MPI/PETSc all show
>> the same results and they agree with beam theory.
>>
>> The output of the modified step-8.cc and step-17.cc are attached along
>> with a plot of peak virtual memory and peak resident memory
>> vs. DOF. The changes between the original distributed step-8 and step-17,
>> with comments and extra newlines excluded (modified file >), are as below:
>>
>> Thanks beforehand
>>
>> Pete Griffin
>>
>>
>> ======================================================================================
>> diff ~/Documents/Zipstore2/dealii-8.4.1-PETSc/examples/step-8/step-8.cc
>> step-8.cc
>> 56c47,48
>> < // This again is C++:
>> ---
>> > #include <deal.II/base/utilities.h>
>>
>> 767a394,402
>> >
>> >         Utilities::System::MemoryStats stats;
>> >         Utilities::System::get_memory_stats(stats);
>> >         std::stringstream Str;
>> >
>> >         Str.str("");
>> >         Str << "   Peak virtual memory: " << stats.VmSize/1024 << " MB,
>> > Peak resident memory: "
>> >                << stats.VmRSS/1024 << " MB" << std::endl;
>> >         std::cout << Str.str();
>>
>> 781c411
>> <       Step8::ElasticProblem<2> elastic_problem_2d;
>> ---
>> >       Step8::ElasticProblem<3> elastic_problem_2d;
>>
>>
>> ======================================================================================
>> diff ~/Documents/Zipstore2/dealii-8.4.1-PETSc/examples/step-17/step-17.cc
>> ../step-17/step-17.cc
>> 84a50
>> > #include <deal.II/base/utilities.h>
>> 1015c355
>> <     for (unsigned int cycle=0; cycle<10; ++cycle)
>> ---
>> >     for (unsigned int cycle=0; cycle<8; ++cycle)
>> 1018d357
>> <
>> 1022c361
>> <             triangulation.refine_global (3);
>> ---
>> >             triangulation.refine_global (2);
>> 1049a383,391
>> >
>> >         Utilities::System::MemoryStats stats;
>> >         Utilities::System::get_memory_stats(stats);
>> >         std::stringstream Str;
>> >
>> >         Str.str("");
>> >         Str << "   Peak virtual memory: " << stats.VmSize/1024 << " MB,
>> > Peak resident memory: "
>> >                << stats.VmRSS/1024 << " MB" << std::endl;
>> >         std::cout << Str.str();
>> 1073,1074c403
>> <
>> <       ElasticProblem<2> elastic_problem;
>> ---
>> >       ElasticProblem<3> elastic_problem;
>>
>>
>> =============================================================================================
>>
> --
> 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.



-- 
Timo Heister
http://www.math.clemson.edu/~heister/

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