22 dbcsr_type, dbcsr_type_antisymmetric, dbcsr_type_no_symmetry, dbcsr_type_symmetric
81#include "./base/base_uses.f90"
87 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_tddfpt2_stda_utils'
110 CHARACTER(len=*),
PARAMETER :: routinen =
'stda_init_matrices'
113 LOGICAL :: do_coulomb
114 TYPE(
cell_type),
POINTER :: cell, cell_ref
120 CALL timeset(routinen, handle)
122 do_coulomb = .NOT. tddfpt_control%rks_triplets
125 CALL setup_gamma(qs_env, stda_kernel, sub_env, work%gamma_exchange)
132 IF (tddfpt_control%stda_control%do_ewald)
THEN
133 NULLIFY (ewald_env, ewald_pw)
137 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
140 CALL get_qs_env(qs_env, cell=cell, cell_ref=cell_ref)
142 cell_periodic=cell%perd)
144 CALL ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section=print_section)
145 work%ewald_env => ewald_env
146 work%ewald_pw => ewald_pw
149 CALL timestop(handle)
161 SUBROUTINE setup_gamma(qs_env, stda_env, sub_env, gamma_matrix, ndim)
166 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: gamma_matrix
167 INTEGER,
INTENT(IN),
OPTIONAL :: ndim
169 CHARACTER(len=*),
PARAMETER :: routinen =
'setup_gamma'
170 REAL(kind=
dp),
PARAMETER :: rsmooth = 1.0_dp
172 INTEGER :: handle, i, iatom, icol, ikind, imat, &
173 irow, jatom, jkind, natom, nmat
174 INTEGER,
DIMENSION(:),
POINTER :: row_blk_sizes
176 REAL(kind=
dp) :: dfcut, dgb, dr, eta, fcut, r, rcut, &
178 REAL(kind=
dp),
DIMENSION(3) :: rij
179 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: dgblock, gblock
182 DIMENSION(:),
POINTER :: nl_iterator
186 CALL timeset(routinen, handle)
189 dbcsr_dist => sub_env%dbcsr_dist
192 n_list => sub_env%sab_orb
194 IF (
PRESENT(ndim))
THEN
199 cpassert(nmat == 1 .OR. nmat == 4)
200 cpassert(.NOT.
ASSOCIATED(gamma_matrix))
203 ALLOCATE (row_blk_sizes(natom))
204 row_blk_sizes(1:natom) = 1
206 ALLOCATE (gamma_matrix(imat)%matrix)
209 CALL dbcsr_create(gamma_matrix(1)%matrix, name=
"gamma", dist=dbcsr_dist, &
210 matrix_type=dbcsr_type_symmetric, row_blk_size=row_blk_sizes, &
211 col_blk_size=row_blk_sizes)
213 CALL dbcsr_create(gamma_matrix(imat)%matrix, name=
"dgamma", dist=dbcsr_dist, &
214 matrix_type=dbcsr_type_antisymmetric, row_blk_size=row_blk_sizes, &
215 col_blk_size=row_blk_sizes)
218 DEALLOCATE (row_blk_sizes)
223 CALL dbcsr_set(gamma_matrix(imat)%matrix, 0.0_dp)
226 NULLIFY (nl_iterator)
230 iatom=iatom, jatom=jatom, r=rij)
232 dr = sqrt(sum(rij(:)**2))
234 eta = (stda_env%kind_param_set(ikind)%kind_param%hardness_param + &
235 stda_env%kind_param_set(jkind)%kind_param%hardness_param)/2.0_dp
237 icol = max(iatom, jatom)
238 irow = min(iatom, jatom)
242 row=irow, col=icol, block=gblock, found=found)
246 rcuta = stda_env%kind_param_set(ikind)%kind_param%rcut
247 rcutb = stda_env%kind_param_set(jkind)%kind_param%rcut
254 gblock(:, :) = gblock(:, :) + eta
255 ELSEIF (dr > rcut)
THEN
258 IF (dr < rcut - rsmooth)
THEN
261 r = dr - (rcut - rsmooth)
263 fcut = -6._dp*x**5 + 15._dp*x**4 - 10._dp*x**3 + 1._dp
265 gblock(:, :) = gblock(:, :) + &
266 fcut*(1._dp/(dr**(stda_env%alpha_param) + eta**(-stda_env%alpha_param))) &
267 **(1._dp/stda_env%alpha_param) - fcut/dr
274 IF (dr < 1.e-6 .OR. dr > rcut)
THEN
278 IF (dr < rcut - rsmooth)
THEN
282 r = dr - (rcut - rsmooth)
284 fcut = -6._dp*x**5 + 15._dp*x**4 - 10._dp*x**3 + 1._dp
285 dfcut = -30._dp*x**4 + 60._dp*x**3 - 30._dp*x**2
286 dfcut = dfcut/rsmooth
288 dgb = dfcut*(1._dp/(dr**(stda_env%alpha_param) + eta**(-stda_env%alpha_param))) &
289 **(1._dp/stda_env%alpha_param)
290 dgb = dgb - dfcut/dr + fcut/dr**2
291 dgb = dgb - fcut*(1._dp/(dr**(stda_env%alpha_param) + eta**(-stda_env%alpha_param))) &
292 **(1._dp/stda_env%alpha_param + 1._dp)*dr**(stda_env%alpha_param - 1._dp)
297 row=irow, col=icol, block=dgblock, found=found)
301 IF (irow == iatom)
THEN
302 dgblock(:, :) = dgblock(:, :) + dgb*rij(i)/dr
304 dgblock(:, :) = dgblock(:, :) - dgb*rij(i)/dr
319 CALL timestop(handle)
335 CHARACTER(len=*),
PARAMETER :: routinen =
'get_lowdin_mo_coefficients'
337 INTEGER :: handle, i, iounit, ispin, j, &
338 max_iter_lanczos, nactive, ndep, nsgf, &
339 nspins, order_lanczos
341 REAL(kind=
dp) :: eps_lanczos, sij, threshold
342 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: slam
343 REAL(kind=
dp),
CONTIGUOUS,
DIMENSION(:, :), &
344 POINTER :: local_data
348 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: matrixkp_s
351 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
354 CALL timeset(routinen, handle)
362 WRITE (iounit,
"(1X,A)")
"", &
363 "-------------------------------------------------------------------------------", &
364 "- Create Matrix SQRT(S) -", &
365 "-------------------------------------------------------------------------------"
368 IF (sub_env%is_split)
THEN
371 CALL get_qs_env(qs_env=qs_env, matrix_s_kp=matrixkp_s)
372 cpassert(
ASSOCIATED(matrixkp_s))
373 cpwarn_if(
SIZE(matrixkp_s, 2) > 1,
"not implemented for k-points.")
374 sm_s => matrixkp_s(1, 1)%matrix
380 threshold = 1.0e-8_dp
382 eps_lanczos = 1.0e-4_dp
383 max_iter_lanczos = 40
385 threshold, order_lanczos, eps_lanczos, max_iter_lanczos, &
389 NULLIFY (qs_kind_set)
390 CALL get_qs_env(qs_env=qs_env, qs_kind_set=qs_kind_set)
394 IF (.NOT. converged)
THEN
396 WRITE (iounit,
"(T3,A)")
"STDA| Newton-Schulz iteration did not converge"
397 WRITE (iounit,
"(T3,A)")
"STDA| Calculate SQRT(S) from diagonalization"
399 CALL get_qs_env(qs_env=qs_env, scf_control=scf_control)
402 para_env=sub_env%para_env, &
403 context=sub_env%blacs_env, &
406 CALL cp_fm_create(matrix=fm_s_half, matrix_struct=fmstruct, name=
"S^(1/2) MATRIX")
407 CALL cp_fm_create(matrix=fm_work1, matrix_struct=fmstruct, name=
"TMP MATRIX")
410 CALL cp_fm_power(fm_s_half, fm_work1, 0.5_dp, scf_control%eps_eigval, ndep)
412 CALL cp_warn(__location__, &
413 "Overlap matrix exhibits linear dependencies. At least some "// &
414 "eigenvalues have been quenched.")
418 IF (iounit > 0)
WRITE (iounit, *)
421 nspins =
SIZE(sub_env%mos_occ)
424 CALL cp_fm_get_info(work%ctransformed(ispin), ncol_global=nactive)
426 work%ctransformed(ispin), nactive, alpha=1.0_dp, beta=0.0_dp)
430 CALL cp_fm_create(matrix=fm_work1, matrix_struct=work%S_eigenvectors%matrix_struct, name=
"TMP MATRIX")
435 ALLOCATE (slam(nsgf, 1))
437 IF (work%S_eigenvalues(i) > 0._dp)
THEN
438 slam(i, 1) = sqrt(work%S_eigenvalues(i))
440 cpabort(
"S matrix not positive definit")
450 DO i = 1,
SIZE(local_data, 2)
451 DO j = 1,
SIZE(local_data, 1)
452 sij = local_data(j, i)
453 IF (sij > 0.0_dp) sij = 1.0_dp/sij
454 local_data(j, i) = sij
459 CALL timestop(handle)
473 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(IN) :: xvec
474 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(INOUT) :: xt
476 CHARACTER(len=*),
PARAMETER :: routinen =
'get_lowdin_x'
478 INTEGER :: handle, ispin, nactive, nspins
480 CALL timeset(routinen, handle)
488 xt(ispin), nactive, alpha=1.0_dp, beta=0.0_dp)
491 CALL timestop(handle)
509 work, is_rks_triplets, X, res)
516 LOGICAL,
INTENT(IN) :: is_rks_triplets
517 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(IN) :: x
518 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(INOUT) :: res
520 CHARACTER(len=*),
PARAMETER :: routinen =
'stda_calculate_kernel'
522 INTEGER :: ewald_type, handle, ia, iatom, ikind, &
523 is, ispin, jatom, jkind, jspin, natom, &
525 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: first_sgf, kind_of, last_sgf
526 INTEGER,
DIMENSION(2) :: nactive, nlim
527 LOGICAL :: calculate_forces, do_coulomb, do_ewald, &
528 do_exchange, use_virial
529 REAL(kind=
dp) :: alpha, bp, dr, eta, gabr, hfx, rbeta, &
531 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: tcharge, tv
532 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: gtcharge
533 REAL(kind=
dp),
DIMENSION(3) :: rij
534 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: gab, pblock
539 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:) :: xtransformed
550 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
553 CALL timeset(routinen, handle)
555 nactive(:) = stda_env%nactive(:)
558 IF (nspins == 2) spinfac = 1.0_dp
560 IF (nspins == 1 .AND. is_rks_triplets)
THEN
565 do_ewald = stda_control%do_ewald
566 do_exchange = stda_control%do_exchange
568 para_env => sub_env%para_env
570 CALL get_qs_env(qs_env, natom=natom, cell=cell, &
571 particle_set=particle_set, qs_kind_set=qs_kind_set)
572 ALLOCATE (first_sgf(natom))
573 ALLOCATE (last_sgf(natom))
574 CALL get_particle_set(particle_set, qs_kind_set, first_sgf=first_sgf, last_sgf=last_sgf)
578 ALLOCATE (xtransformed(nspins))
581 ct => work%ctransformed(ispin)
583 CALL cp_fm_create(matrix=xtransformed(ispin), matrix_struct=fmstruct, name=
"XTRANSFORMED")
587 ALLOCATE (tcharge(natom), gtcharge(natom, 1))
594 ct => work%ctransformed(ispin)
604 ctjspin => work%ctransformed(jspin)
606 CALL cp_fm_get_info(ctjspin, matrix_struct=fmstructjspin, nrow_global=nsgf)
615 DO is = first_sgf(ia), last_sgf(ia)
616 tcharge(ia) = tcharge(ia) + tv(is)
626 tempmat => work%gamma_exchange(1)%matrix
630 gtcharge(iatom, 1) = gtcharge(iatom, 1) + gab(1, 1)*tcharge(jatom)
631 IF (iatom /= jatom)
THEN
632 gtcharge(jatom, 1) = gtcharge(jatom, 1) + gab(1, 1)*tcharge(iatom)
638 ewald_env => work%ewald_env
639 ewald_pw => work%ewald_pw
640 CALL ewald_env_get(ewald_env, alpha=alpha, ewald_type=ewald_type)
643 calculate_forces = .false.
644 n_list => sub_env%sab_orb
647 gtcharge, tcharge, calculate_forces, virial, use_virial)
649 IF (para_env%is_source())
THEN
650 gtcharge(:, 1) = gtcharge(:, 1) - 2._dp*alpha*
oorootpi*tcharge(:)
653 nlim =
get_limit(natom, para_env%num_pe, para_env%mepos)
654 DO iatom = nlim(1), nlim(2)
655 DO jatom = 1, iatom - 1
656 rij = particle_set(iatom)%r - particle_set(jatom)%r
658 dr = sqrt(sum(rij(:)**2))
659 IF (dr > 1.e-6_dp)
THEN
660 gtcharge(iatom, 1) = gtcharge(iatom, 1) + tcharge(jatom)/dr
661 gtcharge(jatom, 1) = gtcharge(jatom, 1) + tcharge(iatom)/dr
666 CALL para_env%sum(gtcharge)
671 DO is = first_sgf(ia), last_sgf(ia)
672 tv(is) = gtcharge(ia, 1)
676 ct => work%ctransformed(ispin)
685 IF (do_exchange)
THEN
687 CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set)
690 tempmat => work%shalf
691 CALL dbcsr_create(pdens, template=tempmat, matrix_type=dbcsr_type_no_symmetry)
693 ct => work%ctransformed(ispin)
696 1.0_dp, keep_sparsity=.false.)
700 bp = stda_env%beta_param
701 hfx = stda_env%hfx_fraction
705 rij = particle_set(iatom)%r - particle_set(jatom)%r
707 dr = sqrt(sum(rij(:)**2))
708 ikind = kind_of(iatom)
709 jkind = kind_of(jatom)
710 eta = (stda_env%kind_param_set(ikind)%kind_param%hardness_param + &
711 stda_env%kind_param_set(jkind)%kind_param%hardness_param)/2.0_dp
713 IF (hfx > 0.0_dp)
THEN
714 gabr = (1._dp/(rbeta + (hfx*eta)**(-bp)))**(1._dp/bp)
739 DEALLOCATE (tcharge, gtcharge)
740 DEALLOCATE (first_sgf, last_sgf)
742 CALL timestop(handle)
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind_set(atomic_kind_set, atom_of_kind, kind_of, natom_of_kind, maxatom, natom, nshell, fist_potential_present, shell_present, shell_adiabatic, shell_check_distance, damping_present)
Get attributes of an atomic kind set.
Handles all functions related to the CELL.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset, transposed)
...
subroutine, public dbcsr_filter(matrix, eps)
...
subroutine, public dbcsr_finalize(matrix)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_release(matrix)
...
subroutine, public dbcsr_add_on_diag(matrix, alpha)
Adds the given scalar to the diagonal of the matrix. Reserves any missing diagonal blocks.
Routines that link DBCSR and CP2K concepts together.
subroutine, public cp_dbcsr_alloc_block_from_nbl(matrix, sab_orb, desymmetrize)
allocate the blocks of a dbcsr based on the neighbor list
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.
subroutine, public cp_dbcsr_plus_fm_fm_t(sparse_matrix, matrix_v, matrix_g, ncol, alpha, keep_sparsity, symmetry_mode)
performs the multiplication sparse_matrix+dense_mat*dens_mat^T if matrix_g is not explicitly given,...
subroutine, public copy_fm_to_dbcsr(fm, matrix, keep_sparsity)
Copy a BLACS matrix to a dbcsr matrix.
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_row_scale(matrixa, scaling)
scales row i of matrix a with scaling(i)
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
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
subroutine, public cp_fm_power(matrix, work, exponent, threshold, n_dependent, verbose, eigvals)
...
subroutine, public choose_eigv_solver(matrix, eigenvectors, eigenvalues, info)
Choose the Eigensolver depending on which library is available ELPA seems to be unstable for small sy...
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_vectorssum(matrix, sum_array, dir)
summing up all the elements along the matrix's i-th index or
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_set_submatrix(fm, new_values, start_row, start_col, n_rows, n_cols, alpha, beta, transpose)
sets a submatrix of a full matrix fm(start_row:start_row+n_rows,start_col:start_col+n_cols) = alpha*o...
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 ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
subroutine, public ewald_env_set(ewald_env, ewald_type, alpha, epsilon, eps_pol, gmax, ns_max, precs, o_spline, para_env, poisson_section, interaction_cutoffs, cell_hmat)
Purpose: Set the EWALD environment.
subroutine, public ewald_env_create(ewald_env, para_env)
allocates and intitializes a ewald_env
subroutine, public read_ewald_section_tb(ewald_env, ewald_section, hmat, silent, pset, cell_periodic)
Purpose: read the EWALD section for TB methods.
subroutine, public ewald_env_get(ewald_env, ewald_type, alpha, eps_pol, epsilon, gmax, ns_max, o_spline, group, para_env, poisson_section, precs, rcut, do_multipoles, max_multipole, do_ipol, max_ipol_iter, interaction_cutoffs, cell_hmat)
Purpose: Get the EWALD environment.
Calculation of Ewald contributions in DFTB.
subroutine, public tb_ewald_overlap(gmcharge, mcharge, alpha, n_list, virial, use_virial)
...
subroutine, public tb_spme_evaluate(ewald_env, ewald_pw, particle_set, box, gmcharge, mcharge, calculate_forces, virial, use_virial)
...
subroutine, public ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section)
creates the structure ewald_pw_type
Routines useful for iterative matrix calculations.
subroutine, public matrix_sqrt_newton_schulz(matrix_sqrt, matrix_sqrt_inv, matrix, threshold, order, eps_lanczos, max_iter_lanczos, symmetrize, converged, iounit)
compute the sqrt of a matrix via the sign function and the corresponding Newton-Schulz iterations the...
Defines the basic variable types.
integer, parameter, public dp
Definition of mathematical constants and functions.
real(kind=dp), parameter, public oorootpi
Interface to the message passing library MPI.
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.
Define the data structure for the particle information.
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.
subroutine, public get_qs_kind_set(qs_kind_set, all_potential_present, tnadd_potential_present, gth_potential_present, sgp_potential_present, paw_atom_present, dft_plus_u_atom_present, maxcgf, maxsgf, maxco, maxco_proj, maxgtops, maxlgto, maxlprj, maxnset, maxsgf_set, ncgf, npgf, nset, nsgf, nshell, maxpol, maxlppl, maxlppnl, maxppnl, nelectron, maxder, max_ngrid_rad, max_sph_harm, maxg_iso_not0, lmax_rho0, basis_rcut, basis_type, total_zeff_corr, npgf_seg, cneo_potential_present, nkind_q, natom_q)
Get attributes of an atomic kind set.
Define the neighbor list data types and the corresponding functionality.
subroutine, public neighbor_list_iterator_create(iterator_set, nl, search, nthread)
Neighbor list iterator functions.
subroutine, public neighbor_list_iterator_release(iterator_set)
...
integer function, public neighbor_list_iterate(iterator_set, mepos)
...
subroutine, public get_iterator_info(iterator_set, mepos, ikind, jkind, nkind, ilist, nlist, inode, nnode, iatom, jatom, r, cell)
...
Simplified Tamm Dancoff approach (sTDA).
Simplified Tamm Dancoff approach (sTDA).
subroutine, public get_lowdin_mo_coefficients(qs_env, sub_env, work)
Calculate Lowdin MO coefficients.
subroutine, public get_lowdin_x(shalf, xvec, xt)
Calculate Lowdin transformed Davidson trial vector X shalf (dbcsr), xvec, xt (fm) are defined in the ...
subroutine, public setup_gamma(qs_env, stda_env, sub_env, gamma_matrix, ndim)
Calculate sTDA exchange-type contributions.
subroutine, public stda_init_matrices(qs_env, stda_kernel, sub_env, work, tddfpt_control)
Calculate sTDA matrices.
subroutine, public stda_calculate_kernel(qs_env, stda_control, stda_env, sub_env, work, is_rks_triplets, x, res)
...Calculate the sTDA kernel contribution by contracting the Lowdin MO coefficients – transition char...
parameters that control an scf iteration
All kind of helpful little routines.
pure integer function, dimension(2), public get_limit(m, n, me)
divide m entries into n parts, return size of part me
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
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...
to build arrays of pointers
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
Parallel (sub)group environment.
Set of temporary ("work") matrices.