tem_init_convergence Subroutine

public subroutine tem_init_convergence(me, tree, varSys, bc_prop, globProc, stencil, nDofs)

Initialize the convergence subtreee

Identify, how many and which elements exist on my local process and are requested from the convergences Empty convergence are removed, so the convergence(:) might be re-allocated

Arguments

Type IntentOptional Attributes Name
type(tem_convergence_type), intent(inout), allocatable :: me(:)

convergence descriptions

type(treelmesh_type), intent(in) :: tree

Global mesh from which the elements are identified and then stored to sub-meshes inside the convergences

type(tem_varSys_type), intent(in) :: varSys

solver-provided variable systems

type(tem_BC_prop_type), intent(in) :: bc_prop

bc property that used to identify elements of certain BCs

type(tem_comm_env_type), intent(in) :: globProc

Process description to use.

type(tem_stencilHeader_type), intent(in), optional :: stencil

stencil used to create subTree of boundary type

integer, intent(in), optional :: nDofs

The number of dofs for each scalar variable of the equation system


Calls

proc~~tem_init_convergence~~CallsGraph proc~tem_init_convergence tem_init_convergence proc~tem_horizontalspacer tem_horizontalSpacer proc~tem_init_convergence->proc~tem_horizontalspacer proc~tem_create_subtree_of tem_create_subTree_of proc~tem_init_convergence->proc~tem_create_subtree_of proc~tem_reduction_spatial_init tem_reduction_spatial_init proc~tem_init_convergence->proc~tem_reduction_spatial_init proc~tem_abort tem_abort proc~tem_init_convergence->proc~tem_abort proc~tem_create_varmap tem_create_varMap proc~tem_init_convergence->proc~tem_create_varmap proc~tem_create_subtree_of->proc~tem_horizontalspacer proc~tem_create_subtree_of->proc~tem_abort proc~tem_shape_initlocal tem_shape_initLocal proc~tem_create_subtree_of->proc~tem_shape_initlocal mpi_comm_split mpi_comm_split proc~tem_create_subtree_of->mpi_comm_split mpi_comm_size mpi_comm_size proc~tem_create_subtree_of->mpi_comm_size mpi_comm_rank mpi_comm_rank proc~tem_create_subtree_of->mpi_comm_rank proc~tem_subtree_from tem_subTree_from proc~tem_create_subtree_of->proc~tem_subtree_from interface~tem_copypropertybits tem_copyPropertyBits proc~tem_create_subtree_of->interface~tem_copypropertybits mpi_allreduce mpi_allreduce proc~tem_create_subtree_of->mpi_allreduce proc~tem_shape2subtree tem_shape2subTree proc~tem_create_subtree_of->proc~tem_shape2subtree interface~init~22 init proc~tem_create_subtree_of->interface~init~22 interface~destroy~22 destroy proc~tem_create_subtree_of->interface~destroy~22 interface~tem_seteffboundingbox tem_setEffBoundingBox proc~tem_create_subtree_of->interface~tem_seteffboundingbox mpi_abort mpi_abort proc~tem_abort->mpi_abort interface~append~25 append proc~tem_create_varmap->interface~append~25 interface~positionofval~5 positionofval proc~tem_create_varmap->interface~positionofval~5 interface~truncate~9 truncate proc~tem_create_varmap->interface~truncate~9 interface~init~9 init proc~tem_create_varmap->interface~init~9

Contents

Source Code


