35 use,
intrinsic :: iso_fortran_env
69 integer(int64) :: type
72 real(real64) :: line_tol
73 real(real64) :: fire_mass
74 integer :: fire_integrator
75 real(real64) :: tolgrad
78 integer :: what2minimize
82 type(ions_t),
pointer :: ions
83 type(hamiltonian_elec_t),
pointer :: hm
84 type(electrons_t),
pointer :: syst
85 class(mesh_t),
pointer :: mesh
86 type(states_elec_t),
pointer :: st
88 integer :: periodic_dim
90 integer :: fixed_atom = 0
92 real(real64),
allocatable :: cell_force(:, :)
93 logical :: symmetrize = .false.
94 real(real64),
allocatable :: initial_length(:)
95 real(real64),
allocatable :: initial_rlattice(:, :)
96 real(real64),
allocatable :: inv_initial_rlattice(:, :)
97 real(real64) :: pressure
99 logical :: poscar_output = .false.
103 type(geom_opt_t),
save :: g_opt
105 integer,
parameter :: &
106 MINWHAT_ENERGY = 1, &
109 integer,
parameter :: &
118 class(*),
intent(inout) :: system
119 logical,
intent(inout) :: from_scratch
125 message(1) =
"CalculationMode = go not implemented for multi-system calculations"
136 type(electrons_t),
target,
intent(inout) :: sys
137 logical,
intent(inout) :: fromscratch
140 real(real64),
allocatable :: coords(:)
141 real(real64) :: energy
143 real(real64),
allocatable :: mass(:)
144 integer :: iatom, imass
145 type(restart_t) :: restart_load
149 if (sys%space%periodic_dim == 1 .or. sys%space%periodic_dim == 2)
then
150 message(1) =
"Geometry optimization not allowed for systems periodic in 1D and 2D, "
151 message(2) =
"as in those cases the total energy is not correct."
156 if (sys%hm%pcm%run_pcm)
then
160 if (sys%kpoints%use_symmetries)
then
167 if (.not. fromscratch)
then
170 call states_elec_load(restart_load, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, ierr)
172 call restart_load%end()
174 message(1) =
"Unable to read wavefunctions: Starting from scratch."
180 call scf_init(g_opt%scfv, sys%namespace, sys%gr, sys%ions, sys%st, sys%mc, sys%hm, sys%space)
183 if (.not. g_opt%scfv%calc_stress)
then
184 message(1) =
"In order to optimize the cell, one needs to set SCFCalculateStress = yes."
189 if (fromscratch)
then
190 call lcao_run(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm, lmm_r = g_opt%scfv%lmm_r)
193 message(1) =
'Info: Setting up Hamiltonian.'
195 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
198 g_opt%symmetrize = sys%kpoints%use_symmetries .or. sys%st%symmetrize_density
201 safe_allocate(coords(1:g_opt%size))
204 if (sys%st%pack_states .and. sys%hm%apply_packed())
call sys%st%pack()
207 select case (g_opt%method)
209 call minimize_multidim_nograd(g_opt%method, g_opt%size, coords, g_opt%step,&
210 g_opt%toldr, g_opt%max_iter, &
214 safe_allocate(mass(1:g_opt%size))
215 mass = g_opt%fire_mass
217 do iatom = 1, sys%ions%natoms
218 if (g_opt%fixed_atom == iatom) cycle
219 if (g_opt%ions%fixed(iatom)) cycle
220 if (g_opt%fire_mass <=
m_zero) mass(imass:imass + 2) = sys%ions%mass(iatom)
225 call minimize_fire(g_opt%size, g_opt%ions%space%dim, coords, g_opt%step, g_opt%tolgrad, &
227 safe_deallocate_a(mass)
230 call minimize_multidim(g_opt%method, g_opt%size, coords, g_opt%step ,&
231 g_opt%line_tol , g_opt%tolgrad, g_opt%toldr, g_opt%max_iter, &
235 if (ierr == 1025)
then
237 message(1) =
"Reached maximum number of iterations allowed by GOMaxIter."
240 message(1) =
"Error occurred during the GSL minimization procedure:"
245 if (sys%st%pack_states .and. sys%hm%apply_packed())
call sys%st%unpack()
249 message(1) =
"Writing final coordinates to min.xyz"
252 call g_opt%ions%write_xyz(
'./min')
254 safe_deallocate_a(coords)
257 call g_opt%scfv%criterion_list%empty()
264 subroutine init_(fromscratch)
265 logical,
intent(inout) :: fromscratch
267 logical :: center, does_exist
268 integer :: iter, iatom, idir
269 character(len=100) :: filename
270 real(real64) :: default_toldr
271 real(real64) :: default_step
276 if (sys%space%is_periodic())
then
304 write(
message(1),
'(a)')
'Input: [GOType = '
305 if (
bitand(g_opt%type, go_ions) /= 0)
then
309 if (len_trim(
message(1)) > 16)
then
315 if (len_trim(
message(1)) > 16)
then
324 message(1) =
"Cell and volume optimization cannot be used simultaneously."
338 message(1) =
"Cell dynamics not supported on GPUs."
344 do iatom = 1, sys%ions%natoms
345 select type(spec=>sys%ions%atom(iatom)%species)
347 write(
message(1),
'(a)')
"Geometry optimization for all-electron potential is not implemented."
357 g_opt%ions => sys%ions
361 g_opt%dim = sys%space%dim
362 g_opt%periodic_dim = sys%space%periodic_dim
366 if (
bitand(g_opt%type, go_ions) /= 0)
then
367 g_opt%size = g_opt%dim * g_opt%ions%natoms
372 g_opt%size = g_opt%size + (g_opt%periodic_dim +1) * g_opt%periodic_dim / 2
373 safe_allocate(g_opt%cell_force(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
378 g_opt%size = g_opt%size + g_opt%periodic_dim
379 safe_allocate(g_opt%cell_force(1:g_opt%periodic_dim, 1:1))
381 safe_allocate(g_opt%initial_length(1:g_opt%periodic_dim))
382 do idir = 1, g_opt%periodic_dim
383 g_opt%initial_length(idir) = norm2(g_opt%ions%latt%rlattice(1:g_opt%periodic_dim, idir))
389 safe_allocate(g_opt%initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
390 g_opt%initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim) &
391 = g_opt%ions%latt%rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
392 safe_allocate(g_opt%inv_initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
393 g_opt%inv_initial_rlattice(:, :) = g_opt%initial_rlattice(:, :)
394 call lalg_inverse(g_opt%periodic_dim, g_opt%inv_initial_rlattice,
'dir')
397 if(g_opt%ions%space%is_periodic())
then
401 assert(g_opt%size > 0)
417 g_opt%size = g_opt%size - g_opt%dim
422 do iatom = 1, g_opt%ions%natoms
423 if (g_opt%ions%fixed(iatom))
then
424 g_opt%size = g_opt%size - g_opt%dim
502 default_toldr = 0.001_real64
504 default_toldr = -
m_one
521 call parse_variable(sys%namespace,
'GOStep', default_step, g_opt%step)
536 call parse_variable(sys%namespace,
'GOLineTol', 0.1_real64, g_opt%line_tol)
546 call parse_variable(sys%namespace,
'GOMaxIter', 200, g_opt%max_iter)
547 if (g_opt%max_iter <= 0)
then
548 message(1) =
"GOMaxIter has to be larger than 0"
585 call parse_variable(sys%namespace,
'GOFireIntegrator', option__gofireintegrator__verlet, g_opt%fire_integrator)
606 call parse_variable(sys%namespace,
'GOObjective', minwhat_energy, g_opt%what2minimize)
667 if (g_opt%ions%natoms /= xyz%n)
then
668 write(
message(1),
'(a,i4,a,i4)')
'I need exactly ', g_opt%ions%natoms,
' constrains, but I found ', xyz%n
672 do iatom = 1, g_opt%ions%natoms
673 where(abs(xyz%atom(iatom)%x) <=
m_epsilon)
674 g_opt%ions%atom(iatom)%c =
m_zero
676 g_opt%ions%atom(iatom)%c =
m_one
683 if (g_opt%fixed_atom > 0)
then
687 do iatom = 1, g_opt%ions%natoms
688 g_opt%ions%atom(iatom)%c =
m_zero
693 call io_rm(
"geom/optimization.log", sys%namespace)
695 call io_rm(
"work-geom.xyz", sys%namespace)
697 if (.not. fromscratch)
then
698 inquire(file =
'./last.xyz', exist = does_exist)
699 if (.not. does_exist) fromscratch = .
true.
702 if (.not. fromscratch)
call g_opt%ions%read_xyz(
'./last')
707 write(filename,
'(a,i4.4,a)')
"geom/go.", iter,
".xyz"
708 inquire(file = trim(filename), exist = does_exist)
710 call io_rm(trim(filename), sys%namespace)
730 call g_opt%scfv%restart_dump%end()
738 safe_deallocate_a(g_opt%cell_force)
749 subroutine calc_point(size, coords, objective, getgrad, df)
750 integer,
intent(in) :: size
751 real(real64),
intent(in) :: coords(size)
752 real(real64),
intent(inout) :: objective
753 integer,
intent(in) :: getgrad
754 real(real64),
intent(inout) :: df(size)
756 integer :: iatom, idir, jdir
757 real(real64),
dimension(g_opt%periodic_dim, g_opt%periodic_dim) :: stress, scaling
762 assert(
size == g_opt%size)
766 if (g_opt%fixed_atom /= 0)
then
767 call g_opt%ions%translate(g_opt%ions%center())
772 call g_opt%ions%fold_atoms_into_cell()
775 do iatom = 1, g_opt%ions%natoms
776 if (.not. g_opt%syst%gr%box%contains_point(g_opt%syst%ions%pos(:, iatom)))
then
777 if (g_opt%syst%space%periodic_dim /= g_opt%syst%space%dim)
then
781 write(
message(1),
'(a,i5,a)')
"Atom ", iatom,
" has moved outside the box during the geometry optimization."
787 call g_opt%ions%write_xyz(
'./work-geom', append = .
true.)
794 g_opt%syst%space, g_opt%syst%hm%psolver, g_opt%syst%hm%kpoints, &
795 g_opt%syst%mc, g_opt%syst%st%qtot, g_opt%ions%latt)
799 g_opt%ions, g_opt%syst%ext_partners, g_opt%st)
800 call density_calc(g_opt%st, g_opt%syst%gr, g_opt%st%rho)
801 call v_ks_calc(g_opt%syst%ks, g_opt%syst%namespace, g_opt%syst%space, g_opt%hm, g_opt%st, &
802 g_opt%ions,g_opt%syst%ext_partners, calc_eigenval = .
true.)
803 call energy_calc_total(g_opt%syst%namespace, g_opt%syst%space, g_opt%hm, g_opt%syst%gr, g_opt%st, g_opt%syst%ext_partners)
806 call scf_run(g_opt%scfv, g_opt%syst%namespace, g_opt%syst%space, g_opt%syst%mc, g_opt%syst%gr, &
807 g_opt%ions, g_opt%syst%ext_partners, &
808 g_opt%st, g_opt%syst%ks, g_opt%hm, outp = g_opt%syst%outp, verbosity =
verb_compact, restart_dump=g_opt%scfv%restart_dump)
811 if (.not. g_opt%syst%space%is_periodic())
then
825 stress = g_opt%syst%st%stress_tensors%total(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
826 do idir = 1, g_opt%periodic_dim
827 stress(idir, idir) = -stress(idir, idir) - g_opt%pressure/g_opt%periodic_dim
829 g_opt%cell_force(:, :) = stress(:, :) * g_opt%ions%latt%rcell_volume
832 scaling = matmul(g_opt%ions%latt%rlattice(1:g_opt%periodic_dim,1:g_opt%periodic_dim), g_opt%inv_initial_rlattice)
834 g_opt%cell_force = matmul(g_opt%cell_force, transpose(scaling))
838 g_opt%cell_force =
m_half * (g_opt%cell_force + transpose(g_opt%cell_force))
843 stress = g_opt%syst%st%stress_tensors%total(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
844 do idir = 1, g_opt%periodic_dim
845 g_opt%cell_force(idir, 1) = -(g_opt%pressure/g_opt%periodic_dim + stress(idir, idir)) * g_opt%ions%latt%rcell_volume
851 do idir = 1, g_opt%periodic_dim
853 jdir = 1, g_opt%periodic_dim)
855 call messages_info(1+g_opt%periodic_dim, namespace=g_opt%ions%namespace, debug_only=.
true.)
857 do idir = 1, ubound(g_opt%cell_force, 2)
859 jdir = 1, g_opt%periodic_dim)
861 call messages_info(1+g_opt%periodic_dim, namespace=g_opt%ions%namespace, debug_only=.
true.)
866 if (getgrad == 1)
call to_grad(g_opt, df)
870 do iatom = 1, g_opt%ions%natoms
871 if (g_opt%ions%fixed(iatom)) cycle
872 objective = objective + sum(g_opt%ions%tot_force(:, iatom)**2)
875 do idir = 1, g_opt%periodic_dim
876 objective = objective + sum(g_opt%cell_force(:, idir)**2)
880 objective = objective + sum(g_opt%cell_force(:,1)**2)
882 objective =
sqrt(objective)
884 objective = g_opt%hm%energy%total
898 real(real64) :: coords(size)
899 real(real64) :: objective
902 real(real64),
allocatable :: df(:)
906 assert(
size == g_opt%size)
909 safe_allocate(df(1:size))
912 call calc_point(
size, coords, objective, getgrad, df)
913 safe_deallocate_a(df)
921 subroutine write_iter_info(geom_iter, size, energy, maxdx, maxdf, coords)
922 integer,
intent(in) :: geom_iter
923 integer,
intent(in) :: size
924 real(real64),
intent(in) :: energy, maxdx, maxdf
925 real(real64),
intent(in) :: coords(size)
927 character(len=256) :: c_geom_iter, title, c_forces_iter
932 write(c_geom_iter,
'(a,i4.4)')
"go.", geom_iter
934 call io_mkdir(
'geom', g_opt%ions%namespace)
935 call g_opt%ions%write_xyz(
'geom/'//trim(c_geom_iter), comment = trim(title))
936 call g_opt%ions%write_xyz(
'./last')
938 if(g_opt%periodic_dim > 0)
then
939 call g_opt%ions%write_xyz(
'geom/'//trim(c_geom_iter), comment =
'Reduced coordinates', reduce_coordinates = .
true.)
941 g_opt%ions%pos, g_opt%ions%atom, g_opt%syst%gr, g_opt%syst%namespace)
944 if (g_opt%syst%outp%what(option__output__forces))
then
945 write(c_forces_iter,
'(a,i4.4)')
"forces.", geom_iter
946 if (
bitand(g_opt%syst%outp%how(option__output__forces), option__outputformat__bild) /= 0)
then
947 call g_opt%ions%write_bild_forces_file(
'forces', trim(c_forces_iter))
950 g_opt%ions%pos, g_opt%ions%atom, g_opt%syst%gr, g_opt%syst%namespace, total_forces=g_opt%ions%tot_force)
955 iunit =
io_open(trim(
'geom/optimization.log'), g_opt%syst%namespace, &
956 action =
'write', position =
'append')
958 if (geom_iter == 1)
then
960 write(iunit,
'(a10,5(5x,a20),a)')
'# iter',
'energy [' // trim(
units_abbrev(
units_out%energy)) //
']', &
965 ' alpha, beta, gamma [degrees]'
967 write(iunit,
'(a10,3(5x,a20))')
'# iter',
'energy [' // trim(
units_abbrev(
units_out%energy)) //
']', &
981 g_opt%ions%latt%alpha, g_opt%ions%latt%beta, g_opt%ions%latt%gamma
994 call messages_write(
"++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .
true.)
1003 if (g_opt%periodic_dim == 0)
then
1014 call messages_write(maxdf, fmt =
"f16.10,1x", print_units = .false., new_line = .
true.)
1018 call messages_write(maxdx, fmt =
"f16.10,1x", print_units = .false., new_line = .
true.)
1021 call messages_write(
"++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .
true.)
1022 call messages_write(
"++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .
true.)
1034 real(real64),
intent(out) :: coords(:)
1036 integer :: iatom, idir, jdir, icoord
1037 real(real64) :: tmp_pos(gopt%dim), strain(g_opt%periodic_dim,g_opt%periodic_dim)
1043 if (
bitand(g_opt%type, go_ions) /= 0)
then
1044 do iatom = 1, gopt%ions%natoms
1045 if (gopt%fixed_atom == iatom) cycle
1046 if (gopt%ions%fixed(iatom)) cycle
1047 tmp_pos = gopt%ions%pos(1:gopt%dim, iatom)
1048 if (gopt%fixed_atom > 0) tmp_pos = tmp_pos - gopt%ions%pos(1:gopt%dim, gopt%fixed_atom)
1049 tmp_pos = gopt%ions%latt%cart_to_red(tmp_pos)
1050 do idir = 1, gopt%dim
1051 coords(icoord) = tmp_pos(idir)
1060 strain = matmul(gopt%ions%latt%rlattice, g_opt%inv_initial_rlattice)
1061 do idir = 1, g_opt%periodic_dim
1062 do jdir = idir, g_opt%periodic_dim
1063 coords(icoord) = strain(idir, jdir)
1071 do idir = 1, g_opt%periodic_dim
1072 coords(icoord) = norm2(gopt%ions%latt%rlattice(1:g_opt%periodic_dim, idir))/g_opt%initial_length(idir)
1083 subroutine to_grad(gopt, grad)
1085 real(real64),
intent(out) :: grad(:)
1087 integer :: iatom, idir, jdir, icoord
1088 real(real64) :: tmp_force(1:gopt%dim)
1094 if (
bitand(g_opt%type, go_ions) /= 0)
then
1095 do iatom = 1, gopt%ions%natoms
1096 if (gopt%fixed_atom == iatom) cycle
1097 if (gopt%ions%fixed(iatom)) cycle
1098 do idir = 1, gopt%dim
1099 if (abs(gopt%ions%atom(iatom)%c(idir)) <=
m_epsilon)
then
1100 tmp_force(idir) = -gopt%ions%tot_force(idir, iatom)
1104 if (gopt%fixed_atom > 0)
then
1105 tmp_force(idir) = tmp_force(idir) + gopt%ions%tot_force(idir, gopt%fixed_atom)
1108 tmp_force = gopt%ions%latt%cart_to_red(tmp_force)
1109 do idir = 1, gopt%dim
1110 grad(icoord) = tmp_force(idir)
1118 do idir = 1, g_opt%periodic_dim
1119 do jdir = idir, g_opt%periodic_dim
1120 grad(icoord) = -g_opt%cell_force(idir, jdir)
1128 do idir = 1, g_opt%periodic_dim
1129 grad(icoord) = -g_opt%cell_force(idir, 1)
1142 real(real64),
intent(in) :: coords(:)
1144 integer :: iatom, idir, jdir, icoord
1145 real(real64) :: tmp_pos(gopt%dim, gopt%ions%natoms), strain(g_opt%periodic_dim,g_opt%periodic_dim)
1153 if (
bitand(g_opt%type, go_ions) /= 0)
then
1154 do iatom = 1, gopt%ions%natoms
1155 if (gopt%fixed_atom == iatom) cycle
1156 if (gopt%ions%fixed(iatom)) cycle
1157 do idir = 1, gopt%dim
1158 tmp_pos(idir, iatom) = coords(icoord)
1163 do iatom = 1, gopt%ions%natoms
1164 tmp_pos(:, iatom) = gopt%ions%latt%cart_to_red(gopt%ions%pos(:, iatom))
1170 do idir = 1, g_opt%periodic_dim
1171 do jdir = idir, g_opt%periodic_dim
1172 strain(idir, jdir) = coords(icoord)
1180 gopt%ions%latt%rlattice = matmul(strain, g_opt%initial_rlattice)
1185 do idir = 1, g_opt%periodic_dim
1186 gopt%ions%latt%rlattice(1:g_opt%periodic_dim, idir) = coords(icoord) &
1187 * gopt%initial_rlattice(1:g_opt%periodic_dim, idir)
1194 call g_opt%syst%gr%symmetrizer%symmetrize_lattice_vectors(g_opt%periodic_dim, g_opt%initial_rlattice, &
1195 gopt%ions%latt%rlattice, gopt%symmetrize)
1196 call gopt%ions%update_lattice_vectors(gopt%ions%latt, gopt%symmetrize)
1200 if (
bitand(g_opt%type, go_ions) /= 0)
then
1202 do iatom = 1, gopt%ions%natoms
1203 if (gopt%fixed_atom == iatom) cycle
1204 if (gopt%ions%fixed(iatom)) cycle
1205 tmp_pos(:, iatom) = gopt%ions%latt%red_to_cart(tmp_pos(:, iatom))
1206 do idir = 1, gopt%dim
1207 if (abs(gopt%ions%atom(iatom)%c(idir)) <=
m_epsilon)
then
1208 gopt%ions%pos(idir, iatom) = tmp_pos(idir, iatom)
1211 if (gopt%fixed_atom > 0)
then
1212 gopt%ions%pos(:, iatom) = gopt%ions%pos(:, iatom) + gopt%ions%pos(:, gopt%fixed_atom)
1216 do iatom = 1, gopt%ions%natoms
1217 gopt%ions%pos(:, iatom) = gopt%ions%latt%red_to_cart(tmp_pos(:, iatom))
1221 if (gopt%symmetrize)
then
1222 call gopt%ions%symmetrize_atomic_coord()
1225 if (
debug%info)
then
1226 call gopt%ions%print_spacegroup()
1235 integer,
intent(in) :: geom_iter
1236 integer,
intent(in) :: size
1237 real(real64),
intent(in) :: energy, maxdx
1238 real(real64),
intent(in) :: coords(size)
subroutine init_(fromscratch)
Define which routines can be seen from the outside.
pure logical function, public accel_is_enabled()
type(debug_t), save, public debug
This module implements a calculator for the density and defines related functions.
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
subroutine, public energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit, full)
This subroutine calculates the total energy of the system. Basically, it adds up the KS eigenvalues,...
subroutine, public forces_set_total_to_zero(ions, force)
subroutine, public geom_opt_run(system, from_scratch)
subroutine calc_point_ng(size, coords, objective)
Same as calc_point, but without the gradients. No intents here is unfortunately required because the ...
subroutine to_grad(gopt, grad)
Transfer data from the forces to the work array for the gradients (grad)
integer, parameter go_cell
integer, parameter minwhat_forces
subroutine write_iter_info_ng(geom_iter, size, energy, maxdx, coords)
Same as write_iter_info, but without the gradients.
subroutine write_iter_info(geom_iter, size, energy, maxdx, maxdf, coords)
Output the information after each iteration of the geometry optimization.
subroutine calc_point(size, coords, objective, getgrad, df)
Note: you might think it would be better to change the arguments with '(size)' below to '(:)'....
subroutine to_coords(gopt, coords)
Transfer the data from the data structures to the work array (coords)
integer, parameter go_volume
subroutine from_coords(gopt, coords)
Transfer the data from the work array (coords) to the actual data structures.
subroutine geom_opt_run_legacy(sys, fromscratch)
real(real64), parameter, public m_zero
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
subroutine, public hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
subroutine, public write_xsf_geometry_file(dir, fname, space, latt, pos, atoms, mesh, namespace, total_forces)
subroutine, public io_close(iunit, grp)
subroutine, public io_rm(fname, namespace)
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
subroutine, public electrons_lattice_vectors_update(namespace, gr, space, psolver, kpoints, mc, qtot, new_latt)
subroutine, public lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r)
This module is intended to contain "only mathematical" functions and procedures.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
subroutine, public messages_new_line()
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
integer, parameter, public minmethod_nmsimplex
integer, parameter, public minmethod_fire
type(mpi_grp_t), public mpi_world
This module implements the basic mulsisystem class, a container system for other systems.
logical function, public parse_is_defined(namespace, name)
integer, parameter, public read_coords_err
for read_coords_info::file_type
subroutine, public read_coords_init(gf)
subroutine, public read_coords_end(gf)
subroutine, public read_coords_read(what, gf, space, namespace)
integer, parameter, public restart_gs
integer, parameter, public restart_type_dump
integer, parameter, public restart_type_load
subroutine, public scf_print_mem_use(namespace)
subroutine, public scf_mix_clear(scf)
integer, parameter, public verb_compact
subroutine, public scf_init(scf, namespace, gr, ions, st, mc, hm, space)
subroutine, public scf_end(scf)
subroutine, public scf_run(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, verbosity, iters_done, restart_dump)
Legacy version of the SCF code.
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, ierr, iter, lr, lowest_missing, label, verbose, skip)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
This module defines the unit system, used for input and output.
type(unit_t), public unit_femtosecond
Time in femtoseconds.
type(unit_t), public unit_amu
Mass in atomic mass units (AKA Dalton).
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
subroutine, public v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
An abstract type for all electron species.
Class describing the electron system.
Container class for lists of system_oct_m::system_t.