<p><b>mhoffman@lanl.gov</b> 2012-05-07 15:41:10 -0600 (Mon, 07 May 2012)</p><p>BRANCH COMMIT - land ice<br>
<br>
Added timers.  <br>
In addition to general timers, I added timers around each LifeV call.  Note that the LifeV timers get reported in the timer table under 'initialize' rather than 'land ice core run' because that is the first place they are called (but you can see from the number of calls that all calls are represented).<br>
</p><hr noshade><pre><font color="gray">Modified: branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_lifev.F
===================================================================
--- branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_lifev.F        2012-05-07 21:28:32 UTC (rev 1875)
+++ branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_lifev.F        2012-05-07 21:41:10 UTC (rev 1876)
@@ -178,6 +178,8 @@
 
    subroutine land_ice_lifev_block_init(block, err)
 
+       use mpas_timer
+
       !-----------------------------------------------------------------
       !
       ! input variables
@@ -260,9 +262,11 @@
 
       !zCell is supposed to be zero for L1L2 and FO solvers
       !nVertLevels should be equal to nVertLevelsSolve (no split the domain in the vertical direction)
+      call mpas_timer_start(&quot;velocity_solver_set_grid_data&quot;)
       call velocity_solver_set_grid_data(nCells, nEdges, nVertices, nVertLevels, nCellsSolve, nEdgesSolve, nVerticesSolve, &amp;
         cellsOnEdge, cellsOnVertex, verticesOnCell, verticesOnEdge, xCell, yCell, zCell, &amp;
         sendCellsArray, recvCellsArray, sendEdgesArray, recvEdgesArray, sendVerticesArray, recvVerticesArray)
+      call mpas_timer_stop(&quot;velocity_solver_set_grid_data&quot;)
 #endif
 
       !these can be deallocated because they have been copied on the c++ side
@@ -282,6 +286,8 @@
 
    subroutine land_ice_lifev_solve(mesh, state, err)
 
+       use mpas_timer
+
       !-----------------------------------------------------------------
       !
       ! input variables
@@ -344,33 +350,55 @@
 #ifdef USE_LIFEV
       ! LifeV calls to be made only when vertex mask changes
       if (anyVertexMaskChanged .eq. 1) then
+          call mpas_timer_start(&quot;velocity_solver_compute_2d_grid&quot;)
           call velocity_solver_compute_2d_grid(vertexMask)
+          call mpas_timer_stop(&quot;velocity_solver_compute_2d_grid&quot;)
 
           select case (config_dycore)
           case ('L1L2')
+              call mpas_timer_start(&quot;velocity_solver_init_L1L2&quot;)
               call velocity_solver_init_L1L2(LayerThicknessFractions)
+              call mpas_timer_stop(&quot;velocity_solver_init_L1L2&quot;)
           case ('FO')
+              call mpas_timer_start(&quot;velocity_solver_extrude_3d_grid&quot;)
               call velocity_solver_extrude_3d_grid(LayerThicknessFractions, lowerSurface, thickness)
+              call mpas_timer_stop(&quot;velocity_solver_extrude_3d_grid&quot;)
+              call mpas_timer_start(&quot;velocity_solver_init_FO&quot;)
               call velocity_solver_init_FO(LayerThicknessFractions)
+              call mpas_timer_stop(&quot;velocity_solver_init_FO&quot;)
           case ('Stokes')
+              call mpas_timer_start(&quot;velocity_solver_extrude_3d_grid&quot;)
               call velocity_solver_extrude_3d_grid(LayerThicknessFractions, lowerSurface, thickness)
+              call mpas_timer_stop(&quot;velocity_solver_extrude_3d_grid&quot;)
+              call mpas_timer_start(&quot;velocity_solver_init_stokes&quot;)
               call velocity_solver_init_stokes(LayerThicknessFractions)
+              call mpas_timer_stop(&quot;velocity_solver_init_stokes&quot;)
           end select
       endif
 
       ! LifeV calls to be made every time step (solve velocity!)
       select case (config_dycore)
       case ('L1L2')
+          call mpas_timer_start(&quot;velocity_solver_solve_L1L2&quot;)
           call velocity_solver_solve_L1L2(lowerSurface, thickness, beta, tracers(index_temperature,:,:), normalVelocity)
+          call mpas_timer_stop(&quot;velocity_solver_solve_L1L2&quot;)
           ! Optional calls to have LifeV output data files
+          call mpas_timer_start(&quot;velocity_solver export&quot;)
           call velocity_solver_export_2d_data(lowerSurface, thickness, beta)
           call velocity_solver_export_L1L2_velocity();
+          call mpas_timer_stop(&quot;velocity_solver export&quot;)
       case ('FO')
+          call mpas_timer_start(&quot;velocity_solver_solve_FO&quot;)
           call velocity_solver_solve_FO(lowerSurface, thickness, beta, tracers(index_temperature,:,:), normalVelocity)
