TABLE OF CONTENTS
- abinit/m_pawcprj
- m_pawcprj/paw_overlap
- m_pawcprj/pawcprj_alloc
- m_pawcprj/pawcprj_axpby
- m_pawcprj/pawcprj_bcast
- m_pawcprj/pawcprj_conjg
- m_pawcprj/pawcprj_copy
- m_pawcprj/pawcprj_free
- m_pawcprj/pawcprj_gather_spin
- m_pawcprj/pawcprj_get
- m_pawcprj/pawcprj_getdim
- m_pawcprj/pawcprj_lincom
- m_pawcprj/pawcprj_mpi_allgather
- m_pawcprj/pawcprj_mpi_exch
- m_pawcprj/pawcprj_mpi_recv
- m_pawcprj/pawcprj_mpi_send
- m_pawcprj/pawcprj_mpi_sum
- m_pawcprj/pawcprj_output
- m_pawcprj/pawcprj_pack
- m_pawcprj/pawcprj_projbd
- m_pawcprj/pawcprj_put
- m_pawcprj/pawcprj_reorder
- m_pawcprj/pawcprj_set_zero
- m_pawcprj/pawcprj_symkn
- m_pawcprj/pawcprj_transpose
- m_pawcprj/pawcprj_type
- m_pawcprj/pawcprj_unpack
- m_pawcprj/pawcprj_zaxpby
abinit/m_pawcprj [ Modules ]
[ Top ] [ abinit ] [ Modules ]
NAME
m_pawcprj
FUNCTION
This module contains functions used to manipulate variables of structured datatype pawcprj_type. pawcprj_type variables are <p_lmn|Cnk> projected quantities, where |p_lmn> are non-local projectors |Cnk> are wave functions
COPYRIGHT
Copyright (C) 2012-2024 ABINIT group (MT,JWZ) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
NOTES
FOR DEVELOPPERS: in order to preserve the portability of libPAW library, please consult ~abinit/src/??_libpaw/libpaw-coding-rules.txt
SOURCE
24 #include "libpaw.h" 25 26 module m_pawcprj 27 28 USE_DEFS 29 USE_MSG_HANDLING 30 USE_MPI_WRAPPERS 31 USE_MEMORY_PROFILING 32 33 use m_pawtab, only : pawtab_type 34 35 implicit none 36 37 private
m_pawcprj/paw_overlap [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
paw_overlap
FUNCTION
Helper function returning the onsite contribution to the overlap between two states.
INPUTS
spinor_comm= (optional) communicator over spinorial components typat(:)=The type of each atom. Pawtab(ntypat)<type(pawtab_type)>=paw tabulated starting data. cprj1,cprj2<pawcprj_type> Projected wave functions <Proj_i|Cnk> with all NL projectors for the left and the right wavefunction,respectively.
OUTPUT
SOURCE
2745 function paw_overlap(cprj1,cprj2,typat,pawtab,spinor_comm) result(onsite) 2746 2747 !Arguments ------------------------------------ 2748 !scalars 2749 integer,intent(in),optional :: spinor_comm 2750 !arrays 2751 integer,intent(in) :: typat(:) 2752 real(dp) :: onsite(2) 2753 type(pawcprj_type),intent(in) :: cprj1(:,:),cprj2(:,:) 2754 type(pawtab_type),intent(in) :: pawtab(:) 2755 2756 !Local variables------------------------------- 2757 !scalars 2758 integer :: iatom,ilmn,itypat,j0lmn,jlmn,klmn,natom,nspinor,isp 2759 real(dp) :: sij 2760 character(len=500) :: msg 2761 !arrays 2762 2763 ! ************************************************************************* 2764 2765 natom=SIZE(typat) 2766 2767 if (SIZE(cprj1,DIM=1)/=SIZE(cprj2,DIM=1) .or. SIZE(cprj1,DIM=1)/=natom) then 2768 write(msg,'(a,3i4)')' Wrong size in typat, cprj1, cprj2 : ',natom,SIZE(cprj1),SIZE(cprj2) 2769 LIBPAW_ERROR(msg) 2770 end if 2771 2772 nspinor = SIZE(cprj1,DIM=2) 2773 2774 onsite=zero 2775 do iatom=1,natom 2776 itypat=typat(iatom) 2777 do jlmn=1,pawtab(itypat)%lmn_size 2778 j0lmn=jlmn*(jlmn-1)/2 2779 do ilmn=1,jlmn 2780 klmn=j0lmn+ilmn 2781 sij=pawtab(itypat)%sij(klmn); if (jlmn==ilmn) sij=sij*half 2782 if (ABS(sij)>tol16) then 2783 do isp=1,nspinor 2784 2785 onsite(1)=onsite(1) + sij*( & 2786 & cprj1(iatom,isp)%cp(1,ilmn) * cprj2(iatom,isp)%cp(1,jlmn) & 2787 & +cprj1(iatom,isp)%cp(2,ilmn) * cprj2(iatom,isp)%cp(2,jlmn) & 2788 & +cprj1(iatom,isp)%cp(1,jlmn) * cprj2(iatom,isp)%cp(1,ilmn) & 2789 & +cprj1(iatom,isp)%cp(2,jlmn) * cprj2(iatom,isp)%cp(2,ilmn) & 2790 & ) 2791 2792 onsite(2)=onsite(2) + sij*( & 2793 & cprj1(iatom,isp)%cp(1,ilmn) * cprj2(iatom,isp)%cp(2,jlmn) & 2794 & -cprj1(iatom,isp)%cp(2,ilmn) * cprj2(iatom,isp)%cp(1,jlmn) & 2795 & +cprj1(iatom,isp)%cp(1,jlmn) * cprj2(iatom,isp)%cp(2,ilmn) & 2796 & -cprj1(iatom,isp)%cp(2,jlmn) * cprj2(iatom,isp)%cp(1,ilmn) & 2797 & ) 2798 end do 2799 end if 2800 end do 2801 end do 2802 end do 2803 2804 if (present(spinor_comm)) then 2805 call xmpi_sum(onsite,spinor_comm,isp) 2806 end if 2807 2808 end function paw_overlap
m_pawcprj/pawcprj_alloc [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_alloc
FUNCTION
Allocation of a cprj datastructure
INPUTS
ncpgr=number of gradients to be allocated nlmn(:)=sizes of cprj%cp
SIDE EFFECTS
cprj(:,:) <type(pawcprj_type)>= cprj datastructure
SOURCE
124 subroutine pawcprj_alloc(cprj,ncpgr,nlmn) 125 126 !Arguments ------------------------------------ 127 !scalars 128 integer,intent(in) :: ncpgr 129 !arrays 130 integer,intent(in) :: nlmn(:) 131 type(pawcprj_type),intent(inout) :: cprj(:,:) 132 133 !Local variables------------------------------- 134 !scalars 135 integer :: ii,jj,n1dim,n2dim,nn 136 character(len=500) :: msg 137 138 ! ************************************************************************* 139 140 n1dim=size(cprj,dim=1);n2dim=size(cprj,dim=2);nn=size(nlmn,dim=1) 141 if (nn/=n1dim) then 142 write(msg,*) 'wrong sizes (pawcprj_alloc)! :',nn,n1dim 143 LIBPAW_ERROR(msg) 144 end if 145 146 do jj=1,n2dim 147 do ii=1,n1dim 148 if (allocated(cprj(ii,jj)%cp)) then 149 LIBPAW_DEALLOCATE(cprj(ii,jj)%cp) 150 end if 151 if (allocated(cprj(ii,jj)%dcp)) then 152 LIBPAW_DEALLOCATE(cprj(ii,jj)%dcp) 153 end if 154 nn=nlmn(ii) 155 cprj(ii,jj)%nlmn=nn 156 LIBPAW_ALLOCATE(cprj(ii,jj)%cp,(2,nn)) 157 cprj(ii,jj)%cp=zero 158 cprj(ii,jj)%ncpgr=ncpgr 159 if (ncpgr>0) then 160 LIBPAW_ALLOCATE(cprj(ii,jj)%dcp,(2,ncpgr,nn)) 161 cprj(ii,jj)%dcp=zero 162 end if 163 end do 164 end do 165 166 end subroutine pawcprj_alloc
m_pawcprj/pawcprj_axpby [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_axpby
FUNCTION
Apply AXPBY (blas-like) operation with 2 cprj datastructures: cprjy(:,:) <- alpha.cprjx(:,:)+beta.cprjy(:,:) alpha and beta are REAL scalars
INPUTS
alpha,beta= alpha,beta REAL factors cprjx(:,:) <type(pawcprj_type)>= input cprjx datastructure
SIDE EFFECTS
cprjy(:,:) <type(pawcprj_type)>= input/output cprjy datastructure
SOURCE
376 subroutine pawcprj_axpby(alpha,beta,cprjx,cprjy) 377 378 !Arguments ------------------------------------ 379 !scalars 380 real(dp),intent(in) :: alpha,beta 381 !arrays 382 type(pawcprj_type),intent(in) :: cprjx(:,:) 383 type(pawcprj_type),intent(inout) :: cprjy(:,:) 384 385 !Local variables------------------------------- 386 !scalars 387 integer :: ii,jj,kk,n1dimx,n1dimy,n2dimx,n2dimy,ncpgrx,ncpgry,nlmn 388 character(len=500) :: msg 389 390 ! ************************************************************************* 391 392 n1dimy=size(cprjy,dim=1);n2dimy=size(cprjy,dim=2);ncpgry=cprjy(1,1)%ncpgr 393 if (abs(alpha)>tol16) then 394 n1dimx=size(cprjx,dim=1);n2dimx=size(cprjx,dim=2);ncpgrx=cprjx(1,1)%ncpgr 395 msg = "" 396 if (n1dimx/=n1dimy) msg = TRIM(msg)//"Error in pawcprj_axpby: n1 wrong sizes !"//ch10 397 if (n2dimx/=n2dimy) msg = TRIM(msg)//"Error in pawcprj_axpby: n2 wrong sizes !"//ch10 398 if (ncpgrx/=ncpgry) msg = TRIM(msg)//"Error in pawcprj_axpby: ncpgr wrong sizes !"//ch10 399 if (LEN_TRIM(msg) > 0) then 400 LIBPAW_ERROR(msg) 401 end if 402 else 403 n1dimx=0;n2dimx=0;ncpgrx=0 404 end if 405 406 if (abs(alpha)<=tol16) then 407 do jj=1,n2dimy 408 do ii=1,n1dimy 409 nlmn=cprjy(ii,jj)%nlmn 410 do kk=1,nlmn 411 cprjy(ii,jj)%cp(1:2,kk)=beta*cprjy(ii,jj)%cp(1:2,kk) 412 end do 413 end do 414 end do 415 if (ncpgry>0) then 416 do jj=1,n2dimy 417 do ii=1,n1dimy 418 nlmn=cprjy(ii,jj)%nlmn 419 do kk=1,nlmn 420 cprjy(ii,jj)%dcp(1:2,1:ncpgry,kk)=beta*cprjy(ii,jj)%dcp(1:2,1:ncpgry,kk) 421 end do 422 end do 423 end do 424 end if 425 else if (abs(beta)<=tol16) then 426 do jj=1,n2dimx 427 do ii=1,n1dimx 428 nlmn=cprjx(ii,jj)%nlmn 429 cprjy(ii,jj)%nlmn=nlmn 430 do kk=1,nlmn 431 cprjy(ii,jj)%cp(1:2,kk)=alpha*cprjx(ii,jj)%cp(1:2,kk) 432 end do 433 end do 434 end do 435 if (ncpgrx>0) then 436 do jj=1,n2dimx 437 do ii=1,n1dimx 438 nlmn=cprjx(ii,jj)%nlmn 439 do kk=1,nlmn 440 cprjy(ii,jj)%dcp(1:2,1:ncpgrx,kk)=alpha*cprjx(ii,jj)%dcp(1:2,1:ncpgrx,kk) 441 end do 442 end do 443 end do 444 end if 445 else ! alpha/=0 and beta/=0 446 do jj=1,n2dimx 447 do ii=1,n1dimx 448 nlmn=cprjx(ii,jj)%nlmn 449 cprjy(ii,jj)%nlmn=nlmn 450 do kk=1,nlmn 451 cprjy(ii,jj)%cp(1:2,kk)=alpha*cprjx(ii,jj)%cp(1:2,kk) & 452 & +beta *cprjy(ii,jj)%cp(1:2,kk) 453 end do 454 end do 455 end do 456 if (ncpgrx>0) then 457 do jj=1,n2dimx 458 do ii=1,n1dimx 459 nlmn=cprjx(ii,jj)%nlmn 460 do kk=1,nlmn 461 cprjy(ii,jj)%dcp(1:2,1:ncpgrx,kk)=alpha*cprjx(ii,jj)%dcp(1:2,1:ncpgrx,kk) & 462 & +beta *cprjy(ii,jj)%dcp(1:2,1:ncpgrx,kk) 463 end do 464 end do 465 end do 466 end if 467 end if 468 469 end subroutine pawcprj_axpby
m_pawcprj/pawcprj_bcast [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_bcast
FUNCTION
Broadcast a pawcprj_type from master to all nodes inside a MPI communicator.
INPUTS
natom=Number of atoms (size of the first dimension of Cprj). n2dim=Size of the second dimension of Cprj. ncpgr=Number of gradients that have to be cast. It is a bit redundant but, it can be used to broad cast only the %cp"s without caring about the gradients. Just set it to 0 but be careful! nlmn(natom)=Number of nlm partial waves for each atom. master=ID of the sending node in spaceComm. spaceComm=MPI Communicator.
OUTPUT
ierr=Error status. Cprj(natom,n2dim)<pawcprj_type>=The datatype to be transmitted by master and received by the others nodes.
SOURCE
2215 subroutine pawcprj_bcast(Cprj,natom,n2dim,nlmn,ncpgr,master,spaceComm,ierr) 2216 2217 !Arguments ------------------------------------ 2218 !scalars 2219 integer,intent(in) :: natom,n2dim,ncpgr,master,spaceComm 2220 integer,intent(out) :: ierr 2221 !arrays 2222 integer,intent(in) :: nlmn(natom) 2223 type(pawcprj_type),intent(inout) :: Cprj(natom,n2dim) 2224 2225 !Local variables------------------------------- 2226 !scalars 2227 integer :: iat,jj,n1dim,nn 2228 integer :: ntotcp,ipck,rank,nprocs 2229 character(len=100) :: msg 2230 !arrays 2231 real(dp),allocatable :: buffer_cp(:,:),buffer_cpgr(:,:,:) 2232 2233 ! ************************************************************************* 2234 2235 ierr=0 2236 nprocs = xmpi_comm_size(spaceComm) 2237 if (nprocs==1) return 2238 2239 rank = xmpi_comm_rank(spaceComm) 2240 2241 nn=size(nlmn,dim=1) 2242 n1dim=size(Cprj,dim=1) 2243 if (nn/=n1dim) then 2244 msg='size mismatch in natom (pawcprj_bcast)!' 2245 LIBPAW_BUG(msg) 2246 end if 2247 2248 ntotcp=n2dim*SUM(nlmn(:)) 2249 2250 LIBPAW_ALLOCATE(buffer_cp,(2,ntotcp)) 2251 if (ncpgr/=0) then 2252 LIBPAW_ALLOCATE(buffer_cpgr,(2,ncpgr,ntotcp)) 2253 end if 2254 2255 !=== Master packs Cprj === 2256 !Write a routine to pack/unpack? 2257 if (rank==master) then 2258 ipck=0 2259 do jj=1,n2dim 2260 do iat=1,natom 2261 nn=nlmn(iat) 2262 buffer_cp(:,ipck+1:ipck+nn)=Cprj(iat,jj)%cp(:,1:nn) 2263 if (ncpgr/=0) buffer_cpgr(:,:,ipck+1:ipck+nn)=Cprj(iat,jj)%dcp(:,:,1:nn) 2264 ipck=ipck+nn 2265 end do 2266 end do 2267 end if 2268 2269 !=== Transmit data === 2270 call xmpi_bcast(buffer_cp,master,spaceComm,ierr) 2271 if (ncpgr/=0) then 2272 call xmpi_bcast(buffer_cpgr,master,spaceComm,ierr) 2273 end if 2274 2275 !=== UnPack the received buffer === 2276 if (rank/=master) then 2277 ipck=0 2278 do jj=1,n2dim 2279 do iat=1,natom 2280 nn=nlmn(iat) 2281 Cprj(iat,jj)%cp(:,1:nn)=buffer_cp(:,ipck+1:ipck+nn) 2282 if (ncpgr/=0) Cprj(iat,jj)%dcp(:,:,1:nn)=buffer_cpgr(:,:,ipck+1:ipck+nn) 2283 ipck=ipck+nn 2284 end do 2285 end do 2286 end if 2287 2288 LIBPAW_DEALLOCATE(buffer_cp) 2289 if (ncpgr/=0) then 2290 LIBPAW_DEALLOCATE(buffer_cpgr) 2291 end if 2292 2293 end subroutine pawcprj_bcast
m_pawcprj/pawcprj_conjg [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_conjg
FUNCTION
conjugate a cprj datastructures: cprj(:,:) <- conjugate(cprj(:,:))
INPUTS
SIDE EFFECTS
cprj(:,:) <type(pawcprj_type)>= input/output cprj datastructure
SOURCE
905 subroutine pawcprj_conjg(cprj) 906 907 !Arguments ------------------------------------ 908 !scalars 909 !arrays 910 type(pawcprj_type),intent(inout) :: cprj(:,:) 911 912 !Local variables------------------------------- 913 !scalars 914 integer :: ii,jj,kk,n1dim,n2dim,ncpgr,nlmn 915 916 ! ************************************************************************* 917 918 919 n1dim=size(cprj,dim=1);n2dim=size(cprj,dim=2);ncpgr=cprj(1,1)%ncpgr 920 921 do jj=1,n2dim 922 do ii=1,n1dim 923 nlmn=cprj(ii,jj)%nlmn 924 do kk=1,nlmn 925 cprj(ii,jj)%cp(2,kk)=-cprj(ii,jj)%cp(2,kk) 926 end do 927 end do 928 end do 929 if (ncpgr>0) then 930 do jj=1,n2dim 931 do ii=1,n1dim 932 nlmn=cprj(ii,jj)%nlmn 933 do kk=1,nlmn 934 cprj(ii,jj)%dcp(2,1:ncpgr,kk)=-cprj(ii,jj)%dcp(2,1:ncpgr,kk) 935 end do 936 end do 937 end do 938 end if 939 940 end subroutine pawcprj_conjg
m_pawcprj/pawcprj_copy [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_copy
FUNCTION
Copy a cprj datastructure into another
INPUTS
icpgr= (optional argument) if present, only component icpgr of input cprj gradient is copied into output cprj Not used if cprj(:,:)%ncpgr<icpgr -1 only copy cp cprj_in(:,:) <type(pawcprj_type)>= input cprj datastructure
OUTPUT
cprj_out(:,:) <type(pawcprj_type)>= output cprj datastructure
NOTES
MG: What about an option to report a pointer to cprj_in?
SOURCE
275 subroutine pawcprj_copy(cprj_in,cprj_out,& 276 & icpgr) ! optional argument 277 278 !Arguments ------------------------------------ 279 !scalars 280 integer,intent(in),optional :: icpgr 281 !arrays 282 type(pawcprj_type),intent(in) :: cprj_in(:,:) 283 type(pawcprj_type),intent(inout) :: cprj_out(:,:) 284 285 !Local variables------------------------------- 286 !scalars 287 integer :: ii,jj,kk,n1dim_in,n1dim_out,n2dim_in,n2dim_out,ncpgr_in,ncpgr_out,nlmn 288 logical :: has_icpgr,copy_dcp 289 character(len=500) :: msg 290 291 ! ************************************************************************* 292 293 n1dim_in=size(cprj_in,dim=1); n1dim_out=size(cprj_out,dim=1) 294 n2dim_in=size(cprj_in,dim=2); n2dim_out=size(cprj_out,dim=2) 295 ncpgr_in=cprj_in(1,1)%ncpgr; ncpgr_out=cprj_out(1,1)%ncpgr 296 297 if (n1dim_in/=n1dim_out) then 298 write(msg,'(a,2(1x,i0))')" Error in pawcprj_copy: n1 wrong sizes ",n1dim_in,n1dim_out 299 LIBPAW_ERROR(msg) 300 end if 301 if (n2dim_in/=n2dim_out) then 302 write(msg,'(a,2(1x,i0))')" Error in pawcprj_copy: n2 wrong sizes ",n2dim_in,n2dim_out 303 LIBPAW_ERROR(msg) 304 end if 305 if (ncpgr_in<ncpgr_out) then 306 write(msg,'(a,2(1x,i0))')" Error in pawcprj_copy: ncpgr wrong sizes ",ncpgr_in,ncpgr_out 307 LIBPAW_ERROR(msg) 308 end if 309 310 !Check if icgr is present and if dcp have to be copy 311 has_icpgr=present(icpgr) 312 copy_dcp = .TRUE. 313 if(has_icpgr)then 314 copy_dcp = icpgr>=0 315 end if 316 317 do jj=1,n2dim_in 318 do ii=1,n1dim_in 319 nlmn=cprj_in(ii,jj)%nlmn 320 cprj_out(ii,jj)%nlmn =nlmn 321 do kk=1,nlmn 322 cprj_out(ii,jj)%cp(1:2,kk)=cprj_in(ii,jj)%cp(1:2,kk) 323 end do 324 end do 325 end do 326 327 if (ncpgr_in>0.and.copy_dcp) then 328 if (has_icpgr) has_icpgr=(ncpgr_out>0.and.icpgr>0.or.icpgr<=ncpgr_in) 329 330 if (has_icpgr) then 331 do jj=1,n2dim_in 332 do ii=1,n1dim_in 333 nlmn=cprj_in(ii,jj)%nlmn 334 do kk=1,nlmn 335 cprj_out(ii,jj)%dcp(1:2,1,kk)=cprj_in(ii,jj)%dcp(1:2,icpgr,kk) 336 end do 337 end do 338 end do 339 else 340 if (ncpgr_out>=ncpgr_in) then 341 do jj=1,n2dim_in 342 do ii=1,n1dim_in 343 nlmn=cprj_in(ii,jj)%nlmn 344 do kk=1,nlmn 345 cprj_out(ii,jj)%dcp(1:2,1:ncpgr_in,kk)=cprj_in(ii,jj)%dcp(1:2,1:ncpgr_in,kk) 346 end do 347 end do 348 end do 349 end if 350 end if 351 end if 352 353 end subroutine pawcprj_copy
m_pawcprj/pawcprj_free [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_free
FUNCTION
Deallocation of a cprj datastructure
SIDE EFFECTS
cprj(:,:) <type(pawcprj_type)>= cprj datastructure
SOURCE
183 subroutine pawcprj_free(cprj) 184 185 !Arguments ------------------------------------ 186 !scalars 187 !arrays 188 type(pawcprj_type),intent(inout) :: cprj(:,:) 189 190 !Local variables------------------------------- 191 !scalars 192 integer :: ii,jj,n1dim,n2dim 193 194 ! ************************************************************************* 195 196 n1dim=size(cprj,dim=1);n2dim=size(cprj,dim=2) 197 198 do jj=1,n2dim 199 do ii=1,n1dim 200 if (allocated(cprj(ii,jj)%cp)) then 201 LIBPAW_DEALLOCATE(cprj(ii,jj)%cp) 202 end if 203 if (allocated(cprj(ii,jj)%dcp)) then 204 LIBPAW_DEALLOCATE(cprj(ii,jj)%dcp) 205 end if 206 end do 207 end do 208 209 end subroutine pawcprj_free
m_pawcprj/pawcprj_gather_spin [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_gather_spin
FUNCTION
INPUTS
cprj(:,:)=the input cprj datastructure n2size=number of cprj datastructures to be gathered (second dim) nspinor : number of spinorial component (on current proc) nspinortot : total number of spinorial component
OUTPUT
cprj_gat(:,:) = the cprj containing all nspinor componants
NOTES
The cprj has been built like the following: loop on nsppol loop on k point loop over band or block of band These quantities were build only if treated by the current proc the inner quantities being nspinor
SOURCE
2571 subroutine pawcprj_gather_spin(cprj,cprj_gat,natom,n2size,nspinor,nspinortot,& 2572 & spaceComm_spin,ierr) 2573 2574 !Arguments ------------------------------------ 2575 !scalars 2576 integer,intent(in) :: natom,nspinor,nspinortot,n2size 2577 integer,intent(in) :: spaceComm_spin 2578 integer,intent(out) :: ierr 2579 !arrays 2580 type(pawcprj_type),intent(in) :: cprj(:,:) 2581 type(pawcprj_type),intent(inout) :: cprj_gat(:,:) 2582 2583 !Local variables------------------------------- 2584 !scalars 2585 integer :: i1,iatom,ibsp,icpgr,ilmn,isp,ispinor,jj,lmndim,n2dim,n2dim_gat,ncpgr 2586 character(len=100) :: msg 2587 !arrays 2588 integer :: nlmn(natom) 2589 real(dp),allocatable :: buffer1(:),buffer2(:) 2590 2591 ! ************************************************************************* 2592 2593 n2dim =size(cprj,dim=2) 2594 n2dim_gat=size(cprj_gat,dim=2) 2595 if (n2dim_gat/=(nspinortot/nspinor)*n2dim) then 2596 msg='wrong dims (pawcprj_gather_spin)!' 2597 LIBPAW_BUG(msg) 2598 end if 2599 2600 do iatom=1,natom 2601 nlmn(iatom)=size(cprj(iatom,1)%cp(1,:)) 2602 end do 2603 ncpgr=cprj(1,1)%ncpgr 2604 lmndim=2*n2size*sum(nlmn(1:natom))*(1+ncpgr) 2605 LIBPAW_ALLOCATE(buffer1,(lmndim)) 2606 LIBPAW_ALLOCATE(buffer2,(lmndim*nspinortot)) 2607 2608 isp=0;ibsp=0 2609 jj=1 2610 do i1=1,n2size 2611 isp=isp+1 2612 do iatom=1,natom 2613 do ilmn=1,nlmn(iatom) 2614 buffer1(jj:jj+1)=cprj(iatom,isp)%cp(1:2,ilmn) 2615 jj=jj+2 2616 end do 2617 if (ncpgr>0) then 2618 do ilmn=1,nlmn(iatom) 2619 do icpgr=1,ncpgr 2620 buffer1(jj:jj+1)=cprj(iatom,isp)%dcp(1:2,icpgr,ilmn) 2621 jj=jj+2 2622 end do 2623 end do 2624 end if 2625 end do 2626 end do 2627 2628 call xmpi_allgather(buffer1,lmndim,buffer2,spaceComm_spin,ierr) 2629 2630 jj=1 2631 do ispinor=1,nspinortot 2632 do i1 =1,n2size 2633 ibsp=(i1-1)*nspinortot + ispinor 2634 do iatom=1,natom 2635 do ilmn=1,nlmn(iatom) 2636 cprj_gat(iatom,ibsp)%cp(1:2,ilmn)=buffer2(jj:jj+1) 2637 jj=jj+2 2638 end do 2639 if (ncpgr>0) then 2640 do ilmn=1,nlmn(iatom) 2641 do icpgr=1,ncpgr 2642 cprj_gat(iatom,ibsp)%dcp(1:2,icpgr,ilmn)=buffer2(jj:jj+1) 2643 jj=jj+2 2644 end do 2645 end do 2646 end if 2647 end do 2648 end do 2649 end do 2650 2651 LIBPAW_DEALLOCATE(buffer1) 2652 LIBPAW_DEALLOCATE(buffer2) 2653 2654 end subroutine pawcprj_gather_spin
m_pawcprj/pawcprj_get [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_get
FUNCTION
Read the cprj_k for a given k-point from memory in cprj or from a temporary file
INPUTS
atind(natom)=index table for atoms (see iorder below) cprj(dimcp,nspinor*mband*mkmem*nsppol)=input cprj (used if mkmem/=0) dimcp=first dimension of cprj_k,cprj arrays (1 or natom) iband1=index of first band in cprj ibg=shift in cprj array to locate current k-point [icpgr]= (optional argument) if present, only component icpgr of input cprj gradient is copied into output cprj Not used if cprj(:,:)%ncpgr<icpgr (mkmem>0) or ncpgr(optional)<icpgr (mkmem=0) ikpt=index of current k-point (only needed for the parallel distribution) iorder=0 if cprj ordering does not change during reading 1 if cprj ordering changes during reading, depending on content of atind array: - if atind=atindx (type-sorted=>unsorted) - if atind=atindx1 (unsorted=>type-sorted) isppol=index of current spin component mband=maximum number of bands mkmem=number of k points which can fit in memory; set to 0 if use disk [mpi_comm]=(optional argument) MPI communicator over (k-pts,bands,spins) Must be used in association with proc_distrb argument natom=number of atoms in cell nband=number of bands to import (usually 1 or nband_k) nband_k=total number of bands for this k-point [ncpgr]=(optional argument) second dimension of cprj%dcp(2,ncpgr,nlmn stored in memory (mkmem>0) or present on disk (mkmem=0)) needed only when optional argument icpgr is present nspinor=number of spinorial components of the wavefunctions (on current proc) nsppol=1 for unpolarized, 2 for spin-polarized [proc_distrb(nkpt,nband,nsppol)]=(optional argument) processor distribution Describe how cprj datastructures are distributed over processors When present, mpicomm argument must be also present uncp=unit number for cprj data (used if mkmem=0)
OUTPUT
cprj_k(dimcp,nspinor*nband) <type(pawcprj_type)>= output cprj datastructure
SOURCE
1151 subroutine pawcprj_get(atind,cprj_k,cprj,dimcp,iband1,ibg,ikpt,iorder,isppol,mband,& 1152 & mkmem,natom,nband,nband_k,nspinor,nsppol,uncp,& 1153 & icpgr,ncpgr,mpicomm,proc_distrb) ! optionals arguments 1154 1155 !Arguments ------------------------------------ 1156 !scalars 1157 integer,intent(in) :: dimcp,iband1,ibg,ikpt,iorder,isppol,mband,mkmem,natom 1158 integer,intent(in) :: nband,nband_k,nspinor,nsppol,uncp 1159 integer,intent(in),optional :: icpgr,mpicomm,ncpgr 1160 !arrays 1161 integer,intent(in) :: atind(natom) 1162 integer,intent(in),optional :: proc_distrb(:,:,:) 1163 type(pawcprj_type),intent(in) :: cprj(dimcp,nspinor*mband*mkmem*nsppol) 1164 type(pawcprj_type),intent(inout) :: cprj_k(dimcp,nspinor*nband) 1165 1166 !Local variables------------------------------- 1167 !scalars 1168 integer :: iatm,iatom,ib,ibsp,icpgr_,isp,ispinor,jband,me,nband0,ncpgr_ 1169 logical :: has_distrb,has_icpgr 1170 character(len=500) :: msg 1171 !arrays 1172 real(dp),allocatable :: tmp(:,:,:) 1173 1174 ! ************************************************************************* 1175 1176 ncpgr_=cprj_k(1,1)%ncpgr;if (present(ncpgr)) ncpgr_=ncpgr 1177 icpgr_=-1;if(present(icpgr)) icpgr_=icpgr 1178 has_icpgr=(icpgr_>0.and.icpgr_<=ncpgr_) 1179 if (present(icpgr).and.(.not.present(ncpgr))) then 1180 msg='ncpgr must be present when icpgr is present (pawcprj_get)!' 1181 LIBPAW_BUG(msg) 1182 end if 1183 if (has_icpgr.and.cprj_k(1,1)%ncpgr<1) then 1184 msg='cprj_k%ncpgr not consistent with icpgr (pawcprj_get)!' 1185 LIBPAW_BUG(msg) 1186 end if 1187 1188 !MPI data 1189 has_distrb=present(proc_distrb) 1190 if (has_distrb) then 1191 if (.not.present(mpicomm)) then 1192 msg='mpicomm must be present when proc_distrb is present (pawcprj_get)!' 1193 LIBPAW_BUG(msg) 1194 end if 1195 me=xmpi_comm_rank(mpicomm) 1196 end if 1197 1198 if (mkmem==0) then 1199 1200 if (iband1==1) then 1201 read(uncp) nband0 1202 if (nband_k/=nband0) then 1203 msg='_PAW file was not created with the right options (pawcprj_get)!' 1204 LIBPAW_BUG(msg) 1205 end if 1206 end if 1207 1208 isp=0;jband=iband1-1 1209 do ib=1,nband 1210 jband=jband+1 1211 if (has_distrb) then 1212 if (abs(proc_distrb(ikpt,jband,isppol)-me)/=0) then 1213 isp=isp+nspinor 1214 cycle 1215 end if 1216 end if 1217 do ispinor=1,nspinor 1218 isp=isp+1 1219 if (iorder==0) then 1220 if (ncpgr_==0) then 1221 do iatom=1,dimcp 1222 read(uncp) cprj_k(iatom,isp)%cp(:,:) 1223 end do 1224 else 1225 if (has_icpgr) then 1226 do iatom=1,dimcp 1227 LIBPAW_ALLOCATE(tmp,(2,ncpgr_,cprj_k(iatom,1)%nlmn)) 1228 read(uncp) cprj_k(iatom,isp)%cp(:,:),tmp(:,:,:) 1229 cprj_k(iatom,isp)%dcp(:,1,:)=tmp(:,icpgr_,:) 1230 LIBPAW_DEALLOCATE(tmp) 1231 end do 1232 else 1233 do iatom=1,dimcp 1234 read(uncp) cprj_k(iatom,isp)%cp(:,:),cprj_k(iatom,isp)%dcp(:,:,:) 1235 end do 1236 end if 1237 end if 1238 else 1239 if (ncpgr_==0) then 1240 do iatom=1,dimcp 1241 iatm=min(atind(iatom),dimcp) 1242 read(uncp) cprj_k(iatm,isp)%cp(:,:) 1243 end do 1244 else 1245 if (has_icpgr) then 1246 do iatom=1,dimcp 1247 iatm=min(atind(iatom),dimcp) 1248 LIBPAW_ALLOCATE(tmp,(2,ncpgr_,cprj_k(iatm,1)%nlmn)) 1249 read(uncp) cprj_k(iatm,isp)%cp(:,:),tmp(:,:,:) 1250 cprj_k(iatm,isp)%dcp(:,1,:)=tmp(:,icpgr_,:) 1251 LIBPAW_DEALLOCATE(tmp) 1252 end do 1253 else 1254 do iatom=1,dimcp 1255 iatm=min(atind(iatom),dimcp) 1256 read(uncp) cprj_k(iatm,isp)%cp(:,:),cprj_k(iatm,isp)%dcp(:,:,:) 1257 end do 1258 end if 1259 end if 1260 end if 1261 end do 1262 end do 1263 1264 else 1265 1266 isp=0;ibsp=ibg+nspinor*(iband1-1);jband=iband1-1 1267 do ib=1,nband 1268 jband=jband+1 1269 if (has_distrb) then 1270 if (abs(proc_distrb(ikpt,jband,isppol)-me)/=0) then 1271 isp=isp+nspinor;ibsp=ibsp+nspinor 1272 cycle 1273 end if 1274 end if 1275 do ispinor=1,nspinor 1276 isp=isp+1;ibsp=ibsp+1 1277 if (iorder==0) then 1278 if (ncpgr_==0) then 1279 do iatom=1,dimcp 1280 cprj_k(iatom,isp)%cp(:,:)=cprj(iatom,ibsp)%cp(:,:) 1281 end do 1282 else 1283 if (has_icpgr) then 1284 do iatom=1,dimcp 1285 cprj_k(iatom,isp)%cp(:,:) =cprj(iatom,ibsp)%cp(:,:) 1286 cprj_k(iatom,isp)%dcp(:,1,:)=cprj(iatom,ibsp)%dcp(:,icpgr_,:) 1287 end do 1288 else 1289 do iatom=1,dimcp 1290 cprj_k(iatom,isp)%cp(:,:) =cprj(iatom,ibsp)%cp(:,:) 1291 cprj_k(iatom,isp)%dcp(:,:,:)=cprj(iatom,ibsp)%dcp(:,:,:) 1292 end do 1293 end if 1294 end if 1295 else 1296 if (ncpgr_==0) then 1297 do iatom=1,dimcp 1298 iatm=min(atind(iatom),dimcp) 1299 cprj_k(iatm,isp)%cp(:,:)=cprj(iatom,ibsp)%cp(:,:) 1300 end do 1301 else 1302 if (has_icpgr) then 1303 do iatom=1,dimcp 1304 iatm=min(atind(iatom),dimcp) 1305 cprj_k(iatm,isp)%cp(:,:) =cprj(iatom,ibsp)%cp(:,:) 1306 cprj_k(iatm,isp)%dcp(:,1,:)=cprj(iatom,ibsp)%dcp(:,icpgr_,:) 1307 end do 1308 else 1309 do iatom=1,dimcp 1310 iatm=min(atind(iatom),dimcp) 1311 cprj_k(iatm,isp)%cp(:,:) =cprj(iatom,ibsp)%cp(:,:) 1312 cprj_k(iatm,isp)%dcp(:,:,:)=cprj(iatom,ibsp)%dcp(:,:,:) 1313 end do 1314 end if 1315 end if 1316 end if 1317 end do 1318 end do 1319 1320 end if 1321 1322 end subroutine pawcprj_get
m_pawcprj/pawcprj_getdim [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_getdim
FUNCTION
Helper function returning the number of lmn components in the <p_{lmn}^i|\psi> for the i-th atom. Used to initialize the dimensioning array that is passed to the pawcprj_alloc routines when the pawcprj_type structure is allocated and initialized.
INPUTS
natom=number of atoms in the unit cell nattyp(ntypat)=number of atoms of each type ntypat=number of atom types typat(natom-= type of each atom Pawtab(ntypat)<pawtab_type>=PAW tabulated starting data. sort_mode(len=*)=String defining the sorting of the atoms in the Cprj arrays. Two modes are possible: -- "O[rdered]", if atoms are sorted by atom type. -- "R[andom]", if atoms are sorted randomly i.e. according the values of typat specified in the input file.
OUTPUT
dimcprj(natom)=Number of nlm elements in the <p_{lmn}^i|\psi> matrix elements for i=1,...,natom.
SOURCE
2684 subroutine pawcprj_getdim(dimcprj,natom,nattyp,ntypat,typat,Pawtab,sort_mode) 2685 2686 !Arguments ------------------------------------ 2687 integer,intent(in) :: natom,ntypat 2688 character(len=*),intent(in) :: sort_mode 2689 !arrays 2690 integer,intent(in) :: nattyp(:),typat(natom) 2691 integer,intent(inout) :: dimcprj(natom) 2692 type(Pawtab_type),intent(in) :: Pawtab(ntypat) 2693 2694 !Local variables------------------------------- 2695 integer :: iatom,itypat 2696 character(len=500) :: msg 2697 2698 ! ************************************************************************* 2699 2700 SELECT CASE (sort_mode(1:1)) 2701 2702 CASE ("o","O") ! Ordered by atom-type 2703 2704 iatom=0 2705 do itypat=1,ntypat 2706 dimcprj(iatom+1:iatom+nattyp(itypat))=Pawtab(itypat)%lmn_size 2707 iatom=iatom+nattyp(itypat) 2708 end do 2709 2710 CASE ("r","R") ! Randomly ordered (typat from input file) 2711 2712 do iatom=1,natom 2713 itypat=typat(iatom) 2714 dimcprj(iatom)=Pawtab(itypat)%lmn_size 2715 end do 2716 2717 CASE DEFAULT 2718 msg='Wrong value for sort_mode: '//TRIM(sort_mode) 2719 LIBPAW_ERROR(msg) 2720 END SELECT 2721 2722 end subroutine pawcprj_getdim
m_pawcprj/pawcprj_lincom [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_lincom
FUNCTION
Compute a LINear COMbination of cprj datastructure: cprj_out(:,:) <--- Sum_i [ alpha_i . cprj_i(:,:) ] alpha_i are COMPLEX scalars
INPUTS
alpha(2,nn)= alpha COMPLEX factors cprj_in(:,:) <type(pawcprj_type)>= input cprj_in datastructure nn= number of cprj involved in the linear combination
OUTPUT
cprj_out(:,:) <type(pawcprj_type)>= output cprj_out datastructure
NOTES
cprj_in and cprj_out must be dimensionned as cprj_in(n1,n2*nn) and cprj_in(n1,n2)
SOURCE
967 subroutine pawcprj_lincom(alpha,cprj_in,cprj_out,nn) 968 969 !Arguments ------------------------------------ 970 !scalars 971 integer,intent(in) :: nn 972 real(dp),intent(in) :: alpha(2,nn) 973 !arrays 974 type(pawcprj_type),intent(in) :: cprj_in(:,:) 975 type(pawcprj_type),intent(inout) :: cprj_out(:,:) 976 977 !Local variables------------------------------- 978 !scalars 979 integer :: ii,in,jj,jn,kk,ll,n1in,n1out,n2in,n2out,ncpgrin,ncpgrout,nlmn 980 real(dp) :: cp1,cp2 981 character(len=500) :: msg 982 983 ! ************************************************************************* 984 985 n1in=size(cprj_in,dim=1);n1out=size(cprj_out,dim=1) 986 n2in=size(cprj_in,dim=2);n2out=size(cprj_out,dim=2) 987 ncpgrin=cprj_in(1,1)%ncpgr;ncpgrout=cprj_out(1,1)%ncpgr 988 989 msg = "" 990 if (n1in/=n1out) msg = TRIM(msg)//"Bug in pawcprj_lincom: n1 wrong sizes!"//ch10 991 if (n2in/=n2out*nn) msg = TRIM(msg)//"Bug in pawcprj_lincom: n2 wrong sizes!"//ch10 992 if (ncpgrin/=ncpgrout) msg = TRIM(msg)//"Bug in pawcprj_lincom: ncpgr wrong sizes!"//ch10 993 if (LEN_TRIM(msg) > 0) then 994 LIBPAW_ERROR(msg) 995 end if 996 997 do jj=1,n2out 998 do ii=1,n1out 999 nlmn=cprj_in(ii,jj)%nlmn 1000 cprj_out(ii,jj)%nlmn=nlmn 1001 cprj_out(ii,jj)%cp(1:2,1:nlmn)=zero 1002 jn=jj 1003 do in=1,nn 1004 do kk=1,nlmn 1005 cp1=cprj_out(ii,jj)%cp(1,kk) & 1006 & +alpha(1,in)*cprj_in(ii,jn)%cp(1,kk)-alpha(2,in)*cprj_in(ii,jn)%cp(2,kk) 1007 cp2=cprj_out(ii,jj)%cp(2,kk) & 1008 & +alpha(1,in)*cprj_in(ii,jn)%cp(2,kk)+alpha(2,in)*cprj_in(ii,jn)%cp(1,kk) 1009 cprj_out(ii,jj)%cp(1,kk)=cp1 1010 cprj_out(ii,jj)%cp(2,kk)=cp2 1011 end do 1012 jn=jn+n2out 1013 end do 1014 end do 1015 end do 1016 1017 if (ncpgrin>0) then 1018 do jj=1,n2out 1019 do ii=1,n1out 1020 nlmn=cprj_in(ii,jj)%nlmn 1021 cprj_out(ii,jj)%dcp(1:2,1:ncpgrin,1:nlmn)=zero 1022 jn=jj 1023 do in=1,nn 1024 do kk=1,nlmn 1025 do ll=1,ncpgrin 1026 cp1=cprj_out(ii,jj)%dcp(1,ll,kk) & 1027 & +alpha(1,in)*cprj_in(ii,jn)%dcp(1,ll,kk) & 1028 & -alpha(2,in)*cprj_in(ii,jn)%dcp(2,ll,kk) 1029 cp2=cprj_out(ii,jj)%dcp(2,ll,kk) & 1030 & +alpha(1,in)*cprj_in(ii,jn)%dcp(2,ll,kk) & 1031 +alpha(2,in)*cprj_in(ii,jn)%dcp(1,ll,kk) 1032 cprj_out(ii,jj)%dcp(1,ll,kk)=cp1 1033 cprj_out(ii,jj)%dcp(2,ll,kk)=cp2 1034 end do 1035 end do 1036 jn=jn+n2out 1037 end do 1038 end do 1039 end do 1040 end if 1041 1042 end subroutine pawcprj_lincom
m_pawcprj/pawcprj_mpi_allgather [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_mpi_allgather
FUNCTION
Perform MPI_ALLGATHER on a pawcprj_type inside a MPI communicator.
INPUTS
cprj_loc= The cprj on the local proc being all-gathered natom=Number of atoms (size of first dimension of cprj_loc). n2dim=Size of the second dimension of cprj_loc. n2std=Stride of n2 dimension if n2std=1, cprj_loc(:,1) is on proc 0, cprj_loc(:,2) is on proc 1, cprj_loc(:,3) is on proc 2, etc. if n2std>1, cprj_loc(:,1:n2std) are on proc 0, cprj_loc(:,n2std+1,2*n2std) are on proc 1, etc. nlmn(natom)=Number of nlm partial waves for each atom. ncpgr = number of gradients in cprj_loc nproc=number of processors being gathered spaceComm=MPI Communicator. [rank_ordered]= optional, default=FALSE TRUE: second dimension of gathered datastructure is rank-ordered FALSE: second dimension of gathered datastructure is not rank-ordered
OUTPUT
cprj_gat=the gathered cprjs ierr=Error status.
SOURCE
2085 subroutine pawcprj_mpi_allgather(cprj_loc,cprj_gat,natom,n2dim,n2std,nlmn,ncpgr,nproc,spaceComm,ierr,& 2086 & rank_ordered) 2087 2088 !Arguments ------------------------------------ 2089 !scalars 2090 integer,intent(in) :: natom,n2dim,n2std,ncpgr,nproc,spaceComm 2091 integer,intent(out) :: ierr 2092 logical,optional,intent(in) :: rank_ordered 2093 !arrays 2094 integer,intent(in) :: nlmn(natom) 2095 type(pawcprj_type),intent(in) :: cprj_loc(:,:) 2096 type(pawcprj_type),intent(inout) :: cprj_gat(:,:) 2097 2098 !Local variables------------------------------- 2099 !scalars 2100 integer :: iat,ii,jj,t2dim,tcpgr,tg2dim,n1dim,nn 2101 integer :: ntotcp,ibuf,ipck,iproc 2102 logical :: rank_ordered_ 2103 character(len=100) :: msg 2104 !arrays 2105 real(dp),allocatable :: buffer_cpgr(:,:,:),buffer_cpgr_all(:,:,:) 2106 2107 ! ************************************************************************* 2108 2109 n1dim=0 2110 t2dim=0 2111 tg2dim=0 2112 tcpgr=0 2113 ierr=0 2114 2115 nn=size(nlmn,dim=1) 2116 n1dim=size(cprj_loc,dim=1) 2117 t2dim=size(cprj_loc,dim=2) 2118 tg2dim=size(cprj_gat,dim=2) 2119 tcpgr=cprj_loc(1,1)%ncpgr 2120 2121 if (nn/=n1dim) then 2122 msg='size mismatch in natom (pawcprj_mpi_allgather)!' 2123 LIBPAW_BUG(msg) 2124 end if 2125 if (t2dim/=n2dim) then 2126 msg='size mismatch in dim=2 (pawcprj_mpi_allgather)!' 2127 LIBPAW_BUG(msg) 2128 end if 2129 if (tg2dim/=n2dim*nproc) then 2130 msg='size mismatch in dim=2 (pawcprj_mpi_allgather)!' 2131 LIBPAW_BUG(msg) 2132 end if 2133 if (tcpgr/=ncpgr) then 2134 msg='size mismatch in ncpgr (pawcprj_mpi_allgather)!' 2135 LIBPAW_BUG(msg) 2136 end if 2137 if (mod(n2dim,n2std)/=0) then 2138 msg='n2std should divide n2dim (pawcprj_mpi_allgather)!' 2139 LIBPAW_BUG(msg) 2140 end if 2141 2142 rank_ordered_=.false.;if(present(rank_ordered)) rank_ordered_=rank_ordered 2143 2144 ntotcp=n2dim*SUM(nlmn(:)) 2145 LIBPAW_ALLOCATE(buffer_cpgr,(2,1+ncpgr,ntotcp)) 2146 LIBPAW_ALLOCATE(buffer_cpgr_all,(2,1+ncpgr,nproc*ntotcp)) 2147 2148 !=== Pack cprj_loc ==== 2149 ipck=0 2150 do jj=1,n2dim 2151 do iat=1,natom 2152 nn=nlmn(iat) 2153 buffer_cpgr(:,1,ipck+1:ipck+nn)=cprj_loc(iat,jj)%cp(:,1:nn) 2154 if (ncpgr/=0) buffer_cpgr(:,2:1+ncpgr,ipck+1:ipck+nn)=cprj_loc(iat,jj)%dcp(:,:,1:nn) 2155 ipck=ipck+nn 2156 end do 2157 end do 2158 2159 !=== allgather data === 2160 call xmpi_allgather(buffer_cpgr,2*(ncpgr+1)*ntotcp,buffer_cpgr_all,spaceComm,ierr) 2161 2162 !=== unpack gathered data into cprj(natom,n2dim*nproc) 2163 !=== second dimension is rank-ordered if rank_ordered_=true 2164 ipck=0 2165 do iproc=1,nproc 2166 do jj=1,n2dim/n2std 2167 do ii=1,n2std 2168 if (rank_ordered_) then 2169 ibuf=(iproc-1)*n2dim+(jj-1)*n2std+ii 2170 else 2171 ibuf=(iproc+(jj-1)*nproc-1)*n2std+ii 2172 end if 2173 do iat=1,natom 2174 nn=nlmn(iat) 2175 cprj_gat(iat,ibuf)%cp(:,1:nn)=buffer_cpgr_all(:,1,ipck+1:ipck+nn) 2176 if (ncpgr/=0) cprj_gat(iat,ibuf)%dcp(:,1:ncpgr,1:nn)=& 2177 & buffer_cpgr_all(:,2:1+ncpgr,ipck+1:ipck+nn) 2178 ipck=ipck+nn 2179 end do 2180 end do 2181 end do 2182 end do 2183 2184 !=== Clean up === 2185 LIBPAW_DEALLOCATE(buffer_cpgr) 2186 LIBPAW_DEALLOCATE(buffer_cpgr_all) 2187 2188 end subroutine pawcprj_mpi_allgather
m_pawcprj/pawcprj_mpi_exch [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_mpi_exch
FUNCTION
Exchange a pawcprj_type between two processors inside a MPI communicator.
INPUTS
natom=Number of atoms (size of first dimension of Cprj_send and Cprj_recv). n2dim=Size of the second dimension. nlmn(natom)=Number of nlm partial waves for each atom. Cprj_send= The datatype to be transmitted. receiver=ID of the receiver in spaceComm. sender=ID of the sender in spaceComm. spaceComm=MPI Communicator. mtag= message tag
OUTPUT
ierr=Error status. Cprj_recv=The datatype copied on proc. receiver.
NOTES
If sender==receiver, Cprj_send is copied into Cprj_recv. It should be easy to avoid this additional copy in the calling routine.
SOURCE
1661 subroutine pawcprj_mpi_exch(natom,n2dim,nlmn,ncpgr,Cprj_send,Cprj_recv,sender,receiver,spaceComm,mtag,ierr) 1662 1663 !Arguments ------------------------------------ 1664 !scalars 1665 integer,intent(in) :: mtag,natom,n2dim,ncpgr 1666 integer,intent(in) :: sender,receiver,spaceComm 1667 integer,intent(out) :: ierr 1668 !arrays 1669 integer,intent(in) :: nlmn(natom) 1670 type(pawcprj_type),intent(in) :: Cprj_send(:,:) 1671 type(pawcprj_type),intent(inout) :: Cprj_recv(:,:) 1672 1673 !Local variables------------------------------- 1674 !scalars 1675 integer :: iat,jj,t2dim,tcpgr,n1dim,nn 1676 integer :: ntotcp,ipck,rank 1677 character(len=500) :: msg 1678 !arrays 1679 real(dp),allocatable :: buffer_cp(:,:),buffer_cpgr(:,:,:) 1680 1681 ! ************************************************************************* 1682 1683 n1dim=0 1684 t2dim=0 1685 tcpgr=0 1686 ierr=0 1687 if (sender==receiver) then 1688 call pawcprj_copy(Cprj_send,Cprj_recv) 1689 return 1690 end if 1691 1692 rank = xmpi_comm_rank(spaceComm) 1693 1694 nn=size(nlmn,dim=1) 1695 if (rank==sender) then 1696 n1dim=size(Cprj_send,dim=1) 1697 t2dim=size(Cprj_send,dim=2) 1698 tcpgr=Cprj_send(1,1)%ncpgr 1699 end if 1700 if (rank==receiver) then 1701 n1dim=size(Cprj_recv,dim=1) 1702 t2dim=size(Cprj_recv,dim=2) 1703 tcpgr=Cprj_recv(1,1)%ncpgr 1704 end if 1705 if (rank/=sender.and.rank/=receiver) then 1706 write(msg,'(a,3i0)') & 1707 & 'rank is not equal to sender or receiver (pawcprj_mpi_exch): ',rank, sender, receiver 1708 LIBPAW_BUG(msg) 1709 end if 1710 1711 ntotcp=n2dim*SUM(nlmn(:)) 1712 1713 LIBPAW_ALLOCATE(buffer_cp,(2,ntotcp)) 1714 if (ncpgr/=0) then 1715 LIBPAW_ALLOCATE(buffer_cpgr,(2,ncpgr,ntotcp)) 1716 end if 1717 1718 !=== Pack Cprj_send === 1719 if (rank==sender) then 1720 ipck=0 1721 do jj=1,n2dim 1722 do iat=1,natom 1723 nn=nlmn(iat) 1724 buffer_cp(:,ipck+1:ipck+nn)=Cprj_send(iat,jj)%cp(:,1:nn) 1725 if (ncpgr/=0) buffer_cpgr(:,:,ipck+1:ipck+nn)=Cprj_send(iat,jj)%dcp(:,:,1:nn) 1726 ipck=ipck+nn 1727 end do 1728 end do 1729 end if 1730 1731 !=== Transmit data === 1732 call xmpi_exch(buffer_cp,2*ntotcp,sender,buffer_cp,receiver,spaceComm,2*mtag,ierr) 1733 if (ncpgr/=0) then 1734 call xmpi_exch(buffer_cpgr,2*ncpgr*ntotcp,sender,buffer_cpgr,receiver,spaceComm,2*mtag+1,ierr) 1735 end if 1736 1737 !=== UnPack buffers into Cprj_recv === 1738 if (rank==receiver) then 1739 ipck=0 1740 do jj=1,n2dim 1741 do iat=1,natom 1742 nn=nlmn(iat) 1743 Cprj_recv(iat,jj)%cp(:,1:nn)=buffer_cp(:,ipck+1:ipck+nn) 1744 if (ncpgr/=0) Cprj_recv(iat,jj)%dcp(:,:,1:nn)=buffer_cpgr(:,:,ipck+1:ipck+nn) 1745 ipck=ipck+nn 1746 end do 1747 end do 1748 end if 1749 1750 LIBPAW_DEALLOCATE(buffer_cp) 1751 if (ncpgr/=0) then 1752 LIBPAW_DEALLOCATE(buffer_cpgr) 1753 end if 1754 1755 end subroutine pawcprj_mpi_exch
m_pawcprj/pawcprj_mpi_recv [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_mpi_recv
FUNCTION
Receive a pawcprj_type inside a MPI communicator.
INPUTS
natom=Number of atoms (size of first dimension of Cprj_in). n2dim=Size of the second dimension. nlmn(natom)=Number of nlm partial waves for each atom. ncpgr = number of gradients in cprj_in sender=ID of the sender in spaceComm. spaceComm=MPI Communicator.
OUTPUT
ierr=Error status. cprj_in=The datatype copied on proc. receiver.
NOTES
Perhaps in general it is more efficient to use pawcprj_mpi_exch but it is convenient for coding to have separate send and receive routines.
SOURCE
1891 subroutine pawcprj_mpi_recv(natom,n2dim,nlmn,ncpgr,cprj_in,sender,spaceComm,ierr) 1892 1893 !Arguments ------------------------------------ 1894 !scalars 1895 integer,intent(in) :: natom,n2dim,ncpgr 1896 integer,intent(in) :: sender,spaceComm 1897 integer,intent(out) :: ierr 1898 !arrays 1899 integer,intent(in) :: nlmn(natom) 1900 type(pawcprj_type),intent(inout) :: cprj_in(:,:) 1901 1902 !Local variables------------------------------- 1903 !scalars 1904 integer :: iat,jj,t2dim,tcpgr,n1dim,nn 1905 integer :: ntotcp,ipck,tag 1906 character(len=100) :: msg 1907 !arrays 1908 real(dp),allocatable :: buffer_cp(:,:),buffer_cpgr(:,:,:) 1909 1910 ! ************************************************************************* 1911 1912 n1dim=0 1913 t2dim=0 1914 tcpgr=0 1915 ierr=0 1916 1917 nn=size(nlmn,dim=1) 1918 n1dim=size(cprj_in,dim=1) 1919 t2dim=size(cprj_in,dim=2) 1920 tcpgr=cprj_in(1,1)%ncpgr 1921 1922 if (nn/=n1dim) then 1923 msg='size mismatch in natom (pawcprj_mpi_recv)!' 1924 LIBPAW_BUG(msg) 1925 end if 1926 if (t2dim/=n2dim) then 1927 msg='size mismatch in dim=2 (pawcprj_mpi_recv)!' 1928 LIBPAW_BUG(msg) 1929 end if 1930 if (tcpgr/=ncpgr) then 1931 msg='size mismatch in ncpgr (pawcprj_mpi_recv)!' 1932 LIBPAW_BUG(msg) 1933 end if 1934 1935 ntotcp=n2dim*SUM(nlmn(:)) 1936 1937 LIBPAW_ALLOCATE(buffer_cp,(2,ntotcp)) 1938 if (ncpgr/=0) then 1939 LIBPAW_ALLOCATE(buffer_cpgr,(2,ncpgr,ntotcp)) 1940 end if 1941 1942 !=== Receive data === 1943 tag = 2*ntotcp 1944 call xmpi_recv(buffer_cp,sender,tag,spaceComm,ierr) 1945 if (ncpgr/=0) then 1946 tag=tag*ncpgr 1947 call xmpi_recv(buffer_cpgr,sender,tag,spaceComm,ierr) 1948 end if 1949 1950 !=== UnPack buffers into cprj_in === 1951 ipck=0 1952 do jj=1,n2dim 1953 do iat=1,natom 1954 nn=nlmn(iat) 1955 cprj_in(iat,jj)%cp(:,1:nn)=buffer_cp(:,ipck+1:ipck+nn) 1956 if (ncpgr/=0) cprj_in(iat,jj)%dcp(:,:,1:nn)=buffer_cpgr(:,:,ipck+1:ipck+nn) 1957 ipck=ipck+nn 1958 end do 1959 end do 1960 1961 !=== Clean up === 1962 LIBPAW_DEALLOCATE(buffer_cp) 1963 if (ncpgr/=0) then 1964 LIBPAW_DEALLOCATE(buffer_cpgr) 1965 end if 1966 1967 end subroutine pawcprj_mpi_recv
m_pawcprj/pawcprj_mpi_send [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_mpi_send
FUNCTION
Send a pawcprj_type inside a MPI communicator.
INPUTS
natom=Number of atoms (size of first dimension of cprj_out). n2dim=Size of the second dimension. nlmn(natom)=Number of nlm partial waves for each atom. ncpgr = number of gradients in cprj_out cprj_out= The datatype to be transmitted. receiver=ID of the receiver in spaceComm. spaceComm=MPI Communicator.
OUTPUT
ierr=Error status.
NOTES
perhaps in general it is more efficient to use pawcprj_mpi_exch but it is convenient for coding to have separate send and recieve routines.
SOURCE
1785 subroutine pawcprj_mpi_send(natom,n2dim,nlmn,ncpgr,cprj_out,receiver,spaceComm,ierr) 1786 1787 !Arguments ------------------------------------ 1788 !scalars 1789 integer,intent(in) :: natom,n2dim,ncpgr 1790 integer,intent(in) :: receiver,spaceComm 1791 integer,intent(out) :: ierr 1792 !arrays 1793 integer,intent(in) :: nlmn(natom) 1794 type(pawcprj_type),intent(in) :: cprj_out(:,:) 1795 1796 !Local variables------------------------------- 1797 !scalars 1798 integer :: iat,jj,t2dim,tcpgr,n1dim,nn 1799 integer :: ntotcp,ipck,tag 1800 character(len=100) :: msg 1801 !arrays 1802 real(dp),allocatable :: buffer_cp(:,:),buffer_cpgr(:,:,:) 1803 1804 ! ************************************************************************* 1805 1806 n1dim=0 1807 t2dim=0 1808 tcpgr=0 1809 ierr=0 1810 1811 nn=size(nlmn,dim=1) 1812 n1dim=size(cprj_out,dim=1) 1813 t2dim=size(cprj_out,dim=2) 1814 tcpgr=cprj_out(1,1)%ncpgr 1815 1816 if (nn/=n1dim) then 1817 msg='size mismatch in natom (pawcprj_mpi_send)!' 1818 LIBPAW_BUG(msg) 1819 end if 1820 if (t2dim/=n2dim) then 1821 msg='size mismatch in dim=2 (pawcprj_mpi_send)!' 1822 LIBPAW_BUG(msg) 1823 end if 1824 if (tcpgr/=ncpgr) then 1825 msg='size mismatch in ncpgr (pawcprj_mpi_send)!' 1826 LIBPAW_BUG(msg) 1827 end if 1828 1829 ntotcp=n2dim*SUM(nlmn(:)) 1830 1831 LIBPAW_ALLOCATE(buffer_cp,(2,ntotcp)) 1832 if (ncpgr/=0) then 1833 LIBPAW_ALLOCATE(buffer_cpgr,(2,ncpgr,ntotcp)) 1834 end if 1835 1836 !=== Pack cprj_out ==== 1837 ipck=0 1838 do jj=1,n2dim 1839 do iat=1,natom 1840 nn=nlmn(iat) 1841 buffer_cp(:,ipck+1:ipck+nn)=cprj_out(iat,jj)%cp(:,1:nn) 1842 if (ncpgr/=0) buffer_cpgr(:,:,ipck+1:ipck+nn)=cprj_out(iat,jj)%dcp(:,:,1:nn) 1843 ipck=ipck+nn 1844 end do 1845 end do 1846 1847 !=== Transmit data === 1848 tag = 2*ntotcp 1849 call xmpi_send(buffer_cp,receiver,tag,spaceComm,ierr) 1850 if (ncpgr/=0) then 1851 tag=tag*ncpgr 1852 call xmpi_send(buffer_cpgr,receiver,tag,spaceComm,ierr) 1853 end if 1854 1855 !=== Clean up === 1856 LIBPAW_DEALLOCATE(buffer_cp) 1857 if (ncpgr/=0) then 1858 LIBPAW_DEALLOCATE(buffer_cpgr) 1859 end if 1860 1861 end subroutine pawcprj_mpi_send
m_pawcprj/pawcprj_mpi_sum [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_mpi_sum
FUNCTION
Perform MPI_SUM on a pawcprj_type inside a MPI communicator.
INPUTS
spaceComm=MPI Communicator.
SIDE EFFECTS
cprj=the cprj datastructure ierr=Error status.
SOURCE
1988 subroutine pawcprj_mpi_sum(cprj,spaceComm,ierr) 1989 1990 !Arguments ------------------------------------ 1991 !scalars 1992 integer,intent(in) :: spaceComm 1993 integer,intent(out) :: ierr 1994 !arrays 1995 type(pawcprj_type),intent(inout) :: cprj(:,:) 1996 1997 !Local variables------------------------------- 1998 !scalars 1999 integer,parameter :: maxBytes=100*1024*1024 ! 100 MBytes 2000 integer :: ii,ipck,jj,ncpgr,nlmn,nn,n1dim,n2dim,n2dim1,n2dim2,sizeBytes,step 2001 logical,parameter :: save_memory=.true. 2002 !arrays 2003 real(dp),allocatable :: buffer_cprj(:,:,:) 2004 2005 ! ************************************************************************* 2006 2007 if (xmpi_comm_size(spaceComm)<2) return 2008 2009 n1dim=size(cprj,1);n2dim=size(cprj,2) 2010 nlmn=sum(cprj(:,:)%nlmn) 2011 ncpgr=maxval(cprj(:,:)%ncpgr) 2012 2013 step=n2dim 2014 if (save_memory) then 2015 sizeBytes=2*(1+ncpgr)*nlmn *8 2016 step=n2dim/max(1,sizeBytes/maxBytes) 2017 if (step==0) step=1 2018 end if 2019 2020 do n2dim1=1,n2dim,step 2021 2022 n2dim2=min(n2dim1+step-1,n2dim) 2023 nlmn=sum(cprj(:,n2dim1:n2dim2)%nlmn) 2024 LIBPAW_ALLOCATE(buffer_cprj,(2,1+ncpgr,nlmn)) 2025 2026 ipck=0 ; buffer_cprj=zero 2027 do jj=n2dim1,n2dim2 2028 do ii=1,n1dim 2029 nn=cprj(ii,jj)%nlmn 2030 buffer_cprj(:,1,ipck+1:ipck+nn)=cprj(ii,jj)%cp(:,1:nn) 2031 if (cprj(ii,jj)%ncpgr/=0) buffer_cprj(:,2:1+ncpgr,ipck+1:ipck+nn)=cprj(ii,jj)%dcp(:,1:ncpgr,1:nn) 2032 ipck=ipck+nn 2033 end do 2034 end do 2035 2036 call xmpi_sum(buffer_cprj,spaceComm,ierr) 2037 2038 ipck=0 2039 do jj=n2dim1,n2dim2 2040 do ii=1,n1dim 2041 nn=cprj(ii,jj)%nlmn 2042 cprj(ii,jj)%cp(:,1:nn)=buffer_cprj(:,1,ipck+1:ipck+nn) 2043 if (cprj(ii,jj)%ncpgr/=0) cprj(ii,jj)%dcp(:,1:ncpgr,1:nn)=buffer_cprj(:,2:1+ncpgr,ipck+1:ipck+nn) 2044 ipck=ipck+nn 2045 end do 2046 end do 2047 2048 LIBPAW_DEALLOCATE(buffer_cprj) 2049 2050 end do 2051 2052 end subroutine pawcprj_mpi_sum
m_pawcprj/pawcprj_output [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_output
FUNCTION
Output a cprj. Useful for debugging.
INPUTS
cprj(:,:) <type(pawcprj_type)>= cprj datastructure prtgrads :: optional, 1 to print gradients also
OUTPUT
SOURCE
1062 subroutine pawcprj_output(cprj,prtgrads) 1063 1064 !Arguments ------------------------------------ 1065 !scalar 1066 integer,optional :: prtgrads 1067 !arrays 1068 type(pawcprj_type),intent(in) :: cprj(:,:) 1069 1070 !Local variables------------------------------- 1071 !scalar 1072 integer :: ii,jj,kk,nlmn,n1dim,n2dim 1073 logical :: gradoutput 1074 1075 ! ************************************************************************* 1076 1077 n1dim=size(cprj,dim=1) 1078 n2dim=size(cprj,dim=2) 1079 gradoutput = .FALSE. 1080 if(present(prtgrads)) then 1081 gradoutput = (prtgrads .EQ. 1) 1082 end if 1083 1084 write(std_out,'(a)')' pawcprj_output ' 1085 1086 do jj=1,n2dim 1087 do ii=1,n1dim 1088 write(std_out,'(a,i4,a,i4)')'atom ',ii,' band*k ',jj 1089 nlmn=cprj(ii,jj)%nlmn 1090 do kk=1,nlmn 1091 write(std_out,'(2f12.8)')cprj(ii,jj)%cp(1,kk),cprj(ii,jj)%cp(2,kk) 1092 if(gradoutput) then 1093 write(std_out,'(6f12.8)')cprj(ii,jj)%dcp(1,1,kk),cprj(ii,jj)%dcp(2,1,kk),& 1094 &cprj(ii,jj)%dcp(1,2,kk),cprj(ii,jj)%dcp(2,2,kk),& 1095 &cprj(ii,jj)%dcp(1,3,kk),cprj(ii,jj)%dcp(2,3,kk) 1096 end if 1097 end do 1098 end do 1099 end do 1100 1101 end subroutine pawcprj_output
m_pawcprj/pawcprj_pack [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_pack
FUNCTION
Pack structured data into a simple buffer
INPUTS
nlmn(natom)=Number of nlm partial waves for each atom. ncpgr = number of gradients in cprj_out cprj= The datatype to be packed.
OUTPUT
buffer = the data packed, dim : (2, n2dim*sum(nlmn)) [buffer_gr] = if present the gradient data packed, dim : (2, ncpgr, n2dim*sum(nlmn))
SOURCE
2831 subroutine pawcprj_pack(nlmn,cprj,buffer,buffer_gr) 2832 2833 !Arguments ------------------------------------ 2834 !scalars 2835 !arrays 2836 integer,intent(in) :: nlmn(:) 2837 type(pawcprj_type),intent(in) :: cprj(:,:) 2838 real(dp),intent(out) :: buffer(:,:) 2839 real(dp),intent(out),optional :: buffer_gr(:,:,:) 2840 2841 !Local variables------------------------------- 2842 !scalars 2843 integer :: natom,n2buffer,ncpgr,n2dim 2844 integer :: iat,jj,n1dim,nn 2845 integer :: ipck 2846 character(len=100) :: msg 2847 !arrays 2848 2849 ! ************************************************************************* 2850 2851 natom=size(nlmn,dim=1) 2852 n2buffer=size(buffer,dim=2) 2853 n1dim=size(cprj,dim=1) 2854 n2dim=size(cprj,dim=2) 2855 2856 if (natom/=n1dim) then 2857 msg='size mismatch in natom (pawcprj_pack)!' 2858 LIBPAW_BUG(msg) 2859 end if 2860 if (n2dim*SUM(nlmn)/=n2buffer) then 2861 msg='size mismatch in dim=2 (pawcprj_pack)!' 2862 LIBPAW_BUG(msg) 2863 end if 2864 ncpgr=0 2865 if (present(buffer_gr)) then 2866 ncpgr=size(buffer_gr,dim=2) 2867 end if 2868 2869 !=== Pack cprj ==== 2870 ipck=0 2871 do jj=1,n2dim 2872 do iat=1,natom 2873 nn=nlmn(iat) 2874 buffer(:,ipck+1:ipck+nn)=cprj(iat,jj)%cp(:,1:nn) 2875 if (ncpgr/=0) then 2876 buffer_gr(:,:,ipck+1:ipck+nn)=cprj(iat,jj)%dcp(:,:,1:nn) 2877 end if 2878 ipck=ipck+nn 2879 end do 2880 end do 2881 2882 end subroutine pawcprj_pack
m_pawcprj/pawcprj_projbd [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_projbd
FUNCTION
Apply ZAXPBY (blas-like) operation with 2 cprj datastructures: cprjy(:,:) <- alpha.cprjx(:,:)+beta.cprjy(:,:) alpha and beta are COMPLEX scalars
INPUTS
alpha(2),beta(2)= alpha,beta COMPLEX factors cprjx(:,:) <type(pawcprj_type)>= input cprjx datastructure
SIDE EFFECTS
cprjy(:,:) <type(pawcprj_type)>= input/output cprjy datastructure
SOURCE
671 subroutine pawcprj_projbd(alpha,cprjx,cprjy) 672 673 !Arguments ------------------------------------ 674 !scalars 675 real(dp),intent(in) :: alpha(:,:) 676 !arrays 677 type(pawcprj_type),intent(in) :: cprjx(:,:) 678 type(pawcprj_type),intent(inout) :: cprjy(:,:) 679 680 !Local variables------------------------------- 681 !scalars 682 integer :: ia,ii,jj,kk,ll,n1dima,n1dimx,n1dimy,n2dimx,n2dimy,n2dima,ncpgrx,ncpgry,nlmn 683 real(dp) :: cp1,cp2,norma 684 character(len=500) :: msg 685 686 ! ************************************************************************* 687 688 n1dimy=size(cprjy,dim=1);n2dimy=size(cprjy,dim=2);ncpgry=cprjy(1,1)%ncpgr 689 n1dimx=size(cprjx,dim=1);n2dimx=size(cprjx,dim=2);ncpgrx=cprjx(1,1)%ncpgr 690 n1dima=size(alpha,dim=1);n2dima=size(alpha,dim=2) 691 msg = "" 692 if (n1dima/=2) msg = TRIM(msg)//"Error in pawcprj_projbd: alpha n1 wrong sizes !"//ch10 693 if (n1dimx/=n1dimy) msg = TRIM(msg)//"Error in pawcprj_projbd: n1 wrong sizes !"//ch10 694 if (n2dimx/=n2dimy*n2dima) msg = TRIM(msg)//"Error in pawcprj_projbd: n2 wrong sizes !"//ch10 695 if (ncpgrx/=ncpgry) msg = TRIM(msg)//"Error in pawcprj_projbd: ncpgr wrong sizes !"//ch10 696 if (LEN_TRIM(msg) > 0) then 697 LIBPAW_ERROR(msg) 698 end if 699 700 do ia=1,n2dima 701 norma=alpha(1,ia)**2+alpha(2,ia)**2 702 if (norma>tol16*tol16) then 703 do jj=1,n2dimy 704 do ii=1,n1dimx 705 nlmn=cprjy(ii,jj)%nlmn 706 cprjy(ii,jj)%nlmn =nlmn 707 do kk=1,nlmn 708 cp1=alpha(1,ia)*cprjx(ii,jj+(ia-1)*n2dimy)%cp(1,kk)-alpha(2,ia)*cprjx(ii,jj+(ia-1)*n2dimy)%cp(2,kk) 709 cp2=alpha(1,ia)*cprjx(ii,jj+(ia-1)*n2dimy)%cp(2,kk)+alpha(2,ia)*cprjx(ii,jj+(ia-1)*n2dimy)%cp(1,kk) 710 cprjy(ii,jj)%cp(1,kk)=cprjy(ii,jj)%cp(1,kk)+cp1 711 cprjy(ii,jj)%cp(2,kk)=cprjy(ii,jj)%cp(2,kk)+cp2 712 end do 713 end do 714 end do 715 if (ncpgrx>0) then 716 do jj=1,n2dimy 717 do ii=1,n1dimx 718 nlmn=cprjy(ii,jj)%nlmn 719 do kk=1,nlmn 720 do ll=1,ncpgrx 721 cp1=alpha(1,ia)*cprjx(ii,jj+(ia-1)*n2dimy)%dcp(1,ll,kk)-alpha(2,ia)*cprjx(ii,jj+(ia-1)*n2dimy)%dcp(2,ll,kk) 722 cp2=alpha(1,ia)*cprjx(ii,jj+(ia-1)*n2dimy)%dcp(2,ll,kk)+alpha(2,ia)*cprjx(ii,jj+(ia-1)*n2dimy)%dcp(1,ll,kk) 723 cprjy(ii,jj)%dcp(1,ll,kk)=cprjy(ii,jj)%dcp(1,ll,kk)+cp1 724 cprjy(ii,jj)%dcp(2,ll,kk)=cprjy(ii,jj)%dcp(2,ll,kk)+cp2 725 end do 726 end do 727 end do 728 end do 729 end if 730 end if 731 end do 732 733 end subroutine pawcprj_projbd
m_pawcprj/pawcprj_put [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_put
FUNCTION
Write cprj_k for a given set of (n,k) into memory in cprj, or into a temporary file
INPUTS
atind(natom)=index table for atoms (see iorder below) cprj_k(dimcp,nspinor*nband) <type(pawcprj_type)>= input cprj datastructure dimcp=first dimension of cprj_k,cprjnk arrays (1 or natom) iband1=index of first band in cprj ibg=shift in cprj array to locate current k-point ikpt=index of current k-point (only needed for the parallel distribution) iorder=0 if cprj ordering does not change during reading 1 if cprj ordering changes during writing, depending on content of atind array: - if atind=atindx (type-sorted->unsorted) - if atind=atindx1 (unsorted->type-sorted) isppol=index of current spin component mband=maximum number of bands mkmem=number of k points which can fit in memory; set to 0 if use disk [mpi_comm]=(optional argument) MPI communicator over (k-pts,bands,spins) Must be used in association with proc_distrb argument [mpi_comm_band]=(optional argument) MPI communicator over bands Must be used in association with proc_distrb argument natom=number of atoms in cell nband=number of bands to export (usually 1, nband_k or nblockbd) nband_k=total number of bands for this k-point nlmn(dimcp)=array of dimensions of cprj_k,cprjnk datastructures nspinor=number of spinorial components of the wavefunctions (on current proc) nsppol=1 for unpolarized, 2 for spin-polarized [proc_distrb(nkpt,nband,nsppol)]=(optional argument) processor distribution Describe how cprj datastructures are distributed over processors When present, mpicomm argument must be also present [to_be_gathered]=(optional argument) TRUE if cprj_k arrays have to be gathered between procs (band-fft parallelism only) uncp=unit number for cprj data (used if mkmem=0)
SIDE EFFECTS
cprj(dimcp,nspinor*mband*mkmem*nsppol)=output cprj (used if mkmem/=0)
SOURCE
1370 subroutine pawcprj_put(atind,cprj_k,cprj,dimcp,iband1,ibg,ikpt,iorder,isppol,mband,& 1371 & mkmem,natom,nband,nband_k,nlmn,nspinor,nsppol,uncp,& 1372 & mpicomm,mpi_comm_band,proc_distrb,to_be_gathered) ! Optional arguments 1373 1374 !Arguments ------------------------------------ 1375 !scalars 1376 integer,intent(in) :: iband1,ibg,ikpt,iorder,isppol,dimcp,mband,mkmem 1377 integer,intent(in) :: natom,nband,nband_k,nspinor,nsppol,uncp 1378 integer,intent(in),optional :: mpicomm,mpi_comm_band 1379 logical,optional,intent(in) :: to_be_gathered 1380 !arrays 1381 integer,intent(in) :: atind(natom),nlmn(dimcp) 1382 integer,intent(in),optional :: proc_distrb(:,:,:) 1383 type(pawcprj_type),intent(inout) :: cprj(dimcp,nspinor*mband*mkmem*nsppol) 1384 type(pawcprj_type),intent(in) :: cprj_k(dimcp,nspinor*nband) 1385 1386 !Local variables------------------------------- 1387 !scalars 1388 integer :: iatm,iatom,iband,ibsp,icpgr,ierr,ii,ilmn,isp,ispinor,jband,jj 1389 integer :: lmndim,me,ncpgr,nproc_band 1390 logical :: has_distrb,to_be_gathered_ 1391 character(len=500) :: msg 1392 !arrays 1393 real(dp),allocatable :: buffer1(:),buffer2(:) 1394 ! ************************************************************************* 1395 1396 ncpgr=cprj_k(1,1)%ncpgr 1397 to_be_gathered_=.false.;if (present(to_be_gathered)) to_be_gathered_=to_be_gathered 1398 1399 !MPI data 1400 nproc_band=1;if (present(mpi_comm_band)) nproc_band=xmpi_comm_size(mpi_comm_band) 1401 has_distrb=present(proc_distrb) 1402 if (has_distrb) then 1403 if (.not.present(mpicomm)) then 1404 msg='mpicomm must be present when proc_distrb is present (pawcprj_put)!' 1405 LIBPAW_BUG(msg) 1406 end if 1407 me=xmpi_comm_rank(mpicomm) 1408 end if 1409 1410 if (nproc_band==1.or.(.not.to_be_gathered_)) then 1411 1412 if (mkmem==0) then 1413 1414 if (iband1==1) write(uncp) nband_k 1415 1416 isp=0;jband=iband1-1 1417 do iband=1,nband 1418 jband=jband+1 1419 if (has_distrb) then 1420 if (abs(proc_distrb(ikpt,jband,isppol)-me)/=0) then 1421 isp=isp+nspinor 1422 cycle 1423 end if 1424 end if 1425 do ispinor=1,nspinor 1426 isp=isp+1 1427 if (iorder==0) then 1428 do iatom=1,dimcp 1429 if (ncpgr==0) then 1430 write(uncp) cprj_k(iatom,isp)%cp(:,:) 1431 else 1432 write(uncp) cprj_k(iatom,isp)%cp(:,:),cprj_k(iatom,isp)%dcp(:,:,:) 1433 end if 1434 end do 1435 else 1436 do iatom=1,dimcp 1437 iatm=min(atind(iatom),dimcp) 1438 if (ncpgr==0) then 1439 write(uncp) cprj_k(iatm,isp)%cp(:,:) 1440 else 1441 write(uncp) cprj_k(iatm,isp)%cp(:,:),cprj_k(iatm,isp)%dcp(:,:,:) 1442 end if 1443 end do 1444 end if 1445 end do 1446 end do 1447 1448 else 1449 1450 isp=0;ibsp=ibg+nspinor*(iband1-1);jband=iband1-1 1451 do iband=1,nband 1452 jband=jband+1 1453 if (has_distrb) then 1454 if (abs(proc_distrb(ikpt,jband,isppol)-me)/=0) then 1455 isp=isp+nspinor;ibsp=ibsp+nspinor 1456 cycle 1457 end if 1458 end if 1459 do ispinor=1,nspinor 1460 isp=isp+1;ibsp=ibsp+1 1461 if (iorder==0) then 1462 do iatom=1,dimcp 1463 cprj(iatom,ibsp)%cp(:,:)=cprj_k(iatom,isp)%cp(:,:) 1464 if (ncpgr>0) cprj(iatom,ibsp)%dcp(:,:,:)=cprj_k(iatom,isp)%dcp(:,:,:) 1465 end do 1466 else 1467 do iatom=1,dimcp 1468 iatm=min(atind(iatom),dimcp) 1469 cprj(iatom,ibsp)%cp(:,:)=cprj_k(iatm,isp)%cp(:,:) 1470 if (ncpgr>0) cprj(iatom,ibsp)%dcp(:,:,:)=cprj_k(iatm,isp)%dcp(:,:,:) 1471 end do 1472 end if 1473 end do 1474 end do 1475 1476 end if 1477 1478 else ! np_band>1 1479 1480 lmndim=2*sum(nlmn(1:dimcp))*(1+ncpgr)*nspinor 1481 LIBPAW_ALLOCATE(buffer1,(lmndim)) 1482 LIBPAW_ALLOCATE(buffer2,(lmndim*nproc_band)) 1483 isp=0;ibsp=ibg+nspinor*(iband1-1) 1484 do iband=1,nband ! must be nblockbd for band-fft parallelism 1485 jj=1 1486 do ispinor=1,nspinor 1487 isp=isp+1 1488 do iatom=1,dimcp 1489 if (iorder==0) then 1490 iatm=iatom 1491 else 1492 iatm=min(atind(iatom),dimcp) 1493 end if 1494 do ilmn=1,nlmn(iatm) 1495 buffer1(jj:jj+1)=cprj_k(iatm,isp)%cp(1:2,ilmn) 1496 jj=jj+2 1497 end do 1498 if (ncpgr>0) then 1499 do ilmn=1,nlmn(iatm) 1500 do icpgr=1,ncpgr 1501 buffer1(jj:jj+1)=cprj_k(iatm,isp)%dcp(1:2,icpgr,ilmn) 1502 jj=jj+2 1503 end do 1504 end do 1505 end if 1506 end do !iatom 1507 end do !ispinor 1508 call xmpi_allgather(buffer1,lmndim,buffer2,mpi_comm_band,ierr) 1509 jj=1 1510 do ii=1,nproc_band 1511 do ispinor=1,nspinor 1512 ibsp=ibsp+1 1513 do iatom=1,dimcp 1514 if (iorder==0) then 1515 iatm=iatom 1516 else 1517 iatm=min(atind(iatom),dimcp) 1518 end if 1519 do ilmn=1,nlmn(iatm) 1520 cprj(iatom,ibsp)%cp(1:2,ilmn)=buffer2(jj:jj+1) 1521 jj=jj+2 1522 end do 1523 if (ncpgr>0) then 1524 do ilmn=1,nlmn(iatm) 1525 do icpgr=1,ncpgr 1526 cprj(iatom,ibsp)%dcp(1:2,icpgr,ilmn)=buffer2(jj:jj+1) 1527 jj=jj+2 1528 end do 1529 end do 1530 end if 1531 end do !iatom 1532 end do !ispinor 1533 end do !ii=1,nproc_band 1534 end do !iband 1535 LIBPAW_DEALLOCATE(buffer1) 1536 LIBPAW_DEALLOCATE(buffer2) 1537 1538 end if ! mode_para=b, nband 1539 1540 end subroutine pawcprj_put
m_pawcprj/pawcprj_reorder [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_reorder
FUNCTION
Change the order of a cprj datastructure From unsorted cprj to atom-sorted cprj (atm_indx=atindx) From atom-sorted cprj to unsorted cprj (atm_indx=atindx1)
INPUTS
atm_indx(natom)=index table for atoms From unsorted cprj to atom-sorted cprj (atm_indx=atindx) From atom-sorted cprj to unsorted cprj (atm_indx=atindx1)
OUTPUT
SIDE EFFECTS
cprj(:,:) <type(pawcprj_type)>= cprj datastructure
SOURCE
1566 subroutine pawcprj_reorder(cprj,atm_indx) 1567 1568 !Arguments ------------------------------------ 1569 !scalars 1570 !arrays 1571 integer,intent(in) :: atm_indx(:) 1572 type(pawcprj_type),intent(inout) :: cprj(:,:) 1573 1574 !Local variables------------------------------- 1575 !scalars 1576 integer :: iexit,ii,jj,kk,n1atindx,n1cprj,n2cprj,ncpgr 1577 character(len=100) :: msg 1578 !arrays 1579 integer,allocatable :: nlmn(:) 1580 type(pawcprj_type),allocatable :: cprj_tmp(:,:) 1581 1582 ! ************************************************************************* 1583 1584 n1cprj=size(cprj,dim=1);n2cprj=size(cprj,dim=2) 1585 n1atindx=size(atm_indx,dim=1) 1586 if (n1cprj==0.or.n2cprj==0.or.n1atindx<=1) return 1587 1588 if (n1cprj/=n1atindx) then 1589 msg='wrong sizes (pawcprj_reorder)!' 1590 LIBPAW_BUG(msg) 1591 end if 1592 1593 !Nothing to do when the atoms are already sorted 1594 iexit=1;ii=0 1595 do while (iexit==1.and.ii<n1atindx) 1596 ii=ii+1 1597 if (atm_indx(ii)/=ii) iexit=0 1598 end do 1599 if (iexit==1) return 1600 1601 LIBPAW_ALLOCATE(nlmn,(n1cprj)) 1602 do ii=1,n1cprj 1603 nlmn(ii)=cprj(ii,1)%nlmn 1604 end do 1605 ncpgr=cprj(1,1)%ncpgr 1606 LIBPAW_DATATYPE_ALLOCATE(cprj_tmp,(n1cprj,n2cprj)) 1607 call pawcprj_alloc(cprj_tmp,ncpgr,nlmn) 1608 call pawcprj_copy(cprj,cprj_tmp) 1609 call pawcprj_free(cprj) 1610 1611 do jj=1,n2cprj 1612 do ii=1,n1cprj 1613 kk=atm_indx(ii) 1614 cprj(kk,jj)%nlmn=nlmn(ii) 1615 cprj(kk,jj)%ncpgr=ncpgr 1616 LIBPAW_ALLOCATE(cprj(kk,jj)%cp,(2,nlmn(ii))) 1617 cprj(kk,jj)%cp(:,:)=cprj_tmp(ii,jj)%cp(:,:) 1618 if (ncpgr>0) then 1619 LIBPAW_ALLOCATE(cprj(kk,jj)%dcp,(2,ncpgr,nlmn(ii))) 1620 cprj(kk,jj)%dcp(:,:,:)=cprj_tmp(ii,jj)%dcp(:,:,:) 1621 end if 1622 end do 1623 end do 1624 1625 call pawcprj_free(cprj_tmp) 1626 LIBPAW_DATATYPE_DEALLOCATE(cprj_tmp) 1627 LIBPAW_DEALLOCATE(nlmn) 1628 1629 end subroutine pawcprj_reorder
m_pawcprj/pawcprj_set_zero [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_set_zero
FUNCTION
Set to zero all arrays in a cprj datastructure
SIDE EFFECTS
cprj(:,:) <type(pawcprj_type)>= cprj datastructure
SOURCE
226 subroutine pawcprj_set_zero(cprj) 227 228 !Arguments ------------------------------------ 229 !scalars 230 !arrays 231 type(pawcprj_type),intent(inout) :: cprj(:,:) 232 233 !Local variables------------------------------- 234 !scalars 235 integer :: ii,jj,n1dim,n2dim 236 237 ! ************************************************************************* 238 239 n1dim=size(cprj,dim=1);n2dim=size(cprj,dim=2) 240 241 do jj=1,n2dim 242 do ii=1,n1dim 243 if (cprj(ii,jj)%nlmn>0) cprj(ii,jj)%cp(:,:)=zero 244 if (cprj(ii,jj)%ncpgr>0) cprj(ii,jj)%dcp(:,:,:)=zero 245 end do 246 end do 247 248 end subroutine pawcprj_set_zero
m_pawcprj/pawcprj_symkn [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_symkn
FUNCTION
compute cprj for a given band and k point based on cprj at a symmetry-related k point.
INPUTS
cprj_ikn (pawcprj_type) :: cprj for a single band and k point, typically a k point in the IBZ cprj_sym(4,nsym,natom) :: 1:3 shift, and 4 final atom, of symmetry isym operating on iatom (S^{-1}(R - t) = r0 + L, see symatm.F90 dimlmn(natom) :: ln dimension of each atom iband :: number of bands to treat, use -1 to treat all nband bands indlmn(6,lmnmax,ntypat) :: n,l,m dimensions for each atom type (see psps type) isym :: symmetry element used in current application itim :: 1 if time reversal also used, 0 else kpt(3) :: kpt vector used lmax :: max l value lmnmax :: max lmn value mband :: maximum number of bands natom :: number of atoms in cell nband :: number of bands in cprj_ikn nspinor :: number of spinors nsym :: total number of symmetry elements ntypat :: number of types of atoms typat(natom) :: type of each atom zarot(2*lmax+1,2*lmax+1,lmax+1,nsym) :: elements of rotation matrix for angular momentum states and symmetry operations. See m_paw_sphharm/setsym_ylm.
OUTPUT
cprj_fkn (pawcprj_type) :: cprj for a single band and k point where the k point is related to the input k point by a symmetry operation
SIDE EFFECTS
NOTES
This routine is based on M. Giantomassi's doctoral dissertation, formula 7.77. It is not clear whether it is implemented correctly for nonsymmorphic symmetries.
SOURCE
780 subroutine pawcprj_symkn(cprj_fkn,cprj_ikn,cprj_sym,dimlmn,iband,indlmn,& 781 & isym,itim,kpt,lmax,lmnmax,mband,natom,nband,nspinor,nsym,ntypat,& 782 & typat,zarot) 783 784 !Arguments--------------------------- 785 !scalars 786 integer,intent(in) :: iband,isym,itim,lmax,lmnmax,mband 787 integer,intent(in) :: natom,nband,nspinor,nsym,ntypat 788 789 !arrays 790 integer,intent(in) :: cprj_sym(4,nsym,natom),dimlmn(natom) 791 integer,intent(in) :: indlmn(6,lmnmax,ntypat),typat(natom) 792 real(dp),intent(in) :: kpt(3) 793 real(dp),intent(in) :: zarot(2*lmax+1,2*lmax+1,lmax+1,nsym) 794 type(pawcprj_type),intent(in) :: cprj_ikn(natom,mband*nspinor) 795 type(pawcprj_type),intent(inout) :: cprj_fkn(natom,mband*nspinor) !vz_i 796 797 !Local variables--------------------------- 798 !scalars 799 integer :: iatm,iatom, ibct, ibnd, ibsp, ibst, icpgr, iin, il, il0, im 800 integer :: ilmn, iln, iln0, ilpm, indexi, ispinor, itypat, jatm,jatom, mm, nlmn 801 real(dp) :: kdotL, phr, phi 802 !arrays 803 real(dp) :: rl(3), t1(2), t2(2) 804 805 ! ************************************************************************* 806 807 if (iband == -1) then 808 ibst = 1 809 ibnd = nband 810 else 811 ibst = iband 812 ibnd = iband 813 end if 814 815 do iatom = 1, natom 816 iatm=iatom 817 itypat = typat(iatom) 818 nlmn = dimlmn(iatm) 819 jatom = cprj_sym(4,isym,iatom) 820 jatm=jatom 821 rl(:) = cprj_sym(1:3,isym,iatom) 822 kdotL = dot_product(rl,kpt) 823 phr = cos(two_pi*kdotL) 824 phi = sin(two_pi*kdotL) 825 826 il0 = -1; iln0 = -1; indexi = 1 827 do ilmn = 1, nlmn 828 829 il = indlmn(1,ilmn,itypat) 830 im = indlmn(2,ilmn,itypat) 831 iin = indlmn(3,ilmn,itypat) 832 iln = indlmn(5,ilmn,itypat) 833 ilpm = 1 + il + im 834 if (iln /= iln0) indexi = indexi + 2*il0 + 1 835 836 do ibct = ibst, ibnd 837 838 do ispinor = 1, nspinor 839 840 ibsp = nspinor*(ibct-1) + ispinor 841 842 t1(:) = zero 843 do mm = 1, 2*il+1 844 t1(1) = t1(1) + zarot(mm,ilpm,il+1,isym)*cprj_ikn(jatm,ibsp)%cp(1,indexi+mm) 845 t1(2) = t1(2) + zarot(mm,ilpm,il+1,isym)*cprj_ikn(jatm,ibsp)%cp(2,indexi+mm) 846 end do 847 t2(1) = t1(1)*phr - t1(2)*phi 848 t2(2) = t1(2)*phr + t1(1)*phi 849 850 if (itim == 1) t2(2) = -t2(2) 851 852 cprj_fkn(iatm,ibsp)%cp(1,ilmn) = t2(1) 853 cprj_fkn(iatm,ibsp)%cp(2,ilmn) = t2(2) 854 855 ! do same transformations for gradients of cprj_ikn 856 ! note that ncpgr = 0 if no gradients present so this loop will not be executed 857 ! in this case 858 859 do icpgr = 1, cprj_ikn(jatom,ibsp)%ncpgr 860 t1(:) = zero 861 862 do mm = 1, 2*il+1 863 t1(1) = t1(1) + zarot(mm,ilpm,il+1,isym)*cprj_ikn(jatm,ibsp)%dcp(1,icpgr,indexi+mm) 864 t1(2) = t1(2) + zarot(mm,ilpm,il+1,isym)*cprj_ikn(jatm,ibsp)%dcp(2,icpgr,indexi+mm) 865 end do 866 867 t2(1) = t1(1)*phr - t1(2)*phi 868 t2(2) = t1(2)*phr + t1(1)*phi 869 870 if (itim == 1) t2(2) = -t2(2) 871 872 cprj_fkn(iatm,ibsp)%dcp(1,icpgr,ilmn) = t2(1) 873 cprj_fkn(iatm,ibsp)%dcp(2,icpgr,ilmn) = t2(2) 874 875 end do ! end loop over ncpgr 876 877 end do ! end loop over nspinor 878 879 end do ! end loop over bands 880 881 il0 = il; iln0 = iln 882 end do ! end loop over ilmn 883 end do ! end loop over atoms 884 885 end subroutine pawcprj_symkn
m_pawcprj/pawcprj_transpose [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_transpose
FUNCTION
Transpose a cprj datastructure FOR A GIVEN (K,SPIN) in order to change the parallel distribution from atom to band (or the contrary). At input, cprj is distributed over bands (or atoms); at output, it is distributed over atoms (or bands)
INPUTS
cprjin(n1indim,n2indim)<pawcprj_type>=the input cprj datastructure cprj_bandpp=number of bands to be treated simultaneoulsy by a processor natom=number of atoms in cell nband=number of bands nspinor=number of spinorial components spaceComm=MPI Communicator.
OUTPUT
cprjout(n1outdim,n2outdim)<pawcprj_type>=the output cprj datastructure with another distribution
NOTES
On the dimensions: To transfer cprj from band distribution to atom distribution, dimensions should be: n1indim =natom n2indim =nband/nproc*nspinor n1outdim=natom/nproc n2outdim=nband*nspinor To transfer cprj from atom distribution to band distribution, dimensions should be: n1indim =natom n2indim =nband/nproc*nspinor n1outdim=natom/nproc n2outdim=nband*nspinor
SOURCE
2329 subroutine pawcprj_transpose(cprjin,cprjout,cprj_bandpp,natom,nband,nspinor,spaceComm) 2330 2331 !Arguments------------------------------------- 2332 !scalars 2333 integer :: cprj_bandpp,natom,nband,nspinor,spaceComm 2334 !arrays 2335 type(pawcprj_type),intent(in) :: cprjin(:,:) 2336 type(pawcprj_type),intent(out) :: cprjout(:,:) 2337 2338 !Local variables------------------------------- 2339 !scalars 2340 integer :: bpp,buf_indx 2341 integer :: iashft,iatom,iatom_max_sd,iatom_max_rc,iatom_1,iatom_2,iatm1_sd,iatm1_rc,iatm2_sd,iatm2_rc 2342 integer :: ib,iband,iband_1,iband_2,iband_shift,iblock_atom,iblock_band,ibshft 2343 integer :: ierr,ip,ispinor,me,nba,nbb,nbnp_sd,nbnp_rc,ncpgr,nlmn,np 2344 integer :: rbufsize,sbufsize,size11,size12,size21,size22,transpose_mode 2345 character(len=100) :: msg 2346 !arrays 2347 integer,allocatable :: cprjsz_atom(:),cprjsz_block(:,:) 2348 integer,allocatable,target :: count_atom(:),count_band(:),displ_atom(:),displ_band(:) 2349 integer,pointer :: scount(:),sdispl(:),rcount(:),rdispl(:) 2350 real(dp),allocatable :: rbuf(:),sbuf(:) 2351 2352 ! ************************************************************************* 2353 2354 !MPI data 2355 me = xmpi_comm_rank(spaceComm) 2356 np = xmpi_comm_size(spaceComm) 2357 2358 !Nothing to do if nprocs=1 2359 if (np==1) then 2360 call pawcprj_copy(cprjin,cprjout) 2361 return 2362 end if 2363 2364 !Compute bloc sizes 2365 bpp=cprj_bandpp 2366 nba=natom/np;if (mod(natom,np)/=0) nba=nba+1 2367 nbb=nband/(np*bpp) 2368 2369 !Check sizes, select direction of transposition 2370 transpose_mode=0 2371 size11=size(cprjin,1);size12=size(cprjin,2) 2372 size21=size(cprjout,1);size22=size(cprjout,2) 2373 if (size11==natom.and.size12==nbb*bpp*nspinor.and.& 2374 & size21==nba.and.size22==nband*nspinor) then 2375 transpose_mode=1 2376 else if (size11==nba.and.size12==nband*nspinor.and.& 2377 & size21==natom.and.size22==nbb*bpp*nspinor) then 2378 else 2379 msg='wrong cprjin/cprjout sizes (pawcprj_transpose)!' 2380 LIBPAW_BUG(msg) 2381 end if 2382 2383 !Compute size of atom bloc (wr to cprj) 2384 LIBPAW_ALLOCATE(cprjsz_atom,(natom)) 2385 LIBPAW_ALLOCATE(cprjsz_block,(np,nba)) 2386 cprjsz_atom=0;cprjsz_block=0 2387 if (transpose_mode==1) then 2388 do iatom=1,natom 2389 cprjsz_atom(iatom)=2*cprjin(iatom,1)%nlmn*(1+cprjin(iatom,1)%ncpgr) 2390 end do 2391 else 2392 do iblock_atom=1,nba 2393 iatom=(iblock_atom-1)*np+1+me 2394 if (iatom<=natom) cprjsz_atom(iatom)=2*cprjin(iblock_atom,1)%nlmn*(1+cprjin(iblock_atom,1)%ncpgr) 2395 end do 2396 call xmpi_sum(cprjsz_atom,spaceComm,ierr) 2397 end if 2398 do iblock_atom=1,nba 2399 iashft=(iblock_atom-1)*np 2400 iatom_1=iashft+1;iatom_2=iashft+np 2401 if (iatom_1>natom) cycle 2402 if (iatom_2>natom) iatom_2=natom 2403 do iatom=iatom_1,iatom_2 2404 cprjsz_block(iatom-iashft,iblock_atom)=cprjsz_atom(iatom)+2 ! +2 for nlmn et ncpgr 2405 end do 2406 end do 2407 LIBPAW_DEALLOCATE(cprjsz_atom) 2408 2409 !Allocations for MPI_ALLTOALL 2410 LIBPAW_ALLOCATE(count_atom,(np)) 2411 LIBPAW_ALLOCATE(displ_atom,(np)) 2412 LIBPAW_ALLOCATE(count_band,(np)) 2413 LIBPAW_ALLOCATE(displ_band,(np)) 2414 2415 !Loop on blocks of bands 2416 do iblock_band=1,nbb !(note: np divides nband) 2417 ibshft=(iblock_band-1)*np*bpp 2418 iband_1=ibshft+1;iband_2=ibshft+np*bpp 2419 if (iband_1>nband.or.iband_2>nband) cycle ! for security 2420 2421 ! Loop on blocks of atoms 2422 do iblock_atom=1,nba 2423 iashft=(iblock_atom-1)*np 2424 iatom_1=iashft+1;iatom_2=iashft+np 2425 if (iatom_1>natom) cycle 2426 if (iatom_2>natom) iatom_2=natom 2427 2428 ! Computation of displacements and sizes of blocks when data are band-distributed 2429 count_band(1)=cprjsz_block(1,iblock_atom)*nspinor*bpp;displ_band(1)=0 2430 do ip=2,np 2431 count_band(ip)=cprjsz_block(ip,iblock_atom)*nspinor*bpp 2432 displ_band(ip)=displ_band(ip-1)+count_band(ip-1) 2433 end do 2434 2435 ! Computation of displacements and sizes of blocks when data are atom-distributed 2436 count_atom(1)=cprjsz_block(1+me,iblock_atom)*bpp*nspinor;displ_atom(1)=0 2437 do ip=2,np 2438 count_atom(ip)=count_atom(1) 2439 displ_atom(ip)=displ_atom(ip-1)+count_atom(ip-1) 2440 end do 2441 2442 ! According to transposition mode, select 2443 ! - displacements and sizes of blocks 2444 ! - shifts in arrays 2445 if (transpose_mode==1) then 2446 scount => count_band ; sdispl => displ_band 2447 rcount => count_atom ; rdispl => displ_atom 2448 nbnp_sd=bpp;nbnp_rc=np*bpp 2449 iatm1_sd=iatom_1;iatm2_sd=iatom_2 2450 iatm1_rc=iblock_atom;iatm2_rc=iatm1_rc 2451 iatom_max_sd=iatom_2;iatom_max_rc=iashft+1+me 2452 else 2453 scount => count_atom ; sdispl => displ_atom 2454 rcount => count_band ; rdispl => displ_band 2455 nbnp_sd=np*bpp;nbnp_rc=bpp 2456 iatm1_sd=iblock_atom;iatm2_sd=iatm1_sd 2457 iatm1_rc=iatom_1;iatm2_rc=iatom_2 2458 iatom_max_sd=iashft+1+me;iatom_max_rc=iatom_2 2459 end if 2460 2461 ! Allocation of buffers 2462 sbufsize=sdispl(np)+scount(np) 2463 rbufsize=rdispl(np)+rcount(np) 2464 LIBPAW_ALLOCATE(sbuf,(sbufsize)) 2465 LIBPAW_ALLOCATE(rbuf,(rbufsize)) 2466 2467 ! Coying of input cprj to buffer for sending 2468 buf_indx=0 2469 iband_shift=(iblock_band-1)*nbnp_sd-1 2470 if (iatom_max_sd<=natom) then 2471 do iatom=iatm1_sd,iatm2_sd 2472 do ib=1,nbnp_sd 2473 iband=(iband_shift+ib)*nspinor 2474 do ispinor=1,nspinor 2475 iband=iband+1 2476 nlmn=cprjin(iatom,iband)%nlmn;ncpgr=cprjin(iatom,iband)%ncpgr 2477 sbuf(buf_indx+1)=dble(nlmn) ;buf_indx=buf_indx+1 2478 sbuf(buf_indx+1)=dble(ncpgr);buf_indx=buf_indx+1 2479 sbuf(buf_indx+1:buf_indx+2*nlmn)=reshape(cprjin(iatom,iband)%cp(1:2,1:nlmn),(/2*nlmn/)) 2480 buf_indx=buf_indx+2*nlmn 2481 if (ncpgr>0) then 2482 sbuf(buf_indx+1:buf_indx+2*ncpgr*nlmn)=reshape(cprjin(iatom,iband)%dcp(1:2,1:ncpgr,1:nlmn),(/2*ncpgr*nlmn/)) 2483 buf_indx=buf_indx+2*ncpgr*nlmn 2484 end if 2485 end do 2486 end do 2487 end do 2488 end if 2489 if (buf_indx/=sbufsize) then 2490 msg='wrong buffer size for sending (pawcprj_transpose)!' 2491 LIBPAW_BUG(msg) 2492 end if 2493 2494 ! Main call to MPI_ALLTOALL 2495 call xmpi_alltoallv(sbuf,scount,sdispl,rbuf,rcount,rdispl,spaceComm,ierr) 2496 2497 ! Retrieving of output cprj for received buffer 2498 buf_indx=0 2499 iband_shift=(iblock_band-1)*nbnp_rc-1 2500 if (iatom_max_rc<=natom) then 2501 do iatom=iatm1_rc,iatm2_rc 2502 do ib=1,nbnp_rc 2503 iband=(iband_shift+ib)*nspinor 2504 do ispinor=1,nspinor 2505 iband=iband+1 2506 nlmn =int(rbuf(buf_indx+1));buf_indx=buf_indx+1 2507 ncpgr=int(rbuf(buf_indx+1));buf_indx=buf_indx+1 2508 cprjout(iatom,iband)%nlmn=nlmn;cprjout(iatom,iband)%ncpgr=ncpgr 2509 cprjout(iatom,iband)%cp(1:2,1:nlmn)=reshape(rbuf(buf_indx+1:buf_indx+2*nlmn),(/2,nlmn/)) 2510 buf_indx=buf_indx+2*nlmn 2511 if (ncpgr>0) then 2512 cprjout(iatom,iband)%dcp(1:2,1:ncpgr,1:nlmn)=reshape(rbuf(buf_indx+1:buf_indx+2*nlmn*ncpgr),(/2,ncpgr,nlmn/)) 2513 buf_indx=buf_indx+2*nlmn*ncpgr 2514 end if 2515 end do 2516 end do 2517 end do 2518 else 2519 cprjout(iatom,iband)%nlmn=0;cprjout(iatom,iband)%ncpgr=0 2520 end if 2521 if (buf_indx/=rbufsize) then 2522 msg='wrong buffer size for receiving (pawcprj_transpose)!' 2523 LIBPAW_BUG(msg) 2524 end if 2525 2526 ! Deallocation of buffers 2527 LIBPAW_DEALLOCATE(sbuf) 2528 LIBPAW_DEALLOCATE(rbuf) 2529 2530 ! End of loops 2531 end do ! do iblock_atom 2532 end do ! do iblock_atom 2533 2534 !Free memory 2535 LIBPAW_DEALLOCATE(count_atom) 2536 LIBPAW_DEALLOCATE(displ_atom) 2537 LIBPAW_DEALLOCATE(count_band) 2538 LIBPAW_DEALLOCATE(displ_band) 2539 LIBPAW_DEALLOCATE(cprjsz_block) 2540 nullify(scount,rcount,sdispl,rdispl) 2541 2542 end subroutine pawcprj_transpose
m_pawcprj/pawcprj_type [ Types ]
[ Top ] [ m_pawcprj ] [ Types ]
NAME
pawcprj_type
FUNCTION
This structured datatype contains <p_lmn|Cnk> projected scalars and derivatives where |p_lmn> are non-local projectors for a given atom |Cnk> is a wave function Used only for PAW calculations.
SOURCE
52 type,public :: pawcprj_type 53 54 !Integer scalars 55 56 integer :: ncpgr=0 57 ! Number of gradients of cp=<p_lmn|Cnk> 58 59 integer :: nlmn=0 60 ! Number of (l,m,n) non-local projectors 61 62 !Real (real(dp)) arrays 63 64 real(dp), allocatable :: cp (:,:) 65 ! cp(2,nlmn) 66 ! <p_lmn|Cnk> projected scalars for a given atom and wave function 67 68 real(dp), allocatable :: dcp (:,:,:) 69 ! dcp(2,ncpgr,nlmn) 70 ! derivatives of <p_lmn|Cnk> projected scalars for a given atom and wave function 71 72 end type pawcprj_type 73 74 !public procedures. 75 public :: pawcprj_alloc ! Allocation 76 public :: pawcprj_free ! Deallocation 77 public :: pawcprj_set_zero ! Set to zero all arrays in a cprj datastructure 78 public :: pawcprj_copy ! Copy a cprj datastructure into another 79 public :: pawcprj_axpby ! cprjy(:,:) <- alpha.cprjx(:,:)+beta.cprjy(:,:) 80 public :: pawcprj_zaxpby ! cprjy(:,:) <- alpha.cprjx(:,:)+beta.cprjy(:,:), alpha and beta are COMPLEX scalars 81 public :: pawcprj_projbd ! cprjy(:,:) <- alpha.cprjx(:,:)+beta.cprjy(:,:), alpha and beta are COMPLEX scalars 82 public :: pawcprj_conjg ! cprj(:,:) <- conjugate(cprj(:,:)) 83 public :: pawcprj_symkn ! construct cprj from that at a symmetry related k point 84 public :: pawcprj_lincom ! Compute a LINear COMbination of cprj datastructure: 85 public :: pawcprj_output ! Output a cprj. Useful for debugging. 86 public :: pawcprj_get ! Read the cprj for a given k-point from memory or from a temporary file 87 public :: pawcprj_put ! Write the cprj for a given set of (n,k) into memory or into a temporary file 88 public :: pawcprj_reorder ! Change the order of a cprj datastructure 89 public :: pawcprj_mpi_allgather ! Perform MPI_ALLGATHER on a pawcprj_type inside a MPI communicator. 90 public :: pawcprj_bcast ! Broadcast a pawcprj_type from master to all nodes inside a MPI communicator. 91 public :: pawcprj_transpose ! Transpose a cprj datastructure FOR A GIVEN (K,SPIN) 92 public :: pawcprj_gather_spin ! Collect spin distributed cprjs. 93 public :: pawcprj_mpi_exch ! Exchange a pawcprj_type between two processors inside a MPI communicator. 94 public :: pawcprj_mpi_send ! Send a pawcprj_type inside a MPI communicator. 95 public :: pawcprj_mpi_recv ! Receive a pawcprj_type inside a MPI communicator. 96 public :: pawcprj_mpi_sum ! Perform MPI_SUM on a pawcprj_type inside a MPI communicator. 97 public :: pawcprj_getdim ! Returns the number of lmn components in the <p_{lmn}^i|\psi> for the i-th atom. 98 public :: paw_overlap ! Compute the onsite contribution to the overlap between two states. 99 public :: pawcprj_pack ! Copy data from a cprj to a simple real buffer 100 public :: pawcprj_unpack ! Copy data from a simple real buffer to a cprj
m_pawcprj/pawcprj_unpack [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_unpack
FUNCTION
Unpack structured data from a simple buffer
INPUTS
nlmn(natom)=Number of nlm partial waves for each atom. ncpgr = number of gradients in cprj_in buffer = the data to be unpacked, dim : (2, n2dim*sum(nlmn)) [buffer_gr] = if present the gradient data to be unpacked, dim : (2, ncpgr, n2dim*sum(nlmn))
OUTPUT
cprj=The datatype unpacked
SOURCE
2905 subroutine pawcprj_unpack(nlmn,cprj,buffer,buffer_gr) 2906 2907 !Arguments ------------------------------------ 2908 !scalars 2909 !arrays 2910 integer,intent(in) :: nlmn(:) 2911 real(dp),intent(in) :: buffer(:,:) 2912 real(dp),intent(in),optional :: buffer_gr(:,:,:) 2913 type(pawcprj_type),intent(inout) :: cprj(:,:) 2914 2915 !Local variables------------------------------- 2916 !scalars 2917 integer :: natom,n2buffer,ncpgr,n2dim 2918 integer :: iat,jj,n1dim,nn 2919 integer :: ipck 2920 character(len=100) :: msg 2921 !arrays 2922 2923 ! ************************************************************************* 2924 2925 natom=size(nlmn,dim=1) 2926 n2buffer=size(buffer,dim=2) 2927 n1dim=size(cprj,dim=1) 2928 n2dim=size(cprj,dim=2) 2929 2930 if (natom/=n1dim) then 2931 msg='size mismatch in natom (pawcprj_unpack)!' 2932 LIBPAW_BUG(msg) 2933 end if 2934 if (n2dim*SUM(nlmn)/=n2buffer) then 2935 msg='size mismatch in dim=2 (pawcprj_unpack)!' 2936 LIBPAW_BUG(msg) 2937 end if 2938 ncpgr=0 2939 if (present(buffer_gr)) then 2940 ncpgr=size(buffer_gr,dim=2) 2941 end if 2942 2943 !=== Unpack buffers into cprj === 2944 ipck=0 2945 do jj=1,n2dim 2946 do iat=1,natom 2947 nn=nlmn(iat) 2948 cprj(iat,jj)%cp(:,1:nn)=buffer(:,ipck+1:ipck+nn) 2949 if (ncpgr/=0) then 2950 cprj(iat,jj)%dcp(:,:,1:nn)=buffer_gr(:,:,ipck+1:ipck+nn) 2951 end if 2952 ipck=ipck+nn 2953 end do 2954 end do 2955 2956 end subroutine pawcprj_unpack 2957 2958 end module m_pawcprj
m_pawcprj/pawcprj_zaxpby [ Functions ]
[ Top ] [ m_pawcprj ] [ Functions ]
NAME
pawcprj_zaxpby
FUNCTION
Apply ZAXPBY (blas-like) operation with 2 cprj datastructures: cprjy(:,:) <- alpha.cprjx(:,:)+beta.cprjy(:,:) alpha and beta are COMPLEX scalars
INPUTS
alpha(2),beta(2)= alpha,beta COMPLEX factors cprjx(:,:) <type(pawcprj_type)>= input cprjx datastructure
SIDE EFFECTS
cprjy(:,:) <type(pawcprj_type)>= input/output cprjy datastructure
SOURCE
492 subroutine pawcprj_zaxpby(alpha,beta,cprjx,cprjy) 493 494 !Arguments ------------------------------------ 495 !scalars 496 real(dp),intent(in) :: alpha(2),beta(2) 497 !arrays 498 type(pawcprj_type),intent(in) :: cprjx(:,:) 499 type(pawcprj_type),intent(inout) :: cprjy(:,:) 500 501 !Local variables------------------------------- 502 !scalars 503 integer :: ii,jj,kk,ll,n1dimx,n1dimy,n2dimx,n2dimy,ncpgrx,ncpgry,nlmn 504 real(dp) :: cp1,cp2,norma,normb 505 character(len=500) :: msg 506 507 ! ************************************************************************* 508 509 norma=alpha(1)**2+alpha(2)**2 510 normb=beta(1) **2+beta(2) **2 511 n1dimy=size(cprjy,dim=1);n2dimy=size(cprjy,dim=2);ncpgry=cprjy(1,1)%ncpgr 512 if (norma>tol16*tol16) then 513 n1dimx=size(cprjx,dim=1);n2dimx=size(cprjx,dim=2);ncpgrx=cprjx(1,1)%ncpgr 514 msg = "" 515 if (n1dimx/=n1dimy) msg = TRIM(msg)//"Error in pawcprj_zaxpby: n1 wrong sizes !"//ch10 516 if (n2dimx/=n2dimy) msg = TRIM(msg)//"Error in pawcprj_zaxpby: n2 wrong sizes !"//ch10 517 if (ncpgrx/=ncpgry) msg = TRIM(msg)//"Error in pawcprj_zaxpby: ncpgr wrong sizes !"//ch10 518 if (LEN_TRIM(msg) > 0) then 519 LIBPAW_ERROR(msg) 520 end if 521 end if 522 523 if (norma<=tol16*tol16) then 524 do jj=1,n2dimy 525 do ii=1,n1dimy 526 nlmn=cprjy(ii,jj)%nlmn 527 do kk=1,nlmn 528 cp1=beta(1)*cprjy(ii,jj)%cp(1,kk)-beta(2)*cprjy(ii,jj)%cp(2,kk) 529 cp2=beta(1)*cprjy(ii,jj)%cp(2,kk)+beta(2)*cprjy(ii,jj)%cp(1,kk) 530 cprjy(ii,jj)%cp(1,kk)=cp1 531 cprjy(ii,jj)%cp(2,kk)=cp2 532 end do 533 end do 534 end do 535 if (ncpgry>0) then 536 do jj=1,n2dimy 537 do ii=1,n1dimy 538 nlmn=cprjy(ii,jj)%nlmn 539 do kk=1,nlmn 540 do ll=1,ncpgry 541 cp1=beta(1)*cprjy(ii,jj)%dcp(1,ll,kk)-beta(2)*cprjy(ii,jj)%dcp(2,ll,kk) 542 cp2=beta(1)*cprjy(ii,jj)%dcp(2,ll,kk)+beta(2)*cprjy(ii,jj)%dcp(1,ll,kk) 543 cprjy(ii,jj)%dcp(1,ll,kk)=cp1 544 cprjy(ii,jj)%dcp(2,ll,kk)=cp2 545 end do 546 end do 547 end do 548 end do 549 end if 550 else if (normb<=tol16*tol16) then 551 do jj=1,n2dimx 552 do ii=1,n1dimx 553 nlmn=cprjx(ii,jj)%nlmn 554 cprjy(ii,jj)%nlmn=nlmn 555 do kk=1,nlmn 556 cprjy(ii,jj)%cp(1,kk)=alpha(1)*cprjx(ii,jj)%cp(1,kk)-alpha(2)*cprjx(ii,jj)%cp(2,kk) 557 cprjy(ii,jj)%cp(2,kk)=alpha(1)*cprjx(ii,jj)%cp(2,kk)+alpha(2)*cprjx(ii,jj)%cp(1,kk) 558 end do 559 end do 560 end do 561 if (ncpgrx>0) then 562 do jj=1,n2dimx 563 do ii=1,n1dimx 564 nlmn=cprjx(ii,jj)%nlmn 565 do kk=1,nlmn 566 cprjy(ii,jj)%dcp(1,1:ncpgrx,kk)=alpha(1)*cprjx(ii,jj)%dcp(1,1:ncpgrx,kk) & 567 & -alpha(2)*cprjx(ii,jj)%dcp(2,1:ncpgrx,kk) 568 cprjy(ii,jj)%dcp(2,1:ncpgrx,kk)=alpha(1)*cprjx(ii,jj)%dcp(2,1:ncpgrx,kk) & 569 & +alpha(2)*cprjx(ii,jj)%dcp(1,1:ncpgrx,kk) 570 end do 571 end do 572 end do 573 end if 574 ! else if (abs(beta(1)-one)<tol16.and.abs(beta(2))<tol16) then 575 ! do jj=1,n2dimx 576 ! do ii=1,n1dimx 577 ! nlmn=cprjx(ii,jj)%nlmn 578 ! cprjy(ii,jj)%nlmn =nlmn 579 ! do kk=1,nlmn 580 ! cp1=cprjy(ii,jj)%cp(1,kk) 581 ! cp2=cprjy(ii,jj)%cp(2,kk) 582 ! cp1=cp1+alpha(1)*cprjx(ii,jj)%cp(1,kk)-alpha(2)*cprjx(ii,jj)%cp(2,kk) 583 ! cp2=cp2+alpha(1)*cprjx(ii,jj)%cp(2,kk)+alpha(2)*cprjx(ii,jj)%cp(1,kk) 584 ! cprjy(ii,jj)%cp(1,kk)=cp1 585 ! cprjy(ii,jj)%cp(2,kk)=cp2 586 ! end do 587 ! end do 588 ! end do 589 ! if (ncpgrx>0) then 590 ! do jj=1,n2dimx 591 ! do ii=1,n1dimx 592 ! nlmn=cprjx(ii,jj)%nlmn 593 ! do kk=1,nlmn 594 ! do ll=1,ncpgrx 595 ! cp1=cprjy(ii,jj)%dcp(1,ll,kk) 596 ! cp2=cprjy(ii,jj)%dcp(2,ll,kk) 597 ! cp1=cp1+alpha(1)*cprjx(ii,jj)%dcp(1,ll,kk)-alpha(2)*cprjx(ii,jj)%dcp(2,ll,kk) 598 ! cp2=cp2+alpha(1)*cprjx(ii,jj)%dcp(2,ll,kk)+alpha(2)*cprjx(ii,jj)%dcp(1,ll,kk) 599 ! cprjy(ii,jj)%dcp(1,ll,kk)=cp1 600 ! cprjy(ii,jj)%dcp(2,ll,kk)=cp2 601 ! end do 602 ! end do 603 ! end do 604 ! end do 605 ! end if 606 else 607 do jj=1,n2dimx 608 do ii=1,n1dimx 609 nlmn=cprjx(ii,jj)%nlmn 610 cprjy(ii,jj)%nlmn =nlmn 611 do kk=1,nlmn 612 cp1=alpha(1)*cprjx(ii,jj)%cp(1,kk)-alpha(2)*cprjx(ii,jj)%cp(2,kk) & 613 & +beta(1) *cprjy(ii,jj)%cp(1,kk)-beta(2) *cprjy(ii,jj)%cp(2,kk) 614 cp2=alpha(1)*cprjx(ii,jj)%cp(2,kk)+alpha(2)*cprjx(ii,jj)%cp(1,kk) & 615 & +beta(1) *cprjy(ii,jj)%cp(2,kk)+beta(2) *cprjy(ii,jj)%cp(1,kk) 616 ! cp1=beta(1) *cprjy(ii,jj)%cp(1,kk)-beta(2) *cprjy(ii,jj)%cp(2,kk) 617 ! cp1=cp1+alpha(1)*cprjx(ii,jj)%cp(1,kk)-alpha(2)*cprjx(ii,jj)%cp(2,kk) 618 ! cp2=beta(1) *cprjy(ii,jj)%cp(2,kk)+beta(2) *cprjy(ii,jj)%cp(1,kk) 619 ! cp2=cp2+alpha(1)*cprjx(ii,jj)%cp(2,kk)+alpha(2)*cprjx(ii,jj)%cp(1,kk) 620 cprjy(ii,jj)%cp(1,kk)=cp1 621 cprjy(ii,jj)%cp(2,kk)=cp2 622 end do 623 end do 624 end do 625 if (ncpgrx>0) then 626 do jj=1,n2dimx 627 do ii=1,n1dimx 628 nlmn=cprjx(ii,jj)%nlmn 629 do kk=1,nlmn 630 do ll=1,ncpgrx 631 cp1=alpha(1)*cprjx(ii,jj)%dcp(1,ll,kk)-alpha(2)*cprjx(ii,jj)%dcp(2,ll,kk) & 632 & +beta(1) *cprjy(ii,jj)%dcp(1,ll,kk)-beta(2) *cprjy(ii,jj)%dcp(2,ll,kk) 633 cp2=alpha(1)*cprjx(ii,jj)%dcp(2,ll,kk)+alpha(2)*cprjx(ii,jj)%dcp(1,ll,kk) & 634 & +beta(1) *cprjy(ii,jj)%dcp(2,ll,kk)+beta(2) *cprjy(ii,jj)%dcp(1,ll,kk) 635 ! cp1=beta(1) *cprjy(ii,jj)%dcp(1,ll,kk)-beta(2) *cprjy(ii,jj)%dcp(2,ll,kk) 636 ! cp1=cp1+alpha(1)*cprjx(ii,jj)%dcp(1,ll,kk)-alpha(2)*cprjx(ii,jj)%dcp(2,ll,kk) 637 ! cp2=beta(1) *cprjy(ii,jj)%dcp(2,ll,kk)+beta(2) *cprjy(ii,jj)%dcp(1,ll,kk) 638 ! cp2=cp2+alpha(1)*cprjx(ii,jj)%dcp(2,ll,kk)+alpha(2)*cprjx(ii,jj)%dcp(1,ll,kk) 639 cprjy(ii,jj)%dcp(1,ll,kk)=cp1 640 cprjy(ii,jj)%dcp(2,ll,kk)=cp2 641 end do 642 end do 643 end do 644 end do 645 end if 646 end if 647 648 end subroutine pawcprj_zaxpby