60#include "./base/base_uses.f90"
65 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_tddfpt2_soc_utils'
70 TYPE dbcsr_soc_package_type
71 TYPE(dbcsr_type),
POINTER :: dbcsr_sg => null()
72 TYPE(dbcsr_type),
POINTER :: dbcsr_tp => null()
73 TYPE(dbcsr_type),
POINTER :: dbcsr_sc => null()
74 TYPE(dbcsr_type),
POINTER :: dbcsr_sf => null()
75 TYPE(dbcsr_type),
POINTER :: dbcsr_prod => null()
76 TYPE(dbcsr_type),
POINTER :: dbcsr_ovlp => null()
77 TYPE(dbcsr_type),
POINTER :: dbcsr_tmp => null()
78 TYPE(dbcsr_type),
POINTER :: dbcsr_work => null()
79 END TYPE dbcsr_soc_package_type
97 CHARACTER(len=*),
PARAMETER :: routinen =
'soc_dipole_operator'
99 INTEGER :: dim_op, handle, i_dim, nao, nspin
100 REAL(kind=
dp),
DIMENSION(3) :: reference_point
103 CALL timeset(routinen, handle)
108 cpabort(
"BERRY DIPOLE FORM NOT IMPLEMENTED FOR SOC")
119 ALLOCATE (soc_env%dipmat_ao(dim_op))
121 ALLOCATE (soc_env%dipmat_ao(i_dim)%matrix)
122 CALL dbcsr_copy(soc_env%dipmat_ao(i_dim)%matrix, &
123 matrix_s(1)%matrix, &
124 name=
"dipole operator matrix")
127 SELECT CASE (tddfpt_control%dipole_form)
133 reference=tddfpt_control%dipole_reference, &
134 ref_point=tddfpt_control%dipole_ref_point)
136 CALL rrc_xyz_ao(op=soc_env%dipmat_ao, qs_env=qs_env, rc=reference_point, order=1, &
137 minimum_image=.false., soft=.false.)
139 CALL length_rep(qs_env, gs_mos, soc_env)
143 CALL p_xyz_ao(soc_env%dipmat_ao, qs_env, minimum_image=.false.)
145 CALL velocity_rep(qs_env, gs_mos, soc_env)
147 cpabort(
"Unimplemented form of the dipole operator")
150 CALL timestop(handle)
160 SUBROUTINE length_rep(qs_env, gs_mos, soc_env)
166 INTEGER :: ideriv, ispin, nao, nderivs, nspins
167 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: nmo_virt
170 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:) :: s_mos_virt
171 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:, :) :: dipole_op_mos_occ
172 TYPE(
cp_fm_type),
POINTER :: dipmat_tmp, wfm_ao_ao
177 CALL get_qs_env(qs_env, matrix_s=matrix_s, blacs_env=blacs_env, para_env=para_env)
182 ALLOCATE (s_mos_virt(nspins), dipole_op_mos_occ(3, nspins), &
183 wfm_ao_ao, nmo_virt(nspins), symm_tmp, dipmat_tmp)
185 CALL cp_fm_struct_create(dip_struct, context=blacs_env, ncol_global=nao, nrow_global=nao, para_env=para_env)
189 DO ideriv = 1, nderivs
190 ALLOCATE (soc_env%dipmat(ideriv)%matrix)
191 CALL dbcsr_create(soc_env%dipmat(ideriv)%matrix, template=symm_tmp, &
192 name=
"contracted operator", matrix_type=
"N")
194 CALL cp_fm_create(dipole_op_mos_occ(ideriv, ispin), matrix_struct=dip_struct)
199 DEALLOCATE (symm_tmp)
202 nmo_virt(ispin) =
SIZE(gs_mos(ispin)%evals_virt)
203 CALL cp_fm_get_info(gs_mos(ispin)%mos_virt, matrix_struct=fm_struct)
208 gs_mos(ispin)%mos_virt, &
210 ncol=nmo_virt(ispin), alpha=1.0_dp, beta=0.0_dp)
212 1.0_dp, s_mos_virt(ispin), gs_mos(ispin)%mos_virt, &
215 DO ideriv = 1, nderivs
219 1.0_dp, wfm_ao_ao, dipmat_tmp, &
220 0.0_dp, dipole_op_mos_occ(ideriv, ispin))
221 CALL copy_fm_to_dbcsr(dipole_op_mos_occ(ideriv, ispin), soc_env%dipmat(ideriv)%matrix)
225 DEALLOCATE (wfm_ao_ao)
231 DO ideriv = 1, nderivs
235 DEALLOCATE (s_mos_virt, dipole_op_mos_occ, nmo_virt, dipmat_tmp)
237 END SUBROUTINE length_rep
245 SUBROUTINE velocity_rep(qs_env, gs_mos, soc_env)
251 INTEGER :: ici, icol, ideriv, irow, ispin, n_act, &
252 n_virt, nao, ncols_local, nderivs, &
254 INTEGER,
DIMENSION(:),
POINTER :: col_indices, row_indices
255 REAL(kind=
dp) :: eval_occ
256 REAL(kind=
dp),
CONTIGUOUS,
DIMENSION(:, :), &
257 POINTER :: local_data_ediff
260 fm_struct, scrm_struct
262 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_s, scrm
267 NULLIFY (scrm, scrm_struct, blacs_env, matrix_s, ao_cvirt_struct, cvirt_ao_struct)
270 ALLOCATE (soc_env%SC(nspins), soc_env%CdS(nspins, nderivs), soc_env%ediff(nspins))
272 CALL get_qs_env(qs_env, ks_env=ks_env, sab_orb=sab_orb, blacs_env=blacs_env, matrix_s=matrix_s)
276 CALL cp_fm_get_info(gs_mos(1)%mos_virt, matrix_struct=ao_cvirt_struct)
279 basis_type_a=
"ORB", basis_type_b=
"ORB", &
285 n_act = gs_mos(ispin)%nmo_active
286 n_virt =
SIZE(gs_mos(ispin)%evals_virt)
288 ncol_global=n_act, context=blacs_env)
290 ncol_global=nao, context=blacs_env)
295 gs_mos(ispin)%mos_virt, &
297 ncol=n_virt, alpha=1.0_dp, beta=0.0_dp)
299 CALL cp_fm_get_info(soc_env%ediff(ispin), nrow_local=nrows_local, ncol_local=ncols_local, &
300 row_indices=row_indices, col_indices=col_indices, local_data=local_data_ediff)
305 DO icol = 1, ncols_local
307 ici = gs_mos(ispin)%index_active(col_indices(icol))
308 eval_occ = gs_mos(ispin)%evals_occ(ici)
310 DO irow = 1, nrows_local
313 local_data_ediff(irow, icol) = 1.0_dp/(gs_mos(ispin)%evals_virt(row_indices(irow)) - eval_occ)
318 DO ideriv = 1, nderivs
319 CALL cp_fm_create(soc_env%CdS(ispin, ideriv), cvirt_ao_struct)
322 CALL parallel_gemm(
'T',
'N', n_virt, nao, nao, 1.0_dp, gs_mos(ispin)%mos_virt, &
323 scrm_fm, 0.0_dp, soc_env%CdS(ispin, ideriv))
334 END SUBROUTINE velocity_rep
347 SUBROUTINE dip_vel_op(soc_env, qs_env, evec_fm, op, ideriv, tp, gs_coeffs, sggs_fm)
350 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(IN) :: evec_fm
352 INTEGER,
INTENT(IN) :: ideriv
353 LOGICAL,
INTENT(IN) :: tp
354 TYPE(
cp_fm_type),
OPTIONAL,
POINTER :: gs_coeffs
355 TYPE(
cp_fm_type),
INTENT(INOUT),
OPTIONAL :: sggs_fm
357 INTEGER :: iex, ispin, n_act, n_virt, nao, nex
361 TYPE(
cp_fm_type) :: cdsc, op_fm, scwcdsc, wcdsc
362 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:, :) :: wcdsc_tmp
366 NULLIFY (virt_occ_struct, virt_occ_struct, op_struct, blacs_env, para_env, coeff)
369 coeff => soc_env%b_coeff
371 coeff => soc_env%a_coeff
375 IF (
PRESENT(gs_coeffs)) sggs = .true.
378 nex =
SIZE(evec_fm, 2)
379 IF (.NOT. sggs)
ALLOCATE (wcdsc_tmp(ispin, nex))
380 CALL get_qs_env(qs_env, blacs_env=blacs_env, para_env=para_env)
381 CALL cp_fm_get_info(soc_env%CdS(ispin, ideriv), ncol_global=nao, nrow_global=n_virt)
385 CALL cp_fm_struct_create(virt_occ_struct, context=blacs_env, para_env=para_env, nrow_global=n_virt, &
387 CALL cp_fm_struct_create(op_struct, context=blacs_env, para_env=para_env, nrow_global=n_act*nex, &
390 CALL cp_fm_struct_create(virt_occ_struct, context=blacs_env, para_env=para_env, nrow_global=n_virt, &
391 ncol_global=n_act*nex)
392 CALL cp_fm_struct_create(op_struct, context=blacs_env, para_env=para_env, nrow_global=n_act*nex, &
393 ncol_global=n_act*nex)
396 CALL cp_fm_create(cdsc, soc_env%ediff(ispin)%matrix_struct)
401 CALL cp_fm_create(wcdsc, soc_env%ediff(ispin)%matrix_struct)
402 CALL parallel_gemm(
'N',
'N', n_virt, n_act, nao, 1.0_dp, soc_env%CdS(ispin, ideriv), &
403 gs_coeffs, 0.0_dp, cdsc)
408 CALL cp_fm_create(wcdsc_tmp(ispin, iex), soc_env%ediff(ispin)%matrix_struct)
409 CALL parallel_gemm(
'N',
'N', n_virt, n_act, nao, 1.0_dp, soc_env%CdS(ispin, ideriv), &
410 evec_fm(ispin, iex), 0.0_dp, cdsc)
418 DEALLOCATE (wcdsc_tmp)
422 CALL parallel_gemm(
'N',
'N', nao, n_act, n_virt, 1.0_dp, soc_env%SC(ispin), wcdsc, 0.0_dp, scwcdsc)
423 CALL parallel_gemm(
'T',
'N', n_act*nex, n_act, nao, 1.0_dp, soc_env%a_coeff, scwcdsc, 0.0_dp, op_fm)
425 CALL parallel_gemm(
'N',
'N', nao, n_act*nex, n_virt, 1.0_dp, soc_env%SC(ispin), wcdsc, 0.0_dp, scwcdsc)
426 CALL parallel_gemm(
'T',
'N', n_act*nex, n_act*nex, nao, 1.0_dp, coeff, scwcdsc, 0.0_dp, op_fm)
451 TYPE(
cp_fm_type),
DIMENSION(:, :),
INTENT(in) :: fm_start
454 CHARACTER(len=*),
PARAMETER :: routinen =
'soc_contract_evect'
456 INTEGER :: handle, ii, jj, nactive, nao, nspins, &
457 nstates, ntmp1, ntmp2
459 CALL timeset(routinen, handle)
461 nstates =
SIZE(fm_start, 2)
462 nspins =
SIZE(fm_start, 1)
468 CALL cp_fm_get_info(fm_start(jj, ii), nrow_global=nao, ncol_global=nactive)
474 1 + nactive*(ii - 1) + (jj - 1)*nao*nstates)
478 CALL timestop(handle)
489 SUBROUTINE test_repetition(vec, new_entry, res, res_int)
490 INTEGER,
DIMENSION(:),
INTENT(IN) :: vec
491 INTEGER,
INTENT(IN) :: new_entry
492 LOGICAL,
INTENT(OUT) :: res
493 INTEGER,
INTENT(OUT),
OPTIONAL :: res_int
498 IF (
PRESENT(res_int)) res_int = -1
501 IF (vec(i) == new_entry)
THEN
503 IF (
PRESENT(res_int)) res_int = i
508 END SUBROUTINE test_repetition
517 INTEGER,
ALLOCATABLE,
DIMENSION(:),
INTENT(OUT) :: sort
519 COMPLEX(dp),
ALLOCATABLE,
DIMENSION(:, :) :: cpl_tmp
520 INTEGER :: i_rep, ii, jj, ntot, tmp
521 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: rep_int
523 REAL(
dp) :: max_dev, max_wfn, wfn_sq
526 ALLOCATE (cpl_tmp(ntot, ntot))
527 ALLOCATE (sort(ntot), rep_int(ntot))
538 wfn_sq = abs(real(cpl_tmp(ii, jj)**2 - aimag(cpl_tmp(ii, jj)**2)))
539 IF (max_wfn <= wfn_sq)
THEN
540 CALL test_repetition(sort, ii, rep, rep_int(ii))
551 IF (rep_int(i_rep) > 0)
THEN
552 max_wfn = abs(real(cpl_tmp(sort(i_rep), jj)**2 - aimag(cpl_tmp(sort(i_rep), jj)**2))) - max_dev
554 wfn_sq = abs(real(cpl_tmp(ii, jj)**2 - aimag(cpl_tmp(ii, jj)**2)))
555 IF ((max_wfn - wfn_sq)/max_wfn <= max_dev)
THEN
556 CALL test_repetition(sort, ii, rep)
557 IF (rep .AND. ii /= i_rep)
THEN
558 sort(jj) = sort(i_rep)
568 DEALLOCATE (cpl_tmp, rep_int)
methods related to the blacs parallel environment
Represents a complex full matrix distributed on many processors.
subroutine, public cp_cfm_get_submatrix(fm, target_m, start_row, start_col, n_rows, n_cols, transpose)
Extract a sub-matrix from the full matrix: op(target_m)(1:n_rows,1:n_cols) = fm(start_row:start_row+n...
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.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_desymmetrize(matrix_a, matrix_b)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_release(matrix)
...
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 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_schur_product(matrix_a, matrix_b, matrix_c)
computes the schur product of two matrices c_ij = a_ij * b_ij
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_to_fm_submat(msource, mtarget, nrow, ncol, s_firstrow, s_firstcol, t_firstrow, t_firstcol)
copy just a part ot the 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
Defines the basic variable types.
integer, parameter, public dp
Interface to the message passing library MPI.
Calculates the moment integrals <a|r^m|b>
subroutine, public get_reference_point(rpoint, drpoint, qs_env, fist_env, reference, ref_point, ifirst, ilast)
...
basic linear algebra operations for full matrixes
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, 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, 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 neighbor list data types and the corresponding functionality.
subroutine, public p_xyz_ao(op, qs_env, minimum_image)
Calculation of the components of the dipole operator in the velocity form The elements of the sparse ...
subroutine, public rrc_xyz_ao(op, qs_env, rc, order, minimum_image, soft)
Calculation of the components of the dipole operator in the length form by taking the relative positi...
Calculation of overlap matrix, its derivatives and forces.
subroutine, public build_overlap_matrix(ks_env, matrix_s, matrixkp_s, matrix_name, nderivative, basis_type_a, basis_type_b, sab_nl, calculate_forces, matrix_p, matrixkp_p)
Calculation of the overlap matrix over Cartesian Gaussian functions.
Utilities absorption spectroscopy using TDDFPT with SOC.
subroutine, public soc_dipole_operator(soc_env, tddfpt_control, qs_env, gs_mos)
Build the atomic dipole operator.
subroutine, public resort_evects(evects_cfm, sort)
Used to find out, which state has which spin-multiplicity.
subroutine, public soc_contract_evect(fm_start, fm_res)
...
subroutine, public dip_vel_op(soc_env, qs_env, evec_fm, op, ideriv, tp, gs_coeffs, sggs_fm)
This routine will construct the dipol operator within velocity representation.
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
Represent a complex full matrix.
keeps the information about the structure of a full matrix
stores all the informations relevant to an mpi environment
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...
Ground state molecular orbitals.