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
146 logical :: known_lower_bound
150 if (sys%space%periodic_dim == 1 .or. sys%space%periodic_dim == 2)
then
151 message(1) =
"Geometry optimization not allowed for systems periodic in 1D and 2D, "
152 message(2) =
"as in those cases the total energy is not correct."
157 if (sys%hm%pcm%run_pcm)
then
161 if (sys%kpoints%use_symmetries)
then
168 if (.not. fromscratch)
then
171 call states_elec_load(restart_load, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
172 fixed_occ=.false., ierr=ierr)
174 call restart_load%end()
176 message(1) =
"Unable to read wavefunctions: Starting from scratch."
182 call scf_init(g_opt%scfv, sys%namespace, sys%gr, sys%ions, sys%st, sys%mc, sys%hm, sys%space)
185 if (.not. g_opt%scfv%calc_stress)
then
186 message(1) =
"In order to optimize the cell, one needs to set SCFCalculateStress = yes."
191 if (fromscratch)
then
192 call lcao_run(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm, &
193 lmm_r = g_opt%scfv%lmm_r, known_lower_bound=known_lower_bound)
197 message(1) =
'Info: Setting up Hamiltonian.'
199 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
203 g_opt%symmetrize = sys%kpoints%use_symmetries .or. sys%st%symmetrize_density
206 safe_allocate(coords(1:g_opt%size))
209 if (sys%st%pack_states .and. sys%hm%apply_packed())
call sys%st%pack()
212 select case (g_opt%method)
214 call minimize_multidim_nograd(g_opt%method, g_opt%size, coords, g_opt%step,&
215 g_opt%toldr, g_opt%max_iter, &
219 safe_allocate(mass(1:g_opt%size))
220 mass = g_opt%fire_mass
222 do iatom = 1, sys%ions%natoms
223 if (g_opt%fixed_atom == iatom) cycle
224 if (g_opt%ions%fixed(iatom)) cycle
225 if (g_opt%fire_mass <=
m_zero) mass(imass:imass + 2) = sys%ions%mass(iatom)
226 imass = imass + g_opt%dim
230 call minimize_fire(g_opt%size, g_opt%ions%space%dim, coords, g_opt%step, g_opt%tolgrad, &
232 safe_deallocate_a(mass)
235 call minimize_multidim(g_opt%method, g_opt%size, coords, g_opt%step ,&
236 g_opt%line_tol , g_opt%tolgrad, g_opt%toldr, g_opt%max_iter, &
240 if (ierr == 1025)
then
242 message(1) =
"Reached maximum number of iterations allowed by GOMaxIter."
245 message(1) =
"Error occurred during the GSL minimization procedure:"
250 if (sys%st%pack_states .and. sys%hm%apply_packed())
call sys%st%unpack()
254 message(1) =
"Writing final coordinates to min.xyz"
257 call g_opt%ions%write_xyz(
'./min')
259 safe_deallocate_a(coords)
262 call g_opt%scfv%criterion_list%empty()
269 subroutine init_(fromscratch)
270 logical,
intent(inout) :: fromscratch
272 logical :: center, does_exist
273 integer :: iter, iatom, idir
274 character(len=100) :: filename
275 real(real64) :: default_toldr
276 real(real64) :: default_step
281 if (sys%space%is_periodic())
then
309 write(
message(1),
'(a)')
'Input: [GOType = '
310 if (
bitand(g_opt%type, go_ions) /= 0)
then
314 if (len_trim(
message(1)) > 16)
then
320 if (len_trim(
message(1)) > 16)
then
329 message(1) =
"Cell and volume optimization cannot be used simultaneously."
343 message(1) =
"Cell dynamics not supported on GPUs."
349 do iatom = 1, sys%ions%natoms
350 select type(spec=>sys%ions%atom(iatom)%species)
352 write(
message(1),
'(a)')
"Geometry optimization for all-electron potential is not implemented."
362 g_opt%ions => sys%ions
366 g_opt%dim = sys%space%dim
367 g_opt%periodic_dim = sys%space%periodic_dim
371 if (
bitand(g_opt%type, go_ions) /= 0)
then
372 g_opt%size = g_opt%dim * g_opt%ions%natoms
377 g_opt%size = g_opt%size + (g_opt%periodic_dim +1) * g_opt%periodic_dim / 2
378 safe_allocate(g_opt%cell_force(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
383 g_opt%size = g_opt%size + g_opt%periodic_dim
384 safe_allocate(g_opt%cell_force(1:g_opt%periodic_dim, 1:1))
386 safe_allocate(g_opt%initial_length(1:g_opt%periodic_dim))
387 do idir = 1, g_opt%periodic_dim
388 g_opt%initial_length(idir) = norm2(g_opt%ions%latt%rlattice(1:g_opt%periodic_dim, idir))
394 safe_allocate(g_opt%initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
395 g_opt%initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim) &
396 = g_opt%ions%latt%rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
397 safe_allocate(g_opt%inv_initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
398 g_opt%inv_initial_rlattice(:, :) = g_opt%initial_rlattice(:, :)
399 call lalg_inverse(g_opt%periodic_dim, g_opt%inv_initial_rlattice,
'dir')
402 if(g_opt%ions%space%is_periodic())
then
418 if (center .and.
bitand(g_opt%type, go_ions) /= 0)
then
420 g_opt%size = g_opt%size - g_opt%dim
425 do iatom = 1, g_opt%ions%natoms
426 if (g_opt%ions%fixed(iatom) .and.
bitand(g_opt%type, go_ions) /= 0)
then
427 g_opt%size = g_opt%size - g_opt%dim
431 assert(g_opt%size > 0)
507 default_toldr = 0.001_real64
509 default_toldr = -
m_one
526 call parse_variable(sys%namespace,
'GOStep', default_step, g_opt%step)
541 call parse_variable(sys%namespace,
'GOLineTol', 0.1_real64, g_opt%line_tol)
551 call parse_variable(sys%namespace,
'GOMaxIter', 200, g_opt%max_iter)
552 if (g_opt%max_iter <= 0)
then
553 message(1) =
"GOMaxIter has to be larger than 0"
590 call parse_variable(sys%namespace,
'GOFireIntegrator', option__gofireintegrator__verlet, g_opt%fire_integrator)
611 call parse_variable(sys%namespace,
'GOObjective', minwhat_energy, g_opt%what2minimize)
672 if (g_opt%ions%natoms /= xyz%n)
then
673 write(
message(1),
'(a,i4,a,i4)')
'I need exactly ', g_opt%ions%natoms,
' constrains, but I found ', xyz%n
677 do iatom = 1, g_opt%ions%natoms
678 where(abs(xyz%atom(iatom)%x) <=
m_epsilon)
679 g_opt%ions%atom(iatom)%c =
m_zero
681 g_opt%ions%atom(iatom)%c =
m_one
688 if (g_opt%fixed_atom > 0)
then
692 do iatom = 1, g_opt%ions%natoms
693 g_opt%ions%atom(iatom)%c =
m_zero
698 call io_rm(
"geom/optimization.log", sys%namespace)
700 call io_rm(
"work-geom.xyz", sys%namespace)
702 if (.not. fromscratch)
then
703 inquire(file =
'./last.xyz', exist = does_exist)
704 if (.not. does_exist) fromscratch = .
true.
707 if (.not. fromscratch)
call g_opt%ions%read_xyz(
'./last')
712 write(filename,
'(a,i4.4,a)')
"geom/go.", iter,
".xyz"
713 inquire(file = trim(filename), exist = does_exist)
715 call io_rm(trim(filename), sys%namespace)
735 call g_opt%scfv%restart_dump%end()
743 safe_deallocate_a(g_opt%cell_force)
754 subroutine calc_point(size, coords, objective, getgrad, df)
755 integer,
intent(in) :: size
756 real(real64),
intent(in) :: coords(size)
757 real(real64),
intent(inout) :: objective
758 integer,
intent(in) :: getgrad
759 real(real64),
intent(inout) :: df(size)
761 integer :: iatom, idir, jdir
762 real(real64),
dimension(g_opt%periodic_dim, g_opt%periodic_dim) :: stress, strain, right_stretch, rotation, sym_stress
767 assert(
size == g_opt%size)
775 if (g_opt%fixed_atom /= 0)
then
776 call g_opt%ions%translate(g_opt%ions%center())
781 call g_opt%ions%fold_atoms_into_cell()
784 do iatom = 1, g_opt%ions%natoms
785 if (.not. g_opt%syst%gr%box%contains_point(g_opt%syst%ions%pos(:, iatom)))
then
786 if (g_opt%syst%space%periodic_dim /= g_opt%syst%space%dim)
then
790 write(
message(1),
'(a,i5,a)')
"Atom ", iatom,
" has moved outside the box during the geometry optimization."
796 call g_opt%ions%write_xyz(
'./work-geom', append = .
true.)
803 g_opt%syst%space, g_opt%syst%hm%psolver, g_opt%syst%hm%kpoints, &
804 g_opt%syst%mc, g_opt%syst%st%qtot, g_opt%ions%latt)
808 g_opt%ions, g_opt%syst%ext_partners, g_opt%st)
809 call density_calc(g_opt%st, g_opt%syst%gr, g_opt%st%rho)
810 call v_ks_calc(g_opt%syst%ks, g_opt%syst%namespace, g_opt%syst%space, g_opt%hm, g_opt%st, &
811 g_opt%ions,g_opt%syst%ext_partners, calc_eigenval = .
true.)
812 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)
815 call scf_run(g_opt%scfv, g_opt%syst%namespace, g_opt%syst%space, g_opt%syst%mc, g_opt%syst%gr, &
816 g_opt%ions, g_opt%syst%ext_partners, &
817 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)
820 if (g_opt%ions%force_total_enforce)
then
836 stress = -g_opt%syst%st%stress_tensors%total(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
840 strain = matmul(g_opt%ions%latt%rlattice(1:g_opt%periodic_dim,1:g_opt%periodic_dim), g_opt%inv_initial_rlattice)
842 call lalg_inverse(g_opt%periodic_dim, right_stretch,
'dir')
843 rotation = matmul(strain, right_stretch)
846 sym_stress = matmul(transpose(rotation), matmul(stress, rotation))
847 sym_stress =
m_half*(sym_stress + transpose(sym_stress))
849 do idir = 1, g_opt%periodic_dim
850 sym_stress(idir, idir) = sym_stress(idir, idir) - g_opt%pressure/g_opt%periodic_dim
852 g_opt%cell_force = sym_stress * g_opt%ions%latt%rcell_volume
854 g_opt%cell_force = matmul(g_opt%cell_force, right_stretch)
859 stress = g_opt%syst%st%stress_tensors%total(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
860 do idir = 1, g_opt%periodic_dim
861 g_opt%cell_force(idir, 1) = -(g_opt%pressure/g_opt%periodic_dim + stress(idir, idir)) * g_opt%ions%latt%rcell_volume
867 do idir = 1, g_opt%periodic_dim
869 jdir = 1, g_opt%periodic_dim)
871 call messages_info(1+g_opt%periodic_dim, namespace=g_opt%ions%namespace, debug_only=.
true.)
873 do idir = 1, ubound(g_opt%cell_force, 2)
875 jdir = 1, g_opt%periodic_dim)
877 call messages_info(1+g_opt%periodic_dim, namespace=g_opt%ions%namespace, debug_only=.
true.)
882 if (getgrad == 1)
call to_grad(g_opt, df)
886 do iatom = 1, g_opt%ions%natoms
887 if (g_opt%ions%fixed(iatom)) cycle
888 objective = objective + sum(g_opt%ions%tot_force(:, iatom)**2)
891 do idir = 1, g_opt%periodic_dim
892 objective = objective + sum(g_opt%cell_force(:, idir)**2)
896 objective = objective + sum(g_opt%cell_force(:,1)**2)
898 objective =
sqrt(objective)
900 objective = g_opt%hm%energy%total
914 real(real64) :: coords(size)
915 real(real64) :: objective
918 real(real64),
allocatable :: df(:)
922 assert(
size == g_opt%size)
925 safe_allocate(df(1:size))
928 call calc_point(
size, coords, objective, getgrad, df)
929 safe_deallocate_a(df)
937 subroutine write_iter_info(geom_iter, size, energy, maxdx, maxdf, coords)
938 integer,
intent(in) :: geom_iter
939 integer,
intent(in) :: size
940 real(real64),
intent(in) :: energy, maxdx, maxdf
941 real(real64),
intent(in) :: coords(size)
943 character(len=256) :: c_geom_iter, title, c_forces_iter
948 write(c_geom_iter,
'(a,i4.4)')
"go.", geom_iter
950 call io_mkdir(
'geom', g_opt%ions%namespace)
951 call g_opt%ions%write_xyz(
'geom/'//trim(c_geom_iter), comment = trim(title))
952 call g_opt%ions%write_xyz(
'./last')
954 if(g_opt%periodic_dim > 0)
then
955 call g_opt%ions%write_xyz(
'geom/'//trim(c_geom_iter), comment =
'Reduced coordinates', reduce_coordinates = .
true.)
957 g_opt%ions%pos, g_opt%ions%atom, g_opt%syst%gr, g_opt%syst%namespace)
960 if (g_opt%syst%outp%what(option__output__forces))
then
961 write(c_forces_iter,
'(a,i4.4)')
"forces.", geom_iter
962 if (
bitand(g_opt%syst%outp%how(option__output__forces), option__outputformat__bild) /= 0)
then
963 call g_opt%ions%write_bild_forces_file(
'forces', trim(c_forces_iter))
966 g_opt%ions%pos, g_opt%ions%atom, g_opt%syst%gr, g_opt%syst%namespace, total_forces=g_opt%ions%tot_force)
970 if (g_opt%syst%st%system_grp%is_root())
then
971 iunit =
io_open(trim(
'geom/optimization.log'), g_opt%syst%namespace, &
972 action =
'write', position =
'append')
974 if (geom_iter == 1)
then
976 write(iunit,
'(a10,5(5x,a20),a)')
'# iter',
'energy [' // trim(
units_abbrev(
units_out%energy)) //
']', &
981 ' alpha, beta, gamma [degrees]'
983 write(iunit,
'(a10,3(5x,a20))')
'# iter',
'energy [' // trim(
units_abbrev(
units_out%energy)) //
']', &
997 g_opt%ions%latt%alpha, g_opt%ions%latt%beta, g_opt%ions%latt%gamma
1010 call messages_write(
"++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .
true.)
1012 call messages_write(
"+++++++++++++++++++++ MINIMIZATION ITER #:")
1019 if (g_opt%periodic_dim == 0)
then
1030 call messages_write(maxdf, fmt =
"f16.10,1x", print_units = .false., new_line = .
true.)
1034 call messages_write(maxdx, fmt =
"f16.10,1x", print_units = .false., new_line = .
true.)
1037 call messages_write(
"++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .
true.)
1038 call messages_write(
"++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .
true.)
1050 real(real64),
intent(out) :: coords(:)
1052 integer :: iatom, idir, jdir, icoord
1053 real(real64) :: tmp_pos(gopt%dim), strain(g_opt%periodic_dim,g_opt%periodic_dim)
1059 if (
bitand(g_opt%type, go_ions) /= 0)
then
1060 do iatom = 1, gopt%ions%natoms
1061 if (gopt%fixed_atom == iatom) cycle
1062 if (gopt%ions%fixed(iatom)) cycle
1063 tmp_pos = gopt%ions%pos(1:gopt%dim, iatom)
1064 if (gopt%fixed_atom > 0) tmp_pos = tmp_pos - gopt%ions%pos(1:gopt%dim, gopt%fixed_atom)
1065 tmp_pos = gopt%ions%latt%cart_to_red(tmp_pos)
1066 do idir = 1, gopt%dim
1067 coords(icoord) = tmp_pos(idir)
1076 strain = matmul(gopt%ions%latt%rlattice, g_opt%inv_initial_rlattice)
1077 do idir = 1, g_opt%periodic_dim
1078 do jdir = idir, g_opt%periodic_dim
1079 coords(icoord) = strain(idir, jdir)
1087 do idir = 1, g_opt%periodic_dim
1088 coords(icoord) = norm2(gopt%ions%latt%rlattice(1:g_opt%periodic_dim, idir))/g_opt%initial_length(idir)
1099 subroutine to_grad(gopt, grad)
1101 real(real64),
intent(out) :: grad(:)
1103 integer :: iatom, idir, jdir, icoord
1104 real(real64) :: tmp_force(1:gopt%dim)
1110 if (
bitand(g_opt%type, go_ions) /= 0)
then
1111 do iatom = 1, gopt%ions%natoms
1112 if (gopt%fixed_atom == iatom) cycle
1113 if (gopt%ions%fixed(iatom)) cycle
1114 do idir = 1, gopt%dim
1115 if (abs(gopt%ions%atom(iatom)%c(idir)) <=
m_epsilon)
then
1116 tmp_force(idir) = -gopt%ions%tot_force(idir, iatom)
1120 if (gopt%fixed_atom > 0)
then
1121 tmp_force(idir) = tmp_force(idir) + gopt%ions%tot_force(idir, gopt%fixed_atom)
1124 tmp_force = gopt%ions%latt%cart_to_red(tmp_force)
1125 do idir = 1, gopt%dim
1126 grad(icoord) = tmp_force(idir)
1134 do idir = 1, g_opt%periodic_dim
1135 do jdir = idir, g_opt%periodic_dim
1136 grad(icoord) = -g_opt%cell_force(idir, jdir)
1144 do idir = 1, g_opt%periodic_dim
1145 grad(icoord) = -g_opt%cell_force(idir, 1)
1158 real(real64),
intent(in) :: coords(:)
1160 integer :: iatom, idir, jdir, icoord
1161 real(real64) :: tmp_pos(gopt%dim, gopt%ions%natoms), strain(g_opt%periodic_dim,g_opt%periodic_dim)
1169 if (
bitand(g_opt%type, go_ions) /= 0)
then
1170 do iatom = 1, gopt%ions%natoms
1171 if (gopt%fixed_atom == iatom) cycle
1172 if (gopt%ions%fixed(iatom)) cycle
1173 do idir = 1, gopt%dim
1174 tmp_pos(idir, iatom) = coords(icoord)
1179 do iatom = 1, gopt%ions%natoms
1180 tmp_pos(:, iatom) = gopt%ions%latt%cart_to_red(gopt%ions%pos(:, iatom))
1186 do idir = 1, g_opt%periodic_dim
1187 do jdir = idir, g_opt%periodic_dim
1188 strain(idir, jdir) = coords(icoord)
1196 gopt%ions%latt%rlattice = matmul(strain, g_opt%initial_rlattice)
1201 do idir = 1, g_opt%periodic_dim
1202 gopt%ions%latt%rlattice(1:g_opt%periodic_dim, idir) = coords(icoord) &
1203 * gopt%initial_rlattice(1:g_opt%periodic_dim, idir)
1210 call g_opt%syst%gr%symmetrizer%symmetrize_lattice_vectors(g_opt%periodic_dim, g_opt%initial_rlattice, &
1211 gopt%ions%latt%rlattice, gopt%symmetrize)
1212 call gopt%ions%update_lattice_vectors(gopt%ions%latt, gopt%symmetrize)
1216 if (
bitand(g_opt%type, go_ions) /= 0)
then
1218 do iatom = 1, gopt%ions%natoms
1219 if (gopt%fixed_atom == iatom) cycle
1220 if (gopt%ions%fixed(iatom)) cycle
1221 tmp_pos(:, iatom) = gopt%ions%latt%red_to_cart(tmp_pos(:, iatom))
1222 do idir = 1, gopt%dim
1223 if (abs(gopt%ions%atom(iatom)%c(idir)) <=
m_epsilon)
then
1224 gopt%ions%pos(idir, iatom) = tmp_pos(idir, iatom)
1227 if (gopt%fixed_atom > 0)
then
1228 gopt%ions%pos(:, iatom) = gopt%ions%pos(:, iatom) + gopt%ions%pos(:, gopt%fixed_atom)
1232 do iatom = 1, gopt%ions%natoms
1233 gopt%ions%pos(:, iatom) = gopt%ions%latt%red_to_cart(tmp_pos(:, iatom))
1237 if (gopt%symmetrize)
then
1238 call gopt%ions%symmetrize_atomic_coord()
1241 if (
debug%info)
then
1242 call gopt%ions%print_spacegroup()
1251 integer,
intent(in) :: geom_iter
1252 integer,
intent(in) :: size
1253 real(real64),
intent(in) :: energy, maxdx
1254 real(real64),
intent(in) :: coords(size)
subroutine init_(fromscratch)
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 ion_dynamics_box_update(namespace, gr, space, new_latt)
real(real64) function, dimension(1:n, 1:n), public lalg_remove_rotation(n, A)
Remove rotation from affine transformation A by computing the polar decomposition and discarding the ...
subroutine, public lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r, known_lower_bound)
System information (time, memory, sysname)
subroutine, public loct_strerror(errno, res)
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
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_set_lower_bound_is_known(scf, known_lower_bound)
Set the flag lower_bound_is_known.
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, fixed_occ, 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.