My Project
Loading...
Searching...
No Matches
cppflow_LCLin.h
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)
2{
3 //Parallelized over framework atoms//
4 //Serialized over adsorbate atoms (only one adsorbate molecule)//
5 size_t ij = blockIdx.x * blockDim.x + threadIdx.x;
6
7 if(ij < totalThreads)
8 {
9 size_t comp = 0;
10
11 const Atoms Component=System[comp];
12 //const size_t typeA = Component.Type[ij];
13 //const size_t MoleculeID = Component.MolID[ij];
14
15 Atoms Adsorbate = Old;
16 if(!UsedForAMove) //Then it is for calculating the total energy//
17 Adsorbate = System[1];
18 for(size_t j = 0; j < Molsize; j++)
19 {
20 if(!ConsiderThisAdsorbateAtom[j]) continue;
21
22 size_t new_j = j + skip;
23 //size_t typeB = Adsorbate.Type[new_j];
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];
28
29 Distances[WriteTo] = rr_dot;
30 //if(ij == 0) printf("ij: %lu, j: %lu, WriteTo: %lu, TypeA(Framework): %lu, TypeB(Adsorbate): %lu, distance: %.5f\n", ij, j, WriteTo, typeA, typeB, rr_dot);
31 }
32 }
33}
34
35static inline size_t MatchInteractionIndex(std::vector<int3>& DNNInteractionList, size_t typeA, size_t typeB)
36{
37 size_t val=0;
38 for(size_t i = 0; i < DNNInteractionList.size(); i++)
39 {
40 if(typeA == DNNInteractionList[i].x && typeB == DNNInteractionList[i].y)
41 {
42 val = i;
43 break;
44 }
45 }
46 return val;
47}
48
49std::vector<std::vector<double>> CalculatePairDistances_CPU(Atoms* System, Boxsize Box, bool* ConsiderThisAdsorbateAtom, std::vector<int3>& DNNInteractionList)
50{
51 std::vector<std::vector<double>>CPU_Distances(DNNInteractionList.size());
52 size_t ads_comp = 1;
53 for(size_t i = 0; i < System[0].size; i++)
54 for(size_t j = 0; j < System[ads_comp].Molsize; j++)
55 {
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);
65 }
66
67 return CPU_Distances;
68}
69
70void PrepareOutliersFiles()
71{
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();
78
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);
84
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();
89
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();
94
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();
99}
100
101static inline cppflow::model load_model(std::string& NAME)
102{
103 setenv("CUDA_VISIBLE_DEVICES", "", 1); //DO NOT use GPU for the tf model//
104 cppflow::model DNNModel(NAME);
105 setenv("CUDA_VISIBLE_DEVICES", "1", 1); //After setting up tf model, set the GPU as visible again//
106 return DNNModel;
107}
108
109static inline void Prepare_Model(Components& SystemComponents, std::string& NAME)
110{
111 cppflow::model TempModel = load_model(NAME);
112 SystemComponents.DNNModel.push_back(TempModel);
113}
114
115void Read_LCLin_Model(Components& SystemComponents)
116{
117 // Try not to let tensorflow occupy the whole GPU//
118 // Serialized config options (example of 30% memory fraction)
119 // Read more to see how to obtain the serialized options
120
121 //Use the model on the CPU//
122
123 size_t ModelID = 0;
124
125 //SystemComponents.DNNModel.emplace_back(std::move(cppflow::model (SystemComponents.ModelName[ModelID])));
126
127 Prepare_Model(SystemComponents, SystemComponents.ModelName[ModelID]);
128
129 /*
130 auto second_output = SystemComponents.DNNModel[0]({{SystemComponents.InputLayer[ModelID], real_input}},{"StatefulPartitionedCall:0"});
131 //auto second_output = SystemComponents.DNNModel[0](real_input);
132 double asdasd = static_cast<double>(second_output[0].get_data<double>()[0]);
133 */
134
135 //Read Min Max//
136 SystemComponents.DNNMinMax = ReadMinMax();
137 printf("MinMax size: %zu\n", SystemComponents.DNNMinMax.size());
138
139 //Prepare for the outlier files//
140 PrepareOutliersFiles();
141
142}
143
144double DNN_Evaluate(Components& SystemComponents, std::vector<double>& Feature)
145{
146
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});
150
151 size_t ModelID = 0;
152
153 double asd = 0.0;
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);
157 return asd;
158}
159
161 bool operator()(const double3& a, const double3& b) const
162 {
163 return a.z < b.z;
164 }
165};
166
167static inline void NormalizeMinMax(std::vector<double>& r_features, std::vector<double2>& MinMax, size_t start)
168{
169
170 for(size_t i = 0; i < r_features.size(); i++)
171 {
172 r_features[i] = (r_features[i] - MinMax[start + i].x) / (MinMax[start + i].y - MinMax[start + i].x);
173 //printf("r_feature: %.5f\n", r_features[i]);
174 }
175}
176
177static inline std::vector<double> Calcuate_Feature_from_Distance(double dist_sq, size_t start, std::vector<double2>& MinMax)
178{
179 //Add Normalization in this function//
180 double dist = sqrt(dist_sq);
181 double one_over_dist = 1.0 / dist;
182 double one_over_distsq = 1.0 / dist_sq;
183
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)});
185 //printf("FeatureSize: %zu, r_featureSize: %zu, start: %zu\n", Features.size(), r_features.size(), start);
186 NormalizeMinMax(r_features, MinMax, start);
187 return r_features;
188}
189
190std::vector<std::vector<double>> CalculatePairDistances_GPU(Simulations& Sim, Components& SystemComponents, std::vector<double>& Features, int DNN_CalcType, size_t NMol)
191{
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);
194 size_t ads_comp = 1;
195
196 double time = omp_get_wtime();
197 switch(DNN_CalcType)
198 {
199 case TOTAL:
200 {
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);
202 break;
203 }
204 case OLD:
205 {
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);
207 break;
208 }
209 case NEW:
210 {
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);
212 break;
213 }
214 case REINSERTION_OLD:
215 {
216 size_t skip = 0;
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);
218 break;
219 }
220 case REINSERTION_NEW:
221 {
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);
224 break;
225 }
226 case DNN_INSERTION: case DNN_DELETION:
227 {
228 size_t skip = 0;
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);
230 break;
231 }
232 }
233 SystemComponents.DNNGPUTime += omp_get_wtime() - time;
234
235 time = omp_get_wtime();
236 //Create 2D vector of distances//
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();
241
242 std::vector<size_t> starts; size_t start_sum = 0;
243 for(size_t i = 0; i < SystemComponents.DNNInteractionList.size(); i++)
244 {
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());
248 }
249
250
251 //#pragma omp parallel for
252 for(size_t i = 0; i < SystemComponents.DNNInteractionList.size(); i++)
253 {
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);
256
257 //Sort the vector//
258 std::sort(std::execution::unseq, Distances[i].begin(), Distances[i].end());
259 }
260 SystemComponents.DNNstdsortTime += omp_get_wtime() - stdsort_time;
261
262 double FeaturizationTime = omp_get_wtime();
263
264 #pragma omp parallel for
265 for(size_t ij = 0; ij < SystemComponents.DNNInteractionList.size() * top; ij++)
266 {
267 size_t i = ij / top;
268 size_t j = ij % top;
269 SORTED_Dist[i][j] = Distances[i][j];
270 //Calculate, Normalize//
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];
275 //Features.insert(std::end(Features), std::begin(r_features), std::end(r_features));
276 }
277
278 SystemComponents.DNNFeaturizationTime += omp_get_wtime() - FeaturizationTime;
279 SystemComponents.DNNSortTime += omp_get_wtime() - time;
280
281 return Distances;
282}
283
284void Prepare_FeatureMatrix(Simulations& Sim, Components& SystemComponents, Atoms* HostSystem, Boxsize Host_Box)
285{
286 //Idea: Write a Kernel to get the distances, write them into pregenerated arrays, then sort using thrust//
287 //First Step, count number of pairs for each type of interactions//
288 size_t ads_comp = 1; //Zhao's note: Assuming the first component is the adsorbate we are interested in//
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]);
294
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]);
298
299 //Hard-coded types for Framework//
300 std::vector<size_t> FrameworkAtomTypesForDNN = {0, 1, 2, 3};
301
302
303 //Ow-Mg (4-0), Ow-O(4-1), Ow-C(4-2), Ow-H(4-3), Hw-Mg(5-0), Hw-O(5-1), Hw-C(5-2), Hw-H(5-3)//
304 //DNNInteractionList is a double3 vector, x: adsorbate-type, y: framework-type, z: Number of framework PseudoAtomsfor this type//
305 for(size_t i = 0; i < AdsorbateAtomTypesForDNN.size(); i++)
306 {
307 size_t ads_type = AdsorbateAtomTypesForDNN[i];
308 for(size_t j = 0; j < FrameworkAtomTypesForDNN.size(); j++)
309 {
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");
313
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]);
316 }
317 }
318 //This determines the number of features for the model//
319 //NTypes: Number of interaction types
320 //Ndistances: top N distances for each type of interaction//
321 //Ninteractions: for each distance, calculate 6 interaction forms//
322 size_t Ntypes = SystemComponents.DNNInteractionList.size(); size_t Ndistances = 9; size_t Ninteractions = 6;
323 SystemComponents.Nfeatures = Ntypes * Ndistances * Ninteractions;
324 //You need an array to determine where to store distances for each pair when calculating the distance on the GPU, this is an array that does the trick.//
325 //Zhao's note: this list will be created only once at the beginning of the simulation//
326 SystemComponents.IndexList.resize(SystemComponents.DNNInteractionList.size());
327
328 size_t count = 0;
329
330 for(size_t i = 0; i < HostSystem[ads_comp].Molsize; i++)
331 printf("ConsiderThisATom? %d\n", SystemComponents.ConsiderThisAdsorbateAtom[i]);
332
333 //Loop over framework atoms and the pseudo atoms in the adsorbate, check the Index of the Interaction//
334 for(size_t i = 0; i < HostSystem[0].size; i++)
335 for(size_t j = 0; j < HostSystem[ads_comp].Molsize; j++)
336 {
337 count++;
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);
342 if(i < 5)
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);
345 }
346 //Then, flatten the 2d vector, and reverse the values that you are recording//
347 //stored value --> index, index --> stored value//
348 std::vector<size_t>FlatIndexList;
349 for(size_t i = 0; i < SystemComponents.IndexList.size(); i++)
350 {
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]);
354 }
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;
359 //Do a test//
360 //Now the elements in the flatIndexList contains the "count" we recorded//
361 //the count will then be used to match the global pair-interaction ids//
362 //Test: Use the value in New as the adsorbate positions//
363 //size_t test_count=0;
364 for(size_t i = 0; i < HostSystem[0].size; i++)
365 for(size_t j = 0; j < HostSystem[ads_comp].Molsize; j++)
366 {
367 if(!SystemComponents.ConsiderThisAdsorbateAtom[j]) continue;
368 if(i < 5)
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]);
370 }
371
372
373 //Now the elements in the flatIndexList contains the "count" we recorded//
374 //the count will then be used to match the global pair-interaction ids//
375 //Test: Use the value in New as the adsorbate positions//
376 cudaMalloc(&SystemComponents.device_InverseIndexList, Listsize * sizeof(size_t));
377 cudaMemcpy(SystemComponents.device_InverseIndexList, InverseIndexList.data(), Listsize * sizeof(size_t), cudaMemcpyHostToDevice);
378
379 //Use Pinned Memory for slightly better mem transfer//
380 printf("Listsize: %zu\n", Listsize);
381 //cudaMallocHost(&SystemComponents.device_Distances, Listsize * sizeof(double));
382 //cudaMalloc(&SystemComponents.device_Distances, Listsize * sizeof(double));
383 printf("Size of the device Distance list: %zu\n", Listsize);
384}
385
386static inline void WriteDistances(Components& SystemComponents, std::vector<std::vector<double>> Distances)
387{
388 //Write to a file for checking//
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();
393
394 std::filesystem::path directoryName = cwd /dirname;
395 std::filesystem::path fileName = cwd /fname;
396 std::filesystem::create_directories(directoryName);
397
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++)
401 {
402 //std::sort(CPU_Distances[i].begin(), CPU_Distances[i].end());
403 for(size_t j = 0; j < Distances[i].size(); j++)
404 {
405 //printf("TypeA(Framework): %.5f TypeB(Adsorbate): %.5f distance: %.5f\n", Distances[i][j].x, Distances[i][j].y, Distances[i][j].z);
406 textrestartFile << Distances[i][j] << "\n";
407 }
408 }
409 textrestartFile.close();
410}
411
412static inline void WriteFeatures(Components& SystemComponents, std::vector<double> Features)
413{
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();
418
419 std::filesystem::path directoryName = cwd /dirname;
420 std::filesystem::path fileName = cwd /fname;
421 std::filesystem::create_directories(directoryName);
422
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();
428}
429
430std::vector<float> Convert_Precision(std::vector<double>& Features)
431{
432 std::vector<float> a; float temp = 0.0;
433 for(size_t i = 0; i < Features.size(); i++)
434 {
435 temp = static_cast<float>(Features[i]); a.push_back(temp);
436 }
437 return a;
438}
439
440double Predict_From_FeatureMatrix_Move(Simulations& Sim, Components& SystemComponents, int DNN_CalcType)
441{
442 if(!SystemComponents.UseDNNforHostGuest || !SystemComponents.UseLCLin) return 0.0;
443 double time;
444 //Get GPU distances//
445 time = omp_get_wtime();
446
447 std::vector<double>Features(SystemComponents.Nfeatures); //feature elements//
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;
451
452 time = omp_get_wtime();
453 size_t ModelID = 0;
454
455 //Use the DNN Model to predict//
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]);
459
460 SystemComponents.DNNPredictTime += omp_get_wtime() - time;
461 //printf("DNN Prediction: %.5f [kJ/mol], %.5f [internal unit]\n", prediction, prediction * SystemComponents.DNNEnergyConversion);
462 return prediction * SystemComponents.DNNEnergyConversion;
463 //Check CPU distances//
464 //std::vector<std::vector<double>> CPU_Distances = CalculatePairDistances_CPU(HostSystem, Host_Box, SystemComponents.ConsiderThisAdsorbateAtom, SystemComponents.DNNInteractionList);
465
466
467 //Write to a file for checking//
468 //WriteFeatureMatrix(SystemComponents, Distances);
469}
470
471double Predict_From_FeatureMatrix_Total(Simulations& Sim, Components& SystemComponents)
472{
473 if(!SystemComponents.UseDNNforHostGuest || !SystemComponents.UseLCLin) return 0.0;
474 size_t ads_comp = 1;
475 size_t NMol = SystemComponents.NumberOfMolecule_for_Component[ads_comp];
476 //Get GPU distances//
477 std::vector<double>AllFeatures;
478 for(size_t iMol = 0; iMol < NMol; iMol++)
479 {
480 std::vector<double>Features(SystemComponents.Nfeatures); //feature elements//
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));
483 if(iMol == 0)
484 {
485 WriteDistances(SystemComponents, Distances);
486 WriteFeatures(SystemComponents, Features);
487 }
488 }
489 std::vector<float> AllFloatFeatures = Convert_Precision(AllFeatures);
490 size_t ModelID = 0;
491
492 //Use the DNN Model to predict//
493 auto real_input = cppflow::tensor(AllFloatFeatures, {NMol, SystemComponents.Nfeatures});
494 auto output = SystemComponents.DNNModel[ModelID]({{SystemComponents.InputLayer[ModelID], real_input}},{"StatefulPartitionedCall:0"});
495
496 double tot = 0.0;
497 for(size_t i = 0; i < NMol; i++)
498 {
499 double prediction = static_cast<double>(output[0].get_data<float>()[i]);
500 tot += prediction;
501 //printf("Molecule %zu, DNN Prediction: %.5f [kJ/mol], %.5f [internal unit]\n", i, prediction, prediction * SystemComponents.DNNEnergyConversion);
502 }
503 return tot * SystemComponents.DNNEnergyConversion;
504}
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