Octopus
partial_charges.F90
Go to the documentation of this file.
1!! Copyright (C) 2015 X. Andrade
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use debug_oct_m
23 use global_oct_m
25 use ions_oct_m
26 use mesh_oct_m
28 use mpi_oct_m
31
32 implicit none
33
34 private
35
36 public :: &
39
40contains
41
42 !----------------------------------------------
43 subroutine partial_charges_calculate(mesh, st, ions, hirshfeld_charges)
44 class(mesh_t), intent(in) :: mesh
45 type(states_elec_t), intent(in) :: st
46 type(ions_t), intent(in) :: ions
47 real(real64), intent(out) :: hirshfeld_charges(:)
48
49 integer :: iatom
50 type(hirshfeld_t) :: hirshfeld
51
53 call profiling_in('PARTIAL_CHARGES')
54
55 call hirshfeld_init(hirshfeld, mesh, ions%namespace, ions%space, ions%latt, ions%atom, &
56 ions%natoms, ions%pos, st%d%spin_channels)
57
58 do iatom = 1, ions%natoms
59 call hirshfeld_charge(hirshfeld, ions%space, iatom, st%rho, hirshfeld_charges(iatom))
60 end do
61
62 call hirshfeld_end(hirshfeld)
63
64 call profiling_out('PARTIAL_CHARGES')
66
67 end subroutine partial_charges_calculate
68
69
71 subroutine partial_charges_compute_and_print_charges(mesh, st, ions, iunit)
72 class(mesh_t), intent(in) :: mesh
73 type(states_elec_t), intent(in) :: st
74 type(ions_t), intent(in) :: ions
75 integer, intent(in) :: iunit
76
77 integer :: iatom
78 real(real64), allocatable :: hirshfeld_charges(:)
79
81
82 safe_allocate(hirshfeld_charges(1:ions%natoms))
83
84 call partial_charges_calculate(mesh, st, ions, hirshfeld_charges)
85
86 if(mpi_world%is_root()) then
87
88 write(iunit,'(a)') 'Partial ionic charges'
89 write(iunit,'(a)') ' Ion Hirshfeld'
90
91 do iatom = 1, ions%natoms
92 write(iunit,'(i4,a10,f16.3)') iatom, trim(ions%atom(iatom)%species%get_label()), hirshfeld_charges(iatom)
93 end do
94
95 end if
96
97 safe_deallocate_a(hirshfeld_charges)
98
101
102
103end module partial_charges_oct_m
104
105!! Local Variables:
106!! mode: f90
107!! coding: utf-8
108!! End:
subroutine, public hirshfeld_init(this, mesh, namespace, space, latt, atom, natoms, pos, nspin)
Definition: hirshfeld.F90:166
subroutine, public hirshfeld_charge(this, space, iatom, density, charge)
Definition: hirshfeld.F90:277
subroutine, public hirshfeld_end(this)
Definition: hirshfeld.F90:258
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:270
subroutine, public partial_charges_compute_and_print_charges(mesh, st, ions, iunit)
Computes and write partial charges to a file.
subroutine, public partial_charges_calculate(mesh, st, ions, hirshfeld_charges)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552