-
Notifications
You must be signed in to change notification settings - Fork 1
/
library_tests.F
1519 lines (1339 loc) · 66.9 KB
/
library_tests.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!--------------------------------------------------------------------------------------------------!
! CP2K: A general program to perform molecular dynamics simulations !
! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
! !
! SPDX-License-Identifier: GPL-2.0-or-later !
!--------------------------------------------------------------------------------------------------!
! **************************************************************************************************
!> \brief Performance tests for basic tasks like matrix multiplies, copy, fft.
!> \par History
!> 30-Nov-2000 (JGH) added input
!> 02-Jan-2001 (JGH) Parallel FFT
!> 28-Feb-2002 (JGH) Clebsch-Gordon Coefficients
!> 06-Jun-2003 (JGH) Real space grid test
!> Eigensolver test (29.08.05,MK)
!> \author JGH 6-NOV-2000
! **************************************************************************************************
MODULE library_tests
USE ai_coulomb_test, ONLY: eri_test
USE cell_methods, ONLY: cell_create,&
init_cell
USE cell_types, ONLY: cell_release,&
cell_type
USE cg_test, ONLY: clebsch_gordon_test
USE cp_blacs_env, ONLY: cp_blacs_env_create,&
cp_blacs_env_release,&
cp_blacs_env_type
USE cp_dbcsr_api, ONLY: dbcsr_reset_randmat_seed,&
dbcsr_run_tests
USE cp_eri_mme_interface, ONLY: cp_eri_mme_perf_acc_test
USE cp_files, ONLY: close_file,&
open_file
USE cp_fm_basic_linalg, ONLY: cp_fm_gemm
USE cp_fm_diag, ONLY: cp_fm_syevd,&
cp_fm_syevx
USE cp_fm_struct, ONLY: cp_fm_struct_create,&
cp_fm_struct_get,&
cp_fm_struct_release,&
cp_fm_struct_type
USE cp_fm_types, ONLY: cp_fm_create,&
cp_fm_pilaenv,&
cp_fm_release,&
cp_fm_set_all,&
cp_fm_set_submatrix,&
cp_fm_to_fm,&
cp_fm_type
USE cp_log_handling, ONLY: cp_get_default_logger,&
cp_logger_type
USE cp_output_handling, ONLY: cp_print_key_finished_output,&
cp_print_key_unit_nr
USE cp_realspace_grid_init, ONLY: init_input_type
USE dbm_tests, ONLY: dbm_run_tests
USE fft_tools, ONLY: BWFFT,&
FFT_RADIX_CLOSEST,&
FWFFT,&
fft3d,&
fft_radix_operations,&
finalize_fft,&
init_fft
USE global_types, ONLY: global_environment_type
USE input_constants, ONLY: do_diag_syevd,&
do_diag_syevx,&
do_mat_random,&
do_mat_read,&
do_pwgrid_ns_fullspace,&
do_pwgrid_ns_halfspace,&
do_pwgrid_spherical
USE input_section_types, ONLY: section_vals_get,&
section_vals_get_subs_vals,&
section_vals_type,&
section_vals_val_get
USE kinds, ONLY: dp
USE machine, ONLY: m_flush,&
m_walltime
USE message_passing, ONLY: mp_para_env_type
USE minimax_exp, ONLY: validate_exp_minimax
USE mp2_grids, ONLY: test_least_square_ft
USE mp_perf_test, ONLY: mpi_perf_test
USE parallel_gemm_api, ONLY: parallel_gemm
USE parallel_rng_types, ONLY: UNIFORM,&
rng_stream_type
USE pw_grid_types, ONLY: FULLSPACE,&
HALFSPACE,&
pw_grid_type
USE pw_grids, ONLY: pw_grid_create,&
pw_grid_release
USE pw_methods, ONLY: pw_transfer,&
pw_zero
USE pw_types, ONLY: pw_c1d_gs_type,&
pw_c3d_rs_type,&
pw_r3d_rs_type
USE realspace_grid_types, ONLY: &
realspace_grid_desc_type, realspace_grid_input_type, realspace_grid_type, rs_grid_create, &
rs_grid_create_descriptor, rs_grid_print, rs_grid_release, rs_grid_release_descriptor, &
rs_grid_zero, transfer_pw2rs, transfer_rs2pw
USE shg_integrals_test, ONLY: shg_integrals_perf_acc_test
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
PUBLIC :: lib_test
INTEGER :: runtest(100)
REAL(KIND=dp) :: max_memory
REAL(KIND=dp), PARAMETER :: threshold = 1.0E-8_dp
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'library_tests'
CONTAINS
! **************************************************************************************************
!> \brief Master routine for tests
!> \param root_section ...
!> \param para_env ...
!> \param globenv ...
!> \par History
!> none
!> \author JGH 6-NOV-2000
! **************************************************************************************************
SUBROUTINE lib_test(root_section, para_env, globenv)
TYPE(section_vals_type), POINTER :: root_section
TYPE(mp_para_env_type), POINTER :: para_env
TYPE(global_environment_type), POINTER :: globenv
CHARACTER(LEN=*), PARAMETER :: routineN = 'lib_test'
INTEGER :: handle, iw
LOGICAL :: explicit
TYPE(cp_logger_type), POINTER :: logger
TYPE(section_vals_type), POINTER :: cp_dbcsr_test_section, cp_fm_gemm_test_section, &
dbm_test_section, eigensolver_section, eri_mme_test_section, pw_transfer_section, &
rs_pw_transfer_section, shg_integrals_test_section
CALL timeset(routineN, handle)
logger => cp_get_default_logger()
iw = cp_print_key_unit_nr(logger, root_section, "TEST%PROGRAM_RUN_INFO", extension=".log")
IF (iw > 0) THEN
WRITE (iw, '(T2,79("*"))')
WRITE (iw, '(A,T31,A,T80,A)') ' *', ' PERFORMANCE TESTS ', '*'
WRITE (iw, '(T2,79("*"))')
END IF
!
CALL test_input(root_section, para_env)
!
IF (runtest(1) /= 0) CALL copy_test(para_env, iw)
!
IF (runtest(2) /= 0) CALL matmul_test(para_env, test_matmul=.TRUE., test_dgemm=.FALSE., iw=iw)
IF (runtest(5) /= 0) CALL matmul_test(para_env, test_matmul=.FALSE., test_dgemm=.TRUE., iw=iw)
!
IF (runtest(3) /= 0) CALL fft_test(para_env, iw, globenv%fftw_plan_type, &
globenv%fftw_wisdom_file_name)
!
IF (runtest(4) /= 0) CALL eri_test(iw)
!
IF (runtest(6) /= 0) CALL clebsch_gordon_test()
!
! runtest 7 has been deleted and can be recycled
!
IF (runtest(8) /= 0) CALL mpi_perf_test(para_env, runtest(8), iw)
!
IF (runtest(10) /= 0) CALL validate_exp_minimax(runtest(10), iw)
!
IF (runtest(11) /= 0) CALL test_least_square_ft(runtest(11), iw)
!
rs_pw_transfer_section => section_vals_get_subs_vals(root_section, "TEST%RS_PW_TRANSFER")
CALL section_vals_get(rs_pw_transfer_section, explicit=explicit)
IF (explicit) THEN
CALL rs_pw_transfer_test(para_env, iw, globenv, rs_pw_transfer_section)
END IF
pw_transfer_section => section_vals_get_subs_vals(root_section, "TEST%PW_TRANSFER")
CALL section_vals_get(pw_transfer_section, explicit=explicit)
IF (explicit) THEN
CALL pw_fft_test(para_env, iw, globenv, pw_transfer_section)
END IF
cp_fm_gemm_test_section => section_vals_get_subs_vals(root_section, "TEST%CP_FM_GEMM")
CALL section_vals_get(cp_fm_gemm_test_section, explicit=explicit)
IF (explicit) THEN
CALL cp_fm_gemm_test(para_env, iw, cp_fm_gemm_test_section)
END IF
eigensolver_section => section_vals_get_subs_vals(root_section, "TEST%EIGENSOLVER")
CALL section_vals_get(eigensolver_section, explicit=explicit)
IF (explicit) THEN
CALL eigensolver_test(para_env, iw, eigensolver_section)
END IF
eri_mme_test_section => section_vals_get_subs_vals(root_section, "TEST%ERI_MME_TEST")
CALL section_vals_get(eri_mme_test_section, explicit=explicit)
IF (explicit) THEN
CALL cp_eri_mme_perf_acc_test(para_env, iw, eri_mme_test_section)
END IF
shg_integrals_test_section => section_vals_get_subs_vals(root_section, "TEST%SHG_INTEGRALS_TEST")
CALL section_vals_get(shg_integrals_test_section, explicit=explicit)
IF (explicit) THEN
CALL shg_integrals_perf_acc_test(iw, shg_integrals_test_section)
END IF
! DBCSR tests
cp_dbcsr_test_section => section_vals_get_subs_vals(root_section, "TEST%CP_DBCSR")
CALL section_vals_get(cp_dbcsr_test_section, explicit=explicit)
IF (explicit) THEN
CALL cp_dbcsr_tests(para_env, iw, cp_dbcsr_test_section)
END IF
! DBM tests
dbm_test_section => section_vals_get_subs_vals(root_section, "TEST%DBM")
CALL section_vals_get(dbm_test_section, explicit=explicit)
IF (explicit) THEN
CALL run_dbm_tests(para_env, iw, dbm_test_section)
END IF
CALL cp_print_key_finished_output(iw, logger, root_section, "TEST%PROGRAM_RUN_INFO")
CALL timestop(handle)
END SUBROUTINE lib_test
! **************************************************************************************************
!> \brief Reads input section &TEST ... &END
!> \param root_section ...
!> \param para_env ...
!> \author JGH 30-NOV-2000
!> \note
!> I---------------------------------------------------------------------------I
!> I SECTION: &TEST ... &END I
!> I I
!> I MEMORY max_memory I
!> I COPY n I
!> I MATMUL n I
!> I FFT n I
!> I ERI n I
!> I PW_FFT n I
!> I Clebsch-Gordon n I
!> I RS_GRIDS n I
!> I MPI n I
!> I RNG n -> Parallel random number generator I
!> I---------------------------------------------------------------------------I
! **************************************************************************************************
SUBROUTINE test_input(root_section, para_env)
TYPE(section_vals_type), POINTER :: root_section
TYPE(mp_para_env_type), POINTER :: para_env
TYPE(section_vals_type), POINTER :: test_section
!
!..defaults
! using this style is not recommended, introduce sections instead (see e.g. cp_fm_gemm)
runtest = 0
test_section => section_vals_get_subs_vals(root_section, "TEST")
CALL section_vals_val_get(test_section, "MEMORY", r_val=max_memory)
CALL section_vals_val_get(test_section, 'COPY', i_val=runtest(1))
CALL section_vals_val_get(test_section, 'MATMUL', i_val=runtest(2))
CALL section_vals_val_get(test_section, 'DGEMM', i_val=runtest(5))
CALL section_vals_val_get(test_section, 'FFT', i_val=runtest(3))
CALL section_vals_val_get(test_section, 'ERI', i_val=runtest(4))
CALL section_vals_val_get(test_section, 'CLEBSCH_GORDON', i_val=runtest(6))
CALL section_vals_val_get(test_section, 'MPI', i_val=runtest(8))
CALL section_vals_val_get(test_section, 'MINIMAX', i_val=runtest(10))
CALL section_vals_val_get(test_section, 'LEAST_SQ_FT', i_val=runtest(11))
CALL para_env%sync()
END SUBROUTINE test_input
! **************************************************************************************************
!> \brief Tests the performance to copy two vectors.
!> \param para_env ...
!> \param iw ...
!> \par History
!> none
!> \author JGH 6-NOV-2000
!> \note
!> The results of these tests allow to determine the size of the cache
!> of the CPU. This can be used to optimize the performance of the
!> FFTSG library.
! **************************************************************************************************
SUBROUTINE copy_test(para_env, iw)
TYPE(mp_para_env_type), POINTER :: para_env
INTEGER :: iw
INTEGER :: i, ierr, j, len, ntim, siz
REAL(KIND=dp) :: perf, t, tend, tstart
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: ca, cb
! test for copy --> Cache size
siz = ABS(runtest(1))
IF (para_env%is_source()) WRITE (iw, '(//,A,/)') " Test of copy ( F95 ) "
DO i = 6, 24
len = 2**i
IF (8.0_dp*REAL(len, KIND=dp) > max_memory*0.5_dp) EXIT
ALLOCATE (ca(len), STAT=ierr)
IF (ierr /= 0) EXIT
ALLOCATE (cb(len), STAT=ierr)
IF (ierr /= 0) EXIT
CALL RANDOM_NUMBER(ca)
ntim = NINT(1.e7_dp/REAL(len, KIND=dp))
ntim = MAX(ntim, 1)
ntim = MIN(ntim, siz*10000)
tstart = m_walltime()
DO j = 1, ntim
cb(:) = ca(:)
ca(1) = REAL(j, KIND=dp)
END DO
tend = m_walltime()
t = tend - tstart + threshold
IF (t > 0.0_dp) THEN
perf = REAL(ntim, KIND=dp)*REAL(len, KIND=dp)*1.e-6_dp/t
ELSE
perf = 0.0_dp
END IF
IF (para_env%is_source()) THEN
WRITE (iw, '(A,i2,i10,A,T59,F14.4,A)') " Copy test: Size = 2^", i, &
len/1024, " Kwords", perf, " Mcopy/s"
END IF
DEALLOCATE (ca)
DEALLOCATE (cb)
END DO
CALL para_env%sync()
END SUBROUTINE copy_test
! **************************************************************************************************
!> \brief Tests the performance of different kinds of matrix matrix multiply
!> kernels for the BLAS and F95 intrinsic matmul.
!> \param para_env ...
!> \param test_matmul ...
!> \param test_dgemm ...
!> \param iw ...
!> \par History
!> none
!> \author JGH 6-NOV-2000
! **************************************************************************************************
SUBROUTINE matmul_test(para_env, test_matmul, test_dgemm, iw)
TYPE(mp_para_env_type), POINTER :: para_env
LOGICAL :: test_matmul, test_dgemm
INTEGER :: iw
INTEGER :: i, ierr, j, len, ntim, siz
REAL(KIND=dp) :: perf, t, tend, tstart, xdum
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: ma, mb, mc
! test for matrix multpies
IF (test_matmul) THEN
siz = ABS(runtest(2))
IF (para_env%is_source()) WRITE (iw, '(//,A,/)') " Test of matmul ( F95 ) "
DO i = 5, siz, 2
len = 2**i + 1
IF (8.0_dp*REAL(len*len, KIND=dp) > max_memory*0.3_dp) EXIT
ALLOCATE (ma(len, len), STAT=ierr)
IF (ierr /= 0) EXIT
ALLOCATE (mb(len, len), STAT=ierr)
IF (ierr /= 0) EXIT
ALLOCATE (mc(len, len), STAT=ierr)
IF (ierr /= 0) EXIT
mc = 0.0_dp
CALL RANDOM_NUMBER(xdum)
ma = xdum
CALL RANDOM_NUMBER(xdum)
mb = xdum
ntim = NINT(1.e8_dp/(2.0_dp*REAL(len, KIND=dp)**3))
ntim = MAX(ntim, 1)
ntim = MIN(ntim, siz*200)
tstart = m_walltime()
DO j = 1, ntim
mc(:, :) = MATMUL(ma, mb)
ma(1, 1) = REAL(j, KIND=dp)
END DO
tend = m_walltime()
t = tend - tstart + threshold
perf = REAL(ntim, KIND=dp)*2.0_dp*REAL(len, KIND=dp)**3*1.e-6_dp/t
IF (para_env%is_source()) THEN
WRITE (iw, '(A,i6,T59,F14.4,A)') &
" Matrix multiply test: c = a * b Size = ", len, perf, " Mflop/s"
END IF
tstart = m_walltime()
DO j = 1, ntim
mc(:, :) = mc + MATMUL(ma, mb)
ma(1, 1) = REAL(j, KIND=dp)
END DO
tend = m_walltime()
t = tend - tstart + threshold
IF (t > 0.0_dp) THEN
perf = REAL(ntim, KIND=dp)*2.0_dp*REAL(len, KIND=dp)**3*1.e-6_dp/t
ELSE
perf = 0.0_dp
END IF
IF (para_env%is_source()) THEN
WRITE (iw, '(A,i6,T59,F14.4,A)') &
" Matrix multiply test: a = a * b Size = ", len, perf, " Mflop/s"
END IF
tstart = m_walltime()
DO j = 1, ntim
mc(:, :) = mc + MATMUL(ma, TRANSPOSE(mb))
ma(1, 1) = REAL(j, KIND=dp)
END DO
tend = m_walltime()
t = tend - tstart + threshold
IF (t > 0.0_dp) THEN
perf = REAL(ntim, KIND=dp)*2.0_dp*REAL(len, KIND=dp)**3*1.e-6_dp/t
ELSE
perf = 0.0_dp
END IF
IF (para_env%is_source()) THEN
WRITE (iw, '(A,i6,T59,F14.4,A)') &
" Matrix multiply test: c = a * b(T) Size = ", len, perf, " Mflop/s"
END IF
tstart = m_walltime()
DO j = 1, ntim
mc(:, :) = mc + MATMUL(TRANSPOSE(ma), mb)
ma(1, 1) = REAL(j, KIND=dp)
END DO
tend = m_walltime()
t = tend - tstart + threshold
IF (t > 0.0_dp) THEN
perf = REAL(ntim, KIND=dp)*2.0_dp*REAL(len, KIND=dp)**3*1.e-6_dp/t
ELSE
perf = 0.0_dp
END IF
IF (para_env%is_source()) THEN
WRITE (iw, '(A,i6,T59,F14.4,A)') &
" Matrix multiply test: c = a(T) * b Size = ", len, perf, " Mflop/s"
END IF
DEALLOCATE (ma)
DEALLOCATE (mb)
DEALLOCATE (mc)
END DO
END IF
! test for matrix multpies
IF (test_dgemm) THEN
siz = ABS(runtest(5))
IF (para_env%is_source()) WRITE (iw, '(//,A,/)') " Test of matmul ( BLAS ) "
DO i = 5, siz, 2
len = 2**i + 1
IF (8.0_dp*REAL(len*len, KIND=dp) > max_memory*0.3_dp) EXIT
ALLOCATE (ma(len, len), STAT=ierr)
IF (ierr /= 0) EXIT
ALLOCATE (mb(len, len), STAT=ierr)
IF (ierr /= 0) EXIT
ALLOCATE (mc(len, len), STAT=ierr)
IF (ierr /= 0) EXIT
mc = 0.0_dp
CALL RANDOM_NUMBER(xdum)
ma = xdum
CALL RANDOM_NUMBER(xdum)
mb = xdum
ntim = NINT(1.e8_dp/(2.0_dp*REAL(len, KIND=dp)**3))
ntim = MAX(ntim, 1)
ntim = MIN(ntim, 1000)
tstart = m_walltime()
DO j = 1, ntim
CALL dgemm("N", "N", len, len, len, 1.0_dp, ma, len, mb, len, 1.0_dp, mc, len)
END DO
tend = m_walltime()
t = tend - tstart + threshold
IF (t > 0.0_dp) THEN
perf = REAL(ntim, KIND=dp)*2.0_dp*REAL(len, KIND=dp)**3*1.e-6_dp/t
ELSE
perf = 0.0_dp
END IF
IF (para_env%is_source()) THEN
WRITE (iw, '(A,i6,T59,F14.4,A)') &
" Matrix multiply test: c = a * b Size = ", len, perf, " Mflop/s"
END IF
tstart = m_walltime()
DO j = 1, ntim
CALL dgemm("N", "N", len, len, len, 1.0_dp, ma, len, mb, len, 1.0_dp, mc, len)
END DO
tend = m_walltime()
t = tend - tstart + threshold
IF (t > 0.0_dp) THEN
perf = REAL(ntim, KIND=dp)*2.0_dp*REAL(len, KIND=dp)**3*1.e-6_dp/t
ELSE
perf = 0.0_dp
END IF
IF (para_env%is_source()) THEN
WRITE (iw, '(A,i6,T59,F14.4,A)') &
" Matrix multiply test: a = a * b Size = ", len, perf, " Mflop/s"
END IF
tstart = m_walltime()
DO j = 1, ntim
CALL dgemm("N", "T", len, len, len, 1.0_dp, ma, len, mb, len, 1.0_dp, mc, len)
END DO
tend = m_walltime()
t = tend - tstart + threshold
IF (t > 0.0_dp) THEN
perf = REAL(ntim, KIND=dp)*2.0_dp*REAL(len, KIND=dp)**3*1.e-6_dp/t
ELSE
perf = 0.0_dp
END IF
IF (para_env%is_source()) THEN
WRITE (iw, '(A,i6,T59,F14.4,A)') &
" Matrix multiply test: c = a * b(T) Size = ", len, perf, " Mflop/s"
END IF
tstart = m_walltime()
DO j = 1, ntim
CALL dgemm("T", "N", len, len, len, 1.0_dp, ma, len, mb, len, 1.0_dp, mc, len)
END DO
tend = m_walltime()
t = tend - tstart + threshold
IF (t > 0.0_dp) THEN
perf = REAL(ntim, KIND=dp)*2.0_dp*REAL(len, KIND=dp)**3*1.e-6_dp/t
ELSE
perf = 0.0_dp
END IF
IF (para_env%is_source()) THEN
WRITE (iw, '(A,i6,T59,F14.4,A)') &
" Matrix multiply test: c = a(T) * b Size = ", len, perf, " Mflop/s"
END IF
DEALLOCATE (ma)
DEALLOCATE (mb)
DEALLOCATE (mc)
END DO
END IF
CALL para_env%sync()
END SUBROUTINE matmul_test
! **************************************************************************************************
!> \brief Tests the performance of all available FFT libraries for 3D FFTs
!> \param para_env ...
!> \param iw ...
!> \param fftw_plan_type ...
!> \param wisdom_file where FFTW3 should look to save/load wisdom
!> \par History
!> none
!> \author JGH 6-NOV-2000
! **************************************************************************************************
SUBROUTINE fft_test(para_env, iw, fftw_plan_type, wisdom_file)
TYPE(mp_para_env_type), POINTER :: para_env
INTEGER :: iw, fftw_plan_type
CHARACTER(LEN=*), INTENT(IN) :: wisdom_file
INTEGER, PARAMETER :: ndate(3) = (/12, 48, 96/)
INTEGER :: iall, ierr, it, j, len, n(3), ntim, &
radix_in, radix_out, siz, stat
COMPLEX(KIND=dp), DIMENSION(4, 4, 4) :: zz
COMPLEX(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: ca, cb, cc
CHARACTER(LEN=7) :: method
REAL(KIND=dp) :: flops, perf, scale, t, tdiff, tend, &
tstart
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: ra
! test for 3d FFT
IF (para_env%is_source()) WRITE (iw, '(//,A,/)') " Test of 3D-FFT "
siz = ABS(runtest(3))
DO iall = 1, 100
SELECT CASE (iall)
CASE DEFAULT
EXIT
CASE (1)
CALL init_fft("FFTSG", alltoall=.FALSE., fftsg_sizes=.TRUE., wisdom_file=wisdom_file, &
pool_limit=10, plan_style=fftw_plan_type)
method = "FFTSG "
CASE (2)
CYCLE
CASE (3)
CALL init_fft("FFTW3", alltoall=.FALSE., fftsg_sizes=.TRUE., wisdom_file=wisdom_file, &
pool_limit=10, plan_style=fftw_plan_type)
method = "FFTW3 "
END SELECT
n = 4
zz = 0.0_dp
CALL fft3d(FWFFT, n, zz, status=stat)
IF (stat == 0) THEN
DO it = 1, 3
radix_in = ndate(it)
CALL fft_radix_operations(radix_in, radix_out, FFT_RADIX_CLOSEST)
len = radix_out
n = len
IF (16.0_dp*REAL(len*len*len, KIND=dp) > max_memory*0.5_dp) EXIT
ALLOCATE (ra(len, len, len), STAT=ierr)
ALLOCATE (ca(len, len, len), STAT=ierr)
CALL RANDOM_NUMBER(ra)
ca(:, :, :) = ra
CALL RANDOM_NUMBER(ra)
ca(:, :, :) = ca + CMPLX(0.0_dp, 1.0_dp, KIND=dp)*ra
flops = REAL(len**3, KIND=dp)*15.0_dp*LOG(REAL(len, KIND=dp))
ntim = NINT(siz*1.e7_dp/flops)
ntim = MAX(ntim, 1)
ntim = MIN(ntim, 200)
scale = 1.0_dp/REAL(len**3, KIND=dp)
tstart = m_walltime()
DO j = 1, ntim
CALL fft3d(FWFFT, n, ca)
CALL fft3d(BWFFT, n, ca)
END DO
tend = m_walltime()
t = tend - tstart + threshold
IF (t > 0.0_dp) THEN
perf = REAL(ntim, KIND=dp)*2.0_dp*flops*1.e-6_dp/t
ELSE
perf = 0.0_dp
END IF
IF (para_env%is_source()) THEN
WRITE (iw, '(T2,A,A,i6,T59,F14.4,A)') &
ADJUSTR(method), " test (in-place) Size = ", len, perf, " Mflop/s"
END IF
DEALLOCATE (ca)
DEALLOCATE (ra)
END DO
IF (para_env%is_source()) WRITE (iw, *)
! test if input data is preserved
len = 24
n = len
ALLOCATE (ra(len, len, len))
ALLOCATE (ca(len, len, len))
ALLOCATE (cb(len, len, len))
ALLOCATE (cc(len, len, len))
CALL RANDOM_NUMBER(ra)
ca(:, :, :) = ra
CALL RANDOM_NUMBER(ra)
ca(:, :, :) = ca + CMPLX(0.0_dp, 1.0_dp, KIND=dp)*ra
cc(:, :, :) = ca
CALL fft3d(FWFFT, n, ca, cb)
tdiff = MAXVAL(ABS(ca - cc))
IF (tdiff > 1.0E-12_dp) THEN
IF (para_env%is_source()) &
WRITE (iw, '(T2,A,A,A)') ADJUSTR(method), " FWFFT ", &
" Input array is changed in out-of-place FFT !"
ELSE
IF (para_env%is_source()) &
WRITE (iw, '(T2,A,A,A)') ADJUSTR(method), " FWFFT ", &
" Input array is not changed in out-of-place FFT !"
END IF
ca(:, :, :) = cc
CALL fft3d(BWFFT, n, ca, cb)
tdiff = MAXVAL(ABS(ca - cc))
IF (tdiff > 1.0E-12_dp) THEN
IF (para_env%is_source()) &
WRITE (iw, '(T2,A,A,A)') ADJUSTR(method), " BWFFT ", &
" Input array is changed in out-of-place FFT !"
ELSE
IF (para_env%is_source()) &
WRITE (iw, '(T2,A,A,A)') ADJUSTR(method), " BWFFT ", &
" Input array is not changed in out-of-place FFT !"
END IF
IF (para_env%is_source()) WRITE (iw, *)
DEALLOCATE (ra)
DEALLOCATE (ca)
DEALLOCATE (cb)
DEALLOCATE (cc)
END IF
CALL finalize_fft(para_env, wisdom_file=wisdom_file)
END DO
END SUBROUTINE fft_test
! **************************************************************************************************
!> \brief test rs_pw_transfer performance
!> \param para_env ...
!> \param iw ...
!> \param globenv ...
!> \param rs_pw_transfer_section ...
!> \author Joost VandeVondele
!> 9.2008 Randomise rs grid [Iain Bethune]
!> (c) The Numerical Algorithms Group (NAG) Ltd, 2008 on behalf of the HECToR project
! **************************************************************************************************
SUBROUTINE rs_pw_transfer_test(para_env, iw, globenv, rs_pw_transfer_section)
TYPE(mp_para_env_type), POINTER :: para_env
INTEGER :: iw
TYPE(global_environment_type), POINTER :: globenv
TYPE(section_vals_type), POINTER :: rs_pw_transfer_section
CHARACTER(LEN=*), PARAMETER :: routineN = 'rs_pw_transfer_test'
INTEGER :: halo_size, handle, i_loop, n_loop, ns_max
INTEGER, DIMENSION(3) :: no, np
INTEGER, DIMENSION(:), POINTER :: i_vals
LOGICAL :: do_rs2pw
REAL(KIND=dp) :: tend, tstart
TYPE(cell_type), POINTER :: box
TYPE(pw_grid_type), POINTER :: grid
TYPE(pw_r3d_rs_type) :: ca
TYPE(realspace_grid_desc_type), POINTER :: rs_desc
TYPE(realspace_grid_input_type) :: input_settings
TYPE(realspace_grid_type) :: rs_grid
TYPE(section_vals_type), POINTER :: rs_grid_section
CALL timeset(routineN, handle)
!..set fft lib
CALL init_fft(globenv%default_fft_library, alltoall=.FALSE., fftsg_sizes=.TRUE., &
pool_limit=globenv%fft_pool_scratch_limit, &
wisdom_file=globenv%fftw_wisdom_file_name, &
plan_style=globenv%fftw_plan_type)
! .. set cell (should otherwise be irrelevant)
NULLIFY (box)
CALL cell_create(box)
box%hmat = RESHAPE((/20.0_dp, 0.0_dp, 0.0_dp, 0.0_dp, 20.0_dp, 0.0_dp, &
0.0_dp, 0.0_dp, 20.0_dp/), (/3, 3/))
CALL init_cell(box)
! .. grid type and pw_grid
NULLIFY (grid)
CALL section_vals_val_get(rs_pw_transfer_section, "GRID", i_vals=i_vals)
np = i_vals
CALL pw_grid_create(grid, para_env, box%hmat, grid_span=FULLSPACE, npts=np, fft_usage=.TRUE., iounit=iw)
no = grid%npts
CALL ca%create(grid)
CALL pw_zero(ca)
! .. rs input setting type
CALL section_vals_val_get(rs_pw_transfer_section, "HALO_SIZE", i_val=halo_size)
rs_grid_section => section_vals_get_subs_vals(rs_pw_transfer_section, "RS_GRID")
ns_max = 2*halo_size + 1
CALL init_input_type(input_settings, ns_max, rs_grid_section, 1, (/-1, -1, -1/))
! .. rs type
NULLIFY (rs_desc)
CALL rs_grid_create_descriptor(rs_desc, pw_grid=grid, input_settings=input_settings)
CALL rs_grid_create(rs_grid, rs_desc)
CALL rs_grid_print(rs_grid, iw)
CALL rs_grid_zero(rs_grid)
! Put random values on the grid, so summation check will pick up errors
CALL RANDOM_NUMBER(rs_grid%r)
CALL section_vals_val_get(rs_pw_transfer_section, "N_loop", i_val=N_loop)
CALL section_vals_val_get(rs_pw_transfer_section, "RS2PW", l_val=do_rs2pw)
! go for the real loops, sync to get max timings
IF (para_env%is_source()) THEN
WRITE (iw, '(T2,A)') ""
WRITE (iw, '(T2,A)') "Timing rs_pw_transfer routine"
WRITE (iw, '(T2,A)') ""
WRITE (iw, '(T2,A)') "iteration time[s]"
END IF
DO i_loop = 1, N_loop
CALL para_env%sync()
tstart = m_walltime()
IF (do_rs2pw) THEN
CALL transfer_rs2pw(rs_grid, ca)
ELSE
CALL transfer_pw2rs(rs_grid, ca)
END IF
CALL para_env%sync()
tend = m_walltime()
IF (para_env%is_source()) THEN
WRITE (iw, '(T2,I9,1X,F12.6)') i_loop, tend - tstart
END IF
END DO
!cleanup
CALL rs_grid_release(rs_grid)
CALL rs_grid_release_descriptor(rs_desc)
CALL ca%release()
CALL pw_grid_release(grid)
CALL cell_release(box)
CALL finalize_fft(para_env, wisdom_file=globenv%fftw_wisdom_file_name)
CALL timestop(handle)
END SUBROUTINE rs_pw_transfer_test
! **************************************************************************************************
!> \brief Tests the performance of PW calls to FFT routines
!> \param para_env ...
!> \param iw ...
!> \param globenv ...
!> \param pw_transfer_section ...
!> \par History
!> JGH 6-Feb-2001 : Test and performance code
!> Made input sensitive [Joost VandeVondele]
!> \author JGH 1-JAN-2001
! **************************************************************************************************
SUBROUTINE pw_fft_test(para_env, iw, globenv, pw_transfer_section)
TYPE(mp_para_env_type), POINTER :: para_env
INTEGER :: iw
TYPE(global_environment_type), POINTER :: globenv
TYPE(section_vals_type), POINTER :: pw_transfer_section
REAL(KIND=dp), PARAMETER :: toler = 1.e-11_dp
INTEGER :: blocked_id, grid_span, i_layout, i_rep, &
ig, ip, itmp, n_loop, n_rep, nn, p, q
INTEGER, ALLOCATABLE, DIMENSION(:, :) :: layouts
INTEGER, DIMENSION(2) :: distribution_layout
INTEGER, DIMENSION(3) :: no, np
INTEGER, DIMENSION(:), POINTER :: i_vals
LOGICAL :: debug, is_fullspace, odd, &
pw_grid_layout_all, spherical
REAL(KIND=dp) :: em, et, flops, gsq, perf, t, t_max, &
t_min, tend, tstart
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: t_end, t_start
TYPE(cell_type), POINTER :: box
TYPE(pw_c1d_gs_type) :: ca, cc
TYPE(pw_c3d_rs_type) :: cb
TYPE(pw_grid_type), POINTER :: grid
!..set fft lib
CALL init_fft(globenv%default_fft_library, alltoall=.FALSE., fftsg_sizes=.TRUE., &
pool_limit=globenv%fft_pool_scratch_limit, &
wisdom_file=globenv%fftw_wisdom_file_name, &
plan_style=globenv%fftw_plan_type)
!..the unit cell (should not really matter, the number of grid points do)
NULLIFY (box, grid)
CALL cell_create(box)
box%hmat = RESHAPE((/10.0_dp, 0.0_dp, 0.0_dp, 0.0_dp, 8.0_dp, 0.0_dp, &
0.0_dp, 0.0_dp, 7.0_dp/), (/3, 3/))
CALL init_cell(box)
CALL section_vals_get(pw_transfer_section, n_repetition=n_rep)
DO i_rep = 1, n_rep
! how often should we do the transfer
CALL section_vals_val_get(pw_transfer_section, "N_loop", i_rep_section=i_rep, i_val=N_loop)
ALLOCATE (t_start(N_loop))
ALLOCATE (t_end(N_loop))
! setup of the grids
CALL section_vals_val_get(pw_transfer_section, "GRID", i_rep_section=i_rep, i_vals=i_vals)
np = i_vals
CALL section_vals_val_get(pw_transfer_section, "PW_GRID_BLOCKED", i_rep_section=i_rep, i_val=blocked_id)
CALL section_vals_val_get(pw_transfer_section, "DEBUG", i_rep_section=i_rep, l_val=debug)
CALL section_vals_val_get(pw_transfer_section, "PW_GRID_LAYOUT_ALL", i_rep_section=i_rep, &
l_val=pw_grid_layout_all)
! prepare to loop over all or a specific layout
IF (pw_grid_layout_all) THEN
! count layouts that fit
itmp = 0
! start from 2, (/1,para_env%num_pe/) is not supported
DO p = 2, para_env%num_pe
q = para_env%num_pe/p
IF (p*q == para_env%num_pe) THEN
itmp = itmp + 1
END IF
END DO
! build list
ALLOCATE (layouts(2, itmp))
itmp = 0
DO p = 2, para_env%num_pe
q = para_env%num_pe/p
IF (p*q == para_env%num_pe) THEN
itmp = itmp + 1
layouts(:, itmp) = (/p, q/)
END IF
END DO
ELSE
CALL section_vals_val_get(pw_transfer_section, "PW_GRID_LAYOUT", i_rep_section=i_rep, i_vals=i_vals)
ALLOCATE (layouts(2, 1))
layouts(:, 1) = i_vals
END IF
DO i_layout = 1, SIZE(layouts, 2)
distribution_layout = layouts(:, i_layout)
CALL section_vals_val_get(pw_transfer_section, "PW_GRID", i_rep_section=i_rep, i_val=itmp)
! from cp_control_utils
SELECT CASE (itmp)
CASE (do_pwgrid_spherical)
spherical = .TRUE.
is_fullspace = .FALSE.
CASE (do_pwgrid_ns_fullspace)
spherical = .FALSE.
is_fullspace = .TRUE.
CASE (do_pwgrid_ns_halfspace)
spherical = .FALSE.
is_fullspace = .FALSE.
END SELECT
! from pw_env_methods
IF (spherical) THEN
grid_span = HALFSPACE
spherical = .TRUE.
odd = .TRUE.
ELSE IF (is_fullspace) THEN
grid_span = FULLSPACE
spherical = .FALSE.
odd = .FALSE.
ELSE
grid_span = HALFSPACE
spherical = .FALSE.
odd = .TRUE.
END IF
! actual setup
CALL pw_grid_create(grid, para_env, box%hmat, grid_span=grid_span, odd=odd, spherical=spherical, &
blocked=blocked_id, npts=np, fft_usage=.TRUE., &
rs_dims=distribution_layout, iounit=iw)
IF (iw > 0) CALL m_flush(iw)
! note that the number of grid points might be different from what the user requested (fft-able needed)
no = grid%npts
CALL ca%create(grid)
CALL cb%create(grid)
CALL cc%create(grid)
! initialize data
CALL pw_zero(ca)
CALL pw_zero(cb)
CALL pw_zero(cc)
nn = SIZE(ca%array)
DO ig = 1, nn
gsq = grid%gsq(ig)
ca%array(ig) = EXP(-gsq)
END DO
flops = PRODUCT(no)*30.0_dp*LOG(REAL(MAXVAL(no), KIND=dp))
tstart = m_walltime()
DO ip = 1, n_loop
CALL para_env%sync()
t_start(ip) = m_walltime()
CALL pw_transfer(ca, cb, debug)
CALL pw_transfer(cb, cc, debug)
CALL para_env%sync()
t_end(ip) = m_walltime()
END DO
tend = m_walltime()
t = tend - tstart + threshold
IF (t > 0.0_dp) THEN
perf = REAL(n_loop, KIND=dp)*2.0_dp*flops*1.e-6_dp/t
ELSE
perf = 0.0_dp
END IF
em = MAXVAL(ABS(ca%array(:) - cc%array(:)))
CALL para_env%max(em)
et = SUM(ABS(ca%array(:) - cc%array(:)))
CALL para_env%sum(et)
t_min = MINVAL(t_end - t_start)
t_max = MAXVAL(t_end - t_start)
IF (para_env%is_source()) THEN
WRITE (iw, *)
WRITE (iw, '(A,T67,E14.6)') " Parallel FFT Tests: Maximal Error ", em
WRITE (iw, '(A,T67,E14.6)') " Parallel FFT Tests: Total Error ", et
WRITE (iw, '(A,T67,F14.0)') &
" Parallel FFT Tests: Performance [Mflops] ", perf
WRITE (iw, '(A,T67,F14.6)') " Best time : ", t_min
WRITE (iw, '(A,T67,F14.6)') " Worst time: ", t_max
IF (iw > 0) CALL m_flush(iw)
END IF
! need debugging ???
IF (em > toler .OR. et > toler) THEN
CPWARN("The FFT results are not accurate ... starting debug pw_transfer")
CALL pw_transfer(ca, cb, .TRUE.)
CALL pw_transfer(cb, cc, .TRUE.)
END IF
! done with these grids
CALL ca%release()
CALL cb%release()
CALL cc%release()
CALL pw_grid_release(grid)
END DO
! local arrays