+          call mpas_timer_stop(&quot;velocity_solver_solve_FO&quot;)
+          call mpas_timer_start(&quot;velocity_solver export&quot;)
           !call velocity_solver_export_2d_data(lowerSurface, thickness, beta)
-          call velocity_solver_export_FO_velocity();
+          call velocity_solver_export_FO_velocity()
+          call mpas_timer_stop(&quot;velocity_solver export&quot;)
       case ('Stokes')
+          call mpas_timer_start(&quot;stokes&quot;)
           ! call stokes
+          call mpas_timer_stop(&quot;stokes&quot;)
       end select
 
 #endif

Modified: branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_mpas_core.F
===================================================================
--- branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_mpas_core.F        2012-05-07 21:28:32 UTC (rev 1875)
+++ branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_mpas_core.F        2012-05-07 21:41:10 UTC (rev 1876)
@@ -47,6 +47,8 @@
 
       call simulation_clock_init(domain, dt, startTimeStamp)
 
+      write(*,*) 'Using ', trim(config_time_integration), ' time integration.'
+
       !call land_ice_temp_init(domain, dt, err)
 
       call land_ice_vel_init(domain, err)
@@ -170,20 +172,25 @@
       call mpas_rbf_interp_initialize(mesh)
       call mpas_init_reconstruct(mesh)
 
+      call mpas_timer_start(&quot;initialize velocity&quot;)
       call land_ice_vel_block_init(block, err_tmp)
+      call mpas_timer_stop(&quot;initialize velocity&quot;)
       err = ior(err, err_tmp)
 
 
       ! Initialize state ==== 
+      call mpas_timer_start(&quot;initial state calculation&quot;)
       call lice_diagnostic_solve(mesh, state, err)
       ! Set the initial flag to 1 so LifeV will calculate its mesh information the first time
       state % anyVertexMaskChanged % scalar = 1
+      call mpas_timer_stop(&quot;initial state calculation&quot;)
 
       ! Calculate an initial (diagnostic) velocity field consistent with the initial thickness field 
       ! \todo: skip this step if a velocity field was supplied in the I.C. input file
+      call mpas_timer_start(&quot;initial velocity calculation&quot;)
       call land_ice_vel_solve(mesh, state, err)
+      call mpas_timer_stop(&quot;initial velocity calculation&quot;)
 
-
  
       if(err == 1) then
           print *, &quot;An error has occurred in mpas_init_block. Aborting...&quot;
@@ -212,10 +219,11 @@
       character(len=32) :: timeStamp
       integer :: ierr
    
-
+      call mpas_timer_start(&quot;land ice core run&quot;)
       currTime = mpas_get_clock_time(clock, MPAS_NOW, ierr)
       call mpas_get_time(curr_time=currTime, dateTimeString=timeStamp, ierr=ierr)         
       write(0,*) 'Initial timestep ', timeStamp
+      write(6,*) 'Initial timestep ', timeStamp
 
       call write_output_frame(output_obj, output_frame, domain)
 
@@ -230,6 +238,7 @@
          currTime = mpas_get_clock_time(clock, MPAS_NOW, ierr)
          call mpas_get_time(curr_time=currTime, dateTimeString=timeStamp, ierr=ierr)         
          write(0,*) 'Doing timestep ', timeStamp
+         write(6,*) 'Doing timestep ', timeStamp
 
          call mpas_timer_start(&quot;time integration&quot;)
          call mpas_timestep(domain, itimestep, dt, timeStamp)
@@ -241,6 +250,7 @@
          !TODO: mpas_get_clock_ringing_alarms is probably faster than multiple mpas_is_alarm_ringing...
 
          if (mpas_is_alarm_ringing(clock, outputAlarmID, ierr=ierr)) then
+            call mpas_timer_start(&quot;write output&quot;)
             call mpas_reset_clock_alarm(clock, outputAlarmID, ierr=ierr)
             ! output_frame will always be &gt; 1 here unless it was reset after the maximum number of frames per outfile was reached
             if(output_frame == 1) then
@@ -248,18 +258,22 @@
                call mpas_output_state_init(output_obj, domain, &quot;OUTPUT&quot;, trim(timeStamp))
             end if
             call write_output_frame(output_obj, output_frame, domain)
+            call mpas_timer_stop(&quot;write output&quot;)
          end if
 
          if (mpas_is_alarm_ringing(clock, restartAlarmID, ierr=ierr)) then
+            call mpas_timer_start(&quot;write restart&quot;)
             call mpas_reset_clock_alarm(clock, restartAlarmID, ierr=ierr)
 
             ! Write one restart time per file
             call mpas_output_state_init(restart_obj, domain, &quot;RESTART&quot;, trim(timeStamp))
             call mpas_output_state_for_domain(restart_obj, domain, 1)
             call mpas_output_state_finalize(restart_obj, domain % dminfo)
