derivedType_memAlloc Program

Uses

  • program~~derivedtype_memalloc~~UsesGraph program~derivedtype_memalloc derivedType_memAlloc module~tem_arrayofarrays_module tem_arrayofarrays_module program~derivedtype_memalloc->module~tem_arrayofarrays_module module~tem_logging_module tem_logging_module program~derivedtype_memalloc->module~tem_logging_module module~tem_debug_module tem_debug_module program~derivedtype_memalloc->module~tem_debug_module module~tem_dyn_array_module tem_dyn_array_module program~derivedtype_memalloc->module~tem_dyn_array_module module~tem_arrayofarrays_module->module~tem_dyn_array_module module~env_module env_module module~tem_arrayofarrays_module->module~env_module module~aot_table_module aot_table_module module~tem_logging_module->module~aot_table_module module~tem_logging_module->module~env_module module~aotus_module aotus_module module~tem_logging_module->module~aotus_module module~tem_debug_module->module~tem_logging_module module~tem_debug_module->module~aot_table_module module~tem_debug_module->module~env_module module~flu_binding flu_binding module~tem_debug_module->module~flu_binding module~tem_tools_module tem_tools_module module~tem_debug_module->module~tem_tools_module module~tem_dyn_array_module->module~env_module module~env_module->module~flu_binding module~env_module->module~aotus_module mpi mpi module~env_module->mpi iso_fortran_env iso_fortran_env module~env_module->iso_fortran_env module~tem_tools_module->module~env_module

Calls

program~~derivedtype_memalloc~~CallsGraph program~derivedtype_memalloc derivedType_memAlloc interface~tem_logging_init tem_logging_init program~derivedtype_memalloc->interface~tem_logging_init interface~init~4 init program~derivedtype_memalloc->interface~init~4 proc~tem_reportstatus tem_reportStatus program~derivedtype_memalloc->proc~tem_reportstatus interface~destroy~4 destroy program~derivedtype_memalloc->interface~destroy~4 interface~append~5 append program~derivedtype_memalloc->interface~append~5 proc~tem_logging_init_logger tem_logging_init_logger interface~tem_logging_init->proc~tem_logging_init_logger proc~tem_logging_init_primary tem_logging_init_primary interface~tem_logging_init->proc~tem_logging_init_primary proc~init_ga_dynlong init_ga_dynlong interface~init~4->proc~init_ga_dynlong proc~my_status_string my_status_string proc~tem_reportstatus->proc~my_status_string interface~tem_logging_isactive tem_logging_isActive proc~tem_reportstatus->interface~tem_logging_isactive proc~destroy_ga_dynlong destroy_ga_dynlong interface~destroy~4->proc~destroy_ga_dynlong proc~append_ga_dynlong_vec append_ga_dynlong_vec interface~append~5->proc~append_ga_dynlong_vec proc~append_ga_dynlong append_ga_dynlong interface~append~5->proc~append_ga_dynlong interface~expand~3 expand proc~append_ga_dynlong_vec->interface~expand~3 proc~tem_connect_tonull tem_connect_toNull proc~tem_logging_init_logger->proc~tem_connect_tonull proc~newunit newunit proc~tem_logging_init_logger->proc~newunit proc~print_self_status print_self_status proc~my_status_string->proc~print_self_status proc~my_status_string->proc~newunit proc~append_ga_dynlong->interface~expand~3 proc~tem_logging_isactive_for tem_logging_isActive_for interface~tem_logging_isactive->proc~tem_logging_isactive_for proc~tem_logging_isactive_primary tem_logging_isActive_primary interface~tem_logging_isactive->proc~tem_logging_isactive_primary proc~tem_logging_init_primary->proc~tem_logging_init_logger proc~tem_connect_tonull->proc~newunit proc~print_self_status->proc~newunit proc~expand_ga_dynlong expand_ga_dynlong interface~expand~3->proc~expand_ga_dynlong

Contents

Source Code


Variables

