Octopus
test.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
21#include "global.h"
22
23module test_oct_m
24 use batch_oct_m
27 use box_oct_m
32 use clock_oct_m
33 use debug_oct_m
38 use global_oct_m
39 use grid_oct_m
44 use iihash_oct_m
46 use iso_c_binding
47 use io_oct_m
49 use, intrinsic :: iso_fortran_env
53 use mesh_oct_m
61 use mpi_oct_m
68 use parser_oct_m
74 use sihash_oct_m
76 use space_oct_m
77 use sphash_oct_m
86 use types_oct_m
87 use v_ks_oct_m
91
92 implicit none
93
95 private
96 integer :: type
97 integer :: repetitions
98 integer :: min_blocksize
99 integer :: max_blocksize
100 end type test_parameters_t
101
102 !Auxiliary quantities needed by the linear solver
103 real(real64) :: shift_aux
104 type(derivatives_t), pointer :: der_aux => null()
105 type(preconditioner_t) :: prec_aux
106
107 public :: test_run
108
109contains
110
111 ! ---------------------------------------------------------
112 subroutine test_run(namespace)
113 type(namespace_t), intent(in) :: namespace
114
115 type(test_parameters_t) :: param
116 integer :: test_mode
117
118 push_sub(test_run)
119
120 call messages_obsolete_variable(namespace, 'WhichTest', 'TestMode')
121
122 !%Variable TestMode
123 !%Type integer
124 !%Default hartree
125 !%Section Calculation Modes::Test
126 !%Description
127 !% Decides what kind of test should be performed.
128 !%Option hartree 1
129 !% Tests the Poisson solvers used to calculate the Hartree potential.
130 !%Option derivatives 2
131 !% Tests and benchmarks the implementation of the finite-difference operators, used to calculate derivatives.
132 !%Option orthogonalization 3
133 !% Tests the implementation of the orthogonalization routines.
134 !%Option interpolation 4
135 !% Test the interpolation routines.
136 !%Option ion_interaction 5
137 !% Tests the ion-ion interaction routines.
138 !%Option projector 6
139 !% Tests the code that applies the nonlocal part of the pseudopotentials
140 !% in case of spin-orbit coupling
141 !%Option dft_u 7
142 !% Tests the DFT+U part of the code for projections on the basis.
143 !%Option hamiltonian_apply 8
144 !% Tests the application of the Hamiltonian, or a part of it
145 !%Option density_calc 9
146 !% Calculation of the density.
147 !%Option exp_apply 10
148 !% Tests the exponential of the Hamiltonian
149 !%Option boundaries 11
150 !% Tests the boundaries conditions
151 !%Option subspace_diag 12
152 !% Tests the subspace diagonalization
153 !%Option batch_ops 13
154 !% Tests the batch operations
155 !%Option clock 18
156 !% Tests for clock
157 !%Option linear_solver 19
158 !% Tests the linear solvers
159 !%Option cgal 20
160 !% Tests for cgal interface
161 !%Option dense_eigensolver 21
162 !% Tests for dense eigensolvers (especially parallel ones)
163 !%Option grid_interpolation 22
164 !% Tests for grid interpolation and multigrid methods.
165 !%Option iihash 23
166 !% Tests for the integer-integer hash table.
167 !%Option sihash 24
168 !% Tests for the string-integer hash table.
169 !%Option sphash 25
170 !% Tests for the string-polymorphic hash table.
171 !%Option mpiwrappers 26
172 !% Tests for the MPI wrappers with large integer displacements.
173 !%Option regridding 27
174 !% Tests the regridding between two different grids.
175 !%Option helmholtz_decomposition 28
176 !% Test for the Helmholtz decomposition subroutines
177 !%Option vecpot_analytical 29
178 !% Tests analytically the vector potential from B field.
179 !%Option current_density 30
180 !% Tests the different contributions to the total electronic current density
181 !%Option mixing_tests 31
182 !% Unit test for the mixing
183 !%Option optimizers 32
184 !% Tests for the optimizers, using standard test functions
185 !%End
186 call parse_variable(namespace, 'TestMode', option__testmode__hartree, test_mode)
188 call messages_obsolete_variable(namespace, 'TestDerivatives', 'TestType')
189 call messages_obsolete_variable(namespace, 'TestOrthogonalization', 'TestType')
191 !%Variable TestType
192 !%Type integer
193 !%Default all
194 !%Section Calculation Modes::Test
195 !%Description
196 !% Decides on what type of values the test should be performed.
197 !%Option real 1
198 !% Test for double-precision real functions.
199 !%Option complex 2
200 !%Option all 3
201 !% Tests for double-precision real and complex functions.
202 !%End
203 call parse_variable(namespace, 'TestType', option__testtype__all, param%type)
204 if (param%type < 1 .or. param%type > 5) then
205 message(1) = "Invalid option for TestType."
206 call messages_fatal(1, only_root_writes = .true., namespace=namespace)
207 end if
208
209 !%Variable TestRepetitions
210 !%Type integer
211 !%Default 1
212 !%Section Calculation Modes::Test
213 !%Description
214 !% This variable controls the behavior of oct-test for performance
215 !% benchmarking purposes. It sets the number of times the
216 !% computational kernel of a test will be executed, in order to
217 !% provide more accurate timings.
218 !%
219 !% Currently this variable is used by the <tt>hartree_test</tt>,
220 !% <tt>derivatives</tt>, and <tt>projector</tt> tests.
221 !%End
222 call parse_variable(namespace, 'TestRepetitions', 1, param%repetitions)
223
224 !%Variable TestMinBlockSize
225 !%Type integer
226 !%Default 1
227 !%Section Calculation Modes::Test
228 !%Description
229 !% Some tests can work with multiple blocksizes, in this case of
230 !% range of blocksizes will be tested. This variable sets the lower
231 !% bound of that range.
232 !%
233 !% Currently this variable is only used by the derivatives test.
234 !%End
235 call parse_variable(namespace, 'TestMinBlockSize', 1, param%min_blocksize)
236
237 !%Variable TestMaxBlockSize
238 !%Type integer
239 !%Default 128
240 !%Section Calculation Modes::Test
241 !%Description
242 !% Some tests can work with multiple blocksizes, in this case of
243 !% range of blocksizes will be tested. This variable sets the lower
244 !% bound of that range.
245 !%
246 !% Currently this variable is only used by the derivatives test.
247 !%End
248 call parse_variable(namespace, 'TestMaxBlockSize', 128, param%max_blocksize)
249
250 call messages_print_with_emphasis(msg="Test mode", namespace=namespace)
251 call messages_print_var_option("TestMode", test_mode, namespace=namespace)
252 call messages_print_var_option("TestType", param%type, namespace=namespace)
253 call messages_print_var_value("TestRepetitions", param%repetitions, namespace=namespace)
254 call messages_print_var_value("TestMinBlockSize", param%min_blocksize, namespace=namespace)
255 call messages_print_var_value("TestMaxBlockSize", param%max_blocksize, namespace=namespace)
256 call messages_print_with_emphasis(namespace=namespace)
257
258 select case (test_mode)
259 case (option__testmode__hartree)
260 call test_hartree(param, namespace)
261 case (option__testmode__derivatives)
262 call test_derivatives(param, namespace)
263 case (option__testmode__orthogonalization)
264 call test_orthogonalization(param, namespace)
265 case (option__testmode__interpolation)
266 call test_interpolation(param, namespace)
267 case (option__testmode__ion_interaction)
268 call test_ion_interaction(namespace)
269 case (option__testmode__projector)
270 call test_projector(param, namespace)
271 case (option__testmode__dft_u)
272 call test_dft_u(param, namespace)
273 case (option__testmode__hamiltonian_apply)
274 call test_hamiltonian(param, namespace)
275 case (option__testmode__density_calc)
276 call test_density_calc(param, namespace)
277 case (option__testmode__exp_apply)
278 call test_exponential(param, namespace)
279 case (option__testmode__boundaries)
280 call test_boundaries(param, namespace)
281 case (option__testmode__subspace_diag)
282 call test_subspace_diagonalization(param, namespace)
283 case (option__testmode__batch_ops)
284 call test_batch_ops(param, namespace)
285 case (option__testmode__clock)
286 call test_clock()
287 case (option__testmode__linear_solver)
288 call test_linear_solver(namespace)
289 case (option__testmode__cgal)
290 call test_cgal()
291 case (option__testmode__dense_eigensolver)
293 case (option__testmode__grid_interpolation)
295 case (option__testmode__iihash)
296 call test_iihash()
297 case (option__testmode__sihash)
298 call test_sihash()
299 case (option__testmode__sphash)
300 call test_sphash(namespace)
301 case (option__testmode__mpiwrappers)
302 call test_mpiwrappers()
303 case (option__testmode__regridding)
304 call test_regridding(namespace)
305 case (option__testmode__helmholtz_decomposition)
306 call test_helmholtz_decomposition(namespace)
307 case (option__testmode__vecpot_analytical)
308 call test_vecpot_analytical(namespace)
309 case (option__testmode__current_density)
310 call test_current_density(namespace)
311 case (option__testmode__mixing_tests)
312 call mix_tests_run()
313 case (option__testmode__optimizers)
314 call test_optimizers(namespace)
315 end select
316
317 pop_sub(test_run)
318 end subroutine test_run
319
320 ! ---------------------------------------------------------
321 subroutine test_hartree(param, namespace)
322 type(test_parameters_t), intent(in) :: param
323 type(namespace_t), intent(in) :: namespace
324
325 type(electrons_t), pointer :: sys
326
327 push_sub(test_hartree)
328
329 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
330
331 sys => electrons_t(namespace, generate_epot=.false.)
332 call sys%init_parallelization(mpi_world)
333 call poisson_test(sys%hm%psolver, sys%space, sys%gr, sys%ions%latt, namespace, param%repetitions)
334 safe_deallocate_p(sys)
335
336 pop_sub(test_hartree)
337 end subroutine test_hartree
338
339 ! ---------------------------------------------------------
340 subroutine test_helmholtz_decomposition(namespace)
341 type(namespace_t), intent(in) :: namespace
342
343 type(electrons_t), pointer :: sys
344 type(helmholtz_decomposition_t) :: helmholtz
345
347
348 ! First of all we have to initialize the grid and the poisson solver
349 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
350
351 sys => electrons_t(namespace, generate_epot=.false.)
352 call sys%init_parallelization(mpi_world)
353
354 call helmholtz%init(namespace, sys%gr, sys%mc, sys%space)
355
356 ! Then we have to initialize the exact fields
357 write(message(1),'(a)') "Helmholtz decomposition: beginning Hertzian dipole test"
358 call messages_info(1, namespace=namespace)
359 call hertzian_dipole_test(helmholtz, sys%gr, sys%namespace, sys%space)
360
361 write(message(1),'(a)') "Helmholtz decomposition: beginning Gaussian test"
362 call messages_info(1, namespace=namespace)
363 call gaussian_test(helmholtz, sys%gr, sys%namespace, sys%space)
364
365 safe_deallocate_p(sys)
366
368 end subroutine test_helmholtz_decomposition
369
370 ! ---------------------------------------------------------
371 subroutine test_linear_solver(namespace)
372 type(namespace_t), intent(in) :: namespace
373
374 type(electrons_t), pointer :: sys
375 real(real64), allocatable :: rho(:), x(:), center(:)
376 real(real64) :: rr, alpha, beta, res
377 integer :: ip, iter
378
379 real(real64), parameter :: threshold = 1.e-7_real64
380
381 push_sub(test_linear_solver)
382
383 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
384
385 sys => electrons_t(namespace, generate_epot=.false.)
386 call sys%init_parallelization(mpi_world)
387
388 ! We need to set up some auxiliary quantities called by the linear solver
389 call mesh_init_mesh_aux(sys%gr)
390 ! Shift of the linear equation
391 shift_aux = 0.25_real64
392 ! Preconditioner used for the QMR algorithm
393 call preconditioner_init(prec_aux, namespace, sys%gr, sys%mc, sys%space)
394 ! Derivative object needed
395 call set_der_aux(sys%gr%der)
396
397 ! Here we put a Gaussian as the right-hand side of the linear solver
398 ! Values are taken from the poisson_test routine
399 alpha = m_four * sys%gr%spacing(1)
400 beta = m_one / (alpha**sys%space%dim * sqrt(m_pi)**sys%space%dim)
401 ! The Gaussian is centered around the origin
402 safe_allocate(center(1:sys%space%dim))
403 center = m_zero
404
405 safe_allocate(rho(1:sys%gr%np))
406 rho = m_zero
407 do ip = 1, sys%gr%np
408 call mesh_r(sys%gr, ip, rr, origin = center(:))
409 rho(ip) = beta*exp(-(rr/alpha)**2)
410 end do
411
412 safe_allocate(x(1:sys%gr%np))
413
414 !Test the CG linear solver
415 x = m_zero
416 iter = 1000
417 call dconjugate_gradients(sys%gr%np, x, rho, laplacian_op, dmf_dotp_aux, iter, res, threshold)
418 write(message(1),'(a,i6,a)') "Info: CG converged with ", iter, " iterations."
419 write(message(2),'(a,e14.6)') "Info: The residue is ", res
420 write(message(3),'(a,e14.6)') "Info: Norm solution CG ", dmf_nrm2(sys%gr, x)
421 call messages_info(3, namespace=namespace)
422
423 call preconditioner_end(prec_aux)
424 safe_deallocate_a(x)
425 safe_deallocate_a(rho)
426 safe_deallocate_p(sys)
427
428 pop_sub(test_linear_solver)
429 contains
430
431 subroutine set_der_aux(der)
432 type(derivatives_t), target, intent(in) :: der
434 der_aux => der
436 end subroutine set_der_aux
437
438 ! ---------------------------------------------------------
440 subroutine laplacian_op(x, hx)
441 real(real64), contiguous,intent(in) :: x(:)
442 real(real64), contiguous,intent(out) :: Hx(:)
443
444 real(real64), allocatable :: tmpx(:)
445
446 assert(associated(mesh_aux))
447
448 safe_allocate(tmpx(1:mesh_aux%np_part))
449 call lalg_copy(mesh_aux%np, x, tmpx)
450 call dderivatives_lapl(der_aux, tmpx, hx)
451 call lalg_scal(mesh_aux%np, -m_one, hx)
452 call lalg_axpy(mesh_aux%np, shift_aux, x, hx)
453 safe_deallocate_a(tmpx)
454
455 end subroutine laplacian_op
456
457 end subroutine test_linear_solver
458
459
460 ! ---------------------------------------------------------
461 subroutine test_projector(param, namespace)
462 type(test_parameters_t), intent(in) :: param
463 type(namespace_t), intent(in) :: namespace
465 type(electrons_t), pointer :: sys
466 type(wfs_elec_t) :: epsib
467 integer :: itime
468 complex(real64), allocatable :: psi(:, :)
469
470 push_sub(test_projector)
471
472 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
473
474 call messages_write('Info: Testing the nonlocal part of the pseudopotential with SOC')
475 call messages_new_line()
476 call messages_new_line()
477 call messages_info(namespace=namespace)
478
479 sys => electrons_t(namespace, generate_epot=.false.)
480 call sys%init_parallelization(mpi_world)
481
482 call states_elec_allocate_wfns(sys%st, sys%gr, wfs_type = type_cmplx)
483 call test_batch_set_gaussian(sys%st%group%psib(1, 1), sys%gr)
484
485 ! Initialize external potential
486 call hamiltonian_elec_epot_generate(sys%hm, sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st)
487
488 call sys%st%group%psib(1, 1)%copy_to(epsib)
489
490 call batch_set_zero(epsib)
491
492 do itime = 1, param%repetitions
493 call zproject_psi_batch(sys%gr, sys%gr%der%boundaries, sys%hm%ep%proj, &
494 sys%hm%ep%natoms, 2, sys%st%group%psib(1, 1), epsib)
495 end do
496
497 safe_allocate(psi(1:sys%gr%np, 1:sys%st%d%dim))
498 do itime = 1, epsib%nst
499 call batch_get_state(epsib, itime, sys%gr%np, psi)
500 write(message(1),'(a,i1,3x, f12.6)') "Norm state ", itime, zmf_nrm2(sys%gr, 2, psi)
501 call messages_info(1, namespace=sys%namespace)
502 end do
503 safe_deallocate_a(psi)
504
505 call epsib%end()
506 call states_elec_deallocate_wfns(sys%st)
507 safe_deallocate_p(sys)
508
509 pop_sub(test_projector)
510 end subroutine test_projector
511
512 ! ---------------------------------------------------------
513 subroutine test_dft_u(param, namespace)
514 type(test_parameters_t), intent(in) :: param
515 type(namespace_t), intent(in) :: namespace
516
517 type(electrons_t), pointer :: sys
518 type(wfs_elec_t) :: epsib, epsib2
519 integer :: itime, ist
520 type(orbitalbasis_t) :: basis
521 real(real64), allocatable :: ddot(:,:,:), dweight(:,:)
522 complex(real64), allocatable :: zdot(:,:,:), zweight(:,:)
523
524 push_sub(test_dft_u)
525
526 call calc_mode_par%unset_parallelization(p_strategy_states)
527 call calc_mode_par%unset_parallelization(p_strategy_kpoints)
528 call calc_mode_par%set_parallelization(p_strategy_domains, default = .true.)
529
530 call messages_write('Info: Testing some DFT+U routines')
531 call messages_new_line()
532 call messages_new_line()
533 call messages_info(namespace=namespace)
534
535 sys => electrons_t(namespace, generate_epot=.false.)
536 call sys%init_parallelization(mpi_world)
537
538 call states_elec_allocate_wfns(sys%st, sys%gr)
539 call test_batch_set_gaussian(sys%st%group%psib(1, 1), sys%gr)
540 if (sys%st%pack_states) call sys%st%pack()
541
542 call sys%st%group%psib(1, 1)%copy_to(epsib2, copy_data = .true.)
543
544 ! We set the phase of the batch if needed
545 if (.not. sys%hm%phase%is_allocated()) then
546 call sys%st%group%psib(1, 1)%copy_to(epsib, copy_data = .true.)
547 else
548 call sys%st%group%psib(1, 1)%copy_to(epsib)
549 call sys%hm%phase%apply_to(sys%gr, sys%gr%np, &
550 .false., epsib, src=sys%st%group%psib(1, 1))
551 epsib2%has_phase = .true.
552 end if
553
554 ! Initialize the orbital basis
555 call orbitalbasis_init(basis, sys%namespace, sys%space%periodic_dim)
556 if (states_are_real(sys%st)) then
557 call dorbitalbasis_build(basis, sys%namespace, sys%ions, sys%gr, sys%st%d%kpt, sys%st%d%dim, &
558 sys%hm%phase%is_allocated(), .false., .false.)
559 safe_allocate(dweight(1:basis%orbsets(1)%norbs, 1:epsib%nst_linear))
560 safe_allocate(ddot(1:sys%st%d%dim, 1:basis%orbsets(1)%norbs, 1:epsib%nst))
561 else
562 call zorbitalbasis_build(basis, sys%namespace, sys%ions, sys%gr, sys%st%d%kpt, sys%st%d%dim, &
563 sys%hm%phase%is_allocated(), .false., .false.)
564 call orbitalset_update_phase(basis%orbsets(1), sys%space%dim, sys%st%d%kpt, sys%kpoints, (sys%st%d%ispin==spin_polarized))
565 safe_allocate(zweight(1:basis%orbsets(1)%norbs, 1:epsib%nst_linear))
566 safe_allocate(zdot(1:sys%st%d%dim, 1:basis%orbsets(1)%norbs, 1:epsib%nst))
567
568 !We set the phase of the orbitals if needed
569 if (sys%hm%phase%is_allocated()) then
570 call orbitalset_update_phase(basis%orbsets(1), sys%space%dim, sys%st%d%kpt, sys%kpoints, &
571 (sys%st%d%ispin==spin_polarized))
572 end if
573 end if
574
575 do itime = 1, param%repetitions
576 call batch_set_zero(epsib2)
577 if (states_are_real(sys%st)) then
578 dweight = m_one
579 ddot = m_zero
580 call dorbitalset_get_coeff_batch(basis%orbsets(1), sys%st%d%dim, sys%st%group%psib(1, 1), ddot)
581 call dorbitalset_add_to_batch(basis%orbsets(1), sys%st%d%dim, epsib2, dweight)
582 else
583 zweight = m_one
584 zdot = m_zero
585 call zorbitalset_get_coeff_batch(basis%orbsets(1), sys%st%d%dim, epsib, zdot)
586 call zorbitalset_add_to_batch(basis%orbsets(1), sys%st%d%dim, epsib2, zweight)
587 end if
588 end do
589
590 if (epsib%is_packed()) then
591 call epsib%do_unpack(force = .true.)
592 end if
593
594 do ist = 1, epsib%nst
595 if (states_are_real(sys%st)) then
596 write(message(1),'(a,i2,3x,e13.6)') "Dotp state ", ist, ddot(1,1,ist)
597 else
598 write(message(1),'(a,i2,2(3x,e13.6))') "Dotp state ", ist, zdot(1,1,ist)
599 end if
600 call messages_info(1)
601 end do
602
603
604 call test_prints_info_batch(sys%st, sys%gr, epsib2)
605
606 safe_deallocate_a(dweight)
607 safe_deallocate_a(zweight)
608 safe_deallocate_a(ddot)
609 safe_deallocate_a(zdot)
610
611 call epsib%end()
612 call epsib2%end()
613 call orbitalbasis_end(basis)
614 call states_elec_deallocate_wfns(sys%st)
615 safe_deallocate_p(sys)
616
617 pop_sub(test_dft_u)
618 end subroutine test_dft_u
619
620 ! ---------------------------------------------------------
621 subroutine test_hamiltonian(param, namespace)
622 type(test_parameters_t), intent(in) :: param
623 type(namespace_t), intent(in) :: namespace
624
625 type(electrons_t), pointer :: sys
626 type(wfs_elec_t) :: hpsib
627 integer :: itime, terms
628
629 push_sub(test_hamiltonian)
630
631 !%Variable TestHamiltonianApply
632 !%Type integer
633 !%Default term_all
634 !%Section Calculation Modes::Test
635 !%Description
636 !% Decides which part of the Hamiltonian is applied.
637 !%Option term_all 0
638 !% Apply the full Hamiltonian.
639 !%Option term_kinetic 1
640 !% Apply only the kinetic operator
641 !%Option term_local_potential 2
642 !% Apply only the local potential.
643 !%Option term_non_local_potential 4
644 !% Apply only the non_local potential.
645 !%End
646 call parse_variable(namespace, 'TestHamiltonianApply', option__testhamiltonianapply__term_all, terms)
647 if (terms == 0) terms = huge(1)
648
649
650 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
651
652 call messages_write('Info: Testing the application of the Hamiltonian')
653 call messages_new_line()
654 call messages_new_line()
655 call messages_info(namespace=namespace)
656
657 sys => electrons_t(namespace, generate_epot=.false.)
658 call sys%init_parallelization(mpi_world)
659
660 call states_elec_allocate_wfns(sys%st, sys%gr)
661 call test_batch_set_gaussian(sys%st%group%psib(1, 1), sys%gr)
662
663 ! Initialize external potential
664 if (sys%st%pack_states .and. sys%hm%apply_packed()) call sys%st%pack()
665 call hamiltonian_elec_epot_generate(sys%hm, sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st)
666 call density_calc(sys%st, sys%gr, sys%st%rho)
667 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, sys%st, sys%ions, sys%ext_partners)
668
669 call boundaries_set(sys%gr%der%boundaries, sys%gr, sys%st%group%psib(1, 1))
670
671 call sys%st%group%psib(1, 1)%copy_to(hpsib)
672
673 if (sys%hm%apply_packed()) then
674 call sys%st%group%psib(1, 1)%do_pack()
675 call hpsib%do_pack(copy = .false.)
676 end if
677
678 do itime = 1, param%repetitions
679 if (states_are_real(sys%st)) then
680 call dhamiltonian_elec_apply_batch(sys%hm, sys%namespace, sys%gr, sys%st%group%psib(1, 1), hpsib, terms = terms, &
681 set_bc = .false.)
682 else
683 call zhamiltonian_elec_apply_batch(sys%hm, sys%namespace, sys%gr, sys%st%group%psib(1, 1), hpsib, terms = terms, &
684 set_bc = .false.)
685 end if
686 end do
687
688 if (hpsib%is_packed()) then
689 call hpsib%do_unpack(force = .true.)
690 end if
691
692 call test_prints_info_batch(sys%st, sys%gr, hpsib)
693
694 call hpsib%end(copy = .false.)
695 call states_elec_deallocate_wfns(sys%st)
696 safe_deallocate_p(sys)
697
698 pop_sub(test_hamiltonian)
699 end subroutine test_hamiltonian
700
701
702 ! ---------------------------------------------------------
703 subroutine test_density_calc(param, namespace)
704 type(test_parameters_t), intent(in) :: param
705 type(namespace_t), intent(in) :: namespace
706
707 type(electrons_t), pointer :: sys
708 integer :: itime
709
710 push_sub(test_density_calc)
711
712 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
713
714 call messages_write('Info: Testing density calculation')
715 call messages_new_line()
716 call messages_new_line()
717 call messages_info(namespace=namespace)
718
719 sys => electrons_t(namespace, generate_epot=.false.)
720 call sys%init_parallelization(mpi_world)
721
722 call states_elec_allocate_wfns(sys%st, sys%gr)
723 call test_batch_set_gaussian(sys%st%group%psib(1, 1), sys%gr)
724 if (sys%st%pack_states) call sys%st%pack()
725
726 do itime = 1, param%repetitions
727 call density_calc(sys%st, sys%gr, sys%st%rho)
728 end do
729
730 write(message(1),'(a,3x, f12.6)') "Norm density ", dmf_nrm2(sys%gr, sys%st%rho(:,1))
731 call messages_info(1, namespace=namespace)
732
733 call states_elec_deallocate_wfns(sys%st)
734 safe_deallocate_p(sys)
735
736 pop_sub(test_density_calc)
737 end subroutine test_density_calc
738
739
740 ! ---------------------------------------------------------
741 subroutine test_boundaries(param, namespace)
742 type(test_parameters_t), intent(in) :: param
743 type(namespace_t), intent(in) :: namespace
744
745 type(electrons_t), pointer :: sys
746 integer :: itime
747
748 push_sub(test_density_calc)
749
750 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
751
752 call messages_write('Info: Testing boundary conditions')
753 call messages_new_line()
754 call messages_new_line()
755 call messages_info(namespace=namespace)
756
757 sys => electrons_t(namespace, generate_epot=.false.)
758 call sys%init_parallelization(mpi_world)
759
760 call states_elec_allocate_wfns(sys%st, sys%gr)
761 call test_batch_set_gaussian(sys%st%group%psib(1, 1), sys%gr)
762 if (sys%st%pack_states) call sys%st%pack()
763
764 do itime = 1, param%repetitions
765 call boundaries_set(sys%gr%der%boundaries, sys%gr, sys%st%group%psib(1, 1))
766 end do
767
768 call test_prints_info_batch(sys%st, sys%gr, sys%st%group%psib(1, 1))
769
770 call states_elec_deallocate_wfns(sys%st)
771 safe_deallocate_p(sys)
772
773 pop_sub(test_density_calc)
774 end subroutine test_boundaries
775
776
777 ! ---------------------------------------------------------
778 subroutine test_exponential(param, namespace)
779 type(test_parameters_t), intent(in) :: param
780 type(namespace_t), intent(in) :: namespace
781
782 type(electrons_t), pointer :: sys
783 type(exponential_t) :: te
784 integer :: itime
785
786 push_sub(test_exponential)
787
788 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
789
790 call messages_write('Info: Testing exponential')
791 call messages_new_line()
792 call messages_new_line()
793 call messages_info(namespace=namespace)
794
795 sys => electrons_t(namespace, generate_epot=.false.)
796 call sys%init_parallelization(mpi_world)
797
798 call states_elec_allocate_wfns(sys%st, sys%gr, wfs_type=type_cmplx)
799 call test_batch_set_gaussian(sys%st%group%psib(1, 1), sys%gr)
800
801 ! Initialize external potential
802 if (sys%st%pack_states .and. sys%hm%apply_packed()) call sys%st%pack()
803 call hamiltonian_elec_epot_generate(sys%hm, sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st)
804 call density_calc(sys%st, sys%gr, sys%st%rho)
805 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, sys%st, sys%ions, sys%ext_partners)
806
807 call exponential_init(te, namespace)
808
809 if (sys%hm%apply_packed()) then
810 call sys%st%group%psib(1, 1)%do_pack()
811 end if
812
813 do itime = 1, param%repetitions
814 call te%apply_batch(sys%namespace, sys%gr, sys%hm, sys%st%group%psib(1, 1), 1e-3_real64)
815 end do
816
817 call test_prints_info_batch(sys%st, sys%gr, sys%st%group%psib(1, 1))
818
819 call states_elec_deallocate_wfns(sys%st)
820 safe_deallocate_p(sys)
821
822 pop_sub(test_exponential)
823 end subroutine test_exponential
824
825
826 ! ---------------------------------------------------------
827 subroutine test_subspace_diagonalization(param, namespace)
828 type(test_parameters_t), intent(in) :: param
829 type(namespace_t), intent(in) :: namespace
830
831 type(electrons_t), pointer :: sys
832 integer :: itime
833 type(subspace_t) :: sdiag
834 real(real64), allocatable :: diff(:)
835
837
838 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
839
840 call messages_write('Info: Testing boundary conditions')
841 call messages_new_line()
842 call messages_new_line()
843 call messages_info(namespace=namespace)
844
845 sys => electrons_t(namespace, generate_epot=.false.)
846 call sys%init_parallelization(mpi_world)
847
848 call states_elec_allocate_wfns(sys%st, sys%gr)
849 call test_batch_set_gaussian(sys%st%group%psib(1, 1), sys%gr)
850
851 if (sys%st%pack_states .and. sys%hm%apply_packed()) call sys%st%pack()
852 call hamiltonian_elec_epot_generate(sys%hm, sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st)
853 call density_calc(sys%st, sys%gr, sys%st%rho)
854 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, sys%st, sys%ions, sys%ext_partners)
855
856 call sdiag%init(sys%namespace, sys%st)
857
858 safe_allocate(diff(1:sys%st%nst))
859
860 do itime = 1, param%repetitions
861 if (states_are_real(sys%st)) then
862 call dsubspace_diag(sdiag, sys%namespace, sys%gr, sys%st, sys%hm, 1, sys%st%eigenval(:, 1), diff)
863 else
864 call zsubspace_diag(sdiag, sys%namespace, sys%gr, sys%st, sys%hm, 1, sys%st%eigenval(:, 1), diff)
865 end if
866 end do
867
868 safe_deallocate_a(diff)
869
870 call test_prints_info_batch(sys%st, sys%gr, sys%st%group%psib(1, 1))
872 call states_elec_deallocate_wfns(sys%st)
873 safe_deallocate_p(sys)
874
876 end subroutine test_subspace_diagonalization
877
878
879 ! ---------------------------------------------------------
880 subroutine test_batch_ops(param, namespace)
881 type(test_parameters_t), intent(in) :: param
882 type(namespace_t), intent(in) :: namespace
883
884 type(electrons_t), pointer :: sys
885 integer :: itime, ops, ops_default, ist, jst, nst
886 type(wfs_elec_t) :: xx, yy
887 real(real64), allocatable :: tmp(:)
888 real(real64), allocatable :: ddotv(:)
889 complex(real64), allocatable :: zdotv(:)
890 real(real64), allocatable :: ddot(:,:)
891 complex(real64), allocatable :: zdot(:,:)
892
893
894 push_sub(test_batch_ops)
895
896 !%Variable TestBatchOps
897 !%Type flag
898 !%Default ops_axpy + ops_scal + ops_nrm2
899 !%Section Calculation Modes::Test
900 !%Description
901 !% Decides which part of the Hamiltonian is applied.
902 !%Option ops_axpy bit(1)
903 !% Tests batch_axpy operation
904 !%Option ops_scal bit(2)
905 !% Tests batch_scal operation
906 !%Option ops_nrm2 bit(3)
907 !% Tests batch_nrm2 operation
908 !%Option ops_dotp_matrix bit(4)
909 !% Tests X(mesh_batch_dotp_matrix)
910 !%Option ops_dotp_self bit(5)
911 !% Tests X(mesh_batch_dotp_self)
912 !%Option ops_dotp_vector bit(6)
913 !% Tests X(mesh_batch_dotp_vector)
914 !%End
915 ops_default = &
916 option__testbatchops__ops_axpy + &
917 option__testbatchops__ops_scal + &
918 option__testbatchops__ops_nrm2 + &
919 option__testbatchops__ops_dotp_matrix + &
920 option__testbatchops__ops_dotp_self + &
921 option__testbatchops__ops_dotp_vector
922
923 call parse_variable(namespace, 'TestBatchOps', ops_default, ops)
924
925 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
926
927 call messages_write('Info: Testing batch operations')
928 call messages_new_line()
929 call messages_new_line()
930 call messages_info(namespace=namespace)
931
932 sys => electrons_t(namespace, generate_epot=.false.)
933 call sys%init_parallelization(mpi_world)
934
935 call states_elec_allocate_wfns(sys%st, sys%gr)
936 call test_batch_set_gaussian(sys%st%group%psib(1, 1), sys%gr)
937 if (sys%st%pack_states) call sys%st%pack()
938
939 if (bitand(ops, option__testbatchops__ops_axpy) /= 0) then
940 message(1) = 'Info: Testing axpy'
941 call messages_info(1, namespace=sys%namespace)
942
943 call sys%st%group%psib(1, 1)%copy_to(xx, copy_data = .true.)
944 call sys%st%group%psib(1, 1)%copy_to(yy, copy_data = .true.)
945
946 do itime = 1, param%repetitions
947 call batch_axpy(sys%gr%np, 0.1_real64, xx, yy)
948 end do
949 call test_prints_info_batch(sys%st, sys%gr, yy, string = "axpy")
950
951 call xx%end()
952 call yy%end()
953 end if
954
955 if (bitand(ops, option__testbatchops__ops_scal) /= 0) then
956 message(1) = 'Info: Testing scal'
957 call messages_info(1, namespace=sys%namespace)
958
959 call sys%st%group%psib(1, 1)%copy_to(xx, copy_data = .true.)
960 call sys%st%group%psib(1, 1)%copy_to(yy, copy_data = .true.)
961
962 do itime = 1, param%repetitions
963 call batch_scal(sys%gr%np, 0.1_real64, yy)
964 end do
965 call test_prints_info_batch(sys%st, sys%gr, yy, string="scal")
966
967 call xx%end()
968 call yy%end()
969 end if
970
971 if (bitand(ops, option__testbatchops__ops_nrm2) /= 0) then
972 message(1) = 'Info: Testing nrm2'
973 call messages_info(1, namespace=sys%namespace)
974
975 call sys%st%group%psib(1, 1)%copy_to(xx, copy_data = .true.)
976 call sys%st%group%psib(1, 1)%copy_to(yy, copy_data = .true.)
977
978 safe_allocate(tmp(1:xx%nst))
979
980 do itime = 1, param%repetitions
981 call mesh_batch_nrm2(sys%gr, yy, tmp)
982 end do
983 do itime = 1, xx%nst
984 write(message(1),'(a,i1,3x,e23.16)') "Nrm2 norm state ", itime, tmp(itime)
985 call messages_info(1, namespace=sys%namespace)
986 end do
987
988 safe_deallocate_a(tmp)
989
990 call xx%end()
991 call yy%end()
992 end if
993
994 if (bitand(ops, option__testbatchops__ops_dotp_matrix) /= 0) then
995
996 message(1) = 'Info: Testing dotp_matrix'
997 call messages_info(1, namespace=sys%namespace)
998
999 call sys%st%group%psib(1, 1)%copy_to(xx, copy_data = .true.)
1000 call sys%st%group%psib(1, 1)%copy_to(yy, copy_data = .true.)
1001
1002 nst = sys%st%group%psib(1, 1)%nst
1003
1004 if (states_are_real(sys%st)) then
1005 safe_allocate(ddot(1:nst, 1:nst))
1006 do itime = 1, param%repetitions
1007 call dmesh_batch_dotp_matrix(sys%gr, xx, yy, ddot)
1008 end do
1009
1010 do ist = 1, nst
1011 do jst = 1, nst
1012 write(message(jst), '(a,2i3,3x,e23.16)') 'Dotp_matrix states', ist, jst, ddot(ist,jst)
1013 end do
1014 call messages_info(nst, namespace=sys%namespace)
1015 end do
1016 safe_deallocate_a(ddot)
1017 else
1018 safe_allocate(zdot(1:nst, 1:nst))
1019 do itime = 1, param%repetitions
1020 call zmesh_batch_dotp_matrix(sys%gr, xx, yy, zdot)
1021 end do
1022
1023 do ist = 1, nst
1024 do jst = 1, nst
1025 write(message(jst), '(a,2i3,3x,2(e23.16,1x))') 'Dotp_matrix states', ist, jst, zdot(ist,jst)
1026 end do
1027 call messages_info(nst, namespace=sys%namespace)
1028 end do
1029 safe_deallocate_a(zdot)
1030 end if
1031
1032 call xx%end()
1033 call yy%end()
1034 end if
1035
1036 if (bitand(ops, option__testbatchops__ops_dotp_vector) /= 0) then
1037
1038 message(1) = 'Info: Testing dotp_vector'
1039 call messages_info(1, namespace=sys%namespace)
1040
1041 call sys%st%group%psib(1, 1)%copy_to(xx, copy_data = .true.)
1042 call sys%st%group%psib(1, 1)%copy_to(yy, copy_data = .true.)
1043
1044 nst = sys%st%group%psib(1, 1)%nst
1045
1046 if (states_are_real(sys%st)) then
1047 safe_allocate(ddotv(1:nst))
1048 do itime = 1, param%repetitions
1049 call dmesh_batch_dotp_vector(sys%gr, xx, yy, ddotv)
1050 end do
1051
1052 do ist = 1, nst
1053 write(message(ist), '(a,i3,3x,e23.16)') 'Dotp_vector state', ist, ddotv(ist)
1054 end do
1055 call messages_info(nst, namespace=sys%namespace)
1056 safe_deallocate_a(ddotv)
1057 else
1058 safe_allocate(zdotv(1:nst))
1059 do itime = 1, param%repetitions
1060 call zmesh_batch_dotp_vector(sys%gr, xx, yy, zdotv)
1061 end do
1062
1063 do ist = 1, nst
1064 write(message(ist), '(a,i3,3x,2(e23.16,1x))') 'Dotp_vector state', ist, zdotv(ist)
1065 end do
1066 call messages_info(nst, namespace=sys%namespace)
1067 safe_deallocate_a(zdotv)
1068 end if
1069
1070 call xx%end()
1071 call yy%end()
1072 end if
1073
1074 if (bitand(ops, option__testbatchops__ops_dotp_self) /= 0) then
1075
1076 message(1) = 'Info: Testing dotp_self'
1077 call messages_info(1, namespace=sys%namespace)
1078
1079 call sys%st%group%psib(1, 1)%copy_to(xx, copy_data = .true.)
1080
1081 nst = sys%st%group%psib(1, 1)%nst
1082
1083 if (states_are_real(sys%st)) then
1084 safe_allocate(ddot(1:nst, 1:nst))
1085 do itime = 1, param%repetitions
1086 call dmesh_batch_dotp_self(sys%gr, xx, ddot)
1087 end do
1088
1089 do ist = 1, nst
1090 do jst = 1, nst
1091 write(message(jst), '(a,2i3,3x,e23.16)') 'Dotp_self states', ist, jst, ddot(ist,jst)
1092 end do
1093 call messages_info(nst, namespace=sys%namespace)
1094 end do
1095 safe_deallocate_a(ddot)
1096 else
1097 safe_allocate(zdot(1:nst, 1:nst))
1098 do itime = 1, param%repetitions
1099 call zmesh_batch_dotp_self(sys%gr, xx, zdot)
1100 end do
1101
1102 do ist = 1, nst
1103 do jst = 1, nst
1104 write(message(jst), '(a,2i3,3x,2(e23.16,1x))') 'Dotp_self states', ist, jst, zdot(ist,jst)
1105 end do
1106 call messages_info(nst, namespace=sys%namespace)
1107 end do
1108 safe_deallocate_a(zdot)
1109 end if
1110
1111 call xx%end()
1112 end if
1113
1114 call states_elec_deallocate_wfns(sys%st)
1115 safe_deallocate_p(sys)
1116
1117 pop_sub(test_batch_ops)
1118 end subroutine test_batch_ops
1119
1120! ---------------------------------------------------------
1121 subroutine test_derivatives(param, namespace)
1122 type(test_parameters_t), intent(in) :: param
1123 type(namespace_t), intent(in) :: namespace
1124
1125 type(electrons_t), pointer :: sys
1126
1127 push_sub(test_derivatives)
1128
1129 sys => electrons_t(namespace, generate_epot=.false.)
1130 call sys%init_parallelization(mpi_world)
1131
1132 message(1) = 'Info: Testing the finite-differences derivatives.'
1133 message(2) = ''
1134 call messages_info(2, namespace=namespace)
1135
1136 if (param%type == option__testtype__all .or. param%type == option__testtype__real) then
1137 call dderivatives_test(sys%gr%der, sys%namespace, param%repetitions, param%min_blocksize, param%max_blocksize)
1138 end if
1139
1140 if (param%type == option__testtype__all .or. param%type == option__testtype__complex) then
1141 call zderivatives_test(sys%gr%der, sys%namespace, param%repetitions, param%min_blocksize, param%max_blocksize)
1142 end if
1143
1144 safe_deallocate_p(sys)
1145
1146 pop_sub(test_derivatives)
1147 end subroutine test_derivatives
1148
1149 ! ---------------------------------------------------------
1150
1151 subroutine test_orthogonalization(param, namespace)
1152 type(test_parameters_t), intent(in) :: param
1153 type(namespace_t), intent(in) :: namespace
1154
1155 type(electrons_t), pointer :: sys
1156 integer :: itime
1157
1158 push_sub(test_orthogonalization)
1159
1160 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
1161 call calc_mode_par%set_scalapack_compat()
1162
1163 sys => electrons_t(namespace, generate_epot=.false.)
1164 call sys%init_parallelization(mpi_world)
1165
1166 message(1) = 'Info: Testing orthogonalization.'
1167 message(2) = ''
1168 call messages_info(2, namespace=namespace)
1169
1170 if (param%type == option__testtype__all .or. param%type == option__testtype__real) then
1171 message(1) = 'Info: Real wave-functions.'
1172 call messages_info(1, namespace=namespace)
1173 do itime = 1, param%repetitions
1174 call dstates_elec_calc_orth_test(sys%st, sys%namespace, sys%gr, sys%kpoints)
1175 end do
1176 end if
1177
1178 if (param%type == option__testtype__all .or. param%type == option__testtype__complex) then
1179 message(1) = 'Info: Complex wave-functions.'
1180 call messages_info(1, namespace=namespace)
1181 do itime = 1, param%repetitions
1182 call zstates_elec_calc_orth_test(sys%st, sys%namespace, sys%gr, sys%kpoints)
1183 end do
1184 end if
1185
1186 safe_deallocate_p(sys)
1187
1188 pop_sub(test_orthogonalization)
1189 end subroutine test_orthogonalization
1190
1191 ! ---------------------------------------------------------
1192
1193 subroutine test_interpolation(param, namespace)
1194 type(test_parameters_t), intent(in) :: param
1195 type(namespace_t), intent(in) :: namespace
1196
1197 type(electrons_t), pointer :: sys
1198
1199 push_sub(test_interpolation)
1200
1201 sys => electrons_t(namespace, generate_epot=.false.)
1202 call sys%init_parallelization(mpi_world)
1203
1204 if (param%type == option__testtype__all .or. param%type == option__testtype__real) then
1205 call messages_write('Info: Testing real interpolation routines')
1206 call messages_new_line()
1207 call messages_new_line()
1208 call messages_info(namespace=namespace)
1209
1210 call dmesh_interpolation_test(sys%gr)
1211 end if
1212
1213 if (param%type == option__testtype__all .or. param%type == option__testtype__complex) then
1215 call messages_write('Info: Testing complex interpolation routines')
1216 call messages_new_line()
1217 call messages_new_line()
1218 call messages_info(namespace=namespace)
1219
1220 call zmesh_interpolation_test(sys%gr)
1221 end if
1222
1223 safe_deallocate_p(sys)
1224
1225 pop_sub(test_interpolation)
1226 end subroutine test_interpolation
1227
1228
1229 ! ---------------------------------------------------------
1230
1231 subroutine test_ion_interaction(namespace)
1232 type(namespace_t), intent(in) :: namespace
1233
1234 type(electrons_t), pointer :: sys
1235
1236 push_sub(test_ion_interaction)
1237
1238 sys => electrons_t(namespace, generate_epot=.false.)
1239 call sys%init_parallelization(mpi_world)
1240
1241 call ion_interaction_test(sys%space, sys%ions%latt, sys%ions%atom, sys%ions%natoms, sys%ions%pos, &
1242 sys%gr%box%bounding_box_l, namespace, sys%mc)
1243
1244 safe_deallocate_p(sys)
1245
1246 pop_sub(test_ion_interaction)
1247 end subroutine test_ion_interaction
1248
1249 ! ---------------------------------------------------------
1250
1251 subroutine test_prints_info_batch(st, gr, psib, string)
1252 type(states_elec_t), intent(in) :: st
1253 type(grid_t), intent(in) :: gr
1254 class(batch_t), intent(inout) :: psib
1255 character(*), optional, intent(in) :: string
1256
1257 integer :: itime
1258 complex(real64), allocatable :: zpsi(:, :)
1259 real(real64), allocatable :: dpsi(:, :)
1260
1261 character(80) :: string_
1262
1263 string_ = optional_default(string, "")
1264
1265 push_sub(test_prints_info_batch)
1266
1267 if (states_are_real(st)) then
1268 safe_allocate(dpsi(1:gr%np, 1:st%d%dim))
1269 else
1270 safe_allocate(zpsi(1:gr%np, 1:st%d%dim))
1271 end if
1272
1273 do itime = 1, psib%nst
1274 if (states_are_real(st)) then
1275 call batch_get_state(psib, itime, gr%np, dpsi)
1276 write(message(1),'(a,i2,3x,e23.16)') "Norm state "//trim(string_)//" ", itime, dmf_nrm2(gr, st%d%dim, dpsi)
1277 else
1278 call batch_get_state(psib, itime, gr%np, zpsi)
1279 write(message(1),'(a,i2,3x,e23.16)') "Norm state "//trim(string_)//" ", itime, zmf_nrm2(gr, st%d%dim, zpsi)
1280 end if
1281 call messages_info(1)
1282 end do
1283
1284 if (states_are_real(st)) then
1285 safe_deallocate_a(dpsi)
1286 else
1287 safe_deallocate_a(zpsi)
1288 end if
1289
1290 pop_sub(test_prints_info_batch)
1291
1292 end subroutine test_prints_info_batch
1293
1294
1295 ! ---------------------------------------------------------
1296 subroutine test_clock()
1297
1298 type(iteration_counter_t) :: clock, other_clock
1299
1300 push_sub(test_clock)
1301
1302 ! Clock operations
1303 write(message(1),'(a)') " Operation Counter Time Global step"
1304 call messages_info(1)
1305
1306 clock = clock_t(time_step=m_two, initial_iteration=2)
1307 call write_clock("initialization")
1308
1309 clock = clock + 5
1310 call write_clock("addition (+5)")
1311
1312 clock = clock - 3
1313 call write_clock("subtraction (-3)")
1314
1315 call clock%reset()
1316 call write_clock("reset (=0)")
1317
1318 call clock%set(clock_t(time_step=m_four, initial_iteration=1))
1319 call write_clock("set (=2)")
1320
1321 ! Clock comparisons:
1322 write(message(1),'(a)') " Clock comparisons:"
1323 call messages_info(1)
1325 other_clock = clock_t(time_step=m_one, initial_iteration=5)
1326
1327 call write_condition_result("4 < 5", clock < other_clock)
1328 call write_condition_result("4 <= 5", clock <= other_clock)
1329 call write_condition_result("4 > 5", clock > other_clock)
1330 call write_condition_result("4 >= 5", clock >= other_clock)
1331 call write_condition_result("4 == 5", clock == other_clock)
1332 call write_condition_result("4 /= 5", clock /= other_clock)
1333
1334 pop_sub(test_clock)
1335 contains
1336 subroutine write_clock(operation)
1337 character(len=*), intent(in) :: operation
1338
1339 write(message(1),'(a17,x,i6,x,f10.1,x,i16)') operation, clock%counter(), clock%value(), clock%global_step()
1340 call messages_info(1)
1341 end subroutine write_clock
1342
1343 subroutine write_condition_result(condition, result)
1344 character(len=*), intent(in) :: condition
1345 logical, intent(in) :: result
1346
1347 write(message(1),'(a10," = ",i1," (",l1,")")') condition, abs(transfer(result, 0)), result
1348 call messages_info(1)
1349 end subroutine write_condition_result
1350
1351 end subroutine test_clock
1352
1353
1354 ! ---------------------------------------------------------
1355 subroutine test_cgal()
1356
1357 type(cgal_polyhedra_t) :: cgal_poly
1358
1359 push_sub(test_cgal)
1360
1361 call cgal_polyhedron_init(cgal_poly, "28-cgal.02-X.off", verbose = .true.)
1362
1363 if (cgal_polyhedron_point_inside(cgal_poly, 30._real64, 10._real64, 30._real64)) then
1364 message(1) = "cgal_polyhedron_point_inside"
1365 call messages_info(1)
1366 end if
1367
1368 call cgal_polyhedron_end(cgal_poly)
1369
1370 pop_sub(test_cgal)
1371 end subroutine test_cgal
1372
1373
1374 ! ---------------------------------------------------------
1375 subroutine test_dense_eigensolver()
1376 integer :: N, ii, jj, N_list(4), i_N
1377 real(real64), allocatable :: matrix(:, :), eigenvectors(:, :), eigenvalues(:), test(:)
1378 real(real64), allocatable :: differences(:)
1379
1380 push_sub(test_dense_eigensolver)
1381
1382 n_list = [15, 32, 100, 500]
1383
1384 do i_n = 1, 4
1385 n = n_list(i_n)
1386 safe_allocate(matrix(1:n, 1:n))
1387 safe_allocate(eigenvectors(1:n, 1:n))
1388 safe_allocate(eigenvalues(1:n))
1389 safe_allocate(test(1:n))
1390 safe_allocate(differences(1:n))
1391
1392
1393 ! set up symmetrix matrix
1394 do jj = 1, n
1395 do ii = 1, n
1396 matrix(ii, jj) = ii * jj
1397 end do
1398 end do
1399
1400 ! parallel solver
1401 eigenvectors(1:n, 1:n) = matrix(1:n, 1:n)
1402 call lalg_eigensolve_parallel(n, eigenvectors, eigenvalues)
1403
1404 do ii = 1, n
1405 test(:) = matmul(matrix, eigenvectors(:, ii)) - eigenvalues(ii) * eigenvectors(:, ii)
1406 differences(ii) = sum(abs(test)) / sum(abs(eigenvectors(:, ii)))
1407 end do
1408 write(message(1), "(A, I3, A, E13.6)") "Parallel solver - N: ", n, &
1409 ", average difference: ", sum(differences)/n
1410 call messages_info(1)
1411
1412 ! serial solver
1413 eigenvectors(1:n, 1:n) = matrix(1:n, 1:n)
1414 call lalg_eigensolve(n, eigenvectors, eigenvalues)
1415
1416 do ii = 1, n
1417 test(:) = matmul(matrix, eigenvectors(:, ii)) - eigenvalues(ii) * eigenvectors(:, ii)
1418 differences(ii) = sum(abs(test)) / sum(abs(eigenvectors(:, ii)))
1419 end do
1420 write(message(1), "(A, I3, A, E13.6)") "Serial solver - N: ", n, &
1421 ", average difference: ", sum(differences)/n
1422 call messages_info(1)
1423
1424 safe_deallocate_a(matrix)
1425 safe_deallocate_a(eigenvectors)
1426 safe_deallocate_a(eigenvalues)
1427 safe_deallocate_a(test)
1428 safe_deallocate_a(differences)
1429 end do
1430
1431 pop_sub(test_dense_eigensolver)
1432 end subroutine test_dense_eigensolver
1433
1434 subroutine test_batch_set_gaussian(psib, mesh)
1435 class(batch_t), intent(inout) :: psib
1436 class(mesh_t), intent(in) :: mesh
1437
1438 real(real64), allocatable :: dff(:)
1439 complex(real64), allocatable :: zff(:)
1440 integer :: ist, ip
1441 real(real64) :: da, db, dc
1442 complex(real64) :: za, zb, zc
1443
1444 push_sub(test_batch_set_gaussian)
1445
1446 ! use a similar function as in the derivatives test
1447 da = m_one/mesh%box%bounding_box_l(1)
1448 db = 10.0_real64
1449 dc = 100.0_real64
1450
1451 if (type_is_complex(psib%type())) then
1452 ! we make things more "complex"
1453 za = da + m_zi*0.01_real64
1454 zb = db*exp(m_zi*0.345_real64)
1455 zc = dc - m_zi*50.0_real64
1456
1457 safe_allocate(zff(1:mesh%np))
1458 do ist = 1, psib%nst_linear
1459 za = za * ist
1460 zb = zb / ist
1461 do ip = 1, mesh%np
1462 zff(ip) = zb*exp(-za*sum(mesh%x(ip, :)**2)) + zc
1463 end do
1464 call batch_set_state(psib, ist, mesh%np, zff)
1465 end do
1466 call zmesh_batch_normalize(mesh, psib)
1467 safe_deallocate_a(zff)
1468 else
1469 safe_allocate(dff(1:mesh%np))
1470 do ist = 1, psib%nst_linear
1471 da = da * ist
1472 db = db / ist
1473 do ip = 1, mesh%np
1474 dff(ip) = db*exp(-da*sum(mesh%x(ip, :)**2)) + dc
1475 end do
1476 call batch_set_state(psib, ist, mesh%np, dff)
1477 end do
1478 call dmesh_batch_normalize(mesh, psib)
1479 safe_deallocate_a(dff)
1480 end if
1481
1483 end subroutine test_batch_set_gaussian
1484
1485 ! ---------------------------------------------------------
1486 subroutine test_grid_interpolation()
1487 type(electrons_t), pointer :: sys
1488 type(multigrid_t) :: mgrid
1489
1490 push_sub(test_grid_interpolation)
1491
1492 sys => electrons_t(global_namespace, generate_epot=.false.)
1493 call sys%init_parallelization(mpi_world)
1494
1495 call multigrid_init(mgrid, global_namespace, sys%space, sys%gr, sys%gr%der, &
1496 sys%gr%stencil, sys%mc, nlevels=3)
1497
1498 call multigrid_test_interpolation(mgrid, sys%space)
1499
1500 call multigrid_end(mgrid)
1501
1502 safe_deallocate_p(sys)
1503
1505 end subroutine test_grid_interpolation
1506
1507 subroutine test_iihash()
1508
1509 type(iihash_t) :: h
1510 integer :: value
1511 logical :: found
1512
1513 call iihash_init(h)
1514 call iihash_insert(h, 1,5)
1515 call iihash_insert(h, 3,7)
1516
1517 value = iihash_lookup(h, 1, found)
1518 write(message(1),*) 'hash[1] :', found, value
1519 call messages_info()
1520
1521 value = iihash_lookup(h, 2, found)
1522 write(message(1),*) 'hash[2] :', found, value
1523 call messages_info()
1524
1525 value = iihash_lookup(h, 3, found)
1526 write(message(1),*) 'hash[3] :', found, value
1528
1529 call iihash_end(h)
1530
1531 end subroutine test_iihash
1532
1533 subroutine test_sihash()
1534 type(sihash_t) :: h
1535 type(sihash_iterator_t) :: it
1536
1537 integer :: counter
1538 integer :: value, sum
1539 logical :: found
1540
1541 call sihash_init(h)
1542 call sihash_insert(h, "one",5)
1543 call sihash_insert(h, "three",7)
1544 call sihash_insert(h, "the answer", 42)
1545
1546 value = sihash_lookup(h, "one", found)
1547 write(message(1),*) 'hash["one"]: ', found, value
1548 call messages_info()
1549
1550 value = sihash_lookup(h, "two", found)
1551 write(message(1),*) 'hash["two"]: ', found, value
1552 call messages_info()
1553
1554 value = sihash_lookup(h, "three", found)
1555 write(message(1),*) 'hash["three"]: ', found, value
1556 call messages_info()
1557
1558 sum = 0
1559 counter = 1
1560 call it%start(h)
1561
1562 do while (it%has_next())
1563 value = it%get_next()
1564 sum = sum + value
1565 write(message(1),'(I3,A,I5)') counter,': hash[...] = ',value
1566 call messages_info()
1567 counter = counter + 1
1568 end do
1569 write(message(1),*) 'counter = ', counter
1570 write(message(2),*) 'sum = ', sum
1571 call messages_info(2)
1572
1573
1574 call sihash_end(h)
1575
1576 end subroutine test_sihash
1577
1578
1579 subroutine test_sphash(namespace)
1580 type(namespace_t), intent(in) :: namespace
1581
1582 type(sphash_t) :: h
1583 type(sphash_iterator_t) :: it
1584
1585 logical :: found
1586
1587
1588 type(iteration_counter_t) :: clock_1
1589 type(iteration_counter_t), allocatable :: clock_2
1590 type(space_t) :: space_1
1591
1592 class(*), pointer :: value
1593
1594 integer :: count_clock, count_space
1595
1596 safe_allocate(clock_2)
1597
1598 clock_1 = clock_t(1.d-5, 1)
1599 clock_2 = clock_t(2.d-5, 2)
1600 space_1 = space_t(namespace)
1601
1602 call sphash_init(h)
1603 call sphash_insert(h, "one", clock_1)
1604 call sphash_insert(h, "two", space_1)
1605 call sphash_insert(h, "three", clock_2, clone=.true.)
1606
1607 value => sphash_lookup(h, "one", found)
1608 select type(value)
1609 type is (iteration_counter_t)
1610 write(message(1),*) 'hash["one"]: ', found, value%counter()
1611 call messages_info()
1612 type is (space_t)
1613 write(message(1),*) 'hash["one"]: ', found, value%short_info()
1614 call messages_info()
1615 class default
1616 write(message(1),*) 'wrong type. found = ', found
1617 call messages_info()
1618 end select
1619
1620 value => sphash_lookup(h, "two", found)
1621 select type(value)
1622 type is (iteration_counter_t)
1623 write(message(1),*) 'hash["two"]: ', found, value%counter()
1624 call messages_info()
1625 type is (space_t)
1626 write(message(1),*) 'hash["two"]: ', found, value%short_info()
1627 call messages_info()
1628 class default
1629 write(message(1),*) 'wrong type. found = ',found
1630 call messages_info()
1631 end select
1632
1633 safe_deallocate_a(clock_2)
1634
1635 value => sphash_lookup(h, "three", found)
1636 select type(value)
1637 type is (iteration_counter_t)
1638 write(message(1),*) 'hash["three"]: ', found, value%counter()
1639 call messages_info()
1640 type is (space_t)
1641 write(message(1),*) 'hash["three"]: ', found, value%short_info()
1642 call messages_info()
1643 class default
1644 write(message(1),*) 'wrong type. found = ',found
1645 call messages_info()
1646 end select
1647
1648 count_clock = 0
1649 count_space = 0
1650
1651 call it%start(h)
1652
1653 do while (it%has_next())
1654 value => it%get_next()
1655 select type(value)
1656 type is (iteration_counter_t)
1657 count_clock = count_clock + 1
1658 type is (space_t)
1659 count_space = count_space + 1
1660 end select
1661 end do
1662
1663 write(message(1), *) 'Count_clock = ', count_clock
1664 write(message(2), *) 'Count_space = ', count_space
1665 call messages_info(2)
1666
1667 call sphash_end(h)
1668
1669 end subroutine test_sphash
1670
1671 ! ---------------------------------------------------------
1672 subroutine test_regridding(namespace)
1673 type(namespace_t), intent(in) :: namespace
1674
1675 type(electrons_t), pointer :: sysA, sysB
1676 type(regridding_t), pointer :: regridding
1677 real(real64), allocatable :: ff_A(:), ff_A_reference(:), ff_B(:), ff_B_reference(:), diff_A(:), diff_B(:)
1678 real(real64) :: norm_ff, norm_diff
1679 integer :: ip, ierr
1680
1681 push_sub(test_regridding)
1682
1683 sysa => electrons_t(namespace_t("A", namespace), generate_epot=.false.)
1684 sysb => electrons_t(namespace_t("B", namespace), generate_epot=.false.)
1685
1686 call calc_mode_par%set_parallelization(p_strategy_states, default=.true.)
1687 call sysa%init_parallelization(mpi_world)
1688 call sysb%init_parallelization(mpi_world)
1689
1690
1691 safe_allocate(ff_a(1:sysa%gr%np))
1692 safe_allocate(ff_a_reference(1:sysa%gr%np))
1693 safe_allocate(diff_a(1:sysa%gr%np))
1694 safe_allocate(ff_b(1:sysb%gr%np))
1695 safe_allocate(ff_b_reference(1:sysb%gr%np))
1696 safe_allocate(diff_b(1:sysb%gr%np))
1697
1698 do ip = 1, sysa%gr%np
1699 ff_a_reference(ip) = values(sysa%gr%x(ip, :))
1700 end do
1701 do ip = 1, sysb%gr%np
1702 ff_b_reference(ip) = values(sysb%gr%x(ip, :))
1703 end do
1704
1705 ! forward mapping A->B
1706 regridding => regridding_t(sysb%gr, sysa%gr, sysa%space, sysa%namespace)
1707 call regridding%do_transfer(ff_b, ff_a_reference)
1708 safe_deallocate_p(regridding)
1709
1710 ! check that mapped function is exactly the same for the points that are available on both grids
1711 do ip = 1, sysb%gr%np
1712 if (abs(ff_b(ip)) <= m_epsilon) then
1713 ff_b_reference(ip) = m_zero
1714 diff_b(ip) = m_zero
1715 else
1716 diff_b(ip) = abs(ff_b_reference(ip) - ff_b(ip))
1717 end if
1718 end do
1719 norm_ff = dmf_nrm2(sysb%gr, ff_b_reference)
1720 norm_diff = dmf_nrm2(sysb%gr, diff_b)
1721
1722 write(message(1),'(a, E14.6)') "Forward: difference of reference to mapped function (rel.): ", &
1723 norm_diff/norm_ff
1724 call messages_info(1, namespace=namespace)
1725
1726 call dio_function_output(io_function_fill_how('AxisX'), ".", "forward_reference", namespace, sysb%space, &
1727 sysb%gr, ff_b_reference, unit_one, ierr)
1728 call dio_function_output(io_function_fill_how('AxisX'), ".", "forward_mapped", namespace, sysb%space, &
1729 sysb%gr, ff_b, unit_one, ierr)
1730 call dio_function_output(io_function_fill_how('AxisX'), ".", "forward_original", namespace, sysa%space, &
1731 sysa%gr, ff_a_reference, unit_one, ierr)
1732
1733 ! backward mapping B->A
1734 regridding => regridding_t(sysa%gr, sysb%gr, sysb%space, sysb%namespace)
1735 call regridding%do_transfer(ff_a, ff_b_reference)
1736 safe_deallocate_p(regridding)
1737
1738 do ip = 1, sysa%gr%np
1739 if (abs(ff_a(ip)) <= m_epsilon) then
1740 ff_a_reference(ip) = m_zero
1741 diff_a(ip) = m_zero
1742 else
1743 diff_a(ip) = abs(ff_a_reference(ip) - ff_a(ip))
1744 end if
1745 end do
1746 norm_ff = dmf_nrm2(sysa%gr, ff_a_reference)
1747 norm_diff = dmf_nrm2(sysa%gr, diff_a)
1748
1749 write(message(1),'(a, E14.6)') "Backward: difference of reference to mapped function (rel.): ", &
1750 norm_diff/norm_ff
1751 call messages_info(1, namespace=namespace)
1752
1753 call dio_function_output(io_function_fill_how('AxisX'), ".", "backward_reference", namespace, sysa%space, &
1754 sysa%gr, ff_a_reference, unit_one, ierr)
1755 call dio_function_output(io_function_fill_how('AxisX'), ".", "backward_mapped", namespace, sysa%space, &
1756 sysa%gr, ff_a, unit_one, ierr)
1757 call dio_function_output(io_function_fill_how('AxisX'), ".", "backward_original", namespace, sysb%space, &
1758 sysb%gr, ff_b_reference, unit_one, ierr)
1759
1760 safe_deallocate_a(ff_a)
1761 safe_deallocate_a(ff_a_reference)
1762 safe_deallocate_a(ff_b)
1763 safe_deallocate_a(ff_b_reference)
1764 safe_deallocate_a(diff_a)
1765 safe_deallocate_a(diff_b)
1766 safe_deallocate_p(sysa)
1767 safe_deallocate_p(sysb)
1768
1769 call messages_info(1, namespace=namespace)
1770
1771 pop_sub(test_regridding)
1772 contains
1773 real(real64) function values(xx)
1774 real(real64), intent(in) :: xx(:)
1775 real(real64) :: xx0(1:size(xx, dim=1))
1776 real(real64), parameter :: aa = m_half
1777 real(real64), parameter :: bb = m_four
1778
1779 ! no push_sub/pop_sub because this function is called often
1780 xx0(:) = m_one
1781 values = bb * exp(-aa*sum((xx-xx0)**2))
1782 end function values
1783 end subroutine test_regridding
1784
1790 subroutine test_vecpot_analytical(namespace)
1791 type(namespace_t), intent(in) :: namespace
1792
1793 class(maxwell_t), pointer :: maxwell_system
1794
1795 real(real64), allocatable :: magnetic_field(:,:)
1796 real(real64), allocatable :: vector_potential_mag(:,:)
1797 real(real64), allocatable :: vector_potential_analytical(:,:)
1798 real(real64), allocatable :: delta(:,:)
1799 real(real64) :: exp_factor
1800 real(real64) :: xx
1801 real(real64) :: yy
1802 real(real64) :: zz
1803 real(real64) :: sigma
1804 integer :: ip, j, ierr, nn
1805 integer(int64) :: out_how
1806 character(len=MAX_PATH_LEN) :: fname, fname2, fname3
1807
1808 out_how = 32
1809 maxwell_system => maxwell_t(namespace)
1810 sigma = maxwell_system%gr%box%bounding_box_l(1)/10_real64 ! this is exponential width
1811 call maxwell_system%init_parallelization(mpi_world)
1812
1813 safe_allocate(magnetic_field(1:maxwell_system%gr%np_part, 1:3))
1814 safe_allocate(vector_potential_mag(1:maxwell_system%gr%np_part, 1:3))
1815 safe_allocate(vector_potential_analytical(1:maxwell_system%gr%np_part, 1:3))
1816 safe_allocate(delta(1:maxwell_system%gr%np, 1:3))
1817
1818 !%Variable TestVectorPotentialType
1819 !%Type integer
1820 !%Default bounded
1821 !%Section Calculation Modes::Test
1822 !%Description
1823 !% Select whether bounded or unbounded type will be used for vector potential tests
1824 !%Option bounded 1
1825 !% Analytical Vector Potential formulation is bounded by spatial gaussian
1826 !%Option unbounded 2
1827 !% Analytical Vector Potential is not bounded
1828 !%End
1829 call parse_variable(namespace, 'TestVectorPotentialType', option__testvectorpotentialtype__bounded, nn)
1830
1831 select case (nn)
1832 case (option__testvectorpotentialtype__bounded) ! bounded input
1833 do ip = 1, maxwell_system%gr%np_part
1834 xx = maxwell_system%gr%x(ip, 1)
1835 yy = maxwell_system%gr%x(ip, 2)
1836 zz = maxwell_system%gr%x(ip, 3)
1837 exp_factor = exp((-xx**2 - yy**2 - zz**2)*1/(2*sigma**2))
1838 magnetic_field(ip, 1) = exp_factor*yy*(1 - (-xx**2 + yy**2)/(3*sigma**2) - zz**2/(3*sigma**2))
1839 magnetic_field(ip, 2) = exp_factor * xx * (1 + (-xx**2 + yy**2)/(3*sigma**2) - zz**2/(3*sigma**2))
1840 magnetic_field(ip, 3) = exp_factor * 2 * xx * yy * zz * 1/(3*sigma**2)
1841
1842 vector_potential_analytical(ip, 1) = m_third * xx * zz * exp_factor
1843 vector_potential_analytical(ip, 2) = - m_third * yy * zz * exp_factor
1844 vector_potential_analytical(ip, 3) = m_third * (-xx**2 + yy**2) * exp_factor
1845 end do
1846 case (option__testvectorpotentialtype__unbounded) ! unbounded input, TODO this unit test requires implementation of BCs for Helmholtz decomposition
1847 do ip = 1, maxwell_system%gr%np_part
1848 magnetic_field(ip, 1) = maxwell_system%gr%x(ip, 2)
1849 magnetic_field(ip, 2) = maxwell_system%gr%x(ip, 1)
1850 magnetic_field(ip, 3) = m_zero
1851
1852 vector_potential_analytical(ip, 1) = m_third * maxwell_system%gr%x(ip, 1) * maxwell_system%gr%x(ip, 3)
1853 vector_potential_analytical(ip, 2) = - m_third * maxwell_system%gr%x(ip, 2) * maxwell_system%gr%x(ip, 3)
1854 vector_potential_analytical(ip, 3) = - m_third * (maxwell_system%gr%x(ip, 1)**2 - maxwell_system%gr%x(ip, 2)**2)
1855 end do
1856 end select
1857 call maxwell_system%helmholtz%get_vector_potential(namespace, vector_potential_mag, magnetic_field)
1858
1859 do j = 1, 3
1860 delta(:,:) = m_zero
1861 do ip = 1, maxwell_system%gr%np
1862 delta(ip,j) = vector_potential_analytical(ip, j) - vector_potential_mag(ip, j)
1863 end do
1864 end do
1865
1866 do j = 1, 3
1867 write(message(j),*) 'j, norm2(delta)', j, norm2(delta(:,j))
1868 end do
1869 call messages_info(3)
1870
1871 write(fname, '(a)') 'deviation_from_analytical_formulation' ! use messages later
1872 call io_function_output_vector(out_how , './' , trim(fname), namespace, maxwell_system%space, maxwell_system%gr, &
1873 delta, unit_one, ierr)
1874 write(fname2, '(a)') 'vector_potential_analytical'
1875 call io_function_output_vector(out_how , './' , trim(fname2), namespace, maxwell_system%space, maxwell_system%gr, &
1876 vector_potential_analytical, unit_one, ierr)
1877 write(fname3, '(a)') 'vector_potential_mag'
1878 call io_function_output_vector(out_how , './' , trim(fname3), namespace, maxwell_system%space, maxwell_system%gr, &
1879 vector_potential_mag, unit_one, ierr)
1880
1881 safe_deallocate_a(magnetic_field)
1882 safe_deallocate_a(vector_potential_mag)
1883 safe_deallocate_a(vector_potential_analytical)
1884 safe_deallocate_a(delta)
1885 safe_deallocate_p(maxwell_system)
1886
1887 end subroutine test_vecpot_analytical
1888
1889 ! ---------------------------------------------------------
1890 subroutine multigrid_test_interpolation(mgrid, space)
1891 type(multigrid_t), intent(in) :: mgrid
1892 class(space_t), intent(in) :: space
1893
1894 real(real64), allocatable :: guess0(:), res0(:), guess1(:)
1895 type(mesh_t), pointer :: mesh0, mesh1
1896 real(real64) :: delta, xx(3,2), alpha, beta, rr
1897 integer :: nn, ip, ierr
1898
1900
1901 message(1) = 'Info: Testing the grid interpolation.'
1902 message(2) = ''
1903 call messages_info(2)
1904
1905 mesh0 => mgrid%level(0)%mesh
1906 mesh1 => mgrid%level(1)%mesh
1907
1908 safe_allocate(guess0(1:mesh0%np_part))
1909 safe_allocate(res0(1:mesh0%np_part))
1910 safe_allocate(guess1(1:mesh1%np_part))
1911
1912 alpha = m_four*mesh0%spacing(1)
1913 beta = m_one / (alpha**space%dim * sqrt(m_pi)**space%dim)
1914
1915 ! Set the centers of the Gaussians by hand
1916 xx(1, 1) = m_one
1917 xx(2, 1) = -m_half
1918 xx(3, 1) = m_two
1919 xx(1, 2) = -m_two
1920 xx(2, 2) = m_zero
1921 xx(3, 2) = -m_one
1922 xx = xx * alpha
1923
1924 ! Density as sum of Gaussians
1925 guess0 = m_zero
1926 do nn = 1, 2
1927 do ip = 1, mesh0%np
1928 call mesh_r(mesh0, ip, rr, origin = xx(:, nn))
1929 guess0(ip) = guess0(ip) + (-1)**nn * beta*exp(-(rr/alpha)**2)
1930 end do
1931 end do
1932
1933 call dio_function_output (io_function_fill_how('AxisX'), ".", "interpolation_target", global_namespace, &
1934 space, mesh0, guess0, unit_one, ierr)
1935 call dio_function_output (io_function_fill_how('AxisZ'), ".", "interpolation_target", global_namespace, &
1936 space, mesh0, guess0, unit_one, ierr)
1937 call dio_function_output (io_function_fill_how('PlaneZ'), ".", "interpolation_target", global_namespace, &
1938 space, mesh0, guess0, unit_one, ierr)
1939
1940 ! We start by testing the interpolation scheme. For this, we generate a function on the fine grid
1941 ! and we inject it on the coarse grid. Then we interpolate it back to the fine grid and we compare
1942 ! This allows for testing the quality of the interpolation scheme
1943
1944 ! move to level 1
1945 call dmultigrid_fine2coarse(mgrid%level(1)%tt, mgrid%level(0)%der, mesh1, guess0, guess1, injection)
1946 ! back to level 0
1947 call dmultigrid_coarse2fine(mgrid%level(1)%tt, mgrid%level(1)%der, mesh0, guess1, res0)
1948
1949 call dio_function_output (io_function_fill_how('AxisX'), ".", "interpolation_result", global_namespace, &
1950 space, mesh0, res0, unit_one, ierr)
1951 call dio_function_output (io_function_fill_how('AxisZ'), ".", "interpolation_result", global_namespace, &
1952 space, mesh0, res0, unit_one, ierr)
1953 call dio_function_output (io_function_fill_how('PlaneZ'), ".", "interpolation_result", global_namespace, &
1954 space, mesh0, res0, unit_one, ierr)
1955
1956 delta = dmf_nrm2(mesh0, guess0(1:mesh0%np)-res0(1:mesh0%np))
1957 write(message(1),'(a,e13.6)') 'Interpolation test (abs.) = ', delta
1958
1959 ! Now we test if the restriction+interpolation combination returns the original result or not
1960 ! This allows to test if restriction and interpolation are adjoint operators or not.
1961
1962 ! move to level 1
1963 call dmultigrid_fine2coarse(mgrid%level(1)%tt, mgrid%level(0)%der, mesh1, guess0, guess1, fullweight)
1964 ! back to level 0
1965 call dmultigrid_coarse2fine(mgrid%level(1)%tt, mgrid%level(1)%der, mesh0, guess1, res0)
1966
1967 call dio_function_output (io_function_fill_how('AxisX'), ".", "restriction_result", global_namespace, &
1968 space, mesh0, res0, unit_one, ierr)
1969 call dio_function_output (io_function_fill_how('AxisZ'), ".", "restriction_result", global_namespace, &
1970 space, mesh0, res0, unit_one, ierr)
1971 call dio_function_output (io_function_fill_how('PlaneZ'), ".", "restriction_result", global_namespace, &
1972 space, mesh0, res0, unit_one, ierr)
1973
1974 delta = dmf_nrm2(mesh0, guess0(1:mesh0%np)-res0(1:mesh0%np))
1975 write(message(2),'(a,e13.6)') 'Restriction test (abs.) = ', delta
1976 call messages_info(2)
1977
1978 safe_deallocate_a(guess0)
1979 safe_deallocate_a(res0)
1980 safe_deallocate_a(guess1)
1981
1984
1985
1987 subroutine test_current_density(namespace)
1988 type(namespace_t), intent(in) :: namespace
1989 type(electrons_t), pointer :: sys
1990
1991 type(current_t) :: current
1992 character(len=MAX_PATH_LEN) :: fname
1993 integer :: ierr, ip, idir
1994 integer(int64) :: out_how
1995 real(real64), allocatable :: current_para_ref(:,:), current_dia_ref(:,:), current_mag_ref(:,:), delta(:)
1996 real(real64) :: xx(3), rr, a0, mag_curr, sin_thet, sin_phi, cos_phi, vec_pot_slope
1997 complex(real64) :: alpha
1998
1999 sys => electrons_t(namespace, generate_epot=.false.)
2000 call sys%init_parallelization(mpi_world)
2001
2002 alpha = (0.0_real64, 0.5_real64)
2003 a0 = m_one
2004 vec_pot_slope = 0.4_real64 ! units of B field
2005
2006 call states_elec_allocate_wfns(sys%st, sys%gr, wfs_type = type_cmplx)
2007 call set_hydrogen_states(sys%st%group%psib(1, 1), sys%gr, namespace, alpha, a0)
2008
2009 call current_init(current, namespace)
2010
2011 call hamiltonian_elec_epot_generate(sys%hm, sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st)
2012
2013 safe_allocate(sys%hm%hm_base%vector_potential(1:3, 1:sys%gr%np))
2014
2015 sys%hm%hm_base%vector_potential = m_zero
2016 do ip = 1, sys%gr%np
2017 xx = sys%gr%x(ip, 1:3)
2018 sys%hm%hm_base%vector_potential(2, ip) = vec_pot_slope * xx(1) / p_c ! vector potential is here devided by c_0
2019 end do
2020
2021 call states_elec_allocate_current(sys%st, sys%space, sys%gr)
2022 call density_calc(sys%st, sys%gr, sys%st%rho)
2023 ! paramagnetic + diamagnetic current densities
2024 call current_calculate(current, namespace, sys%gr, sys%hm, sys%space, sys%st)
2025
2026 safe_allocate(current_para_ref(1:sys%gr%np,1:3))
2027 safe_allocate(current_dia_ref(1:sys%gr%np,1:3))
2028 safe_allocate(current_mag_ref(1:sys%gr%np,1:3))
2029 safe_allocate(delta(1:sys%gr%np))
2030
2031 ! analytic paramagnetic
2032 current_para_ref(:,:) = m_zero
2033 do ip = 1, sys%gr%np
2034 call mesh_r(sys%gr, ip, rr)
2035 xx = sys%gr%x(ip, 1:3)
2036 if (rr > r_small) then
2037 current_para_ref(ip,1:3) = - (psi_2s(rr, a0) * dr_psi_1s(rr, a0) - &
2038 psi_1s(rr, a0) * dr_psi_2s(rr, a0) ) * imag(alpha) / (1 + abs(alpha)**2) * xx(1:3) / rr
2039 end if
2040 end do
2041
2042 write(fname, '(a)') 'current_para'
2043 out_how = io_function_fill_how("PlaneZ")
2044 call io_function_output_vector(out_how , './' , trim(fname), namespace, sys%space, sys%gr, &
2045 sys%st%current_para(:,:,1), unit_one, ierr)
2046
2047 write(fname, '(a)') 'current_para-ref'
2048 out_how = io_function_fill_how("PlaneZ")
2049 call io_function_output_vector(out_how , './' , trim(fname), namespace, sys%space, sys%gr, &
2050 current_para_ref(:,:), unit_one, ierr)
2051
2052 do idir = 1, 3
2053 delta = m_zero
2054 delta(:) = current_para_ref(1:sys%gr%np, idir) - sys%st%current_para(1:sys%gr%np, idir, 1)
2055 write(message(idir),*) 'idir =',idir,', norm2(delta paramagnetic)',norm2(delta)
2056 end do
2057 call messages_info(3)
2058
2059 ! analytic diamagnetic
2060 current_dia_ref(:,:) = m_zero
2061 do ip = 1, sys%gr%np
2062 call mesh_r(sys%gr, ip, rr)
2063 current_dia_ref(ip,1:3) = - sys%hm%hm_base%vector_potential(1:3,ip) *&
2064 m_one / (1 + abs(alpha)**2) * abs(lc_hydrogen_state(rr, alpha, a0))**2
2065 end do
2066
2067 write(fname, '(a)') 'current_dia'
2068 out_how = io_function_fill_how("PlaneZ")
2069 call io_function_output_vector(out_how , './' , trim(fname), namespace, sys%space, sys%gr, &
2070 sys%st%current_dia(:,:,1), unit_one, ierr)
2071
2072 write(fname, '(a)') 'current_dia-ref'
2073 out_how = io_function_fill_how("PlaneZ")
2074 call io_function_output_vector(out_how , './' , trim(fname), namespace, sys%space, sys%gr, &
2075 current_dia_ref(:,:), unit_one, ierr)
2076
2077 do idir = 1, 3
2078 delta = m_zero
2079 delta(:) = current_dia_ref(1:sys%gr%np, idir) - sys%st%current_dia(1:sys%gr%np, idir, 1)
2080 write(message(idir),*) 'idir =',idir,', norm2(delta diamagnetic)',norm2(delta)
2081 end do
2082 call messages_info(3)
2083
2084 ! magnetization current
2085 call current_calculate_mag(sys%gr%der, sys%st)
2086
2087 ! analytic magnetization
2088 current_mag_ref(:,:) = m_zero
2089 do ip = 1, sys%gr%np
2090 call mesh_r(sys%gr, ip, rr)
2091 xx = sys%gr%x(ip, 1:3)
2092 if (norm2(xx(1:2)) > r_small .and. rr > r_small) then
2093 sin_thet = norm2(xx(1:2)) / rr
2094 sin_phi = xx(2) / norm2(xx(1:2))
2095 cos_phi = xx(1) / norm2(xx(1:2))
2096 mag_curr = m_two * (psi_1s(rr,a0)*dr_psi_1s(rr,a0) + real(alpha)*(psi_1s(rr,a0)*dr_psi_2s(rr,a0) &
2097 + dr_psi_1s(rr,a0)*psi_2s(rr,a0)) + abs(alpha)**2 * psi_2s(rr,a0)*dr_psi_2s(rr,a0))
2098 ! minus signs are reversed because of the charge of the electron
2099 current_mag_ref(ip,1) = m_half * mag_curr * sin_thet * sin_phi / (1+abs(alpha)**2)
2100 current_mag_ref(ip,2) = -m_half * mag_curr * sin_thet * cos_phi / (1+abs(alpha)**2)
2101 endif
2102 end do
2103
2104 write(fname, '(a)') 'current_mag'
2105 out_how = io_function_fill_how("PlaneZ")
2106 call io_function_output_vector(out_how , './' , trim(fname), namespace, sys%space, sys%gr, &
2107 sys%st%current_mag(:,:,1), unit_one, ierr)
2108
2109 write(fname, '(a)') 'current_mag-ref'
2110 out_how = io_function_fill_how("PlaneZ")
2111 call io_function_output_vector(out_how , './' , trim(fname), namespace, sys%space, sys%gr, &
2112 current_mag_ref(:,:), unit_one, ierr)
2113
2114 do idir = 1, 3
2115 delta = m_zero
2116 delta(:) = current_mag_ref(1:sys%gr%np, idir) - sys%st%current_mag(1:sys%gr%np, idir, 1)
2117 write(message(idir),*) 'idir =',idir,', norm2(delta magnetization)',norm2(delta)
2118 end do
2119 call messages_info(3)
2120
2121 safe_deallocate_a(current_para_ref)
2122 safe_deallocate_a(current_dia_ref)
2123 safe_deallocate_a(current_mag_ref)
2124 safe_deallocate_a(delta)
2125 safe_deallocate_p(sys)
2126
2127 end subroutine test_current_density
2128
2129
2130 subroutine set_hydrogen_states(psib, mesh, namespace, alpha, a0)
2131 class(batch_t), intent(inout) :: psib
2132 class(mesh_t), intent(in) :: mesh
2133 type(namespace_t), intent(in) :: namespace
2134 complex(real64), intent(in) :: alpha
2135 real(real64), intent(in) :: a0
2136
2137 complex(real64), allocatable :: zff(:)
2138 integer :: ip
2139 real(real64) :: rr
2140
2141 push_sub(set_hydrogen_states)
2142
2143 ! psi = phi_{1s} + \alpha \phi_{2s}
2144 safe_allocate(zff(1:mesh%np))
2145 if (type_is_complex(psib%type())) then
2146 do ip = 1, mesh%np
2147 call mesh_r(mesh, ip, rr)
2148 zff(ip) = m_one / sqrt(1 + abs(alpha)**2) * lc_hydrogen_state(rr, alpha, a0)
2149 call batch_set_state(psib, 1, mesh%np, zff)
2150 end do
2151 safe_deallocate_a(zff)
2152 else
2153 write(message(1),*) "States should be complex for the linear combination of hydrogenic states to work"
2154 call messages_info(1, namespace=namespace)
2155 end if
2156
2157 pop_sub(set_hydrogen_states)
2158 end subroutine set_hydrogen_states
2159
2160 complex(real64) function lc_hydrogen_state(rr, alpha, a0)
2161 real(real64), intent(in) :: a0, rr
2162 complex(real64), intent(in) :: alpha
2163
2164 lc_hydrogen_state = psi_1s(rr, a0) + alpha * psi_2s(rr, a0)
2165 end function lc_hydrogen_state
2166
2167 ! phi_{1s} = R_{10}(r)Y_{00}(\theta, \phi) = 1/(2\sqrt(\pi)) * 2 a0^{-3/2} * \exp(-r/a0)
2168 real(real64) function psi_1s(rr, a0)
2169 real(real64), intent(in) :: a0, rr
2170
2171 psi_1s = a0**(-m_three/m_two) * exp(-rr / a0) / sqrt(m_pi)
2172 end function psi_1s
2173
2174 ! phi_{2s} = R_{20}(r)Y_{00}(\theta, \phi) = 1/(2\sqrt(\pi)) * \sqrt(2) a0^{-3/2} (1 - r/(2a0)) \exp(-r/(2a0))
2175 real(real64) function psi_2s(rr, a0)
2176 real(real64), intent(in) :: a0, rr
2177
2178 psi_2s = sqrt(m_two) * a0**(-m_three/m_two) * (m_one - rr/(m_two * a0)) * exp(-rr/(m_two * a0)) &
2179 / (m_two * sqrt(m_pi))
2180 end function psi_2s
2181
2182 real(real64) function dr_psi_1s(rr, a0)
2183 real(real64), intent(in) :: a0, rr
2184
2185 dr_psi_1s = -(m_one / a0) * psi_1s(rr, a0)
2186 end function dr_psi_1s
2187
2188 real(real64) function dr_psi_2s(rr, a0)
2189 real(real64), intent(in) :: a0, rr
2190
2191 dr_psi_2s = -(m_half / a0) * psi_2s(rr, a0) - &
2192 a0**(-m_five/m_two) * exp(-rr/(m_two * a0)) / (sqrt(m_two) * m_two * sqrt(m_pi))
2193 end function dr_psi_2s
2194
2195end module test_oct_m
2196
2197!! Local Variables:
2198!! mode: f90
2199!! coding: utf-8
2200!! End:
batchified version of the BLAS axpy routine:
Definition: batch_ops.F90:154
scale a batch by a constant or vector
Definition: batch_ops.F90:162
There are several ways how to call batch_set_state and batch_get_state:
Definition: batch_ops.F90:201
constant times a vector plus a vector
Definition: lalg_basic.F90:170
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:185
scales a vector by a constant
Definition: lalg_basic.F90:156
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:180
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
Definition: batch_ops.F90:242
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:122
Module, implementing a factory for boxes.
This module handles the calculation mode.
integer, parameter, public p_strategy_kpoints
parallelization in k-points
integer, parameter, public p_strategy_domains
parallelization in domains
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
subroutine, public cgal_polyhedron_init(cgal_poly, fname, verbose)
logical function, public cgal_polyhedron_point_inside(cgal_poly, xq, yq, zq)
subroutine, public cgal_polyhedron_end(cgal_poly)
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:608
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_test(this, namespace, repetitions, min_blocksize, max_blocksize)
unit test for derivatives
subroutine, public dderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
subroutine, public zderivatives_test(this, namespace, repetitions, min_blocksize, max_blocksize)
unit test for derivatives
integer, parameter, public spin_polarized
subroutine, public exponential_init(te, namespace, full_batch)
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_four
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public zhamiltonian_elec_apply_batch(hm, namespace, mesh, psib, hpsib, terms, set_bc)
subroutine, public hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
subroutine, public dhamiltonian_elec_apply_batch(hm, namespace, mesh, psib, hpsib, terms, set_bc)
The Helmholtz decomposition is intended to contain "only mathematical" functions and procedures to co...
Test suit for the Helmholtz decomposition module.
subroutine, public gaussian_test(helmholtz, sys_grid, namespace, space)
subroutine, public hertzian_dipole_test(helmholtz, sys_grid, namespace, space)
This module implements a simple hash table for non-negative integer keys and integer values.
Definition: iihash.F90:125
subroutine, public iihash_end(h)
Free a hash table.
Definition: iihash.F90:184
subroutine, public iihash_insert(h, key, val)
Insert a (key, val) pair into the hash table h.
Definition: iihash.F90:206
integer function, public iihash_lookup(h, key, found)
Look up a value in the hash table h. If found is present, it indicates if key could be found in the t...
Definition: iihash.F90:231
subroutine, public iihash_init(h)
Initialize a hash table h.
Definition: iihash.F90:161
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
integer(int64) function, public io_function_fill_how(where)
Use this function to quickly plot functions for debugging purposes: call dio_function_output(io_funct...
Definition: io.F90:114
subroutine, public ion_interaction_test(space, latt, atom, natoms, pos, lsize, namespace, mc)
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:116
subroutine, public dmesh_batch_dotp_matrix(mesh, aa, bb, dot, reduce)
Calculate the overlap matrix of two batches.
Definition: mesh_batch.F90:271
subroutine, public dmesh_batch_dotp_self(mesh, aa, dot, reduce)
calculate the overlap matrix of a batch with itself
Definition: mesh_batch.F90:517
subroutine, public mesh_batch_nrm2(mesh, aa, nrm2, reduce)
Calculate the norms (norm2, not the square!) of a batch of mesh functions.
Definition: mesh_batch.F90:177
subroutine, public zmesh_batch_dotp_vector(mesh, aa, bb, dot, reduce, cproduct)
calculate the vector of dot-products of mesh functions between two batches
subroutine, public zmesh_batch_dotp_matrix(mesh, aa, bb, dot, reduce)
Calculate the overlap matrix of two batches.
subroutine, public zmesh_batch_normalize(mesh, psib, norm)
Normalize a batch.
subroutine, public zmesh_batch_dotp_self(mesh, aa, dot, reduce)
calculate the overlap matrix of a batch with itself
subroutine, public dmesh_batch_normalize(mesh, psib, norm)
Normalize a batch.
subroutine, public dmesh_batch_dotp_vector(mesh, aa, bb, dot, reduce, cproduct)
calculate the vector of dot-products of mesh functions between two batches
Definition: mesh_batch.F90:632
This module defines various routines, operating on mesh functions.
class(mesh_t), pointer, public mesh_aux
Globally-scoped pointer to the mesh instance.
real(real64) function, public dmf_dotp_aux(f1, f2)
dot product between two vectors (mesh functions)
subroutine, public mesh_init_mesh_aux(mesh)
Initialise a pointer to the grid/mesh, that is globally exposed, such that low level mesh operations ...
subroutine, public zmesh_interpolation_test(mesh)
subroutine, public dmesh_interpolation_test(mesh)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:336
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:930
character(len=512), private msg
Definition: messages.F90:165
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
subroutine, public messages_new_line()
Definition: messages.F90:1146
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_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:624
This modules takes care of testing optimizers using standard test functions.
subroutine, public test_optimizers(namespace)
Unit tests for different optimizers.
This module implements unit tests for the mixing methods.
Definition: mix_tests.F90:115
subroutine, public mix_tests_run()
Definition: mix_tests.F90:137
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
subroutine, public test_mpiwrappers
Definition: mpi_test.F90:128
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
subroutine, public multigrid_end(mgrid)
Definition: multigrid.F90:502
subroutine, public multigrid_init(mgrid, namespace, space, mesh, der, stencil, mc, nlevels)
Definition: multigrid.F90:185
type(namespace_t), public global_namespace
Definition: namespace.F90:132
subroutine, public orbitalbasis_end(this)
subroutine, public zorbitalbasis_build(this, namespace, ions, mesh, kpt, ndim, skip_s_orb, use_all_orb, verbose)
This routine is an interface for constructing the orbital basis.
subroutine, public dorbitalbasis_build(this, namespace, ions, mesh, kpt, ndim, skip_s_orb, use_all_orb, verbose)
This routine is an interface for constructing the orbital basis.
subroutine, public orbitalbasis_init(this, namespace, periodic_dim)
subroutine, public dorbitalset_get_coeff_batch(os, ndim, psib, dot)
Definition: orbitalset.F90:595
subroutine, public dorbitalset_add_to_batch(os, ndim, psib, weight)
Definition: orbitalset.F90:739
subroutine, public zorbitalset_add_to_batch(os, ndim, psib, weight)
subroutine, public orbitalset_update_phase(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
Build the phase correction to the global phase in case the orbital crosses the border of the simulato...
Definition: orbitalset.F90:279
subroutine, public zorbitalset_get_coeff_batch(os, ndim, psib, dot)
subroutine, public poisson_test(this, space, mesh, latt, namespace, repetitions)
This routine checks the Hartree solver selected in the input file by calculating numerically and anal...
Definition: poisson.F90:1104
subroutine, public preconditioner_end(this)
subroutine, public preconditioner_init(this, namespace, gr, mc, space)
subroutine, public zproject_psi_batch(mesh, bnd, pj, npj, dim, psib, ppsib)
To optimize the application of the non-local operator in parallel, the projectors are applied in step...
Definition: projector.F90:1258
Implementation details for regridding.
Definition: regridding.F90:170
This module implements a simple hash table for string valued keys and integer values using the C++ ST...
Definition: sihash.F90:118
subroutine, public sihash_insert(h, key, val)
Insert a (key, val) pair into the hash table h.
Definition: sihash.F90:203
subroutine, public sihash_init(h)
Initialize a hash table h with size entries. Since we use separate chaining, the number of entries in...
Definition: sihash.F90:162
integer function, public sihash_lookup(h, key, found)
Look up a value in the hash table h. If found is present, it indicates if key could be found in the t...
Definition: sihash.F90:229
subroutine, public sihash_end(h)
Free a hash table.
Definition: sihash.F90:183
This module is intended to contain "only mathematical" functions and procedures.
Definition: solvers.F90:115
This module implements a simple hash table for string valued keys and integer values using the C++ ST...
Definition: sphash.F90:118
subroutine, public sphash_init(h)
Initialize a hash table h with size entries. Since we use separate chaining, the number of entries in...
Definition: sphash.F90:222
subroutine, public sphash_insert(h, key, val, clone)
Insert a (key, val) pair into the hash table h. If clone=.true., the object will be copied.
Definition: sphash.F90:291
subroutine, public sphash_end(h)
Free a hash table.
Definition: sphash.F90:248
class(*) function, pointer, public sphash_lookup(h, key, found)
Look up a value in the hash table h. If found is present, it indicates if key could be found in the t...
Definition: sphash.F90:323
pure logical function, public states_are_real(st)
subroutine, public zstates_elec_calc_orth_test(st, namespace, mesh, kpoints)
subroutine, public dstates_elec_calc_orth_test(st, namespace, mesh, kpoints)
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public dsubspace_diag(this, namespace, mesh, st, hm, ik, eigenval, diff, nonortho)
Diagonalises the Hamiltonian in the subspace defined by the states.
Definition: subspace.F90:304
subroutine, public zsubspace_diag(this, namespace, mesh, st, hm, ik, eigenval, diff, nonortho)
Diagonalises the Hamiltonian in the subspace defined by the states.
Definition: subspace.F90:915
This module implements a unit-test like runmode for Octopus.
Definition: test.F90:116
real(real64) function psi_2s(rr, a0)
Definition: test.F90:2269
subroutine test_sihash()
Definition: test.F90:1627
subroutine test_cgal()
Definition: test.F90:1449
subroutine test_helmholtz_decomposition(namespace)
Definition: test.F90:434
subroutine test_hartree(param, namespace)
Definition: test.F90:415
subroutine test_derivatives(param, namespace)
Definition: test.F90:1215
subroutine test_subspace_diagonalization(param, namespace)
Definition: test.F90:921
subroutine, public test_run(namespace)
Definition: test.F90:206
subroutine test_density_calc(param, namespace)
Definition: test.F90:797
subroutine test_iihash()
Definition: test.F90:1601
real(real64) function psi_1s(rr, a0)
Definition: test.F90:2262
subroutine test_current_density(namespace)
Here we test the different contributions to the total electronic current density.
Definition: test.F90:2081
subroutine test_batch_set_gaussian(psib, mesh)
Definition: test.F90:1528
subroutine test_regridding(namespace)
Definition: test.F90:1766
subroutine test_sphash(namespace)
Definition: test.F90:1673
subroutine test_hamiltonian(param, namespace)
Definition: test.F90:715
subroutine test_dense_eigensolver()
Definition: test.F90:1469
subroutine test_interpolation(param, namespace)
Definition: test.F90:1287
subroutine multigrid_test_interpolation(mgrid, space)
Definition: test.F90:1984
subroutine test_ion_interaction(namespace)
Definition: test.F90:1325
subroutine test_boundaries(param, namespace)
Definition: test.F90:835
subroutine test_orthogonalization(param, namespace)
Definition: test.F90:1245
real(real64) function dr_psi_2s(rr, a0)
Definition: test.F90:2282
subroutine test_batch_ops(param, namespace)
Definition: test.F90:974
complex(real64) function lc_hydrogen_state(rr, alpha, a0)
Definition: test.F90:2254
subroutine test_projector(param, namespace)
Definition: test.F90:555
subroutine test_dft_u(param, namespace)
Definition: test.F90:607
subroutine test_prints_info_batch(st, gr, psib, string)
Definition: test.F90:1345
real(real64) function dr_psi_1s(rr, a0)
Definition: test.F90:2276
subroutine test_linear_solver(namespace)
Definition: test.F90:465
subroutine test_grid_interpolation()
Definition: test.F90:1580
subroutine set_hydrogen_states(psib, mesh, namespace, alpha, a0)
Definition: test.F90:2224
subroutine test_exponential(param, namespace)
Definition: test.F90:872
subroutine test_vecpot_analytical(namespace)
Here, analytical formulation for vector potential and B field are used. Ref: Sangita Sen and Erik I....
Definition: test.F90:1884
subroutine test_clock()
Definition: test.F90:1390
type(type_t), public type_cmplx
Definition: types.F90:134
logical pure function, public type_is_complex(this)
Definition: types.F90:180
This module defines the unit system, used for input and output.
type(unit_t), public unit_one
some special units required for particular quantities
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Definition: v_ks.F90:732
Class defining batches of mesh functions.
Definition: batch.F90:159
class representing derivatives
Class describing the electron system.
Definition: electrons.F90:218
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
This class implements the iteration counter used by the multisystem algorithms. As any iteration coun...
Describes mesh distribution to nodes.
Definition: mesh.F90:186
contains the information of the meshes and provides the transfer functions
Definition: regridding.F90:197
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:138
int true(void)
subroutine write_clock(operation)
Definition: test.F90:1430
subroutine laplacian_op(x, hx)
Computes .
Definition: test.F90:534
subroutine set_der_aux(der)
Definition: test.F90:525
subroutine write_condition_result(condition, result)
Definition: test.F90:1437
real(real64) function values(xx)
Definition: test.F90:1867