clone_gdata Subroutine

private subroutine clone_gdata(this, gdata_buf, tmp_fvals_size, fvals_size, tmp_apvals_size, apvals_size, nsigns)

Type Bound

gdata_io_t

Arguments

Type IntentOptional Attributes Name
class(gdata_io_t), intent(inout) :: this
integer(kind=hsize_t), intent(inout), allocatable :: gdata_buf(:,:)
integer, intent(in) :: tmp_fvals_size
integer, intent(in) :: fvals_size
integer, intent(in) :: tmp_apvals_size
integer, intent(in) :: apvals_size
integer(kind=int64), intent(in) :: nsigns

Contents

Source Code


Source Code

    subroutine clone_gdata(this, gdata_buf, tmp_fvals_size, fvals_size, tmp_apvals_size, &
                           apvals_size, nsigns)
        ! expand the global det data:
        ! clone the fvals from tmp_fvals_size to fvals_size, leaving the rest of
        ! the data as it is
        ! Input: gdata_buf - On input, the read-in gdata sized with tmp_fvals_size
        !                    on return, the re-sized gdata sized with fvals_size
        !        tmp_fvals_size - size of the acc. rate part of gdata on input
        !        fvals_size - size of the acc. rate part of gdata on return
        !        nsigns - number of determinants affected
        class(gdata_io_t), intent(inout) :: this
        integer(hsize_t), allocatable, intent(inout) :: gdata_buf(:, :)
        integer, intent(in) :: tmp_fvals_size, fvals_size, tmp_apvals_size, apvals_size
        integer(int64), intent(in) :: nsigns

        integer(hsize_t), allocatable :: tmp_fvals_acc(:, :), tmp_fvals_tot(:, :), tmp_mr(:, :)
        integer(hsize_t), allocatable :: tmp_apvals_sum(:, :), tmp_apvals_iter(:, :)
        integer :: max_ratio_size, gdata_size
        logical :: t_aas, t_sb, t_ap
        integer :: tmp_tot_start, tot_start, tmp_acc_size, acc_size
        integer :: tmp_iter_start, iter_start, tmp_sum_size, sum_size

        ! if size of data is above 0, the option has been set and memory is reserved
        t_aas = this%fvals_end - this%fvals_start + 1 > 0
        max_ratio_size = this%max_ratio_end - this%max_ratio_start + 1
        t_sb = max_ratio_size > 0
        t_ap = this%apvals_end - this%apvals_start + 1 > 0

        ! Only the sizes fvals and apvals depend on lenof_sign
        if (t_aas .or. t_ap) then
            ! can only resize buffers with correct size
            if (this%entry_size() == size(gdata_buf, dim=1)) then
                if (t_aas) then
                    ! copy the fvals to a temporary - one for tot. and one for acc. spawns
                    ! each half the total size
                    ! we do this to clone each of them independently
                    tmp_acc_size = tmp_fvals_size.div.2
                    tmp_tot_start = this%fvals_start + tmp_acc_size
                    allocate(tmp_fvals_acc(tmp_acc_size, nsigns))
                    allocate(tmp_fvals_tot(tmp_acc_size, nsigns))
                    tmp_fvals_acc(:, :) = gdata_buf(this%fvals_start:tmp_tot_start - 1, :)
                    tmp_fvals_tot(:, :) = gdata_buf(tmp_tot_start:this%fvals_end, :)

                    acc_size = fvals_size.div.2
                    ! clone the content of the temporary - clone the first
                    ! and second half seperately
                    call clone_signs(tmp_fvals_acc, tmp_acc_size, acc_size, nsigns)
                    ! tot and acc have to have the same size
                    call clone_signs(tmp_fvals_tot, tmp_acc_size, acc_size, nsigns)
                end if
                if (t_sb) then
                    allocate(tmp_mr(max_ratio_size, nsigns))
                    tmp_mr(:, :) = gdata_buf(this%max_ratio_start:this%max_ratio_end, :)
                end if
                if (t_ap) then
                    ! copy the apvals (pops sum and pops iter) to a temporary -
                    ! only the size of pops sum dependes on lenof_sign
                    tmp_sum_size = tmp_apvals_size - 1
                    tmp_iter_start = this%apvals_start + tmp_sum_size
                    allocate(tmp_apvals_sum(tmp_sum_size, nsigns))
                    allocate(tmp_apvals_iter(1, nsigns))
                    tmp_apvals_sum(:, :) = gdata_buf(this%apvals_start:tmp_iter_start - 1, :)
                    tmp_apvals_iter(:, :) = gdata_buf(tmp_iter_start:this%apvals_end, :)

                    sum_size = apvals_size - 1
                    ! clone the content of the temporary pops sum
                    call clone_signs(tmp_apvals_sum, tmp_sum_size, sum_size, nsigns)
                end if
                deallocate(gdata_buf)
                ! adjust the gdata offsets of this io handler
                call this%init_gdata_io(t_aas, t_sb, t_ap, fvals_size, max_ratio_size, apvals_size)

                ! resize the buffer - with the new gdata_size
                gdata_size = this%entry_size()
                allocate(gdata_buf(gdata_size, nsigns))

                if (t_aas) then
                    tot_start = this%fvals_start + acc_size
                    ! fill in the resized data
                    gdata_buf(this%fvals_start:tot_start - 1, :) = tmp_fvals_acc(:, :)
                    gdata_buf(tot_start:this%fvals_end, :) = tmp_fvals_tot(:, :)
                    deallocate(tmp_fvals_tot)
                    deallocate(tmp_fvals_acc)
                end if
                if (t_sb) then
                    gdata_buf(this%max_ratio_start:this%max_ratio_end, :) = tmp_mr(:, :)
                    deallocate(tmp_mr)
                end if
                if (t_ap) then
                    iter_start = this%apvals_start + sum_size
                    ! fill in the resized data
                    gdata_buf(this%apvals_start:iter_start - 1, :) = tmp_apvals_sum(:, :)
                    gdata_buf(iter_start:this%apvals_end, :) = tmp_apvals_iter(:, :)
                    deallocate(tmp_apvals_sum)
                    deallocate(tmp_apvals_iter)

                end if
            else
                write(stderr, *) "WARNING: Dimension mismatch in clone_gdata. No data read"
            end if
        end if
    end subroutine clone_gdata