From fa43416f5c8fc0d922a15d419ffab2da07e43f29 Mon Sep 17 00:00:00 2001 From: GregorySchwing <39970712+GregorySchwing@users.noreply.github.com> Date: Mon, 15 Mar 2021 16:05:02 -0400 Subject: [PATCH] Revert "Revert "MultiParticle Brownian like motion"" --- CMake/FileLists.cmake | 3 + lib/BasicTypes.h | 12 + src/BoxDimensions.h | 5 + src/ConfigSetup.cpp | 35 +- src/ConfigSetup.h | 4 +- src/ConsoleOutput.cpp | 13 + src/GOMCEventsProfileDef.h | 5 + src/GPU/CalculateMinImageCUDAKernel.cuh | 41 ++ src/GPU/ConstantDefinitionsCUDAKernel.cu | 1 + src/GPU/TransformParticlesCUDAKernel.cu | 520 +++++++++++++++++++++- src/GPU/TransformParticlesCUDAKernel.cuh | 99 +++++ src/GPU/VariablesCUDA.cuh | 2 + src/MoleculeLookup.cpp | 20 +- src/MoleculeLookup.h | 4 +- src/MoveConst.h | 86 ++-- src/MoveSettings.cpp | 10 +- src/PRNG.h | 26 +- src/Random123Wrapper.cpp | 58 +++ src/Random123Wrapper.h | 46 +- src/StaticVals.cpp | 5 +- src/System.cpp | 9 +- src/moves/MultiParticle.h | 32 +- src/moves/MultiParticleBrownianMotion.h | 530 +++++++++++++++++++++++ 23 files changed, 1447 insertions(+), 119 deletions(-) create mode 100644 src/Random123Wrapper.cpp create mode 100644 src/moves/MultiParticleBrownianMotion.h diff --git a/CMake/FileLists.cmake b/CMake/FileLists.cmake index 6f45571c6..5e7bb14a2 100644 --- a/CMake/FileLists.cmake +++ b/CMake/FileLists.cmake @@ -42,6 +42,7 @@ set(sources src/PDBOutput.cpp src/PRNGSetup.cpp src/PSFOutput.cpp + src/Random123Wrapper.cpp src/Reader.cpp src/Simulation.cpp src/StaticVals.cpp @@ -131,6 +132,7 @@ set(headers src/PRNG.h src/PRNGSetup.h src/PSFOutput.h + src/Random123Wrapper.h src/Reader.h src/SeedReader.h src/Setup.h @@ -174,6 +176,7 @@ set(headers src/moves/MoleculeTransfer.h src/moves/MoveBase.h src/moves/MultiParticle.h + src/moves/MultiParticleBrownianMotion.h src/moves/Regrowth.h src/moves/Rotation.h src/moves/TargetedSwap.h diff --git a/lib/BasicTypes.h b/lib/BasicTypes.h index a8a770cd6..1a7f594c8 100644 --- a/lib/BasicTypes.h +++ b/lib/BasicTypes.h @@ -123,6 +123,18 @@ struct XYZ { return true; return false; } + bool operator<(XYZ const& rhs) + { + if(x < rhs.x && y < rhs.y && z < rhs.z) + return true; + return false; + } + bool operator>(XYZ const& rhs) + { + if(x > rhs.x || y > rhs.y || z > rhs.z) + return true; + return false; + } XYZ& operator+=(XYZ const& rhs) { x += rhs.x; diff --git a/src/BoxDimensions.h b/src/BoxDimensions.h index ecf0261c6..8dbbb854b 100644 --- a/src/BoxDimensions.h +++ b/src/BoxDimensions.h @@ -49,6 +49,11 @@ class BoxDimensions { return axis.Get(b); } + + XYZ GetHalfAxis(const uint b) const + { + return halfAx.Get(b); + } double GetTotVolume(const uint b1, const uint b2) const; diff --git a/src/ConfigSetup.cpp b/src/ConfigSetup.cpp index bbf099198..c3550ac3e 100644 --- a/src/ConfigSetup.cpp +++ b/src/ConfigSetup.cpp @@ -76,6 +76,7 @@ ConfigSetup::ConfigSetup(void) sys.moves.intraSwap = DBL_MAX; sys.moves.multiParticleEnabled = false; sys.moves.multiParticle = DBL_MAX; + sys.moves.multiParticleBrownian = DBL_MAX; sys.moves.regrowth = DBL_MAX; sys.moves.crankShaft = DBL_MAX; sys.moves.intraMemc = DBL_MAX; @@ -678,6 +679,14 @@ void ConfigSetup::Init(const char *fileName, MultiSim const*const& multisim) printf("%-40s %-4.4f \n", "Info: Multi-Particle move frequency", sys.moves.multiParticle); + } else if(CheckString(line[0], "MultiParticleBrownianFreq")) { + sys.moves.multiParticleBrownian = stringtod(line[1]); + if(sys.moves.multiParticleBrownian > 0.00) { + sys.moves.multiParticleEnabled = true; + } + printf("%-40s %-4.4f \n", + "Info: Multi-Particle Brownian move frequency", + sys.moves.multiParticleBrownian); } else if(CheckString(line[0], "IntraSwapFreq")) { sys.moves.intraSwap = stringtod(line[1]); printf("%-40s %-4.4f \n", "Info: Intra-Swap move frequency", @@ -1251,6 +1260,13 @@ void ConfigSetup::fillDefaults(void) sys.moves.multiParticle); } + if(sys.moves.multiParticleBrownian == DBL_MAX) { + sys.moves.multiParticleBrownian = 0.000; + printf("%-40s %-4.4f \n", + "Default: Multi-Particle Brownian move frequency", + sys.moves.multiParticleBrownian); + } + if(sys.moves.intraMemc == DBL_MAX) { sys.moves.intraMemc = 0.0; printf("%-40s %-4.4f \n", "Default: Intra-MEMC move frequency", @@ -1506,6 +1522,12 @@ void ConfigSetup::verifyInputs(void) } } #endif + if(abs(sys.moves.multiParticle) > 0.0000001 && + abs(sys.moves.multiParticleBrownian) > 0.0000001) { + std::cout << "Error: Both multi-Particle and multi-Particle Brownian! " << + " cannot be used at the same time!" << std::endl; + exit(EXIT_FAILURE); + } if(!sys.elect.enable && sys.elect.oneFourScale != DBL_MAX) { printf("Warning: 1-4 Electrostatic scaling set, but will be ignored.\n"); @@ -1660,8 +1682,8 @@ void ConfigSetup::verifyInputs(void) if(std::abs(sys.moves.displace + sys.moves.rotate + sys.moves.transfer + sys.moves.intraSwap + sys.moves.volume + sys.moves.regrowth + sys.moves.memc + sys.moves.intraMemc + sys.moves.crankShaft + - sys.moves.multiParticle + sys.moves.cfcmc + - sys.moves.targetedSwap - 1.0) > 0.001) { + sys.moves.multiParticle + sys.moves.multiParticleBrownian + + sys.moves.cfcmc + sys.moves.targetedSwap - 1.0) > 0.001) { std::cout << "Error: Sum of move frequencies is not equal to one!\n"; exit(EXIT_FAILURE); } @@ -1673,7 +1695,8 @@ void ConfigSetup::verifyInputs(void) if(std::abs(sys.moves.displace + sys.moves.rotate + sys.moves.intraSwap + sys.moves.volume + sys.moves.regrowth + sys.moves.intraMemc + - sys.moves.crankShaft + sys.moves.multiParticle - 1.0) > 0.001) { + sys.moves.crankShaft + sys.moves.multiParticle + + sys.moves.multiParticleBrownian - 1.0) > 0.001) { std::cout << "Error: Sum of move frequencies is not equal to one!\n"; exit(EXIT_FAILURE); } @@ -1686,15 +1709,15 @@ void ConfigSetup::verifyInputs(void) if(std::abs(sys.moves.displace + sys.moves.rotate + sys.moves.intraSwap + sys.moves.transfer + sys.moves.regrowth + sys.moves.memc + sys.moves.intraMemc + sys.moves.crankShaft + - sys.moves.multiParticle + sys.moves.cfcmc + - sys.moves.targetedSwap - 1.0) > 0.001) { + sys.moves.multiParticle + sys.moves.multiParticleBrownian + + sys.moves.cfcmc + sys.moves.targetedSwap - 1.0) > 0.001) { std::cout << "Error: Sum of move frequencies is not equal to one!!\n"; exit(EXIT_FAILURE); } #else if(std::abs(sys.moves.displace + sys.moves.rotate + sys.moves.intraSwap + sys.moves.regrowth + sys.moves.intraMemc + sys.moves.crankShaft + - sys.moves.multiParticle - 1.0) > 0.001) { + sys.moves.multiParticle + sys.moves.multiParticleBrownian - 1.0) > 0.001) { std::cout << "Error: Sum of move frequencies is not equal to one!!\n"; exit(EXIT_FAILURE); } diff --git a/src/ConfigSetup.h b/src/ConfigSetup.h index 947e33bb6..fd1280a53 100644 --- a/src/ConfigSetup.h +++ b/src/ConfigSetup.h @@ -183,8 +183,8 @@ struct Step { //Holds the percentage of each kind of move for this ensemble. struct MovePercents { double displace, rotate, intraSwap, intraMemc, regrowth, crankShaft, - multiParticle; - bool multiParticleEnabled; + multiParticle, multiParticleBrownian; + bool multiParticleEnabled; // for both multiparticle and multiparticleBrownian #ifdef VARIABLE_VOLUME double volume; #endif diff --git a/src/ConsoleOutput.cpp b/src/ConsoleOutput.cpp index d5f7899d2..2e3d951df 100644 --- a/src/ConsoleOutput.cpp +++ b/src/ConsoleOutput.cpp @@ -121,6 +121,13 @@ void ConsoleOutput::PrintMove(const uint box, const ulong step) const printElement(var->GetAcceptPercent(box, sub), elementWidth); } + if(var->Performed(mv::MULTIPARTICLE_BM)) { + sub = mv::MULTIPARTICLE_BM; + printElement(var->GetTries(box, sub), elementWidth); + printElement(var->GetAccepted(box, sub), elementWidth); + printElement(var->GetAcceptPercent(box, sub), elementWidth); + } + if(var->Performed(mv::INTRA_SWAP)) { sub = mv::INTRA_SWAP; printElement(var->GetTries(box, sub), elementWidth); @@ -371,6 +378,12 @@ void ConsoleOutput::PrintMoveTitle() printElement("MPACCEPT%", elementWidth); } + if(var->Performed(mv::MULTIPARTICLE_BM)) { + printElement("MULTIPARTICLEBM", elementWidth); + printElement("MPBMACCEPT", elementWidth); + printElement("MPBMACCEPT%", elementWidth); + } + if(var->Performed(mv::INTRA_SWAP)) { printElement("INTRASWAP", elementWidth); printElement("INTACCEPT", elementWidth); diff --git a/src/GOMCEventsProfileDef.h b/src/GOMCEventsProfileDef.h index f9dc81db7..a057c8992 100644 --- a/src/GOMCEventsProfileDef.h +++ b/src/GOMCEventsProfileDef.h @@ -20,6 +20,7 @@ GOMC_PROFILE_EVENT(FREE_ENERGY_OUTPUT, "free_energy_file") GOMC_PROFILE_EVENT(DISPLACE, "translate_move") GOMC_PROFILE_EVENT(ROTATE, "rotate_move") GOMC_PROFILE_EVENT(MULTIPARTICLE, "multiParticle_move") +GOMC_PROFILE_EVENT(MULTIPARTICLE_BM, "multiParticleBM_move") GOMC_PROFILE_EVENT(INTRA_SWAP, "intraSwap_move") GOMC_PROFILE_EVENT(INTRA_MEMC, "intraMEMC_move") GOMC_PROFILE_EVENT(CRANKSHAFT, "crankshaft_move") @@ -33,6 +34,7 @@ GOMC_PROFILE_EVENT(CFCMC, "CFCMC_move") GOMC_PROFILE_EVENT(PREP_DISPLACE, "prepare_translate_move") GOMC_PROFILE_EVENT(PREP_ROTATE, "prepare_rotate_move") GOMC_PROFILE_EVENT(PREP_MULTIPARTICLE, "prepare_multiParticle_move") +GOMC_PROFILE_EVENT(PREP_MULTIPARTICLE_BM, "prepare_multiParticleBM_move") GOMC_PROFILE_EVENT(PREP_INTRA_SWAP, "prepare_intraSwap_move") GOMC_PROFILE_EVENT(PREP_INTRA_MEMC, "prepare_intraMEMC_move") GOMC_PROFILE_EVENT(PREP_CRANKSHAFT, "prepare_crankshaft_move") @@ -46,6 +48,7 @@ GOMC_PROFILE_EVENT(PREP_CFCMC, "prepare_CFCMC_move") GOMC_PROFILE_EVENT(TRANS_DISPLACE, "transform_translate_move") GOMC_PROFILE_EVENT(TRANS_ROTATE, "transform_rotate_move") GOMC_PROFILE_EVENT(TRANS_MULTIPARTICLE, "transform_multiParticle_move") +GOMC_PROFILE_EVENT(TRANS_MULTIPARTICLE_BM, "transform_multiParticleBM_move") GOMC_PROFILE_EVENT(TRANS_INTRA_SWAP, "transform_intraSwap_move") GOMC_PROFILE_EVENT(TRANS_INTRA_MEMC, "transform_intraMEMC_move") GOMC_PROFILE_EVENT(TRANS_CRANKSHAFT, "transform_crankshaft_move") @@ -59,6 +62,7 @@ GOMC_PROFILE_EVENT(TRANS_CFCMC, "transform_CFCMC_move") GOMC_PROFILE_EVENT(CALC_EN_DISPLACE, "cal_energy_translate_move") GOMC_PROFILE_EVENT(CALC_EN_ROTATE, "cal_energy_rotate_move") GOMC_PROFILE_EVENT(CALC_EN_MULTIPARTICLE, "cal_energy_multiParticle_move") +GOMC_PROFILE_EVENT(CALC_EN_MULTIPARTICLE_BM, "cal_energy_multiParticleBM_move") GOMC_PROFILE_EVENT(CALC_EN_INTRA_SWAP, "cal_energy_intraSwap_move") GOMC_PROFILE_EVENT(CALC_EN_INTRA_MEMC, "cal_energy_intraMEMC_move") GOMC_PROFILE_EVENT(CALC_EN_CRANKSHAFT, "cal_energy_crankshaft_move") @@ -72,6 +76,7 @@ GOMC_PROFILE_EVENT(CALC_EN_CFCMC, "cal_energy_CFCMC_move") GOMC_PROFILE_EVENT(ACC_DISPLACE, "accept_translate_move") GOMC_PROFILE_EVENT(ACC_ROTATE, "accept_rotate_move") GOMC_PROFILE_EVENT(ACC_MULTIPARTICLE, "accept_multiParticle_move") +GOMC_PROFILE_EVENT(ACC_MULTIPARTICLE_BM, "accept_multiParticleBM_move") GOMC_PROFILE_EVENT(ACC_INTRA_SWAP, "accept_intraSwap_move") GOMC_PROFILE_EVENT(ACC_INTRA_MEMC, "accept_intraMEMC_move") GOMC_PROFILE_EVENT(ACC_CRANKSHAFT, "accept_crankshaft_move") diff --git a/src/GPU/CalculateMinImageCUDAKernel.cuh b/src/GPU/CalculateMinImageCUDAKernel.cuh index f8e14fe93..7ca74252a 100644 --- a/src/GPU/CalculateMinImageCUDAKernel.cuh +++ b/src/GPU/CalculateMinImageCUDAKernel.cuh @@ -38,6 +38,47 @@ __device__ inline void TransformUnSlantGPU(double3 & dist, dist.z = slant.x * gpu_Invcell_z[0] + slant.y * gpu_Invcell_z[1] + slant.z * gpu_Invcell_z[2]; } + +__device__ inline void WrapPBC(double &v, double &ax) +{ + if(v >= ax) + v -= ax; + else if(v < 0) + v += ax; +} + +__device__ inline void WrapPBC_f3(double3 &v, double3 &ax) +{ + WrapPBC(v.x, ax.x); + WrapPBC(v.y, ax.y); + WrapPBC(v.z, ax.z); +} + +__device__ inline void UnwrapPBC( + double &v, + double &ref, + double &ax, + double &halfax) +{ + if(abs(ref - v) > halfax) { + if(ref < halfax) + v -= ax; + else + v += ax; + } +} + +__device__ inline void UnwrapPBC_f3( + double3 &v, + double3 &ref, + double3 &ax, + double3 &halfax) +{ + UnwrapPBC(v.x, ref.x, ax.x, halfax.x); + UnwrapPBC(v.y, ref.y, ax.y, halfax.y); + UnwrapPBC(v.z, ref.z, ax.z, halfax.z); +} + __device__ inline double MinImageSignedGPU(double raw, double ax, double halfAx) { if (raw > halfAx) diff --git a/src/GPU/ConstantDefinitionsCUDAKernel.cu b/src/GPU/ConstantDefinitionsCUDAKernel.cu index 822d86deb..0a22f4bc6 100644 --- a/src/GPU/ConstantDefinitionsCUDAKernel.cu +++ b/src/GPU/ConstantDefinitionsCUDAKernel.cu @@ -373,6 +373,7 @@ void DestroyCUDAVars(VariablesCUDA *vars) CUFREE(vars->gpu_cellVector); CUFREE(vars->gpu_mapParticleToCell); CUFREE(vars->gpu_nonOrth); + CUFREE(vars->gpu_startAtomIdx); for(uint b = 0; b < BOX_TOTAL; b++) { CUFREE(vars->gpu_cell_x[b]); CUFREE(vars->gpu_cell_y[b]); diff --git a/src/GPU/TransformParticlesCUDAKernel.cu b/src/GPU/TransformParticlesCUDAKernel.cu index d928c5d5e..a2ccf3068 100644 --- a/src/GPU/TransformParticlesCUDAKernel.cu +++ b/src/GPU/TransformParticlesCUDAKernel.cu @@ -6,7 +6,9 @@ along with this program, also can be found at . ********************************************************************************/ #ifdef GOMC_CUDA #include "TransformParticlesCUDAKernel.cuh" +#include "CalculateMinImageCUDAKernel.cuh" #include "CUDAMemoryManager.cuh" +#include "Random123/boxmuller.hpp" #define MIN_FORCE 1E-12 #define MAX_FORCE 30 @@ -23,24 +25,19 @@ __device__ inline double randomGPU(unsigned int counter, unsigned int step, unsi return (double)r[0] / UINT_MAX; } -__device__ inline double WrapPBC(double &v, double ax) +__device__ inline double randomGaussianGPU(unsigned int counter, unsigned int step, + unsigned int seed, double mean, double stdDev) { - if(v >= ax) - v -= ax; - else if(v < 0) - v += ax; - return v; -} - -__device__ inline double UnwrapPBC(double &v, double ref, double ax, double halfax) -{ - if(abs(ref - v) > halfax) { - if(ref < halfax) - v -= ax; - else - v += ax; - } - return v; + RNG::ctr_type c = {{}}; + RNG::ukey_type uk = {{}}; + uk[0] = step; + uk[1] = seed; + RNG::key_type k = uk; + c[0] = counter; + RNG::ctr_type r = philox4x32(c, k); + float2 normalf2 = r123::boxmuller(r[0], r[1]); + double shiftedVal = mean + (double)(normalf2.x) * stdDev; + return shiftedVal; } __device__ inline void ApplyRotation(double &x, double &y, double &z, @@ -53,6 +50,9 @@ __device__ inline void ApplyRotation(double &x, double &y, double &z, double axisy = roty * (1.0 / rotLen); double axisz = rotz * (1.0 / rotLen); double matrix[3][3], cross[3][3], tensor[3][3]; + double halfAxx = axisx * 0.5; + double halfAxy = axisy * 0.5; + double halfAxz = axisz * 0.5; // build cross cross[0][0] = 0.0; @@ -95,9 +95,9 @@ __device__ inline void ApplyRotation(double &x, double &y, double &z, } // unwrap molecule - UnwrapPBC(x, comx, axx, axx / 2.0); - UnwrapPBC(y, comy, axy, axy / 2.0); - UnwrapPBC(z, comz, axz, axz / 2.0); + UnwrapPBC(x, comx, axx, halfAxx); + UnwrapPBC(y, comy, axy, halfAxy); + UnwrapPBC(z, comz, axz, halfAxz); // move particle to zero x -= comx; @@ -461,4 +461,484 @@ __global__ void RotateParticlesKernel(unsigned int numberOfMolecules, rotx, roty, rotz, xAxes, yAxes, zAxes); } +// CUDA implementation of MultiParticle Brownian transformation + +void BrownianMotionRotateParticlesGPU( + VariablesCUDA *vars, + const std::vector &moleculeInvolved, + XYZArray &mTorque, + XYZArray &newMolPos, + XYZArray &newCOMs, + XYZArray &r_k, + const XYZ &boxAxes, + const double BETA, + const double r_max, + unsigned int step, + unsigned int seed, + const int box, + bool isOrthogonal) +{ + int atomCount = newMolPos.Count(); + int molCount = newCOMs.Count(); + int molCountInBox = moleculeInvolved.size(); + int *gpu_moleculeInvolved; + // Each block would handle one molecule + int threadsPerBlock = 64; + int blocksPerGrid = molCountInBox; + + CUMALLOC((void **) &gpu_moleculeInvolved, molCountInBox * sizeof(int)); + + cudaMemcpy(vars->gpu_mTorquex, mTorque.x, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_mTorquey, mTorque.y, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_mTorquez, mTorque.z, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_x, newMolPos.x, atomCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_y, newMolPos.y, atomCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_z, newMolPos.z, atomCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_comx, newCOMs.x, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_comy, newCOMs.y, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_comz, newCOMs.z, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(gpu_moleculeInvolved, &moleculeInvolved[0], molCountInBox * sizeof(int), cudaMemcpyHostToDevice); + + double3 axis = make_double3(boxAxes.x, boxAxes.y, boxAxes.z); + double3 halfAx = make_double3(boxAxes.x / 2.0, boxAxes.y / 2.0, boxAxes.z / 2.0); + + if(isOrthogonal) { + BrownianMotionRotateKernel<<< blocksPerGrid, threadsPerBlock>>>( + vars->gpu_startAtomIdx, + vars->gpu_x, + vars->gpu_y, + vars->gpu_z, + vars->gpu_mTorquex, + vars->gpu_mTorquey, + vars->gpu_mTorquez, + vars->gpu_comx, + vars->gpu_comy, + vars->gpu_comz, + vars->gpu_r_k_x, + vars->gpu_r_k_y, + vars->gpu_r_k_z, + gpu_moleculeInvolved, + vars->gpu_cell_x[box], + vars->gpu_cell_y[box], + vars->gpu_cell_z[box], + vars->gpu_Invcell_x[box], + vars->gpu_Invcell_y[box], + vars->gpu_Invcell_z[box], + axis, + halfAx, + atomCount, + r_max, + step, + seed, + BETA); + } else { + BrownianMotionRotateKernel<<< blocksPerGrid, threadsPerBlock>>>( + vars->gpu_startAtomIdx, + vars->gpu_x, + vars->gpu_y, + vars->gpu_z, + vars->gpu_mTorquex, + vars->gpu_mTorquey, + vars->gpu_mTorquez, + vars->gpu_comx, + vars->gpu_comy, + vars->gpu_comz, + vars->gpu_r_k_x, + vars->gpu_r_k_y, + vars->gpu_r_k_z, + gpu_moleculeInvolved, + vars->gpu_cell_x[box], + vars->gpu_cell_y[box], + vars->gpu_cell_z[box], + vars->gpu_Invcell_x[box], + vars->gpu_Invcell_y[box], + vars->gpu_Invcell_z[box], + axis, + halfAx, + atomCount, + r_max, + step, + seed, + BETA); + } + + cudaDeviceSynchronize(); + checkLastErrorCUDA(__FILE__, __LINE__); + + cudaMemcpy(newMolPos.x, vars->gpu_x, atomCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(newMolPos.y, vars->gpu_y, atomCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(newMolPos.z, vars->gpu_z, atomCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(r_k.x, vars->gpu_r_k_x, molCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(r_k.y, vars->gpu_r_k_y, molCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(r_k.z, vars->gpu_r_k_z, molCount * sizeof(double), cudaMemcpyDeviceToHost); + CUFREE(gpu_moleculeInvolved); + checkLastErrorCUDA(__FILE__, __LINE__); +} + +template +__global__ void BrownianMotionRotateKernel( + int *startAtomIdx, + double *gpu_x, + double *gpu_y, + double *gpu_z, + double *molTorquex, + double *molTorquey, + double *molTorquez, + double *gpu_comx, + double *gpu_comy, + double *gpu_comz, + double *gpu_r_k_x, + double *gpu_r_k_y, + double *gpu_r_k_z, + int *moleculeInvolved, + double *gpu_cell_x, + double *gpu_cell_y, + double *gpu_cell_z, + double *gpu_Invcell_x, + double *gpu_Invcell_y, + double *gpu_Invcell_z, + double3 axis, + double3 halfAx, + int atomCount, + double r_max, + unsigned int step, + unsigned int seed, + double BETA) +{ + //Each grid take cares of one molecule + int molIndex = moleculeInvolved[blockIdx.x]; + int startIdx = startAtomIdx[molIndex]; + int endIdx = startAtomIdx[molIndex + 1]; + int atomIdx; + + __shared__ double matrix[3][3]; + __shared__ double3 com; + + // thread 0 will setup the matrix and update the gpu_r_k + if(threadIdx.x == 0) { + com = make_double3(gpu_comx[molIndex], gpu_comy[molIndex], gpu_comz[molIndex]); + // This section calculates the amount of rotation + double stdDev = sqrt(2.0 * r_max); + double btm_x = molTorquex[molIndex] * BETA * r_max; + double btm_y = molTorquey[molIndex] * BETA * r_max; + double btm_z = molTorquez[molIndex] * BETA * r_max; + + double rot_x = btm_x + randomGaussianGPU(molIndex * 3, step, seed, 0.0, stdDev); + double rot_y = btm_y + randomGaussianGPU(molIndex * 3 + 1, step, seed, 0.0, stdDev); + double rot_z = btm_z + randomGaussianGPU(molIndex * 3 + 2, step, seed, 0.0, stdDev); + // update the trial torque + gpu_r_k_x[molIndex] = rot_x; + gpu_r_k_y[molIndex] = rot_y; + gpu_r_k_z[molIndex] = rot_z; + // build rotation matrix + double cross[3][3], tensor[3][3]; + double rotLen = sqrt(rot_x * rot_x + rot_y * rot_y + rot_z * rot_z); + double axisx = rot_x * (1.0 / rotLen); + double axisy = rot_y * (1.0 / rotLen); + double axisz = rot_z * (1.0 / rotLen); + // build cross + cross[0][0] = 0.0; cross[0][1] = -axisz; cross[0][2] = axisy; + cross[1][0] = axisz; cross[1][1] = 0.0; cross[1][2] = -axisx; + cross[2][0] = -axisy; cross[2][1] = axisx; cross[2][2] = 0.0; + // build tensor + int i, j; + for(i = 0; i < 3; ++i) { + tensor[0][i] = axisx; + tensor[1][i] = axisy; + tensor[2][i] = axisz; + } + for(i = 0; i < 3; ++i) { + tensor[i][0] *= axisx; + tensor[i][1] *= axisy; + tensor[i][2] *= axisz; + } + // build matrix + double s, c; + sincos(rotLen, &s, &c); + for(i = 0; i < 3; ++i) { + for(j = 0; j < 3; ++j) { + matrix[i][j] = 0.0; + } + matrix[i][i] = c; + } + for(i = 0; i < 3; ++i) { + for(j = 0; j < 3; ++j) { + matrix[i][j] += s * cross[i][j] + (1 - c) * tensor[i][j]; + } + } + } + + __syncthreads(); + // use strid of blockDim.x, which is 64 + // each thread handles one atom rotation + for(atomIdx = startIdx + threadIdx.x; atomIdx < endIdx; atomIdx += blockDim.x) { + double3 coor = make_double3(gpu_x[atomIdx], gpu_y[atomIdx], gpu_z[atomIdx]); + // unwrap molecule + if(isOrthogonal) { + UnwrapPBC_f3(coor, com, axis, halfAx); + } else { + double3 unSlant = make_double3(0.0, 0.0, 0.0); + TransformUnSlantGPU(unSlant, coor, gpu_Invcell_x, gpu_Invcell_y, gpu_Invcell_z); + UnwrapPBC_f3(unSlant, com, axis, halfAx); + TransformSlantGPU(coor, unSlant, gpu_cell_x, gpu_cell_y, gpu_cell_z); + } + + // move COM of molecule to zero + coor.x -= com.x; + coor.y -= com.y; + coor.z -= com.z; + // rotate + double newx = matrix[0][0] * coor.x + matrix[0][1] * coor.y + matrix[0][2] * coor.z; + double newy = matrix[1][0] * coor.x + matrix[1][1] * coor.y + matrix[1][2] * coor.z; + double newz = matrix[2][0] * coor.x + matrix[2][1] * coor.y + matrix[2][2] * coor.z; + + // move back to com + coor.x = newx + com.x; + coor.y = newy + com.y; + coor.z = newz + com.z; + + // wrap again + if(isOrthogonal) { + WrapPBC_f3(coor, axis); + } else { + double3 unSlant = make_double3(0.0, 0.0, 0.0); + TransformUnSlantGPU(unSlant, coor, gpu_Invcell_x, gpu_Invcell_y, gpu_Invcell_z); + WrapPBC_f3(unSlant, axis); + TransformSlantGPU(coor, unSlant, gpu_cell_x, gpu_cell_y, gpu_cell_z); + } + // update the new position + gpu_x[atomIdx] = coor.x; + gpu_y[atomIdx] = coor.y; + gpu_z[atomIdx] = coor.z; + } +} + +void BrownianMotionTranslateParticlesGPU( + VariablesCUDA *vars, + const std::vector &moleculeInvolved, + XYZArray &mForce, + XYZArray &mForceRec, + XYZArray &newMolPos, + XYZArray &newCOMs, + XYZArray &t_k, + const XYZ &boxAxes, + const double BETA, + const double t_max, + unsigned int step, + unsigned int seed, + const int box, + bool isOrthogonal) +{ + int atomCount = newMolPos.Count(); + int molCount = newCOMs.Count(); + int molCountInBox = moleculeInvolved.size(); + int *gpu_moleculeInvolved; + // Each block would handle one molecule + int threadsPerBlock = 64; + int blocksPerGrid = molCountInBox; + + CUMALLOC((void **) &gpu_moleculeInvolved, molCountInBox * sizeof(int)); + + cudaMemcpy(vars->gpu_mForcex, mForce.x, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_mForcey, mForce.y, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_mForcez, mForce.z, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_mForceRecx, mForceRec.x, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_mForceRecy, mForceRec.y, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_mForceRecz, mForceRec.z, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_x, newMolPos.x, atomCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_y, newMolPos.y, atomCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_z, newMolPos.z, atomCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_comx, newCOMs.x, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_comy, newCOMs.y, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(vars->gpu_comz, newCOMs.z, molCount * sizeof(double), cudaMemcpyHostToDevice); + cudaMemcpy(gpu_moleculeInvolved, &moleculeInvolved[0], molCountInBox * sizeof(int), cudaMemcpyHostToDevice); + + double3 axis = make_double3(boxAxes.x, boxAxes.y, boxAxes.z); + double3 halfAx = make_double3(boxAxes.x / 2.0, boxAxes.y / 2.0, boxAxes.z / 2.0); + + if(isOrthogonal) { + BrownianMotionTranslateKernel<<< blocksPerGrid, threadsPerBlock>>>( + vars->gpu_startAtomIdx, + vars->gpu_x, + vars->gpu_y, + vars->gpu_z, + vars->gpu_mForcex, + vars->gpu_mForcey, + vars->gpu_mForcez, + vars->gpu_mForceRecx, + vars->gpu_mForceRecy, + vars->gpu_mForceRecz, + vars->gpu_comx, + vars->gpu_comy, + vars->gpu_comz, + vars->gpu_t_k_x, + vars->gpu_t_k_y, + vars->gpu_t_k_z, + gpu_moleculeInvolved, + vars->gpu_cell_x[box], + vars->gpu_cell_y[box], + vars->gpu_cell_z[box], + vars->gpu_Invcell_x[box], + vars->gpu_Invcell_y[box], + vars->gpu_Invcell_z[box], + axis, + halfAx, + atomCount, + t_max, + step, + seed, + BETA); + } else { + BrownianMotionTranslateKernel<<< blocksPerGrid, threadsPerBlock>>>( + vars->gpu_startAtomIdx, + vars->gpu_x, + vars->gpu_y, + vars->gpu_z, + vars->gpu_mForcex, + vars->gpu_mForcey, + vars->gpu_mForcez, + vars->gpu_mForceRecx, + vars->gpu_mForceRecy, + vars->gpu_mForceRecz, + vars->gpu_comx, + vars->gpu_comy, + vars->gpu_comz, + vars->gpu_t_k_x, + vars->gpu_t_k_y, + vars->gpu_t_k_z, + gpu_moleculeInvolved, + vars->gpu_cell_x[box], + vars->gpu_cell_y[box], + vars->gpu_cell_z[box], + vars->gpu_Invcell_x[box], + vars->gpu_Invcell_y[box], + vars->gpu_Invcell_z[box], + axis, + halfAx, + atomCount, + t_max, + step, + seed, + BETA); + } + + cudaDeviceSynchronize(); + checkLastErrorCUDA(__FILE__, __LINE__); + + cudaMemcpy(newMolPos.x, vars->gpu_x, atomCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(newMolPos.y, vars->gpu_y, atomCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(newMolPos.z, vars->gpu_z, atomCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(newCOMs.x, vars->gpu_comx, molCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(newCOMs.y, vars->gpu_comy, molCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(newCOMs.z, vars->gpu_comz, molCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(t_k.x, vars->gpu_t_k_x, molCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(t_k.y, vars->gpu_t_k_y, molCount * sizeof(double), cudaMemcpyDeviceToHost); + cudaMemcpy(t_k.z, vars->gpu_t_k_z, molCount * sizeof(double), cudaMemcpyDeviceToHost); + CUFREE(gpu_moleculeInvolved); + checkLastErrorCUDA(__FILE__, __LINE__); +} + + +template +__global__ void BrownianMotionTranslateKernel( + int *startAtomIdx, + double *gpu_x, + double *gpu_y, + double *gpu_z, + double *molForcex, + double *molForcey, + double *molForcez, + double *molForceRecx, + double *molForceRecy, + double *molForceRecz, + double *gpu_comx, + double *gpu_comy, + double *gpu_comz, + double *gpu_t_k_x, + double *gpu_t_k_y, + double *gpu_t_k_z, + int *moleculeInvolved, + double *gpu_cell_x, + double *gpu_cell_y, + double *gpu_cell_z, + double *gpu_Invcell_x, + double *gpu_Invcell_y, + double *gpu_Invcell_z, + double3 axis, + double3 halfAx, + int atomCount, + double t_max, + unsigned int step, + unsigned int seed, + double BETA) +{ + //Each grid take cares of one molecule + int molIndex = moleculeInvolved[blockIdx.x]; + int startIdx = startAtomIdx[molIndex]; + int endIdx = startAtomIdx[molIndex + 1]; + int atomIdx; + + __shared__ double3 shift; + + // thread 0 will calculate the shift vector and update COM and gpu_t_k + if(threadIdx.x == 0) { + double3 com = make_double3(gpu_comx[molIndex], gpu_comy[molIndex], gpu_comz[molIndex]); + // This section calculates the amount of shift + double stdDev = sqrt(2.0 * t_max); + double bfm_x = (molForcex[molIndex] + molForceRecx[molIndex]) * BETA * t_max; + double bfm_y = (molForcey[molIndex] + molForceRecy[molIndex]) * BETA * t_max; + double bfm_z = (molForcez[molIndex] + molForceRecz[molIndex]) * BETA * t_max; + + shift.x = bfm_x + randomGaussianGPU(molIndex * 3, step, seed, 0.0, stdDev); + shift.y = bfm_y + randomGaussianGPU(molIndex * 3 + 1, step, seed, 0.0, stdDev); + shift.z = bfm_z + randomGaussianGPU(molIndex * 3 + 2, step, seed, 0.0, stdDev); + // update the trial translate + gpu_t_k_x[molIndex] = shift.x; + gpu_t_k_y[molIndex] = shift.y; + gpu_t_k_z[molIndex] = shift.z; + // shift COM + com.x += shift.x; + com.y += shift.y; + com.z += shift.z; + // wrap COM + if(isOrthogonal) { + WrapPBC_f3(com, axis); + } else { + double3 unSlant = make_double3(0.0, 0.0, 0.0); + TransformUnSlantGPU(unSlant, com, gpu_Invcell_x, gpu_Invcell_y, gpu_Invcell_z); + WrapPBC_f3(unSlant, axis); + TransformSlantGPU(com, unSlant, gpu_cell_x, gpu_cell_y, gpu_cell_z); + } + //update COM + gpu_comx[molIndex] = com.x; + gpu_comy[molIndex] = com.y; + gpu_comz[molIndex] = com.z; + } + + __syncthreads(); + // use strid of blockDim.x, which is 64 + // each thread handles one atom translation + for(atomIdx = startIdx + threadIdx.x; atomIdx < endIdx; atomIdx += blockDim.x) { + double3 coor = make_double3(gpu_x[atomIdx], gpu_y[atomIdx], gpu_z[atomIdx]); + + // translate the atom + coor.x += shift.x; + coor.y += shift.y; + coor.z += shift.z; + // wrap coordinate + if(isOrthogonal) { + WrapPBC_f3(coor, axis); + } else { + double3 unSlant = make_double3(0.0, 0.0, 0.0); + TransformUnSlantGPU(unSlant, coor, gpu_Invcell_x, gpu_Invcell_y, gpu_Invcell_z); + WrapPBC_f3(unSlant, axis); + TransformSlantGPU(coor, unSlant, gpu_cell_x, gpu_cell_y, gpu_cell_z); + } + // update the new position + gpu_x[atomIdx] = coor.x; + gpu_y[atomIdx] = coor.y; + gpu_z[atomIdx] = coor.z; + } +} + #endif diff --git a/src/GPU/TransformParticlesCUDAKernel.cuh b/src/GPU/TransformParticlesCUDAKernel.cuh index 0bba24300..ac2e9d67d 100644 --- a/src/GPU/TransformParticlesCUDAKernel.cuh +++ b/src/GPU/TransformParticlesCUDAKernel.cuh @@ -104,4 +104,103 @@ __global__ void RotateParticlesKernel(unsigned int numberOfMolecules, double *gpu_r_k_y, double *gpu_r_k_z, int8_t *gpu_isMoleculeInvolved); + +// Brownian Motion multiparticle +void BrownianMotionRotateParticlesGPU( + VariablesCUDA *vars, + const std::vector &moleculeInvolved, + XYZArray &mTorque, + XYZArray &newMolPos, + XYZArray &newCOMs, + XYZArray &r_k, + const XYZ &boxAxes, + const double BETA, + const double r_max, + unsigned int step, + unsigned int seed, + const int box, + bool isOrthogonal); + + +void BrownianMotionTranslateParticlesGPU( + VariablesCUDA *vars, + const std::vector &moleculeInvolved, + XYZArray &mForce, + XYZArray &mForceRec, + XYZArray &newMolPos, + XYZArray &newCOMs, + XYZArray &t_k, + const XYZ &boxAxes, + const double BETA, + const double t_max, + unsigned int step, + unsigned int seed, + const int box, + bool isOrthogonal); + + +template +__global__ void BrownianMotionRotateKernel( + int *startAtomIdx, + double *gpu_x, + double *gpu_y, + double *gpu_z, + double *molTorquex, + double *molTorquey, + double *molTorquez, + double *gpu_comx, + double *gpu_comy, + double *gpu_comz, + double *gpu_r_k_x, + double *gpu_r_k_y, + double *gpu_r_k_z, + int *moleculeInvolved, + double *gpu_cell_x, + double *gpu_cell_y, + double *gpu_cell_z, + double *gpu_Invcell_x, + double *gpu_Invcell_y, + double *gpu_Invcell_z, + double3 axis, + double3 halfAx, + int atomCount, + double r_max, + unsigned int step, + unsigned int seed, + double BETA); + + +template +__global__ void BrownianMotionTranslateKernel( + int *startAtomIdx, + double *gpu_x, + double *gpu_y, + double *gpu_z, + double *molForcex, + double *molForcey, + double *molForcez, + double *molForceRecx, + double *molForceRecy, + double *molForceRecz, + double *gpu_comx, + double *gpu_comy, + double *gpu_comz, + double *gpu_t_k_x, + double *gpu_t_k_y, + double *gpu_t_k_z, + int *moleculeInvolved, + double *gpu_cell_x, + double *gpu_cell_y, + double *gpu_cell_z, + double *gpu_Invcell_x, + double *gpu_Invcell_y, + double *gpu_Invcell_z, + double3 axis, + double3 halfAx, + int atomCount, + double t_max, + unsigned int step, + unsigned int seed, + double BETA); + #endif diff --git a/src/GPU/VariablesCUDA.cuh b/src/GPU/VariablesCUDA.cuh index a601d6f65..4c4568124 100644 --- a/src/GPU/VariablesCUDA.cuh +++ b/src/GPU/VariablesCUDA.cuh @@ -72,6 +72,7 @@ public: gpu_mForcex = NULL; gpu_mForcey = NULL; gpu_mForcez = NULL; + gpu_startAtomIdx = NULL; // setting lambda values to null gpu_molIndex = NULL; @@ -86,6 +87,7 @@ public: int *gpu_VDW_Kind; int *gpu_isMartini; int *gpu_count; + int *gpu_startAtomIdx; //start atom index of the molecule double *gpu_rCut; double *gpu_rCutCoulomb; double *gpu_rCutLow; diff --git a/src/MoleculeLookup.cpp b/src/MoleculeLookup.cpp index dedcd1e6a..c7764e2bf 100644 --- a/src/MoleculeLookup.cpp +++ b/src/MoleculeLookup.cpp @@ -12,9 +12,16 @@ along with this program, also can be found at . #include #include #include +#ifdef GOMC_CUDA +#include +#include +#include "CUDAMemoryManager.cuh" +#include "VariablesCUDA.cuh" +#endif void MoleculeLookup::Init(const Molecules& mols, - const pdb_setup::Atoms& atomData) + const pdb_setup::Atoms& atomData, + Forcefield &ff) { numKinds = mols.GetKindsCount(); @@ -97,6 +104,17 @@ void MoleculeLookup::Init(const Molecules& mols, } } boxAndKindStart[numKinds * BOX_TOTAL] = mols.count; + +// allocate and set gpu variables +#ifdef GOMC_CUDA + VariablesCUDA *cudaVars = ff.particles->getCUDAVars(); + int numMol = mols.count + 1; + // allocate memory to store molecule start atom index + CUMALLOC((void**) &cudaVars->gpu_startAtomIdx, numMol * sizeof(int)); + // copy start atom index + cudaMemcpy(cudaVars->gpu_startAtomIdx, mols.start, numMol * sizeof(int), cudaMemcpyHostToDevice); +#endif + } uint MoleculeLookup::NumInBox(const uint box) const diff --git a/src/MoleculeLookup.h b/src/MoleculeLookup.h index f98f29c7c..f9085c61e 100644 --- a/src/MoleculeLookup.h +++ b/src/MoleculeLookup.h @@ -10,6 +10,7 @@ along with this program, also can be found at . #include "EnsemblePreprocessor.h" //for BOX_TOTAL #include "BasicTypes.h" //For uint #include "algorithm" +#include "Forcefield.h" #include class CheckpointOutput; @@ -44,7 +45,8 @@ class MoleculeLookup } //Initialize this object to be consistent with Molecules mols - void Init(Molecules const& mols, const pdb_setup::Atoms& atomData); + void Init(Molecules const& mols, const pdb_setup::Atoms& atomData, + Forcefield &ff); uint GetNumKind(void) const { diff --git a/src/MoveConst.h b/src/MoveConst.h index d4a3d5331..2dc1d56e7 100644 --- a/src/MoveConst.h +++ b/src/MoveConst.h @@ -29,24 +29,14 @@ const uint GEMC_NPT = 1; const uint DISPLACE = 0; const uint ROTATE = 1; const uint MULTIPARTICLE = 2; +const uint MULTIPARTICLE_BM = 3; #if ENSEMBLE == NVT -const uint INTRA_SWAP = 3; -const uint REGROWTH = 4; -const uint INTRA_MEMC = 5; -const uint CRANKSHAFT = 6; -const uint MOVE_KINDS_TOTAL = 7; +const uint INTRA_SWAP = 4; +const uint REGROWTH = 5; +const uint INTRA_MEMC = 6; +const uint CRANKSHAFT = 7; +const uint MOVE_KINDS_TOTAL = 8; #elif ENSEMBLE == GCMC -const uint INTRA_SWAP = 3; -const uint REGROWTH = 4; -const uint INTRA_MEMC = 5; -const uint CRANKSHAFT = 6; -const uint MEMC = 7; -const uint MOL_TRANSFER = 8; -const uint CFCMC = 9; -const uint TARGETED_SWAP = 10; -const uint MOVE_KINDS_TOTAL = 11; -#elif ENSEMBLE == GEMC -const uint VOL_TRANSFER = 3; const uint INTRA_SWAP = 4; const uint REGROWTH = 5; const uint INTRA_MEMC = 6; @@ -56,13 +46,24 @@ const uint MOL_TRANSFER = 9; const uint CFCMC = 10; const uint TARGETED_SWAP = 11; const uint MOVE_KINDS_TOTAL = 12; +#elif ENSEMBLE == GEMC +const uint VOL_TRANSFER = 4; +const uint INTRA_SWAP = 5; +const uint REGROWTH = 6; +const uint INTRA_MEMC = 7; +const uint CRANKSHAFT = 8; +const uint MEMC = 9; +const uint MOL_TRANSFER = 10; +const uint CFCMC = 11; +const uint TARGETED_SWAP = 12; +const uint MOVE_KINDS_TOTAL = 13; #elif ENSEMBLE == NPT -const uint VOL_TRANSFER = 3; -const uint INTRA_SWAP = 4; -const uint REGROWTH = 5; -const uint INTRA_MEMC = 6; -const uint CRANKSHAFT = 7; -const uint MOVE_KINDS_TOTAL = 8; +const uint VOL_TRANSFER = 4; +const uint INTRA_SWAP = 5; +const uint REGROWTH = 6; +const uint INTRA_MEMC = 7; +const uint CRANKSHAFT = 8; +const uint MOVE_KINDS_TOTAL = 9; #endif const uint BOX0 = 0; @@ -70,32 +71,35 @@ const uint BOX1 = 1; ////////////////////////////////////////////////////////// -//NVT : 1. Disp (box 0) 2. Rotate (box 0) 3. MultiParticle (box 0) -// 4. IntraSwap (box 0) 5. Regrowth (box 0) 6. IntraMEMC (box 0) -// 7. CrankShaft (box 0) +//NVT : 1. Disp (box 0) 2. Rotate (box 0) 3. MultiParticle (box 0) +// 4. MultiParticle_BM(box0) 5. IntraSwap (box 0) 6. Regrowth (box 0) +// 7. IntraMEMC (box 0) 8. CrankShaft (box 0) // //GCMC: 1. Disp (box 0) 2. Rotate (box 0) 3. MultiParticle (box 0) -// 4. IntraSwap (box 0) 5. Regrowth (box 0) 6. IntraMEMC (box 0) -// 7. CrankShaft (box 0) 8. MEMC (box 0) 9. Deletion (box 0) -// 10. Insertion (box 0) 11. CFCMC (box 0) 12. CFCMC (box 1) -// 13. TargetedSwap (box 0) 14. TargetedSwap (box 1) +// 4. MultiParticle_BM(box0) +// 5. IntraSwap (box 0) 6. Regrowth (box 0) 7. IntraMEMC (box 0) +// 8. CrankShaft (box 0) 9. MEMC (box 0) 10. Deletion (box 0) +// 11. Insertion (box 0) 12. CFCMC (box 0) 13. CFCMC (box 1) +// 14. TargetedSwap (box 0) 15. TargetedSwap (box 1) // //GEMC: 1. Disp (box 0) 2. Disp (box 1) // 3. MultiParticle (box 0) 4. MultiParticle (box 1) -// 5. Rotate (box 0) 6. Rotate (box 1) -// 7. Vol. (b0->b1) 8. Vol. (b1->b0) -// 9. IntraSwap (box 0) 10. IntraSwap (box 1) -// 11. Regrowth (box 0) 12. Regrowth (box 1) -// 13. IntraMEMC (box 0) 14. IntraMEMC (box 1) -// 15. CrankShaft (box 0) 16. CrankShaft (box 1) -// 17. MEMC (b0->b1) 18. MEMC (b1->b0) -// 19. Mol Trans (b0->b1) 20. Mol Trans (b1->b0) -// 21. CFCMC (b0->b1) 22. CFCMC (b1->b0) -// 23. TargetedSwap (b0->b1) 24. TargetedSwap (b1->b0) +// 5. MultiParticle_BM(box0) 6. MultiParticle_BM(box1) +// 7. Rotate (box 0) 8. Rotate (box 1) +// 9. Vol. (b0->b1) 10. Vol. (b1->b0) +// 11. IntraSwap (box 0) 12. IntraSwap (box 1) +// 13. Regrowth (box 0) 14. Regrowth (box 1) +// 15. IntraMEMC (box 0) 16. IntraMEMC (box 1) +// 17. CrankShaft (box 0) 18. CrankShaft (box 1) +// 19. MEMC (b0->b1) 20. MEMC (b1->b0) +// 21. Mol Trans (b0->b1) 22. Mol Trans (b1->b0) +// 23. CFCMC (b0->b1) 24. CFCMC (b1->b0) +// 25. TargetedSwap (b0->b1) 26. TargetedSwap (b1->b0) // //NPT : 1. Disp (box 0) 2. Rotate (box 0) 3. MultiParticle (box 0) -// 4. Vol. (box 0) 5. IntraSwap (box 0) 6. Regrowth (box 0) -// 7. IntraMEMC (box 0) 8. CrankShaft (box 0) +// 4. MultiParticle_BM (box0) +// 5. Vol. (box 0) 6. IntraSwap (box 0) 7. Regrowth (box 0) +// 8. IntraMEMC (box 0) 9. CrankShaft (box 0) //AUTO REJECTION OR ACCEPTANCE FLAGS diff --git a/src/MoveSettings.cpp b/src/MoveSettings.cpp index e6a4d63fa..58057937b 100644 --- a/src/MoveSettings.cpp +++ b/src/MoveSettings.cpp @@ -93,7 +93,7 @@ void MoveSettings::Update(const uint move, const bool isAccepted, tempAccepted[box][move][kind]++; accepted[box][move][kind]++; - if(move != mv::MULTIPARTICLE) { + if(move != mv::MULTIPARTICLE || move != mv::MULTIPARTICLE_BM) { SetSingleMoveAccepted(box); #if ENSEMBLE == GEMC //GEMC has multiple boxes and this move changed both boxes, so we @@ -109,7 +109,7 @@ void MoveSettings::Update(const uint move, const bool isAccepted, } } - if(move == mv::MULTIPARTICLE) + if(move == mv::MULTIPARTICLE || move == mv::MULTIPARTICLE_BM) UnsetSingleMoveAccepted(box); acceptPercent[box][move][kind] = (double)(accepted[box][move][kind]) / @@ -117,7 +117,7 @@ void MoveSettings::Update(const uint move, const bool isAccepted, //for any move that we don't care about kind of molecule, it should be included //in the if condition - if (move == mv::INTRA_MEMC || move == mv::MULTIPARTICLE + if (move == mv::INTRA_MEMC || move == mv::MULTIPARTICLE || move == mv::MULTIPARTICLE_BM #if ENSEMBLE == GEMC || ENSEMBLE == GCMC || move == mv::MEMC #endif @@ -255,7 +255,7 @@ uint MoveSettings::GetAcceptTot(const uint box, const uint move) const sum += accepted[box][move][k]; } - if(move == mv::INTRA_MEMC || move == mv::MULTIPARTICLE + if(move == mv::INTRA_MEMC || move == mv::MULTIPARTICLE || move == mv::MULTIPARTICLE_BM #if ENSEMBLE == GEMC || ENSEMBLE == GCMC || move == mv::MEMC #endif @@ -284,7 +284,7 @@ uint MoveSettings::GetTrialTot(const uint box, const uint move) const sum += tries[box][move][k]; } - if(move == mv::INTRA_MEMC || move == mv::MULTIPARTICLE + if(move == mv::INTRA_MEMC || move == mv::MULTIPARTICLE || move == mv::MULTIPARTICLE_BM #if ENSEMBLE == GEMC || ENSEMBLE == GCMC || move == mv::MEMC #endif diff --git a/src/PRNG.h b/src/PRNG.h index 57729f165..227cd75ed 100644 --- a/src/PRNG.h +++ b/src/PRNG.h @@ -31,7 +31,7 @@ along with this program, also can be found at . class PRNG { public: - PRNG(MoleculeLookup & molLook) : molLookRef(molLook) {} + PRNG(MoleculeLookup & molLook) : molLookRef(molLook), hasSecondGaussian(false) {} ~PRNG(void) { delete gen; @@ -517,6 +517,28 @@ class PRNG return PickMol(m, mk, b, subDraw, boxDiv); } + //draw a uniform distribution with average mean + // and standard deviation stdDev + double Gaussian(double mean, double stdDev) + { + if(hasSecondGaussian) { + hasSecondGaussian = false; + return (mean + secondGaussian * stdDev); + } else { + double r, v1, v2, factor; + do { + v1 = Sym(1.0); + v2 = Sym(1.0); + r = v1*v1 + v2*v2; + } while (r >= 1.0 || r < 1.523e-8); + + factor = sqrt(-2.0 * log(r) / r); + hasSecondGaussian = true; + secondGaussian = v1 * factor; + return (mean + v2 * factor * stdDev); + } + } + MTRand * GetGenerator() { return gen; @@ -525,6 +547,8 @@ class PRNG private: MTRand * gen; MoleculeLookup & molLookRef; + bool hasSecondGaussian; + double secondGaussian; }; //Saves the current state of the PRNG as ./filename diff --git a/src/Random123Wrapper.cpp b/src/Random123Wrapper.cpp new file mode 100644 index 000000000..e3ffb65ce --- /dev/null +++ b/src/Random123Wrapper.cpp @@ -0,0 +1,58 @@ +#include "Random123Wrapper.h" +#include "Random123/boxmuller.hpp" +#include + +Random123Wrapper::Random123Wrapper() +{ + c = {{}}; + uk = {{}}; +} + +void Random123Wrapper::SetStep(unsigned int step) +{ + uk[0] = step; +} +void Random123Wrapper::SetRandomSeed(unsigned int seedValue) +{ + uk[1] = seedValue; +} + +double Random123Wrapper::GetRandomNumber(unsigned int counter) +{ + c[0] = counter; + RNG::key_type k = uk; + RNG::ctr_type r = rng(c, k); + double r01 = r[0]; + r01 /= UINT_MAX; + return r01; +} + +unsigned int Random123Wrapper::GetStep() +{ + return uk[0]; +} + +unsigned int Random123Wrapper::GetSeedValue() +{ + return uk[1]; +} + +double Random123Wrapper::operator() (unsigned int counter) +{ + return GetRandomNumber(counter); +} + +double Random123Wrapper::GetGaussian(unsigned int counter) +{ + c[0] = counter; + RNG::key_type k = uk; + RNG::ctr_type r = rng(c, k); + r123::float2 normalf2 = r123::boxmuller(r[0], r[1]); + return double(normalf2.x); +} + +double Random123Wrapper::GetGaussianNumber(unsigned int counter, double mean, double stdDev) +{ + double gNum = this->GetGaussian(counter); + return (mean + gNum * stdDev); +} diff --git a/src/Random123Wrapper.h b/src/Random123Wrapper.h index f1dc29687..730a919c1 100644 --- a/src/Random123Wrapper.h +++ b/src/Random123Wrapper.h @@ -6,41 +6,19 @@ typedef r123::Philox4x32 RNG; class Random123Wrapper { public: - Random123Wrapper() - { - c = {{}}; - uk = {{}}; - } + Random123Wrapper(); - void SetStep(unsigned int step) - { - uk[0] = step; - } - void SetRandomSeed(unsigned int seedValue) - { - uk[1] = seedValue; - } - double GetRandomNumber(unsigned int counter) - { - c[0] = counter; - RNG::key_type k = uk; - RNG::ctr_type r = rng(c, k); - double r01 = r[0]; - r01 /= UINT_MAX; - return r01; - } - unsigned int GetStep() - { - return uk[0]; - } - unsigned int GetSeedValue() - { - return uk[1]; - } - double operator() (unsigned int counter) - { - return GetRandomNumber(counter); - } + void SetStep(unsigned int step); + void SetRandomSeed(unsigned int seedValue); + + double GetRandomNumber(unsigned int counter); + + unsigned int GetStep(); + unsigned int GetSeedValue(); + double operator() (unsigned int counter); + + double GetGaussian(unsigned int counter); + double GetGaussianNumber(unsigned int counter, double mean, double stdDev); private: RNG::ctr_type c; diff --git a/src/StaticVals.cpp b/src/StaticVals.cpp index dd83deee6..187d82a12 100644 --- a/src/StaticVals.cpp +++ b/src/StaticVals.cpp @@ -18,7 +18,7 @@ void StaticVals::Init(Setup & set, System& sys) forcefield.Init(set); mol.Init(set, forcefield, sys); #ifndef VARIABLE_PARTICLE_NUMBER - molLookup.Init(mol, set.pdb.atoms); + molLookup.Init(mol, set.pdb.atoms, forcefield); #endif InitMovePercents(set.config.sys.moves); #if ENSEMBLE == GEMC || ENSEMBLE == NPT @@ -45,6 +45,9 @@ void StaticVals::InitMovePercents(config_setup::MovePercents const& perc) case mv::MULTIPARTICLE: movePerc[m] = perc.multiParticle; break; + case mv::MULTIPARTICLE_BM: + movePerc[m] = perc.multiParticleBrownian; + break; case mv::ROTATE: movePerc[m] = perc.rotate; break; diff --git a/src/System.cpp b/src/System.cpp index 747fedca5..e9b1dbc6c 100644 --- a/src/System.cpp +++ b/src/System.cpp @@ -23,6 +23,7 @@ along with this program, also can be found at . #include "MoleculeTransfer.h" #include "IntraSwap.h" #include "MultiParticle.h" +#include "MultiParticleBrownianMotion.h" #include "Regrowth.h" #include "MoleculeExchange1.h" #include "MoleculeExchange2.h" @@ -68,6 +69,8 @@ System::~System() delete calcEwald; delete moves[mv::DISPLACE]; delete moves[mv::ROTATE]; + delete moves[mv::MULTIPARTICLE]; + delete moves[mv::MULTIPARTICLE_BM]; delete moves[mv::INTRA_SWAP]; delete moves[mv::REGROWTH]; delete moves[mv::INTRA_MEMC]; @@ -97,7 +100,7 @@ void System::Init(Setup & set, ulong & startStep) prngParallelTemp->Init(set.prngParallelTemp.prngMaker.prng); #endif #ifdef VARIABLE_PARTICLE_NUMBER - molLookup.Init(statV.mol, set.pdb.atoms); + molLookup.Init(statV.mol, set.pdb.atoms, statV.forcefield); #endif moveSettings.Init(statV, set.pdb.remarks, molLookupRef.GetNumKind()); @@ -165,6 +168,7 @@ void System::InitMoves(Setup const& set) { moves[mv::DISPLACE] = new Translate(*this, statV); moves[mv::MULTIPARTICLE] = new MultiParticle(*this, statV); + moves[mv::MULTIPARTICLE_BM] = new MultiParticleBrownian(*this, statV); moves[mv::ROTATE] = new Rotate(*this, statV); moves[mv::INTRA_SWAP] = new IntraSwap(*this, statV); moves[mv::REGROWTH] = new Regrowth(*this, statV); @@ -258,7 +262,7 @@ void System::RecalculateTrajectory(Setup &set, uint frameNum) set.pdb.Init(set.config.in.restart, set.config.in.files.pdb.name, frameNum); statV.InitOver(set, *this); #ifdef VARIABLE_PARTICLE_NUMBER - molLookup.Init(statV.mol, set.pdb.atoms); + molLookup.Init(statV.mol, set.pdb.atoms, statV.forcefield); #endif coordinates.InitFromPDB(set.pdb.atoms); com.CalcCOM(); @@ -350,6 +354,7 @@ void System::PrintTime() printf("%-36s %10.4f sec.\n", "Displacement:", moveTime[mv::DISPLACE]); printf("%-36s %10.4f sec.\n", "Rotation:", moveTime[mv::ROTATE]); printf("%-36s %10.4f sec.\n", "MultiParticle:", moveTime[mv::MULTIPARTICLE]); + printf("%-36s %10.4f sec.\n", "MultiParticle-Brownian:", moveTime[mv::MULTIPARTICLE_BM]); printf("%-36s %10.4f sec.\n", "Intra-Swap:", moveTime[mv::INTRA_SWAP]); printf("%-36s %10.4f sec.\n", "Regrowth:", moveTime[mv::REGROWTH]); printf("%-36s %10.4f sec.\n", "Intra-MEMC:", moveTime[mv::INTRA_MEMC]); diff --git a/src/moves/MultiParticle.h b/src/moves/MultiParticle.h index 23257e597..f546a5ead 100644 --- a/src/moves/MultiParticle.h +++ b/src/moves/MultiParticle.h @@ -47,6 +47,7 @@ class MultiParticle : public MoveBase Coordinates newMolsPos; COM newCOMs; int moveType; + bool allTranslate; std::vector moleculeIndex; const MoleculeLookup& molLookup; #ifdef GOMC_CUDA @@ -90,6 +91,16 @@ inline MultiParticle::MultiParticle(System &sys, StaticVals const &statV) : initMol[b] = false; } + // Check to see if we have only monoatomic molecule or not + allTranslate = false; + int numAtomsPerKind = 0; + for (int k = 0; k < molLookup.GetNumKind(); ++k) { + numAtomsPerKind += molRef.NumAtoms(k); + } + // If we have only one atom in each kind, it means all molecule + // in the system is monoatomic + allTranslate = (numAtomsPerKind == molLookup.GetNumKind()); + #ifdef GOMC_CUDA cudaVars = sys.statV.forcefield.particles->getCUDAVars(); @@ -156,17 +167,21 @@ inline uint MultiParticle::Prep(const double subDraw, const double movPerc) prng.PickBox(bPick, subDraw, movPerc); #endif - moveType = prng.randIntExc(mp::MPTOTALTYPES); + // In each step, we perform either: + // 1- All displacement move. + // 2- All rotation move. + if(allTranslate) { + moveType = mp::MPDISPLACE; + } else { + moveType = prng.randIntExc(mp::MPTOTALTYPES); + } + SetMolInBox(bPick); if (moleculeIndex.size() == 0) { std::cout << "Warning: MultiParticle move can't move any molecules, Skipping...\n"; state = mv::fail_state::NO_MOL_OF_KIND_IN_BOX; return state; } - uint length = molRef.GetKind(moleculeIndex[0]).NumAtoms(); - if(length == 1) { - moveType = mp::MPDISPLACE; - } //We don't use forces for non-MP moves, so we need to calculate them for the //current system if any other moves, besides other MP moves, have been accepted. @@ -517,6 +532,13 @@ inline void MultiParticle::RotateForceBiased(uint molIndex) inline void MultiParticle::TranslateForceBiased(uint molIndex) { XYZ shift = t_k.Get(molIndex); + if(shift > boxDimRef.GetHalfAxis(bPick)) { + std::cout << "Error: Trial Displacement exceed half of the box length in Multiparticle \n" + << " move!\n"; + std::cout << " Trial transformation vector: " << shift << std::endl; + std::cout << " Box Dimension: " << boxDimRef.GetAxis(bPick) << std::endl << std::endl; + exit(EXIT_FAILURE); + } XYZ newcom = newCOMs.Get(molIndex); uint stop, start, len; diff --git a/src/moves/MultiParticleBrownianMotion.h b/src/moves/MultiParticleBrownianMotion.h new file mode 100644 index 000000000..f766e7063 --- /dev/null +++ b/src/moves/MultiParticleBrownianMotion.h @@ -0,0 +1,530 @@ +/******************************************************************************* +GPU OPTIMIZED MONTE CARLO (GOMC) 2.50 +Copyright (C) 2018 GOMC Group +A copy of the GNU General Public License can be found in the COPYRIGHT.txt +along with this program, also can be found at . +********************************************************************************/ +#ifndef MULTIPARTICLEBROWNIANMOTION_H +#define MULTIPARTICLEBROWNIANMOTION_H + +#include "MoveBase.h" +#include "System.h" +#include "StaticVals.h" +#include +#include "Random123Wrapper.h" +#ifdef GOMC_CUDA +#include "TransformParticlesCUDAKernel.cuh" +#include "VariablesCUDA.cuh" +#endif + +class MultiParticleBrownian : public MoveBase +{ +public: + MultiParticleBrownian(System &sys, StaticVals const& statV); + + virtual uint Prep(const double subDraw, const double movPerc); + virtual void CalcEn(); + virtual uint Transform(); + virtual void Accept(const uint rejectState, const uint step); + virtual void PrintAcceptKind(); + //used in CFCMC for initialization + void PrepCFCMC(const uint box); + +private: + uint bPick; + bool initMol; + SystemPotential sysPotNew; + XYZArray molTorqueRef; + XYZArray molTorqueNew; + XYZArray atomForceRecNew; + XYZArray molForceRecNew; + XYZArray t_k; + XYZArray r_k; + Coordinates newMolsPos; + COM newCOMs; + int moveType; + std::vector moleculeIndex; + const MoleculeLookup& molLookup; + Random123Wrapper &r123wrapper; + bool allTranslate; +#ifdef GOMC_CUDA + VariablesCUDA *cudaVars; + bool isOrthogonal; +#endif + + double GetCoeff(); + void CalculateTrialDistRot(); + void RotateForceBiased(uint molIndex); + void TranslateForceBiased(uint molIndex); + void SetMolInBox(uint box); + XYZ CalcRandomTransform(XYZ const &lb, double const max, uint molIndex); + double CalculateWRatio(XYZ const &lb_new, XYZ const &lb_old, XYZ const &k, + double max4); +}; + +inline MultiParticleBrownian::MultiParticleBrownian(System &sys, StaticVals const &statV) : + MoveBase(sys, statV), + newMolsPos(sys.boxDimRef, newCOMs, sys.molLookupRef, sys.prng, statV.mol), + newCOMs(sys.boxDimRef, newMolsPos, sys.molLookupRef, statV.mol), + molLookup(sys.molLookup), r123wrapper(sys.r123wrapper) +{ + molTorqueNew.Init(sys.com.Count()); + molTorqueRef.Init(sys.com.Count()); + atomForceRecNew.Init(sys.coordinates.Count()); + molForceRecNew.Init(sys.com.Count()); + + t_k.Init(sys.com.Count()); + r_k.Init(sys.com.Count()); + newMolsPos.Init(sys.coordinates.Count()); + newCOMs.Init(sys.com.Count()); + + initMol = false; + + // Check to see if we have only monoatomic molecule or not + allTranslate = false; + int numAtomsPerKind = 0; + for (int k = 0; k < molLookup.GetNumKind(); ++k) { + numAtomsPerKind += molRef.NumAtoms(k); + } + // If we have only one atom in each kind, it means all molecule + // in the system is monoatomic + allTranslate = (numAtomsPerKind == molLookup.GetNumKind()); + +#ifdef GOMC_CUDA + cudaVars = sys.statV.forcefield.particles->getCUDAVars(); + isOrthogonal = statV.isOrthogonal; +#endif +} + +inline void MultiParticleBrownian::PrintAcceptKind() +{ + printf("%-37s", "% Accepted MultiParticle-Brownian "); + for(uint b = 0; b < BOX_TOTAL; b++) { + printf("%10.5f ", 100.0 * moveSetRef.GetAccept(b, mv::MULTIPARTICLE_BM)); + } + std::cout << std::endl; +} + + +inline void MultiParticleBrownian::SetMolInBox(uint box) +{ + // NEED to check if atom is not fixed! +#if ENSEMBLE == GCMC || ENSEMBLE == GEMC + // need to be initialized for every move since number of atom is changing + moleculeIndex.clear(); + MoleculeLookup::box_iterator thisMol = molLookup.BoxBegin(box); + MoleculeLookup::box_iterator end = molLookup.BoxEnd(box); + while(thisMol != end) { + //Make sure this molecule is not fixed in its position + if(!molLookup.IsFix(*thisMol)) { + moleculeIndex.push_back(*thisMol); + } + thisMol++; + } +#else + // box would be always 0 in NVT or NPT ensemble, initialize it once + if(!initMol) { + moleculeIndex.clear(); + MoleculeLookup::box_iterator thisMol = molLookup.BoxBegin(box); + MoleculeLookup::box_iterator end = molLookup.BoxEnd(box); + while(thisMol != end) { + //Make sure this molecule is not fixed in its position + if(!molLookup.IsFix(*thisMol)) { + moleculeIndex.push_back(*thisMol); + } + thisMol++; + } + initMol = true; + } +#endif +} + +inline uint MultiParticleBrownian::Prep(const double subDraw, const double movPerc) +{ + GOMC_EVENT_START(1, GomcProfileEvent::PREP_MULTIPARTICLE_BM); + uint state = mv::fail_state::NO_FAIL; +#if ENSEMBLE == GCMC + bPick = mv::BOX0; +#else + prng.PickBox(bPick, subDraw, movPerc); +#endif + + // In each step, we perform either: + // 1- All displacement move. + // 2- All rotation move. + if(allTranslate) { + moveType = mp::MPDISPLACE; + } else { + moveType = prng.randIntExc(mp::MPTOTALTYPES); + } + + SetMolInBox(bPick); + if (moleculeIndex.size() == 0) { + std::cout << "Warning: MultiParticleBrownian move can't move any molecules, Skipping...\n"; + state = mv::fail_state::NO_MOL_OF_KIND_IN_BOX; + return state; + } + + //We don't use forces for non-MP moves, so we need to calculate them for the + //current system if any other moves, besides other MP moves, have been accepted. + //Or, if this is the first MP move, which is handled with the same flag. + if(moveSetRef.GetSingleMoveAccepted(bPick)) { + GOMC_EVENT_START(1, GomcProfileEvent::CALC_EN_MULTIPARTICLE_BM); + //Copy ref reciprocal terms to new for calculation with old positions + calcEwald->CopyRecip(bPick); + + //Calculate force for long range electrostatic using old position + calcEwald->BoxForceReciprocal(coordCurrRef, atomForceRecRef, molForceRecRef, bPick); + + //calculate short range energy and force for old positions + calcEnRef.BoxForce(sysPotRef, coordCurrRef, atomForceRef, molForceRef, + boxDimRef, bPick); + + if(moveType == mp::MPROTATE) { + //Calculate Torque for old positions + calcEnRef.CalculateTorque(moleculeIndex, coordCurrRef, comCurrRef, + atomForceRef, atomForceRecRef, molTorqueRef, bPick); + } + + sysPotRef.Total(); + GOMC_EVENT_STOP(1, GomcProfileEvent::CALC_EN_MULTIPARTICLE_BM); + } + coordCurrRef.CopyRange(newMolsPos, 0, 0, coordCurrRef.Count()); + comCurrRef.CopyRange(newCOMs, 0, 0, comCurrRef.Count()); + GOMC_EVENT_STOP(1, GomcProfileEvent::PREP_MULTIPARTICLE_BM); + return state; +} + +inline void MultiParticleBrownian::PrepCFCMC(const uint box) +{ + bPick = box; + moveType = prng.randIntExc(mp::MPTOTALTYPES); + SetMolInBox(bPick); + uint length = molRef.GetKind(moleculeIndex[0]).NumAtoms(); + if(length == 1) { + moveType = mp::MPDISPLACE; + } + + //We don't use forces for non-MP moves, so we need to calculate them for the + //current system if any other moves, besides other MP moves, have been accepted. + //Or, if this is the first MP move, which is handled with the same flag. + if(moveSetRef.GetSingleMoveAccepted(bPick)) { + //Copy ref reciprocal terms to new for calculation with old positions + calcEwald->CopyRecip(bPick); + + //Calculate long range electrostatic force for old positions + calcEwald->BoxForceReciprocal(coordCurrRef, atomForceRecRef, molForceRecRef, bPick); + + //Calculate short range energy and force for old positions + calcEnRef.BoxForce(sysPotRef, coordCurrRef, atomForceRef, molForceRef, + boxDimRef, bPick); + + if(moveType == mp::MPROTATE) { + //Calculate Torque for old positions + calcEnRef.CalculateTorque(moleculeIndex, coordCurrRef, comCurrRef, + atomForceRef, atomForceRecRef, molTorqueRef, bPick); + } + + sysPotRef.Total(); + } + coordCurrRef.CopyRange(newMolsPos, 0, 0, coordCurrRef.Count()); + comCurrRef.CopyRange(newCOMs, 0, 0, comCurrRef.Count()); +} + +inline uint MultiParticleBrownian::Transform() +{ + GOMC_EVENT_START(1, GomcProfileEvent::TRANS_MULTIPARTICLE_BM); + // Based on the reference force decided whether to displace or rotate each + // individual particle. + uint state = mv::fail_state::NO_FAIL; + +#ifdef GOMC_CUDA + // This kernel will calculate translation/rotation amount + shifting/rotating + if(moveType == mp::MPROTATE) { + double r_max = moveSetRef.GetRMAX(bPick); + BrownianMotionRotateParticlesGPU( + cudaVars, + moleculeIndex, + molTorqueRef, + newMolsPos, + newCOMs, + r_k, + boxDimRef.GetAxis(bPick), + BETA, + r_max, + r123wrapper.GetStep(), + r123wrapper.GetSeedValue(), + bPick, + isOrthogonal); + + } else { + double t_max = moveSetRef.GetTMAX(bPick); + BrownianMotionTranslateParticlesGPU( + cudaVars, + moleculeIndex, + molForceRef, + molForceRecRef, + newMolsPos, + newCOMs, + t_k, + boxDimRef.GetAxis(bPick), + BETA, + t_max, + r123wrapper.GetStep(), + r123wrapper.GetSeedValue(), + bPick, + isOrthogonal); + + } +#else + // Calculate trial translate and rotate + // move particles according to force and torque and store them in the new pos + CalculateTrialDistRot(); +#endif + GOMC_EVENT_STOP(1, GomcProfileEvent::TRANS_MULTIPARTICLE_BM); + return state; +} + +inline void MultiParticleBrownian::CalcEn() +{ + GOMC_EVENT_START(1, GomcProfileEvent::CALC_EN_MULTIPARTICLE_BM); + // Calculate the new force and energy and we will compare that to the + // reference values in Accept() function + cellList.GridAll(boxDimRef, newMolsPos, molLookup); + + //back up cached fourier term + calcEwald->backupMolCache(); + + //setup reciprocate vectors for new positions + calcEwald->BoxReciprocalSetup(bPick, newMolsPos); + + sysPotNew = sysPotRef; + //calculate short range energy and force + sysPotNew = calcEnRef.BoxForce(sysPotNew, newMolsPos, atomForceNew, + molForceNew, boxDimRef, bPick); + //calculate long range of new electrostatic energy + sysPotNew.boxEnergy[bPick].recip = calcEwald->BoxReciprocal(bPick, false); + //Calculate long range of new electrostatic force + calcEwald->BoxForceReciprocal(newMolsPos, atomForceRecNew, molForceRecNew, + bPick); + + if(moveType == mp::MPROTATE) { + //Calculate Torque for new positions + calcEnRef.CalculateTorque(moleculeIndex, newMolsPos, newCOMs, atomForceNew, + atomForceRecNew, molTorqueNew, bPick); + } + sysPotNew.Total(); + GOMC_EVENT_STOP(1, GomcProfileEvent::CALC_EN_MULTIPARTICLE); +} + +inline double MultiParticleBrownian::CalculateWRatio(XYZ const &lb_new, XYZ const &lb_old, + XYZ const &k, double max4) +{ + double w_ratio = 0.0; + XYZ old_var = lb_old - k; + XYZ new_var = lb_new + k; + + //Note: we could factor max4 and multiply at the end, but + // for the move, where we translate and rotate all molecules, + // this method would not work. Hence, I did not factor it. + // its actually is w_ratio += -1.0* but we simplify it + w_ratio -= (new_var.LengthSq() / max4); + // its actually is w_ratio -= -1.0* but we simplify it + w_ratio += (old_var.LengthSq() / max4); + + return w_ratio; +} + +inline double MultiParticleBrownian::GetCoeff() +{ + // calculate (w_new->old/w_old->new) and return it. + double w_ratio = 0.0; + double r_max = moveSetRef.GetRMAX(bPick); + double t_max = moveSetRef.GetTMAX(bPick); + double r_max4 = 4.0 * r_max; + double t_max4 = 4.0 * t_max; + + if(moveType == mp::MPROTATE) {// rotate, +#ifdef _OPENMP + #pragma omp parallel for default(none) shared(moleculeIndex, r_max, t_max, r_max4, t_max4) reduction(+:w_ratio) +#endif + for(int m = 0; m < moleculeIndex.size(); m++) { + uint molNumber = moleculeIndex[m]; + // rotate, bt_ = BETA * force * maxForce + XYZ bt_old = molTorqueRef.Get(molNumber) * BETA * r_max; + XYZ bt_new = molTorqueNew.Get(molNumber) * BETA * r_max; + w_ratio += CalculateWRatio(bt_new, bt_old, r_k.Get(molNumber), r_max4); + } + } else {// displace +#ifdef _OPENMP + #pragma omp parallel for default(none) shared(moleculeIndex, r_max, t_max, r_max4, t_max4) reduction(+:w_ratio) +#endif + for(int m = 0; m < moleculeIndex.size(); m++) { + uint molNumber = moleculeIndex[m]; + // bf_ = BETA * torque * maxTorque + XYZ bf_old = (molForceRef.Get(molNumber) + molForceRecRef.Get(molNumber)) * + BETA * t_max; + XYZ bf_new = (molForceNew.Get(molNumber) + molForceRecNew.Get(molNumber)) * + BETA * t_max; + w_ratio += CalculateWRatio(bf_new, bf_old, t_k.Get(molNumber), t_max4); + } + } + + return w_ratio; +} + +inline void MultiParticleBrownian::Accept(const uint rejectState, const uint step) +{ + GOMC_EVENT_START(1, GomcProfileEvent::ACC_MULTIPARTICLE_BM); + // Here we compare the values of reference and trial and decide whether to + // accept or reject the move + double MPCoeff = GetCoeff(); + double accept = exp(-BETA * (sysPotNew.Total() - sysPotRef.Total()) + MPCoeff); + bool result = (rejectState == mv::fail_state::NO_FAIL) && prng() < accept; + if(result) { + sysPotRef = sysPotNew; + swap(coordCurrRef, newMolsPos); + swap(comCurrRef, newCOMs); + swap(molForceRef, molForceNew); + swap(atomForceRef, atomForceNew); + swap(molForceRecRef, molForceRecNew); + swap(atomForceRecRef, atomForceRecNew); + swap(molTorqueRef, molTorqueNew); + //update reciprocate value + calcEwald->UpdateRecip(bPick); + } else { + cellList.GridAll(boxDimRef, coordCurrRef, molLookup); + calcEwald->exgMolCache(); + } + + moveSetRef.UpdateMoveSettingMultiParticle(bPick, result, moveType); + moveSetRef.Update(mv::MULTIPARTICLE_BM, result, bPick); + GOMC_EVENT_STOP(1, GomcProfileEvent::ACC_MULTIPARTICLE_BM); +} + +inline XYZ MultiParticleBrownian::CalcRandomTransform(XYZ const &lb, double const max, uint molIndex) +{ + XYZ lbmax = lb * max; + XYZ num; + //variance is 2A according to the paper, so stdDev is sqrt(variance) + double stdDev = sqrt(2.0 * max); + + num.x = lbmax.x + r123wrapper.GetGaussianNumber(molIndex * 3 + 0, 0.0, stdDev); + num.y = lbmax.y + r123wrapper.GetGaussianNumber(molIndex * 3 + 1, 0.0, stdDev); + num.z = lbmax.z + r123wrapper.GetGaussianNumber(molIndex * 3 + 2, 0.0, stdDev); + + + if(num.Length() >= boxDimRef.axis.Min(bPick)) { + std::cout << "Trial Displacement exceeds half of the box length in Brownian Motion MultiParticle move.\n"; + std::cout << "Trial transform: " << num; + exit(EXIT_FAILURE); + } else if (!std::isfinite(num.Length())) { + std::cout << "Error: Trial transform is not a finite number in Brownian Motion Multiparticle move.\n"; + std::cout << " Trial transform: " << num; + exit(EXIT_FAILURE); + } + // We can possible bound them + return num; +} + +inline void MultiParticleBrownian::CalculateTrialDistRot() +{ + uint molIndex; + double r_max = moveSetRef.GetRMAX(bPick); + double t_max = moveSetRef.GetTMAX(bPick); + if(moveType == mp::MPROTATE) { // rotate + double *x = r_k.x; + double *y = r_k.y; + double *z = r_k.z; +#ifdef _OPENMP + #pragma omp parallel for default(none) shared(moleculeIndex, r_max, x, y, z) +#endif + for(int m = 0; m < moleculeIndex.size(); m++) { + uint molIndex = moleculeIndex[m]; + XYZ bt = molTorqueRef.Get(molIndex) * BETA; + XYZ val = CalcRandomTransform(bt, r_max, molIndex); + x[molIndex] = val.x; + y[molIndex] = val.y; + z[molIndex] = val.z; + RotateForceBiased(molIndex); + } + } else { // displace + double *x = t_k.x; + double *y = t_k.y; + double *z = t_k.z; +#ifdef _OPENMP + #pragma omp parallel for default(none) shared(moleculeIndex, t_max, x, y, z) +#endif + for(int m = 0; m < moleculeIndex.size(); m++) { + uint molIndex = moleculeIndex[m]; + XYZ bf = (molForceRef.Get(molIndex) + molForceRecRef.Get(molIndex)) * BETA; + XYZ val = CalcRandomTransform(bf, t_max, molIndex); + x[molIndex] = val.x; + y[molIndex] = val.y; + z[molIndex] = val.z; + TranslateForceBiased(molIndex); + } + } +} + +inline void MultiParticleBrownian::RotateForceBiased(uint molIndex) +{ + XYZ rot = r_k.Get(molIndex); + double rotLen = rot.Length(); + RotationMatrix matrix; + + matrix = RotationMatrix::FromAxisAngle(rotLen, rot.Normalize()); + + XYZ center = newCOMs.Get(molIndex); + uint start, stop, len; + molRef.GetRange(start, stop, len, molIndex); + + // Copy the range into temporary array + XYZArray temp(len); + newMolsPos.CopyRange(temp, start, 0, len); + boxDimRef.UnwrapPBC(temp, bPick, center); + + // Do Rotation + for(uint p = 0; p < len; p++) { + temp.Add(p, -center); + temp.Set(p, matrix.Apply(temp[p])); + temp.Add(p, center); + } + boxDimRef.WrapPBC(temp, bPick); + // Copy back the result + temp.CopyRange(newMolsPos, 0, start, len); +} + +inline void MultiParticleBrownian::TranslateForceBiased(uint molIndex) +{ + XYZ shift = t_k.Get(molIndex); + if(shift > boxDimRef.GetHalfAxis(bPick)) { + std::cout << "Error: Trial Displacement exceed half of the box length in Multiparticle \n" + << " Brownian Motion move!\n"; + std::cout << " Trial transformation vector: " << shift << std::endl; + std::cout << " Box Dimension: " << boxDimRef.GetAxis(bPick) << std::endl << std::endl; + std::cout << "This might be due to bad initial configuration, where atom of the molecules \n" + << "are too close to each other or overlaps. Please equilibrate your system using \n" + << "rigid body translation or rotation MC moves, before using the Multiparticle \n" + << "Brownian Motion move. \n\n"; + exit(EXIT_FAILURE); + } + + XYZ newcom = newCOMs.Get(molIndex); + uint stop, start, len; + molRef.GetRange(start, stop, len, molIndex); + // Copy the range into temporary array + XYZArray temp(len); + newMolsPos.CopyRange(temp, start, 0, len); + //Shift the coordinate and COM + temp.AddAll(shift); + newcom += shift; + //rewrapping + boxDimRef.WrapPBC(temp, bPick); + newcom = boxDimRef.WrapPBC(newcom, bPick); + //set the new coordinate + temp.CopyRange(newMolsPos, 0, start, len); + newCOMs.Set(molIndex, newcom); +} + +#endif