RMPCDMD
setup_simple_colloids.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 md
17  use iso_c_binding
18  use omp_lib
19  implicit none
20 
21  type(cell_system_t) :: solvent_cells
22  type(particle_system_t) :: solvent
23  type(particle_system_t) :: colloids
24  type(neighbor_list_t) :: neigh
25  type(lj_params_t) :: solvent_colloid_lj
26  type(lj_params_t) :: colloid_lj
27 
28  integer, parameter :: rho = 10
29  double precision, parameter :: sigma = 1
30  double precision, parameter :: sigma_cut = sigma*2**(1.d0/6.d0)
31  integer :: error
32 
33  double precision :: epsilon
34  double precision :: mass
35  double precision :: so_max, co_max
36 
37  double precision :: e1, e2
38  double precision :: tau, dt
39  integer :: N_MD_steps
40  integer :: n_extra_sorting
41  integer :: n_threads
42 
43  type(threefry_rng_t), allocatable :: state(:)
44 
45  double precision :: skin
46 
47  integer :: i, L(3), N
48  integer :: j, k
49 
50  n_threads = omp_get_max_threads()
51  allocate(state(n_threads))
52  call threefry_rng_init(state, 719287321987291_c_long)
53 
54  call h5open_f(error)
55 
56  call colloids% init_from_file('input_data.h5', 'colloids', h5md_linear, 4)
57  write(*, *) colloids% pos
58  colloids% species = 1
59  colloids% vel = 0
60 
61  l = 10
62  n = rho*l(1)*l(2)*l(3) - int(rho*4*3.142/3 * sigma**3*colloids%Nmax)
63  if (n <= 0) error stop 'Not enough volume available for solvent'
64 
65  epsilon = 1
66 
67  call solvent_colloid_lj% init( reshape( [ epsilon ], [1, 1] ), &
68  reshape( [ sigma ], [1, 1] ), reshape( [ sigma_cut ], [1, 1] ) )
69 
70  epsilon = 1
71 
72  call colloid_lj% init( reshape( [ epsilon ], [1, 1] ), &
73  reshape( [ 2*sigma ], [1, 1] ), reshape( [ 2*sigma_cut ], [1, 1] ) )
74 
75  mass = 3000.d0 / (l(1)*l(2)*l(3)) * sigma**3 * 4 * 3.141/3
76 
77  call solvent% init(n)
78 
79  call random_number(solvent% vel(:, :))
80  solvent% vel = (solvent% vel - 0.5d0)*sqrt(6.d0*2)
81  solvent% vel = solvent% vel - spread(sum(solvent% vel, dim=2)/solvent% Nmax, 2, solvent% Nmax)
82  solvent% force = 0
83  solvent% species = 1
84 
85  call solvent_cells%init(l, 1.d0)
86 
87  call solvent% random_placement(solvent_cells% edges, colloids, solvent_colloid_lj, state=state(1))
88 
89  call solvent% sort(solvent_cells)
90 
91  call neigh% init(colloids% Nmax, 500)
92 
93  tau = 0.1d0
94  n_md_steps = 100
95  dt = tau / n_md_steps
96 
97  skin = 0.8
98  n_extra_sorting = 0
99 
100  call neigh% make_stencil(solvent_cells, sigma_cut+skin)
101  call neigh% update_list(colloids, solvent, sigma_cut+skin, solvent_cells)
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  write(*,*) ''
109  write(*,*) ' e co so | e co co | kin co | kin so | total | temp |'
110  write(*,*) ''
111 
112  solvent% pos_old = solvent% pos
113  colloids% pos_old = colloids% pos
114  do i = 1, 2000
115  md_loop: do j = 1, n_md_steps
116  call md_pos(solvent, dt)
117  do k = 1, colloids% Nmax
118  colloids% pos(:,k) = colloids% pos(:,k) + dt * colloids% vel(:,k) + dt**2 * colloids% force(:,k) / (2 * mass)
119  end do
120  so_max = solvent% maximum_displacement()
121  co_max = colloids% maximum_displacement()
122 
123  if ( (co_max >= skin*0.1) .or. (so_max >= skin*0.9) ) then
124  call apply_pbc(solvent, solvent_cells% edges)
125  call apply_pbc(colloids, solvent_cells% edges)
126  call solvent% sort(solvent_cells)
127  call neigh% update_list(colloids, solvent, sigma_cut+skin, solvent_cells)
128  solvent% pos_old = solvent% pos
129  colloids% pos_old = colloids% pos
130  n_extra_sorting = n_extra_sorting + 1
131  end if
132 
133  call switch(solvent% force, solvent% force_old)
134  call switch(colloids% force, colloids% force_old)
135  solvent% force = 0
136  colloids% force = 0
137  e1 = compute_force(colloids, solvent, neigh, solvent_cells% edges, solvent_colloid_lj)
138  e2 = compute_force_n2(colloids, solvent_cells% edges, colloid_lj)
139 
140  call md_vel(solvent, dt)
141  colloids% vel = colloids% vel + dt * ( colloids% force + colloids% force_old ) / (2 * mass)
142 
143  write(15,*) colloids% pos + colloids% image * spread(solvent_cells% edges, dim=2, ncopies=colloids% Nmax)
144 
145  end do md_loop
146 
147  call apply_pbc(solvent, solvent_cells% edges)
148  call apply_pbc(colloids, solvent_cells% edges)
149  call solvent% sort(solvent_cells)
150  call neigh% update_list(colloids, solvent, sigma_cut + skin, solvent_cells)
151  solvent% pos_old = solvent% pos
152  colloids% pos_old = colloids% pos
153 
154  call simple_mpcd_step(solvent, solvent_cells, state)
155 
156  if (modulo(i,10)==0) &
157  write(*,'(6f16.3)') e1, e2, mass*sum(colloids% vel**2)/2, sum(solvent% vel**2)/2, &
158  e1+e2+mass*sum(colloids% vel**2)/2+sum(solvent% vel**2)/2, &
159  compute_temperature(solvent, solvent_cells)
160 
161  end do
162 
163  write(*,*) 'n extra sorting', n_extra_sorting
164 
165  call h5close_f(error)
166 
167 end program setup_simple_colloids
Derived type and routines for neighbor listing.
subroutine, public md_vel(particles, dt)
Definition: md.f90:220
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)
subroutine, public apply_pbc(particles, edges)
Definition: md.f90:203
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
subroutine, public md_pos(particles, dt)
Definition: md.f90:62
program setup_simple_colloids
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
Lennard-Jones potential definition.
Definition: interaction.f90:10