-
Notifications
You must be signed in to change notification settings - Fork 12
BDHI FCM
BDHI_FCM is a BDHI (Brownian Hydrodynamics) submodule[4]. Source is at Integrator/BDHI/BDHI_FCM.cuh
As this is a BDHI module. BDHI_FCM computes the terms M·F and B·dW in the differential equation:
dR = K·R·dt + M·F·dt + sqrt(2Tdt)· B·dW
Where:
R: Particle positions
K: Shear matrix (FCM does not implement this in PBC)
M: FCM Mobility tensor
T: Temperature
F: Forces acting on the particles.
dW: Collection of Weinner processes
B: B·B^T = M
Using periodic boundary conditions.
The mobility, M, is computed according to the Force Coupling Method (FCM) tensor[1].
Use as any other BDHI method. You have an example in examples/FCM.cu.
BDHI::FCM has the following parameters:
...
BDHI::FCM::Parameters par;
par.temperature = 1.0;
par.viscosity = 1.0;
par.hydrodynamicRadius = 1.0;
par.dt = 1e-3;
par.box = box;
par.tolerance = 1e-3;
auto bdhi = make_shared<BDHI::EulerMaruyama<BDHI::FCM>>(pd, pg, sys, par);
...
The algorithm behind FCM is very similar to the PSE one but without the near field part. The particle forcing is spreaded to a grid using Gaussian kernels and then the fluid velocity (including fluctuations) is solved in Fourier space. Finally the fluid velocity is interpolated back to the particles using the same Gaussian kernels.
Using the operator terminology in [2]:
M·F = σ·St·FFTi·B·FFTf·(S·F+f)
-σ: The volume of a grid cell
-S: An operator that spreads each element of a vector to a regular grid using a gaussian kernel.
-F: Forces acting on the particles
-f: Force density acting on the fluid
-FFT: Fast fourier transform operator.
-B: 1/(\eta k^2)·(I-k\otimes k) in Fourier space. Multiplies fluid forcing by the FCM kernel, resulting in fluid velocity, see eq.59 in [1].
Related functions:
FFT: cufftExecR2C (forward), cufftC2R(inverse)
S: FCM_ns::particles2Grid (S), FCM_ns::grid2Particles (σ·St)
B: FCM_ns::forceFourier2vel
sqrt(M)·dW: The stochastic contribution is computed in fourier space along M·F as:
M·F + sqrt(M)·dW = σ·St·FFTi·B·FFTf·(S·F+f)+ √σ·St·FFTi·√B·dW =
= σ·St·FFTi( B·FFTf·(S·F+f) + 1/√σ·√B·dW)
Only one St·FFTi is needed, the stochastic term is added as a velocity in fourier space. dW is a gaussian random vector of complex numbers, special care must be taken to ensure the correct conjugacy properties needed for the FFT. See FCM_ns::fourierBrownianNoise
[1] Fluctuating force-coupling method for simulations of colloidal suspensions. Eric E. Keaveny. 2014.
[2] Rapid Sampling of Stochastic Displacements in Brownian Dynamics Simulations - https://arxiv.org/pdf/1611.09322.pdf
[3] Spectral accuracy in fast Ewald-based methods for particle simulations
[4] https://github.com/RaulPPelaez/UAMMD/wiki/Brownian-Hydrodynamics
-
-
1. PairForces
2. NbodyForces
3. ExternalForces
4. BondedForces
5. AngularBondedForces
6. TorsionalBondedForces
7. Poisson (Electrostatics) -
-
MD (Molecular Dynamics)
1. VerletNVT
2. VerletNVE - BD Brownian Dynamics
-
BDHI Brownian Dynamics with Hydrodynamic Interactions
1. EulerMaruyama
1.1 BDHI_Cholesky Brownian displacements through Cholesky factorization.
1.2 BDHI_Lanczos Brownian displacements through Lanczos algorithm.
1.3 BDHI_PSE Positively Split Edwald.
1.4 BDHI_FCM Force Coupling Method. - DPD Dissipative Particle Dynamics
- SPH Smoothed Particle Hydrodynamics
-
Hydrodynamics
1. ICM Inertial Coupling Method.
2. FIB Fluctuating Immerse Boundary.
3. Quasi2D Quasi2D hydrodynamics
-
MD (Molecular Dynamics)
-
- 1. Neighbour Lists
-
1. Programming Tools
2. Utils
-
1. Transverser
2. Functor
3. Potential
-
1. Particle Data
2. Particle Group
3. System
4. Parameter updatable
-
1. Tabulated Function
2. Postprocessing tools
3. InputFile
4. Tests
5. Allocator
6. Temporary memory
7. Immersed Boundary (IBM)
-
1. NBody
2. Neighbour Lists
3. Python wrappers
- 1. Superpunto