Source Code

  subroutine tem_init_convergence( me, tree, varSys, bc_prop, globProc,&
    &                              stencil, nDofs )
    ! -------------------------------------------------------------------- !
    !> convergence descriptions
    type(tem_convergence_type),intent(inout), allocatable :: me(:)
    !> Global mesh from which the elements are identified and then stored to
    !! sub-meshes inside the convergences
    type(treelmesh_type), intent(in) :: tree
    !> bc property that used to identify elements of certain BCs
    type( tem_bc_prop_type ), intent(in) :: bc_prop
    !> solver-provided variable systems
    type(tem_varSys_type), intent(in) :: varSys
    !> Process description to use.
    type(tem_comm_env_type), intent(in) :: globProc
    !> stencil used to create subTree of boundary type
    type(tem_stencilHeader_type), optional, intent(in) :: stencil
    !> The number of dofs for each scalar variable of the equation system
    integer, intent(in), optional :: nDofs
    ! -------------------------------------------------------------------- !
    integer :: iConv, nConv, nVars
    integer :: nChunks, chunkSize, nElems, maxComponents, nPoints
    ! temporary convergence array
    type( tem_convergence_type ),allocatable :: tempConv(:)
    ! -------------------------------------------------------------------- !

    call tem_horizontalSpacer(fUnit=logUnit(1))
    write(logUnit(1),*) 'Setting up the convergence infrastructure'
    nConv = 0


    ! Allocate the temporary convergence
    allocate(tempConv(size(me)))

    do iConv = 1, size(me)
      write(logUnit(10),"(A,I0)") 'Initializing convergence object ', iConv

      ! map variables
      ! create convergence variable position in the global varSys
      call tem_create_varMap( varname = me(iConv)%header%varname, &
        &                     varSys  = varSys,                   &
        &                     varMap  = me(iConv)%varMap          )

      ! @todo KM: If variable not found in varSys then skip that variable
      ! reduction and condition info while copying convergence
      !
      ! Terminate if number of variables to check for convergence does not
      ! match with variables found in varMap
      if (me(iConv)%varMap%varPos%nVals /= me(iConv)%header%nRequestedVars) then
        write(logUnit(1),*) 'Error: Mapping Convergences variables'
        write(logUnit(1),*) 'Variables defined in convergence '// &
          &                 'table ', iConv, ' are not found in varSys'
        call tem_abort()
      end if

      ! -----------------------------------------------------------------------
      ! identify convergence elements
      ! -----------------------------------------------------------------------
      call tem_create_subTree_of( inTree    = tree,                           &
        &                         bc_prop   = bc_prop,                        &
        &                         stencil   = stencil,                        &
        &                         subTree   = me( iConv )%subTree,            &
        &                         storePnts = me (iConv )%header%useGetPoint, &
        &                         inShape   = me( iConv )%header%geometry     )

      ! get rid of the empty convergence in order
      ! to avoid empty writes to disk
      if ( me(iConv)%subTree%useGlobalMesh ) then

        ! set convergence communicator
        me(iConv)%proc = globProc

        nConv = nConv + 1
        tempConv( nConv ) = me( iConv )

      else if ( me(iConv)%subTree%nElems > 0 ) then

        nConv = nConv + 1 ! Increase number of log

        ! set convergence communicator
        me(iConv)%proc%comm      = me(iConv)%subTree%global%comm
        me(iConv)%proc%rank      = me(iConv)%subTree%global%myPart
        me(iConv)%proc%comm_size = me(iConv)%subTree%global%nParts
        me(iConv)%proc%root      = 0

        ! Copy all entries from the log derived type to the temporary one.
        ! this might not work on all compilers!!
        ! This assignment is realized by operator overloader Copy_convergence
        tempConv( nConv ) = me( iConv )
      end if ! useGlobalMesh
    end do  ! nConv

    deallocate(me)
    allocate( me(nConv) )

    do iConv = 1, nConv
      ! Copy the stuff from the temporary track
      me(iConv) = tempConv(iConv)

      ! number of variables in varMap
      nVars = me(iConv)%varMap%varPos%nVals

      ! nDofs is valid only for get_element
      if (me(iConv)%header%useGetPoint) then
        me(iConv)%nDofs = 1
      else
        if (present(nDofs)) then
          ! out_config%nDofs is set to -1 if unspecied
          ! in the config file. In this case all the dof's
          ! should be dumped
          if (me(iConv)%header%nDofs < 0) then
            me(iConv)%nDofs = nDofs
          else
            ! Otherwise the number of dofs dumped should
            ! be what's specified in the config
            me(iConv)%nDofs = me(iConv)%header%nDofs
          end if
        else
          me(iConv)%nDofs = 1
        end if
      end if

      if (me(iConv)%subTree%useGlobalMesh) then
        nElems = tree%nElems
        nPoints = tree%nElems
      else
        nElems = me(iConv)%subTree%nElems
        nPoints = me(iConv)%subTree%nPoints
      end if

      ! max nComponent in current convergence variables
      maxComponents = maxval(varSys%method%val(me(iConv)%varMap               &
        &                                     %varPos%val(:nVars))%nComponents)

      if (me(iConv)%header%useGetPoint) then
        chunkSize = min(io_buffer_size/maxComponents, nPoints)
      else
        chunkSize = min(io_buffer_size/(maxComponents*me(iConv)%nDofs), nElems)
      end if
      if ( (nElems > 0) .and. (chunkSize == 0) ) then
        write(logUnit(0),*)'Error in init_convergence: '
        write(logUnit(0),*) 'The chosen io_buffer_size of ', io_buffer_size
        write(logUnit(0),*) 'is too small for evaluating ', maxComponents
        write(logUnit(0),*) 'scalar values'
        write(logUnit(0),*) 'Please increase the io_buffer_size to at &
          & least ', real(maxComponents*me(iConv)%nDofs) / real(131072), ' MB!'
        call tem_abort()
      end if

      nChunks = 0
      if (chunkSize>0) then
        if (me(iConv)%header%useGetPoint) then
          nChunks = ceiling(real(nPoints, kind=rk)/real(chunkSize, kind=rk))
        else
          nChunks = ceiling(real(nElems, kind=rk)/real(chunkSize, kind=rk))
        end if
      else
        nChunks = 0
      end if

      me(iConv)%nChunks = nChunks
      me(iConv)%chunkSize = chunkSize
      me(iConv)%nChecks = 0

      ! Initialize reduction
      call tem_reduction_spatial_init(                                 &
        &      me                = me(iConv)%redSpatial,               &
        &      redSpatial_config = me(iConv)%header%redSpatial_config, &
        &      varSys            = varSys,                             &
        &      varPos            = me(iConv)%varMap%varPos%val(:nVars) )

      ! Allocate some arrays
      allocate( me(iConv)%lastState( me(iConv)%header%nLastVals, &
        &                           me(iConv)%varMap%nScalars ) )
      me(iConv)%lastState = huge( me(iConv)%lastState(1,1) )          &
        &                 / real( me(iConv)%header%nLastVals, kind=rk )

    end do

    deallocate(tempConv)

    call tem_horizontalSpacer(fUnit=logUnit(1))

  end subroutine tem_init_convergence