(git:374b731)
Loading...
Searching...
No Matches
sirius_interface.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8!***************************************************************************************************
9!> \brief Interface to the SIRIUS Library
10!> \par History
11!> 07.2018 initial create
12!> \author JHU
13!***************************************************************************************************
14#if defined(__SIRIUS)
16 USE iso_c_binding, ONLY: c_double,&
17 c_int,&
18 c_loc
22 USE atom_upf, ONLY: atom_upfpot_type
26 USE cell_types, ONLY: cell_type,&
31 USE cp_output_handling, ONLY: cp_p_file,&
45 USE kinds, ONLY: default_string_length,&
46 dp
47 USE mathconstants, ONLY: fourpi,&
48 gamma1
51 USE physcon, ONLY: massunit
60 USE qs_kind_types, ONLY: get_qs_kind,&
64 USE sirius, ONLY: &
65 sirius_integer_array_type, sirius_integer_type, sirius_logical_array_type, &
66 sirius_logical_type, sirius_number_array_type, sirius_number_type, &
67 sirius_string_array_type, sirius_string_type, sirius_add_atom, sirius_add_atom_type, &
68 sirius_add_atom_type_radial_function, sirius_add_xc_functional, sirius_context_handler, &
69 sirius_create_context, sirius_create_ground_state, sirius_create_kset_from_grid, &
70 sirius_finalize, sirius_find_ground_state, sirius_get_band_energies, &
71 sirius_get_band_occupancies, sirius_get_energy, sirius_get_forces, &
72 sirius_get_kpoint_properties, sirius_get_num_kpoints, sirius_get_parameters, &
73 sirius_get_stress_tensor, sirius_ground_state_handler, sirius_import_parameters, &
74 sirius_initialize, sirius_initialize_context, sirius_kpoint_set_handler, &
75 sirius_option_get_info, sirius_option_get_section_length, sirius_option_set, &
76 sirius_set_atom_position, sirius_set_atom_type_dion, sirius_set_atom_type_hubbard, &
77 sirius_set_atom_type_radial_grid, sirius_set_lattice_vectors, sirius_set_mpi_grid_dims, &
78 sirius_update_ground_state
79#include "./base/base_uses.f90"
80
81 IMPLICIT NONE
82
83 PRIVATE
84
85! *** Global parameters ***
86
87 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'sirius_interface'
88
89! *** Public subroutines ***
90
93
94CONTAINS
95
96!***************************************************************************************************
97!> \brief ...
98!> \param
99!> \par History
100!> 07.2018 start the Sirius library
101!> \author JHU
102! **************************************************************************************************
103 SUBROUTINE cp_sirius_init()
104 CALL sirius_initialize(.false.)
105 END SUBROUTINE cp_sirius_init
106
107!***************************************************************************************************
108!> \brief ...
109!> \param
110!> \par History
111!> 07.2018 stop the Sirius library
112!> \author JHU
113! **************************************************************************************************
114 SUBROUTINE cp_sirius_finalize()
115 CALL sirius_finalize(.false., .false., .false.)
116 END SUBROUTINE cp_sirius_finalize
117
118!***************************************************************************************************
119!> \brief ...
120!> \param pwdft_env ...
121!> \param
122!> \par History
123!> 07.2018 Create the Sirius environment
124!> \author JHU
125! **************************************************************************************************
126 SUBROUTINE cp_sirius_create_env(pwdft_env)
127 TYPE(pwdft_environment_type), POINTER :: pwdft_env
128#if defined(__SIRIUS)
129
130 CHARACTER(len=2) :: element_symbol
131 CHARACTER(len=default_string_length) :: label
132 INTEGER :: i, iatom, ibeta, ifun, ikind, iwf, j, l, &
133 n, natom, nbeta, nkind, nmesh, &
134 num_mag_dims, sirius_mpi_comm, vdw_func, nu, lu
135 INTEGER, DIMENSION(:), POINTER :: mpi_grid_dims
136 INTEGER(KIND=C_INT), DIMENSION(3) :: k_grid, k_shift
137 INTEGER, DIMENSION(:), POINTER :: kk
138 LOGICAL :: up, use_ref_cell
139 LOGICAL(4) :: use_so, use_symmetry, dft_plus_u_atom
140 REAL(KIND=c_double), ALLOCATABLE, DIMENSION(:) :: fun
141 REAL(KIND=c_double), ALLOCATABLE, DIMENSION(:, :) :: dion
142 REAL(KIND=c_double), DIMENSION(3) :: a1, a2, a3, v1, v2
143 REAL(KIND=dp) :: al, angle1, angle2, cval, focc, &
144 magnetization, mass, pf, rl, zeff, alpha_u, beta_u, &
145 j0_u, j_u, u_u, occ_u, u_minus_j
146 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: beta, corden, ef, fe, locpot, rc, rp
147 REAL(KIND=dp), DIMENSION(3) :: vr, vs
148 REAL(KIND=dp), DIMENSION(:), POINTER :: density
149 REAL(KIND=dp), DIMENSION(:, :), POINTER :: wavefunction, wfninfo
150 TYPE(atom_gthpot_type), POINTER :: gth_atompot
151 TYPE(atom_upfpot_type), POINTER :: upf_pot
152 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
153 TYPE(atomic_kind_type), POINTER :: atomic_kind
154 TYPE(cell_type), POINTER :: my_cell
155 TYPE(mp_para_env_type), POINTER :: para_env
156 TYPE(grid_atom_type), POINTER :: atom_grid
157 TYPE(gth_potential_type), POINTER :: gth_potential
158 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
159 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
160 TYPE(qs_subsys_type), POINTER :: qs_subsys
161 TYPE(section_vals_type), POINTER :: pwdft_section, pwdft_sub_section, &
162 xc_fun, xc_section
163 TYPE(sirius_context_handler) :: sctx
164 TYPE(sirius_ground_state_handler) :: gs_handler
165 TYPE(sirius_kpoint_set_handler) :: ks_handler
166
167 cpassert(ASSOCIATED(pwdft_env))
168
169 ! create context of simulation
170 CALL pwdft_env_get(pwdft_env, para_env=para_env)
171 sirius_mpi_comm = para_env%get_handle()
172 CALL sirius_create_context(sirius_mpi_comm, sctx)
173
174! the "fun" starts.
175
176 CALL pwdft_env_get(pwdft_env=pwdft_env, pwdft_input=pwdft_section, xc_input=xc_section)
177
178 ! cp2k should *have* a function that return all xc_functionals. Doing
179 ! manually is prone to errors
180
181 IF (ASSOCIATED(xc_section)) THEN
182 ifun = 0
183 DO
184 ifun = ifun + 1
185 xc_fun => section_vals_get_subs_vals2(xc_section, i_section=ifun)
186 IF (.NOT. ASSOCIATED(xc_fun)) EXIT
187 ! Here, we do not have to check whether the functional name starts with XC_
188 ! because we only allow the shorter form w/o XC_
189 CALL sirius_add_xc_functional(sctx, "XC_"//trim(xc_fun%section%name))
190 END DO
191 END IF
192
193 ! import control section
194 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "control")
195 IF (ASSOCIATED(pwdft_sub_section)) THEN
196 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "control")
197 CALL section_vals_val_get(pwdft_sub_section, "mpi_grid_dims", i_vals=mpi_grid_dims)
198 END IF
199
200! import parameters section
201 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "parameters")
202
203 IF (ASSOCIATED(pwdft_sub_section)) THEN
204 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "parameters")
205 CALL section_vals_val_get(pwdft_sub_section, "ngridk", i_vals=kk)
206 k_grid(1) = kk(1)
207 k_grid(2) = kk(2)
208 k_grid(3) = kk(3)
209
210 CALL section_vals_val_get(pwdft_sub_section, "shiftk", i_vals=kk)
211 k_shift(1) = kk(1)
212 k_shift(2) = kk(2)
213 k_shift(3) = kk(3)
214 CALL section_vals_val_get(pwdft_sub_section, "num_mag_dims", i_val=num_mag_dims)
215 CALL section_vals_val_get(pwdft_sub_section, "use_symmetry", l_val=use_symmetry)
216 CALL section_vals_val_get(pwdft_sub_section, "so_correction", l_val=use_so)
217
218! now check if van der walls corrections are needed
219 vdw_func = -1
220#ifdef __LIBVDWXC
221 CALL section_vals_val_get(pwdft_sub_section, "vdw_functional", i_val=vdw_func)
222 SELECT CASE (vdw_func)
223 CASE (sirius_func_vdwdf)
224 CALL sirius_add_xc_functional(sctx, "XC_FUNC_VDWDF")
225 CASE (sirius_func_vdwdf2)
226 CALL sirius_add_xc_functional(sctx, "XC_FUNC_VDWDF2")
228 CALL sirius_add_xc_functional(sctx, "XC_FUNC_VDWDF2")
229 CASE default
230 END SELECT
231#endif
232
233 END IF
234
235! import mixer section
236 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "mixer")
237 IF (ASSOCIATED(pwdft_sub_section)) THEN
238 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "mixer")
239 END IF
240
241! import settings section
242 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "settings")
243 IF (ASSOCIATED(pwdft_sub_section)) THEN
244 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "settings")
245 END IF
246
247 ! import solver section
248 pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "iterative_solver")
249 IF (ASSOCIATED(pwdft_sub_section)) THEN
250 CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "iterative_solver")
251 END IF
252
253 !
254 ! uncomment these lines when nlcg is officially supported
255 !
256
257 ! import nlcg section
258 ! pwdft_sub_section => section_vals_get_subs_vals(pwdft_section, "nlcg")
259 ! IF (ASSOCIATED(pwdft_sub_section)) THEN
260 ! CALL cp_sirius_fill_in_section(sctx, pwdft_sub_section, "nlcg")
261 ! ENDIF
262
263 !CALL sirius_dump_runtime_setup(sctx, "runtime.json")
264 CALL sirius_import_parameters(sctx, '{}')
265
266! lattice vectors of the unit cell should be in [a.u.] (length is in [a.u.])
267 CALL pwdft_env_get(pwdft_env=pwdft_env, qs_subsys=qs_subsys)
268 CALL qs_subsys_get(qs_subsys, cell=my_cell, use_ref_cell=use_ref_cell)
269 a1(:) = my_cell%hmat(:, 1)
270 a2(:) = my_cell%hmat(:, 2)
271 a3(:) = my_cell%hmat(:, 3)
272 CALL sirius_set_lattice_vectors(sctx, a1(1), a2(1), a3(1))
273
274 IF (use_ref_cell) THEN
275 cpwarn("SIRIUS| The specified CELL_REF will be ignored for PW_DFT calculations")
276 END IF
277
278! set up the atomic type definitions
279 CALL qs_subsys_get(qs_subsys, &
280 atomic_kind_set=atomic_kind_set, &
281 qs_kind_set=qs_kind_set, &
282 particle_set=particle_set)
283 nkind = SIZE(atomic_kind_set)
284 DO ikind = 1, nkind
285 CALL get_atomic_kind(atomic_kind_set(ikind), &
286 name=label, element_symbol=element_symbol, mass=mass)
287 CALL get_qs_kind(qs_kind_set(ikind), zeff=zeff)
288 NULLIFY (upf_pot, gth_potential)
289 CALL get_qs_kind(qs_kind_set(ikind), upf_potential=upf_pot, gth_potential=gth_potential)
290
291 IF (ASSOCIATED(upf_pot)) THEN
292 CALL sirius_add_atom_type(sctx, label, fname=upf_pot%filename, &
293 symbol=element_symbol, &
294 mass=real(mass/massunit, kind=c_double))
295
296 ELSEIF (ASSOCIATED(gth_potential)) THEN
297!
298 NULLIFY (atom_grid)
299 CALL allocate_grid_atom(atom_grid)
300 nmesh = 929
301 atom_grid%nr = nmesh
302 CALL create_grid_atom(atom_grid, nmesh, 1, 1, 0, do_gapw_log)
303 ALLOCATE (rp(nmesh), fun(nmesh))
304 IF (atom_grid%rad(1) < atom_grid%rad(nmesh)) THEN
305 up = .true.
306 ELSE
307 up = .false.
308 END IF
309 IF (up) THEN
310 rp(1:nmesh) = atom_grid%rad(1:nmesh)
311 ELSE
312 DO i = 1, nmesh
313 rp(i) = atom_grid%rad(nmesh - i + 1)
314 END DO
315 END IF
316! add new atom type
317 CALL sirius_add_atom_type(sctx, label, &
318 zn=nint(zeff + 0.001d0), &
319 symbol=element_symbol, &
320 mass=real(mass/massunit, kind=c_double), &
321 spin_orbit=.false.)
322!
323 ALLOCATE (gth_atompot)
324 CALL gth_potential_conversion(gth_potential, gth_atompot)
325! set radial grid
326 fun(1:nmesh) = rp(1:nmesh)
327 CALL sirius_set_atom_type_radial_grid(sctx, label, nmesh, fun(1))
328! set beta-projectors
329 ALLOCATE (ef(nmesh), beta(nmesh))
330 ibeta = 0
331 DO l = 0, 3
332 IF (gth_atompot%nl(l) == 0) cycle
333 rl = gth_atompot%rcnl(l)
334! we need to multiply by r so that data transferred to sirius are r \beta(r) not beta(r)
335 ef(1:nmesh) = exp(-0.5_dp*rp(1:nmesh)*rp(1:nmesh)/(rl*rl))
336 DO i = 1, gth_atompot%nl(l)
337 pf = rl**(l + 0.5_dp*(4._dp*i - 1._dp))
338 j = l + 2*i - 1
339 pf = sqrt(2._dp)/(pf*sqrt(gamma1(j)))
340 beta(:) = pf*rp**(l + 2*i - 2)*ef
341 ibeta = ibeta + 1
342 fun(1:nmesh) = beta(1:nmesh)*rp(1:nmesh)
343 CALL sirius_add_atom_type_radial_function(sctx, label, &
344 "beta", fun(1), nmesh, l=l)
345 END DO
346 END DO
347 DEALLOCATE (ef, beta)
348 nbeta = ibeta
349
350! nonlocal PP matrix elements
351 ALLOCATE (dion(nbeta, nbeta))
352 dion = 0.0_dp
353 DO l = 0, 3
354 IF (gth_atompot%nl(l) == 0) cycle
355 ibeta = sum(gth_atompot%nl(0:l - 1)) + 1
356 i = ibeta + gth_atompot%nl(l) - 1
357 dion(ibeta:i, ibeta:i) = gth_atompot%hnl(1:gth_atompot%nl(l), 1:gth_atompot%nl(l), l)
358 END DO
359 CALL sirius_set_atom_type_dion(sctx, label, nbeta, dion(1, 1))
360 DEALLOCATE (dion)
361
362! set non-linear core correction
363 IF (gth_atompot%nlcc) THEN
364 ALLOCATE (corden(nmesh), fe(nmesh), rc(nmesh))
365 corden(:) = 0.0_dp
366 n = gth_atompot%nexp_nlcc
367 DO i = 1, n
368 al = gth_atompot%alpha_nlcc(i)
369 rc(:) = rp(:)/al
370 fe(:) = exp(-0.5_dp*rc(:)*rc(:))
371 DO j = 1, gth_atompot%nct_nlcc(i)
372 cval = gth_atompot%cval_nlcc(j, i)
373 corden(:) = corden(:) + fe(:)*rc(:)**(2*j - 2)*cval
374 END DO
375 END DO
376 fun(1:nmesh) = corden(1:nmesh)*rp(1:nmesh)
377 CALL sirius_add_atom_type_radial_function(sctx, label, "ps_rho_core", &
378 fun(1), nmesh)
379 DEALLOCATE (corden, fe, rc)
380 END IF
381
382! local potential
383 ALLOCATE (locpot(nmesh))
384 locpot(:) = 0.0_dp
385 CALL atom_local_potential(locpot, gth_atompot, rp)
386 fun(1:nmesh) = locpot(1:nmesh)
387 CALL sirius_add_atom_type_radial_function(sctx, label, "vloc", &
388 fun(1), nmesh)
389 DEALLOCATE (locpot)
390!
391 NULLIFY (density, wavefunction, wfninfo)
392 CALL calculate_atomic_orbitals(atomic_kind_set(ikind), qs_kind_set(ikind), &
393 density=density, wavefunction=wavefunction, &
394 wfninfo=wfninfo, agrid=atom_grid)
395
396! set the atomic radial functions
397 DO iwf = 1, SIZE(wavefunction, 2)
398 focc = wfninfo(1, iwf)
399 l = nint(wfninfo(2, iwf))
400! we can not easily get the principal quantum number
401 nu = -1
402 IF (up) THEN
403 fun(1:nmesh) = wavefunction(1:nmesh, iwf)*rp(i)
404 ELSE
405 DO i = 1, nmesh
406 fun(i) = wavefunction(nmesh - i + 1, iwf)*rp(i)
407 END DO
408 END IF
409 CALL sirius_add_atom_type_radial_function(sctx, &
410 label, "ps_atomic_wf", &
411 fun(1), nmesh, l=l, occ=real(focc, kind=c_double), n=nu)
412 END DO
413
414! set total charge density of a free atom (to compute initial rho(r))
415 IF (up) THEN
416 fun(1:nmesh) = fourpi*density(1:nmesh)*atom_grid%rad(1:nmesh)**2
417 ELSE
418 DO i = 1, nmesh
419 fun(i) = fourpi*density(nmesh - i + 1)*atom_grid%rad(nmesh - i + 1)**2
420 END DO
421 END IF
422 CALL sirius_add_atom_type_radial_function(sctx, label, "ps_rho_total", &
423 fun(1), nmesh)
424
425 IF (ASSOCIATED(density)) DEALLOCATE (density)
426 IF (ASSOCIATED(wavefunction)) DEALLOCATE (wavefunction)
427 IF (ASSOCIATED(wfninfo)) DEALLOCATE (wfninfo)
428
429 CALL deallocate_grid_atom(atom_grid)
430 DEALLOCATE (rp, fun)
431 DEALLOCATE (gth_atompot)
432!
433 ELSE
434 CALL cp_abort(__location__, &
435 "CP2K/SIRIUS: atomic kind needs UPF or GTH potential definition")
436 END IF
437
438 CALL get_qs_kind(qs_kind_set(ikind), &
439 dft_plus_u_atom=dft_plus_u_atom, &
440 l_of_dft_plus_u=lu, &
441 n_of_dft_plus_u=nu, &
442 u_minus_j_target=u_minus_j, &
443 u_of_dft_plus_u=u_u, &
444 j_of_dft_plus_u=j_u, &
445 alpha_of_dft_plus_u=alpha_u, &
446 beta_of_dft_plus_u=beta_u, &
447 j0_of_dft_plus_u=j0_u, &
448 occupation_of_dft_plus_u=occ_u)
449
450 IF (dft_plus_u_atom) THEN
451 IF (nu < 1) THEN
452 cpabort("CP2K/SIRIUS (hubbard): principal quantum number not specified")
453 END IF
454
455 IF (lu < 0) THEN
456 cpabort("CP2K/SIRIUS (hubbard): l can not be negative.")
457 END IF
458
459 IF (occ_u < 0.0) THEN
460 cpabort("CP2K/SIRIUS (hubbard): the occupation number can not be negative.")
461 END IF
462
463 IF (abs(u_minus_j) < 1e-8) THEN
464 CALL sirius_set_atom_type_hubbard(sctx, label, lu, nu, &
465 occ_u, u_u, j_u, alpha_u, beta_u, j0_u)
466 ELSE
467 CALL sirius_set_atom_type_hubbard(sctx, label, lu, nu, &
468 occ_u, u_minus_j, 0.0_dp, alpha_u, beta_u, j0_u)
469 END IF
470 END IF
471
472 END DO
473
474! add atoms to the unit cell
475! WARNING: sirius accepts only fractional coordinates;
476 natom = SIZE(particle_set)
477 DO iatom = 1, natom
478 vr(1:3) = particle_set(iatom)%r(1:3)
479 CALL real_to_scaled(vs, vr, my_cell)
480 atomic_kind => particle_set(iatom)%atomic_kind
481 ikind = atomic_kind%kind_number
482 CALL get_atomic_kind(atomic_kind, name=label)
483 CALL get_qs_kind(qs_kind_set(ikind), zeff=zeff, magnetization=magnetization)
484! angle of magnetization might come from input Atom x y z mx my mz
485! or as an angle?
486! Answer : SIRIUS only accept the magnetization as mx, my, mz
487 IF (num_mag_dims .EQ. 3) THEN
488 angle1 = 0.0_dp
489 angle2 = 0.0_dp
490 v1(1) = magnetization*sin(angle1)*cos(angle2)
491 v1(2) = magnetization*sin(angle1)*sin(angle2)
492 v1(3) = magnetization*cos(angle1)
493 ELSE
494 v1 = 0._dp
495 v1(3) = magnetization
496 END IF
497 v2(1:3) = vs(1:3)
498 CALL sirius_add_atom(sctx, label, v2(1), v1(1))
499 END DO
500
501 CALL sirius_set_mpi_grid_dims(sctx, 2, mpi_grid_dims)
502
503! initialize global variables/indices/arrays/etc. of the simulation
504 CALL sirius_initialize_context(sctx)
505
506 ! strictly speaking the parameter use_symmetry is initialized at the
507 ! beginning but it does no harm to do it that way
508 IF (use_symmetry) THEN
509 CALL sirius_create_kset_from_grid(sctx, k_grid(1), k_shift(1), use_symmetry=.true., kset_handler=ks_handler)
510 ELSE
511 CALL sirius_create_kset_from_grid(sctx, k_grid(1), k_shift(1), use_symmetry=.false., kset_handler=ks_handler)
512 END IF
513! create ground-state class
514 CALL sirius_create_ground_state(ks_handler, gs_handler)
515
516 CALL pwdft_env_set(pwdft_env, sctx=sctx, gs_handler=gs_handler, ks_handler=ks_handler)
517#endif
518 END SUBROUTINE cp_sirius_create_env
519
520!***************************************************************************************************
521!> \brief ...
522!> \param pwdft_env ...
523!> \param
524!> \par History
525!> 07.2018 Update the Sirius environment
526!> \author JHU
527! **************************************************************************************************
528 SUBROUTINE cp_sirius_update_context(pwdft_env)
529 TYPE(pwdft_environment_type), POINTER :: pwdft_env
530
531 INTEGER :: iatom, natom
532 REAL(KIND=c_double), DIMENSION(3) :: a1, a2, a3, v2
533 REAL(KIND=dp), DIMENSION(3) :: vr, vs
534 TYPE(cell_type), POINTER :: my_cell
535 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
536 TYPE(qs_subsys_type), POINTER :: qs_subsys
537 TYPE(sirius_context_handler) :: sctx
538 TYPE(sirius_ground_state_handler) :: gs_handler
539
540 cpassert(ASSOCIATED(pwdft_env))
541 CALL pwdft_env_get(pwdft_env, sctx=sctx, gs_handler=gs_handler)
542
543! get current positions and lattice vectors
544 CALL pwdft_env_get(pwdft_env=pwdft_env, qs_subsys=qs_subsys)
545
546! lattice vectors of the unit cell should be in [a.u.] (length is in [a.u.])
547 CALL qs_subsys_get(qs_subsys, cell=my_cell)
548 a1(:) = my_cell%hmat(:, 1)
549 a2(:) = my_cell%hmat(:, 2)
550 a3(:) = my_cell%hmat(:, 3)
551 CALL sirius_set_lattice_vectors(sctx, a1(1), a2(1), a3(1))
552
553! new atomic positions
554 CALL qs_subsys_get(qs_subsys, particle_set=particle_set)
555 natom = SIZE(particle_set)
556 DO iatom = 1, natom
557 vr(1:3) = particle_set(iatom)%r(1:3)
558 CALL real_to_scaled(vs, vr, my_cell)
559 v2(1:3) = vs(1:3)
560 CALL sirius_set_atom_position(sctx, iatom, v2(1))
561 END DO
562
563! update ground-state class
564 CALL sirius_update_ground_state(gs_handler)
565
566 CALL pwdft_env_set(pwdft_env, sctx=sctx, gs_handler=gs_handler)
567
568 END SUBROUTINE cp_sirius_update_context
569
570! **************************************************************************************************
571!> \brief ...
572!> \param sctx ...
573!> \param section ...
574!> \param section_name ...
575! **************************************************************************************************
576 SUBROUTINE cp_sirius_fill_in_section(sctx, section, section_name)
577 TYPE(sirius_context_handler), INTENT(INOUT) :: sctx
578 TYPE(section_vals_type), POINTER :: section
579 CHARACTER(*), INTENT(in) :: section_name
580
581 CHARACTER(len=256), TARGET :: option_name
582 CHARACTER(len=4096) :: description, usage
583 CHARACTER(len=80), DIMENSION(:), POINTER :: tmp
584 CHARACTER(len=80), TARGET :: str
585 INTEGER :: ctype, elem, ic, j
586 INTEGER, DIMENSION(:), POINTER :: ivals
587 INTEGER, TARGET :: enum_length, ival, length, &
588 num_possible_values, number_of_options
589 LOGICAL :: explicit
590 LOGICAL, DIMENSION(:), POINTER :: lvals
591 LOGICAL, TARGET :: found, lval
592 REAL(kind=dp), DIMENSION(:), POINTER :: rvals
593 REAL(kind=dp), TARGET :: rval
594
595 NULLIFY (rvals)
596 NULLIFY (ivals)
597 CALL sirius_option_get_section_length(section_name, number_of_options)
598
599 DO elem = 1, number_of_options
600 option_name = ''
601 CALL sirius_option_get_info(section_name, &
602 elem, &
603 option_name, &
604 256, &
605 ctype, &
606 num_possible_values, &
607 enum_length, &
608 description, &
609 4096, &
610 usage, &
611 4096)
612 IF ((option_name /= 'memory_usage') .AND. (option_name /= 'xc_functionals') .AND. (option_name /= 'vk')) THEN
613 CALL section_vals_val_get(section, option_name, explicit=found)
614 IF (found) THEN
615 SELECT CASE (ctype)
616 CASE (sirius_integer_type)
617 CALL section_vals_val_get(section, option_name, i_val=ival)
618 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(ival))
619 CASE (sirius_number_type)
620 CALL section_vals_val_get(section, option_name, r_val=rval)
621 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(rval))
622 CASE (sirius_logical_type)
623 CALL section_vals_val_get(section, option_name, l_val=lval)
624 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(lval))
625 CASE (sirius_string_type) ! string nightmare
626 str = ''
627 CALL section_vals_val_get(section, option_name, explicit=explicit, c_val=str)
628 str = trim(adjustl(str))
629 DO j = 1, len(str)
630 ic = ichar(str(j:j))
631 IF (ic >= 65 .AND. ic < 90) str(j:j) = char(ic + 32)
632 END DO
633
634 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(str), max_length=len_trim(str))
635 CASE (sirius_integer_array_type)
636 CALL section_vals_val_get(section, option_name, i_vals=ivals)
637 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(ivals(1)), &
638 max_length=num_possible_values)
639 CASE (sirius_number_array_type)
640 CALL section_vals_val_get(section, option_name, r_vals=rvals)
641 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(rvals(1)), &
642 max_length=num_possible_values)
643 CASE (sirius_logical_array_type)
644 CALL section_vals_val_get(section, option_name, l_vals=lvals)
645 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(lvals(1)), &
646 max_length=num_possible_values)
647 CASE (sirius_string_array_type)
648 CALL section_vals_val_get(section, option_name, explicit=explicit, n_rep_val=length)
649 DO j = 1, length
650 str = ''
651 CALL section_vals_val_get(section, option_name, i_rep_val=j, explicit=explicit, c_vals=tmp)
652 str = trim(adjustl(tmp(j)))
653 CALL sirius_option_set(sctx, section_name, option_name, ctype, c_loc(str), &
654 max_length=len_trim(str), append=.true.)
655 END DO
656 CASE DEFAULT
657 END SELECT
658 END IF
659 END IF
660 END DO
661 END SUBROUTINE cp_sirius_fill_in_section
662
663!***************************************************************************************************
664!> \brief ...
665!> \param pwdft_env ...
666!> \param calculate_forces ...
667!> \param calculate_stress_tensor ...
668!> \param
669!> \par History
670!> 07.2018 start the Sirius library
671!> \author JHU
672! **************************************************************************************************
673 SUBROUTINE cp_sirius_energy_force(pwdft_env, calculate_forces, calculate_stress_tensor)
674 TYPE(pwdft_environment_type), INTENT(INOUT), &
675 POINTER :: pwdft_env
676 LOGICAL, INTENT(IN) :: calculate_forces, calculate_stress_tensor
677
678 INTEGER :: n1, n2
679 LOGICAL :: do_print
680 REAL(KIND=c_double) :: etotal
681 REAL(KIND=c_double), ALLOCATABLE, DIMENSION(:, :) :: cforces
682 REAL(KIND=c_double), DIMENSION(3, 3) :: cstress
683 REAL(KIND=dp), DIMENSION(3, 3) :: stress
684 REAL(KIND=dp), DIMENSION(:, :), POINTER :: forces
685 TYPE(pwdft_energy_type), POINTER :: energy
686 TYPE(section_vals_type), POINTER :: print_section, pwdft_input
687 TYPE(sirius_ground_state_handler) :: gs_handler
688
689 cpassert(ASSOCIATED(pwdft_env))
690
691 CALL pwdft_env_get(pwdft_env=pwdft_env, gs_handler=gs_handler)
692 CALL sirius_find_ground_state(gs_handler)
693 CALL pwdft_env_get(pwdft_env=pwdft_env, energy=energy)
694 etotal = 0.0_c_double
695
696 CALL sirius_get_energy(gs_handler, 'band-gap', etotal)
697 energy%band_gap = etotal
698
699 etotal = 0.0_c_double
700 CALL sirius_get_energy(gs_handler, 'total', etotal)
701 energy%etotal = etotal
702
703 ! extract entropy (TS returned by sirius is always negative, sign
704 ! convention in QE)
705 etotal = 0.0_c_double
706 CALL sirius_get_energy(gs_handler, 'demet', etotal)
707 energy%entropy = -etotal
708
709 IF (calculate_forces) THEN
710 CALL pwdft_env_get(pwdft_env=pwdft_env, forces=forces)
711 n1 = SIZE(forces, 1)
712 n2 = SIZE(forces, 2)
713
714 ALLOCATE (cforces(n2, n1))
715 cforces = 0.0_c_double
716 CALL sirius_get_forces(gs_handler, 'total', cforces)
717 ! Sirius computes the forces but cp2k use the gradient everywhere
718 ! so a minus sign is needed.
719 ! note also that sirius and cp2k store the forces transpose to each other
720 ! sirius : forces(coordinates, atoms)
721 ! cp2k : forces(atoms, coordinates)
722 forces = -transpose(cforces(:, :))
723 DEALLOCATE (cforces)
724 END IF
725
726 IF (calculate_stress_tensor) THEN
727 cstress = 0.0_c_double
728 CALL sirius_get_stress_tensor(gs_handler, 'total', cstress)
729 stress(1:3, 1:3) = cstress(1:3, 1:3)
730 CALL pwdft_env_set(pwdft_env=pwdft_env, stress=stress)
731 END IF
732
733 CALL pwdft_env_get(pwdft_env=pwdft_env, pwdft_input=pwdft_input)
734 print_section => section_vals_get_subs_vals(pwdft_input, "PRINT")
735 CALL section_vals_get(print_section, explicit=do_print)
736 IF (do_print) THEN
737 CALL cp_sirius_print_results(pwdft_env, print_section)
738 END IF
739
740 END SUBROUTINE cp_sirius_energy_force
741
742!***************************************************************************************************
743!> \brief ...
744!> \param pwdft_env ...
745!> \param print_section ...
746!> \param
747!> \par History
748!> 12.2019 init
749!> \author JHU
750! **************************************************************************************************
751 SUBROUTINE cp_sirius_print_results(pwdft_env, print_section)
752 TYPE(pwdft_environment_type), INTENT(INOUT), &
753 POINTER :: pwdft_env
754 TYPE(section_vals_type), POINTER :: print_section
755
756 CHARACTER(LEN=default_string_length) :: my_act, my_pos
757 INTEGER :: i, ik, iounit, ispn, iterstep, iv, iw, &
758 nbands, nhist, nkpts, nspins
759 INTEGER(KIND=C_INT) :: cint
760 LOGICAL :: append, dos, ionode
761 REAL(KIND=c_double) :: creal
762 REAL(KIND=c_double), ALLOCATABLE, DIMENSION(:) :: slist
763 REAL(KIND=dp) :: de, e_fermi(2), emax, emin, eval
764 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: wkpt
765 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: ehist, hist, occval
766 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: energies, occupations
767 TYPE(cp_logger_type), POINTER :: logger
768 TYPE(sirius_context_handler) :: sctx
769 TYPE(sirius_ground_state_handler) :: gs_handler
770 TYPE(sirius_kpoint_set_handler) :: ks_handler
771
772 NULLIFY (logger)
773 logger => cp_get_default_logger()
774 ionode = logger%para_env%is_source()
775 iounit = cp_logger_get_default_io_unit(logger)
776
777 ! Density of States
778 dos = btest(cp_print_key_should_output(logger%iter_info, print_section, "DOS"), cp_p_file)
779 IF (dos) THEN
780 CALL pwdft_env_get(pwdft_env, ks_handler=ks_handler)
781 CALL pwdft_env_get(pwdft_env, gs_handler=gs_handler)
782 CALL pwdft_env_get(pwdft_env, sctx=sctx)
783
784 CALL section_vals_val_get(print_section, "DOS%DELTA_E", r_val=de)
785 CALL section_vals_val_get(print_section, "DOS%APPEND", l_val=append)
786
787 CALL sirius_get_num_kpoints(ks_handler, cint)
788 nkpts = cint
789 CALL sirius_get_parameters(sctx, num_bands=cint)
790 nbands = cint
791 CALL sirius_get_parameters(sctx, num_spins=cint)
792 nspins = cint
793 e_fermi(:) = 0.0_dp
794 ALLOCATE (energies(nbands, nspins, nkpts))
795 energies = 0.0_dp
796 ALLOCATE (occupations(nbands, nspins, nkpts))
797 occupations = 0.0_dp
798 ALLOCATE (wkpt(nkpts))
799 ALLOCATE (slist(nbands))
800 DO ik = 1, nkpts
801 CALL sirius_get_kpoint_properties(ks_handler, ik, creal)
802 wkpt(ik) = creal
803 END DO
804 DO ik = 1, nkpts
805 DO ispn = 1, nspins
806 CALL sirius_get_band_energies(ks_handler, ik, ispn, slist)
807 energies(1:nbands, ispn, ik) = slist(1:nbands)
808 CALL sirius_get_band_occupancies(ks_handler, ik, ispn, slist)
809 occupations(1:nbands, ispn, ik) = slist(1:nbands)
810 END DO
811 END DO
812 emin = minval(energies)
813 emax = maxval(energies)
814 nhist = nint((emax - emin)/de) + 1
815 ALLOCATE (hist(nhist, nspins), occval(nhist, nspins), ehist(nhist, nspins))
816 hist = 0.0_dp
817 occval = 0.0_dp
818 ehist = 0.0_dp
819
820 DO ik = 1, nkpts
821 DO ispn = 1, nspins
822 DO i = 1, nbands
823 eval = energies(i, ispn, ik) - emin
824 iv = nint(eval/de) + 1
825 cpassert((iv > 0) .AND. (iv <= nhist))
826 hist(iv, ispn) = hist(iv, ispn) + wkpt(ik)
827 occval(iv, ispn) = occval(iv, ispn) + wkpt(ik)*occupations(i, ispn, ik)
828 END DO
829 END DO
830 END DO
831 hist = hist/real(nbands, kind=dp)
832 DO i = 1, nhist
833 ehist(i, 1:nspins) = emin + (i - 1)*de
834 END DO
835
836 iterstep = logger%iter_info%iteration(logger%iter_info%n_rlevel)
837 my_act = "WRITE"
838 IF (append .AND. iterstep > 1) THEN
839 my_pos = "APPEND"
840 ELSE
841 my_pos = "REWIND"
842 END IF
843
844 iw = cp_print_key_unit_nr(logger, print_section, "DOS", &
845 extension=".dos", file_position=my_pos, file_action=my_act, &
846 file_form="FORMATTED")
847 IF (iw > 0) THEN
848 IF (nspins == 2) THEN
849 WRITE (unit=iw, fmt="(T2,A,I0,A,2F12.6)") &
850 "# DOS at iteration step i = ", iterstep, ", E_Fermi[a.u.] = ", e_fermi(1:2)
851 WRITE (unit=iw, fmt="(T2,A, A)") " Energy[a.u.] Alpha_Density Occupation", &
852 " Beta_Density Occupation"
853 ELSE
854 WRITE (unit=iw, fmt="(T2,A,I0,A,F12.6)") &
855 "# DOS at iteration step i = ", iterstep, ", E_Fermi[a.u.] = ", e_fermi(1)
856 WRITE (unit=iw, fmt="(T2,A)") " Energy[a.u.] Density Occupation"
857 END IF
858 DO i = 1, nhist
859 eval = emin + (i - 1)*de
860 IF (nspins == 2) THEN
861 WRITE (unit=iw, fmt="(F15.8,4F15.4)") eval, hist(i, 1), occval(i, 1), &
862 hist(i, 2), occval(i, 2)
863 ELSE
864 WRITE (unit=iw, fmt="(F15.8,2F15.4)") eval, hist(i, 1), occval(i, 1)
865 END IF
866 END DO
867 END IF
868 CALL cp_print_key_finished_output(iw, logger, print_section, "DOS")
869
870 DEALLOCATE (energies, occupations, wkpt, slist)
871 DEALLOCATE (hist, occval, ehist)
872 END IF
873 END SUBROUTINE cp_sirius_print_results
874
875END MODULE sirius_interface
876
877#else
878
879!***************************************************************************************************
880!> \brief Empty implementation in case SIRIUS is not compiled in.
881!***************************************************************************************************
884#include "./base/base_uses.f90"
885
886 IMPLICIT NONE
887 PRIVATE
888
891
892CONTAINS
893
894! **************************************************************************************************
895!> \brief Empty implementation in case SIRIUS is not compiled in.
896! **************************************************************************************************
897 SUBROUTINE cp_sirius_init()
898 END SUBROUTINE cp_sirius_init
899
900! **************************************************************************************************
901!> \brief Empty implementation in case SIRIUS is not compiled in.
902! **************************************************************************************************
904 END SUBROUTINE cp_sirius_finalize
905
906! **************************************************************************************************
907!> \brief Empty implementation in case SIRIUS is not compiled in.
908!> \param pwdft_env ...
909! **************************************************************************************************
910 SUBROUTINE cp_sirius_create_env(pwdft_env)
911 TYPE(pwdft_environment_type), POINTER :: pwdft_env
912
913 mark_used(pwdft_env)
914 cpabort("Sirius library is missing")
915 END SUBROUTINE cp_sirius_create_env
916
917! **************************************************************************************************
918!> \brief Empty implementation in case SIRIUS is not compiled in.
919!> \param pwdft_env ...
920!> \param calculate_forces ...
921!> \param calculate_stress ...
922! **************************************************************************************************
923 SUBROUTINE cp_sirius_energy_force(pwdft_env, calculate_forces, calculate_stress)
924 TYPE(pwdft_environment_type), POINTER :: pwdft_env
925 LOGICAL :: calculate_forces, calculate_stress
926
927 mark_used(pwdft_env)
928 mark_used(calculate_forces)
929 mark_used(calculate_stress)
930 cpabort("Sirius library is missing")
931 END SUBROUTINE cp_sirius_energy_force
932
933! **************************************************************************************************
934!> \brief Empty implementation in case SIRIUS is not compiled in.
935!> \param pwdft_env ...
936! **************************************************************************************************
937 SUBROUTINE cp_sirius_update_context(pwdft_env)
938 TYPE(pwdft_environment_type), POINTER :: pwdft_env
939
940 mark_used(pwdft_env)
941 cpabort("Sirius library is missing")
942 END SUBROUTINE cp_sirius_update_context
943
944END MODULE sirius_interface
945
946#endif
calculate the orbitals for a given atomic kind type
subroutine, public calculate_atomic_orbitals(atomic_kind, qs_kind, agrid, iunit, pmat, fmat, density, wavefunction, wfninfo, confine, xc_section, nocc)
...
subroutine, public gth_potential_conversion(gth_potential, gth_atompot)
...
Define the atom type and its sub types.
Definition atom_types.F:15
Routines that process Quantum Espresso UPF files.
Definition atom_upf.F:14
Some basic routines for atomic calculations.
Definition atom_utils.F:15
pure subroutine, public atom_local_potential(locpot, gthpot, rr)
...
Definition atom_utils.F:800
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.
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
Definition cell_types.F:486
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
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...
Definition of the atomic potential types.
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_gapw_log
integer, parameter, public sirius_func_vdwdfcx
integer, parameter, public sirius_func_vdwdf2
integer, parameter, public sirius_func_vdwdf
objects that represent the structure of input sections and the data contained in an input section
type(section_vals_type) function, pointer, public section_vals_get_subs_vals2(section_vals, i_section, i_rep_section)
returns the values of the n-th non default subsection (null if no such section exists (not so many no...
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
Definition of mathematical constants and functions.
real(kind=dp), dimension(0:maxfac), parameter, public gamma1
real(kind=dp), parameter, public fourpi
Interface to the message passing library MPI.
Define the data structure for the particle information.
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public massunit
Definition physcon.F:141
The type definitions for the PWDFT environment.
subroutine, public pwdft_env_get(pwdft_env, pwdft_input, force_env_input, xc_input, cp_subsys, qs_subsys, para_env, energy, forces, stress, sctx, gs_handler, ks_handler)
Returns various attributes of the pwdft environment.
subroutine, public pwdft_env_set(pwdft_env, pwdft_input, force_env_input, xc_input, qs_subsys, cp_subsys, para_env, energy, forces, stress, sctx, gs_handler, ks_handler)
Sets various attributes of the pwdft environment.
subroutine, public deallocate_grid_atom(grid_atom)
Deallocate a Gaussian-type orbital (GTO) basis set data set.
subroutine, public allocate_grid_atom(grid_atom)
Initialize components of the grid_atom_type structure.
subroutine, public create_grid_atom(grid_atom, nr, na, llmax, ll, quadrature)
...
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, 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_r3d_rs_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_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
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)
...
Interface to the SIRIUS Library.
subroutine, public cp_sirius_update_context(pwdft_env)
Empty implementation in case SIRIUS is not compiled in.
subroutine, public cp_sirius_init()
Empty implementation in case SIRIUS is not compiled in.
subroutine, public cp_sirius_energy_force(pwdft_env, calculate_forces, calculate_stress)
Empty implementation in case SIRIUS is not compiled in.
subroutine, public cp_sirius_finalize()
Empty implementation in case SIRIUS is not compiled in.
subroutine, public cp_sirius_create_env(pwdft_env)
Empty implementation in case SIRIUS is not compiled in.
Provides all information about a pseudopotential.
Definition atom_types.F:98
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
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
Provides all information about a quickstep kind.