tem_load_convergenceHeader Subroutine

private subroutine tem_load_convergenceHeader(me, conf, sub_handle)

Read the convergence variables from convergence subtables defined in configuration from the main lua file

If convergence is just a single table with single convergence entry then load only one convergence log exists with one or more variables using tem_load_convergenceHeader_single. Else if convergence is table of many log then allocate log and load each log type using tem_load_convergenceHeader_single Setup the values for the convergence entities

Arguments

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

list of the convergence entities to create

type(flu_State) :: conf

handle of the lua config file

integer, intent(in) :: sub_handle

table sub-handle for the convergence table


Calls

proc~~tem_load_convergenceheader~~CallsGraph proc~tem_load_convergenceheader tem_load_convergenceHeader interface~tem_load_shape tem_load_shape proc~tem_load_convergenceheader->interface~tem_load_shape proc~upper_to_lower upper_to_lower proc~tem_load_convergenceheader->proc~upper_to_lower proc~tem_load_condition tem_load_condition proc~tem_load_convergenceheader->proc~tem_load_condition proc~tem_load_reduction_spatial tem_load_reduction_spatial proc~tem_load_convergenceheader->proc~tem_load_reduction_spatial proc~aot_get_val aot_get_val proc~tem_load_convergenceheader->proc~aot_get_val proc~tem_timecontrol_dump tem_timeControl_dump proc~tem_load_convergenceheader->proc~tem_timecontrol_dump proc~tem_abort tem_abort proc~tem_load_convergenceheader->proc~tem_abort proc~tem_timecontrol_load tem_timeControl_load proc~tem_load_convergenceheader->proc~tem_timecontrol_load proc~aot_table_close aot_table_close proc~tem_load_convergenceheader->proc~aot_table_close proc~tem_load_shapes tem_load_shapes interface~tem_load_shape->proc~tem_load_shapes proc~tem_load_shape_single tem_load_shape_single interface~tem_load_shape->proc~tem_load_shape_single proc~tem_load_condition->proc~aot_table_close proc~aot_table_open aot_table_open proc~tem_load_condition->proc~aot_table_open proc~tem_load_cond_single tem_load_cond_single proc~tem_load_condition->proc~tem_load_cond_single proc~aot_table_length aot_table_length proc~tem_load_condition->proc~aot_table_length proc~tem_load_reduction_spatial->proc~aot_table_close proc~tem_load_reduction_spatial->proc~aot_table_open proc~tem_load_reduction_spatial->proc~aot_table_length proc~tem_load_reduction_single tem_load_reduction_single proc~tem_load_reduction_spatial->proc~tem_load_reduction_single proc~tem_time_dump tem_time_dump proc~tem_timecontrol_dump->proc~tem_time_dump mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~tem_timecontrol_load->proc~aot_get_val proc~tem_timecontrol_load->proc~aot_table_close proc~tem_timecontrol_load->proc~aot_table_open proc~tem_time_never tem_time_never proc~tem_timecontrol_load->proc~tem_time_never proc~tem_time_default_zero tem_time_default_zero proc~tem_timecontrol_load->proc~tem_time_default_zero proc~tem_time_needs_reduce tem_time_needs_reduce proc~tem_timecontrol_load->proc~tem_time_needs_reduce proc~tem_time_load tem_time_load proc~tem_timecontrol_load->proc~tem_time_load

Called by

proc~~tem_load_convergenceheader~~CalledByGraph proc~tem_load_convergenceheader tem_load_convergenceHeader proc~tem_convergence_load tem_convergence_load proc~tem_convergence_load->proc~tem_load_convergenceheader proc~tem_abortcriteria_load tem_abortCriteria_load proc~tem_abortcriteria_load->proc~tem_convergence_load proc~tem_simcontrol_load tem_simControl_load proc~tem_simcontrol_load->proc~tem_abortcriteria_load proc~tem_load_general tem_load_general proc~tem_load_general->proc~tem_simcontrol_load

Contents


