(git:ed6f26b)
Loading...
Searching...
No Matches
library_tests.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Performance tests for basic tasks like matrix multiplies, copy, fft.
10!> \par History
11!> 30-Nov-2000 (JGH) added input
12!> 02-Jan-2001 (JGH) Parallel FFT
13!> 28-Feb-2002 (JGH) Clebsch-Gordon Coefficients
14!> 06-Jun-2003 (JGH) Real space grid test
15!> Eigensolver test (29.08.05,MK)
16!> \author JGH 6-NOV-2000
17! **************************************************************************************************
19
20 USE ai_coulomb_test, ONLY: eri_test
21 USE cell_methods, ONLY: cell_create,&
23 USE cell_types, ONLY: cell_release,&
29 USE cp_dbcsr_api, ONLY: dbcsr_reset_randmat_seed,&
30 dbcsr_run_tests
32 USE cp_files, ONLY: close_file,&
35 USE cp_fm_diag, ONLY: cp_fm_syevd,&
41 USE cp_fm_types, ONLY: cp_fm_create,&
53 USE dbm_tests, ONLY: dbm_run_tests
54 USE fft_tools, ONLY: bwfft,&
56 fwfft,&
57 fft3d,&
73 USE kinds, ONLY: dp
74 USE machine, ONLY: m_flush,&
81 USE parallel_rng_types, ONLY: uniform,&
83 USE pw_grid_types, ONLY: fullspace,&
84 halfspace,&
86 USE pw_grids, ONLY: pw_grid_create,&
88 USE pw_methods, ONLY: pw_transfer,&
90 USE pw_types, ONLY: pw_c1d_gs_type,&
93 USE realspace_grid_types, ONLY: &
98#include "./base/base_uses.f90"
99
100 IMPLICIT NONE
101
102 PRIVATE
103 PUBLIC :: lib_test
104
105 INTEGER :: runtest(100)
106 REAL(KIND=dp) :: max_memory
107 REAL(KIND=dp), PARAMETER :: threshold = 1.0e-8_dp
108 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'library_tests'
109
110CONTAINS
111
112! **************************************************************************************************
113!> \brief Master routine for tests
114!> \param root_section ...
115!> \param para_env ...
116!> \param globenv ...
117!> \par History
118!> none
119!> \author JGH 6-NOV-2000
120! **************************************************************************************************
121 SUBROUTINE lib_test(root_section, para_env, globenv)
122
123 TYPE(section_vals_type), POINTER :: root_section
124 TYPE(mp_para_env_type), POINTER :: para_env
125 TYPE(global_environment_type), POINTER :: globenv
126
127 CHARACTER(LEN=*), PARAMETER :: routinen = 'lib_test'
128
129 INTEGER :: handle, iw
130 LOGICAL :: explicit
131 TYPE(cp_logger_type), POINTER :: logger
132 TYPE(section_vals_type), POINTER :: cp_dbcsr_test_section, cp_fm_gemm_test_section, &
133 dbm_test_section, eigensolver_section, eri_mme_test_section, pw_transfer_section, &
134 rs_pw_transfer_section, shg_integrals_test_section
135
136 CALL timeset(routinen, handle)
137
138 logger => cp_get_default_logger()
139 iw = cp_print_key_unit_nr(logger, root_section, "TEST%PROGRAM_RUN_INFO", extension=".log")
140
141 IF (iw > 0) THEN
142 WRITE (iw, '(T2,79("*"))')
143 WRITE (iw, '(A,T31,A,T80,A)') ' *', ' PERFORMANCE TESTS ', '*'
144 WRITE (iw, '(T2,79("*"))')
145 END IF
146 !
147 CALL test_input(root_section, para_env)
148 !
149 IF (runtest(1) /= 0) CALL copy_test(para_env, iw)
150 !
151 IF (runtest(2) /= 0) CALL matmul_test(para_env, test_matmul=.true., test_dgemm=.false., iw=iw)
152 IF (runtest(5) /= 0) CALL matmul_test(para_env, test_matmul=.false., test_dgemm=.true., iw=iw)
153 !
154 IF (runtest(3) /= 0) CALL fft_test(para_env, iw, globenv%fftw_plan_type, &
155 globenv%fftw_wisdom_file_name)
156 !
157 IF (runtest(4) /= 0) CALL eri_test(iw)
158 !
159 IF (runtest(6) /= 0) CALL clebsch_gordon_test()
160 !
161 ! runtest 7 has been deleted and can be recycled
162 !
163 IF (runtest(8) /= 0) CALL mpi_perf_test(para_env, runtest(8), iw)
164 !
165 IF (runtest(10) /= 0) CALL validate_exp_minimax(runtest(10), iw)
166 !
167 IF (runtest(11) /= 0) CALL test_least_square_ft(runtest(11), iw)
168 !
169
170 rs_pw_transfer_section => section_vals_get_subs_vals(root_section, "TEST%RS_PW_TRANSFER")
171 CALL section_vals_get(rs_pw_transfer_section, explicit=explicit)
172 IF (explicit) THEN
173 CALL rs_pw_transfer_test(para_env, iw, globenv, rs_pw_transfer_section)
174 END IF
175
176 pw_transfer_section => section_vals_get_subs_vals(root_section, "TEST%PW_TRANSFER")
177 CALL section_vals_get(pw_transfer_section, explicit=explicit)
178 IF (explicit) THEN
179 CALL pw_fft_test(para_env, iw, globenv, pw_transfer_section)
180 END IF
181
182 cp_fm_gemm_test_section => section_vals_get_subs_vals(root_section, "TEST%CP_FM_GEMM")
183 CALL section_vals_get(cp_fm_gemm_test_section, explicit=explicit)
184 IF (explicit) THEN
185 CALL cp_fm_gemm_test(para_env, iw, cp_fm_gemm_test_section)
186 END IF
187
188 eigensolver_section => section_vals_get_subs_vals(root_section, "TEST%EIGENSOLVER")
189 CALL section_vals_get(eigensolver_section, explicit=explicit)
190 IF (explicit) THEN
191 CALL eigensolver_test(para_env, iw, eigensolver_section)
192 END IF
193
194 eri_mme_test_section => section_vals_get_subs_vals(root_section, "TEST%ERI_MME_TEST")
195 CALL section_vals_get(eri_mme_test_section, explicit=explicit)
196 IF (explicit) THEN
197 CALL cp_eri_mme_perf_acc_test(para_env, iw, eri_mme_test_section)
198 END IF
199
200 shg_integrals_test_section => section_vals_get_subs_vals(root_section, "TEST%SHG_INTEGRALS_TEST")
201 CALL section_vals_get(shg_integrals_test_section, explicit=explicit)
202 IF (explicit) THEN
203 CALL shg_integrals_perf_acc_test(iw, shg_integrals_test_section)
204 END IF
205
206 ! DBCSR tests
207 cp_dbcsr_test_section => section_vals_get_subs_vals(root_section, "TEST%CP_DBCSR")
208 CALL section_vals_get(cp_dbcsr_test_section, explicit=explicit)
209 IF (explicit) THEN
210 CALL cp_dbcsr_tests(para_env, iw, cp_dbcsr_test_section)
211 END IF
212
213 ! DBM tests
214 dbm_test_section => section_vals_get_subs_vals(root_section, "TEST%DBM")
215 CALL section_vals_get(dbm_test_section, explicit=explicit)
216 IF (explicit) THEN
217 CALL run_dbm_tests(para_env, iw, dbm_test_section)
218 END IF
219
220 CALL cp_print_key_finished_output(iw, logger, root_section, "TEST%PROGRAM_RUN_INFO")
221
222 CALL timestop(handle)
223
224 END SUBROUTINE lib_test
225
226! **************************************************************************************************
227!> \brief Reads input section &TEST ... &END
228!> \param root_section ...
229!> \param para_env ...
230!> \author JGH 30-NOV-2000
231!> \note
232!> I---------------------------------------------------------------------------I
233!> I SECTION: &TEST ... &END I
234!> I I
235!> I MEMORY max_memory I
236!> I COPY n I
237!> I MATMUL n I
238!> I FFT n I
239!> I ERI n I
240!> I PW_FFT n I
241!> I Clebsch-Gordon n I
242!> I RS_GRIDS n I
243!> I MPI n I
244!> I RNG n -> Parallel random number generator I
245!> I---------------------------------------------------------------------------I
246! **************************************************************************************************
247 SUBROUTINE test_input(root_section, para_env)
248 TYPE(section_vals_type), POINTER :: root_section
249 TYPE(mp_para_env_type), POINTER :: para_env
250
251 TYPE(section_vals_type), POINTER :: test_section
252
253!
254!..defaults
255! using this style is not recommended, introduce sections instead (see e.g. cp_fm_gemm)
256
257 runtest = 0
258 test_section => section_vals_get_subs_vals(root_section, "TEST")
259 CALL section_vals_val_get(test_section, "MEMORY", r_val=max_memory)
260 CALL section_vals_val_get(test_section, 'COPY', i_val=runtest(1))
261 CALL section_vals_val_get(test_section, 'MATMUL', i_val=runtest(2))
262 CALL section_vals_val_get(test_section, 'DGEMM', i_val=runtest(5))
263 CALL section_vals_val_get(test_section, 'FFT', i_val=runtest(3))
264 CALL section_vals_val_get(test_section, 'ERI', i_val=runtest(4))
265 CALL section_vals_val_get(test_section, 'CLEBSCH_GORDON', i_val=runtest(6))
266 CALL section_vals_val_get(test_section, 'MPI', i_val=runtest(8))
267 CALL section_vals_val_get(test_section, 'MINIMAX', i_val=runtest(10))
268 CALL section_vals_val_get(test_section, 'LEAST_SQ_FT', i_val=runtest(11))
269
270 CALL para_env%sync()
271 END SUBROUTINE test_input
272
273! **************************************************************************************************
274!> \brief Tests the performance to copy two vectors.
275!> \param para_env ...
276!> \param iw ...
277!> \par History
278!> none
279!> \author JGH 6-NOV-2000
280!> \note
281!> The results of these tests allow to determine the size of the cache
282!> of the CPU. This can be used to optimize the performance of the
283!> FFTSG library.
284! **************************************************************************************************
285 SUBROUTINE copy_test(para_env, iw)
286 TYPE(mp_para_env_type), POINTER :: para_env
287 INTEGER :: iw
288
289 INTEGER :: i, ierr, j, len, ntim, siz
290 REAL(kind=dp) :: perf, t, tend, tstart
291 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: ca, cb
292
293! test for copy --> Cache size
294
295 siz = abs(runtest(1))
296 IF (para_env%is_source()) WRITE (iw, '(//,A,/)') " Test of copy ( F95 ) "
297 DO i = 6, 24
298 len = 2**i
299 IF (8.0_dp*real(len, kind=dp) > max_memory*0.5_dp) EXIT
300 ALLOCATE (ca(len), stat=ierr)
301 IF (ierr /= 0) EXIT
302 ALLOCATE (cb(len), stat=ierr)
303 IF (ierr /= 0) EXIT
304
305 CALL random_number(ca)
306 ntim = nint(1.e7_dp/real(len, kind=dp))
307 ntim = max(ntim, 1)
308 ntim = min(ntim, siz*10000)
309
310 tstart = m_walltime()
311 DO j = 1, ntim
312 cb(:) = ca(:)
313 ca(1) = real(j, kind=dp)
314 END DO
315 tend = m_walltime()
316 t = tend - tstart + threshold
317 IF (t > 0.0_dp) THEN
318 perf = real(ntim, kind=dp)*real(len, kind=dp)*1.e-6_dp/t
319 ELSE
320 perf = 0.0_dp
321 END IF
322
323 IF (para_env%is_source()) THEN
324 WRITE (iw, '(A,i2,i10,A,T59,F14.4,A)') " Copy test: Size = 2^", i, &
325 len/1024, " Kwords", perf, " Mcopy/s"
326 END IF
327
328 DEALLOCATE (ca)
329 DEALLOCATE (cb)
330 END DO
331 CALL para_env%sync()
332 END SUBROUTINE copy_test
333
334! **************************************************************************************************
335!> \brief Tests the performance of different kinds of matrix matrix multiply
336!> kernels for the BLAS and F95 intrinsic matmul.
337!> \param para_env ...
338!> \param test_matmul ...
339!> \param test_dgemm ...
340!> \param iw ...
341!> \par History
342!> none
343!> \author JGH 6-NOV-2000
344! **************************************************************************************************
345 SUBROUTINE matmul_test(para_env, test_matmul, test_dgemm, iw)
346 TYPE(mp_para_env_type), POINTER :: para_env
347 LOGICAL :: test_matmul, test_dgemm
348 INTEGER :: iw
349
350 INTEGER :: i, ierr, j, len, ntim, siz
351 REAL(kind=dp) :: perf, t, tend, tstart, xdum
352 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: ma, mb, mc
353
354! test for matrix multpies
355
356 IF (test_matmul) THEN
357 siz = abs(runtest(2))
358 IF (para_env%is_source()) WRITE (iw, '(//,A,/)') " Test of matmul ( F95 ) "
359 DO i = 5, siz, 2
360 len = 2**i + 1
361 IF (8.0_dp*real(len*len, kind=dp) > max_memory*0.3_dp) EXIT
362 ALLOCATE (ma(len, len), stat=ierr)
363 IF (ierr /= 0) EXIT
364 ALLOCATE (mb(len, len), stat=ierr)
365 IF (ierr /= 0) EXIT
366 ALLOCATE (mc(len, len), stat=ierr)
367 IF (ierr /= 0) EXIT
368 mc = 0.0_dp
369
370 CALL random_number(xdum)
371 ma = xdum
372 CALL random_number(xdum)
373 mb = xdum
374 ntim = nint(1.e8_dp/(2.0_dp*real(len, kind=dp)**3))
375 ntim = max(ntim, 1)
376 ntim = min(ntim, siz*200)
377 tstart = m_walltime()
378 DO j = 1, ntim
379 mc(:, :) = matmul(ma, mb)
380 ma(1, 1) = real(j, kind=dp)
381 END DO
382 tend = m_walltime()
383 t = tend - tstart + threshold
384 perf = real(ntim, kind=dp)*2.0_dp*real(len, kind=dp)**3*1.e-6_dp/t
385 IF (para_env%is_source()) THEN
386 WRITE (iw, '(A,i6,T59,F14.4,A)') &
387 " Matrix multiply test: c = a * b Size = ", len, perf, " Mflop/s"
388 END IF
389 tstart = m_walltime()
390 DO j = 1, ntim
391 mc(:, :) = mc + matmul(ma, mb)
392 ma(1, 1) = real(j, kind=dp)
393 END DO
394 tend = m_walltime()
395 t = tend - tstart + threshold
396 IF (t > 0.0_dp) THEN
397 perf = real(ntim, kind=dp)*2.0_dp*real(len, kind=dp)**3*1.e-6_dp/t
398 ELSE
399 perf = 0.0_dp
400 END IF
401
402 IF (para_env%is_source()) THEN
403 WRITE (iw, '(A,i6,T59,F14.4,A)') &
404 " Matrix multiply test: a = a * b Size = ", len, perf, " Mflop/s"
405 END IF
406
407 tstart = m_walltime()
408 DO j = 1, ntim
409 mc(:, :) = mc + matmul(ma, transpose(mb))
410 ma(1, 1) = real(j, kind=dp)
411 END DO
412 tend = m_walltime()
413 t = tend - tstart + threshold
414 IF (t > 0.0_dp) THEN
415 perf = real(ntim, kind=dp)*2.0_dp*real(len, kind=dp)**3*1.e-6_dp/t
416 ELSE
417 perf = 0.0_dp
418 END IF
419
420 IF (para_env%is_source()) THEN
421 WRITE (iw, '(A,i6,T59,F14.4,A)') &
422 " Matrix multiply test: c = a * b(T) Size = ", len, perf, " Mflop/s"
423 END IF
424
425 tstart = m_walltime()
426 DO j = 1, ntim
427 mc(:, :) = mc + matmul(transpose(ma), mb)
428 ma(1, 1) = real(j, kind=dp)
429 END DO
430 tend = m_walltime()
431 t = tend - tstart + threshold
432 IF (t > 0.0_dp) THEN
433 perf = real(ntim, kind=dp)*2.0_dp*real(len, kind=dp)**3*1.e-6_dp/t
434 ELSE
435 perf = 0.0_dp
436 END IF
437
438 IF (para_env%is_source()) THEN
439 WRITE (iw, '(A,i6,T59,F14.4,A)') &
440 " Matrix multiply test: c = a(T) * b Size = ", len, perf, " Mflop/s"
441 END IF
442
443 DEALLOCATE (ma)
444 DEALLOCATE (mb)
445 DEALLOCATE (mc)
446 END DO
447 END IF
448
449 ! test for matrix multpies
450 IF (test_dgemm) THEN
451 siz = abs(runtest(5))
452 IF (para_env%is_source()) WRITE (iw, '(//,A,/)') " Test of matmul ( BLAS ) "
453 DO i = 5, siz, 2
454 len = 2**i + 1
455 IF (8.0_dp*real(len*len, kind=dp) > max_memory*0.3_dp) EXIT
456 ALLOCATE (ma(len, len), stat=ierr)
457 IF (ierr /= 0) EXIT
458 ALLOCATE (mb(len, len), stat=ierr)
459 IF (ierr /= 0) EXIT
460 ALLOCATE (mc(len, len), stat=ierr)
461 IF (ierr /= 0) EXIT
462 mc = 0.0_dp
463
464 CALL random_number(xdum)
465 ma = xdum
466 CALL random_number(xdum)
467 mb = xdum
468 ntim = nint(1.e8_dp/(2.0_dp*real(len, kind=dp)**3))
469 ntim = max(ntim, 1)
470 ntim = min(ntim, 1000)
471
472 tstart = m_walltime()
473 DO j = 1, ntim
474 CALL dgemm("N", "N", len, len, len, 1.0_dp, ma, len, mb, len, 1.0_dp, mc, len)
475 END DO
476 tend = m_walltime()
477 t = tend - tstart + threshold
478 IF (t > 0.0_dp) THEN
479 perf = real(ntim, kind=dp)*2.0_dp*real(len, kind=dp)**3*1.e-6_dp/t
480 ELSE
481 perf = 0.0_dp
482 END IF
483
484 IF (para_env%is_source()) THEN
485 WRITE (iw, '(A,i6,T59,F14.4,A)') &
486 " Matrix multiply test: c = a * b Size = ", len, perf, " Mflop/s"
487 END IF
488
489 tstart = m_walltime()
490 DO j = 1, ntim
491 CALL dgemm("N", "N", len, len, len, 1.0_dp, ma, len, mb, len, 1.0_dp, mc, len)
492 END DO
493 tend = m_walltime()
494 t = tend - tstart + threshold
495 IF (t > 0.0_dp) THEN
496 perf = real(ntim, kind=dp)*2.0_dp*real(len, kind=dp)**3*1.e-6_dp/t
497 ELSE
498 perf = 0.0_dp
499 END IF
500
501 IF (para_env%is_source()) THEN
502 WRITE (iw, '(A,i6,T59,F14.4,A)') &
503 " Matrix multiply test: a = a * b Size = ", len, perf, " Mflop/s"
504 END IF
505
506 tstart = m_walltime()
507 DO j = 1, ntim
508 CALL dgemm("N", "T", len, len, len, 1.0_dp, ma, len, mb, len, 1.0_dp, mc, len)
509 END DO
510 tend = m_walltime()
511 t = tend - tstart + threshold
512 IF (t > 0.0_dp) THEN
513 perf = real(ntim, kind=dp)*2.0_dp*real(len, kind=dp)**3*1.e-6_dp/t
514 ELSE
515 perf = 0.0_dp
516 END IF
517
518 IF (para_env%is_source()) THEN
519 WRITE (iw, '(A,i6,T59,F14.4,A)') &
520 " Matrix multiply test: c = a * b(T) Size = ", len, perf, " Mflop/s"
521 END IF
522
523 tstart = m_walltime()
524 DO j = 1, ntim
525 CALL dgemm("T", "N", len, len, len, 1.0_dp, ma, len, mb, len, 1.0_dp, mc, len)
526 END DO
527 tend = m_walltime()
528 t = tend - tstart + threshold
529 IF (t > 0.0_dp) THEN
530 perf = real(ntim, kind=dp)*2.0_dp*real(len, kind=dp)**3*1.e-6_dp/t
531 ELSE
532 perf = 0.0_dp
533 END IF
534
535 IF (para_env%is_source()) THEN
536 WRITE (iw, '(A,i6,T59,F14.4,A)') &
537 " Matrix multiply test: c = a(T) * b Size = ", len, perf, " Mflop/s"
538 END IF
539
540 DEALLOCATE (ma)
541 DEALLOCATE (mb)
542 DEALLOCATE (mc)
543 END DO
544 END IF
545
546 CALL para_env%sync()
547
548 END SUBROUTINE matmul_test
549
550! **************************************************************************************************
551!> \brief Tests the performance of all available FFT libraries for 3D FFTs
552!> \param para_env ...
553!> \param iw ...
554!> \param fftw_plan_type ...
555!> \param wisdom_file where FFTW3 should look to save/load wisdom
556!> \par History
557!> none
558!> \author JGH 6-NOV-2000
559! **************************************************************************************************
560 SUBROUTINE fft_test(para_env, iw, fftw_plan_type, wisdom_file)
561
562 TYPE(mp_para_env_type), POINTER :: para_env
563 INTEGER :: iw, fftw_plan_type
564 CHARACTER(LEN=*), INTENT(IN) :: wisdom_file
565
566 INTEGER, PARAMETER :: ndate(3) = (/12, 48, 96/)
567
568 INTEGER :: iall, ierr, it, j, len, n(3), ntim, &
569 radix_in, radix_out, siz, stat
570 COMPLEX(KIND=dp), DIMENSION(4, 4, 4) :: zz
571 COMPLEX(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: ca, cb, cc
572 CHARACTER(LEN=7) :: method
573 REAL(kind=dp) :: flops, perf, scale, t, tdiff, tend, &
574 tstart
575 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: ra
576
577! test for 3d FFT
578
579 IF (para_env%is_source()) WRITE (iw, '(//,A,/)') " Test of 3D-FFT "
580 siz = abs(runtest(3))
581
582 DO iall = 1, 100
583 SELECT CASE (iall)
584 CASE DEFAULT
585 EXIT
586 CASE (1)
587 CALL init_fft("FFTSG", alltoall=.false., fftsg_sizes=.true., wisdom_file=wisdom_file, &
588 pool_limit=10, plan_style=fftw_plan_type)
589 method = "FFTSG "
590 CASE (2)
591 cycle
592 CASE (3)
593 CALL init_fft("FFTW3", alltoall=.false., fftsg_sizes=.true., wisdom_file=wisdom_file, &
594 pool_limit=10, plan_style=fftw_plan_type)
595 method = "FFTW3 "
596 END SELECT
597 n = 4
598 zz = 0.0_dp
599 CALL fft3d(fwfft, n, zz, status=stat)
600 IF (stat == 0) THEN
601 DO it = 1, 3
602 radix_in = ndate(it)
603 CALL fft_radix_operations(radix_in, radix_out, fft_radix_closest)
604 len = radix_out
605 n = len
606 IF (16.0_dp*real(len*len*len, kind=dp) > max_memory*0.5_dp) EXIT
607 ALLOCATE (ra(len, len, len), stat=ierr)
608 ALLOCATE (ca(len, len, len), stat=ierr)
609 CALL random_number(ra)
610 ca(:, :, :) = ra
611 CALL random_number(ra)
612 ca(:, :, :) = ca + cmplx(0.0_dp, 1.0_dp, kind=dp)*ra
613 flops = real(len**3, kind=dp)*15.0_dp*log(real(len, kind=dp))
614 ntim = nint(siz*1.e7_dp/flops)
615 ntim = max(ntim, 1)
616 ntim = min(ntim, 200)
617 scale = 1.0_dp/real(len**3, kind=dp)
618 tstart = m_walltime()
619 DO j = 1, ntim
620 CALL fft3d(fwfft, n, ca)
621 CALL fft3d(bwfft, n, ca)
622 END DO
623 tend = m_walltime()
624 t = tend - tstart + threshold
625 IF (t > 0.0_dp) THEN
626 perf = real(ntim, kind=dp)*2.0_dp*flops*1.e-6_dp/t
627 ELSE
628 perf = 0.0_dp
629 END IF
630
631 IF (para_env%is_source()) THEN
632 WRITE (iw, '(T2,A,A,i6,T59,F14.4,A)') &
633 adjustr(method), " test (in-place) Size = ", len, perf, " Mflop/s"
634 END IF
635 DEALLOCATE (ca)
636 DEALLOCATE (ra)
637 END DO
638 IF (para_env%is_source()) WRITE (iw, *)
639 ! test if input data is preserved
640 len = 24
641 n = len
642 ALLOCATE (ra(len, len, len))
643 ALLOCATE (ca(len, len, len))
644 ALLOCATE (cb(len, len, len))
645 ALLOCATE (cc(len, len, len))
646 CALL random_number(ra)
647 ca(:, :, :) = ra
648 CALL random_number(ra)
649 ca(:, :, :) = ca + cmplx(0.0_dp, 1.0_dp, kind=dp)*ra
650 cc(:, :, :) = ca
651 CALL fft3d(fwfft, n, ca, cb)
652 tdiff = maxval(abs(ca - cc))
653 IF (tdiff > 1.0e-12_dp) THEN
654 IF (para_env%is_source()) &
655 WRITE (iw, '(T2,A,A,A)') adjustr(method), " FWFFT ", &
656 " Input array is changed in out-of-place FFT !"
657 ELSE
658 IF (para_env%is_source()) &
659 WRITE (iw, '(T2,A,A,A)') adjustr(method), " FWFFT ", &
660 " Input array is not changed in out-of-place FFT !"
661 END IF
662 ca(:, :, :) = cc
663 CALL fft3d(bwfft, n, ca, cb)
664 tdiff = maxval(abs(ca - cc))
665 IF (tdiff > 1.0e-12_dp) THEN
666 IF (para_env%is_source()) &
667 WRITE (iw, '(T2,A,A,A)') adjustr(method), " BWFFT ", &
668 " Input array is changed in out-of-place FFT !"
669 ELSE
670 IF (para_env%is_source()) &
671 WRITE (iw, '(T2,A,A,A)') adjustr(method), " BWFFT ", &
672 " Input array is not changed in out-of-place FFT !"
673 END IF
674 IF (para_env%is_source()) WRITE (iw, *)
675
676 DEALLOCATE (ra)
677 DEALLOCATE (ca)
678 DEALLOCATE (cb)
679 DEALLOCATE (cc)
680 END IF
681 CALL finalize_fft(para_env, wisdom_file=wisdom_file)
682 END DO
683
684 END SUBROUTINE fft_test
685
686! **************************************************************************************************
687!> \brief test rs_pw_transfer performance
688!> \param para_env ...
689!> \param iw ...
690!> \param globenv ...
691!> \param rs_pw_transfer_section ...
692!> \author Joost VandeVondele
693!> 9.2008 Randomise rs grid [Iain Bethune]
694!> (c) The Numerical Algorithms Group (NAG) Ltd, 2008 on behalf of the HECToR project
695! **************************************************************************************************
696 SUBROUTINE rs_pw_transfer_test(para_env, iw, globenv, rs_pw_transfer_section)
697
698 TYPE(mp_para_env_type), POINTER :: para_env
699 INTEGER :: iw
700 TYPE(global_environment_type), POINTER :: globenv
701 TYPE(section_vals_type), POINTER :: rs_pw_transfer_section
702
703 CHARACTER(LEN=*), PARAMETER :: routinen = 'rs_pw_transfer_test'
704
705 INTEGER :: halo_size, handle, i_loop, n_loop, ns_max
706 INTEGER, DIMENSION(3) :: no, np
707 INTEGER, DIMENSION(:), POINTER :: i_vals
708 LOGICAL :: do_rs2pw
709 REAL(kind=dp) :: tend, tstart
710 TYPE(cell_type), POINTER :: box
711 TYPE(pw_grid_type), POINTER :: grid
712 TYPE(pw_r3d_rs_type) :: ca
713 TYPE(realspace_grid_desc_type), POINTER :: rs_desc
714 TYPE(realspace_grid_input_type) :: input_settings
715 TYPE(realspace_grid_type) :: rs_grid
716 TYPE(section_vals_type), POINTER :: rs_grid_section
717
718 CALL timeset(routinen, handle)
719
720 !..set fft lib
721 CALL init_fft(globenv%default_fft_library, alltoall=.false., fftsg_sizes=.true., &
722 pool_limit=globenv%fft_pool_scratch_limit, &
723 wisdom_file=globenv%fftw_wisdom_file_name, &
724 plan_style=globenv%fftw_plan_type)
725
726 ! .. set cell (should otherwise be irrelevant)
727 NULLIFY (box)
728 CALL cell_create(box)
729 box%hmat = reshape((/20.0_dp, 0.0_dp, 0.0_dp, 0.0_dp, 20.0_dp, 0.0_dp, &
730 0.0_dp, 0.0_dp, 20.0_dp/), (/3, 3/))
731 CALL init_cell(box)
732
733 ! .. grid type and pw_grid
734 NULLIFY (grid)
735 CALL section_vals_val_get(rs_pw_transfer_section, "GRID", i_vals=i_vals)
736 np = i_vals
737 CALL pw_grid_create(grid, para_env, box%hmat, grid_span=fullspace, npts=np, fft_usage=.true., iounit=iw)
738 no = grid%npts
739
740 CALL ca%create(grid)
741 CALL pw_zero(ca)
742
743 ! .. rs input setting type
744 CALL section_vals_val_get(rs_pw_transfer_section, "HALO_SIZE", i_val=halo_size)
745 rs_grid_section => section_vals_get_subs_vals(rs_pw_transfer_section, "RS_GRID")
746 ns_max = 2*halo_size + 1
747 CALL init_input_type(input_settings, ns_max, rs_grid_section, 1, (/-1, -1, -1/))
748
749 ! .. rs type
750 NULLIFY (rs_desc)
751 CALL rs_grid_create_descriptor(rs_desc, pw_grid=grid, input_settings=input_settings)
752 CALL rs_grid_create(rs_grid, rs_desc)
753 CALL rs_grid_print(rs_grid, iw)
754 CALL rs_grid_zero(rs_grid)
755
756 ! Put random values on the grid, so summation check will pick up errors
757 CALL random_number(rs_grid%r)
758
759 CALL section_vals_val_get(rs_pw_transfer_section, "N_loop", i_val=n_loop)
760 CALL section_vals_val_get(rs_pw_transfer_section, "RS2PW", l_val=do_rs2pw)
761
762 ! go for the real loops, sync to get max timings
763 IF (para_env%is_source()) THEN
764 WRITE (iw, '(T2,A)') ""
765 WRITE (iw, '(T2,A)') "Timing rs_pw_transfer routine"
766 WRITE (iw, '(T2,A)') ""
767 WRITE (iw, '(T2,A)') "iteration time[s]"
768 END IF
769 DO i_loop = 1, n_loop
770 CALL para_env%sync()
771 tstart = m_walltime()
772 IF (do_rs2pw) THEN
773 CALL transfer_rs2pw(rs_grid, ca)
774 ELSE
775 CALL transfer_pw2rs(rs_grid, ca)
776 END IF
777 CALL para_env%sync()
778 tend = m_walltime()
779 IF (para_env%is_source()) THEN
780 WRITE (iw, '(T2,I9,1X,F12.6)') i_loop, tend - tstart
781 END IF
782 END DO
783
784 !cleanup
785 CALL rs_grid_release(rs_grid)
786 CALL rs_grid_release_descriptor(rs_desc)
787 CALL ca%release()
788 CALL pw_grid_release(grid)
789 CALL cell_release(box)
790 CALL finalize_fft(para_env, wisdom_file=globenv%fftw_wisdom_file_name)
791
792 CALL timestop(handle)
793
794 END SUBROUTINE rs_pw_transfer_test
795
796! **************************************************************************************************
797!> \brief Tests the performance of PW calls to FFT routines
798!> \param para_env ...
799!> \param iw ...
800!> \param globenv ...
801!> \param pw_transfer_section ...
802!> \par History
803!> JGH 6-Feb-2001 : Test and performance code
804!> Made input sensitive [Joost VandeVondele]
805!> \author JGH 1-JAN-2001
806! **************************************************************************************************
807 SUBROUTINE pw_fft_test(para_env, iw, globenv, pw_transfer_section)
808
809 TYPE(mp_para_env_type), POINTER :: para_env
810 INTEGER :: iw
811 TYPE(global_environment_type), POINTER :: globenv
812 TYPE(section_vals_type), POINTER :: pw_transfer_section
813
814 REAL(kind=dp), PARAMETER :: toler = 1.e-11_dp
815
816 INTEGER :: blocked_id, grid_span, i_layout, i_rep, &
817 ig, ip, itmp, n_loop, n_rep, nn, p, q
818 INTEGER, ALLOCATABLE, DIMENSION(:, :) :: layouts
819 INTEGER, DIMENSION(2) :: distribution_layout
820 INTEGER, DIMENSION(3) :: no, np
821 INTEGER, DIMENSION(:), POINTER :: i_vals
822 LOGICAL :: debug, is_fullspace, odd, &
823 pw_grid_layout_all, spherical
824 REAL(kind=dp) :: em, et, flops, gsq, perf, t, t_max, &
825 t_min, tend, tstart
826 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: t_end, t_start
827 TYPE(cell_type), POINTER :: box
828 TYPE(pw_c1d_gs_type) :: ca, cc
829 TYPE(pw_c3d_rs_type) :: cb
830 TYPE(pw_grid_type), POINTER :: grid
831
832!..set fft lib
833
834 CALL init_fft(globenv%default_fft_library, alltoall=.false., fftsg_sizes=.true., &
835 pool_limit=globenv%fft_pool_scratch_limit, &
836 wisdom_file=globenv%fftw_wisdom_file_name, &
837 plan_style=globenv%fftw_plan_type)
838
839 !..the unit cell (should not really matter, the number of grid points do)
840 NULLIFY (box, grid)
841 CALL cell_create(box)
842 box%hmat = reshape((/10.0_dp, 0.0_dp, 0.0_dp, 0.0_dp, 8.0_dp, 0.0_dp, &
843 0.0_dp, 0.0_dp, 7.0_dp/), (/3, 3/))
844 CALL init_cell(box)
845
846 CALL section_vals_get(pw_transfer_section, n_repetition=n_rep)
847 DO i_rep = 1, n_rep
848
849 ! how often should we do the transfer
850 CALL section_vals_val_get(pw_transfer_section, "N_loop", i_rep_section=i_rep, i_val=n_loop)
851 ALLOCATE (t_start(n_loop))
852 ALLOCATE (t_end(n_loop))
853
854 ! setup of the grids
855 CALL section_vals_val_get(pw_transfer_section, "GRID", i_rep_section=i_rep, i_vals=i_vals)
856 np = i_vals
857
858 CALL section_vals_val_get(pw_transfer_section, "PW_GRID_BLOCKED", i_rep_section=i_rep, i_val=blocked_id)
859 CALL section_vals_val_get(pw_transfer_section, "DEBUG", i_rep_section=i_rep, l_val=debug)
860
861 CALL section_vals_val_get(pw_transfer_section, "PW_GRID_LAYOUT_ALL", i_rep_section=i_rep, &
862 l_val=pw_grid_layout_all)
863
864 ! prepare to loop over all or a specific layout
865 IF (pw_grid_layout_all) THEN
866 ! count layouts that fit
867 itmp = 0
868 ! start from 2, (/1,para_env%num_pe/) is not supported
869 DO p = 2, para_env%num_pe
870 q = para_env%num_pe/p
871 IF (p*q == para_env%num_pe) THEN
872 itmp = itmp + 1
873 END IF
874 END DO
875 ! build list
876 ALLOCATE (layouts(2, itmp))
877 itmp = 0
878 DO p = 2, para_env%num_pe
879 q = para_env%num_pe/p
880 IF (p*q == para_env%num_pe) THEN
881 itmp = itmp + 1
882 layouts(:, itmp) = (/p, q/)
883 END IF
884 END DO
885 ELSE
886 CALL section_vals_val_get(pw_transfer_section, "PW_GRID_LAYOUT", i_rep_section=i_rep, i_vals=i_vals)
887 ALLOCATE (layouts(2, 1))
888 layouts(:, 1) = i_vals
889 END IF
890
891 DO i_layout = 1, SIZE(layouts, 2)
892
893 distribution_layout = layouts(:, i_layout)
894
895 CALL section_vals_val_get(pw_transfer_section, "PW_GRID", i_rep_section=i_rep, i_val=itmp)
896
897 ! from cp_control_utils
898 SELECT CASE (itmp)
900 spherical = .true.
901 is_fullspace = .false.
903 spherical = .false.
904 is_fullspace = .true.
906 spherical = .false.
907 is_fullspace = .false.
908 END SELECT
909
910 ! from pw_env_methods
911 IF (spherical) THEN
912 grid_span = halfspace
913 spherical = .true.
914 odd = .true.
915 ELSE IF (is_fullspace) THEN
916 grid_span = fullspace
917 spherical = .false.
918 odd = .false.
919 ELSE
920 grid_span = halfspace
921 spherical = .false.
922 odd = .true.
923 END IF
924
925 ! actual setup
926 CALL pw_grid_create(grid, para_env, box%hmat, grid_span=grid_span, odd=odd, spherical=spherical, &
927 blocked=blocked_id, npts=np, fft_usage=.true., &
928 rs_dims=distribution_layout, iounit=iw)
929
930 IF (iw > 0) CALL m_flush(iw)
931
932 ! note that the number of grid points might be different from what the user requested (fft-able needed)
933 no = grid%npts
934
935 CALL ca%create(grid)
936 CALL cb%create(grid)
937 CALL cc%create(grid)
938
939 ! initialize data
940 CALL pw_zero(ca)
941 CALL pw_zero(cb)
942 CALL pw_zero(cc)
943 nn = SIZE(ca%array)
944 DO ig = 1, nn
945 gsq = grid%gsq(ig)
946 ca%array(ig) = exp(-gsq)
947 END DO
948
949 flops = product(no)*30.0_dp*log(real(maxval(no), kind=dp))
950 tstart = m_walltime()
951 DO ip = 1, n_loop
952 CALL para_env%sync()
953 t_start(ip) = m_walltime()
954 CALL pw_transfer(ca, cb, debug)
955 CALL pw_transfer(cb, cc, debug)
956 CALL para_env%sync()
957 t_end(ip) = m_walltime()
958 END DO
959 tend = m_walltime()
960 t = tend - tstart + threshold
961 IF (t > 0.0_dp) THEN
962 perf = real(n_loop, kind=dp)*2.0_dp*flops*1.e-6_dp/t
963 ELSE
964 perf = 0.0_dp
965 END IF
966
967 em = maxval(abs(ca%array(:) - cc%array(:)))
968 CALL para_env%max(em)
969 et = sum(abs(ca%array(:) - cc%array(:)))
970 CALL para_env%sum(et)
971 t_min = minval(t_end - t_start)
972 t_max = maxval(t_end - t_start)
973
974 IF (para_env%is_source()) THEN
975 WRITE (iw, *)
976 WRITE (iw, '(A,T67,E14.6)') " Parallel FFT Tests: Maximal Error ", em
977 WRITE (iw, '(A,T67,E14.6)') " Parallel FFT Tests: Total Error ", et
978 WRITE (iw, '(A,T67,F14.0)') &
979 " Parallel FFT Tests: Performance [Mflops] ", perf
980 WRITE (iw, '(A,T67,F14.6)') " Best time : ", t_min
981 WRITE (iw, '(A,T67,F14.6)') " Worst time: ", t_max
982 IF (iw > 0) CALL m_flush(iw)
983 END IF
984
985 ! need debugging ???
986 IF (em > toler .OR. et > toler) THEN
987 cpwarn("The FFT results are not accurate ... starting debug pw_transfer")
988 CALL pw_transfer(ca, cb, .true.)
989 CALL pw_transfer(cb, cc, .true.)
990 END IF
991
992 ! done with these grids
993 CALL ca%release()
994 CALL cb%release()
995 CALL cc%release()
996 CALL pw_grid_release(grid)
997
998 END DO
999
1000 ! local arrays
1001 DEALLOCATE (layouts)
1002 DEALLOCATE (t_start)
1003 DEALLOCATE (t_end)
1004
1005 END DO
1006
1007 ! cleanup
1008 CALL cell_release(box)
1009 CALL finalize_fft(para_env, wisdom_file=globenv%fftw_wisdom_file_name)
1010
1011 END SUBROUTINE pw_fft_test
1012
1013! **************************************************************************************************
1014!> \brief Tests the eigensolver library routines
1015!> \param para_env ...
1016!> \param iw ...
1017!> \param eigensolver_section ...
1018!> \par History
1019!> JGH 6-Feb-2001 : Test and performance code
1020!> \author JGH 1-JAN-2001
1021! **************************************************************************************************
1022 SUBROUTINE eigensolver_test(para_env, iw, eigensolver_section)
1023
1024 TYPE(mp_para_env_type), POINTER :: para_env
1025 INTEGER :: iw
1026 TYPE(section_vals_type), POINTER :: eigensolver_section
1027
1028 INTEGER :: diag_method, i, i_loop, i_rep, &
1029 init_method, j, n, n_loop, n_rep, &
1030 neig, unit_number
1031 REAL(kind=dp) :: t1, t2
1032 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: eigenvalues
1033 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: buffer
1034 TYPE(cp_blacs_env_type), POINTER :: blacs_env
1035 TYPE(cp_fm_struct_type), POINTER :: fmstruct
1036 TYPE(cp_fm_type) :: eigenvectors, matrix, work
1037 TYPE(rng_stream_type), ALLOCATABLE :: rng_stream
1038
1039 IF (iw > 0) THEN
1040 WRITE (unit=iw, fmt="(/,/,T2,A,/)") "EIGENSOLVER TEST"
1041 END IF
1042
1043 ! create blacs env corresponding to para_env
1044 NULLIFY (blacs_env)
1045 CALL cp_blacs_env_create(blacs_env=blacs_env, &
1046 para_env=para_env)
1047
1048 ! loop over all tests
1049 CALL section_vals_get(eigensolver_section, n_repetition=n_rep)
1050 DO i_rep = 1, n_rep
1051
1052 ! parse section
1053 CALL section_vals_val_get(eigensolver_section, "N", i_rep_section=i_rep, i_val=n)
1054 CALL section_vals_val_get(eigensolver_section, "EIGENVALUES", i_rep_section=i_rep, i_val=neig)
1055 CALL section_vals_val_get(eigensolver_section, "DIAG_METHOD", i_rep_section=i_rep, i_val=diag_method)
1056 CALL section_vals_val_get(eigensolver_section, "INIT_METHOD", i_rep_section=i_rep, i_val=init_method)
1057 CALL section_vals_val_get(eigensolver_section, "N_loop", i_rep_section=i_rep, i_val=n_loop)
1058
1059 ! proper number of eigs
1060 IF (neig < 0) neig = n
1061 neig = min(neig, n)
1062
1063 ! report
1064 IF (iw > 0) THEN
1065 WRITE (iw, *) "Matrix size", n
1066 WRITE (iw, *) "Number of eigenvalues", neig
1067 WRITE (iw, *) "Timing loops", n_loop
1068 SELECT CASE (diag_method)
1069 CASE (do_diag_syevd)
1070 WRITE (iw, *) "Diag using syevd"
1071 CASE (do_diag_syevx)
1072 WRITE (iw, *) "Diag using syevx"
1073 CASE DEFAULT
1074 ! stop
1075 END SELECT
1076
1077 SELECT CASE (init_method)
1078 CASE (do_mat_random)
1079 WRITE (iw, *) "using random matrix"
1080 CASE (do_mat_read)
1081 WRITE (iw, *) "reading from file"
1082 CASE DEFAULT
1083 ! stop
1084 END SELECT
1085 END IF
1086
1087 ! create matrix struct type
1088 NULLIFY (fmstruct)
1089 CALL cp_fm_struct_create(fmstruct=fmstruct, &
1090 para_env=para_env, &
1091 context=blacs_env, &
1092 nrow_global=n, &
1093 ncol_global=n)
1094
1095 ! create all needed matrices, and buffers for the eigenvalues
1096 CALL cp_fm_create(matrix=matrix, &
1097 matrix_struct=fmstruct, &
1098 name="MATRIX")
1099 CALL cp_fm_set_all(matrix, 0.0_dp)
1100
1101 CALL cp_fm_create(matrix=eigenvectors, &
1102 matrix_struct=fmstruct, &
1103 name="EIGENVECTORS")
1104 CALL cp_fm_set_all(eigenvectors, 0.0_dp)
1105
1106 CALL cp_fm_create(matrix=work, &
1107 matrix_struct=fmstruct, &
1108 name="WORK")
1109 CALL cp_fm_set_all(matrix, 0.0_dp)
1110
1111 ALLOCATE (eigenvalues(n))
1112 eigenvalues = 0.0_dp
1113 ALLOCATE (buffer(1, n))
1114
1115 ! generate initial matrix, either by reading a file, or using random numbers
1116 IF (para_env%is_source()) THEN
1117 SELECT CASE (init_method)
1118 CASE (do_mat_random)
1119 rng_stream = rng_stream_type( &
1120 name="rng_stream", &
1121 distribution_type=uniform, &
1122 extended_precision=.true.)
1123 CASE (do_mat_read)
1124 CALL open_file(file_name="MATRIX", &
1125 file_action="READ", &
1126 file_form="FORMATTED", &
1127 file_status="OLD", &
1128 unit_number=unit_number)
1129 END SELECT
1130 END IF
1131
1132 DO i = 1, n
1133 IF (para_env%is_source()) THEN
1134 SELECT CASE (init_method)
1135 CASE (do_mat_random)
1136 DO j = i, n
1137 buffer(1, j) = rng_stream%next() - 0.5_dp
1138 END DO
1139 !MK activate/modify for a diagonal dominant symmetric matrix:
1140 !MK buffer(1,i) = 10.0_dp*buffer(1,i)
1141 CASE (do_mat_read)
1142 READ (unit=unit_number, fmt=*) buffer(1, 1:n)
1143 END SELECT
1144 END IF
1145 CALL para_env%bcast(buffer)
1146 SELECT CASE (init_method)
1147 CASE (do_mat_random)
1148 CALL cp_fm_set_submatrix(fm=matrix, &
1149 new_values=buffer, &
1150 start_row=i, &
1151 start_col=i, &
1152 n_rows=1, &
1153 n_cols=n - i + 1, &
1154 alpha=1.0_dp, &
1155 beta=0.0_dp, &
1156 transpose=.false.)
1157 CALL cp_fm_set_submatrix(fm=matrix, &
1158 new_values=buffer, &
1159 start_row=i, &
1160 start_col=i, &
1161 n_rows=n - i + 1, &
1162 n_cols=1, &
1163 alpha=1.0_dp, &
1164 beta=0.0_dp, &
1165 transpose=.true.)
1166 CASE (do_mat_read)
1167 CALL cp_fm_set_submatrix(fm=matrix, &
1168 new_values=buffer, &
1169 start_row=i, &
1170 start_col=1, &
1171 n_rows=1, &
1172 n_cols=n, &
1173 alpha=1.0_dp, &
1174 beta=0.0_dp, &
1175 transpose=.false.)
1176 END SELECT
1177 END DO
1178
1179 DEALLOCATE (buffer)
1180
1181 IF (para_env%is_source()) THEN
1182 SELECT CASE (init_method)
1183 CASE (do_mat_read)
1184 CALL close_file(unit_number=unit_number)
1185 END SELECT
1186 END IF
1187
1188 DO i_loop = 1, n_loop
1189 eigenvalues = 0.0_dp
1190 CALL cp_fm_set_all(eigenvectors, 0.0_dp)
1191 CALL cp_fm_to_fm(source=matrix, &
1192 destination=work)
1193
1194 ! DONE, now testing
1195 t1 = m_walltime()
1196 SELECT CASE (diag_method)
1197 CASE (do_diag_syevd)
1198 CALL cp_fm_syevd(matrix=work, &
1199 eigenvectors=eigenvectors, &
1200 eigenvalues=eigenvalues)
1201 CASE (do_diag_syevx)
1202 CALL cp_fm_syevx(matrix=work, &
1203 eigenvectors=eigenvectors, &
1204 eigenvalues=eigenvalues, &
1205 neig=neig, &
1206 work_syevx=1.0_dp)
1207 END SELECT
1208 t2 = m_walltime()
1209 IF (iw > 0) WRITE (iw, *) "Timing for loop ", i_loop, " : ", t2 - t1
1210 END DO
1211
1212 IF (iw > 0) THEN
1213 WRITE (iw, *) "Eigenvalues: "
1214 WRITE (unit=iw, fmt="(T3,5F14.6)") eigenvalues(1:neig)
1215 WRITE (unit=iw, fmt="(T3,A4,F16.6)") "Sum:", sum(eigenvalues(1:neig))
1216 WRITE (iw, *) ""
1217 END IF
1218
1219 ! Clean up
1220 DEALLOCATE (eigenvalues)
1221 CALL cp_fm_release(matrix=work)
1222 CALL cp_fm_release(matrix=eigenvectors)
1223 CALL cp_fm_release(matrix=matrix)
1224 CALL cp_fm_struct_release(fmstruct=fmstruct)
1225
1226 END DO
1227
1228 CALL cp_blacs_env_release(blacs_env=blacs_env)
1229
1230 END SUBROUTINE eigensolver_test
1231
1232! **************************************************************************************************
1233!> \brief Tests the parallel matrix multiply
1234!> \param para_env ...
1235!> \param iw ...
1236!> \param cp_fm_gemm_test_section ...
1237! **************************************************************************************************
1238 SUBROUTINE cp_fm_gemm_test(para_env, iw, cp_fm_gemm_test_section)
1239
1240 TYPE(mp_para_env_type), POINTER :: para_env
1241 INTEGER :: iw
1242 TYPE(section_vals_type), POINTER :: cp_fm_gemm_test_section
1243
1244 CHARACTER(LEN=1) :: transa, transb
1245 INTEGER :: i_loop, i_rep, k, m, n, n_loop, n_rep, ncol_block, ncol_block_actual, &
1246 ncol_global, np, nrow_block, nrow_block_actual, nrow_global
1247 INTEGER, DIMENSION(:), POINTER :: grid_2d
1248 LOGICAL :: force_blocksize, row_major, transa_p, &
1249 transb_p
1250 REAL(kind=dp) :: t1, t2, t3, t4
1251 TYPE(cp_blacs_env_type), POINTER :: blacs_env
1252 TYPE(cp_fm_struct_type), POINTER :: fmstruct_a, fmstruct_b, fmstruct_c
1253 TYPE(cp_fm_type) :: matrix_a, matrix_b, matrix_c
1254
1255 CALL section_vals_get(cp_fm_gemm_test_section, n_repetition=n_rep)
1256 DO i_rep = 1, n_rep
1257
1258 ! how often should we do the multiply
1259 CALL section_vals_val_get(cp_fm_gemm_test_section, "N_loop", i_rep_section=i_rep, i_val=n_loop)
1260
1261 ! matrices def.
1262 CALL section_vals_val_get(cp_fm_gemm_test_section, "K", i_rep_section=i_rep, i_val=k)
1263 CALL section_vals_val_get(cp_fm_gemm_test_section, "N", i_rep_section=i_rep, i_val=n)
1264 CALL section_vals_val_get(cp_fm_gemm_test_section, "M", i_rep_section=i_rep, i_val=m)
1265 CALL section_vals_val_get(cp_fm_gemm_test_section, "transa", i_rep_section=i_rep, l_val=transa_p)
1266 CALL section_vals_val_get(cp_fm_gemm_test_section, "transb", i_rep_section=i_rep, l_val=transb_p)
1267 CALL section_vals_val_get(cp_fm_gemm_test_section, "nrow_block", i_rep_section=i_rep, i_val=nrow_block)
1268 CALL section_vals_val_get(cp_fm_gemm_test_section, "ncol_block", i_rep_section=i_rep, i_val=ncol_block)
1269 CALL section_vals_val_get(cp_fm_gemm_test_section, "ROW_MAJOR", i_rep_section=i_rep, l_val=row_major)
1270 CALL section_vals_val_get(cp_fm_gemm_test_section, "GRID_2D", i_rep_section=i_rep, i_vals=grid_2d)
1271 CALL section_vals_val_get(cp_fm_gemm_test_section, "FORCE_BLOCKSIZE", i_rep_section=i_rep, l_val=force_blocksize)
1272 transa = "N"
1273 transb = "N"
1274 IF (transa_p) transa = "T"
1275 IF (transb_p) transb = "T"
1276
1277 IF (iw > 0) THEN
1278 WRITE (iw, '(T2,A)') "----------- TESTING PARALLEL MATRIX MULTIPLY -------------"
1279 WRITE (iw, '(T2,A)', advance="NO") "C = "
1280 IF (transa_p) THEN
1281 WRITE (iw, '(A)', advance="NO") "TRANSPOSE(A) x"
1282 ELSE
1283 WRITE (iw, '(A)', advance="NO") "A x "
1284 END IF
1285 IF (transb_p) THEN
1286 WRITE (iw, '(A)') "TRANSPOSE(B) "
1287 ELSE
1288 WRITE (iw, '(A)') "B "
1289 END IF
1290 WRITE (iw, '(T2,A,T50,I5,A,I5)') 'requested block size', nrow_block, ' by ', ncol_block
1291 WRITE (iw, '(T2,A,T50,I5)') 'number of repetitions of cp_fm_gemm ', n_loop
1292 WRITE (iw, '(T2,A,T50,L5)') 'Row Major', row_major
1293 WRITE (iw, '(T2,A,T50,2I7)') 'GRID_2D ', grid_2d
1294 WRITE (iw, '(T2,A,T50,L5)') 'Force blocksize ', force_blocksize
1295 ! check the return value of pilaenv, too small values limit the performance (assuming pdgemm is the vanilla variant)
1296 np = cp_fm_pilaenv(0, 'D')
1297 IF (np > 0) THEN
1298 WRITE (iw, '(T2,A,T50,I5)') 'PILAENV blocksize', np
1299 END IF
1300 END IF
1301
1302 NULLIFY (blacs_env)
1303 CALL cp_blacs_env_create(blacs_env=blacs_env, &
1304 para_env=para_env, &
1305 row_major=row_major, &
1306 grid_2d=grid_2d)
1307
1308 NULLIFY (fmstruct_a)
1309 IF (transa_p) THEN
1310 nrow_global = m; ncol_global = k
1311 ELSE
1312 nrow_global = k; ncol_global = m
1313 END IF
1314 CALL cp_fm_struct_create(fmstruct=fmstruct_a, para_env=para_env, context=blacs_env, &
1315 nrow_global=nrow_global, ncol_global=ncol_global, &
1316 nrow_block=nrow_block, ncol_block=ncol_block, force_block=force_blocksize)
1317 CALL cp_fm_struct_get(fmstruct_a, nrow_block=nrow_block_actual, ncol_block=ncol_block_actual)
1318 IF (iw > 0) WRITE (iw, '(T2,A,I9,A,I9,A,I5,A,I5)') 'matrix A ', nrow_global, " by ", ncol_global, &
1319 ' using blocks of ', nrow_block_actual, ' by ', ncol_block_actual
1320
1321 IF (transb_p) THEN
1322 nrow_global = n; ncol_global = m
1323 ELSE
1324 nrow_global = m; ncol_global = n
1325 END IF
1326 NULLIFY (fmstruct_b)
1327 CALL cp_fm_struct_create(fmstruct=fmstruct_b, para_env=para_env, context=blacs_env, &
1328 nrow_global=nrow_global, ncol_global=ncol_global, &
1329 nrow_block=nrow_block, ncol_block=ncol_block, force_block=force_blocksize)
1330 CALL cp_fm_struct_get(fmstruct_b, nrow_block=nrow_block_actual, ncol_block=ncol_block_actual)
1331 IF (iw > 0) WRITE (iw, '(T2,A,I9,A,I9,A,I5,A,I5)') 'matrix B ', nrow_global, " by ", ncol_global, &
1332 ' using blocks of ', nrow_block_actual, ' by ', ncol_block_actual
1333
1334 NULLIFY (fmstruct_c)
1335 nrow_global = k
1336 ncol_global = n
1337 CALL cp_fm_struct_create(fmstruct=fmstruct_c, para_env=para_env, context=blacs_env, &
1338 nrow_global=nrow_global, ncol_global=ncol_global, &
1339 nrow_block=nrow_block, ncol_block=ncol_block, force_block=force_blocksize)
1340 CALL cp_fm_struct_get(fmstruct_c, nrow_block=nrow_block_actual, ncol_block=ncol_block_actual)
1341 IF (iw > 0) WRITE (iw, '(T2,A,I9,A,I9,A,I5,A,I5)') 'matrix C ', nrow_global, " by ", ncol_global, &
1342 ' using blocks of ', nrow_block_actual, ' by ', ncol_block_actual
1343
1344 CALL cp_fm_create(matrix=matrix_a, matrix_struct=fmstruct_a, name="MATRIX A")
1345 CALL cp_fm_create(matrix=matrix_b, matrix_struct=fmstruct_b, name="MATRIX B")
1346 CALL cp_fm_create(matrix=matrix_c, matrix_struct=fmstruct_c, name="MATRIX C")
1347
1348 CALL random_number(matrix_a%local_data)
1349 CALL random_number(matrix_b%local_data)
1350 CALL random_number(matrix_c%local_data)
1351
1352 IF (iw > 0) CALL m_flush(iw)
1353
1354 t1 = m_walltime()
1355 DO i_loop = 1, n_loop
1356 t3 = m_walltime()
1357 CALL parallel_gemm(transa, transb, k, n, m, 1.0_dp, matrix_a, matrix_b, 0.0_dp, matrix_c)
1358 t4 = m_walltime()
1359 IF (iw > 0) THEN
1360 WRITE (iw, '(T2,A,T50,F12.6)') "cp_fm_gemm timing: ", (t4 - t3)
1361 CALL m_flush(iw)
1362 END IF
1363 END DO
1364 t2 = m_walltime()
1365
1366 IF (iw > 0) THEN
1367 WRITE (iw, '(T2,A,T50,F12.6)') "average cp_fm_gemm timing: ", (t2 - t1)/n_loop
1368 IF (t2 > t1) THEN
1369 WRITE (iw, '(T2,A,T50,F12.6)') "cp_fm_gemm Gflops per MPI task: ", &
1370 2*real(m, kind=dp)*real(n, kind=dp)*real(k, kind=dp)*n_loop/max(0.001_dp, t2 - t1)/1.0e9_dp/para_env%num_pe
1371 END IF
1372 END IF
1373
1374 CALL cp_fm_release(matrix=matrix_a)
1375 CALL cp_fm_release(matrix=matrix_b)
1376 CALL cp_fm_release(matrix=matrix_c)
1377 CALL cp_fm_struct_release(fmstruct=fmstruct_a)
1378 CALL cp_fm_struct_release(fmstruct=fmstruct_b)
1379 CALL cp_fm_struct_release(fmstruct=fmstruct_c)
1380 CALL cp_blacs_env_release(blacs_env=blacs_env)
1381
1382 END DO
1383
1384 END SUBROUTINE cp_fm_gemm_test
1385
1386! **************************************************************************************************
1387!> \brief Tests the DBCSR interface.
1388!> \param para_env ...
1389!> \param iw ...
1390!> \param input_section ...
1391! **************************************************************************************************
1392 SUBROUTINE cp_dbcsr_tests(para_env, iw, input_section)
1393
1394 TYPE(mp_para_env_type), POINTER :: para_env
1395 INTEGER :: iw
1396 TYPE(section_vals_type), POINTER :: input_section
1397
1398 CHARACTER, DIMENSION(3) :: types
1399 INTEGER :: data_type, i_rep, k, m, n, n_loop, &
1400 n_rep, test_type
1401 INTEGER, DIMENSION(:), POINTER :: bs_k, bs_m, bs_n, nproc
1402 LOGICAL :: always_checksum, retain_sparsity, &
1403 transa_p, transb_p
1404 REAL(kind=dp) :: alpha, beta, filter_eps, s_a, s_b, s_c
1405
1406! ---------------------------------------------------------------------------
1407
1408 NULLIFY (bs_m, bs_n, bs_k)
1409 CALL section_vals_get(input_section, n_repetition=n_rep)
1410 CALL dbcsr_reset_randmat_seed()
1411 DO i_rep = 1, n_rep
1412 ! how often should we do the multiply
1413 CALL section_vals_val_get(input_section, "N_loop", i_rep_section=i_rep, i_val=n_loop)
1414
1415 ! matrices def.
1416 CALL section_vals_val_get(input_section, "TEST_TYPE", i_rep_section=i_rep, i_val=test_type)
1417 CALL section_vals_val_get(input_section, "DATA_TYPE", i_rep_section=i_rep, i_val=data_type)
1418 CALL section_vals_val_get(input_section, "K", i_rep_section=i_rep, i_val=k)
1419 CALL section_vals_val_get(input_section, "N", i_rep_section=i_rep, i_val=n)
1420 CALL section_vals_val_get(input_section, "M", i_rep_section=i_rep, i_val=m)
1421 CALL section_vals_val_get(input_section, "transa", i_rep_section=i_rep, l_val=transa_p)
1422 CALL section_vals_val_get(input_section, "transb", i_rep_section=i_rep, l_val=transb_p)
1423 CALL section_vals_val_get(input_section, "bs_m", i_rep_section=i_rep, &
1424 i_vals=bs_m)
1425 CALL section_vals_val_get(input_section, "bs_n", i_rep_section=i_rep, &
1426 i_vals=bs_n)
1427 CALL section_vals_val_get(input_section, "bs_k", i_rep_section=i_rep, &
1428 i_vals=bs_k)
1429 CALL section_vals_val_get(input_section, "keepsparse", i_rep_section=i_rep, l_val=retain_sparsity)
1430 CALL section_vals_val_get(input_section, "asparsity", i_rep_section=i_rep, r_val=s_a)
1431 CALL section_vals_val_get(input_section, "bsparsity", i_rep_section=i_rep, r_val=s_b)
1432 CALL section_vals_val_get(input_section, "csparsity", i_rep_section=i_rep, r_val=s_c)
1433 CALL section_vals_val_get(input_section, "alpha", i_rep_section=i_rep, r_val=alpha)
1434 CALL section_vals_val_get(input_section, "beta", i_rep_section=i_rep, r_val=beta)
1435 CALL section_vals_val_get(input_section, "nproc", i_rep_section=i_rep, &
1436 i_vals=nproc)
1437 CALL section_vals_val_get(input_section, "atype", i_rep_section=i_rep, &
1438 c_val=types(1))
1439 CALL section_vals_val_get(input_section, "btype", i_rep_section=i_rep, &
1440 c_val=types(2))
1441 CALL section_vals_val_get(input_section, "ctype", i_rep_section=i_rep, &
1442 c_val=types(3))
1443 CALL section_vals_val_get(input_section, "filter_eps", &
1444 i_rep_section=i_rep, r_val=filter_eps)
1445 CALL section_vals_val_get(input_section, "ALWAYS_CHECKSUM", i_rep_section=i_rep, l_val=always_checksum)
1446
1447 CALL dbcsr_run_tests(para_env%get_handle(), iw, nproc, &
1448 (/m, n, k/), &
1449 (/transa_p, transb_p/), &
1450 bs_m, bs_n, bs_k, &
1451 (/s_a, s_b, s_c/), &
1452 alpha, beta, &
1453 data_type=data_type, &
1454 test_type=test_type, &
1455 n_loops=n_loop, eps=filter_eps, retain_sparsity=retain_sparsity, &
1456 always_checksum=always_checksum)
1457 END DO
1458 END SUBROUTINE cp_dbcsr_tests
1459
1460! **************************************************************************************************
1461!> \brief Tests the DBM library.
1462!> \param para_env ...
1463!> \param iw ...
1464!> \param input_section ...
1465! **************************************************************************************************
1466 SUBROUTINE run_dbm_tests(para_env, iw, input_section)
1467
1468 TYPE(mp_para_env_type), POINTER :: para_env
1469 INTEGER :: iw
1470 TYPE(section_vals_type), POINTER :: input_section
1471
1472 INTEGER :: i_rep, k, m, n, n_loop, n_rep
1473 INTEGER, DIMENSION(:), POINTER :: bs_k, bs_m, bs_n
1474 LOGICAL :: always_checksum, retain_sparsity, &
1475 transa_p, transb_p
1476 REAL(kind=dp) :: alpha, beta, filter_eps, s_a, s_b, s_c
1477
1478! ---------------------------------------------------------------------------
1479
1480 NULLIFY (bs_m, bs_n, bs_k)
1481 CALL section_vals_get(input_section, n_repetition=n_rep)
1482 CALL dbcsr_reset_randmat_seed()
1483 DO i_rep = 1, n_rep
1484 CALL section_vals_val_get(input_section, "N_loop", i_rep_section=i_rep, i_val=n_loop)
1485 CALL section_vals_val_get(input_section, "K", i_rep_section=i_rep, i_val=k)
1486 CALL section_vals_val_get(input_section, "N", i_rep_section=i_rep, i_val=n)
1487 CALL section_vals_val_get(input_section, "M", i_rep_section=i_rep, i_val=m)
1488 CALL section_vals_val_get(input_section, "transa", i_rep_section=i_rep, l_val=transa_p)
1489 CALL section_vals_val_get(input_section, "transb", i_rep_section=i_rep, l_val=transb_p)
1490 CALL section_vals_val_get(input_section, "bs_m", i_rep_section=i_rep, i_vals=bs_m)
1491 CALL section_vals_val_get(input_section, "bs_n", i_rep_section=i_rep, i_vals=bs_n)
1492 CALL section_vals_val_get(input_section, "bs_k", i_rep_section=i_rep, i_vals=bs_k)
1493 CALL section_vals_val_get(input_section, "keepsparse", i_rep_section=i_rep, l_val=retain_sparsity)
1494 CALL section_vals_val_get(input_section, "asparsity", i_rep_section=i_rep, r_val=s_a)
1495 CALL section_vals_val_get(input_section, "bsparsity", i_rep_section=i_rep, r_val=s_b)
1496 CALL section_vals_val_get(input_section, "csparsity", i_rep_section=i_rep, r_val=s_c)
1497 CALL section_vals_val_get(input_section, "alpha", i_rep_section=i_rep, r_val=alpha)
1498 CALL section_vals_val_get(input_section, "beta", i_rep_section=i_rep, r_val=beta)
1499 CALL section_vals_val_get(input_section, "filter_eps", i_rep_section=i_rep, r_val=filter_eps)
1500 CALL section_vals_val_get(input_section, "ALWAYS_CHECKSUM", i_rep_section=i_rep, l_val=always_checksum)
1501
1502 CALL dbm_run_tests(mp_group=para_env, &
1503 io_unit=iw, &
1504 matrix_sizes=(/m, n, k/), &
1505 trs=(/transa_p, transb_p/), &
1506 bs_m=bs_m, &
1507 bs_n=bs_n, &
1508 bs_k=bs_k, &
1509 sparsities=(/s_a, s_b, s_c/), &
1510 alpha=alpha, &
1511 beta=beta, &
1512 n_loops=n_loop, &
1513 eps=filter_eps, &
1514 retain_sparsity=retain_sparsity, &
1515 always_checksum=always_checksum)
1516 END DO
1517 END SUBROUTINE run_dbm_tests
1518
1519END MODULE library_tests
static void dgemm(const char transa, const char transb, const int m, const int n, const int k, const double alpha, const double *a, const int lda, const double *b, const int ldb, const double beta, double *c, const int ldc)
Convenient wrapper to hide Fortran nature of dgemm_, swapping a and b.
Test of Electron Repulsion Routines (ERI)
real(kind=dp), parameter threshold
subroutine, public eri_test(iw)
...
Handles all functions related to the CELL.
subroutine, public init_cell(cell, hmat, periodic)
Initialise/readjust a simulation cell after hmat has been changed.
subroutine, public cell_create(cell, hmat, periodic, tag)
allocates and initializes a cell
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public cell_release(cell)
releases the given cell (see doc/ReferenceCounting.html)
Definition cell_types.F:559
Test of Clebsch-Gordon Coefficients.
Definition cg_test.F:14
subroutine, public clebsch_gordon_test()
...
Definition cg_test.F:42
methods related to the blacs parallel environment
subroutine, public cp_blacs_env_release(blacs_env)
releases the given blacs_env
subroutine, public cp_blacs_env_create(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
allocates and initializes a type that represent a blacs context
Interface to Minimax-Ewald method for periodic ERI's to be used in CP2K.
subroutine, public cp_eri_mme_perf_acc_test(para_env, iw, eri_mme_test_section)
...
Utility routines to open and close files. Tracking of preconnections.
Definition cp_files.F:16
subroutine, public open_file(file_name, file_status, file_form, file_action, file_position, file_pad, unit_number, debug, skip_get_unit_number, file_access)
Opens the requested file using a free unit number.
Definition cp_files.F:308
subroutine, public close_file(unit_number, file_status, keep_preconnection)
Close an open file given by its logical unit number. Optionally, keep the file and unit preconnected.
Definition cp_files.F:119
basic linear algebra operations for full matrices
subroutine, public cp_fm_gemm(transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, matrix_c, a_first_col, a_first_row, b_first_col, b_first_row, c_first_col, c_first_row)
computes matrix_c = beta * matrix_c + alpha * ( matrix_a ** transa ) * ( matrix_b ** transb )
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
Definition cp_fm_diag.F:17
subroutine, public cp_fm_syevd(matrix, eigenvectors, eigenvalues, info)
Computes all eigenvalues and vectors of a real symmetric matrix significantly faster than syevx,...
Definition cp_fm_diag.F:434
subroutine, public cp_fm_syevx(matrix, eigenvectors, eigenvalues, neig, work_syevx)
compute eigenvalues and optionally eigenvectors of a real symmetric matrix using scalapack....
Definition cp_fm_diag.F:678
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_get(fmstruct, para_env, context, descriptor, ncol_block, nrow_block, nrow_global, ncol_global, first_p_pos, row_indices, col_indices, nrow_local, ncol_local, nrow_locals, ncol_locals, local_leading_dimension)
returns the values of various attributes of the matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_set_submatrix(fm, new_values, start_row, start_col, n_rows, n_cols, alpha, beta, transpose)
sets a submatrix of a full matrix fm(start_row:start_row+n_rows,start_col:start_col+n_cols) = alpha*o...
subroutine, public cp_fm_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
integer function, public cp_fm_pilaenv(ictxt, prec)
...
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp)
creates a new full matrix with the given structure
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
subroutine, public init_input_type(input_settings, nsmax, rs_grid_section, ilevel, higher_grid_layout)
parses an input section to assign the proper values to the input type
subroutine, public dbm_run_tests(mp_group, io_unit, matrix_sizes, trs, bs_m, bs_n, bs_k, sparsities, alpha, beta, n_loops, eps, retain_sparsity, always_checksum)
Tests the DBM library.
Definition dbm_tests.F:54
subroutine, public fft_radix_operations(radix_in, radix_out, operation)
Determine the allowed lengths of FFT's '''.
Definition fft_tools.F:239
integer, parameter, public bwfft
Definition fft_tools.F:146
integer, parameter, public fft_radix_closest
Definition fft_tools.F:147
subroutine, public init_fft(fftlib, alltoall, fftsg_sizes, pool_limit, wisdom_file, plan_style)
...
Definition fft_tools.F:186
integer, parameter, public fwfft
Definition fft_tools.F:146
subroutine, public finalize_fft(para_env, wisdom_file)
does whatever is needed to finalize the current fft setup
Definition fft_tools.F:216
Define type storing the global information of a run. Keep the amount of stored data small....
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_diag_syevd
integer, parameter, public do_mat_read
integer, parameter, public do_pwgrid_ns_fullspace
integer, parameter, public do_diag_syevx
integer, parameter, public do_pwgrid_spherical
integer, parameter, public do_mat_random
integer, parameter, public do_pwgrid_ns_halfspace
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Performance tests for basic tasks like matrix multiplies, copy, fft.
subroutine, public lib_test(root_section, para_env, globenv)
Master routine for tests.
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:130
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Definition machine.F:147
Interface to the message passing library MPI.
Routines to calculate the minimax coefficients in order to approximate 1/x as a sum over exponential ...
Definition minimax_exp.F:29
subroutine, public validate_exp_minimax(n_r, iw)
Unit test checking that numerical error of minimax approximations generated using any k15 or k53 coef...
Routines to calculate frequency and time grids (integration points and weights) for correlation metho...
Definition mp2_grids.F:14
subroutine, public test_least_square_ft(nr, iw)
test the singular value decomposition for the computation of integration weights for the Fourier tran...
Definition mp2_grids.F:1167
Interface to the message passing library MPI.
subroutine, public mpi_perf_test(comm, npow, output_unit)
Tests the MPI library.
basic linear algebra operations for full matrixes
Parallel (pseudo)random number generator (RNG) for multiple streams and substreams of random numbers.
subroutine advance(self, e, c)
Advance the state by n steps, i.e. jump n steps forward, if n > 0, or backward if n < 0.
integer, parameter, public uniform
integer, parameter, public halfspace
integer, parameter, public fullspace
This module defines the grid data type and some basic operations on it.
Definition pw_grids.F:36
subroutine, public pw_grid_release(pw_grid)
releases the given pw grid
Definition pw_grids.F:2163
subroutine, public rs_grid_print(rs, iounit)
Print information on grids to output.
subroutine, public rs_grid_create(rs, desc)
...
subroutine, public rs_grid_create_descriptor(desc, pw_grid, input_settings, border_points)
Determine the setup of real space grids - this is divided up into the creation of a descriptor and th...
subroutine, public transfer_pw2rs(rs, pw)
...
subroutine, public rs_grid_release_descriptor(rs_desc)
releases the given rs grid descriptor (see doc/ReferenceCounting.html)
subroutine, public transfer_rs2pw(rs, pw)
...
subroutine, public rs_grid_release(rs_grid)
releases the given rs grid (see doc/ReferenceCounting.html)
subroutine, public rs_grid_zero(rs)
Initialize grid to zero.
Calculates 2-center integrals for different r12 operators comparing the Solid harmonic Gaussian integ...
subroutine, public shg_integrals_perf_acc_test(iw, shg_integrals_test_section)
Unit test for performance and accuracy of the SHG integrals.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
keeps the information about the structure of a full matrix
represent a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
contains the initially parsed file and the initial parallel environment
stores all the informations relevant to an mpi environment