Octopus
target_gstransformation.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
24 use global_oct_m
25 use grid_oct_m
27 use io_oct_m
28 use ions_oct_m
30 use, intrinsic :: iso_fortran_env
35 use output_oct_m
37 use parser_oct_m
40 use space_oct_m
44 use td_oct_m
45
46 implicit none
47
48 private
49 public :: &
54
55
56contains
57
58 ! ----------------------------------------------------------------------
60 subroutine target_init_gstransformation(gr, namespace, space, tg, td, restart, kpoints)
61 type(grid_t), intent(in) :: gr
62 type(namespace_t), intent(in) :: namespace
63 class(space_t), intent(in) :: space
64 type(target_t), intent(inout) :: tg
65 type(td_t), intent(in) :: td
66 type(restart_t), intent(inout) :: restart
67 type(kpoints_t), intent(in) :: kpoints
68
70
71 message(1) = 'Info: Using Superposition of States for TargetOperator'
72 call messages_info(1, namespace=namespace)
73
74 tg%move_ions = ion_dynamics_ions_move(td%ions_dyn)
75 tg%dt = td%dt
76
77 !%Variable OCTTargetTransformStates
78 !%Type block
79 !%Default no
80 !%Section Calculation Modes::Optimal Control
81 !%Description
82 !% If <tt>OCTTargetOperator = oct_tg_gstransformation</tt>, you must specify a
83 !% <tt>OCTTargetTransformStates</tt> block, in order to specify which linear
84 !% combination of the states present in <tt>restart/gs</tt> is used to
85 !% create the target state.
86 !%
87 !% The syntax is the same as the <tt>TransformStates</tt> block.
88 !%End
89 call states_elec_transform(tg%st, namespace, space, restart, gr, kpoints, prefix = "OCTTarget")
90
91 if (.not. parse_is_defined(namespace, 'OCTTargetTransformStates')) then
92 message(1) = 'If "OCTTargetOperator = oct_tg_superposition", then you must'
93 message(2) = 'supply an "OCTTargetTransformStates" block to create the superposition.'
94 call messages_fatal(2, namespace=namespace)
95 end if
96 call density_calc(tg%st, gr, tg%st%rho)
97
99 end subroutine target_init_gstransformation
100
101
102 ! ----------------------------------------------------------------------
103 subroutine target_output_gstransformation(tg, namespace, space, gr, dir, ions, hm, outp)
104 type(target_t), intent(in) :: tg
105 type(namespace_t), intent(in) :: namespace
106 class(space_t), intent(in) :: space
107 type(grid_t), intent(in) :: gr
108 character(len=*), intent(in) :: dir
109 type(ions_t), intent(in) :: ions
110 type(hamiltonian_elec_t), intent(in) :: hm
111 type(output_t), intent(in) :: outp
113
114 call io_mkdir(trim(dir), namespace)
115 call output_states(outp, namespace, space, trim(dir), tg%st, gr, ions, hm, -1)
116
118 end subroutine target_output_gstransformation
119 ! ----------------------------------------------------------------------
120
121
122
123 ! ----------------------------------------------------------------------
125 real(real64) function target_j1_gstransformation(tg, gr, psi) result(j1)
126 type(target_t), intent(in) :: tg
127 type(grid_t), intent(in) :: gr
128 type(states_elec_t), intent(in) :: psi
129
130 integer :: ik, ist
131
132 complex(real64), allocatable :: zpsi(:, :), zst(:, :)
133
134 push_sub(target_j1_gstransformation)
135
136 safe_allocate(zpsi(1:gr%np, 1:tg%st%d%dim))
137 safe_allocate(zst(1:gr%np, 1:tg%st%d%dim))
138
139 j1 = m_zero
140 do ik = 1, psi%nik
141 do ist = psi%st_start, psi%st_end
142
143 call states_elec_get_state(psi, gr, ist, ik, zpsi)
144 call states_elec_get_state(tg%st, gr, ist, ik, zst)
145
146 j1 = j1 + abs(zmf_dotp(gr, psi%d%dim, zpsi, zst))**2
147
148 end do
149 end do
150
151 safe_deallocate_a(zpsi)
152 safe_deallocate_a(zst)
154 pop_sub(target_j1_gstransformation)
155 end function target_j1_gstransformation
156
157
158 ! ----------------------------------------------------------------------
160 subroutine target_chi_gstransformation(tg, gr, psi_in, chi_out)
161 type(target_t), intent(in) :: tg
162 type(grid_t), intent(in) :: gr
163 type(states_elec_t), intent(in) :: psi_in
164 type(states_elec_t), intent(inout) :: chi_out
165
166 integer :: ik, ist
167 complex(real64) :: olap
168 complex(real64), allocatable :: zpsi(:, :), zst(:, :), zchi(:, :)
169
171
172 safe_allocate(zpsi(1:gr%np, 1:tg%st%d%dim))
173 safe_allocate(zst(1:gr%np, 1:tg%st%d%dim))
174 safe_allocate(zchi(1:gr%np, 1:tg%st%d%dim))
175
176 do ik = 1, psi_in%nik
177 do ist = psi_in%st_start, psi_in%st_end
178
179 call states_elec_get_state(psi_in, gr, ist, ik, zpsi)
180 call states_elec_get_state(tg%st, gr, ist, ik, zst)
181
182 olap = zmf_dotp(gr, zst(:, 1), zpsi(:, 1))
183 zchi(1:gr%np, 1:tg%st%d%dim) = olap*zst(1:gr%np, 1:tg%st%d%dim)
184
185 call states_elec_set_state(chi_out, gr, ist, ik, zchi)
186
187 end do
188 end do
189
190 safe_deallocate_a(zpsi)
191 safe_deallocate_a(zst)
192 safe_deallocate_a(zchi)
193
195 end subroutine target_chi_gstransformation
198
199!! Local Variables:
200!! mode: f90
201!! coding: utf-8
202!! End:
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:609
real(real64), parameter, public m_zero
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
Definition: io.F90:114
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
logical pure function, public ion_dynamics_ions_move(this)
This module defines various routines, operating on mesh functions.
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:414
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
this module contains the low-level part of the output system
Definition: output_low.F90:115
this module contains the output system
Definition: output.F90:115
subroutine, public output_states(outp, namespace, space, dir, st, gr, ions, hm, iter)
Definition: output.F90:1888
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_transform(st, namespace, space, restart, mesh, kpoints, prefix)
real(real64) function, public target_j1_gstransformation(tg, gr, psi)
subroutine, public target_output_gstransformation(tg, namespace, space, gr, dir, ions, hm, outp)
subroutine, public target_chi_gstransformation(tg, gr, psi_in, chi_out)
subroutine, public target_init_gstransformation(gr, namespace, space, tg, td, restart, kpoints)
Definition: td.F90:114