RMPCDMD
setup_single_colloid.f90
Go to the documentation of this file.
1 ! This file is part of RMPCDMD
2 ! Copyright (c) 2015-2016 Pierre de Buyl and contributors
3 ! License: BSD 3-clause (see file LICENSE)
4 
6  use common
7  use cell_system
9  use hilbert
10  use neighbor_list
11  use hdf5
12  use h5md_module
13  use interaction
14  use threefry_module
15  use mpcd
16  use iso_c_binding
17  use omp_lib
18  implicit none
19 
20  type(cell_system_t) :: solvent_cells
21  type(particle_system_t) :: solvent
22  type(particle_system_t) :: colloids
23  type(neighbor_list_t) :: neigh
24  type(lj_params_t) :: solvent_colloid_lj
25  type(lj_params_t) :: colloid_lj
26 
27  integer, parameter :: rho = 10
28  integer :: N
29  integer :: error
30 
31  double precision :: epsilon, sigma, sigma_cut
32  double precision :: mass
33  double precision :: so_max, co_max
34 
35  double precision :: e1, e2
36  double precision :: tau, dt
37  integer :: N_MD_steps
38 
39  integer :: i, L(3)
40  integer :: jump(3)
41  integer :: j, k
42  type(threefry_rng_t), allocatable :: state(:)
43  integer :: n_threads
44 
45  n_threads = omp_get_max_threads()
46  allocate(state(n_threads))
47  call threefry_rng_init(state, 4811119176716995551_c_int64_t)
48 
49  call h5open_f(error)
50 
51  l = [20, 20, 20]
52  n = rho *l(1)*l(2)*l(3)
53 
54  epsilon = 1
55  sigma = 1
56  sigma_cut = sigma*2**(1.d0/6.d0)
57 
58  call solvent_colloid_lj% init( reshape( [ epsilon ], [1, 1] ), &
59  reshape( [ sigma ], [1, 1] ), reshape( [ sigma_cut ], [1, 1] ) )
60 
61  epsilon = 1
62  sigma = 1
63  sigma_cut = sigma*2**(1.d0/6.d0)
64 
65  call colloid_lj% init( reshape( [ epsilon ], [1, 1] ), &
66  reshape( [ sigma ], [1, 1] ), reshape( [ sigma_cut ], [1, 1] ) )
67 
68  mass = rho * sigma**3 * 4 * 3.141/3
69  write(*,*) 'mass =', mass
70 
71  call solvent% init(n)
72 
73  call colloids% init(1)
74 
75  open(15,file ='data56.txt')
76 
77  write(*, *) colloids% pos
78  colloids% species = 1
79  colloids% vel = 0
80 
81  call random_number(solvent% vel(:, :))
82  solvent% vel = (solvent% vel - 0.5d0)*sqrt(6.d0*2)
83  solvent% vel = solvent% vel - spread(sum(solvent% vel, dim=2)/solvent% Nmax, 2, solvent% Nmax)
84  solvent% force = 0
85  solvent% species = 1
86 
87  call solvent_cells%init(l, 1.d0)
88  colloids% pos(:,1) = solvent_cells% edges/2.0
89 
90  call solvent% random_placement(solvent_cells% edges, colloids, solvent_colloid_lj, state(1))
91 
92  call solvent% sort(solvent_cells)
93 
94  call neigh% init(colloids% Nmax, 300)
95  call neigh% make_stencil(solvent_cells, 1.5d0)
96 
97  call neigh% update_list(colloids, solvent, 1.5d0, solvent_cells)
98 
99  tau = 0.1d0
100  n_md_steps = 100
101  dt = tau / n_md_steps
102 
103  e1 = compute_force(colloids, solvent, neigh, solvent_cells% edges, solvent_colloid_lj)
104  e2 = compute_force_n2(colloids, solvent_cells% edges, colloid_lj)
105  solvent% force_old = solvent% force
106  colloids% force_old = colloids% force
107 
108  do i = 1, 100
109  md_equil: do j = 1, n_md_steps
110  solvent% pos_old = solvent% pos + dt * solvent% vel + dt**2 * solvent% force / 2
111  !$omp parallel do private(k)
112  do k = 1, solvent% Nmax
113  solvent% pos(:,k) = modulo(solvent% pos_old(:,k), solvent_cells% edges)
114  end do
115  colloids% pos_old = colloids% pos + dt * colloids% vel + dt**2 * colloids% force / (2 * mass)
116  do k = 1, colloids% Nmax
117  jump = floor(colloids% pos_old(:,k) / solvent_cells% edges)
118  colloids% image(:,k) = colloids% image(:,k) + jump
119  colloids% pos(:,k) = colloids% pos_old(:,k) - jump*solvent_cells% edges
120  end do
121 
122  call switch(solvent% force, solvent% force_old)
123  call switch(colloids% force, colloids% force_old)
124  solvent% force = 0
125  colloids% force = 0
126  e1 = compute_force(colloids, solvent, neigh, solvent_cells% edges, solvent_colloid_lj)
127  e2 = compute_force_n2(colloids, solvent_cells% edges, colloid_lj)
128 
129  !$omp parallel do private(k)
130  do k = 1, solvent% Nmax
131  solvent% vel(:,k) = solvent% vel(:,k) + dt * ( solvent% force(:,k) + solvent% force_old(:,k) ) / 2
132  end do
133  colloids% vel = colloids% vel + dt * ( colloids% force + colloids% force_old ) / (2 * mass)
134 
135  end do md_equil
136 
137  call solvent% sort(solvent_cells)
138  call neigh% update_list(colloids, solvent, 1.5d0, solvent_cells)
139  call simple_mpcd_step(solvent, solvent_cells, state)
140 
141  end do
142 
143  write(*,*) ''
144  write(*,*) ' e co so | e co co | kin co | kin so | total | temp |'
145  write(*,*) ''
146 
147  do i = 1, 750
148  so_max = 0
149  co_max = 0
150  md: do j = 1, n_md_steps
151  solvent% pos_old = solvent% pos + dt * solvent% vel + dt**2 * solvent% force / 2
152  !$omp parallel do private(k)
153  do k = 1, solvent% Nmax
154  solvent% pos(:,k) = modulo(solvent% pos_old(:,k), solvent_cells% edges)
155  end do
156  colloids% pos_old = colloids% pos + dt * colloids% vel + dt**2 * colloids% force / (2 * mass)
157  do k = 1, colloids% Nmax
158  jump = floor(colloids% pos_old(:,k) / solvent_cells% edges)
159  colloids% image(:,k) = colloids% image(:,k) + jump
160  colloids% pos(:,k) = colloids% pos_old(:,k) - jump*solvent_cells% edges
161  end do
162  so_max = max(maxval(sqrt(sum((solvent% pos - solvent% pos_old)**2, dim=1))), so_max)
163  co_max = max(maxval(sqrt(sum((colloids% pos - colloids% pos_old)**2, dim=1))), co_max)
164 
165  call switch(solvent% force, solvent% force_old)
166  call switch(colloids% force, colloids% force_old)
167  solvent% force = 0
168  colloids% force = 0
169  e1 = compute_force(colloids, solvent, neigh, solvent_cells% edges, solvent_colloid_lj)
170  e2 = compute_force_n2(colloids, solvent_cells% edges, colloid_lj)
171 
172  !$omp parallel do private(k)
173  do k = 1, solvent% Nmax
174  solvent% vel(:,k) = solvent% vel(:,k) + dt * ( solvent% force(:,k) + solvent% force_old(:,k) ) / 2
175  end do
176  colloids% vel = colloids% vel + dt * ( colloids% force + colloids% force_old ) / (2 * mass)
177 
178  end do md
179 
180  write(15,*) colloids% pos + colloids% image * spread(solvent_cells% edges, dim=2, ncopies=colloids% Nmax)
181 
182  call solvent% sort(solvent_cells)
183  call neigh% update_list(colloids, solvent, 1.5d0, solvent_cells)
184 
185  call simple_mpcd_step(solvent, solvent_cells, state)
186 
187  write(*,'(6f16.3)') e1, e2, mass*sum(colloids% vel**2)/2, sum(solvent% vel**2)/2, &
188  e1+e2+mass*sum(colloids% vel**2)/2+sum(solvent% vel**2)/2, &
189  compute_temperature(solvent, solvent_cells)
190 
191  end do
192 
193 
194  call h5close_f(error)
195 
196 end program setup_simple_colloid
Derived type and routines for neighbor listing.
double precision function, public compute_temperature(particles, cells, tz)
Compute the temperature of a mpcd solvent.
Definition: mpcd.f90:343
Data for particles.
subroutine, public simple_mpcd_step(particles, cells, state, alpha, thermostat, T, hydro)
Perform a collision.
Definition: mpcd.f90:72
double precision function, public compute_force(ps1, ps2, n_list, L, lj_params)
double precision function, public compute_force_n2(ps, L, lj_params)
Compute compact Hilbert indices.
Definition: hilbert.f90:9
Utility routines.
Definition: common.f90:12
integer, parameter dim
Definition: hilbert.f90:15
Routines to perform MPCD dynamics.
Definition: mpcd.f90:18
Routines to perform Molecular Dynamics integration.
Definition: md.f90:13
Spatial cells.
Definition: cell_system.f90:11
program setup_simple_colloid
Lennard-Jones potential definition.
Definition: interaction.f90:10