51 integer :: interpolation_steps
52 real(real64),
allocatable,
public :: v_old(:, :, :)
53 real(real64),
allocatable,
public :: vtau_old(:, :, :)
54 logical :: mgga_with_exc
62 type(potential_interpolation_t),
intent(inout) :: vkso
63 type(potential_interpolation_t),
intent(in) :: vksi
67 safe_allocate_source_a(vkso%v_old, vksi%v_old)
68 safe_allocate_source_a(vkso%vtau_old, vksi%vtau_old)
69 vkso%mgga_with_exc = vksi%mgga_with_exc
70 vkso%interpolation_steps = vksi%interpolation_steps
78 type(potential_interpolation_t),
intent(inout) :: potential_interpolation
79 integer,
intent(in) :: np, nspin
80 logical,
intent(in) :: mgga_with_exc
81 integer,
optional,
intent(in) :: order
87 safe_allocate(potential_interpolation%v_old(1:np, 1:nspin, 0:potential_interpolation%interpolation_steps))
88 potential_interpolation%v_old(:, :, :) =
m_zero
90 potential_interpolation%mgga_with_exc = mgga_with_exc
92 if (potential_interpolation%mgga_with_exc)
then
94 safe_allocate(potential_interpolation%vtau_old(1:np, 1:nspin, 0:potential_interpolation%interpolation_steps))
95 potential_interpolation%vtau_old(:, :, :) =
m_zero
105 type(potential_interpolation_t),
intent(inout) :: potential_interpolation
109 assert(
allocated(potential_interpolation%v_old))
110 safe_deallocate_a(potential_interpolation%v_old)
111 safe_deallocate_a(potential_interpolation%vtau_old)
119 type(potential_interpolation_t),
intent(inout) :: potential_interpolation
120 integer,
intent(in) :: np, nspin
121 real(real64),
intent(in) :: vhxc(:, :)
122 real(real64),
optional,
intent(in) :: vtau(:, :)
124 integer :: i, ispin, ip
128 do i = 1, potential_interpolation%interpolation_steps
131 potential_interpolation%v_old(ip, ispin, i) = vhxc(ip, ispin)
136 if (potential_interpolation%mgga_with_exc)
then
137 assert(
present(vtau))
138 do i = 1, potential_interpolation%interpolation_steps
141 potential_interpolation%vtau_old(ip, ispin, i) = vtau(ip, ispin)
154 integer,
intent(in) :: np, nspin
155 real(real64),
intent(in) :: time, dt
156 real(real64),
contiguous,
intent(in) :: vhxc(:, :)
157 real(real64),
contiguous,
optional,
intent(in) :: vtau(:, :)
159 real(real64),
allocatable :: times(:)
164 safe_allocate(times(1:potential_interpolation%interpolation_steps))
165 do j = 1, potential_interpolation%interpolation_steps
166 times(j) = time - j*dt
169 do j = potential_interpolation%interpolation_steps, 2, -1
170 call lalg_copy(np, nspin, potential_interpolation%v_old(:, :, j-1), potential_interpolation%v_old(:, :, j))
172 call lalg_copy(np, nspin, vhxc(:, :), potential_interpolation%v_old(:, :, 1))
173 call interpolate( times, potential_interpolation%v_old(:, :, 1:potential_interpolation%interpolation_steps), &
174 time, potential_interpolation%v_old(:, :, 0))
176 if (potential_interpolation%mgga_with_exc)
then
177 assert(
present(vtau))
178 do j = potential_interpolation%interpolation_steps, 2, -1
179 call lalg_copy(np, nspin, potential_interpolation%vtau_old(:, :, j-1), potential_interpolation%vtau_old(:, :, j))
181 call lalg_copy(np, nspin, vtau(:, :), potential_interpolation%vtau_old(:, :, 1))
182 call interpolate( times, potential_interpolation%vtau_old(:, :, 1:potential_interpolation%interpolation_steps), &
183 time, potential_interpolation%vtau_old(:, :, 0))
186 safe_deallocate_a(times)
194 integer,
intent(in) :: np, nspin
195 integer,
intent(in) :: i
196 real(real64),
contiguous,
intent(inout) :: vhxc(:, :)
197 real(real64),
contiguous,
optional,
intent(in) :: vtau(:, :)
201 call lalg_copy(np, nspin, vhxc, potential_interpolation%v_old(:, :, i))
203 if (potential_interpolation%mgga_with_exc)
then
204 assert(
present(vtau))
205 call lalg_copy(np, nspin, vtau, potential_interpolation%vtau_old(:, :, i))
215 integer,
intent(in) :: np, nspin
216 integer,
intent(in) :: i
217 real(real64),
contiguous,
intent(inout) :: vhxc(:, :)
218 real(real64),
optional,
contiguous,
intent(inout) :: vtau(:, :)
222 call lalg_copy(np, nspin, potential_interpolation%v_old(:, :, i), vhxc)
224 if (potential_interpolation%mgga_with_exc)
then
225 assert(
present(vtau))
226 call lalg_copy(np, nspin, potential_interpolation%vtau_old(:, :, i), vtau)
237 integer,
intent(in) :: order
238 real(real64),
intent(in) :: time, dt, t
239 real(real64),
contiguous,
intent(inout) :: vhxc(:, :)
240 real(real64),
contiguous,
optional,
intent(inout) :: vtau(:, :)
243 real(real64),
allocatable :: times(:)
247 assert(order <= potential_interpolation%interpolation_steps)
249 safe_allocate(times(1:potential_interpolation%interpolation_steps))
250 do j = 1, potential_interpolation%interpolation_steps
251 times(j) = time - (j-1)*dt
254 call interpolate( times(1:order), potential_interpolation%v_old(:, :, 0:order-1), t, vhxc(:, :))
256 if (potential_interpolation%mgga_with_exc)
then
257 assert(
present(vtau))
258 call interpolate( times(1:order), potential_interpolation%vtau_old(:, :, 0:order-1), t, vtau(:, :))
261 safe_deallocate_a(times)
270 class(
space_t),
intent(in) :: space
271 class(
mesh_t),
intent(in) :: mesh
272 integer,
intent(in) :: nspin
273 integer,
intent(out) :: err2
275 integer :: ii, is, err
276 character(len=256) :: filename
282 do ii = 1, potential_interpolation%interpolation_steps - 1
284 write(filename,
'(a6,i2.2,i3.3)')
'vprev_', ii, is
286 potential_interpolation%v_old(1:mesh%np, is, ii), err)
289 if (err /= 0) err2 = err2 + 1
294 if (potential_interpolation%mgga_with_exc)
then
295 do ii = 1, potential_interpolation%interpolation_steps - 1
297 write(filename,
'(a9,i2.2,i3.3)')
'vtauprev_', ii, is
299 potential_interpolation%vtau_old(1:mesh%np, is, ii), err)
302 if (err /= 0) err2 = err2 + 1
315 class(
space_t),
intent(in) :: space
317 class(
mesh_t),
intent(in) :: mesh
318 integer,
intent(in) :: nspin
319 integer,
intent(out) :: err2
321 integer :: ii, is, err
322 character(len=256) :: filename
327 do ii = 1, potential_interpolation%interpolation_steps - 1
329 write(filename,
'(a,i2.2,i3.3)')
'vprev_', ii, is
331 potential_interpolation%v_old(1:mesh%np, is, ii), err)
334 message(1) =
"Unable to read VKS restart file '" // trim(filename) //
"'"
340 if (potential_interpolation%mgga_with_exc)
then
341 do ii = 1, potential_interpolation%interpolation_steps - 1
343 write(filename,
'(a,i2.2,i3.3)')
'vtauprev_', ii, is
345 potential_interpolation%vtau_old(1:mesh%np, is, ii), err)
348 message(1) =
"Unable to read VKS restart file '" // trim(filename) //
"'"
Copies a vector x, to a vector y.
This is the common interface to a simple-minded polynomical interpolation procedure (simple use of th...
real(real64), parameter, public m_zero
This module is intended to contain "only mathematical" functions and procedures.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_warning(no_lines, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public potential_interpolation_copy(vkso, vksi)
subroutine, public potential_interpolation_end(potential_interpolation)
subroutine, public potential_interpolation_interpolate(potential_interpolation, order, time, dt, t, vhxc, vtau)
subroutine, public potential_interpolation_set(potential_interpolation, np, nspin, i, vhxc, vtau)
subroutine, public potential_interpolation_new(potential_interpolation, np, nspin, time, dt, vhxc, vtau)
subroutine, public potential_interpolation_run_zero_iter(potential_interpolation, np, nspin, vhxc, vtau)
subroutine, public potential_interpolation_get(potential_interpolation, np, nspin, i, vhxc, vtau)
subroutine, public potential_interpolation_load(potential_interpolation, namespace, space, restart, mesh, nspin, err2)
subroutine, public potential_interpolation_dump(potential_interpolation, space, restart, mesh, nspin, err2)
subroutine, public potential_interpolation_init(potential_interpolation, np, nspin, mgga_with_exc, order)
Describes mesh distribution to nodes.