40 use,
intrinsic :: iso_fortran_env
65 class(box_t),
pointer :: box
66 character(len=500) :: clist
67 character(len=15) :: lab
69 logical,
allocatable :: mesh_mask(:)
70 logical,
allocatable :: ions_mask(:)
71 type(local_write_t) :: writ
74 type(electrons_t),
pointer :: sys
75 integer,
parameter :: BADER = 9
102 safe_deallocate_p(sys)
118 integer :: err, iter, l_start, l_end, l_step, id
119 integer :: ia, n_spec_def, read_data, iunit, ispec
120 integer ::
length, folder_index
121 real(real64) :: default_dt, dt
122 character(len=MAX_PATH_LEN) :: filename, folder, folder_default, radiifile
123 character(len=MAX_PATH_LEN) :: in_folder, restart_folder, basename, frmt, ldrestart_folder
124 character(len=15) :: lab
125 logical :: iterate, ldupdate, ldoverwrite, ldrestart
131 message(1) =
'Info: Creating local domains'
135 write(folder_default,
'(a)')
'restart/gs/'
146 if (index(folder,
'/', back = .
true.) /= len_trim(folder))
then
147 folder = trim(folder)//
'/'
153 write(
message(1),
'(a)')
'Input: TDTimeStep must be positive.'
154 write(
message(2),
'(a)')
'Input: TDTimeStep reset to 0. Check input file.'
167 if (basename ==
" ") basename =
""
169 length = len_trim(basename)
172 basename = trim(basename(1:
length-4))
204 if (trim(radiifile) /=
"default")
then
206 if (n_spec_def > 0) n_spec_def = n_spec_def - 1
209 do ia = 1, sys%ions%nspecies
214 write(
message(1),
'(a,a)')
'Redefining def_rsize from file:', trim(radiifile)
218 default_file:
do ispec = 1, n_spec_def
220 if (trim(lab) == trim(sys%ions%species(ia)%s%get_label()))
then
221 select type(spec=>sys%ions%species(ia)%s)
230 write(
message(1),
'(a,a)')
'Octopus could not open then file:', trim(radiifile)
247 write(folder_default,
'(a)')
'ld.general'
259 if (index(ldrestart_folder,
'/', .
true.) /= len_trim(ldrestart_folder))
then
260 write(ldrestart_folder,
'(a,a1)') trim(ldrestart_folder),
'/'
301 message(1) =
'Info: Computing local multipoles'
305 call local_init(sys%space, sys%gr, sys%ions, nd, loc_domains)
308 if (any(loc_domains(:)%dshape == bader))
then
322 dir=trim(ldrestart_folder), mesh = sys%gr)
332 in_folder = folder(1:len_trim(folder)-1)
333 folder_index = index(in_folder,
'/', .
true.)
334 restart_folder = in_folder(1:folder_index)
336 restart_folder = folder
339 dir=trim(restart_folder), mesh = sys%gr)
342 do iter = l_start, l_end, l_step
345 write(in_folder,
'(a,i0.7,a)') folder(folder_index+1:len_trim(folder)-1),iter,
"/"
346 write(filename,
'(a,a,a)') trim(in_folder), trim(basename)
349 if (l_start /= l_end)
then
353 write(frmt,
'(a,i0,a)')
"(a,i0.",10-len_trim(basename),
")"
354 write(filename, fmt=trim(frmt)) trim(basename), iter
357 write(filename,
'(a,a,a,a)') trim(in_folder),
"/", trim(basename)
369 write(
message(1),*)
'While reading density: "', trim(filename),
'", error code:', err
374 if ((iter == l_start .and. .not. ldrestart) .or. ldupdate)
then
381 loc_domains(id)%mesh_mask, sys%gr, sys%st, sys%hm, sys%ks, sys%ions, &
382 sys%ext_partners, kick, iter, l_start, ldoverwrite)
393 safe_deallocate_a(loc_domains)
395 message(1) =
'Info: Exiting local domains'
406 subroutine local_init(space, mesh, ions, nd, loc_domains)
407 class(
space_t),
intent(in) :: space
408 class(
mesh_t),
intent(in) :: mesh
409 type(
ions_t),
intent(in) :: ions
410 integer,
intent(out) :: nd
457 safe_allocate(loc_domains(1:nd))
460 safe_allocate(loc_domains(id)%mesh_mask(1:mesh%np))
461 safe_allocate(loc_domains(id)%ions_mask(1:ions%natoms))
479 class(
box_t),
pointer :: box
483 if (domain%dshape /= bader)
then
485 safe_deallocate_p(box)
488 safe_deallocate_a(domain%mesh_mask)
489 safe_deallocate_a(domain%ions_mask)
497 class(
space_t),
intent(in) :: space
498 type(
ions_t),
intent(in) :: ions
500 integer,
intent(in) :: row
504 real(real64) :: radius, center(space%dim), half_length(space%dim)
511 select case (domain%dshape)
518 do ia = 1, ions%natoms
520 call minimum_box%add_box(
box_sphere_t(space%dim, ions%pos(:, ia), radius, namespace))
523 domain%box => minimum_box
526 do ia = 1, ions%natoms
527 if (domain%box%contains_point(ions%pos(:, ia)) .and. .not.
loct_isinstringlist(ia, domain%clist))
then
529 if (ic <= 20)
write(
message(ic),
'(a,a,I0,a,a)')
'Atom: ',trim(ions%atom(ia)%species%get_label()), ia, &
530 ' is inside the union box BUT not in list: ', trim(domain%clist)
536 message(1) =
' THIS COULD GIVE INCORRECT RESULTS '
539 message(1) =
' AT LEAST 19 ATOMS ARE NOT PRESENT IN THE LIST'
551 domain%box =>
box_sphere_t(space%dim, center, radius, namespace)
563 domain%box =>
box_cylinder_t(space%dim, center, ions%latt%rlattice, radius, half_length(1)*
m_two, namespace)
589 class(
mesh_t),
intent(in) :: mesh
591 integer,
intent(in) :: nd
594 integer :: id, ip, iunit, ierr
595 character(len=MAX_PATH_LEN) :: dirname
596 character(len=140),
allocatable :: lines(:)
597 real(real64),
allocatable :: ff(:,:)
602 dirname =
"./restart/ld"
603 write(
message(1),
'(a,a)')
'Info: Writing restart info to ', trim(dirname)
608 safe_allocate(lines(1:nd+2))
609 write(lines(1),
'(a,1x,i5)')
'Number of local domains =', nd
610 write(lines(2),
'(a3,1x,a15,1x,a5,1x)')
'#id',
'label',
'shape'
612 write(lines(id+2),
'(i3,1x,a15,1x,i5)') id, trim(loc_domains(id)%lab), loc_domains(id)%dshape
617 safe_deallocate_a(lines)
619 safe_allocate(ff(1:mesh%np, 1))
623 if (loc_domains(id)%mesh_mask(ip)) ff(ip, 1) = ff(ip, 1) + 2**real(id, real64)
627 safe_deallocate_a(ff)
636 class(
space_t),
intent(in) :: space
637 class(
mesh_t),
intent(in) :: mesh
638 type(
ions_t),
intent(in) :: ions
639 integer,
intent(in) :: nd
643 integer :: id, ip, iunit, ierr
644 character(len=31) :: line(1), tmp
645 real(real64),
allocatable :: mask(:)
649 message(1) =
'Info: Reading mesh points inside each local domain'
652 safe_allocate(mask(1:mesh%np))
655 iunit =
restart_open(restart,
"ldomains.info", status=
'old')
657 read(line(1),
'(a25,1x,i5)') tmp, ierr
663 loc_domains(id)%mesh_mask = .false.
665 if (
bitand(int(mask(ip)), 2**id) /= 0) loc_domains(id)%mesh_mask(ip) = .
true.
669 call local_ions_mask(loc_domains(id)%mesh_mask, ions, mesh, loc_domains(id)%ions_mask)
672 safe_deallocate_a(mask)
680 class(
mesh_t),
intent(in) :: mesh
681 type(
ions_t),
intent(in) :: ions
682 integer,
intent(in) :: nd
685 real(real64),
intent(in) :: ff(:)
687 integer(int64) :: how
688 integer :: id, ip, ierr
690 character(len=MAX_PATH_LEN) :: filename
691 real(real64),
allocatable :: dble_domain_map(:), domain_mesh(:)
696 message(1) =
'Info: Assigning mesh points inside each local domain'
700 if (any(loc_domains(:)%dshape == bader))
then
701 if (mesh%parallel_in_domains)
then
702 write(
message(1),
'(a)')
'Bader volumes can only be computed in serial'
710 select case (loc_domains(id)%dshape)
718 call local_ions_mask(loc_domains(id)%mesh_mask, ions, mesh, loc_domains(id)%ions_mask)
721 if (any(loc_domains(:)%dshape == bader))
then
732 safe_allocate(dble_domain_map(1:mesh%np))
733 safe_allocate(domain_mesh(1:mesh%np))
739 if (loc_domains(id)%mesh_mask(ip))
then
740 dble_domain_map(ip) = real(id, real64)
741 domain_mesh(ip) = domain_mesh(ip) + dble_domain_map(ip)
745 write(filename,
'(a,a)')
'domain.', trim(loc_domains(id)%lab)
747 pos=ions%pos, atoms=ions%atom)
751 pos=ions%pos, atoms=ions%atom)
753 safe_deallocate_a(dble_domain_map)
754 safe_deallocate_a(domain_mesh)
761 subroutine create_basins(namespace, space, mesh, ions, ff, basins)
763 class(
space_t),
intent(in) :: space
764 class(
mesh_t),
intent(in) :: mesh
765 type(
ions_t),
intent(in) :: ions
766 real(real64),
intent(in) :: ff(:)
767 type(
basins_t),
intent(inout) :: basins
769 logical :: use_atomic_radii
770 integer :: ip, ierr, ia, is, rankmin
771 integer(int64) :: how
772 real(real64) :: bader_threshold, dmin, dd
773 integer,
allocatable :: ion_map(:)
774 real(real64),
allocatable :: ff2(:,:), ffs(:)
786 call parse_variable(namespace,
'LDBaderThreshold', 0.01_real64, bader_threshold)
789 safe_allocate(ff2(1:mesh%np,1))
790 ff2(1:mesh%np,1) = ff(1:mesh%np)
793 safe_allocate(ffs(1:mesh%np))
794 do ia = 1, ions%natoms
796 ions%pos(:, ia), mesh, ffs)
797 do is = 1, ubound(ff2, dim=2)
798 ff2(1:mesh%np, is) = ff2(1:mesh%np, is) - ffs(1:mesh%np)
801 safe_deallocate_a(ffs)
804 call basins_analyze(basins, namespace, mesh, ff2(:,1), ff2, bader_threshold)
812 call dio_function_output(how,
'local.general',
'basinsmap', namespace, space, mesh, real(basins%map(1:mesh%np), real64) , &
813 unit_one, ierr, pos=ions%pos, atoms=ions%atom)
814 call dio_function_output(how,
'local.general',
'dens_ff2', namespace, space, mesh, ff2(:,1), &
815 unit_one, ierr, pos=ions%pos, atoms=ions%atom)
817 safe_deallocate_a(ff2)
828 if (use_atomic_radii)
then
829 safe_allocate(ion_map(1:ions%natoms))
831 do ia = 1, ions%natoms
837 if (all(ion_map(:) /= basins%map(ip)))
then
839 do ia = 1, ions%natoms
840 dd = norm2(mesh%x(ip, :) - ions%pos(:, ia))
841 if (dd <= ions%atom(ia)%species%get_vdw_radius()) basins%map(ip) = ion_map(ia)
846 safe_deallocate_a(ion_map)
855 class(
mesh_t),
intent(in) :: mesh
859 domain%mesh_mask = domain%box%contains_points(mesh%np, mesh%x)
867 type(
basins_t),
intent(in) :: basins
868 class(
mesh_t),
intent(in) :: mesh
869 type(
ions_t),
intent(in) :: ions
871 integer :: ia, ib, ip, ix, n_basins, rankmin
873 integer,
allocatable :: domain_map(:)
879 do ia = 1, ions%natoms
883 safe_allocate(domain_map(1:n_basins))
888 do ia = 1, ions%natoms
892 domain_map(ib) = basins%map(ix)
898 domain%mesh_mask(ip) = any(domain_map(1:n_basins) == basins%map(ip))
901 safe_deallocate_a(domain_map)
909 logical,
intent(in) :: mesh_mask(:)
910 type(
ions_t),
intent(in) :: ions
911 class(
mesh_t),
intent(in) :: mesh
912 logical,
intent(out) :: ions_mask(:)
915 integer,
allocatable :: mask_tmp(:)
921 safe_allocate(mask_tmp(1:ions%natoms))
924 do ia = 1, ions%natoms
926 if (rankmin /= mesh%mpi_grp%rank) cycle
927 if (mesh_mask(ix)) mask_tmp(ia) = 1
930 if (mesh%parallel_in_domains)
then
931 call mesh%allreduce(mask_tmp, ions%natoms)
933 ions_mask = mask_tmp == 1
935 safe_deallocate_a(mask_tmp)
subroutine create_basins(namespace, space, mesh, ions, ff, basins)
program oct_local_multipoles
subroutine local_restart_read(space, mesh, ions, nd, loc_domains, restart)
subroutine local_ions_mask(mesh_mask, ions, mesh, ions_mask)
subroutine local_restart_write(namespace, space, mesh, mc, nd, loc_domains)
Write restart files for local domains.
subroutine local_inside_domain(space, mesh, ions, nd, loc_domains, namespace, ff)
subroutine bader_domain_create_mask(domain, basins, mesh, ions)
subroutine local_read_from_block(domain, space, ions, blk, row, namespace)
subroutine local_domains()
Reads a global grid and density and make local domains This is a high-level interface that reads the ...
subroutine local_init(space, mesh, ions, nd, loc_domains)
Initialize local_domain_t variable, allocating variable and reading parameters from input file.
subroutine local_end(domain)
Ending local_domain_t variable, allocating variable and reading parameters from input file.
subroutine box_domain_create_mask(domain, mesh)
subroutine, public basins_init(this, namespace, mesh)
subroutine, public basins_analyze(this, namespace, mesh, f, rho, threshold)
subroutine, public basins_end(this)
Module, implementing a factory for boxes.
integer, parameter, public sphere
integer, parameter, public minimum
integer, parameter, public parallelepiped
integer, parameter, public cylinder
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
type(debug_t), save, public debug
real(real64), parameter, public m_two
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
real(real64), parameter, public m_zero
subroutine, public global_init(communicator)
Initialise Octopus.
integer, parameter, public length
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
subroutine, public io_init(defaults)
If the argument defaults is present and set to true, then the routine will not try to read anything f...
subroutine, public io_close(iunit, grp)
subroutine, public io_end()
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
subroutine, public kick_end(kick)
subroutine, public kick_init(kick, namespace, space, kpoints, nspin)
subroutine, public local_write_init(writ, namespace, lab, iter, dt)
subroutine, public local_write_end(writ)
subroutine, public local_write_iter(writ, namespace, space, lab, ions_mask, mesh_mask, mesh, st, hm, ks, ions, ext_partners, kick, iter, l_start, ldoverwrite)
logical function, public loct_isinstringlist(a, s)
This module defines the meshes, which are used in Octopus.
integer function, public mesh_nearest_point(mesh, pos, dmin, rankmin)
Returns the index of the point which is nearest to a given vector position pos.
subroutine, public messages_end()
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
character(len=512), private msg
subroutine, public messages_init(output_dir)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public print_date(str)
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
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)
type(mpi_grp_t), public mpi_world
This module handles the communicators for the various parallelization strategies.
type(namespace_t), public global_namespace
subroutine, public parser_init()
Initialise the Octopus parser.
subroutine, public parser_end()
End the Octopus parser.
integer function, public parse_block(namespace, name, blk, check_varinfo_)
subroutine, public profiling_end(namespace)
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
subroutine, public read_from_default_file(iunit, read_data, spec)
subroutine, public restart_module_init(namespace)
subroutine, public restart_read(restart, iunit, lines, nlines, ierr)
integer, parameter, public restart_undefined
subroutine, public restart_close(restart, iunit)
Close a file previously opened with restart_open.
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
integer, parameter, public restart_type_dump
subroutine, public restart_write(restart, iunit, lines, nlines, ierr)
integer function, public restart_open(restart, filename, status, position, silent)
Open file "filename" found inside the current restart directory. Depending on the type of restart,...
integer, parameter, public restart_type_load
subroutine, public restart_end(restart)
subroutine, public species_get_long_range_density(species, namespace, space, latt, pos, mesh, rho, sphere_inout, nlr_x)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
subroutine, public unit_system_init(namespace)
type(unit_system_t), public units_inp
the units systems for reading and writing
type(unit_t), public unit_one
some special units required for particular quantities
This module is intended to contain simple general-purpose utility functions and procedures.
subroutine, public print_header()
This subroutine prints the logo followed by information about the compilation and the system....
Class implementing a cylinder box. The cylinder axis is always along the first direction defined by t...
class to tell whether a point is inside or outside
Class implementing a parallelepiped box. Currently this is restricted to a rectangular cuboid (all th...
Class implementing a spherical box.
Class implementing a box that is an union other boxes.
Class describing the electron system.
Describes mesh distribution to nodes.
Stores all communicators and groups.