Line data Source code
1 : module m_distrib_vx
2 : use m_types
3 : #ifdef CPP_MPI
4 : use mpi
5 : #endif
6 : use m_types_mpimat
7 : use m_glob_tofrom_loc
8 : contains
9 12 : subroutine distrib_vx(fi, fmpi, nococonv, vx_loc, vx_tmp, hybdat)
10 : implicit none
11 : type(t_fleurinput), intent(in) :: fi
12 : type(t_nococonv), intent(in) :: nococonv
13 : TYPE(t_mpi), INTENT(IN) :: fmpi
14 : integer, intent(in) :: vx_loc(:,:)
15 : type(t_mat), intent(inout) :: vx_tmp(:,:)
16 : type(t_hybdat), intent(inout) :: hybdat
17 :
18 :
19 : integer :: jsp, nk
20 :
21 12 : call timestart("distrib_vx")
22 :
23 12 : if(.not. allocated(hybdat%v_x)) then
24 6 : IF (fmpi%n_size == 1) THEN
25 0 : ALLOCATE (t_mat::hybdat%v_x(fi%kpts%nkpt, fi%input%jspins))
26 : ELSE
27 56 : ALLOCATE (t_mpimat::hybdat%v_x(fi%kpts%nkpt, fi%input%jspins))
28 : END IF
29 : endif
30 :
31 28 : do jsp = 1,fi%input%jspins
32 76 : do nk = 1,fi%kpts%nkpt
33 64 : call distrib_single_vx(fi, fmpi, jsp, nk, vx_loc(nk,jsp), vx_tmp(nk,jsp), hybdat, nococonv=nococonv)
34 : enddo
35 : enddo
36 :
37 12 : call timestop("distrib_vx")
38 12 : end subroutine distrib_vx
39 :
40 48 : subroutine distrib_single_vx(fi, fmpi, jsp, nk, vx_loc, vx_tmp, hybdat, dims, nococonv)
41 : implicit none
42 : type(t_fleurinput), intent(in) :: fi
43 : TYPE(t_mpi), INTENT(IN) :: fmpi
44 : integer, intent(in) :: vx_loc, jsp, nk
45 : type(t_mat), intent(inout) :: vx_tmp
46 : type(t_hybdat), intent(inout) :: hybdat
47 : integer, optional, intent(in) :: dims(2)
48 : type(t_nococonv), intent(in), optional :: nococonv
49 :
50 48 : type(t_lapw) :: lapw
51 : integer :: mat_sz, recver, ierr
52 :
53 48 : call timestart("distrib_vx_single")
54 :
55 48 : if(present(nococonv) .and. (.not. present(dims))) then
56 48 : CALL lapw%init(fi, nococonv, nk)
57 48 : mat_sz = lapw%nv(jsp) + fi%atoms%nlotot
58 0 : elseif((.not. present(nococonv)) .and. present(dims)) then
59 0 : mat_sz = dims(1)
60 : else
61 0 : call juDFT_error("I need either nococonv or dims")
62 : endif
63 :
64 96 : if(any(nk == fmpi%k_list)) then
65 48 : CALL hybdat%v_x(nk, jsp)%init(fi%sym%invs, mat_sz, mat_sz, fmpi%sub_comm, .false.)
66 : else
67 0 : call hybdat%v_x(nk, jsp)%init(fi%sym%invs, 1,1, fmpi%sub_comm, .false.)
68 : endif
69 :
70 48 : call copy_vx_to_distr(fmpi, vx_loc, nk, mat_sz, vx_tmp, hybdat%v_x(nk, jsp))
71 48 : call vx_tmp%free()
72 :
73 48 : call timestop("distrib_vx_single")
74 48 : end subroutine distrib_single_vx
75 :
76 48 : subroutine copy_vx_to_distr(fmpi, sender, nk, mat_sz, vx_den, vx_distr)
77 : implicit none
78 : type(t_mpi), intent(in) :: fmpi
79 : integer, intent(in) :: sender, nk, mat_sz
80 : type(t_mat), intent(in) :: vx_den
81 : class(t_mat), intent(inout) :: vx_distr
82 :
83 : integer :: i_loc, pe_i, i, ierr, recver
84 :
85 48 : pe_i = 0
86 48 : i_loc = 1
87 7256 : do i = 1, mat_sz
88 7208 : call glob_to_loc(fmpi, i, pe_i, i_loc)
89 14424 : if(any(fmpi%k_list == nk) .and. pe_i == fmpi%n_rank) then
90 3604 : recver = fmpi%irank
91 : else
92 3604 : recver = -1
93 : endif
94 : #ifdef CPP_MPI
95 7208 : call MPI_Allreduce(MPI_IN_PLACE, recver, 1, MPI_INTEGER, MPI_MAX, fmpi%mpi_comm, ierr)
96 : #endif
97 :
98 7256 : if(fmpi%irank == sender .and. sender == recver) then
99 1806 : if(vx_den%l_real) then
100 231144 : vx_distr%data_r(:,i_loc) = vx_den%data_r(:,i)
101 : else
102 98302 : vx_distr%data_c(:,i_loc) = vx_den%data_c(:,i)
103 : endif
104 : #ifdef CPP_MPI
105 5402 : elseif(fmpi%irank == sender) then
106 1798 : if(vx_den%l_real) then
107 1258 : call MPI_Send(vx_den%data_r(:,i), vx_den%matsize1, MPI_DOUBLE_PRECISION, recver, i, fmpi%mpi_comm, ierr)
108 : else
109 540 : call MPI_Send(vx_den%data_c(:,i), vx_den%matsize1, MPI_DOUBLE_COMPLEX, recver, i, fmpi%mpi_comm, ierr)
110 : endif
111 3604 : elseif(fmpi%irank == recver) then
112 1798 : if(vx_distr%l_real) then
113 1258 : call MPI_Recv(vx_distr%data_r(:,i_loc), vx_distr%matsize1, MPI_DOUBLE_PRECISION, sender, i, fmpi%mpi_comm, MPI_STATUS_IGNORE, ierr)
114 : else
115 540 : call MPI_Recv(vx_distr%data_c(:,i_loc), vx_distr%matsize1, MPI_DOUBLE_COMPLEX, sender, i, fmpi%mpi_comm, MPI_STATUS_IGNORE, ierr)
116 : endif
117 : #endif
118 : endif
119 : enddo
120 48 : end subroutine copy_vx_to_distr
121 : end module m_distrib_vx
|