+            call mpas_timer_stop(&quot;write restart&quot;)
          end if
 
       end do
+      call mpas_timer_stop(&quot;land ice core run&quot;)
 
    end subroutine mpas_core_run
    

Modified: branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_time_integration.F
===================================================================
--- branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_time_integration.F        2012-05-07 21:28:32 UTC (rev 1875)
+++ branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_time_integration.F        2012-05-07 21:41:10 UTC (rev 1876)
@@ -47,7 +47,7 @@
       integer :: err
       
 
-      write(*,*) 'Using ', trim(config_time_integration), ' time integration.'
+      !write(*,*) 'Using ', trim(config_time_integration), ' time integration.'
       select case (config_time_integration)
       case ('ForwardEuler')
          call land_ice_time_integrator_forwardeuler(domain, dt)

Modified: branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_time_integration_forwardeuler.F
===================================================================
--- branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_time_integration_forwardeuler.F        2012-05-07 21:28:32 UTC (rev 1875)
+++ branches/land_ice_projects/implement_core/src/core_land_ice/mpas_land_ice_time_integration_forwardeuler.F        2012-05-07 21:41:10 UTC (rev 1876)
@@ -18,6 +18,7 @@
    use mpas_configure
    use mpas_constants
    use mpas_dmpar
+   use mpas_timer
    use land_ice_vel, only: land_ice_vel_solve
    use land_ice_tendency
 
@@ -101,6 +102,7 @@
 
 
 ! === Calculate Tendencies ========================
+         call mpas_timer_start(&quot;calculate tendencies&quot;)
          ! Calculate thickness tendency using state at time n
          call lice_tend_h(mesh, normalVelocityOld, layerThicknessOld, thicknessOld, thickness_tend, dt, err)
 
@@ -116,37 +118,38 @@
 
          ! (Calculate tracer tendencies)
          !call ()
+         call mpas_timer_stop(&quot;calculate tendencies&quot;)
 
 
-
 ! === Compute new state for prognostic variables ==================================
 ! (once implicit column physics are added (i.e. temp diffusion), these calculations will need to be adjusted to apply to the new values as needed)
+         call mpas_timer_start(&quot;calc. new prognostic vars&quot;)
          ! Update thickness (and tracers eventually)
          thicknessNew = thicknessOld + thickness_tend * dt / SecondsInYear        
 
          !Optionally print some information about the new state
          !print *, 'thickness_tend maxval:', maxval(thickness_tend(1:mesh % nCellsSolve))       
          !print *, 'thicknessOld maxval:', maxval(thicknessOld(1:mesh % nCellsSolve))
-         print *, 'thicknessNew maxval:', maxval(thicknessNew(1:mesh % nCellsSolve))
+         print *, '  thicknessNew maxval:', maxval(thicknessNew(1:mesh % nCellsSolve))
          mask = 0
          where (thicknessNew .lt. 0.0)
             mask = 1
             thicknessNew = 0.0
          end where
          if (sum(mask) .gt. 0) then
-                 print *,'Cells with negative thickness (set to 0):',sum(mask)
+                 print *,'  Cells with negative thickness (set to 0):',sum(mask)
          endif
 
          mask = 0
          where (thicknessNew .gt. 0.0)
             mask = 1
          end where
-         print *,'Cells with nonzero thickness:', sum(mask)
+         print *,'  Cells with nonzero thickness:', sum(mask)
 
+         call mpas_timer_stop(&quot;calc. new prognostic vars&quot;)
 
-
 ! === Calculate diagnostic variables for new state =====================
-
+         call mpas_timer_start(&quot;calc. diagnostic vars except vel&quot;)
          call lice_diagnostic_solve(mesh, stateNew, err)  ! perhaps velocity solve should move in here.
 
          ! Determine if the vertex mask changed during this time step for this block (needed for LifeV)
@@ -160,12 +163,15 @@
          ! Determine if any blocks on this processor had a change to the vertex mask
          procVertexMaskChanged = max(procVertexMaskChanged, blockVertexMaskChanged)
 
+         call mpas_timer_stop(&quot;calc. diagnostic vars except vel&quot;)
+
          block =&gt; block % next
       end do
        
       ! Determine if the vertex mask has changed on any processor (need to exit the block loop to do so)
       call mpas_dmpar_max_int(dminfo, procVertexMaskChanged, anyVertexMaskChanged)
 
+      call mpas_timer_start(&quot;velocity solve&quot;)
       block =&gt; domain % blocklist
       do while (associated(block))
 
@@ -184,8 +190,8 @@
 
          block =&gt; block % next
       end do
+      call mpas_timer_stop(&quot;velocity solve&quot;)
 
-
 ! === Cleanup &amp; Misc. =============================
 
 

</font>
</pre>