41#include "./base/base_uses.f90"
46 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'molden_utils'
47 LOGICAL,
PARAMETER :: debug_this_module = .false.
49 INTEGER,
PARAMETER :: molden_lmax = 4
50 INTEGER,
PARAMETER :: molden_ncomax = (molden_lmax + 1)*(molden_lmax + 2)/2
70 CHARACTER(LEN=*),
PARAMETER :: routinen =
'write_mos_molden'
71 CHARACTER(LEN=molden_lmax+1),
PARAMETER :: angmom =
"spdfg"
73 CHARACTER(LEN=15) :: fmtstr1, fmtstr2
74 CHARACTER(LEN=2) :: element_symbol
75 INTEGER :: gto_kind, handle, i, iatom, icgf, icol, ikind, ipgf, irow, irow_in, iset, isgf, &
76 ishell, ispin, iw, lshell, ncgf, ncol_global, ndigits, nrow_global, nset, nsgf, z
77 INTEGER,
DIMENSION(:),
POINTER :: npgf, nshell
78 INTEGER,
DIMENSION(:, :),
POINTER :: l
79 INTEGER,
DIMENSION(molden_ncomax, 0:molden_lmax) :: orbmap
81 REAL(kind=
dp) :: expzet, prefac
82 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: cmatrix, smatrix
83 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: zet
84 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: gcc
88 CALL timeset(routinen, handle)
94 extension=
".molden", file_status=
'REPLACE')
99 ndigits = min(max(3, ndigits), 30)
100 WRITE (unit=fmtstr1, fmt=
'("(I6,1X,ES",I0,".",I0,")")') ndigits + 7, ndigits
101 WRITE (unit=fmtstr2, fmt=
'("((T51,2F",I0,".",I0,"))")') ndigits + 10, ndigits
105 IF (mos(1)%use_mo_coeff_b)
THEN
108 DO ispin = 1,
SIZE(mos)
109 IF (.NOT.
ASSOCIATED(mos(ispin)%mo_coeff_b))
THEN
118 WRITE (iw,
'(T2,A)')
"[Molden Format]"
119 WRITE (iw,
'(T2,A)')
"[Atoms] AU"
120 DO i = 1,
SIZE(particle_set)
122 element_symbol=element_symbol)
125 WRITE (iw,
'(T2,A2,I8,I8,3X,3(F12.6,3X))') &
126 element_symbol, i, z, particle_set(i)%r(:)
129 WRITE (iw,
'(T2,A)')
"[GTO]"
131 DO i = 1,
SIZE(particle_set)
132 CALL get_atomic_kind(atomic_kind=particle_set(i)%atomic_kind, kind_number=ikind, &
133 element_symbol=element_symbol)
134 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set)
135 IF (
ASSOCIATED(orb_basis_set))
THEN
136 WRITE (iw,
'(T2,I8,I8)') i, 0
146 DO ishell = 1, nshell(iset)
147 lshell = l(ishell, iset)
148 IF (lshell <= molden_lmax)
THEN
149 WRITE (unit=iw, fmt=
'(T25,A2,4X,I4,4X,F4.2)') &
150 angmom(lshell + 1:lshell + 1), npgf(iset), 1.0_dp
154 prefac = 2_dp**lshell*(2/
pi)**0.75_dp
155 expzet = 0.25_dp*(2*lshell + 3.0_dp)
156 WRITE (unit=iw, fmt=fmtstr2) &
157 (zet(ipgf, iset), gcc(ipgf, ishell, iset)/(prefac*zet(ipgf, iset)**expzet), &
161 CALL cp_warn(__location__, &
162 "MOLDEN format does not support Gaussian orbitals with l > 4.")
169 WRITE (iw,
'(A4)')
" "
176 WRITE (iw,
'(T2,A)')
"[5D7F]"
177 WRITE (iw,
'(T2,A)')
"[9G]"
180 WRITE (iw,
'(T2,A)')
"[MO]"
209 orbmap = reshape((/1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
210 1, 2, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
211 1, 4, 6, 2, 3, 5, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
212 1, 7, 10, 4, 2, 3, 6, 9, 8, 5, 0, 0, 0, 0, 0, &
213 1, 11, 15, 2, 3, 7, 12, 10, 14, 4, 6, 13, 5, 8, 9/), &
214 (/molden_ncomax, molden_lmax + 1/))
219 orbmap = reshape((/1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
220 3, 1, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
221 3, 4, 2, 5, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
222 4, 5, 3, 6, 2, 7, 1, 0, 0, 0, 0, 0, 0, 0, 0, &
223 5, 6, 4, 7, 3, 8, 2, 9, 1, 0, 0, 0, 0, 0, 0/), &
224 (/molden_ncomax, molden_lmax + 1/))
228 DO ispin = 1,
SIZE(mos)
230 nrow_global=nrow_global, &
231 ncol_global=ncol_global)
232 ALLOCATE (smatrix(nrow_global, ncol_global))
239 ALLOCATE (cmatrix(ncgf, ncgf))
247 DO iatom = 1,
SIZE(particle_set)
248 NULLIFY (orb_basis_set)
249 CALL get_atomic_kind(particle_set(iatom)%atomic_kind, kind_number=ikind)
251 basis_set=orb_basis_set)
252 IF (
ASSOCIATED(orb_basis_set))
THEN
258 DO ishell = 1, nshell(iset)
259 lshell = l(ishell, iset)
260 CALL dgemm(
"T",
"N",
nco(lshell), mos(ispin)%nmo,
nso(lshell), 1.0_dp, &
262 smatrix(isgf, 1), nsgf, 0.0_dp, &
263 cmatrix(icgf, 1), ncgf)
264 icgf = icgf +
nco(lshell)
265 isgf = isgf +
nso(lshell)
272 DO icol = 1, mos(ispin)%nmo
281 WRITE (iw,
'(A,ES20.10)')
'Ene=', mos(ispin)%eigenvalues(icol)
283 WRITE (iw,
'(A)')
'Spin= Alpha'
285 WRITE (iw,
'(A)')
'Spin= Beta'
287 WRITE (iw,
'(A,F12.7)')
'Occup=', mos(ispin)%occupation_numbers(icol)
289 DO iatom = 1,
SIZE(particle_set)
290 NULLIFY (orb_basis_set)
292 element_symbol=element_symbol, kind_number=ikind)
294 basis_set=orb_basis_set)
295 IF (
ASSOCIATED(orb_basis_set))
THEN
307 DO ishell = 1, nshell(iset)
308 lshell = l(ishell, iset)
310 IF (lshell <= molden_lmax)
THEN
311 CALL print_coeffs(iw, fmtstr1, ndigits, irow_in, orbmap(:, lshell), &
312 cmatrix(irow:irow +
nco(lshell) - 1, icol))
313 irow_in = irow_in +
nco(lshell)
316 irow = irow +
nco(lshell)
325 DO ishell = 1, nshell(iset)
326 lshell = l(ishell, iset)
328 IF (lshell <= molden_lmax)
THEN
329 CALL print_coeffs(iw, fmtstr1, ndigits, irow_in, orbmap(:, lshell), &
330 smatrix(irow:irow +
nso(lshell) - 1, icol))
331 irow_in = irow_in +
nso(lshell)
334 irow = irow +
nso(lshell)
344 IF (
ALLOCATED(cmatrix))
DEALLOCATE (cmatrix)
345 IF (
ALLOCATED(smatrix))
DEALLOCATE (smatrix)
352 CALL timestop(handle)
365 SUBROUTINE print_coeffs(iw, fmtstr1, ndigits, irow_in, orbmap, mo_coeff)
366 INTEGER,
INTENT(in) :: iw
367 CHARACTER(LEN=*),
INTENT(in) :: fmtstr1
368 INTEGER,
INTENT(in) :: ndigits, irow_in
369 INTEGER,
DIMENSION(molden_ncomax),
INTENT(in) :: orbmap
370 REAL(kind=
dp),
DIMENSION(:),
INTENT(in) :: mo_coeff
376 IF (abs(mo_coeff(orbmap(
orbital))) >= 10.0_dp**(-ndigits))
THEN
377 WRITE (iw, fmtstr1) irow_in +
orbital - 1, mo_coeff(orbmap(
orbital))
382 END SUBROUTINE print_coeffs
398 dump_only_positive, logger, list)
402 REAL(kind=
dp),
DIMENSION(:) :: freq
403 REAL(kind=
dp),
DIMENSION(:, :) :: eigen_vec
404 REAL(kind=
dp),
DIMENSION(:),
POINTER :: intensities
405 LOGICAL,
INTENT(in) :: calc_intens, dump_only_positive
407 INTEGER,
DIMENSION(:),
OPTIONAL,
POINTER ::
list
409 CHARACTER(len=*),
PARAMETER :: routinen =
'write_vibrations_molden'
411 CHARACTER(LEN=2) :: element_symbol
412 INTEGER :: handle, i, iw, j, k, l, z
413 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: my_list
414 REAL(kind=
dp) :: fint
416 CALL timeset(routinen, handle)
419 extension=
".mol", file_status=
'REPLACE')
422 cpassert(mod(
SIZE(eigen_vec, 1), 3) == 0)
423 cpassert(
SIZE(freq, 1) ==
SIZE(eigen_vec, 2))
424 ALLOCATE (my_list(
SIZE(particles)))
427 IF (
PRESENT(
list))
THEN
433 cpassert(
SIZE(particles) ==
SIZE(eigen_vec, 1)/3)
434 DO i = 1,
SIZE(my_list)
438 WRITE (iw,
'(T2,A)')
"[Molden Format]"
439 WRITE (iw,
'(T2,A)')
"[Atoms] AU"
440 DO i = 1,
SIZE(particles)
442 element_symbol=element_symbol)
445 WRITE (iw,
'(T2,A2,I8,I8,3X,3(F12.6,3X))') &
446 element_symbol, i, z, particles(i)%r(:)
449 WRITE (iw,
'(T2,A)')
"[FREQ]"
450 DO i = 1,
SIZE(freq, 1)
451 IF ((.NOT. dump_only_positive) .OR. (freq(i) >= 0._dp))
WRITE (iw,
'(T5,F12.6)') freq(i)
453 WRITE (iw,
'(T2,A)')
"[FR-COORD]"
454 DO i = 1,
SIZE(particles)
456 element_symbol=element_symbol)
457 WRITE (iw,
'(T2,A2,3X,3(F12.6,3X))') &
458 element_symbol, particles(i)%r(:)
460 WRITE (iw,
'(T2,A)')
"[FR-NORM-COORD]"
462 DO i = 1,
SIZE(eigen_vec, 2)
463 IF ((.NOT. dump_only_positive) .OR. (freq(i) >= 0._dp))
THEN
465 WRITE (iw,
'(T2,A,1X,I6)')
"vibration", l
466 DO j = 1,
SIZE(particles)
467 IF (my_list(j) .NE. 0)
THEN
468 k = (my_list(j) - 1)*3
469 WRITE (iw,
'(T2,3(F12.6,3X))') eigen_vec(k + 1, i), eigen_vec(k + 2, i), eigen_vec(k + 3, i)
471 WRITE (iw,
'(T2,3(F12.6,3X))') 0.0_dp, 0.0_dp, 0.0_dp
476 IF (calc_intens)
THEN
479 WRITE (iw,
'(T2,A)')
"[INT]"
480 DO i = 1,
SIZE(intensities)
481 IF ((.NOT. dump_only_positive) .OR. (freq(i) >= 0._dp))
WRITE (iw,
'(3X,F18.6)') fint*intensities(i)**2
488 CALL timestop(handle)
static void dgemm(const char transa, const char transb, const int m, const int n, const int k, const double alpha, const double *a, const int lda, const double *b, const int ldb, const double beta, double *c, const int ldc)
Convenient wrapper to hide Fortran nature of dgemm_, swapping a and b.
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
subroutine, public get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, npgf_sum, nshell_sum, maxder, short_kind_radius, npgf_seg_sum)
...
DBCSR operations in CP2K.
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
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_get_submatrix(fm, target_m, start_row, start_col, n_rows, n_cols, transpose)
gets a submatrix of a full matrix op(target_m)(1:n_rows,1:n_cols) =fm(start_row:start_row+n_rows,...
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...
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...
Defines the basic variable types.
integer, parameter, public dp
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
Functions handling the MOLDEN format. Split from mode_selective.
subroutine, public write_mos_molden(mos, qs_kind_set, particle_set, print_section)
Write out the MOs in molden format for visualisation.
subroutine, public write_vibrations_molden(input, particles, freq, eigen_vec, intensities, calc_intens, dump_only_positive, logger, list)
writes the output for vibrational analysis in MOLDEN format
Provides Cartesian and spherical orbital pointers and indices.
integer, dimension(:), allocatable, public nco
integer, dimension(:), allocatable, public nso
Define the data structure for the particle information.
Periodic Table related data definitions.
subroutine, public get_ptable_info(symbol, number, amass, ielement, covalent_radius, metallic_radius, vdw_radius, found)
Pass information about the kind given the element symbol.
Definition of physical constants:
real(kind=dp), parameter, public massunit
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
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)
Get attributes of an atomic kind set.
Definition and initialisation of the mo data type.
type of a logger, at the moment it contains just a print level starting at which level it should be l...
Orbital angular momentum.
Provides all information about a quickstep kind.