RMPCDMD
test_correlator_0.f90
Go to the documentation of this file.
1 ! This file is part of RMPCDMD
2 ! Copyright (c) 2016 Pierre de Buyl and contributors
3 ! License: BSD 3-clause (see file LICENSE)
4 
6  use correlator
7  use tester
8  implicit none
9 
10  type(tester_t) :: test
11 
12  double precision, parameter :: dt = 0.2
13  double precision, parameter :: gamma = atan(1.d0)*4
14  double precision, parameter :: v0 = 0.03
15  integer, parameter :: l = 8
16  integer, parameter :: b = 4
17  integer, parameter :: N=8*l**(b+1)
18  type(correlator_t) :: acf1, acf2, acf3, acf4
19 
20  integer :: i, j, i_block
21  double precision :: t
22 
23  double precision :: w
24  double precision :: x(3)
25  double precision :: y(4)
26  double precision :: z(2, 5)
27 
28  double precision :: comp(l)
29 
30  call test%init(tolerance64=1d-8)
31 
32  call acf1%init(l, b)
33  call acf2%init(l, b, dim=size(x))
34  call acf3%init(l, b, n=size(y))
35  call acf4%init(l, b, dim=size(z, dim=1), n=size(z, dim=2))
36 
37  do i = 0, n-1
38  t = i*dt
39  w = t*v0
40  x = [-1, 0, 1]
41  y = [-1, 0, 1, 2]
42  do j = 1, size(z, dim=2)
43  z(1, j) = cos(gamma*i)
44  z(2, j) = 0
45  end do
46  call acf1%add(i, correlate_block_distsq, x=w)
47  call acf2%add(i, correlate_block_dot, xvec=x)
48  call acf3%add(i, correlate_block_dot, x_n=y)
49  call acf4%add(i, correlate_block_dot, xvec_n=z)
50  end do
51 
52  do i_block = 1, b
53  comp = [ ((i*dt*v0*l**(i_block-1))**2, i=0, l-1) ]
54  call test%assert_close(acf1%correlation(1,1,:,i_block) / acf1%count(i_block), comp)
55  end do
56 
57  do i_block = 1, b
58  comp = 1
59  call test%assert_close(acf2%correlation(1,1,:,i_block) / acf2%count(i_block), comp)
60  comp = 0
61  call test%assert_close(acf2%correlation(2,1,:,i_block) / acf2%count(i_block), comp)
62  comp = 1
63  call test%assert_close(acf2%correlation(3,1,:,i_block) / acf2%count(i_block), comp)
64  end do
65 
66  do i_block = 1, b
67  comp = 1
68  call test%assert_close(acf2%correlation(1,1,:,i_block) / acf2%count(i_block), comp)
69  comp = 0
70  call test%assert_close(acf2%correlation(2,1,:,i_block) / acf2%count(i_block), comp)
71  comp = 1
72  call test%assert_close(acf2%correlation(3,1,:,i_block) / acf2%count(i_block), comp)
73  end do
74 
75  do i_block = 1, b
76  comp = 1
77  call test%assert_close(acf3%correlation(1,1,:,i_block) / acf3%count(i_block), comp)
78  comp = 0
79  call test%assert_close(acf3%correlation(1,2,:,i_block) / acf3%count(i_block), comp)
80  comp = 1
81  call test%assert_close(acf3%correlation(1,3,:,i_block) / acf3%count(i_block), comp)
82  comp = 4
83  call test%assert_close(acf3%correlation(1,4,:,i_block) / acf3%count(i_block), comp)
84  end do
85 
86  do i_block = 1, b
87  if (i_block==1) then
88  comp = [ ((-1)**(modulo(i+2, 2)), i=0, l-1) ]
89  else
90  comp = 1
91  end if
92  call test%assert_close(acf4%correlation(1,1,:,i_block) / acf4%count(i_block), comp)
93  write(15,*) acf4%correlation(1,1,:,i_block) / acf4%count(i_block)
94  end do
95 
96  call test%print()
97 
98 end program test_correlator_0
subroutine, public correlate_block_dot(s, c, idx, l, n, dim)
Correlation routine for dot product.
Definition: correlator.f90:148
Block correlators.
Definition: correlator.f90:16
program test_correlator_0
subroutine, public correlate_block_distsq(s, c, idx, l, n, dim)
Correlation routine for squared distance.
Definition: correlator.f90:168