param_chkpt_dist Subroutine

public subroutine param_chkpt_dist()

Uses

  • proc~~param_chkpt_dist~~UsesGraph proc~param_chkpt_dist param_chkpt_dist module~w90_constants w90_constants proc~param_chkpt_dist->module~w90_constants module~w90_comms w90_comms proc~param_chkpt_dist->module~w90_comms module~w90_io w90_io proc~param_chkpt_dist->module~w90_io module~w90_comms->module~w90_constants module~w90_comms->module~w90_io module~w90_io->module~w90_constants

Distribute the chk files

Arguments

None

Calls

proc~~param_chkpt_dist~~CallsGraph proc~param_chkpt_dist param_chkpt_dist interface~comms_bcast comms_bcast proc~param_chkpt_dist->interface~comms_bcast proc~comms_bcast_char comms_bcast_char interface~comms_bcast->proc~comms_bcast_char proc~comms_bcast_cmplx comms_bcast_cmplx interface~comms_bcast->proc~comms_bcast_cmplx proc~comms_bcast_int comms_bcast_int interface~comms_bcast->proc~comms_bcast_int proc~comms_bcast_real comms_bcast_real interface~comms_bcast->proc~comms_bcast_real proc~comms_bcast_logical comms_bcast_logical interface~comms_bcast->proc~comms_bcast_logical

Called by

proc~~param_chkpt_dist~~CalledByGraph proc~param_chkpt_dist param_chkpt_dist program~wannier wannier program~wannier->proc~param_chkpt_dist

Contents

Source Code


Source Code

  subroutine param_chkpt_dist
    !===========================================================!
    !                                                           !
    !! Distribute the chk files
    !                                                           !
    !===========================================================!

    use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
    use w90_io, only: io_error, io_file_unit, &
      io_date, io_time, io_stopwatch
    use w90_comms, only: on_root, comms_bcast

    implicit none

    integer :: ierr, loop_kpt, m, i, j

    call comms_bcast(checkpoint, len(checkpoint))

    if (.not. on_root .and. .not. allocated(u_matrix)) then
      allocate (u_matrix(num_wann, num_wann, num_kpts), stat=ierr)
      if (ierr /= 0) &
        call io_error('Error allocating u_matrix in param_chkpt_dist')
    endif
    call comms_bcast(u_matrix(1, 1, 1), num_wann*num_wann*num_kpts)

!    if (.not.on_root .and. .not.allocated(m_matrix)) then
!       allocate(m_matrix(num_wann,num_wann,nntot,num_kpts),stat=ierr)
!       if (ierr/=0)&
!            call io_error('Error allocating m_matrix in param_chkpt_dist')
!    endif
!    call comms_bcast(m_matrix(1,1,1,1),num_wann*num_wann*nntot*num_kpts)

    call comms_bcast(have_disentangled, 1)

    if (have_disentangled) then
      if (.not. on_root) then

        if (.not. allocated(u_matrix_opt)) then
          allocate (u_matrix_opt(num_bands, num_wann, num_kpts), stat=ierr)
          if (ierr /= 0) &
            call io_error('Error allocating u_matrix_opt in param_chkpt_dist')
        endif

        if (.not. allocated(lwindow)) then
          allocate (lwindow(num_bands, num_kpts), stat=ierr)
          if (ierr /= 0) &
            call io_error('Error allocating lwindow in param_chkpt_dist')
        endif

        if (.not. allocated(ndimwin)) then
          allocate (ndimwin(num_kpts), stat=ierr)
          if (ierr /= 0) &
            call io_error('Error allocating ndimwin in param_chkpt_dist')
        endif

      end if

      call comms_bcast(u_matrix_opt(1, 1, 1), num_bands*num_wann*num_kpts)
      call comms_bcast(lwindow(1, 1), num_bands*num_kpts)
      call comms_bcast(ndimwin(1), num_kpts)
      call comms_bcast(omega_invariant, 1)
    end if
    call comms_bcast(wannier_centres(1, 1), 3*num_wann)
    call comms_bcast(wannier_spreads(1), num_wann)

  end subroutine param_chkpt_dist