Type AttributesNameInitial
type(outerDT_type), allocatable:: simpleArr(:)
type(outerDTarr_type), allocatable:: nestedArr(:)
type(outerDTpnt_type), allocatable:: nestedPnt(:)
type(grw_dynlongarray_type) :: array
type(dyn_longarray_type) :: tDynLongArray
type(tem_debug_type) :: debug
integer :: QQ
integer :: nElems
integer :: iElem
integer :: iSize
integer :: iArraySize
integer :: performTest
integer :: iInner
integer, allocatable:: elemSize(:)
integer, allocatable:: arraySize(:)
logical :: error
character(len=64) :: buffer

Derived Types

type :: innerDT_type

Components

TypeVisibility AttributesNameInitial
integer, public :: test

type :: outerDT_type

Components

TypeVisibility AttributesNameInitial
type(innerDT_type), public :: inner

Components

TypeVisibility AttributesNameInitial
type(innerDT_type), public, allocatable:: inner(:)

Components

TypeVisibility AttributesNameInitial
type(innerDT_type), public, pointer:: inner(:)

Source Code

program derivedType_memAlloc

  ! include treelm modules
  use tem_debug_module,         only: tem_reportStatus, tem_debug_type
  use tem_dyn_array_module,     only: dyn_longArray_type, init, append, destroy
  use tem_arrayofarrays_module, only: grw_dynlongArray_type, init, append, &
    &                                 destroy
  use tem_logging_module,       only: tem_logging_init

  implicit none
  type innerDT_type
    integer :: test
  end type innerDT_type

  type outerDT_type
    type( innerDT_type ) :: inner
  end type outerDT_type
  type outerDTarr_type
    type( innerDT_type ), allocatable :: inner(:)
  end type outerDTarr_type
  type outerDTpnt_type
    type( innerDT_type ), pointer :: inner(:)
  end type outerDTpnt_type
  !----------------------------------------------------------------------------
  type(outerDT_type), allocatable :: simpleArr(:)
  type(outerDTarr_type), allocatable :: nestedArr(:)
  type(outerDTpnt_type), allocatable :: nestedPnt(:)
  type(grw_dynlongArray_type) :: array
  type(dyn_longArray_type) :: tDynLongArray
  type(tem_debug_type) :: debug
  integer :: QQ, nElems, iElem, iSize
  integer :: iArraySize, performTest, iInner
  integer, allocatable :: elemSize(:), arraySize(:)
  logical :: error
  character(len=64) :: buffer
  !----------------------------------------------------------------------------
  performTest = 4
  error = .false.

  ! initialize the debug type
  call tem_logging_init( me    = debug%logger,                                 &
    &                    level = 10,                                           &
    &                    rank  = 0,                                            &
    &                    filename = 'tem_dT_memTrace.res' )

  ! Get initial memory demand
  call tem_reportStatus( debug = debug, level = 1, &
                         text = 'starting overhead', string = 'VmRSS')
  allocate(elemSize(3))
  elemSize = [ 10000, 100000, 1000000 ]
  allocate(arraySize(1))
  arraySize = [ 10 ]

  select case( performTest )
  case default
    write(*,*) 'error, no valid test chosen.'
    write(*,*) 'choose: 1 ... grw_dynlongArray_test'
    write(*,*) '        2 ... simple derived type test'
    write(*,*) '        3 ... nested derived type test'

  case( 4 ) ! nested derived type test
    do iArraySize = 1, size( arraySize )
      QQ = arraySize( iArraySize )
      write(*,*) 'Starting the memory analysis for internal arr size of ', QQ
      do iSize = 1, size(elemSize)
        nElems = elemSize( iSize )
        write(*,*) 'Starting the memory analysis for nested pointer derived type arrays'
        write(buffer,'(a,i10,a,i3)') 'starting to allocate with nElems ', nElems
        call tem_reportStatus( debug = debug, level = 1, &
                               text = buffer, string = 'VmRSS')
        allocate( nestedPnt( nElems ))
        do iElem = 1, nElems
          allocate( nestedPnt( iElem )%inner( QQ ))
          do iInner = 1, QQ
            nestedPnt( iElem )%inner( iInner )%test = 1
          enddo
        enddo
        write(buffer,'(a,i10,a,i3)') 'finished allocation    ', nElems
        call tem_reportStatus( debug = debug, level = 1, &
                               text = buffer, string = 'VmRSS')

        do iElem = 1, nElems
          deallocate( nestedPnt( iElem )%inner )
        end do
        deallocate( nestedPnt )
        write(buffer,'(a,i10,a,i3)') 'after deallocation'
        call tem_reportStatus( debug = debug, level = 1, &
                                 text = buffer, string = 'VmRSS')
      end do
    end do
  case( 3 ) ! nested derived type test
    do iArraySize = 1, size( arraySize )
      QQ = arraySize( iArraySize )
      write(*,*) 'Starting the memory analysis for internal arr size of ', QQ
      do iSize = 1, size(elemSize)
        nElems = elemSize( iSize )
        write(*,*) 'Starting the memory analysis for nested allocatable derived type arrays'
        write(buffer,'(a,i10,a,i3)') 'starting to allocate with nElems ', nElems
        call tem_reportStatus( debug = debug, level = 1, &
                               text = buffer, string = 'VmRSS')
        allocate( nestedArr( nElems ))
        do iElem = 1, nElems
          allocate( nestedArr( iElem )%inner( QQ ))
          do iInner = 1, QQ
            nestedArr( iElem )%inner( iInner )%test = 1
          enddo
        enddo
        write(buffer,'(a,i10,a,i3)') 'finished allocation    ', nElems
        call tem_reportStatus( debug = debug, level = 1, &
                               text = buffer, string = 'VmRSS')

        do iElem = 1, nElems
          deallocate( nestedArr( iElem )%inner )
        end do
        deallocate( nestedArr )
        write(buffer,'(a,i10,a,i3)') 'after deallocation'
        call tem_reportStatus( debug = debug, level = 1, &
                                 text = buffer, string = 'VmRSS')
      end do
    end do
  case( 2 ) ! simple derived type test
    do iSize = 1, size(elemSize)
      nElems = elemSize( iSize )
      write(*,*) 'Starting the memory analysis for simple nested derived type arrays'
      write(buffer,'(a,i10,a,i3)') 'starting to allocate with nElems ', nElems
      call tem_reportStatus( debug = debug, level = 1, &
                             text = buffer, string = 'VmRSS')
      allocate( simpleArr( nElems ))
      do iElem = 1, nElems
        simpleArr( iElem )%inner%test = 1
      enddo
      write(buffer,'(a,i10,a,i3)') 'finished allocation    ', nElems
      call tem_reportStatus( debug = debug, level = 1, &
                             text = buffer, string = 'VmRSS')

      deallocate( simpleArr )
      write(buffer,'(a,i10,a,i3)') 'after deallocation'
      call tem_reportStatus( debug = debug, level = 1, &
                               text = buffer, string = 'VmRSS')

    end do
  case( 1 ) ! grw_dynlongArray test
    do iArraySize = 1, size( arraySize )
      QQ = arraySize( iArraySize )
      write(*,*) 'Starting the memory analysis for internal arr size of ', QQ

      do iSize = 1, size(elemSize)
        nElems = elemSize( iSize )
        write(*,*) 'analyzing with nElems ', nElems
        call init( me = array )
        write(buffer,'(a,i10,a,i3)') 'starting to add nElems ', nElems, ' QQ', QQ
        call tem_reportStatus( debug = debug, level = 1, &
                               text = buffer, string = 'VmRSS')

        ! Create a list of dummy entries
        do iElem = 1, nElems
          !write(*,*) 'element number ', iElem
          ! Append a dummy element
          call append( me = array, val = tDynLongArray )
          array%val( iElem )%nVals = 0
        end do

        write(buffer,'(a,i10,a,i3)') 'finished adding nElems ', nElems, ' QQ', QQ
        call tem_reportStatus( debug = debug, level = 1, &
                               text = buffer, string = 'VmRSS')

        call destroy( me = array )
        write(buffer,'(a,i10,a,i3)') 'after destruction'
        call tem_reportStatus( debug = debug, level = 1, &
                               text = buffer, string = 'VmRSS')
      end do
      call destroy( me = tDynLongArray )
    end do
  end select

  close( debug%logger%funit(0) )

  if( error ) then
    write(*,*) 'FAILED'
    stop -1
  else
    write(*,*) 'PASSED'
    stop 0
  end if

end program derivedType_memAlloc