Source Code

  subroutine tem_load_convergenceHeader(me, conf, sub_handle)
    ! --------------------------------------------------------------------------
    !> list of the convergence entities to create
    type( tem_convergence_type ), intent(out) :: me
    !> handle of the lua config file
    type( flu_state ) :: conf
    !> table sub-handle for the convergence table
    integer, intent(in) :: sub_handle
    ! --------------------------------------------------------------------------
    integer :: iError            ! error flag handle
    integer, allocatable :: vError(:)
    character(len=labelLen) :: norm
    ! --------------------------------------------------------------------------

    call aot_get_val( val       = me%header%varname, &
      &               ErrCode   = vError,            &
      &               maxLength = 100,               &
      &               L         = conf,              &
      &               thandle   = sub_handle,        &
      &               key       = "variable"         )

    if ( any(btest(vError, aoterr_Fatal)) ) then
      write(logUnit(1),*) 'FATAL Error occured, while retrieving'
      write(logUnit(1),*) 'list of variables to use in convergence'
      call tem_abort()
    end if

    me%header%nRequestedVars = size(me%header%varName)

    ! load time control to output convergence
    call tem_timeControl_load( conf           = conf,                          &
      &                        parent         = sub_handle,                    &
      &                        me             = me%header%timeControl )
    call tem_timeControl_dump(me%header%timeControl, logUnit(3))

    ! load convergence object shapes like point, line, plane
    call tem_load_shape( conf = conf, parent = sub_handle,                     &
                         me = me%header%geometry )

    if( size( me%header%geometry) < 1) then
      write(logUnit(1),*)'The geometrical objects for the convergence are'//  &
        &            ' not defined correctly.'
      call tem_abort()
    end if

    ! load reductions
    call tem_load_reduction_spatial( conf   = conf,                       &
      &                              parent = sub_handle,                 &
      &                   redSpatial_config = me%header%redSpatial_config )

    if( me%header%redSpatial_config%active ) then
      ! Check if the number of reductions correspond to the number of variables
      ! in the system
      if( size( me%header%redSpatial_config%reduceType ) &
        & /= me%header%nRequestedVars ) then
        write(logUnit(1),*) 'Error: In convergence.'
        write(logUnit(1),*) 'The number of defined reductions does not '// &
          &            'correspond to the '
        write(logUnit(1),*)'number of variables in the system. '
        call tem_abort()
      end if
    else
      write(logUnit(1),*) 'Error: No spatial reduction defined.'
      write(logUnit(1),*) 'NOTE: Convergence requires spatial reduction '//&
        &                 'for each variable'
      call tem_abort()
    end if

    ! get the kind of the convergence norm
    call aot_get_val( L       = conf,                                          &
      &               thandle = sub_handle,                                    &
      &               val     = norm,                                          &
      &               ErrCode = iError,                                        &
      &               key     = 'norm',                                        &
      &               default = 'simple' )
    norm = adjustl(norm)
    norm = upper_to_lower(norm)
    me%header%norm = trim(norm)
    select case( trim( norm ))
    case( 'simple' )
      me%norm_kind = norm_simple
      ! Only need one last value to compare against in the simple case.
      me%header%nLastVals = 1
    case( 'average')
      me%norm_kind = norm_average
      ! Get number of last values to check for convergence
      call aot_get_val( L       = conf,                &
        &               thandle = sub_handle,          &
        &               val     = me%header%nLastVals, &
        &               ErrCode = iError,              &
        &               key     = 'nvals',             &
        &               default = 1                    )
    case default
      write(logUnit(1),*) 'Error: Unknown convergence norm'
      write(logUnit(1),*) 'Solution: Supported norms '
      write(logUnit(1),*) '        * simple'
      write(logUnit(1),*) '        * average'
      call tem_abort
    end select


    ! type of convergence error: relative or absolute
    call aot_get_val( L       = conf,                    &
      &               thandle = sub_handle,              &
      &               val     = me%header%absoluteError, &
      &               ErrCode = iError,                  &
      &               key     = 'absolute',              &
      &               default = .false.                  )

    ! To decide whether to use get_point or get_element
    call aot_get_val( L       = conf,                  &
      &               thandle = sub_handle,            &
      &               val     = me%header%useGetPoint, &
      &               ErrCode = iError,                &
      &               default = .false.,               &
      &               key     = 'use_get_point'        )

    ! Get the number of Dofs to be written in the output
    ! The default is set to -1. If the dofs are not specified,
    ! all the dofs should be dumped
    call aot_get_val( L       = conf,            &
      &               thandle = sub_handle,      &
      &               val     = me%header%nDofs, &
      &               ErrCode = iError,          &
      &               default = -1,              &
      &               key     = 'ndofs'          )

    ! load condition for convergence
    call tem_load_condition( me      = me%header%cond, &
      &                      conf    = conf,           &
      &                      parent  = sub_handle      )

    me%header%nConditions = size(me%header%cond)
    ! check if there is condition for each variable
    if (me%header%nConditions /= me%header%nRequestedVars) then
      write(logUnit(1),*) 'Error: Nr. of conditions \= Nr. of variables '
      write(logUnit(1),"(2(A,I0))") 'nCond: ', me%header%nConditions, &
        &                           'nVars: ', me%header%nRequestedVars
      call tem_abort()
    end if

    write(logUnit(1),"(A,I0)") ' loaded convergence with nConditions=', &
      &                 me%header%nConditions
    write(logUnit(1),*) '    Norm : '//trim(norm)
    if( me%header%absoluteError ) then
      write(logUnit(1),*)'    absolute error'
    else
      write(logUnit(1),*)'    relative error'
    end if
    write(logUnit(1),"(A,I0)")'    nVal : ', me%header%nLastVals
    write(logUnit(7),*) '  Use get_point: ', me%header%useGetPoint
    call aot_table_close(L=conf, thandle=sub_handle )

  end subroutine tem_load_convergenceHeader