34 use,
intrinsic :: iso_fortran_env
68 integer(int64) :: type
71 real(real64) :: line_tol
72 real(real64) :: fire_mass
73 integer :: fire_integrator
74 real(real64) :: tolgrad
77 integer :: what2minimize
81 type(ions_t),
pointer :: ions
82 type(hamiltonian_elec_t),
pointer :: hm
83 type(electrons_t),
pointer :: syst
84 class(mesh_t),
pointer :: mesh
85 type(states_elec_t),
pointer :: st
87 integer :: periodic_dim
89 integer :: fixed_atom = 0
91 real(real64),
allocatable :: cell_force(:, :)
92 logical :: symmetrize = .false.
93 real(real64),
allocatable :: initial_length(:)
94 real(real64),
allocatable :: initial_rlattice(:, :)
95 real(real64),
allocatable :: inv_initial_rlattice(:, :)
96 real(real64) :: pressure
99 type(geom_opt_t),
save :: g_opt
101 integer,
parameter :: &
102 MINWHAT_ENERGY = 1, &
105 integer,
parameter :: &
114 class(*),
intent(inout) :: system
115 logical,
intent(inout) :: from_scratch
121 message(1) =
"CalculationMode = go not implemented for multi-system calculations"
132 type(electrons_t),
target,
intent(inout) :: sys
133 logical,
intent(inout) :: fromscratch
136 real(real64),
allocatable :: coords(:)
137 real(real64) :: energy
139 real(real64),
allocatable :: mass(:)
140 integer :: iatom, imass
141 type(restart_t) :: restart_load
145 if (sys%space%periodic_dim == 1 .or. sys%space%periodic_dim == 2)
then
146 message(1) =
"Geometry optimization not allowed for systems periodic in 1D and 2D, "
147 message(2) =
"as in those cases the total energy is not correct."
152 if (sys%hm%pcm%run_pcm)
then
156 if (sys%kpoints%use_symmetries)
then
160 call init_(fromscratch)
163 if (.not. fromscratch)
then
166 call states_elec_load(restart_load, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
167 fixed_occ=.false., ierr=ierr)
171 message(1) =
"Unable to read wavefunctions: Starting from scratch."
177 call scf_init(g_opt%scfv, sys%namespace, sys%gr, sys%ions, sys%st, sys%mc, sys%hm, sys%space)
180 if (.not. g_opt%scfv%calc_stress)
then
181 message(1) =
"In order to optimize the cell, one needs to set SCFCalculateStress = yes."
186 if (fromscratch)
then
187 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)
190 message(1) =
'Info: Setting up Hamiltonian.'
192 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
195 g_opt%symmetrize = sys%kpoints%use_symmetries .or. sys%st%symmetrize_density
198 safe_allocate(coords(1:g_opt%size))
201 if (sys%st%pack_states .and. sys%hm%apply_packed())
call sys%st%pack()
204 select case (g_opt%method)
206 call minimize_multidim_nograd(g_opt%method, g_opt%size, coords, real(g_opt%step, real64) ,&
207 real(g_opt%toldr, real64) , g_opt%max_iter, &
211 safe_allocate(mass(1:g_opt%size))
212 mass = g_opt%fire_mass
214 do iatom = 1, sys%ions%natoms
215 if (g_opt%fixed_atom == iatom) cycle
216 if (g_opt%ions%fixed(iatom)) cycle
217 if (g_opt%fire_mass <=
m_zero) mass(imass:imass + 2) = sys%ions%mass(iatom)
218 imass = imass + g_opt%dim
222 call minimize_fire(g_opt%size, g_opt%ions%space%dim, coords, real(g_opt%step, real64) , real(g_opt%tolgrad, real64) , &
224 safe_deallocate_a(mass)
227 call minimize_multidim(g_opt%method, g_opt%size, coords, real(g_opt%step, real64) ,&
228 real(g_opt%line_tol, real64) ,
real(g_opt%tolgrad, real64) ,
real(g_opt%toldr, real64) , g_opt%max_iter, &
229 calc_point, write_iter_info, energy, ierr)
232 if (ierr == 1025)
then
234 message(1) =
"Reached maximum number of iterations allowed by GOMaxIter."
237 message(1) =
"Error occurred during the GSL minimization procedure:"
242 if (sys%st%pack_states .and. sys%hm%apply_packed())
call sys%st%unpack()
246 message(1) =
"Writing final coordinates to min.xyz"
249 call g_opt%ions%write_xyz(
'./min')
251 safe_deallocate_a(coords)
254 call g_opt%scfv%criterion_list%empty()
261 subroutine init_(fromscratch)
262 logical,
intent(inout) :: fromscratch
264 logical :: center, does_exist
265 integer :: iter, iatom, idir
266 character(len=100) :: filename
267 real(real64) :: default_toldr
268 real(real64) :: default_step
273 if (sys%space%is_periodic())
then
301 write(
message(1),
'(a)')
'Input: [GOType = '
302 if (
bitand(g_opt%type, go_ions) /= 0)
then
306 if (len_trim(
message(1)) > 16)
then
312 if (len_trim(
message(1)) > 16)
then
321 message(1) =
"Cell and volume optimization cannot be used simultaneously."
335 message(1) =
"Cell dynamics not supported on GPUs."
341 do iatom = 1, sys%ions%natoms
342 select type(spec=>sys%ions%atom(iatom)%species)
344 write(
message(1),
'(a)')
"Geometry optimization for all-electron potential is not implemented."
354 g_opt%ions => sys%ions
358 g_opt%dim = sys%space%dim
359 g_opt%periodic_dim = sys%space%periodic_dim
363 if (
bitand(g_opt%type, go_ions) /= 0)
then
364 g_opt%size = g_opt%dim * g_opt%ions%natoms
369 g_opt%size = g_opt%size + (g_opt%periodic_dim +1) * g_opt%periodic_dim / 2
370 safe_allocate(g_opt%cell_force(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
375 g_opt%size = g_opt%size + g_opt%periodic_dim
376 safe_allocate(g_opt%cell_force(1:g_opt%periodic_dim, 1:1))
378 safe_allocate(g_opt%initial_length(1:g_opt%periodic_dim))
379 do idir = 1, g_opt%periodic_dim
380 g_opt%initial_length(idir) = norm2(g_opt%ions%latt%rlattice(1:g_opt%periodic_dim, idir))
386 safe_allocate(g_opt%initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
387 g_opt%initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim) &
388 = g_opt%ions%latt%rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
389 safe_allocate(g_opt%inv_initial_rlattice(1:g_opt%periodic_dim, 1:g_opt%periodic_dim))
390 g_opt%inv_initial_rlattice(:, :) = g_opt%initial_rlattice(:, :)
391 call lalg_inverse(g_opt%periodic_dim, g_opt%inv_initial_rlattice,
'dir')
394 if(g_opt%ions%space%is_periodic())
then
410 if (center .and.
bitand(g_opt%type, go_ions) /= 0)
then
412 g_opt%size = g_opt%size - g_opt%dim
417 do iatom = 1, g_opt%ions%natoms
418 if (g_opt%ions%fixed(iatom) .and.
bitand(g_opt%type, go_ions) /= 0)
then
419 g_opt%size = g_opt%size - g_opt%dim
423 assert(g_opt%size > 0)
499 default_toldr = 0.001_real64
501 default_toldr = -
m_one
518 call parse_variable(sys%namespace,
'GOStep', default_step, g_opt%step)
533 call parse_variable(sys%namespace,
'GOLineTol', 0.1_real64, g_opt%line_tol)
543 call parse_variable(sys%namespace,
'GOMaxIter', 200, g_opt%max_iter)
544 if (g_opt%max_iter <= 0)
then
545 message(1) =
"GOMaxIter has to be larger than 0"
580 call parse_variable(sys%namespace,
'GOFireIntegrator', option__gofireintegrator__verlet, g_opt%fire_integrator)
601 call parse_variable(sys%namespace,
'GOObjective', minwhat_energy, g_opt%what2minimize)
662 if (g_opt%ions%natoms /= xyz%n)
then
663 write(
message(1),
'(a,i4,a,i4)')
'I need exactly ', g_opt%ions%natoms,
' constrains, but I found ', xyz%n
667 do iatom = 1, g_opt%ions%natoms
668 where(abs(xyz%atom(iatom)%x) <=
m_epsilon)
669 g_opt%ions%atom(iatom)%c =
m_zero
671 g_opt%ions%atom(iatom)%c =
m_one
677 if (g_opt%fixed_atom > 0)
then
681 do iatom = 1, g_opt%ions%natoms
682 g_opt%ions%atom(iatom)%c =
m_zero
686 call io_rm(
"geom/optimization.log", sys%namespace)
688 call io_rm(
"work-geom.xyz", sys%namespace)
690 if (.not. fromscratch)
then
691 inquire(file =
'./last.xyz', exist = does_exist)
692 if (.not. does_exist) fromscratch = .
true.
695 if (.not. fromscratch)
call g_opt%ions%read_xyz(
'./last')
700 write(filename,
'(a,i4.4,a)')
"geom/go.", iter,
".xyz"
701 inquire(file = trim(filename), exist = does_exist)
703 call io_rm(trim(filename), sys%namespace)
731 safe_deallocate_a(g_opt%cell_force)
742 subroutine calc_point(size, coords, objective, getgrad, df)
743 integer,
intent(in) :: size
744 real(real64),
intent(in) :: coords(size)
745 real(real64),
intent(inout) :: objective
746 integer,
intent(in) :: getgrad
747 real(real64),
intent(inout) :: df(size)
749 integer :: iatom, idir, jdir
750 real(real64),
dimension(g_opt%periodic_dim, g_opt%periodic_dim) :: stress, strain, right_stretch, rotation, sym_stress
755 assert(
size == g_opt%size)
763 if (g_opt%fixed_atom /= 0)
then
764 call g_opt%ions%translate(g_opt%ions%center())
768 call g_opt%ions%fold_atoms_into_cell()
771 do iatom = 1, g_opt%ions%natoms
772 if (.not. g_opt%syst%gr%box%contains_point(g_opt%syst%ions%pos(:, iatom)))
then
773 if (g_opt%syst%space%periodic_dim /= g_opt%syst%space%dim)
then
777 write(
message(1),
'(a,i5,a)')
"Atom ", iatom,
" has moved outside the box during the geometry optimization."
783 call g_opt%ions%write_xyz(
'./work-geom', append = .
true.)
790 g_opt%syst%space, g_opt%syst%hm%psolver, g_opt%syst%hm%kpoints, &
791 g_opt%syst%mc, g_opt%syst%st%qtot, g_opt%ions%latt)
795 g_opt%ions, g_opt%syst%ext_partners, g_opt%st)
796 call density_calc(g_opt%st, g_opt%syst%gr, g_opt%st%rho)
797 call v_ks_calc(g_opt%syst%ks, g_opt%syst%namespace, g_opt%syst%space, g_opt%hm, g_opt%st, &
798 g_opt%ions,g_opt%syst%ext_partners, calc_eigenval = .
true.)
799 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)
802 call scf_run(g_opt%scfv, g_opt%syst%namespace, g_opt%syst%space, g_opt%syst%mc, g_opt%syst%gr, &
803 g_opt%ions, g_opt%syst%ext_partners, &
804 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)
818 stress = -g_opt%syst%st%stress_tensors%total(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
822 strain = matmul(g_opt%ions%latt%rlattice(1:g_opt%periodic_dim,1:g_opt%periodic_dim), g_opt%inv_initial_rlattice)
824 call lalg_inverse(g_opt%periodic_dim, right_stretch,
'dir')
825 rotation = matmul(strain, right_stretch)
828 sym_stress = matmul(transpose(rotation), matmul(stress, rotation))
829 sym_stress =
m_half*(sym_stress + transpose(sym_stress))
831 do idir = 1, g_opt%periodic_dim
832 sym_stress(idir, idir) = sym_stress(idir, idir) - g_opt%pressure/g_opt%periodic_dim
834 g_opt%cell_force = sym_stress * g_opt%ions%latt%rcell_volume
836 g_opt%cell_force = matmul(g_opt%cell_force, right_stretch)
841 stress = g_opt%syst%st%stress_tensors%total(1:g_opt%periodic_dim, 1:g_opt%periodic_dim)
842 do idir = 1, g_opt%periodic_dim
843 g_opt%cell_force(idir, 1) = -(g_opt%pressure/g_opt%periodic_dim + stress(idir, idir)) * g_opt%ions%latt%rcell_volume
849 do idir = 1, g_opt%periodic_dim
851 jdir = 1, g_opt%periodic_dim)
853 call messages_info(1+g_opt%periodic_dim, namespace=g_opt%ions%namespace, debug_only=.
true.)
855 do idir = 1, ubound(g_opt%cell_force, 2)
857 jdir = 1, g_opt%periodic_dim)
859 call messages_info(1+g_opt%periodic_dim, namespace=g_opt%ions%namespace, debug_only=.
true.)
864 if (getgrad == 1)
call to_grad(g_opt, df)
868 do iatom = 1, g_opt%ions%natoms
869 if (g_opt%ions%fixed(iatom)) cycle
870 objective = objective + sum(g_opt%ions%tot_force(:, iatom)**2)
873 do idir = 1, g_opt%periodic_dim
874 objective = objective + sum(g_opt%cell_force(:, idir)**2)
878 objective = objective + sum(g_opt%cell_force(:,1)**2)
880 objective =
sqrt(objective)
882 objective = g_opt%hm%energy%total
896 real(real64) :: coords(size)
897 real(real64) :: objective
900 real(real64),
allocatable :: df(:)
904 assert(
size == g_opt%size)
907 safe_allocate(df(1:size))
910 call calc_point(
size, coords, objective, getgrad, df)
911 safe_deallocate_a(df)
919 subroutine write_iter_info(geom_iter, size, energy, maxdx, maxdf, coords)
920 integer,
intent(in) :: geom_iter
921 integer,
intent(in) :: size
922 real(real64),
intent(in) :: energy, maxdx, maxdf
923 real(real64),
intent(in) :: coords(size)
925 character(len=256) :: c_geom_iter, title, c_forces_iter
930 write(c_geom_iter,
'(a,i4.4)')
"go.", geom_iter
932 call io_mkdir(
'geom', g_opt%ions%namespace)
933 call g_opt%ions%write_xyz(
'geom/'//trim(c_geom_iter), comment = trim(title))
934 call g_opt%ions%write_xyz(
'./last')
936 if(g_opt%periodic_dim > 0)
then
937 call g_opt%ions%write_xyz(
'geom/'//trim(c_geom_iter), comment =
'Reduced coordinates', reduce_coordinates = .
true.)
939 g_opt%ions%pos, g_opt%ions%atom, g_opt%syst%gr, g_opt%syst%namespace)
942 if (g_opt%syst%outp%what(option__output__forces))
then
943 write(c_forces_iter,
'(a,i4.4)')
"forces.", geom_iter
944 if (
bitand(g_opt%syst%outp%how(option__output__forces), option__outputformat__bild) /= 0)
then
945 call g_opt%ions%write_bild_forces_file(
'forces', trim(c_forces_iter))
948 g_opt%ions%pos, g_opt%ions%atom, g_opt%syst%gr, g_opt%syst%namespace, total_forces=g_opt%ions%tot_force)
953 iunit =
io_open(trim(
'geom/optimization.log'), g_opt%syst%namespace, &
954 action =
'write', position =
'append')
956 if (geom_iter == 1)
then
958 write(iunit,
'(a10,5(5x,a20),a)')
'# iter',
'energy [' // trim(
units_abbrev(
units_out%energy)) //
']', &
963 ' alpha, beta, gamma [degrees]'
965 write(iunit,
'(a10,3(5x,a20))')
'# iter',
'energy [' // trim(
units_abbrev(
units_out%energy)) //
']', &
979 g_opt%ions%latt%alpha, g_opt%ions%latt%beta, g_opt%ions%latt%gamma
992 call messages_write(
"++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .
true.)
1001 if (g_opt%periodic_dim == 0)
then
1016 call messages_write(maxdx, fmt =
"f16.10,1x", print_units = .false., new_line = .
true.)
1019 call messages_write(
"++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .
true.)
1020 call messages_write(
"++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++", new_line = .
true.)
1032 real(real64),
intent(out) :: coords(:)
1034 integer :: iatom, idir, jdir, icoord
1035 real(real64) :: tmp_pos(gopt%dim), strain(g_opt%periodic_dim,g_opt%periodic_dim)
1041 if (
bitand(g_opt%type, go_ions) /= 0)
then
1042 do iatom = 1, gopt%ions%natoms
1043 if (gopt%fixed_atom == iatom) cycle
1044 if (gopt%ions%fixed(iatom)) cycle
1045 tmp_pos = gopt%ions%pos(1:gopt%dim, iatom)
1046 if (gopt%fixed_atom > 0) tmp_pos = tmp_pos - gopt%ions%pos(1:gopt%dim, gopt%fixed_atom)
1047 tmp_pos = gopt%ions%latt%cart_to_red(tmp_pos)
1048 do idir = 1, gopt%dim
1049 coords(icoord) = tmp_pos(idir)
1058 strain = matmul(gopt%ions%latt%rlattice, g_opt%inv_initial_rlattice)
1059 do idir = 1, g_opt%periodic_dim
1060 do jdir = idir, g_opt%periodic_dim
1061 coords(icoord) = strain(idir, jdir)
1069 do idir = 1, g_opt%periodic_dim
1070 coords(icoord) = norm2(gopt%ions%latt%rlattice(1:g_opt%periodic_dim, idir))/g_opt%initial_length(idir)
1081 subroutine to_grad(gopt, grad)
1083 real(real64),
intent(out) :: grad(:)
1085 integer :: iatom, idir, jdir, icoord
1086 real(real64) :: tmp_force(1:gopt%dim)
1092 if (
bitand(g_opt%type, go_ions) /= 0)
then
1093 do iatom = 1, gopt%ions%natoms
1094 if (gopt%fixed_atom == iatom) cycle
1095 if (gopt%ions%fixed(iatom)) cycle
1096 do idir = 1, gopt%dim
1097 if (abs(gopt%ions%atom(iatom)%c(idir)) <=
m_epsilon)
then
1098 tmp_force(idir) = -gopt%ions%tot_force(idir, iatom)
1102 if (gopt%fixed_atom > 0)
then
1103 tmp_force(idir) = tmp_force(idir) + gopt%ions%tot_force(idir, gopt%fixed_atom)
1106 tmp_force = gopt%ions%latt%cart_to_red(tmp_force)
1107 do idir = 1, gopt%dim
1108 grad(icoord) = tmp_force(idir)
1116 do idir = 1, g_opt%periodic_dim
1117 do jdir = idir, g_opt%periodic_dim
1118 grad(icoord) = -g_opt%cell_force(idir, jdir)
1126 do idir = 1, g_opt%periodic_dim
1127 grad(icoord) = -g_opt%cell_force(idir, 1)
1140 real(real64),
intent(in) :: coords(:)
1142 integer :: iatom, idir, jdir, icoord
1143 real(real64) :: tmp_pos(gopt%dim, gopt%ions%natoms), strain(g_opt%periodic_dim,g_opt%periodic_dim)
1151 if (
bitand(g_opt%type, go_ions) /= 0)
then
1152 do iatom = 1, gopt%ions%natoms
1153 if (gopt%fixed_atom == iatom) cycle
1154 if (gopt%ions%fixed(iatom)) cycle
1155 do idir = 1, gopt%dim
1156 tmp_pos(idir, iatom) = coords(icoord)
1161 do iatom = 1, gopt%ions%natoms
1162 tmp_pos(:, iatom) = gopt%ions%latt%cart_to_red(gopt%ions%pos(:, iatom))
1168 do idir = 1, g_opt%periodic_dim
1169 do jdir = idir, g_opt%periodic_dim
1170 strain(idir, jdir) = coords(icoord)
1178 gopt%ions%latt%rlattice = matmul(strain, g_opt%initial_rlattice)
1183 do idir = 1, g_opt%periodic_dim
1184 gopt%ions%latt%rlattice(1:g_opt%periodic_dim, idir) = coords(icoord) &
1185 * gopt%initial_rlattice(1:g_opt%periodic_dim, idir)
1192 call g_opt%syst%gr%symmetrizer%symmetrize_lattice_vectors(g_opt%periodic_dim, g_opt%initial_rlattice, &
1193 gopt%ions%latt%rlattice, gopt%symmetrize)
1194 call gopt%ions%update_lattice_vectors(gopt%ions%latt, gopt%symmetrize)
1198 if (
bitand(g_opt%type, go_ions) /= 0)
then
1200 do iatom = 1, gopt%ions%natoms
1201 if (gopt%fixed_atom == iatom) cycle
1202 if (gopt%ions%fixed(iatom)) cycle
1203 tmp_pos(:, iatom) = gopt%ions%latt%red_to_cart(tmp_pos(:, iatom))
1204 do idir = 1, gopt%dim
1205 if (abs(gopt%ions%atom(iatom)%c(idir)) <=
m_epsilon)
then
1206 gopt%ions%pos(idir, iatom) = tmp_pos(idir, iatom)
1209 if (gopt%fixed_atom > 0)
then
1210 gopt%ions%pos(:, iatom) = gopt%ions%pos(:, iatom) + gopt%ions%pos(:, gopt%fixed_atom)
1214 do iatom = 1, gopt%ions%natoms
1215 gopt%ions%pos(:, iatom) = gopt%ions%latt%red_to_cart(tmp_pos(:, iatom))
1219 if (gopt%symmetrize)
then
1220 call gopt%ions%symmetrize_atomic_coord()
1223 if (
debug%info)
then
1224 call gopt%ions%print_spacegroup()
1233 integer,
intent(in) :: geom_iter
1234 integer,
intent(in) :: size
1235 real(real64),
intent(in) :: energy, maxdx
1236 real(real64),
intent(in) :: coords(size)
subroutine init_(fromscratch)
Define which routines can be seen from the outside.
subroutine minimize_multidim(method, dim, x, step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum, ierr)
subroutine minimize_fire(dim, space_dim, x, step, tolgrad, maxiter, f, write_iter_info, en, ierr, mass, integrator)
Implementation of the Fast Inertial Relaxation Engine (FIRE)
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 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)
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 minmethod_nmsimplex
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
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
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
integer, parameter, public restart_type_dump
integer, parameter, public restart_type_load
subroutine, public restart_end(restart)
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.