summer school - mpi (scorep cube plugins heatmap)

Issue #62 new
jg piccinali repo owner created an issue

Setup

Instrument (operators.f90)

module operators
...
#ifdef _SCOREP
#include "scorep/SCOREP_User.inc"
#endif
...
contains

subroutine diffusion
#ifdef _SCOREP
    SCOREP_USER_REGION_DEFINE( my_region_handle_i )
#endif
...
#ifdef _SCOREP
    SCOREP_USER_REGION_BEGIN( my_region_handle_i, "operators.interiors",
    SCOREP_USER_REGION_TYPE_DYNAMIC )
#endif
...
#ifdef _SCOREP
    SCOREP_USER_REGION_END( my_region_handle_i )
#endif
...
end subroutine diffusion
end module operators

Compile

  • make clean
  • make FTN="scorep --user ftn" FFLAGS="-D_SCOREP -cpp"
  • mv main GNU+sc142+user

Run

  • aprun -n 24 -N 8 -d 1 -j 1 ./GNU+sc142+user 1024 1024 50 0.001
========================================
                       Welcome to mini-stencil!
 MPI : pid           24
 mesh ::         1024 *        1024     dx =   9.7751710563898087E-004
 time ::           50 time steps from 0 ..    1.0000000000000000E-003
========================================
--------------------------------------------------------------------------------
 simulation took    1.0000000000000000       seconds
                 5031  conjugate gradient iterations   5031.       per second
                  226  nonlinear newton iterations

Report

  • module rm ddt # libQT conflicts => Cube will fail to load plugins
  • cube */profile.cubex

heatmap1.png

heatmap2.png

heatmap3.png

Comments (1)

  1. jg piccinali reporter

    Loop work estimates (CCE compiler only)

    • make FFLAGS="-eZ -homp -hprofile_generate"
    • pat_build -w CRAY
    • aprun -n 24 -N 8 -d 1 -j 1 ./CRAY+pat 1024 1024 50 0.001
    • pat_report CRAY+pat+24358-12t > xf
    Table 2:  Inclusive and Exclusive Time in Loops (from -hprofile_generate)
      Loop |     Loop |     Time |   Loop |    Loop |  Loop |  Loop |Function=/.LOOP[.]
      Incl |     Incl |    (Loop |    Hit |   Trips | Trips | Trips | PE=HIDE
     Time% |     Time |    Adj.) |        |     Avg |   Min |   Max |
    |-----------------------------------------------------------------------------
    | 97.5% | 3.983024 | 0.000027 |      1 |    50.0 |    50 |    50 |diffusion_serial_.LOOP.3.li.133
    => do timestep=1,nt
    | 97.4% | 3.981499 | 0.000535 |     50 |    50.0 |    50 |    50 |diffusion_serial_.LOOP.4.li.138
    => do it=1,50
    | 90.5% | 3.696327 | 0.022371 |    179 |    28.3 |     8 |    91 |ss_cg$linalg_.LOOP.1.li.308
    linalg.f90 => do iter=1, maxiters
    | 21.1% | 0.863581 | 0.028348 |   5651 |   168.7 |   168 |   172 |diffusion$operators_.LOOP.1.li.118
    | 20.4% | 0.835234 | 0.835234 | 953135 |   254.0 |   254 |   254 |diffusion$operators_.LOOP.2.li.119
    | 19.2% | 0.785677 | 0.785677 |  10307 | 43690.7 | 43520 | 44544 |ss_axpy$linalg_.LOOP.1.li.121
    | 19.2% | 0.784759 | 0.784759 |   9949 | 43690.7 | 43520 | 44544 |ss_lcomb$linalg_.LOOP.1.li.212
    | 17.9% | 0.731072 | 0.731072 |   5064 | 43690.7 | 43520 | 44544 |ss_scaled_diff$linalg_.LOOP.1.li.166
    |  5.5% | 0.226664 | 0.226664 |  10307 | 43690.7 | 43520 | 44544 |ss_dot$linalg_.LOOP.1.li.49
    |  1.4% | 0.058559 | 0.058559 |    408 | 43690.7 | 43520 | 44544 |ss_copy$linalg_.LOOP.1.li.232
    
  2. Log in to comment