Octopus
modelmb_density_matrix.F90
Go to the documentation of this file.
1!! Copyright (C) 2009 N. Helbig and M. Verstraete
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 batch_oct_m
23 use comm_oct_m
24 use debug_oct_m
25 use global_oct_m
27 use io_oct_m
28 use index_oct_m
29 use, intrinsic :: iso_fortran_env
31 use mesh_oct_m
35 use mpi_oct_m
38 use parser_oct_m
40 use space_oct_m
42
43 implicit none
44
45 private
46
47 public :: &
53
55 private
56 integer :: ndensmat_to_calculate
57 character(len=200) :: dirname
58 character(80), allocatable :: labels(:)
59 integer, allocatable :: particle_kept(:)
60 integer, allocatable :: nnatorb_prt(:)
61 end type modelmb_denmat_t
62
63contains
64
65 subroutine modelmb_density_matrix_init(dir, namespace, st, denmat)
66 character(len=*), intent(in) :: dir
67 type(namespace_t), intent(in) :: namespace
68 type(states_elec_t), intent(in) :: st
69 type(modelmb_denmat_t), intent(out) :: denmat
70
71 integer :: ncols, ipart
72 type(block_t) :: blk
73
75
76 !%Variable DensitytoCalc
77 !%Type block
78 !%Section States::ModelMB
79 !%Description
80 !% Choice of which particle density (event. matrices) will be calculated and output, in the
81 !% modelmb particles scheme.
82 !%
83 !% <tt>%DensitytoCalc
84 !% <br>&nbsp;&nbsp; "proton" | 1 | 10
85 !% <br>&nbsp;&nbsp; "electron" | 2 | 15
86 !% <br>%</tt>
87 !%
88 !% would ask octopus to calculate the density matrix corresponding to the 1st
89 !% particle (whose coordinates correspond to dimensions 1 to ndim_modelmb),
90 !% which is an proton, then that corresponding to the 2nd particle
91 !% (electron with dimensions ndim_modelmb+1 to 2*ndim_modelmb), printing
92 !% 10 natural orbitals for the first and 15 for the second.
93 !%
94 !% <tt>%DensitytoCalc
95 !% <br>&nbsp;&nbsp; "proton" | 1 | -1
96 !% <br>&nbsp;&nbsp; "electron" | 2 | -1
97 !% <br>%</tt>
98 !%
99 !% would ask octopus to print out just the densities for particles 1 and 2
100 !% without any density matrix output.
101 !%
102 !%End
103
104 call messages_obsolete_variable(namespace, 'DensityMatrixtoCalc', 'DensitytoCalc')
105 call messages_obsolete_variable(namespace, 'DensitiestoCalc', 'DensitytoCalc')
106
107 if (parse_block(namespace, 'DensitytoCalc', blk) /= 0) then
108 message(1) = 'To print out density (matrices), you must specify the DensitytoCalc block in input'
109 call messages_fatal(1, namespace=namespace)
110 end if
111
112 ncols = parse_block_cols(blk, 0)
113 if (ncols /= 3) then
114 call messages_input_error(namespace, "DensitytoCalc")
115 end if
116 denmat%ndensmat_to_calculate=parse_block_n(blk)
117 if (denmat%ndensmat_to_calculate < 0 .or. &
118 denmat%ndensmat_to_calculate > st%modelmbparticles%nparticle) then
119 call messages_input_error(namespace, "DensitytoCalc")
120 end if
121
122 safe_allocate(denmat%labels(1:denmat%ndensmat_to_calculate))
123 safe_allocate(denmat%particle_kept(1:denmat%ndensmat_to_calculate))
124 safe_allocate(denmat%nnatorb_prt(1:denmat%ndensmat_to_calculate))
125
126 do ipart=1,denmat%ndensmat_to_calculate
127 call parse_block_string(blk, ipart-1, 0, denmat%labels(ipart))
128 call parse_block_integer(blk, ipart-1, 1, denmat%particle_kept(ipart))
129 call parse_block_integer(blk, ipart-1, 2, denmat%nnatorb_prt(ipart))
130
131 write (message(1),'(a,a)') 'labels_densmat = ', denmat%labels(ipart)
132 write (message(2),'(a,i6)') 'particle_kept_densmat = ', denmat%particle_kept(ipart)
133 write (message(3),'(a,i6)') 'nnatorb_prt_densmat = ', denmat%nnatorb_prt(ipart)
134 call messages_info(3, namespace=namespace)
135 end do
136 call parse_block_end(blk)
137 ! END reading in of input var block DensitytoCalc
138
139 denmat%dirname = trim(dir)
140
142 end subroutine modelmb_density_matrix_init
143
144 ! ---------------------------------------------------------
145 subroutine modelmb_density_matrix_end(this)
146 type(modelmb_denmat_t), intent(inout) :: this
150 safe_deallocate_a(this%labels)
151 safe_deallocate_a(this%particle_kept)
152 safe_deallocate_a(this%nnatorb_prt)
155 end subroutine modelmb_density_matrix_end
156
157
158#include "undef.F90"
159#include "real.F90"
160#include "modelmb_density_matrix_inc.F90"
161#include "undef.F90"
162
163#include "complex.F90"
164#include "modelmb_density_matrix_inc.F90"
165#include "undef.F90"
166
168
169
170!! Local Variables:
171!! mode: f90
172!! coding: utf-8
173!! End:
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements the index, used for the mesh points.
Definition: index.F90:122
Definition: io.F90:114
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:116
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
general module for modelmb particles (eg 4 electrons in 1D equiv to 1 in 4D). Also calculate differen...
subroutine, public zmodelmb_density_matrix_write(space, mesh, st, wf, mm, denmat, namespace)
subroutine, public dmodelmb_density_matrix_write(space, mesh, st, wf, mm, denmat, namespace)
subroutine, public modelmb_density_matrix_init(dir, namespace, st, denmat)
subroutine, public modelmb_density_matrix_end(this)
Some general things and nomenclature:
Definition: par_vec.F90:171
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618