114#include "./base/base_uses.f90"
122 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_loc_methods'
163 zij_fm_set, para_env, cell, weights, ispin, print_loc_section, &
164 restricted, nextra, nmo, vectors_2, guess_mos)
166 INTEGER,
INTENT(IN) :: method
169 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: op_sm_set
170 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(IN) :: zij_fm_set
173 REAL(
dp),
DIMENSION(:) :: weights
174 INTEGER,
INTENT(IN) :: ispin
176 INTEGER :: restricted
177 INTEGER,
INTENT(IN),
OPTIONAL :: nextra, nmo
178 TYPE(
cp_fm_type),
INTENT(IN),
OPTIONAL :: vectors_2, guess_mos
180 CHARACTER(len=*),
PARAMETER :: routinen =
'optimize_loc_berry'
182 INTEGER :: handle, idim1, idim2, max_iter, nao, &
183 nmoloc, out_each, output_unit, sweeps
184 LOGICAL :: converged, crazy_use_diag, &
185 do_jacobi_refinement, my_do_mixed
186 REAL(
dp) :: crazy_scale, eps_localization, &
187 max_crazy_angle, start_time, &
189 TYPE(
cp_cfm_type),
ALLOCATABLE,
DIMENSION(:, :) :: my_zij_cfm_set
193 NULLIFY (complex_vectors)
195 CALL timeset(routinen, handle)
198 extension=
".locInfo")
203 CALL zij_matrix(vectors, op_sm_set, zij_fm_set)
205 max_iter = qs_loc_env%localized_wfn_control%max_iter
206 max_crazy_angle = qs_loc_env%localized_wfn_control%max_crazy_angle
207 crazy_use_diag = qs_loc_env%localized_wfn_control%crazy_use_diag
208 crazy_scale = qs_loc_env%localized_wfn_control%crazy_scale
209 eps_localization = qs_loc_env%localized_wfn_control%eps_localization
210 out_each = qs_loc_env%localized_wfn_control%out_each
211 target_time = qs_loc_env%target_time
212 start_time = qs_loc_env%start_time
213 do_jacobi_refinement = qs_loc_env%localized_wfn_control%jacobi_refinement
214 my_do_mixed = qs_loc_env%localized_wfn_control%do_mixed
217 ispin, print_loc_section, zij=zij_fm_set, only_initial_out=.true.)
223 CALL jacobi_rotations(weights, zij_fm_set, vectors, para_env, max_iter=max_iter, &
224 eps_localization=eps_localization, sweeps=sweeps, &
225 out_each=out_each, target_time=target_time, start_time=start_time, &
226 restricted=restricted)
228 IF (my_do_mixed)
THEN
230 IF (
PRESENT(guess_mos))
THEN
232 eps_localization, sweeps, out_each, nextra, &
233 qs_loc_env%localized_wfn_control%do_cg_po, &
234 nmo=nmo, vectors_2=vectors_2, mos_guess=guess_mos)
237 eps_localization, sweeps, out_each, nextra, &
238 qs_loc_env%localized_wfn_control%do_cg_po, &
239 nmo=nmo, vectors_2=vectors_2)
243 eps_localization, sweeps, out_each, 0, &
244 qs_loc_env%localized_wfn_control%do_cg_po)
247 cpabort(
"GAPO works only with STATES MIXED")
253 CALL zij_matrix(vectors, op_sm_set, zij_fm_set)
254 IF (do_jacobi_refinement)
THEN
257 ispin, print_loc_section, zij=zij_fm_set, only_initial_out=.true.)
258 CALL jacobi_rotations(weights, zij_fm_set, vectors, para_env, max_iter=max_iter, &
259 eps_localization=eps_localization, sweeps=sweeps, &
260 out_each=out_each, target_time=target_time, start_time=start_time, &
261 restricted=restricted)
264 CALL crazy_rotations(weights, zij_fm_set, vectors, max_iter=max_iter, max_crazy_angle=max_crazy_angle, &
265 crazy_scale=crazy_scale, crazy_use_diag=crazy_use_diag, &
266 eps_localization=eps_localization, iterations=sweeps, converged=converged)
268 IF (.NOT. converged)
THEN
269 IF (qs_loc_env%localized_wfn_control%jacobi_fallback)
THEN
270 IF (output_unit > 0)
WRITE (output_unit,
"(T4,A,I6,/,T4,A)") &
271 " Crazy Wannier localization not converged after ", sweeps, &
272 " iterations, switching to jacobi rotations"
273 CALL jacobi_rotations(weights, zij_fm_set, vectors, para_env, max_iter=max_iter, &
274 eps_localization=eps_localization, sweeps=sweeps, &
275 out_each=out_each, target_time=target_time, start_time=start_time, &
276 restricted=restricted)
278 IF (output_unit > 0)
WRITE (output_unit,
"(T4,A,I6,/,T4,A)") &
279 " Crazy Wannier localization not converged after ", sweeps, &
280 " iterations, and jacobi_fallback switched off "
284 CALL direct_mini(weights, zij_fm_set, vectors, max_iter=max_iter, &
285 eps_localization=eps_localization, iterations=sweeps)
287 IF (.NOT. cell%orthorhombic)
THEN
288 cpabort(
"Non-orthorhombic cell with the selected method NYI")
292 CALL zij_matrix(vectors, op_sm_set, zij_fm_set)
295 ALLOCATE (my_zij_cfm_set(
SIZE(zij_fm_set, 1),
SIZE(zij_fm_set, 2)))
296 DO idim1 = 1,
SIZE(zij_fm_set, 1)
297 DO idim2 = 1,
SIZE(zij_fm_set, 2)
298 CALL cp_cfm_create(my_zij_cfm_set(idim1, idim2), zij_fm_set(idim1, idim2)%matrix_struct)
299 CALL cp_fm_to_cfm(msourcer=zij_fm_set(idim1, idim2), mtarget=my_zij_cfm_set(idim1, idim2))
302 ALLOCATE (complex_vectors)
304 CALL cp_fm_to_cfm(msourcer=vectors, mtarget=complex_vectors)
306 CALL cardoso_souloumiac(weights, my_zij_cfm_set, max_iter, eps_localization, sweeps, &
307 out_each, complex_vectors)
311 DEALLOCATE (complex_vectors)
312 DO idim1 = 1,
SIZE(my_zij_cfm_set, 1)
313 DO idim2 = 1,
SIZE(my_zij_cfm_set, 2)
317 DEALLOCATE (my_zij_cfm_set)
320 IF (output_unit > 0)
THEN
321 WRITE (output_unit,
'(T4,A,I6,A)')
" No MOS localization applied "
324 cpabort(
"Unknown localization method")
326 IF (output_unit > 0)
THEN
327 IF (sweeps <= max_iter)
WRITE (output_unit,
'(T4,A,I3,A,I6,A)')
" Localization for spin ", ispin, &
328 " converged in ", sweeps,
" iterations"
332 ispin, print_loc_section, zij=zij_fm_set)
336 CALL timestop(handle)
354 ispin, print_loc_section)
356 INTEGER,
INTENT(IN) :: method
359 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(IN) :: zij_fm_set
360 INTEGER,
INTENT(IN) :: ispin
363 CHARACTER(len=*),
PARAMETER :: routinen =
'optimize_loc_pipek'
365 INTEGER :: handle, iatom, idim1, idim2, isgf, ldz, &
366 max_iter, nao, natom, ncol, nmoloc, &
367 out_each, output_unit, sweeps
368 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: first_sgf, last_sgf, nsgf
370 TYPE(
cp_cfm_type),
POINTER :: complex_vectors, my_zij_cfm_set(:, :)
378 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
380 NULLIFY (mols, ov_fm)
382 CALL timeset(routinen, handle)
385 extension=
".locInfo")
387 NULLIFY (particle_set)
395 CALL get_qs_env(qs_env=qs_env, matrix_s=matrix_s, &
396 particle_set=particle_set, qs_kind_set=qs_kind_set, molecule_set=mols)
399 max_iter = qs_loc_env%localized_wfn_control%max_iter
400 eps = qs_loc_env%localized_wfn_control%eps_localization
401 out_each = qs_loc_env%localized_wfn_control%out_each
403 natom =
SIZE(particle_set, 1)
404 ALLOCATE (first_sgf(natom))
405 ALLOCATE (last_sgf(natom))
406 ALLOCATE (nsgf(natom))
410 first_sgf=first_sgf, last_sgf=last_sgf, nsgf=nsgf)
413 CALL mpools_get(qs_env%mpools, ao_ao_fm_pools=ao_ao_fm_pools)
422 isgf = first_sgf(iatom)
426 CALL parallel_gemm(
'N',
'N', nao, nmoloc, nao, 1.0_dp, ov_fm, vectors, 0.0_dp, opvec, &
427 a_first_col=1, a_first_row=1, b_first_col=1, b_first_row=1)
429 CALL parallel_gemm(
'T',
'N', nmoloc, nmoloc, ncol, 0.5_dp, vectors, opvec, &
430 0.0_dp, zij_fm_set(iatom, 1), &
431 a_first_col=1, a_first_row=isgf, b_first_col=1, b_first_row=isgf)
433 CALL parallel_gemm(
'N',
'N', nao, nmoloc, ncol, 1.0_dp, ov_fm, vectors, 0.0_dp, opvec, &
434 a_first_col=isgf, a_first_row=1, b_first_col=1, b_first_row=isgf)
436 CALL parallel_gemm(
'T',
'N', nmoloc, nmoloc, nao, 0.5_dp, vectors, opvec, &
437 1.0_dp, zij_fm_set(iatom, 1), &
438 a_first_col=1, a_first_row=1, b_first_col=1, b_first_row=1)
447 cpabort(
"GAPO and Pipek not implemented.")
449 cpabort(
"Crazy and Pipek not implemented.")
451 cpabort(
"L1 norm and Pipek not implemented.")
453 cpabort(
"Direct and Pipek not implemented.")
455 ALLOCATE (my_zij_cfm_set(
SIZE(zij_fm_set, 1),
SIZE(zij_fm_set, 2)))
456 DO idim1 = 1,
SIZE(zij_fm_set, 1)
457 DO idim2 = 1,
SIZE(zij_fm_set, 2)
458 CALL cp_cfm_create(my_zij_cfm_set(idim1, idim2), zij_fm_set(idim1, idim2)%matrix_struct)
459 CALL cp_fm_to_cfm(msourcer=zij_fm_set(idim1, idim2), mtarget=my_zij_cfm_set(idim1, idim2))
462 ALLOCATE (complex_vectors)
464 CALL cp_fm_to_cfm(msourcer=vectors, mtarget=complex_vectors)
470 DEALLOCATE (complex_vectors)
471 DO idim1 = 1,
SIZE(my_zij_cfm_set, 1)
472 DO idim2 = 1,
SIZE(my_zij_cfm_set, 2)
476 DEALLOCATE (my_zij_cfm_set)
479 IF (output_unit > 0)
WRITE (output_unit,
'(A,I6,A)')
" No MOS localization applied "
481 cpabort(
"Unknown localization method")
484 IF (output_unit > 0)
WRITE (output_unit,
'(A,I3,A,I6,A)')
" Localization for spin ", ispin, &
485 " converged in ", sweeps,
" iterations"
487 CALL centers_spreads_pipek(qs_loc_env, zij_fm_set, particle_set, ispin, &
490 DEALLOCATE (first_sgf, last_sgf, nsgf)
497 CALL timestop(handle)
513 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(IN) :: zij_fm_set
517 INTEGER :: iatom, istate, jstate, natom, nstate
518 REAL(kind=
dp) :: aij, bij, ct, mii, mij, mjj, ratio, st, &
526 tolerance = 1.0e10_dp
528 natom =
SIZE(zij_fm_set, 1)
530 DO WHILE (tolerance >= 1.0e-4_dp)
532 DO istate = 1, nstate
533 DO jstate = istate + 1, nstate
540 aij = aij + mij*(mii - mjj)
541 bij = bij + mij*mij - 0.25_dp*(mii - mjj)*(mii - mjj)
543 IF (abs(bij) > 1.e-10_dp)
THEN
545 theta = 0.25_dp*atan(ratio)
552 IF (theta >
pi*0.5_dp)
THEN
553 theta = theta -
pi*0.25_dp
554 ELSE IF (theta < -
pi*0.5_dp)
THEN
555 theta = theta +
pi*0.25_dp
569 CALL check_tolerance_real(zij_fm_set, tolerance)
585 SUBROUTINE check_tolerance_real(zij_fm_set, tolerance)
587 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(IN) :: zij_fm_set
588 REAL(
dp),
INTENT(OUT) :: tolerance
590 INTEGER :: iatom, istate, jstate, natom, &
591 ncol_local, nrow_global, nrow_local
592 INTEGER,
DIMENSION(:),
POINTER :: col_indices, row_indices
593 REAL(
dp) :: grad_ij, zii, zij, zjj
594 REAL(
dp),
DIMENSION(:, :),
POINTER :: diag
597 CALL cp_fm_create(force, zij_fm_set(1, 1)%matrix_struct)
600 NULLIFY (diag, col_indices, row_indices)
601 natom =
SIZE(zij_fm_set, 1)
603 ncol_local=ncol_local, nrow_global=nrow_global, &
604 row_indices=row_indices, col_indices=col_indices)
605 ALLOCATE (diag(nrow_global, natom))
608 DO istate = 1, nrow_global
609 CALL cp_fm_get_element(zij_fm_set(iatom, 1), istate, istate, diag(istate, iatom))
613 DO istate = 1, nrow_local
614 DO jstate = 1, ncol_local
617 zii = diag(row_indices(istate), iatom)
618 zjj = diag(col_indices(jstate), iatom)
619 zij = zij_fm_set(iatom, 1)%local_data(istate, jstate)
620 grad_ij = grad_ij + 4.0_dp*zij*(zjj - zii)
622 force%local_data(istate, jstate) = grad_ij
631 END SUBROUTINE check_tolerance_real
648 print_loc_section, zij, c_zij, only_initial_out)
650 INTEGER,
INTENT(IN) :: nmoloc
652 REAL(
dp),
DIMENSION(:) :: weights
653 INTEGER,
INTENT(IN) :: ispin
655 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(IN), &
658 LOGICAL,
INTENT(IN),
OPTIONAL :: only_initial_out
660 CHARACTER(len=default_path_length) :: file_tmp, iter
661 COMPLEX(KIND=dp) :: z
662 INTEGER :: idir, istate, jdir, nstates, &
663 output_unit, unit_out_s
664 LOGICAL :: my_only_init, zij_is_complex
665 REAL(
dp) :: avg_spread_ii, spread_i, spread_ii, &
666 sum_spread_i, sum_spread_ii
667 REAL(
dp),
DIMENSION(3) :: c, c2, cpbc
668 REAL(
dp),
DIMENSION(:, :),
POINTER :: centers
669 REAL(kind=
dp) :: imagpart, realpart
673 NULLIFY (centers, logger, print_key)
675 my_only_init = .false.
676 IF (
PRESENT(only_initial_out)) my_only_init = only_initial_out
677 IF (
PRESENT(zij))
THEN
678 zij_is_complex = .false.
679 cpassert(.NOT.
PRESENT(c_zij))
682 zij_is_complex = .true.
683 cpassert(
PRESENT(c_zij))
687 file_tmp = trim(qs_loc_env%tag_mo)//
"_spreads_s"//trim(adjustl(
cp_to_string(ispin)))
689 extension=
".locInfo")
691 middle_name=file_tmp, extension=
".data")
693 IF (unit_out_s > 0 .AND. .NOT. my_only_init)
WRITE (unit_out_s,
'(i6,/,A)') nmoloc, trim(iter)
695 cpassert(nstates >= nmoloc)
697 centers => qs_loc_env%localized_wfn_control%centers_set(ispin)%array
698 cpassert(
ASSOCIATED(centers))
699 cpassert(
SIZE(centers, 2) == nmoloc)
700 sum_spread_i = 0.0_dp
701 sum_spread_ii = 0.0_dp
702 avg_spread_ii = 0.0_dp
703 DO istate = 1, nmoloc
708 IF (.NOT. zij_is_complex)
THEN
709 DO jdir = 1,
SIZE(zij, 2)
712 z = cmplx(realpart, imagpart,
dp)
713 spread_i = spread_i - weights(jdir)* &
714 log(realpart*realpart + imagpart*imagpart)/
twopi/
twopi
715 spread_ii = spread_ii + weights(jdir)* &
716 (1.0_dp - (realpart*realpart + imagpart*imagpart))/
twopi/
twopi
719 c(idir) = c(idir) + &
720 (cell%hmat(idir, jdir)/
twopi)*aimag(log(z))
725 DO jdir = 1,
SIZE(c_zij, 2)
727 realpart = real(z,
dp)
729 spread_i = spread_i - weights(jdir)* &
730 log(realpart*realpart + imagpart*imagpart)/
twopi/
twopi
731 spread_ii = spread_ii + weights(jdir)* &
732 (1.0_dp - (realpart*realpart + imagpart*imagpart))/
twopi/
twopi
735 c(idir) = c(idir) + &
736 (cell%hmat(idir, jdir)/
twopi)*aimag(log(z))
742 centers(1:3, istate) = cpbc(1:3)
743 centers(4, istate) = spread_i
744 centers(5, istate) = spread_ii
745 sum_spread_i = sum_spread_i + centers(4, istate)
746 sum_spread_ii = sum_spread_ii + centers(5, istate)
747 IF (unit_out_s > 0 .AND. .NOT. my_only_init)
WRITE (unit_out_s,
'(I6,2F16.8)') istate, centers(4:5, istate)
749 avg_spread_ii = sum_spread_ii/real(nmoloc, kind=
dp)
753 IF (.NOT. my_only_init)
CALL print_wannier_centers(qs_loc_env, print_key, centers, logger, ispin)
755 IF (output_unit > 0)
THEN
756 WRITE (output_unit,
'(T4, A, 2x, 2A26,/,T23, A28)')
" Spread Functional ",
"sum_in -w_i ln(|z_in|^2)", &
757 "sum_in w_i(1-|z_in|^2)",
"sum_in w_i(1-|z_in|^2)/n"
758 IF (my_only_init)
THEN
759 WRITE (output_unit,
'(T4,A,T38,2F20.10,/,T38,F20.10)')
" Initial Spread (Berry) : ", &
760 sum_spread_i, sum_spread_ii, avg_spread_ii
762 WRITE (output_unit,
'(T4,A,T38,2F20.10,/,T38,F20.10)')
" Total Spread (Berry) : ", &
763 sum_spread_i, sum_spread_ii, avg_spread_ii
767 IF (unit_out_s > 0 .AND. .NOT. my_only_init)
WRITE (unit_out_s,
'(A,2F16.10)')
" Total ", sum_spread_i, sum_spread_ii
785 SUBROUTINE centers_spreads_pipek(qs_loc_env, zij_fm_set, particle_set, ispin, &
789 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(IN) :: zij_fm_set
791 INTEGER,
INTENT(IN) :: ispin
794 CHARACTER(len=default_path_length) :: file_tmp, iter
795 INTEGER :: iatom, istate, natom, nstate, unit_out_s
796 INTEGER,
DIMENSION(:),
POINTER :: atom_of_state
798 REAL(
dp),
ALLOCATABLE,
DIMENSION(:) :: qii, ziimax
799 REAL(
dp),
ALLOCATABLE,
DIMENSION(:, :) :: diag
800 REAL(
dp),
DIMENSION(:, :),
POINTER :: centers
804 NULLIFY (centers, logger, print_key)
808 natom =
SIZE(zij_fm_set, 1)
810 centers => qs_loc_env%localized_wfn_control%centers_set(ispin)%array
811 cpassert(
ASSOCIATED(centers))
812 cpassert(
SIZE(centers, 2) == nstate)
814 file_tmp = trim(qs_loc_env%tag_mo)//
"_spreads_s"//trim(adjustl(
cp_to_string(ispin)))
816 middle_name=file_tmp, extension=
".data", log_filename=.false.)
818 IF (unit_out_s > 0)
WRITE (unit_out_s,
'(i6,/,A)') trim(iter)
820 ALLOCATE (atom_of_state(nstate))
823 ALLOCATE (diag(nstate, natom))
827 DO istate = 1, nstate
828 CALL cp_fm_get_element(zij_fm_set(iatom, 1), istate, istate, diag(istate, iatom))
832 ALLOCATE (qii(nstate), ziimax(nstate))
837 DO istate = 1, nstate
838 qii(istate) = qii(istate) + diag(istate, iatom)*diag(istate, iatom)
839 IF (abs(diag(istate, iatom)) > ziimax(istate))
THEN
840 ziimax(istate) = abs(diag(istate, iatom))
841 atom_of_state(istate) = iatom
846 DO istate = 1, nstate
847 iatom = atom_of_state(istate)
848 r(1:3) = particle_set(iatom)%r(1:3)
849 centers(1:3, istate) = r(1:3)
850 centers(4, istate) = 1.0_dp/qii(istate)
851 IF (unit_out_s > 0)
WRITE (unit_out_s,
'(I6,F16.8)') istate,
angstrom*centers(4, istate)
856 CALL print_wannier_centers(qs_loc_env, print_key, centers, logger, ispin)
860 DEALLOCATE (qii, ziimax, atom_of_state, diag)
862 END SUBROUTINE centers_spreads_pipek
884 print_key, root, ispin, idir, state0, file_position)
888 INTEGER,
INTENT(IN) :: nstates
889 INTEGER,
DIMENSION(:),
POINTER :: state_list
890 REAL(
dp),
DIMENSION(:, :),
POINTER :: centers
892 CHARACTER(LEN=*) :: root
893 INTEGER,
INTENT(IN),
OPTIONAL :: ispin, idir
894 INTEGER,
OPTIONAL :: state0
895 CHARACTER(LEN=default_string_length),
OPTIONAL :: file_position
897 CHARACTER(len=*),
PARAMETER :: routinen =
'qs_print_cubes'
898 CHARACTER,
DIMENSION(3),
PARAMETER :: labels = [
'x',
'y',
'z']
901 CHARACTER(LEN=default_path_length) :: file_tmp, filename, my_pos
902 CHARACTER(LEN=default_string_length) :: title
903 INTEGER :: handle, ia, ie, istate, ivector, &
904 my_ispin, my_state0, unit_out_c
905 LOGICAL :: add_idir, add_spin, mpi_io
916 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
919 CALL timeset(routinen, handle)
922 CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, subsys=subsys)
925 CALL pw_env_get(pw_env, auxbas_pw_pool=auxbas_pw_pool)
927 CALL auxbas_pw_pool%create_pw(wf_r)
928 CALL auxbas_pw_pool%create_pw(wf_g)
931 IF (
PRESENT(state0)) my_state0 = state0
935 IF (
PRESENT(ispin))
THEN
940 IF (
PRESENT(idir))
THEN
946 IF (
PRESENT(file_position))
THEN
947 my_pos = file_position
950 DO istate = 1, nstates
951 ivector = state_list(istate) - my_state0
952 CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set, qs_kind_set=qs_kind_set, cell=cell, &
953 dft_control=dft_control, particle_set=particle_set, pw_env=pw_env)
956 qs_kind_set, cell, dft_control, particle_set, pw_env)
959 ivector = state_list(istate)
962 filename = root(ia:ie)//
"_"//label//
"_w"//trim(adjustl(
cp_to_string(ivector)))
964 filename = root(ia:ie)//
"_w"//trim(adjustl(
cp_to_string(ivector)))
969 filename = file_tmp(ia:ie)//
"_s"//trim(adjustl(
cp_to_string(ispin)))
975 extension=
".cube", file_position=my_pos, log_filename=.false., &
977 IF (
SIZE(centers, 1) == 6)
THEN
978 WRITE (title,
'(A7,I5.5,A2,I1.1,A1,6F10.4)')
"WFN ", ivector,
"_s", my_ispin,
" ", &
981 WRITE (title,
'(A7,I5.5,A2,I1.1,A1,3F10.4)')
"WFN ", ivector,
"_s", my_ispin,
" ", &
985 particles=particles, &
990 CALL auxbas_pw_pool%give_back_pw(wf_r)
991 CALL auxbas_pw_pool%give_back_pw(wf_g)
992 CALL timestop(handle)
1009 SUBROUTINE qs_print_density_cubes(qs_env, nstates, state_list, centers, print_key, root, &
1010 mo_coeff, complex_mo_coeff, ispin, file_position)
1013 INTEGER,
INTENT(IN) :: nstates
1014 INTEGER,
DIMENSION(:),
POINTER :: state_list
1015 REAL(
dp),
DIMENSION(:, :),
POINTER :: centers
1017 CHARACTER(LEN=*) :: root
1018 TYPE(
cp_fm_type),
OPTIONAL,
POINTER :: mo_coeff
1019 TYPE(
cp_cfm_type),
OPTIONAL,
POINTER :: complex_mo_coeff
1020 INTEGER,
INTENT(IN),
OPTIONAL :: ispin
1021 CHARACTER(LEN=default_string_length),
OPTIONAL :: file_position
1025 IF (
PRESENT(mo_coeff))
THEN
1026 cpassert(.NOT.
PRESENT(complex_mo_coeff))
1028 squaredmat%local_data = mo_coeff%local_data**2.0
1030 cpassert(
PRESENT(complex_mo_coeff))
1031 CALL cp_fm_create(squaredmat, complex_mo_coeff%matrix_struct)
1032 squaredmat%local_data = real(conjg(complex_mo_coeff%local_data) &
1033 *complex_mo_coeff%local_data,
dp)
1035 CALL qs_print_cubes(qs_env, squaredmat, nstates, state_list, centers, print_key, root, &
1036 ispin=ispin, file_position=file_position)
1039 END SUBROUTINE qs_print_density_cubes
1049 SUBROUTINE print_wannier_centers(qs_loc_env, print_key, center, logger, ispin)
1052 REAL(kind=
dp),
INTENT(IN) :: center(:, :)
1054 INTEGER,
INTENT(IN) :: ispin
1056 CHARACTER(default_string_length) :: iter, unit_str
1057 CHARACTER(LEN=default_string_length) :: my_ext, my_form
1058 INTEGER :: iunit, l, nstates
1059 LOGICAL :: first_time, init_traj
1060 REAL(kind=
dp) :: unit_conv
1062 nstates =
SIZE(center, 2)
1063 my_form =
"formatted"
1072 IF (first_time .OR. (.NOT. qs_loc_env%first_time))
THEN
1074 middle_name=trim(qs_loc_env%tag_mo)//
"_centers_s"//trim(adjustl(
cp_to_string(ispin))), &
1075 log_filename=.false., on_file=.true., is_new_file=init_traj)
1083 CALL print_wannier_traj(qs_loc_env, print_key, center, iunit, init_traj, unit_conv)
1087 WRITE (iunit,
'(i6,/,a)') nstates, trim(iter)
1089 WRITE (iunit,
'(A,4F16.8)')
"X", unit_conv*center(1:4, l)
1096 END SUBROUTINE print_wannier_centers
1107 SUBROUTINE print_wannier_traj(qs_loc_env, print_key, center, iunit, init_traj, unit_conv)
1110 REAL(kind=
dp),
INTENT(IN) :: center(:, :)
1111 INTEGER,
INTENT(IN) :: iunit
1112 LOGICAL,
INTENT(IN) :: init_traj
1113 REAL(kind=
dp),
INTENT(IN) :: unit_conv
1115 CHARACTER(len=default_string_length) :: iter, remark1, remark2, title
1116 INTEGER :: i, iskip, natom, ntot, outformat
1122 NULLIFY (particle_set, atomic_kind_set, atomic_kind, logger)
1124 natom =
SIZE(qs_loc_env%particle_set)
1125 ntot = natom +
SIZE(center, 2)
1127 ALLOCATE (atomic_kind_set(1))
1128 atomic_kind => atomic_kind_set(1)
1130 name=
"X", element_symbol=
"X", mass=0.0_dp)
1133 particle_set(i)%atomic_kind => qs_loc_env%particle_set(i)%atomic_kind
1134 particle_set(i)%r =
pbc(qs_loc_env%particle_set(i)%r, qs_loc_env%cell)
1137 DO i = natom + 1, ntot
1138 particle_set(i)%atomic_kind => atomic_kind
1139 particle_set(i)%r =
pbc(center(1:3, i - natom), qs_loc_env%cell)
1145 SELECT CASE (outformat)
1156 WRITE (iunit)
"CORD", 0, -1, iskip, &
1157 0, 0, 0, 0, 0, 0, real(0, kind=
sp), 1, 0, 0, 0, 0, 0, 0, 0, 0, 24
1158 remark1 =
"REMARK FILETYPE CORD DCD GENERATED BY CP2K"
1159 remark2 =
"REMARK Support new DCD format with cell information"
1160 WRITE (iunit) 2, remark1, remark2
1166 WRITE (unit=title, fmt=
"(A)")
" Particles+Wannier centers. Iteration:"//trim(iter)
1171 unit_conv=unit_conv)
1175 END SUBROUTINE print_wannier_traj
1193 INTEGER,
INTENT(IN) :: ispin
1195 INTEGER,
PARAMETER :: norder = 2
1196 LOGICAL,
PARAMETER :: debug_this_subroutine = .false.
1198 CHARACTER(len=default_path_length) :: file_tmp
1199 INTEGER :: imoment, istate, ncol_global, nm, &
1200 nmoloc, nrow_global, output_unit, &
1202 REAL(
dp),
DIMENSION(:, :),
POINTER :: centers
1203 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: moment_set
1204 REAL(kind=
dp),
DIMENSION(3) :: rcc
1205 REAL(kind=
dp),
DIMENSION(6) :: cov
1208 TYPE(
cp_fm_type),
DIMENSION(:),
POINTER :: moloc_coeff
1211 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_s, moments
1217 extension=
".locInfo")
1219 IF (output_unit > 0)
THEN
1220 WRITE (output_unit,
'(/,T2,A)') &
1221 '!-----------------------------------------------------------------------------!'
1222 WRITE (output_unit,
'(T12,A)')
"Computing second moments of Wannier functions..."
1223 WRITE (output_unit,
'(T2,A)') &
1224 '!-----------------------------------------------------------------------------!'
1227 file_tmp =
"_r_loc_cov_matrix_"//trim(adjustl(
cp_to_string(ispin)))
1229 middle_name=file_tmp, extension=
".data")
1231 CALL get_qs_loc_env(qs_loc_env=qs_loc_env, moloc_coeff=moloc_coeff)
1232 centers => qs_loc_env%localized_wfn_control%centers_set(ispin)%array
1240 para_env => qs_loc_env%para_env
1243 ALLOCATE (moment_set(nm, nmoloc))
1246 mo_localized => moloc_coeff(ispin)
1248 DO istate = 1, nmoloc
1249 rcc = centers(1:3, istate)
1252 ALLOCATE (moments(nm))
1254 ALLOCATE (moments(imoment)%matrix)
1255 CALL dbcsr_copy(moments(imoment)%matrix, matrix_s(ispin)%matrix,
'MOM MAT')
1256 CALL dbcsr_set(moments(imoment)%matrix, 0.0_dp)
1262 CALL cp_fm_get_info(mo_localized, ncol_global=ncol_global, nrow_global=nrow_global)
1264 para_env=mo_localized%matrix_struct%para_env, &
1265 context=mo_localized%matrix_struct%context)
1267 CALL cp_fm_create(omvector, fm_struct, name=
"omvector")
1272 CALL cp_fm_to_fm(mo_localized, mvector, 1, istate, 1)
1278 moment_set(imoment, istate) = sum(momv%local_data)
1288 DEALLOCATE (moments)
1291 CALL para_env%sum(moment_set)
1293 IF (output_unit_sm > 0)
THEN
1294 WRITE (unit=output_unit_sm, fmt=
'(A,T13,A,I1)')
"#", &
1295 "SECOND MOMENTS OF WANNIER CENTERS FOR SPIN ", ispin
1296 WRITE (unit=output_unit_sm, fmt=
'(A,T29,A2,5(14x,A2))')
"#",
"XX",
"XY",
"XZ",
"YY",
"YZ",
"ZZ"
1297 DO istate = 1, nmoloc
1299 cov(1) = moment_set(4, istate) - moment_set(1, istate)*moment_set(1, istate)
1300 cov(2) = moment_set(5, istate) - moment_set(1, istate)*moment_set(2, istate)
1301 cov(3) = moment_set(6, istate) - moment_set(1, istate)*moment_set(3, istate)
1302 cov(4) = moment_set(7, istate) - moment_set(2, istate)*moment_set(2, istate)
1303 cov(5) = moment_set(8, istate) - moment_set(2, istate)*moment_set(3, istate)
1304 cov(6) = moment_set(9, istate) - moment_set(3, istate)*moment_set(3, istate)
1305 WRITE (unit=output_unit_sm, fmt=
'(T4,A,I6,6(T20,6F16.8))')
"Center:", istate, cov(1:6)
1306 IF (debug_this_subroutine)
THEN
1307 WRITE (unit=output_unit_sm, fmt=
'(T4,A,I5,9(T20,9F12.6))')
"Center:", istate, moment_set(1:, istate)
1314 DEALLOCATE (moment_set)
1336 INTEGER,
INTENT(IN) :: ispin
1338 INTEGER,
PARAMETER :: taylor_order = 1
1339 LOGICAL,
PARAMETER :: debug_this_subroutine = .false.
1341 CHARACTER(len=default_path_length) :: file_tmp
1342 COMPLEX(KIND=dp) :: z
1343 INTEGER :: icov, imoment, istate, ncol_global, nm, &
1344 nmoloc, norder, nrow_global, &
1345 output_unit, output_unit_sm
1346 REAL(
dp),
DIMENSION(:, :),
POINTER :: centers
1347 REAL(kind=
dp) :: imagpart, lx, ly, lz, realpart
1348 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: moment_set
1349 REAL(kind=
dp),
DIMENSION(3) :: rcc
1350 REAL(kind=
dp),
DIMENSION(6) :: cov
1354 TYPE(
cp_fm_type),
DIMENSION(:),
POINTER :: moloc_coeff
1357 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_s, moments
1365 extension=
".locInfo")
1367 IF (output_unit > 0)
THEN
1368 WRITE (output_unit,
'(/,T2,A)') &
1369 '!-----------------------------------------------------------------------------!'
1370 WRITE (output_unit,
'(T12,A)')
"Computing second moments of Wannier functions..."
1371 WRITE (output_unit,
'(T2,A)') &
1372 '!-----------------------------------------------------------------------------!'
1375 file_tmp =
"_r_berry_cov_matrix_"//trim(adjustl(
cp_to_string(ispin)))
1377 middle_name=file_tmp, extension=
".data")
1379 norder = 2*(2*taylor_order + 1)
1380 nm = (6 + 11*norder + 6*norder**2 + norder**3)/6 - 1
1383 CALL get_qs_loc_env(qs_loc_env=qs_loc_env, moloc_coeff=moloc_coeff, cell=cell)
1384 lx = cell%hmat(1, 1)
1385 ly = cell%hmat(2, 2)
1386 lz = cell%hmat(3, 3)
1388 centers => qs_loc_env%localized_wfn_control%centers_set(ispin)%array
1390 para_env => qs_loc_env%para_env
1396 ALLOCATE (moment_set(nm, nmoloc))
1399 mo_localized => moloc_coeff(ispin)
1401 DO istate = 1, nmoloc
1402 rcc = centers(1:3, istate)
1405 ALLOCATE (moments(nm))
1407 ALLOCATE (moments(imoment)%matrix)
1408 CALL dbcsr_copy(moments(imoment)%matrix, matrix_s(ispin)%matrix,
'MOM MAT')
1409 CALL dbcsr_set(moments(imoment)%matrix, 0.0_dp)
1415 CALL cp_fm_get_info(mo_localized, ncol_global=ncol_global, nrow_global=nrow_global)
1417 para_env=mo_localized%matrix_struct%para_env, &
1418 context=mo_localized%matrix_struct%context)
1420 CALL cp_fm_create(omvector, fm_struct, name=
"omvector")
1425 CALL cp_fm_to_fm(mo_localized, mvector, 1, istate, 1)
1430 moment_set(imoment, istate) = sum(momv%local_data)
1440 DEALLOCATE (moments)
1443 CALL para_env%sum(moment_set)
1445 IF (output_unit_sm > 0)
THEN
1446 WRITE (unit=output_unit_sm, fmt=
'(A,T13,A,I1)')
"#", &
1447 "SECOND MOMENTS OF WANNIER CENTERS FOR SPIN ", ispin
1448 WRITE (unit=output_unit_sm, fmt=
'(A,T29,A2,5(14x,A2))')
"#",
"XX",
"XY",
"XZ",
"YY",
"YZ",
"ZZ"
1449 DO istate = 1, nmoloc
1454 z = cmplx(realpart, imagpart,
dp)
1459 realpart = 1.0 - 0.5*
twopi*
twopi/lx/lx/lx/lx &
1460 *moment_set(20, istate)
1461 imagpart =
twopi/lx/lx*moment_set(4, istate) &
1469 z = cmplx(realpart, imagpart,
dp)
1470 cov(1) = lx*lx/
twopi*aimag(log(z)) - lx*lx/
twopi/
twopi*moment_set(1, istate)*moment_set(1, istate)
1474 realpart = 1.0 - 0.5*
twopi*
twopi/lx/ly/lx/ly &
1475 *moment_set(23, istate)
1476 imagpart =
twopi/lx/ly*moment_set(5, istate) &
1484 z = cmplx(realpart, imagpart,
dp)
1485 cov(2) = lx*ly/
twopi*aimag(log(z)) - lx*ly/
twopi/
twopi*moment_set(1, istate)*moment_set(2, istate)
1489 realpart = 1.0 - 0.5*
twopi*
twopi/lx/lz/lx/lz &
1490 *moment_set(25, istate)
1491 imagpart =
twopi/lx/lz*moment_set(6, istate) &
1499 z = cmplx(realpart, imagpart,
dp)
1500 cov(3) = lx*lz/
twopi*aimag(log(z)) - lx*lz/
twopi/
twopi*moment_set(1, istate)*moment_set(3, istate)
1504 realpart = 1.0 - 0.5*
twopi*
twopi/ly/ly/ly/ly &
1505 *moment_set(30, istate)
1506 imagpart =
twopi/ly/ly*moment_set(7, istate) &
1514 z = cmplx(realpart, imagpart,
dp)
1515 cov(4) = ly*ly/
twopi*aimag(log(z)) - ly*ly/
twopi/
twopi*moment_set(2, istate)*moment_set(2, istate)
1519 realpart = 1.0 - 0.5*
twopi*
twopi/ly/lz/ly/lz &
1520 *moment_set(32, istate)
1521 imagpart =
twopi/ly/lz*moment_set(8, istate) &
1529 z = cmplx(realpart, imagpart,
dp)
1530 cov(5) = ly*lz/
twopi*aimag(log(z)) - ly*lz/
twopi/
twopi*moment_set(2, istate)*moment_set(3, istate)
1534 realpart = 1.0 - 0.5*
twopi*
twopi/lz/lz/lz/lz &
1535 *moment_set(34, istate)
1536 imagpart =
twopi/lz/lz*moment_set(9, istate) &
1544 z = cmplx(realpart, imagpart,
dp)
1545 cov(6) = lz*lz/
twopi*aimag(log(z)) - lz*lz/
twopi/
twopi*moment_set(3, istate)*moment_set(3, istate)
1549 WRITE (unit=output_unit_sm, fmt=
'(T4,A,I6,6(T20,6F16.8))')
"Center:", istate, cov(1:6)
1550 IF (debug_this_subroutine)
THEN
1551 WRITE (unit=output_unit_sm, fmt=
'(T4,A,I5,9(T20,9F12.6))')
"Center:", istate, moment_set(1:, istate)
1558 DEALLOCATE (moment_set)
Define the atomic kind types and their sub types.
subroutine, public set_atomic_kind(atomic_kind, element_symbol, name, mass, kind_number, natom, atom_list, fist_potential, shell, shell_active, damping)
Set the components of an atomic kind data set.
subroutine, public deallocate_atomic_kind_set(atomic_kind_set)
Destructor routine for a set of atomic kinds.
Handles all functions related to the CELL.
Represents a complex full matrix distributed on many processors.
subroutine, public cp_cfm_get_element(matrix, irow_global, icol_global, alpha)
Get the matrix element by its global index.
subroutine, public cp_cfm_release(matrix)
Releases a full matrix.
subroutine, public cp_fm_to_cfm(msourcer, msourcei, mtarget)
Construct a complex full matrix by taking its real and imaginary parts from two separate real-value f...
subroutine, public cp_cfm_create(matrix, matrix_struct, name, nrow, ncol, set_zero)
Creates a new full matrix with the given structure.
subroutine, public cp_cfm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, matrix_struct, para_env)
Returns information about a full matrix.
subroutine, public cp_cfm_to_fm(msource, mtargetr, mtargeti)
Copy real and imaginary parts of a complex full matrix into separate real-value full matrices.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_deallocate_matrix(matrix)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_set(matrix, alpha)
...
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_sm_fm_multiply(matrix, fm_in, fm_out, ncol, alpha, beta)
multiply a dbcsr with a fm matrix
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_rot_rows(matrix, irow, jrow, cs, sn)
Applies a planar rotation defined by cs and sn to the i'th and j'th rows.
subroutine, public cp_fm_rot_cols(matrix, icol, jcol, cs, sn)
Applies a planar rotation defined by cs and sn to the i'th and j'th columnns.
subroutine, public cp_fm_schur_product(matrix_a, matrix_b, matrix_c)
computes the schur product of two matrices c_ij = a_ij * b_ij
pool for for elements that are retained and released
subroutine, public fm_pool_create_fm(pool, element, name)
returns an element, allocating it if none is in the pool
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
subroutine, public cp_fm_get_element(matrix, irow_global, icol_global, alpha, local)
returns an element of a fm this value is valid on every cpu using this call is expensive
subroutine, public cp_fm_maxabsval(matrix, a_max, ir_max, ic_max)
find the maximum absolute value of the matrix element maxval(abs(matrix))
subroutine, public cp_fm_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
character(len=default_string_length) function, public cp_iter_string(iter_info, print_key, for_file)
returns the iteration string, a string that is useful to create unique filenames (once you trim it)
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
integer, parameter, public cp_p_file
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
A wrapper around pw_to_cube() which accepts particle_list_type.
subroutine, public cp_pw_to_cube(pw, unit_nr, title, particles, zeff, stride, max_file_size_mb, zero_tails, silent, mpi_io)
...
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
integer, parameter, public default_path_length
integer, parameter, public sp
Machine interface based on Fortran 2003 and POSIX.
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
real(kind=dp), parameter, public twopi
Interface to the message passing library MPI.
Define the data structure for the molecule information.
Output Utilities for MOTION_SECTION.
subroutine, public get_output_format(section, path, my_form, my_ext)
Info on the unit to be opened to dump MD informations.
Provides Cartesian and spherical orbital pointers and indices.
integer, dimension(:), allocatable, public ncoset
basic linear algebra operations for full matrixes
represent a simple array based list of the given type
Define methods related to particle_type.
subroutine, public get_particle_set(particle_set, qs_kind_set, first_sgf, last_sgf, nsgf, nmao, basis, ncgf)
Get the components of a particle set.
subroutine, public write_particle_coordinates(particle_set, iunit, output_format, content, title, cell, array, unit_conv, charge_occup, charge_beta, charge_extended, print_kind)
Should be able to write a few formats e.g. xmol, and some binary format (dcd) some format can be used...
Define the data structure for the particle information.
subroutine, public deallocate_particle_set(particle_set)
Deallocate a particle set.
subroutine, public allocate_particle_set(particle_set, nparticle)
Allocate a particle set.
Definition of physical constants:
real(kind=dp), parameter, public angstrom
container for various plainwaves related things
subroutine, public pw_env_get(pw_env, pw_pools, cube_info, gridlevel_info, auxbas_pw_pool, auxbas_grid, auxbas_rs_desc, auxbas_rs_grid, rs_descs, rs_grids, xc_pw_pool, vdw_pw_pool, poisson_env, interp_section)
returns the various attributes of the pw env
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Calculate the plane wave density by collocating the primitive Gaussian functions (pgf).
subroutine, public calculate_wavefunction(mo_vectors, ivector, rho, rho_gspace, atomic_kind_set, qs_kind_set, cell, dft_control, particle_set, pw_env, basis_type)
maps a given wavefunction on the grid
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, mimic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
Driver for the localization that should be general for all the methods available and all the definiti...
subroutine, public optimize_loc_pipek(qs_env, method, qs_loc_env, vectors, zij_fm_set, ispin, print_loc_section)
...
subroutine, public centers_second_moments_berry(qs_env, qs_loc_env, print_loc_section, ispin)
Compute the second moments of the centers using a periodic quadrupole operator.
subroutine, public centers_second_moments_loc(qs_env, qs_loc_env, print_loc_section, ispin)
Compute the second moments of the centers using the local (non-periodic) pos operators.
subroutine, public optimize_loc_berry(method, qs_loc_env, vectors, op_sm_set, zij_fm_set, para_env, cell, weights, ispin, print_loc_section, restricted, nextra, nmo, vectors_2, guess_mos)
Calculate and optimize the spread functional as calculated from the selected mos and by the definitio...
subroutine, public qs_print_cubes(qs_env, mo_coeff, nstates, state_list, centers, print_key, root, ispin, idir, state0, file_position)
write the cube files for a set of selected states
subroutine, public centers_spreads_berry(qs_loc_env, nmoloc, cell, weights, ispin, print_loc_section, zij, c_zij, only_initial_out)
...
subroutine, public jacobi_rotation_pipek(zij_fm_set, vectors, sweeps)
2by2 rotation for the pipek operator in this case the z_ij numbers are reals
New version of the module for the localization of the molecular orbitals This should be able to use d...
subroutine, public get_qs_loc_env(qs_loc_env, cell, local_molecules, localized_wfn_control, moloc_coeff, op_sm_set, op_fm_set, para_env, particle_set, weights, dim_op)
...
Localization methods such as 2x2 Jacobi rotations Steepest Decents Conjugate Gradient.
subroutine, public approx_l1_norm_sd(c, iterations, eps, converged, sweeps)
...
subroutine, public rotate_orbitals(rmat, vectors)
...
subroutine, public cardoso_souloumiac(weights, zij, max_iter, eps_localization, sweeps, out_each, vectors)
Achieves minimisation of the spread functional by simultaneous diagonalisation with Jacobi rotations ...
subroutine, public jacobi_rotations(weights, zij, vectors, para_env, max_iter, eps_localization, sweeps, out_each, target_time, start_time, restricted)
wrapper for the jacobi routines, should be removed if jacobi_rot_para can deal with serial para_envs.
subroutine, public zij_matrix(vectors, op_sm_set, zij_fm_set)
...
subroutine, public direct_mini(weights, zij, vectors, max_iter, eps_localization, iterations)
use the exponential parametrization as described in to perform a direct mini Gerd Berghold et al....
subroutine, public crazy_rotations(weights, zij, vectors, max_iter, max_crazy_angle, crazy_scale, crazy_use_diag, eps_localization, iterations, converged)
yet another crazy try, computes the angles needed to rotate the orbitals first and rotates them all a...
subroutine, public scdm_qrfact(vectors)
...
subroutine, public cardoso_souloumiac_pipek(zij, vec, sweeps, max_iter, eps, out_each)
Pipek-Mezey version of the Cardoso-Souloumiac PADE algorithm for complex-valued matrices.
subroutine, public jacobi_cg_edf_ls(para_env, weights, zij, vectors, max_iter, eps_localization, iter, out_each, nextra, do_cg, nmo, vectors_2, mos_guess)
combine jacobi rotations (serial) and conjugate gradient with golden section line search for partiall...
wrapper for the pools of matrixes
subroutine, public mpools_get(mpools, ao_mo_fm_pools, ao_ao_fm_pools, mo_mo_fm_pools, ao_mosub_fm_pools, mosub_mosub_fm_pools, maxao_maxmo_fm_pool, maxao_maxao_fm_pool, maxmo_maxmo_fm_pool)
returns various attributes of the mpools (notably the pools contained in it)
Calculates the moment integrals <a|r^m|b> and <a|r x d/dr|b>
subroutine, public build_local_moment_matrix(qs_env, moments, nmoments, ref_point, ref_points, basis_type)
...
types that represent a quickstep subsys
subroutine, public qs_subsys_get(subsys, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell, cell_ref, use_ref_cell, energy, force, qs_kind_set, cp_subsys, nelectron_total, nelectron_spin)
...
Utilities for string manipulations.
elemental subroutine, public xstring(string, ia, ib)
...
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Just to build arrays of pointers to matrices.
Represent a complex full matrix.
to create arrays of pools
keeps the information about the structure of a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
stores all the informations relevant to an mpi environment
represent a list of objects
contained for different pw related things
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Provides all information about a quickstep kind.
contains all the info needed by quickstep to calculate the spread of a selected set of orbitals and i...