(git:374b731)
Loading...
Searching...
No Matches
preconditioner.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 computes preconditioners, and implements methods to apply them
10!> currently used in qs_ot
11!> \par History
12!> - [UB] 2009-05-13 Adding stable approximate inverse (full and sparse)
13!> \author Joost VandeVondele (09.2002)
14! **************************************************************************************************
22 USE cp_fm_types, ONLY: cp_fm_create,&
27 USE dbcsr_api, ONLY: dbcsr_p_type,&
28 dbcsr_type
37 USE kinds, ONLY: default_string_length,&
38 dp
53 USE qs_mo_types, ONLY: get_mo_set,&
56#include "./base/base_uses.f90"
57
58 IMPLICIT NONE
59
60 PRIVATE
61
62 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'preconditioner'
63
66
67! The public interface for apply preconditioner, the routines can be found in preconditioner_apply.F
69 MODULE PROCEDURE apply_preconditioner_dbcsr
70 MODULE PROCEDURE apply_preconditioner_fm
71 END INTERFACE
72
73! **************************************************************************************************
74
75CONTAINS
76
77! **************************************************************************************************
78
79! creates a preconditioner for the system (H-energy_homo S)
80! this preconditioner is (must be) symmetric positive definite.
81! currently uses a atom-block-diagonal form
82! each block will be ....
83! might overwrite matrix_h, matrix_t
84
85! **************************************************************************************************
86!> \brief ...
87!> \param preconditioner_env ...
88!> \param precon_type ...
89!> \param solver_type ...
90!> \param matrix_h ...
91!> \param matrix_s ...
92!> \param matrix_t ...
93!> \param mo_set ...
94!> \param energy_gap ...
95!> \param convert_precond_to_dbcsr ...
96!> \param chol_type ...
97!> \par History
98!> 09.2014 removed some unused or unfinished methods
99!> removed sparse preconditioners and the
100!> sparse approximate inverse at rev 14341 [Florian Schiffmann]
101! **************************************************************************************************
102 SUBROUTINE make_preconditioner(preconditioner_env, precon_type, solver_type, matrix_h, matrix_s, &
103 matrix_t, mo_set, energy_gap, convert_precond_to_dbcsr, chol_type)
104
105 TYPE(preconditioner_type) :: preconditioner_env
106 INTEGER, INTENT(IN) :: precon_type, solver_type
107 TYPE(dbcsr_type), POINTER :: matrix_h
108 TYPE(dbcsr_type), OPTIONAL, POINTER :: matrix_s, matrix_t
109 TYPE(mo_set_type), INTENT(IN) :: mo_set
110 REAL(kind=dp) :: energy_gap
111 LOGICAL, INTENT(IN), OPTIONAL :: convert_precond_to_dbcsr
112 INTEGER, INTENT(IN), OPTIONAL :: chol_type
113
114 CHARACTER(len=*), PARAMETER :: routinen = 'make_preconditioner'
115
116 INTEGER :: handle, k, my_solver_type, nao, nhomo
117 LOGICAL :: my_convert_precond_to_dbcsr, &
118 needs_full_spectrum, needs_homo, &
119 use_mo_coeff_b
120 REAL(kind=dp) :: energy_homo
121 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: eigenvalues_ot
122 TYPE(cp_fm_struct_type), POINTER :: fm_struct
123 TYPE(cp_fm_type) :: mo_occ
124 TYPE(cp_fm_type), POINTER :: mo_coeff
125 TYPE(dbcsr_type), POINTER :: mo_coeff_b
126
127 CALL timeset(routinen, handle)
128
129 CALL get_mo_set(mo_set=mo_set, mo_coeff=mo_coeff, mo_coeff_b=mo_coeff_b, homo=nhomo)
130 use_mo_coeff_b = mo_set%use_mo_coeff_b
131 CALL cp_fm_get_info(mo_coeff, ncol_global=k, nrow_global=nao)
132
133 ! Starting some matrix mess, check where to store the result in preconditioner_env, fm or dbcsr_matrix
134 my_convert_precond_to_dbcsr = .false.
135 IF (PRESENT(convert_precond_to_dbcsr)) my_convert_precond_to_dbcsr = convert_precond_to_dbcsr
136
137 ! Thanks to the mess with the matrices we need to make sure in this case that the
138 ! Previous inverse is properly stored as a sparse matrix, fm gets deallocated here
139 ! if it wasn't anyway
140 IF (preconditioner_env%solver == ot_precond_solver_update) &
141 CALL transfer_fm_to_dbcsr(preconditioner_env%fm, preconditioner_env%dbcsr_matrix, matrix_h)
142
143 needs_full_spectrum = .false.
144 needs_homo = .false.
145
146 SELECT CASE (precon_type)
148 needs_full_spectrum = .true.
149 ! both of them need the coefficients as fm's, more matrix mess
150 IF (use_mo_coeff_b) THEN
151 CALL copy_dbcsr_to_fm(mo_coeff_b, mo_coeff)
152 END IF
154 needs_homo = .true.
155 ! XXXX to be removed if homo estimate only is implemented
156 needs_full_spectrum = .true.
158 ! these should be happy without an estimate for the homo energy
159 ! preconditioning can not depend on an absolute eigenvalue, only on eigenvalue differences
160 CASE DEFAULT
161 cpabort("The preconditioner is unknown ...")
162 END SELECT
163
164 ALLOCATE (eigenvalues_ot(k))
165 energy_homo = 0.0_dp
166 IF (needs_full_spectrum) THEN
167 ! XXXXXXXXXXXXXXXX do not touch the initial MOs, could be harmful for either
168 ! the case of non-equivalent MOs but also for the derivate
169 ! we could already have all eigenvalues e.g. full_all and we could skip this
170 ! to be optimised later.
171 ! one flaw is that not all SCF methods (i.e. that go over mo_derivs directly)
172 ! have a 'valid' matrix_h... (we even don't know what evals are in that case)
173 IF (use_mo_coeff_b) THEN
174 CALL calculate_subspace_eigenvalues(mo_coeff_b, matrix_h, &
175 eigenvalues_ot, do_rotation=.false., &
176 para_env=mo_coeff%matrix_struct%para_env, &
177 blacs_env=mo_coeff%matrix_struct%context)
178 ELSE
179 CALL calculate_subspace_eigenvalues(mo_coeff, matrix_h, &
180 eigenvalues_ot, do_rotation=.false.)
181 END IF
182 IF (k > 0) THEN
183 cpassert(nhomo > 0 .AND. nhomo <= k)
184 energy_homo = eigenvalues_ot(nhomo)
185 END IF
186 ELSE
187 IF (needs_homo) THEN
188 cpabort("Not yet implemented")
189 END IF
190 END IF
191
192 ! After all bits and pieces of checking and initialization, here comes the
193 ! part where the preconditioner matrix gets created and solved.
194 ! This will give the matrices for later use
195 my_solver_type = solver_type
196 preconditioner_env%in_use = precon_type
197 preconditioner_env%cholesky_use = cholesky_reduce
198 IF (PRESENT(chol_type)) preconditioner_env%cholesky_use = chol_type
199 preconditioner_env%in_use = precon_type
200 IF (nhomo == k) THEN
201 CALL make_preconditioner_matrix(preconditioner_env, matrix_h, matrix_s, matrix_t, mo_coeff, &
202 energy_homo, eigenvalues_ot, energy_gap, my_solver_type)
203 ELSE
204 CALL cp_fm_struct_create(fm_struct, nrow_global=nao, ncol_global=nhomo, &
205 context=preconditioner_env%ctxt, &
206 para_env=preconditioner_env%para_env)
207 CALL cp_fm_create(mo_occ, fm_struct)
208 CALL cp_fm_to_fm(mo_coeff, mo_occ, nhomo)
209 CALL cp_fm_struct_release(fm_struct)
210 !
211 CALL make_preconditioner_matrix(preconditioner_env, matrix_h, matrix_s, matrix_t, mo_occ, &
212 energy_homo, eigenvalues_ot(1:nhomo), energy_gap, my_solver_type)
213 !
214 CALL cp_fm_release(mo_occ)
215 END IF
216
217 CALL solve_preconditioner(my_solver_type, preconditioner_env, matrix_s, matrix_h)
218
219 ! Here comes more matrix mess, make sure to output the correct matrix format,
220 ! A bit pointless to convert the cholesky factorized version as it doesn't work in
221 ! dbcsr form and will crash later,...
222 IF (my_convert_precond_to_dbcsr) THEN
223 CALL transfer_fm_to_dbcsr(preconditioner_env%fm, preconditioner_env%dbcsr_matrix, matrix_h)
224 ELSE
225 CALL transfer_dbcsr_to_fm(preconditioner_env%dbcsr_matrix, preconditioner_env%fm, &
226 preconditioner_env%para_env, preconditioner_env%ctxt)
227 END IF
228
229 DEALLOCATE (eigenvalues_ot)
230
231 CALL timestop(handle)
232
233 END SUBROUTINE make_preconditioner
234
235! **************************************************************************************************
236!> \brief Allows for a restart of the preconditioner
237!> depending on the method it purges all arrays or keeps them
238!> \param qs_env ...
239!> \param preconditioner ...
240!> \param prec_type ...
241!> \param nspins ...
242! **************************************************************************************************
243 SUBROUTINE restart_preconditioner(qs_env, preconditioner, prec_type, nspins)
244
245 TYPE(qs_environment_type), POINTER :: qs_env
246 TYPE(preconditioner_p_type), DIMENSION(:), POINTER :: preconditioner
247 INTEGER, INTENT(IN) :: prec_type, nspins
248
249 INTEGER :: ispin
250 TYPE(cp_blacs_env_type), POINTER :: blacs_env
251 TYPE(mp_para_env_type), POINTER :: para_env
252
253 NULLIFY (para_env, blacs_env)
254 CALL get_qs_env(qs_env, para_env=para_env, blacs_env=blacs_env)
255
256 IF (ASSOCIATED(preconditioner)) THEN
257 SELECT CASE (prec_type)
258 CASE (ot_precond_full_all, ot_precond_full_single) ! these depend on the ks matrix
259 DO ispin = 1, SIZE(preconditioner)
260 CALL destroy_preconditioner(preconditioner(ispin)%preconditioner)
261 DEALLOCATE (preconditioner(ispin)%preconditioner)
262 END DO
263 DEALLOCATE (preconditioner)
265 ! do nothing
266 CASE DEFAULT
267 cpabort("")
268 END SELECT
269 END IF
270
271 ! add an OT preconditioner if none is present
272 IF (.NOT. ASSOCIATED(preconditioner)) THEN
273 SELECT CASE (prec_type)
275 ALLOCATE (preconditioner(nspins))
276 CASE DEFAULT
277 ALLOCATE (preconditioner(1))
278 END SELECT
279 DO ispin = 1, SIZE(preconditioner)
280 ALLOCATE (preconditioner(ispin)%preconditioner)
281 CALL init_preconditioner(preconditioner(ispin)%preconditioner, &
282 para_env=para_env, &
283 blacs_env=blacs_env)
284 END DO
285 END IF
286
287 END SUBROUTINE restart_preconditioner
288
289! **************************************************************************************************
290!> \brief ...
291!> \param qs_env ...
292!> \param mos ...
293!> \param matrix_ks ...
294!> \param matrix_s ...
295!> \param ot_preconditioner ...
296!> \param prec_type ...
297!> \param solver_type ...
298!> \param energy_gap ...
299!> \param nspins ...
300!> \param has_unit_metric ...
301!> \param convert_to_dbcsr ...
302!> \param chol_type ...
303!> \param full_mo_set ...
304! **************************************************************************************************
305 SUBROUTINE prepare_preconditioner(qs_env, mos, matrix_ks, matrix_s, &
306 ot_preconditioner, prec_type, solver_type, &
307 energy_gap, nspins, has_unit_metric, &
308 convert_to_dbcsr, chol_type, full_mo_set)
309
310 TYPE(qs_environment_type), POINTER :: qs_env
311 TYPE(mo_set_type), DIMENSION(:), INTENT(INOUT) :: mos
312 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_ks, matrix_s
313 TYPE(preconditioner_p_type), DIMENSION(:), POINTER :: ot_preconditioner
314 INTEGER, INTENT(IN) :: prec_type, solver_type
315 REAL(dp), INTENT(IN) :: energy_gap
316 INTEGER, INTENT(IN) :: nspins
317 LOGICAL, INTENT(IN), OPTIONAL :: has_unit_metric, convert_to_dbcsr
318 INTEGER, INTENT(IN), OPTIONAL :: chol_type
319 LOGICAL, INTENT(IN), OPTIONAL :: full_mo_set
320
321 CHARACTER(LEN=*), PARAMETER :: routinen = 'prepare_preconditioner'
322
323 CHARACTER(LEN=default_string_length) :: msg
324 INTEGER :: handle, icall, ispin, n_loops
325 INTEGER, DIMENSION(5) :: nocc, norb
326 LOGICAL :: do_co_rotate, my_convert_to_dbcsr, &
327 my_full_mo_set, my_has_unit_metric, &
328 use_mo_coeff_b
329 TYPE(cp_blacs_env_type), POINTER :: blacs_env
330 TYPE(cp_fm_type), POINTER :: mo_coeff
331 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: kinetic
332 TYPE(dbcsr_type), POINTER :: matrix_t, mo_coeff_b
333 TYPE(dft_control_type), POINTER :: dft_control
334 TYPE(mp_para_env_type), POINTER :: para_env
335
336 CALL timeset(routinen, handle)
337 NULLIFY (matrix_t, mo_coeff_b, mo_coeff, kinetic, dft_control, para_env, blacs_env)
338 my_has_unit_metric = .false.
339 IF (PRESENT(has_unit_metric)) my_has_unit_metric = has_unit_metric
340 my_convert_to_dbcsr = .true.
341 IF (PRESENT(convert_to_dbcsr)) my_convert_to_dbcsr = convert_to_dbcsr
342 my_full_mo_set = .false.
343 IF (PRESENT(full_mo_set)) my_full_mo_set = full_mo_set
344
345 CALL get_qs_env(qs_env, &
346 dft_control=dft_control, &
347 para_env=para_env, &
348 blacs_env=blacs_env)
349
350 IF (dft_control%qs_control%semi_empirical .OR. dft_control%qs_control%dftb .OR. &
351 dft_control%qs_control%xtb) THEN
352 IF (prec_type == ot_precond_full_kinetic) THEN
353 msg = "Full_kinetic not available for semi-empirical methods"
354 cpabort(trim(msg))
355 END IF
356 matrix_t => matrix_s(1)%matrix
357 ELSE
358 cpassert(.NOT. my_has_unit_metric)
359 CALL get_qs_env(qs_env, kinetic=kinetic)
360 matrix_t => kinetic(1)%matrix
361 END IF
362
363 ! use full set of MOs or just occupied MOs
364 nocc = 0
365 norb = 0
366 IF (my_full_mo_set) THEN
367 DO ispin = 1, nspins
368 CALL get_mo_set(mo_set=mos(ispin), homo=nocc(ispin), nmo=norb(ispin))
369 CALL set_mo_set(mo_set=mos(ispin), homo=norb(ispin))
370 END DO
371 END IF
372 !determines how often make preconditioner is called, spin dependent methods have to be called twice
373 n_loops = 1
374 IF (prec_type == ot_precond_full_single_inverse) n_loops = nspins
375 ! check whether we need the ev and rotate the MOs
376 SELECT CASE (prec_type)
378 ! if one of these preconditioners is used every spin needs to call make_preconditioner
379 n_loops = nspins
380
381 do_co_rotate = ASSOCIATED(qs_env%mo_derivs)
382 DO ispin = 1, nspins
383 CALL get_mo_set(mo_set=mos(ispin), mo_coeff_b=mo_coeff_b, mo_coeff=mo_coeff)
384 use_mo_coeff_b = mos(ispin)%use_mo_coeff_b
385 IF (use_mo_coeff_b .AND. do_co_rotate) THEN
386 CALL calculate_subspace_eigenvalues(mo_coeff_b, matrix_ks(ispin)%matrix, &
387 do_rotation=.true., &
388 co_rotate=qs_env%mo_derivs(ispin)%matrix, &
389 para_env=para_env, &
390 blacs_env=blacs_env)
391 ELSEIF (use_mo_coeff_b) THEN
392 CALL calculate_subspace_eigenvalues(mo_coeff_b, matrix_ks(ispin)%matrix, &
393 do_rotation=.true., &
394 para_env=para_env, &
395 blacs_env=blacs_env)
396 ELSE
397 CALL calculate_subspace_eigenvalues(mo_coeff, matrix_ks(ispin)%matrix, &
398 do_rotation=.true.)
399 END IF
400 END DO
401 CASE DEFAULT
402 ! No need to rotate the MOs
403 END SELECT
404
405 ! check whether we have a preconditioner
406 SELECT CASE (prec_type)
407 CASE (ot_precond_none)
408 DO ispin = 1, SIZE(ot_preconditioner)
409 ot_preconditioner(ispin)%preconditioner%in_use = 0
410 END DO
411 CASE DEFAULT
412 DO icall = 1, n_loops
413 IF (my_has_unit_metric) THEN
414 CALL make_preconditioner(ot_preconditioner(icall)%preconditioner, &
415 prec_type, &
416 solver_type, &
417 matrix_h=matrix_ks(icall)%matrix, &
418 mo_set=mos(icall), &
419 energy_gap=energy_gap, &
420 convert_precond_to_dbcsr=my_convert_to_dbcsr)
421 ELSE
422 CALL make_preconditioner(ot_preconditioner(icall)%preconditioner, &
423 prec_type, &
424 solver_type, &
425 matrix_h=matrix_ks(icall)%matrix, &
426 matrix_s=matrix_s(1)%matrix, &
427 matrix_t=matrix_t, &
428 mo_set=mos(icall), &
429 energy_gap=energy_gap, &
430 convert_precond_to_dbcsr=my_convert_to_dbcsr, chol_type=chol_type)
431 END IF
432 END DO
433 END SELECT
434
435 ! reset homo values
436 IF (my_full_mo_set) THEN
437 DO ispin = 1, nspins
438 CALL set_mo_set(mo_set=mos(ispin), homo=nocc(ispin))
439 END DO
440 END IF
441
442 CALL timestop(handle)
443
444 END SUBROUTINE prepare_preconditioner
445
446END MODULE preconditioner
447
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
DBCSR operations in CP2K.
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
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
Definition cp_fm_types.F:15
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)
creates a new full matrix with the given structure
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public ot_precond_full_kinetic
integer, parameter, public cholesky_reduce
integer, parameter, public ot_precond_full_single
integer, parameter, public ot_precond_none
integer, parameter, public ot_precond_full_single_inverse
integer, parameter, public ot_precond_s_inverse
integer, parameter, public ot_precond_solver_update
integer, parameter, public ot_precond_full_all
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
Interface to the message passing library MPI.
computes preconditioners, and implements methods to apply them currently used in qs_ot
subroutine, public apply_preconditioner_fm(preconditioner_env, matrix_in, matrix_out)
applies a previously created preconditioner to a full matrix
subroutine, public apply_preconditioner_dbcsr(preconditioner_env, matrix_in, matrix_out)
...
computes preconditioners, and implements methods to apply them currently used in qs_ot
subroutine, public make_preconditioner_matrix(preconditioner_env, matrix_h, matrix_s, matrix_t, mo_coeff, energy_homo, eigenvalues_ot, energy_gap, my_solver_type)
...
solves the preconditioner, contains to utility function for fm<->dbcsr transfers, should be moved soo...
subroutine, public transfer_dbcsr_to_fm(dbcsr_matrix, fm_matrix, para_env, context)
transfers a dbcsr to a full matrix
subroutine, public solve_preconditioner(my_solver_type, preconditioner_env, matrix_s, matrix_h)
...
subroutine, public transfer_fm_to_dbcsr(fm_matrix, dbcsr_matrix, template_mat)
transfers a full matrix to a dbcsr
types of preconditioners
subroutine, public init_preconditioner(preconditioner_env, para_env, blacs_env)
...
subroutine, public destroy_preconditioner(preconditioner_env)
...
computes preconditioners, and implements methods to apply them currently used in qs_ot
subroutine, public restart_preconditioner(qs_env, preconditioner, prec_type, nspins)
Allows for a restart of the preconditioner depending on the method it purges all arrays or keeps them...
subroutine, public prepare_preconditioner(qs_env, mos, matrix_ks, matrix_s, ot_preconditioner, prec_type, solver_type, energy_gap, nspins, has_unit_metric, convert_to_dbcsr, chol_type, full_mo_set)
...
subroutine, public make_preconditioner(preconditioner_env, precon_type, solver_type, matrix_h, matrix_s, matrix_t, mo_set, energy_gap, convert_precond_to_dbcsr, chol_type)
...
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_nonbond, sab_almo, sab_kp, sab_kp_nosym, 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, ecoul_1c, rho0_s_rs, rho0_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, 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, rhs)
Get the QUICKSTEP environment.
collects routines that perform operations directly related to MOs
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
subroutine, public set_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, uniform_occupation, kts, mu, flexible_electron_count)
Set the components of a MO set data structure.
subroutine, public get_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, mo_coeff, mo_coeff_b, uniform_occupation, kts, mu, flexible_electron_count)
Get the components of a MO set data structure.
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
keeps the information about the structure of a full matrix
represent a full matrix
stores all the informations relevant to an mpi environment