TABLE OF CONTENTS
ABINIT/align_u_matrices [ Functions ]
NAME
align_u_matrices
FUNCTION
INPUTS
OUTPUT
PARENTS
xcart2deloc
CHILDREN
SOURCE
282 subroutine align_u_matrices(natom,ninternal,u_matrix,u_matrix_old,s_matrix,f_eigs) 283 284 use defs_basis 285 use m_profiling_abi 286 287 !This section has been created automatically by the script Abilint (TD). 288 !Do not modify the following lines by hand. 289 #undef ABI_FUNC 290 #define ABI_FUNC 'align_u_matrices' 291 !End of the abilint section 292 293 implicit none 294 295 !Arguments ------------------------------------ 296 !scalars 297 integer,intent(in) :: ninternal,natom 298 !arrays 299 real(dp),intent(in) :: u_matrix_old(ninternal,3*(natom-1)) 300 real(dp),intent(inout) :: f_eigs(3*natom) 301 real(dp),intent(inout) :: s_matrix(3*natom,3*natom) 302 real(dp),intent(inout) :: u_matrix(ninternal,3*(natom-1)) 303 304 !Local variables ------------------------------ 305 !scalars 306 integer :: ii,iint1,imax 307 real(dp) :: ss 308 !arrays 309 integer :: eigv_flag(3*(natom-1)),eigv_ind(3*(natom-1)) 310 real(dp) :: tmps(3*natom,3*natom) 311 real(dp) :: tmpu(ninternal,3*(natom-1)) 312 real(dp) :: tmpf(3*natom) 313 314 !****************************************************************** 315 316 eigv_flag(:) = 0 317 eigv_ind(:) = 0 318 319 !just permit a change in sign 320 do iint1=1,3*(natom-1) 321 ss = zero 322 do ii=1,ninternal 323 ss = ss + u_matrix_old(ii,iint1)*u_matrix(ii,iint1) 324 end do 325 if (ss < -tol12) then 326 imax = -iint1 327 else 328 imax = iint1 329 end if 330 eigv_ind(iint1) = imax 331 eigv_flag(abs(imax)) = 1 332 end do 333 334 tmpu(:,:) = u_matrix 335 tmps(:,:) = s_matrix 336 tmpf(:) = f_eigs 337 !exchange eigenvectors... 338 do iint1=1,3*(natom-1) 339 ss = one 340 if (eigv_ind(iint1) < 0) ss = -one 341 342 imax = abs(eigv_ind(iint1)) 343 344 tmpu(:,imax) = ss*u_matrix(:,iint1) 345 346 tmps(:,imax+3) = ss*s_matrix(:,iint1+3) 347 348 tmpf(imax+3) = f_eigs(iint1+3) 349 end do 350 351 u_matrix(:,:) = tmpu(:,:) 352 s_matrix(:,:) = tmps(:,:) 353 f_eigs(:) = tmpf(:) 354 355 end subroutine align_u_matrices
ABINIT/calc_btinv_matrix [ Functions ]
NAME
calc_btinv_matrix
FUNCTION
COPYRIGHT
Copyright (C) 2008-2018 ABINIT group (MJV). This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
OUTPUT
PARENTS
xcart2deloc
NOTES
bt_inv_matrix is inverse transpose of the delocalized coordinate B matrix. b_matrix is the primitive internal B matrix
CHILDREN
SOURCE
173 subroutine calc_btinv_matrix(b_matrix,natom,ninternal,bt_inv_matrix,u_matrix) 174 175 use defs_basis 176 use m_profiling_abi 177 use m_errors 178 use m_linalg_interfaces 179 180 !This section has been created automatically by the script Abilint (TD). 181 !Do not modify the following lines by hand. 182 #undef ABI_FUNC 183 #define ABI_FUNC 'calc_btinv_matrix' 184 use interfaces_45_geomoptim, except_this_one => calc_btinv_matrix 185 !End of the abilint section 186 187 implicit none 188 189 !Arguments ------------------------------------ 190 integer,intent(in) :: ninternal,natom 191 real(dp),intent(in) :: b_matrix(ninternal,3*natom) 192 real(dp),intent(out) :: bt_inv_matrix(3*(natom-1),3*natom) 193 real(dp),intent(inout) :: u_matrix(ninternal,3*(natom-1)) 194 195 !Local variables ------------------------------------ 196 !scalars 197 integer :: ii,info,lwork 198 !arrays 199 real(dp) :: f_eigs(3*natom),f_matrix(3*natom,3*natom) 200 real(dp) :: s_matrix(3*natom,3*natom) 201 real(dp) :: s_red(3*natom,3*(natom-1)) 202 real(dp) :: u_matrix_old(ninternal,3*(natom-1)) 203 real(dp),allocatable :: work(:) 204 205 !****************************************************************** 206 207 !f matrix = B^{T} B 208 call dgemm('T','N',3*natom,3*natom,ninternal,one,& 209 & b_matrix,ninternal,b_matrix,ninternal,zero,f_matrix,3*natom) 210 211 lwork = max(1,3*3*natom-1) 212 ABI_ALLOCATE(work,(lwork)) 213 s_matrix(:,:) = f_matrix(:,:) 214 215 call dsyev('V','L',3*natom,s_matrix,3*natom,& 216 & f_eigs,work,lwork,info) 217 218 ABI_DEALLOCATE(work) 219 220 if (abs(f_eigs(1)) + abs(f_eigs(2)) + abs(f_eigs(3)) > tol10 ) then 221 write(std_out,*) 'Error: 3 lowest eigenvalues are not zero' 222 write(std_out,*) ' internal coordinates do NOT span the full degrees of freedom !' 223 write(std_out,'(6E16.6)') f_eigs 224 MSG_ERROR("Aborting now") 225 end if 226 if ( abs(f_eigs(4)) < tol10 ) then 227 write(std_out,*) 'Error: fourth eigenvalue is zero' 228 write(std_out,*) ' internal coordinates do NOT span the full degrees of freedom !' 229 write(std_out,'(6E16.6)') f_eigs 230 MSG_ERROR("Aborting now") 231 end if 232 233 !calculate U matrix from U = B * S_red * lambda^{-1/2} 234 do ii=1,3*(natom-1) 235 s_red(:,ii) = s_matrix(:,ii+3)/sqrt(f_eigs(ii+3)) 236 end do 237 238 u_matrix_old(:,:) = u_matrix(:,:) 239 240 call dgemm('N','N',ninternal,3*(natom-1),3*natom,one,& 241 & b_matrix,ninternal,s_red,3*natom,zero,u_matrix,ninternal) 242 243 244 !align eigenvectors, to preserve a form of continuity in convergences 245 !!!! eigenvalues are no longer in increasing order!!! but only s_red is reordered 246 !so that btinv is correct. 247 call align_u_matrices(natom,ninternal,u_matrix,u_matrix_old,s_matrix,f_eigs) 248 249 !calculate B_deloc^{-1} matrix for transformation of forces to deloc coord. 250 !(B^{T}_deloc)^{-1} = (B_deloc B^{T}_deloc)^{-1} B_deloc = lambda^{-3/2} S^{T} F 251 != ( S lambda^{3/2} )^{T} F 252 253 !! DEFINITION 254 !! real(dp),intent(out) :: bt_inv_matrix(3*(natom-1),3*natom) 255 256 !even better: B_deloc^{-1} = lambda^{-1/2} S^{T} 257 do ii=1,3*(natom-1) 258 ! s_red(:,ii) = s_matrix(:,ii+3)*sqrt(f_eigs(ii+3)) 259 bt_inv_matrix(ii,:) = s_matrix(:,ii+3)/sqrt(f_eigs(ii+3)) 260 end do 261 262 end subroutine calc_btinv_matrix
ABINIT/xcart2deloc [ Functions ]
NAME
xcart2deloc
FUNCTION
Calculate values of delocalized coordinates as a function of cartesian ones. First primitive internals, then B matrix, then F, then U then delocalized internals.
COPYRIGHT
Copyright (C) 2003-2018 ABINIT group (MVer) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt.
INPUTS
deloc <type(delocint)>=Important variables for | pred_delocint | | nang = Number of angles | nbond = Number of bonds | ncart = Number of cartesian directions | (used for constraints) | ndihed = Number of dihedrals | nrshift = Dimension of rshift | ninternal= Number of internal coordinates | ninternal=nbond+nang+ndihed+ncart | | angs(2,3,nang) = Indexes to characterize angles | bonds(2,2,nbond)= For a bond between iatom and jatom | bonds(1,1,nbond) = iatom | bonds(2,1,nbond) = icenter | bonds(1,2,nbond) = jatom | bonds(2,2,nbond) = irshift | carts(2,ncart) = Index of total primitive internal, | and atom (carts(2,:)) | dihedrals(2,4,ndihed)= Indexes to characterize dihedrals | | rshift(3,nrshift)= Shift in xred that must be done to find | all neighbors of a given atom within a | given number of neighboring shells natom = Number of atoms rprimd(3,3) = Dimensional real space primitive translations (bohr) xcart(3,natom) = Cartesian coordinates of atoms (bohr)
OUTPUT
bt_inv_matrix(3*(natom-1),3*natom) = Inverse of B^{T} matrix deloc_int(3*(natom-1)) = Delocalized internal coordinates prim_int(ninternal) = Primitive internal coordinates
SIDE EFFECTS
u_matrix(ninternal,3*(natom-1)) = Eigenvectors of BB^T matrix
NOTES
PARENTS
deloc2xcart,pred_delocint,xfh_recover_deloc
CHILDREN
SOURCE
67 #if defined HAVE_CONFIG_H 68 #include "config.h" 69 #endif 70 71 #include "abi_common.h" 72 73 subroutine xcart2deloc(deloc,natom,rprimd,xcart,& 74 & bt_inv_matrix,u_matrix,deloc_int,prim_int) 75 76 use defs_basis 77 use m_errors 78 use m_profiling_abi 79 use m_abimover 80 use m_linalg_interfaces 81 82 !This section has been created automatically by the script Abilint (TD). 83 !Do not modify the following lines by hand. 84 #undef ABI_FUNC 85 #define ABI_FUNC 'xcart2deloc' 86 use interfaces_45_geomoptim, except_this_one => xcart2deloc 87 !End of the abilint section 88 89 implicit none 90 91 !Arguments ------------------------------------ 92 !scalars 93 integer,intent(in) :: natom 94 type(delocint),intent(inout) :: deloc 95 !arrays 96 real(dp),intent(in) :: rprimd(3,3),xcart(3,natom) 97 real(dp),intent(inout) :: u_matrix(deloc%ninternal,3*(natom-1)) 98 real(dp),intent(out) :: bt_inv_matrix(3*(natom-1),3*natom) 99 real(dp),intent(out) :: deloc_int(3*(natom-1)) 100 real(dp),intent(out) :: prim_int(deloc%ninternal) 101 102 !Local variables------------------------------- 103 !scalars 104 integer :: ii 105 logical :: DEBUG=.FALSE. 106 !arrays 107 real(dp) :: b_matrix(deloc%ninternal,3*natom) 108 109 ! ****************************************************************** 110 111 call calc_prim_int(deloc,natom,rprimd,xcart,prim_int) 112 if (DEBUG)then 113 write(std_out,*) 'Primitive Internals' 114 do ii=1,deloc%ninternal 115 write(std_out,*) prim_int(ii) 116 end do 117 end if 118 119 call calc_b_matrix(deloc,natom,rprimd,xcart,b_matrix) 120 if (DEBUG)then 121 write(std_out,*) 'B Matrix' 122 do ii=1,deloc%ninternal 123 write(std_out,*) b_matrix(:,ii) 124 end do 125 end if 126 127 call calc_btinv_matrix(b_matrix,natom,deloc%ninternal,& 128 & bt_inv_matrix,u_matrix) 129 if (DEBUG)then 130 write(std_out,*) 'BT Inverse Matrix' 131 do ii=1,3*natom 132 write(std_out,*) bt_inv_matrix(:,ii) 133 end do 134 end if 135 136 !calculate value of delocalized internals 137 138 call dgemv('T',deloc%ninternal,3*(natom-1),one,& 139 & u_matrix,deloc%ninternal,prim_int,1,zero,deloc_int,1) 140 141 end subroutine xcart2deloc