40 use,
intrinsic :: iso_fortran_env
66 class(box_t),
pointer :: box
67 character(len=500) :: clist
68 character(len=15) :: lab
70 logical,
allocatable :: mesh_mask(:)
71 logical,
allocatable :: ions_mask(:)
72 type(local_write_t) :: writ
75 type(electrons_t),
pointer :: sys
76 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, ld_need_hm
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
212 if (iunit /= -1)
then
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
320 message(1) =
'Info: Local domains output includes local potential and/or local energy.'
324 sys%ext_partners, sys%st, time=
m_zero)
331 dir=trim(ldrestart_folder), mesh = sys%gr)
333 call restart_ld%end()
341 in_folder = folder(1:len_trim(folder)-1)
342 folder_index = index(in_folder,
'/', .
true.)
343 restart_folder = in_folder(1:folder_index)
345 restart_folder = folder
348 dir=trim(restart_folder), mesh = sys%gr)
351 do iter = l_start, l_end, l_step
354 write(in_folder,
'(a,i0.7,a)') folder(folder_index+1:len_trim(folder)-1),iter,
"/"
355 write(filename,
'(a,a,a)') trim(in_folder), trim(basename)
358 if (l_start /= l_end)
then
362 write(frmt,
'(a,i0,a)')
"(a,i0.",10-len_trim(basename),
")"
363 write(filename, fmt=trim(frmt)) trim(basename), iter
366 write(filename,
'(a,a,a,a)') trim(in_folder),
"/", trim(basename)
375 call restart%read_mesh_function(trim(filename), sys%gr, sys%st%rho(:,1), err)
378 write(
message(1),*)
'While reading density: "', trim(filename),
'", error code:', err
383 if ((iter == l_start .and. .not. ldrestart) .or. ldupdate)
then
390 loc_domains(id)%mesh_mask, sys%gr, sys%st, sys%hm, sys%ks, sys%ions, &
391 sys%ext_partners, kick, iter, l_start, ldoverwrite)
402 safe_deallocate_a(loc_domains)
404 message(1) =
'Info: Exiting local domains'
415 subroutine local_init(space, mesh, ions, nd, loc_domains)
416 class(
space_t),
intent(in) :: space
417 class(
mesh_t),
intent(in) :: mesh
418 type(
ions_t),
intent(in) :: ions
419 integer,
intent(out) :: nd
466 safe_allocate(loc_domains(1:nd))
469 safe_allocate(loc_domains(id)%mesh_mask(1:mesh%np))
470 safe_allocate(loc_domains(id)%ions_mask(1:ions%natoms))
488 class(
box_t),
pointer :: box
492 if (domain%dshape /= bader)
then
494 safe_deallocate_p(box)
497 safe_deallocate_a(domain%mesh_mask)
498 safe_deallocate_a(domain%ions_mask)
506 class(
space_t),
intent(in) :: space
507 type(
ions_t),
intent(in) :: ions
508 type(
block_t),
intent(in) :: blk
509 integer,
intent(in) :: row
513 real(real64) :: radius, center(space%dim), half_length(space%dim)
520 select case (domain%dshape)
527 do ia = 1, ions%natoms
529 call minimum_box%add_box(
box_sphere_t(space%dim, ions%pos(:, ia), radius, namespace))
532 domain%box => minimum_box
535 do ia = 1, ions%natoms
536 if (domain%box%contains_point(ions%pos(:, ia)) .and. .not.
loct_isinstringlist(ia, domain%clist))
then
538 if (ic <= 20)
write(
message(ic),
'(a,a,I0,a,a)')
'Atom: ',trim(ions%atom(ia)%species%get_label()), ia, &
539 ' is inside the union box BUT not in list: ', trim(domain%clist)
545 message(1) =
' THIS COULD GIVE INCORRECT RESULTS '
548 message(1) =
' AT LEAST 19 ATOMS ARE NOT PRESENT IN THE LIST'
560 domain%box =>
box_sphere_t(space%dim, center, radius, namespace)
572 domain%box =>
box_cylinder_t(space%dim, center, ions%latt%rlattice, radius, half_length(1)*
m_two, namespace)
597 class(
space_t),
intent(in) :: space
598 class(
mesh_t),
intent(in) :: mesh
600 integer,
intent(in) :: nd
603 integer :: id, ip, iunit, ierr
604 character(len=MAX_PATH_LEN) :: dirname
605 character(len=140),
allocatable :: lines(:)
606 real(real64),
allocatable :: ff(:,:)
611 dirname =
"./restart/ld"
612 write(
message(1),
'(a,a)')
'Info: Writing restart info to ', trim(dirname)
617 safe_allocate(lines(1:nd+2))
618 write(lines(1),
'(a,1x,i5)')
'Number of local domains =', nd
619 write(lines(2),
'(a3,1x,a15,1x,a5,1x)')
'#id',
'label',
'shape'
621 write(lines(id+2),
'(i3,1x,a15,1x,i5)') id, trim(loc_domains(id)%lab), loc_domains(id)%dshape
623 iunit = restart%open(
"ldomains.info")
624 call restart%write(iunit, lines, nd+2, ierr)
626 safe_deallocate_a(lines)
628 safe_allocate(ff(1:mesh%np, 1))
632 if (loc_domains(id)%mesh_mask(ip)) ff(ip, 1) = ff(ip, 1) + 2**real(id, real64)
635 call restart%write_mesh_function(
"ldomains", mesh, ff(1:mesh%np, 1), ierr)
636 safe_deallocate_a(ff)
645 class(
space_t),
intent(in) :: space
646 class(
mesh_t),
intent(in) :: mesh
647 type(
ions_t),
intent(in) :: ions
648 integer,
intent(in) :: nd
652 integer :: id, ip, iunit, ierr
653 character(len=31) :: line(1), tmp
654 real(real64),
allocatable :: mask(:)
658 message(1) =
'Info: Reading mesh points inside each local domain'
661 safe_allocate(mask(1:mesh%np))
664 iunit = restart%open(
"ldomains.info", status=
'old')
665 call restart%read(iunit, line, 1, ierr)
666 read(line(1),
'(a25,1x,i5)') tmp, ierr
667 call restart%close(iunit)
669 call restart%read_mesh_function(
"ldomains", mesh, mask, ierr)
672 loc_domains(id)%mesh_mask = .false.
674 if (
bitand(int(mask(ip)), 2**id) /= 0) loc_domains(id)%mesh_mask(ip) = .
true.
678 call local_ions_mask(loc_domains(id)%mesh_mask, ions, mesh, loc_domains(id)%ions_mask)
681 safe_deallocate_a(mask)
688 class(
space_t),
intent(in) :: space
689 class(
mesh_t),
intent(in) :: mesh
691 integer,
intent(in) :: nd
694 real(real64),
intent(in) :: ff(:)
696 integer(int64) :: how
697 integer :: id, ip, ierr
699 character(len=MAX_PATH_LEN) :: filename
700 real(real64),
allocatable :: dble_domain_map(:), domain_mesh(:)
705 message(1) =
'Info: Assigning mesh points inside each local domain'
709 if (any(loc_domains(:)%dshape == bader))
then
710 if (mesh%parallel_in_domains)
then
711 write(
message(1),
'(a)')
'Bader volumes can only be computed in serial'
719 select case (loc_domains(id)%dshape)
727 call local_ions_mask(loc_domains(id)%mesh_mask, ions, mesh, loc_domains(id)%ions_mask)
730 if (any(loc_domains(:)%dshape == bader))
then
741 safe_allocate(dble_domain_map(1:mesh%np))
742 safe_allocate(domain_mesh(1:mesh%np))
748 if (loc_domains(id)%mesh_mask(ip))
then
749 dble_domain_map(ip) = real(id, real64)
750 domain_mesh(ip) = domain_mesh(ip) + dble_domain_map(ip)
754 write(filename,
'(a,a)')
'domain.', trim(loc_domains(id)%lab)
756 pos=ions%pos, atoms=ions%atom)
760 pos=ions%pos, atoms=ions%atom)
762 safe_deallocate_a(dble_domain_map)
763 safe_deallocate_a(domain_mesh)
770 subroutine create_basins(namespace, space, mesh, ions, ff, basins)
772 class(
space_t),
intent(in) :: space
773 class(
mesh_t),
intent(in) :: mesh
774 type(
ions_t),
intent(in) :: ions
775 real(real64),
intent(in) :: ff(:)
776 type(
basins_t),
intent(inout) :: basins
778 logical :: use_atomic_radii
779 integer :: ip, ierr, ia, is, rankmin
780 integer(int64) :: how
781 real(real64) :: bader_threshold, dmin, dd
782 integer,
allocatable :: ion_map(:)
783 real(real64),
allocatable :: ff2(:,:), ffs(:)
795 call parse_variable(namespace,
'LDBaderThreshold', 0.01_real64, bader_threshold)
798 safe_allocate(ff2(1:mesh%np,1))
799 ff2(1:mesh%np,1) = ff(1:mesh%np)
802 safe_allocate(ffs(1:mesh%np))
803 do ia = 1, ions%natoms
805 ions%pos(:, ia), mesh, ffs)
806 do is = 1, ubound(ff2, dim=2)
807 ff2(1:mesh%np, is) = ff2(1:mesh%np, is) - ffs(1:mesh%np)
810 safe_deallocate_a(ffs)
813 call basins_analyze(basins, namespace, mesh, ff2(:,1), ff2, bader_threshold)
821 call dio_function_output(how,
'local.general',
'basinsmap', namespace, space, mesh, real(basins%map(1:mesh%np), real64) , &
822 unit_one, ierr, pos=ions%pos, atoms=ions%atom)
823 call dio_function_output(how,
'local.general',
'dens_ff2', namespace, space, mesh, ff2(:,1), &
824 unit_one, ierr, pos=ions%pos, atoms=ions%atom)
826 safe_deallocate_a(ff2)
837 if (use_atomic_radii)
then
838 safe_allocate(ion_map(1:ions%natoms))
840 do ia = 1, ions%natoms
846 if (all(ion_map(:) /= basins%map(ip)))
then
848 do ia = 1, ions%natoms
849 dd = norm2(mesh%x(:, ip) - ions%pos(:, ia))
850 if (dd <= ions%atom(ia)%species%get_vdw_radius()) basins%map(ip) = ion_map(ia)
855 safe_deallocate_a(ion_map)
864 class(
mesh_t),
intent(in) :: mesh
868 domain%mesh_mask = domain%box%contains_points(mesh%np, mesh%x_t)
876 type(
basins_t),
intent(in) :: basins
877 class(
mesh_t),
intent(in) :: mesh
878 type(
ions_t),
intent(in) :: ions
880 integer :: ia, ib, ip, ix, n_basins, rankmin
882 integer,
allocatable :: domain_map(:)
888 do ia = 1, ions%natoms
892 safe_allocate(domain_map(1:n_basins))
897 do ia = 1, ions%natoms
901 domain_map(ib) = basins%map(ix)
907 domain%mesh_mask(ip) = any(domain_map(1:n_basins) == basins%map(ip))
910 safe_deallocate_a(domain_map)
918 logical,
intent(in) :: mesh_mask(:)
919 type(
ions_t),
intent(in) :: ions
920 class(
mesh_t),
intent(in) :: mesh
921 logical,
intent(out) :: ions_mask(:)
924 integer,
allocatable :: mask_tmp(:)
930 safe_allocate(mask_tmp(1:ions%natoms))
933 do ia = 1, ions%natoms
935 if (rankmin /= mesh%mpi_grp%rank) cycle
936 if (mesh_mask(ix)) mask_tmp(ia) = 1
939 call mesh%allreduce(mask_tmp, ions%natoms)
940 ions_mask = mask_tmp == 1
942 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 hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
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)
logical function, public local_write_check_hm(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)
System information (time, memory, sysname)
logical function, public loct_isinstringlist(a, s)
subroutine, public loct_progress_bar(a, maxcount)
A wrapper around the progress bar, such that it can be silenced without needing to dress the call wit...
character(kind=c_char, len=1) function, dimension(len_trim(f_string)+1), private string_f_to_c(f_string)
convert a Fortran string to a C string
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)
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)
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 parse_block_string(blk, l, c, res, convert_to_c)
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)
integer, parameter, public restart_undefined
integer, parameter, public restart_type_dump
integer, parameter, public restart_type_load
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.