1__global__
void CalculatePairDistances(
Boxsize Box,
Atoms* System,
Atoms Old,
bool* ConsiderThisAdsorbateAtom,
size_t Molsize,
size_t skip,
size_t* indexList,
double* Distances,
size_t totalThreads,
bool UsedForAMove)
5 size_t ij = blockIdx.x * blockDim.x + threadIdx.x;
11 const Atoms Component=System[comp];
15 Atoms Adsorbate = Old;
17 Adsorbate = System[1];
18 for(
size_t j = 0; j < Molsize; j++)
20 if(!ConsiderThisAdsorbateAtom[j])
continue;
22 size_t new_j = j + skip;
24 double3 posvec = Component.pos[ij] - Adsorbate.pos[new_j];
25 PBC(posvec, Box.Cell, Box.InverseCell, Box.Cubic);
26 const double rr_dot = dot(posvec, posvec);
27 size_t WriteTo = indexList[ij * Molsize + j];
29 Distances[WriteTo] = rr_dot;
35static inline size_t MatchInteractionIndex(std::vector<int3>& DNNInteractionList,
size_t typeA,
size_t typeB)
38 for(
size_t i = 0; i < DNNInteractionList.size(); i++)
40 if(typeA == DNNInteractionList[i].x && typeB == DNNInteractionList[i].y)
49std::vector<std::vector<double>> CalculatePairDistances_CPU(
Atoms* System,
Boxsize Box,
bool* ConsiderThisAdsorbateAtom, std::vector<int3>& DNNInteractionList)
51 std::vector<std::vector<double>>CPU_Distances(DNNInteractionList.size());
53 for(
size_t i = 0; i < System[0].size; i++)
54 for(
size_t j = 0; j < System[ads_comp].Molsize; j++)
56 if(!ConsiderThisAdsorbateAtom[j])
continue;
57 size_t typeAdsorbate = System[ads_comp].Type[j];
58 size_t typeFramework = System[0].Type[i];
59 size_t InteractionIndex = MatchInteractionIndex(DNNInteractionList, typeAdsorbate, typeFramework);
60 double3 posvec = System[0].pos[i] - System[ads_comp].pos[j];
61 if(i < 3) printf(
"adsorbate xyz: %.5f %.5f %.5f\n", System[ads_comp].pos[j].x, System[ads_comp].pos[j].y, System[ads_comp].pos[j].z);
62 PBC(posvec, Box.Cell, Box.InverseCell, Box.Cubic);
63 const double rr_dot = dot(posvec, posvec);
64 CPU_Distances[InteractionIndex].push_back(rr_dot);
70void PrepareOutliersFiles()
72 std::ofstream textrestartFile{};
73 std::string dirname=
"DNN/";
74 std::string Ifname = dirname +
"/" +
"Outliers_INSERTION.data";
75 std::string Dfname = dirname +
"/" +
"Outliers_DELETION.data";
76 std::string TRname = dirname +
"/" +
"Outliers_SINGLE_PARTICLE.data";
77 std::filesystem::path cwd = std::filesystem::current_path();
79 std::filesystem::path directoryName = cwd /dirname;
80 std::filesystem::path IfileName = cwd /Ifname;
81 std::filesystem::path DfileName = cwd /Dfname;
82 std::filesystem::path TRfileName = cwd /TRname;
83 std::filesystem::create_directories(directoryName);
85 textrestartFile = std::ofstream(TRfileName, std::ios::out);
86 textrestartFile <<
"THIS FILE IS RECORDING THE DELTA ENERGIES BETWEEN THE NEW AND OLD STATES (New - Old)" <<
"\n";
87 textrestartFile <<
"x y z Type Move HostGuest_DNN Host_Guest_Classical_MINUS_Host_Guest_Classical" <<
"\n";
88 textrestartFile.close();
90 textrestartFile = std::ofstream(IfileName, std::ios::out);
91 textrestartFile <<
"THIS FILE IS RECORDING THE ENERGIES RELATED TO EITHER THE NEW/INSERTION" <<
"\n";
92 textrestartFile <<
"x y z Type Move HostGuest_DNN Host_Guest_Classical_MINUS_Host_Guest_Classical" <<
"\n";
93 textrestartFile.close();
95 textrestartFile = std::ofstream(DfileName, std::ios::out);
96 textrestartFile <<
"THIS FILE IS RECORDING THE ENERGIES RELATED TO OLD/DELETION (TAKE THE OPPOSITE)" <<
"\n";
97 textrestartFile <<
"x y z Type Move HostGuest_DNN Host_Guest_Classical_MINUS_Host_Guest_Classical" <<
"\n";
98 textrestartFile.close();
101static inline cppflow::model load_model(std::string& NAME)
103 setenv(
"CUDA_VISIBLE_DEVICES",
"", 1);
104 cppflow::model DNNModel(NAME);
105 setenv(
"CUDA_VISIBLE_DEVICES",
"1", 1);
109static inline void Prepare_Model(
Components& SystemComponents, std::string& NAME)
111 cppflow::model TempModel = load_model(NAME);
112 SystemComponents.DNNModel.push_back(TempModel);
115void Read_LCLin_Model(
Components& SystemComponents)
127 Prepare_Model(SystemComponents, SystemComponents.ModelName[ModelID]);
136 SystemComponents.DNNMinMax = ReadMinMax();
137 printf(
"MinMax size: %zu\n", SystemComponents.DNNMinMax.size());
140 PrepareOutliersFiles();
144double DNN_Evaluate(
Components& SystemComponents, std::vector<double>& Feature)
147 std::vector<double>Vector(SystemComponents.Nfeatures, 0.0);
148 for(
size_t i = 0; i < Vector.size(); i++) Vector[i] =
static_cast<double>(Feature[i]);
149 auto real_input = cppflow::tensor(Vector, {1, SystemComponents.Nfeatures});
154 auto output = SystemComponents.DNNModel[ModelID]({{SystemComponents.InputLayer[ModelID], real_input}},{
"StatefulPartitionedCall:0"});
155 asd =
static_cast<double>(output[0].get_data<
double>()[0]);
156 printf(
"REAL DNN Prediction = %.5f\n", asd);
161 bool operator()(
const double3& a,
const double3& b)
const
167static inline void NormalizeMinMax(std::vector<double>& r_features, std::vector<double2>& MinMax,
size_t start)
170 for(
size_t i = 0; i < r_features.size(); i++)
172 r_features[i] = (r_features[i] - MinMax[start + i].x) / (MinMax[start + i].y - MinMax[start + i].x);
177static inline std::vector<double> Calcuate_Feature_from_Distance(
double dist_sq,
size_t start, std::vector<double2>& MinMax)
180 double dist = sqrt(dist_sq);
181 double one_over_dist = 1.0 / dist;
182 double one_over_distsq = 1.0 / dist_sq;
184 std::vector<double>r_features({std::exp(-dist), one_over_dist, std::pow(one_over_distsq, 2), std::pow(one_over_distsq, 3), std::pow(one_over_distsq, 4), std::pow(one_over_distsq, 5)});
186 NormalizeMinMax(r_features, MinMax, start);
190std::vector<std::vector<double>> CalculatePairDistances_GPU(
Simulations& Sim,
Components& SystemComponents, std::vector<double>& Features,
int DNN_CalcType,
size_t NMol)
192 size_t NFrameworkAtoms = SystemComponents.Moleculesize[0] * SystemComponents.NumberOfMolecule_for_Component[0];
193 size_t Nblock = 0;
size_t Nthread = 0; Setup_threadblock(NFrameworkAtoms, &Nblock, &Nthread);
196 double time = omp_get_wtime();
201 CalculatePairDistances<<<Nblock, Nthread>>>(Sim.Box, Sim.d_a, Sim.New, SystemComponents.ConsiderThisAdsorbateAtom, SystemComponents.Moleculesize[ads_comp], NMol * SystemComponents.Moleculesize[ads_comp], SystemComponents.device_InverseIndexList, Sim.Blocksum, NFrameworkAtoms,
false);
206 CalculatePairDistances<<<Nblock, Nthread>>>(Sim.Box, Sim.d_a, Sim.Old, SystemComponents.ConsiderThisAdsorbateAtom, SystemComponents.Moleculesize[ads_comp], 0, SystemComponents.device_InverseIndexList, Sim.Blocksum, NFrameworkAtoms,
true);
211 CalculatePairDistances<<<Nblock, Nthread>>>(Sim.Box, Sim.d_a, Sim.New, SystemComponents.ConsiderThisAdsorbateAtom, SystemComponents.Moleculesize[ads_comp], 0, SystemComponents.device_InverseIndexList, Sim.Blocksum, NFrameworkAtoms,
true);
214 case REINSERTION_OLD:
217 CalculatePairDistances<<<Nblock, Nthread>>>(Sim.Box, Sim.d_a, Sim.Old, SystemComponents.ConsiderThisAdsorbateAtom, SystemComponents.Moleculesize[ads_comp], skip, SystemComponents.device_InverseIndexList, Sim.Blocksum, NFrameworkAtoms,
true);
220 case REINSERTION_NEW:
222 size_t skip = SystemComponents.Moleculesize[ads_comp];
223 CalculatePairDistances<<<Nblock, Nthread>>>(Sim.Box, Sim.d_a, Sim.Old, SystemComponents.ConsiderThisAdsorbateAtom, SystemComponents.Moleculesize[ads_comp], skip, SystemComponents.device_InverseIndexList, Sim.Blocksum, NFrameworkAtoms,
true);
226 case DNN_INSERTION:
case DNN_DELETION:
229 CalculatePairDistances<<<Nblock, Nthread>>>(Sim.Box, Sim.d_a, Sim.Old, SystemComponents.ConsiderThisAdsorbateAtom, SystemComponents.Moleculesize[ads_comp], skip, SystemComponents.device_InverseIndexList, Sim.Blocksum, NFrameworkAtoms,
true);
233 SystemComponents.DNNGPUTime += omp_get_wtime() - time;
235 time = omp_get_wtime();
237 size_t start = 0;
size_t end = 0;
size_t top = 9;
size_t n_features = 6;
double stdsort_time;
238 std::vector<std::vector<double>>Distances(SystemComponents.DNNInteractionList.size());
239 std::vector<std::vector<double>>SORTED_Dist(SystemComponents.DNNInteractionList.size(), std::vector<double>(top));
240 stdsort_time = omp_get_wtime();
242 std::vector<size_t> starts;
size_t start_sum = 0;
243 for(
size_t i = 0; i < SystemComponents.DNNInteractionList.size(); i++)
245 starts.push_back(start_sum);
246 start_sum += SystemComponents.IndexList[i].size();
247 if(SystemComponents.CURRENTCYCLE == 10) printf(
"number of distances [%zu]: %zu\n", i, SystemComponents.IndexList[i].size());
252 for(
size_t i = 0; i < SystemComponents.DNNInteractionList.size(); i++)
254 Distances[i].resize(SystemComponents.IndexList[i].size()); end += SystemComponents.IndexList[i].size();
255 cudaMemcpyAsync(Distances[i].data(), &Sim.Blocksum[starts[i]], SystemComponents.IndexList[i].size()*
sizeof(
double), cudaMemcpyDeviceToHost);
258 std::sort(std::execution::unseq, Distances[i].begin(), Distances[i].end());
260 SystemComponents.DNNstdsortTime += omp_get_wtime() - stdsort_time;
262 double FeaturizationTime = omp_get_wtime();
264 #pragma omp parallel for
265 for(
size_t ij = 0; ij < SystemComponents.DNNInteractionList.size() * top; ij++)
269 SORTED_Dist[i][j] = Distances[i][j];
271 start = (i * top + j) * n_features;
272 std::vector<double> r_feature = Calcuate_Feature_from_Distance(SORTED_Dist[i][j], start, SystemComponents.DNNMinMax);
273 for(
size_t k = 0; k < n_features; k++)
274 Features[(i * top + j) * n_features + k] = r_feature[k];
278 SystemComponents.DNNFeaturizationTime += omp_get_wtime() - FeaturizationTime;
279 SystemComponents.DNNSortTime += omp_get_wtime() - time;
289 printf(
"-------------- Preparing DNN Interaction Types --------------\n");
290 std::vector<size_t>AdsorbateAtomTypesForDNN;
291 for(
size_t i = 0; i < HostSystem[ads_comp].Molsize; i++)
292 if(SystemComponents.ConsiderThisAdsorbateAtom[i])
293 AdsorbateAtomTypesForDNN.push_back(HostSystem[ads_comp].Type[i]);
295 std::sort(AdsorbateAtomTypesForDNN.begin(), AdsorbateAtomTypesForDNN.end());
296 AdsorbateAtomTypesForDNN.erase(std::unique(AdsorbateAtomTypesForDNN.begin(), AdsorbateAtomTypesForDNN.end()), AdsorbateAtomTypesForDNN.end());
297 for(
size_t i = 0; i < AdsorbateAtomTypesForDNN.size(); i++) printf(
"AdsorbateDNN Types %zu\n", AdsorbateAtomTypesForDNN[i]);
300 std::vector<size_t> FrameworkAtomTypesForDNN = {0, 1, 2, 3};
305 for(
size_t i = 0; i < AdsorbateAtomTypesForDNN.size(); i++)
307 size_t ads_type = AdsorbateAtomTypesForDNN[i];
308 for(
size_t j = 0; j < FrameworkAtomTypesForDNN.size(); j++)
310 size_t framework_type = FrameworkAtomTypesForDNN[j];
311 if(ads_type == framework_type)
312 throw std::runtime_error(
"Adsorbate Type equals Framework Type for Pseudo Atoms, Weird");
314 SystemComponents.DNNInteractionList.push_back({(int)ads_type, (
int)framework_type, (int)SystemComponents.NumberOfPseudoAtoms[framework_type]});
315 printf(
"TypeA-B [%zu-%zu], Number: %zu\n", ads_type, framework_type, SystemComponents.NumberOfPseudoAtoms[framework_type]);
322 size_t Ntypes = SystemComponents.DNNInteractionList.size();
size_t Ndistances = 9;
size_t Ninteractions = 6;
323 SystemComponents.Nfeatures = Ntypes * Ndistances * Ninteractions;
326 SystemComponents.IndexList.resize(SystemComponents.DNNInteractionList.size());
330 for(
size_t i = 0; i < HostSystem[ads_comp].Molsize; i++)
331 printf(
"ConsiderThisATom? %d\n", SystemComponents.ConsiderThisAdsorbateAtom[i]);
334 for(
size_t i = 0; i < HostSystem[0].size; i++)
335 for(
size_t j = 0; j < HostSystem[ads_comp].Molsize; j++)
338 if(!SystemComponents.ConsiderThisAdsorbateAtom[j])
continue;
339 size_t typeAdsorbate = HostSystem[ads_comp].Type[j];
340 size_t typeFramework = HostSystem[0].Type[i];
341 size_t InteractionIndex = MatchInteractionIndex(SystemComponents.DNNInteractionList, typeAdsorbate, typeFramework);
343 printf(
"count: %zu, TypeA-B: [%zu-%zu], InteractionIndex: %zu\n", i * HostSystem[ads_comp].Molsize + j, typeAdsorbate, typeFramework, InteractionIndex);
344 SystemComponents.IndexList[InteractionIndex].push_back(i * HostSystem[ads_comp].Molsize + j);
348 std::vector<size_t>FlatIndexList;
349 for(
size_t i = 0; i < SystemComponents.IndexList.size(); i++)
351 printf(
"Interaction [%zu], Amount [%zu]\n", i, SystemComponents.IndexList[i].size());
352 for(
size_t j = 0; j < SystemComponents.IndexList[i].size(); j++)
353 FlatIndexList.push_back(SystemComponents.IndexList[i][j]);
355 size_t Listsize = HostSystem[0].size * HostSystem[ads_comp].Molsize;
356 std::vector<size_t>InverseIndexList(Listsize);
357 for(
size_t i = 0; i < FlatIndexList.size(); i++)
358 InverseIndexList[FlatIndexList[i]] = i;
364 for(
size_t i = 0; i < HostSystem[0].size; i++)
365 for(
size_t j = 0; j < HostSystem[ads_comp].Molsize; j++)
367 if(!SystemComponents.ConsiderThisAdsorbateAtom[j])
continue;
369 printf(
"test_count: %zu, TypeA-B: [%zu-%zu], Where it is stored: %zu\n", i * HostSystem[ads_comp].Molsize + j, HostSystem[ads_comp].Type[j], HostSystem[0].Type[i], InverseIndexList[i * HostSystem[ads_comp].Molsize + j]);
376 cudaMalloc(&SystemComponents.device_InverseIndexList, Listsize *
sizeof(
size_t));
377 cudaMemcpy(SystemComponents.device_InverseIndexList, InverseIndexList.data(), Listsize *
sizeof(
size_t), cudaMemcpyHostToDevice);
380 printf(
"Listsize: %zu\n", Listsize);
383 printf(
"Size of the device Distance list: %zu\n", Listsize);
386static inline void WriteDistances(
Components& SystemComponents, std::vector<std::vector<double>> Distances)
389 std::ofstream textrestartFile{};
390 std::string dirname=
"DNN/";
391 std::string fname = dirname +
"/" +
"Distances_squared.data";
392 std::filesystem::path cwd = std::filesystem::current_path();
394 std::filesystem::path directoryName = cwd /dirname;
395 std::filesystem::path fileName = cwd /fname;
396 std::filesystem::create_directories(directoryName);
398 textrestartFile = std::ofstream(fileName, std::ios::out);
399 textrestartFile <<
"TypeA(Framework) TypeB(Adsorbate) distance (GPU)" <<
"\n";
400 for(
size_t i = 0; i < SystemComponents.DNNInteractionList.size(); i++)
403 for(
size_t j = 0; j < Distances[i].size(); j++)
406 textrestartFile << Distances[i][j] <<
"\n";
409 textrestartFile.close();
412static inline void WriteFeatures(
Components& SystemComponents, std::vector<double> Features)
414 std::ofstream textrestartFile{};
415 std::string dirname=
"DNN/";
416 std::string fname = dirname +
"/" +
"Features.data";
417 std::filesystem::path cwd = std::filesystem::current_path();
419 std::filesystem::path directoryName = cwd /dirname;
420 std::filesystem::path fileName = cwd /fname;
421 std::filesystem::create_directories(directoryName);
423 textrestartFile = std::ofstream(fileName, std::ios::out);
424 textrestartFile <<
"Features" <<
"\n";
425 for(
size_t i = 0; i < Features.size(); i++)
426 textrestartFile << Features[i] <<
"\n";
427 textrestartFile.close();
430std::vector<float> Convert_Precision(std::vector<double>& Features)
432 std::vector<float> a;
float temp = 0.0;
433 for(
size_t i = 0; i < Features.size(); i++)
435 temp =
static_cast<float>(Features[i]); a.push_back(temp);
440double Predict_From_FeatureMatrix_Move(
Simulations& Sim,
Components& SystemComponents,
int DNN_CalcType)
442 if(!SystemComponents.UseDNNforHostGuest || !SystemComponents.UseLCLin)
return 0.0;
445 time = omp_get_wtime();
447 std::vector<double>Features(SystemComponents.Nfeatures);
448 std::vector<std::vector<double>> Distances = CalculatePairDistances_GPU(Sim, SystemComponents, Features, DNN_CalcType, 1);
449 std::vector<float> Float_Features = Convert_Precision(Features);
450 SystemComponents.DNNFeatureTime += omp_get_wtime() - time;
452 time = omp_get_wtime();
456 auto real_input = cppflow::tensor(Float_Features, {1, SystemComponents.Nfeatures});
457 auto output = SystemComponents.DNNModel[ModelID]({{SystemComponents.InputLayer[ModelID], real_input}},{
"StatefulPartitionedCall:0"});
458 double prediction =
static_cast<double>(output[0].get_data<
float>()[0]);
460 SystemComponents.DNNPredictTime += omp_get_wtime() - time;
462 return prediction * SystemComponents.DNNEnergyConversion;
473 if(!SystemComponents.UseDNNforHostGuest || !SystemComponents.UseLCLin)
return 0.0;
475 size_t NMol = SystemComponents.NumberOfMolecule_for_Component[ads_comp];
477 std::vector<double>AllFeatures;
478 for(
size_t iMol = 0; iMol < NMol; iMol++)
480 std::vector<double>Features(SystemComponents.Nfeatures);
481 std::vector<std::vector<double>> Distances = CalculatePairDistances_GPU(Sim, SystemComponents, Features, TOTAL, iMol);
482 AllFeatures.insert(std::end(AllFeatures), std::begin(Features), std::end(Features));
485 WriteDistances(SystemComponents, Distances);
486 WriteFeatures(SystemComponents, Features);
489 std::vector<float> AllFloatFeatures = Convert_Precision(AllFeatures);
493 auto real_input = cppflow::tensor(AllFloatFeatures, {NMol, SystemComponents.Nfeatures});
494 auto output = SystemComponents.DNNModel[ModelID]({{SystemComponents.InputLayer[ModelID], real_input}},{
"StatefulPartitionedCall:0"});
497 for(
size_t i = 0; i < NMol; i++)
499 double prediction =
static_cast<double>(output[0].get_data<
float>()[i]);
503 return tot * SystemComponents.DNNEnergyConversion;
Definition data_struct.h:746
Definition data_struct.h:819
Definition cppflow_LCLin.h:160
Definition data_struct.h:843
Definition data_struct.h:1027