master_sort_and_group Subroutine

private subroutine master_sort_and_group(Array, Array_groups, Array_size, subgroup_info)

Uses

  • proc~~master_sort_and_group~~UsesGraph proc~master_sort_and_group master_sort_and_group module~w90_constants w90_constants proc~master_sort_and_group->module~w90_constants module~w90_io w90_io proc~master_sort_and_group->module~w90_io module~w90_parameters w90_parameters proc~master_sort_and_group->module~w90_parameters module~w90_hamiltonian w90_hamiltonian proc~master_sort_and_group->module~w90_hamiltonian module~w90_io->module~w90_constants module~w90_parameters->module~w90_constants module~w90_parameters->module~w90_io module~w90_hamiltonian->module~w90_constants module~w90_comms w90_comms module~w90_hamiltonian->module~w90_comms module~w90_comms->module~w90_constants module~w90_comms->module~w90_io

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(inout), dimension(2, Array_size):: Array
integer, intent(in), dimension(:):: Array_groups
integer, intent(in) :: Array_size
integer, intent(out), allocatable, dimension(:, :):: subgroup_info

Calls

proc~~master_sort_and_group~~CallsGraph proc~master_sort_and_group master_sort_and_group proc~sort sort proc~master_sort_and_group->proc~sort proc~group group proc~master_sort_and_group->proc~group

Called by

proc~~master_sort_and_group~~CalledByGraph proc~master_sort_and_group master_sort_and_group proc~tran_lcr_2c2_sort tran_lcr_2c2_sort proc~tran_lcr_2c2_sort->proc~master_sort_and_group proc~tran_main tran_main proc~tran_main->proc~tran_lcr_2c2_sort program~wannier wannier program~wannier->proc~tran_main proc~wannier_run wannier_run proc~wannier_run->proc~tran_main

Contents

Source Code


