Octopus
pcm_potential.F90
Go to the documentation of this file.
1!! Copyright (C) 2014 Alain Delgado Gran, Carlo Andrea Rozzi, Stefano Corni, Gabriel Gil
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
23 use comm_oct_m
24 use debug_oct_m
26 use global_oct_m
27 use grid_oct_m
29 use index_oct_m
31 use io_oct_m
32 use ions_oct_m
33 use kick_oct_m
34 use lasers_oct_m
36 use mesh_oct_m
38 use mpi_oct_m
41 use parser_oct_m
42 use pcm_oct_m
47 use space_oct_m
49
50 ! to output debug info
51 use unit_oct_m
55
56 implicit none
57
58 private
59
60 public :: &
62
63contains
64
65 ! -----------------------------------------------------------------------------
67 subroutine pcm_hartree_potential(pcm, space, mesh, psolver, ext_partners, vhartree, density, pcm_corr, &
68 kick, time)
69 type(pcm_t), intent(inout) :: pcm
70 class(space_t), intent(in) :: space
71 class(mesh_t), intent(in) :: mesh
72 type(poisson_t), intent(inout) :: psolver
73 type(partner_list_t),intent(in) :: ext_partners
74 real(real64), intent(in) :: vhartree(:)
75 real(real64), intent(in) :: density(:)
76 real(real64), intent(out) :: pcm_corr
77 type(kick_t), optional, intent(in) :: kick
78 real(real64), optional, intent(in) :: time
79
80 real(real64), allocatable :: potx(:)
81 complex(real64), allocatable :: kick_eval(:)
82 real(real64), allocatable :: kick_real(:)
83 integer :: ii
84
85 logical :: kick_time
86 type(lasers_t), pointer :: lasers
87
88 push_sub(pcm_hartree_potential)
89
90 if (.not. pcm%run_pcm .or. .not. pcm_update(pcm)) then
91 pcm_corr = m_zero
93 return
94 end if
95
97 if (pcm%solute) then
98 call pcm_calc_pot_rs(pcm, mesh, psolver, v_h = vhartree, time_present = present(time))
99 end if
100
105 if (pcm%localf .and. present(time)) then
106 lasers => list_get_lasers(ext_partners)
107 if (associated(lasers) .and. present(kick)) then
108 safe_allocate(potx(1:mesh%np_part))
109 safe_allocate(kick_eval(1:mesh%np_part))
110 safe_allocate(kick_real(1:mesh%np_part))
111 potx = m_zero
112 kick_eval = m_zero
113 do ii = 1, lasers%no_lasers
114 call laser_potential(lasers%lasers(ii), mesh, potx, time)
115 end do
116 kick_real = m_zero
117 kick_time = ((pcm%iter-1)*pcm%dt <= kick%time) .and. (pcm%iter*pcm%dt > kick%time)
118 if (kick_time) then
119 call kick_function_get(space, mesh, kick, kick_eval, 1, to_interpolate = .true.)
120 kick_eval = kick%delta_strength * kick_eval
121 kick_real = real(kick_eval, real64)
122 end if
123 call pcm_calc_pot_rs(pcm, mesh, psolver, v_ext = potx, kick = -kick_real, &
124 time_present = present(time), kick_time = kick_time)
125 safe_deallocate_a(potx)
126 safe_deallocate_a(kick_eval)
127 safe_deallocate_a(kick_real)
128 else if (associated(lasers) .and. .not.present(kick)) then
129 safe_allocate(potx(1:mesh%np_part))
130 potx = m_zero
131 do ii = 1, lasers%no_lasers
132 call laser_potential(lasers%lasers(ii), mesh, potx, time)
133 end do
134 call pcm_calc_pot_rs(pcm, mesh, psolver, v_ext = potx, time_present = present(time))
135 safe_deallocate_a(potx)
136 else if (.not.associated(lasers) .and. present(kick)) then
137 safe_allocate(kick_eval(1:mesh%np_part))
138 safe_allocate(kick_real(1:mesh%np_part))
139 kick_eval = m_zero
140 kick_real = m_zero
141 kick_time =((pcm%iter-1)*pcm%dt <= kick%time) .and. (pcm%iter*pcm%dt > kick%time)
142 if (kick_time) then
143 call kick_function_get(space, mesh, kick, kick_eval, 1, to_interpolate = .true.)
144 kick_eval = kick%delta_strength * kick_eval
145 kick_real = real(kick_eval, real64)
146 end if
147 call pcm_calc_pot_rs(pcm, mesh, psolver, kick = -kick_real, &
148 time_present = present(time), kick_time = kick_time)
149 safe_deallocate_a(kick_eval)
150 safe_deallocate_a(kick_real)
151 end if
152
153 ! Calculating the PCM term renormalizing the sum of the single-particle energies
154 ! to keep the idea of pcm_corr... but it will be added later on
155 pcm_corr = dmf_dotp( mesh, density, pcm%v_e_rs + pcm%v_n_rs + pcm%v_ext_rs)
156 else
157 ! Calculating the PCM term renormalizing the sum of the single-particle energies
158 pcm_corr = dmf_dotp( mesh, density, pcm%v_e_rs + pcm%v_n_rs)
159 end if
161 if (debug%info) then
162 call messages_write(' PCM potential updated')
163 call messages_new_line()
164 call messages_write(' PCM update iteration counter: ')
165 call messages_write(pcm%iter)
166 call messages_info()
167 end if
168
169 pop_sub(pcm_hartree_potential)
170
171 end subroutine pcm_hartree_potential
172
173end module pcm_potential_oct_m
174
175!! Local Variables:
176!! mode: f90
177!! coding: utf-8
178!! End:
type(debug_t), save, public debug
Definition: debug.F90:154
type(lasers_t) function, pointer, public list_get_lasers(partners)
real(real64), parameter, public m_zero
Definition: global.F90:187
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module implements the index, used for the mesh points.
Definition: index.F90:122
This module defines classes and functions for interaction partners.
Definition: io.F90:114
subroutine, public kick_function_get(space, mesh, kick, kick_function, iq, to_interpolate)
Definition: kick.F90:993
subroutine, public laser_potential(laser, mesh, pot, time)
Definition: lasers.F90:1039
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
subroutine, public messages_new_line()
Definition: messages.F90:1146
Some general things and nomenclature:
Definition: par_vec.F90:171
subroutine, public pcm_calc_pot_rs(pcm, mesh, psolver, ions, v_h, v_ext, kick, time_present, kick_time)
Definition: pcm.F90:1215
logical function, public pcm_update(this)
Update pcm potential.
Definition: pcm.F90:3169
subroutine, public pcm_hartree_potential(pcm, space, mesh, psolver, ext_partners, vhartree, density, pcm_corr, kick, time)
PCM reaction field due to the electronic density.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
int true(void)