Source Code

  subroutine master_sort_and_group(Array, Array_groups, Array_size, subgroup_info)
    !=============================================================!
    ! General sorting and grouping subroutine which takes Array,  !
    ! an ordered in conduction direction array of wannier function!
    ! indexes and positions, and returns the ordered (and grouped)!
    ! indexes and positions after considering the other two       !
    ! directions. Sub group info is also return for later checks. !
    !=============================================================!

    use w90_constants, only: dp
    use w90_io, only: io_error, stdout, io_stopwatch
    use w90_parameters, only: iprint, timing_level
    use w90_hamiltonian, only: wannier_centres_translated

    implicit none

    integer, intent(in), dimension(:)                 :: Array_groups
    integer, intent(in)                              :: Array_size

    integer, intent(out), allocatable, dimension(:, :)  :: subgroup_info

    real(dp), intent(inout), dimension(2, Array_size)  :: Array

    integer                                         :: i, j, k, Array_num_groups, increment, ierr, &
                                                       subgroup_increment, group_num_subgroups
    integer, allocatable, dimension(:)                :: group_subgroups

    real(dp), allocatable, dimension(:, :)             :: group_array, sorted_group_array, &
                                                          subgroup_array, sorted_subgroup_array
    character(30)                                   :: fmt_2

    if (timing_level > 2) call io_stopwatch('tran: lcr_2c2_sort: master_sort', 1)

    allocate (subgroup_info(size(Array_groups), maxval(Array_groups)), stat=ierr)
    if (ierr /= 0) call io_error('Error in allocating subgroup_info in master_sort_and_group')
    subgroup_info = 0
    !
    !Number of groups inside the principal layer
    !
    Array_num_groups = size(Array_groups)

    !
    !Convenient variable which will be amended later. Used to appropriately extract the group array from the Array
    !
    increment = 1
    !
    !Loop over groups inside Array
    !
    do j = 1, Array_num_groups
      allocate (group_array(2, Array_groups(j)), stat=ierr)
      if (ierr /= 0) call io_error('Error in allocating group_array in master_sort_and_group')
      allocate (sorted_group_array(2, Array_groups(j)), stat=ierr)
      if (ierr /= 0) call io_error('Error in allocating sorted_group_array in master_sort_and_group')
      !
      !Extract the group from the Array
      !
      group_array = Array(:, increment:increment + Array_groups(j) - 1)
      !
      !Updating group_array to contain coord(2)
      !
      do k = 1, Array_groups(j)
        group_array(2, k) = wannier_centres_translated(coord(2), int(group_array(1, k)))
      enddo

      call sort(group_array, sorted_group_array)
      call group(sorted_group_array, group_subgroups)

      group_num_subgroups = size(group_subgroups)

      if (iprint .ge. 4) then
        !
        !Printing subgroup breakdown
        !
        write (fmt_2, '(i5)') group_num_subgroups
        fmt_2 = adjustl(fmt_2)
        fmt_2 = '(a7,i3,a1,i5,a2,'//trim(fmt_2)//'i4,a1)'
        write (stdout, fmt_2) ' Group ', j, ' ', group_num_subgroups, ' (', (group_subgroups(i), i=1, group_num_subgroups), ')'
      endif
      !
      ! filling up subgroup_info
      !
      do k = 1, group_num_subgroups
        subgroup_info(j, k) = group_subgroups(k)
      enddo
      !
      !Convenient variable which will be amended later. Used to appropriately extract the subgroup array from the group_array
      !
      subgroup_increment = 1
      !
      !Loop over subgroups inside group
      !
      do k = 1, group_num_subgroups
        allocate (subgroup_array(2, group_subgroups(k)), stat=ierr)
        if (ierr /= 0) call io_error('Error in allocating subgroup_array in master_sort_and_group')
        allocate (sorted_subgroup_array(2, group_subgroups(k)), stat=ierr)
        if (ierr /= 0) call io_error('Error in allocating sorted_subgroup_array in master_sort_and_group')
        !
        !Extract the subgroup from the group
        !
        subgroup_array = sorted_group_array(:, subgroup_increment:subgroup_increment + group_subgroups(k) - 1)
        !
        !Updating subgroup_array to contain coord(3)
        !
        do i = 1, group_subgroups(k)
          subgroup_array(2, i) = wannier_centres_translated(coord(3), int(subgroup_array(1, i)))
        enddo

        call sort(subgroup_array, sorted_subgroup_array)
        !
        !Update sorted_group array with the sorted subgroup array
        !
        sorted_group_array(:, subgroup_increment:subgroup_increment + group_subgroups(k) - 1) = sorted_subgroup_array
        !
        !Update the subgroup_increment
        !
        subgroup_increment = subgroup_increment + group_subgroups(k)
        deallocate (sorted_subgroup_array, stat=ierr)
        if (ierr /= 0) call io_error('Error deallocating sorted_subgroup_array in master_sort_and_group')
        deallocate (subgroup_array, stat=ierr)
        if (ierr /= 0) call io_error('Error deallocating subgroup_array in master_sort_and_group')
      enddo
      !
      !Update Array with the sorted group array
      !
      Array(:, increment:increment + Array_groups(j) - 1) = sorted_group_array
      !
      !Update the group increment
      !
      increment = increment + Array_groups(j)
      deallocate (group_array, stat=ierr)
      if (ierr /= 0) call io_error('Error deallocating group_array in master_sort_and_group')
      deallocate (sorted_group_array, stat=ierr)
      if (ierr /= 0) call io_error('Error deallocating sorted_group_array in master_sort_and_group')
      deallocate (group_subgroups, stat=ierr)
      if (ierr /= 0) call io_error('Error deallocating group_subgroups in master_sort_and_group')
    enddo

    if (timing_level > 2) call io_stopwatch('tran: lcr_2c2_sort: master_sort', 2)

    return

  end subroutine master_sort_and_group