OpenFPM_pdata  4.1.0
Project that contain the implementation of distributed structures
 
Loading...
Searching...
No Matches
main.cu
Go to the documentation of this file.
1//
2// Created by jstark on 2023-04-24.
3//
29
30
31
32#include <iostream>
33#include <typeinfo>
34
35#include "CSVReader/CSVReader.hpp"
36
37// Include redistancing files
41#include "RemoveLines.hpp" // For removing thin (diagonal or straight) lines
42
44#include "level_set/redistancing_Sussman/AnalyticalSDF.hpp"
45#include "FiniteDifference/FD_simple.hpp"
46#include "Decomposition/Distribution/BoxDistribution.hpp"
47#include "timer.hpp"
48
49// Help functions for simulating diffusion in the image-based geometry
50#include "include/DiffusionSpace_sparseGrid.hpp"
51#include "include/HelpFunctions_diffusion.hpp"
52
53
54
57
58// input
59const std::string path_to_redistancing_result =
60"/INPUT_PATH/benchmarks/CaCO3/sussman_redistancing/build/output_sussman_maxIter6e3_CaCO3_fluidPhase_531x531x531/";
61
62
63const std::string redistancing_filename = "grid_CaCO3_post_redistancing.hdf5";
64const std::string path_to_size = path_to_redistancing_result;
65
66// output
67const std::string output_name = "output_inhomogDiffusion_CaCO3_fluid_phase";
70// Property indices full grid
71const size_t PHI_FULL = 0;
72
74const openfpm::vector<std::string> prop_names_full = {"Phi_Sussman_Out"};
75
76
77// Property indices sparse grid
78const size_t PHI_PHASE = 0;
79const size_t CONC_N = 1;
80const size_t CONC_NPLUS1 = 2;
81const size_t DIFFUSION_COEFFICIENT = 3;
82
83
85const openfpm::vector<std::string> prop_names_sparse = {"PHI_PHASE",
86 "CONC_N",
87 "CONC_NPLUS1",
88 "DIFFUSION_COEFFICIENT"};
89
90// Space indices
91constexpr size_t x = 0, y = 1, z = 2;
92
93// Parameters for grid size
94const size_t dims = 3;
98
112int main(int argc, char* argv[])
113{
114 // Initialize library.
115 openfpm_init(&argc, &argv);
116 auto & v_cl = create_vcluster();
117
118 timer t_total;
119 timer t_iterative_diffusion_total;
120 timer t_iteration_wct;
121 timer t_GPU;
122
123 t_total.start();
124
127 // Set current working directory, define output paths and create folders where output will be saved
128 std::string cwd = get_cwd();
129 const std::string path_output = cwd + "/" + output_name + "/";
131
132 if(v_cl.rank()==0) std::cout << "Redistancing result will be loaded from " << path_to_redistancing_result << redistancing_filename << std::endl;
133
134 if(v_cl.rank()==0) std::cout << "Outputs will be saved to " << path_output << std::endl;
135
137
139
149 // Create full grid and load redistancing result
151 // Read size of sussman redistancing output grid from csv files
154
155 size_t m, n;
156 read_csv_to_vector(path_to_size + "/N.csv", v_sz, m, n);
157 read_csv_to_vector(path_to_size + "/L.csv", v_L, m, n);
158
159 const float Lx_low = 0.0;
160 const float Lx_up = v_L.get(x);
161 const float Ly_low = 0.0;
162 const float Ly_up = v_L.get(y);
163 const float Lz_low = 0.0;
164 const float Lz_up = v_L.get(z);
165
167 size_t sz[dims] = {v_sz.get(x), v_sz.get(y), v_sz.get(z)};
168 Box<dims, float> box({Lx_low, Ly_low, Lz_low}, {Lx_up, Ly_up, Lz_up});
169 Ghost<dims, long int> ghost(1);
170
171 // Defining the decomposition and the input grid type
174 grid_in_type g_dist(sz, box, ghost);
175
176 g_dist.load(path_to_redistancing_result + "/" + redistancing_filename); // Load SDF
178
180
193 // Get sparse grid of diffusion domain
195 // Create sparse grid
196 std::cout << "Rank " << v_cl.rank() << " starts creating sparse grid." << std::endl;
197 typedef sgrid_dist_id_gpu<dims, float, props_sparse, CudaMemory, Dec> sparse_grid_type;
198 sparse_grid_type g_sparse(sz, box, ghost);
199 g_sparse.setPropNames(prop_names_sparse);
200
201 // Lower and upper bound for SDF indicating fluid phase
202 const float d_low = 0.0, d_up = Lx_up;
203
204 // Alternatively: load sparse grid of diffusion domain from HDF5 file
205// g_sparse.load(path_to_redistancing_result + "/" + redistancing_filename);
206
207 // Obtain sparse grid
208 get_diffusion_domain_sparse_grid<PHI_FULL, PHI_PHASE>(g_dist, g_sparse, d_low, d_up);
209 std::cout << "Rank " << v_cl.rank() << " finished creating sparse grid." << std::endl;
210
212
214
227 // Initialize properties on sparse grid
229 // Diffusion coefficient
230 const float min_porosity = 0.0;
231 const float max_porosity = 1.0;
232
233 float D_molecular = 1.0;
234 float D_min = D_molecular * min_porosity;
235 float D_max = D_molecular * max_porosity;
236
237 float min_sdf_calcite = get_min_val<PHI_PHASE>(g_sparse);
238
239 float x_stretch = 4.0 * g_sparse.spacing(x);
240 float x_shift = min_sdf_calcite * x_stretch;
241 float y_min = D_min;
242 float y_max = D_max - D_min;
243
244 const float c0 = 10.0;
245
246 auto dom = g_sparse.getDomainIterator();
247 while(dom.isNext())
248 {
249 auto key = dom.get();
250
251 Point<dims, float> coords = g_sparse.getPos(key);
252 if(coords.get(x) <= 0.05 * Lx_up)
253 {
254 g_sparse.template insertFlush<CONC_N>(key) = c0;
255 }
256 else g_sparse.template insertFlush<CONC_N>(key) = 0.0;
257
258 // Compute diffusion coefficient and store on grid
259 g_sparse.template insertFlush<DIFFUSION_COEFFICIENT>(key) =
260 get_smooth_sigmoidal(
261 g_sparse.template getProp<PHI_PHASE>(key) - min_sdf_calcite,
262 x_shift,
263 x_stretch,
264 y_min,
265 y_max);
266 ++dom;
267 }
268
271
285 // Defining the functor that will be passed to the GPU to simulate reaction-diffusion
286 // GetCpBlockType<typename SGridGpu, unsigned int prp, unsigned int stencil_size>
287 typedef typename GetCpBlockType<decltype(g_sparse),0,1>::type CpBlockType;
288
289 const float dx = g_sparse.spacing(x), dy = g_sparse.spacing(y), dz = g_sparse.spacing(z);
290 // Stability criterion for FTCS
291 const float dt = diffusion_time_step(g_sparse, D_max);
292
293 std::cout << "dx = " << dx << "dy = " << dy << "dz = " << dz << std::endl;
294 std::cout << "dt = " << dt << std::endl;
295
296
297 auto func_inhomogDiffusion = [dx, dy, dz, dt, d_low] __device__ (
298 float & u_out, // concentration output
299 float & D_out, // diffusion coefficient output (dummy)
300 float & phi_out, // sdf of domain output (dummy)
301 CpBlockType & u, // concentration input
302 CpBlockType & D, // diffusion coefficient input (to get also the half step diffusion coefficients)
303 CpBlockType & phi, // sdf of domain (for boundary conditions)
304 auto & block,
305 int offset,
306 int i,
307 int j,
308 int k
309 )
310 {
312 // Stencil
313 // Concentration
314 float u_c = u(i, j, k);
315 float u_px = u(i+1, j, k);
316 float u_mx = u(i-1, j, k);
317 float u_py = u(i, j+1, k);
318 float u_my = u(i, j-1, k);
319 float u_pz = u(i, j, k+1);
320 float u_mz = u(i, j, k-1);
321
322 // Signed distance function
323 float phi_c = phi(i, j, k);
324 float phi_px = phi(i+1, j, k);
325 float phi_mx = phi(i-1, j, k);
326 float phi_py = phi(i, j+1, k);
327 float phi_my = phi(i, j-1, k);
328 float phi_pz = phi(i, j, k+1);
329 float phi_mz = phi(i, j, k-1);
330
331 // Diffusion coefficient
332 float D_c = D(i, j, k);
333 float D_px = D(i+1, j, k);
334 float D_mx = D(i-1, j, k);
335 float D_py = D(i, j+1, k);
336 float D_my = D(i, j-1, k);
337 float D_pz = D(i, j, k+1);
338 float D_mz = D(i, j, k-1);
339
340 // Impose no-flux boundary conditions
341 if (phi_px <= d_low + std::numeric_limits<float>::epsilon()) {u_px = u_c; D_px = D_c;}
342 if (phi_mx <= d_low + std::numeric_limits<float>::epsilon()) {u_mx = u_c; D_mx = D_c;}
343 if (phi_py <= d_low + std::numeric_limits<float>::epsilon()) {u_py = u_c; D_py = D_c;}
344 if (phi_my <= d_low + std::numeric_limits<float>::epsilon()) {u_my = u_c; D_my = D_c;}
345 if (phi_pz <= d_low + std::numeric_limits<float>::epsilon()) {u_pz = u_c; D_pz = D_c;}
346 if (phi_mz <= d_low + std::numeric_limits<float>::epsilon()) {u_mz = u_c; D_mz = D_c;}
347
348 // Interpolate diffusion constant of half-step neighboring points
349 float D_half_px = (D_c + D_px)/2.0;
350 float D_half_mx = (D_c + D_mx)/2.0;
351 float D_half_py = (D_c + D_py)/2.0;
352 float D_half_my = (D_c + D_my)/2.0;
353 float D_half_pz = (D_c + D_pz)/2.0;
354 float D_half_mz = (D_c + D_mz)/2.0;
355
356 // Compute concentration of next time point
357 u_out = u_c + dt *
358 (1/(dx*dx) * (D_half_px * (u_px - u_c) - D_half_mx * (u_c - u_mx)) +
359 1/(dy*dy) * (D_half_py * (u_py - u_c) - D_half_my * (u_c - u_my)) +
360 1/(dz*dz) * (D_half_pz * (u_pz - u_c) - D_half_mz * (u_c - u_mz)));
361
362 // dummy outs
363 D_out = D_c;
364 phi_out = phi_c;
365 };
366
367
370
380 // Create text file to which wall clock time for iteration incl. ghost communication will be saved
381 std::string path_time_filewct = path_output + "/time_wct_incl_ghostCommunication" + std::to_string(v_cl.rank()) + ".txt";
382 std::ofstream out_wct_file(path_time_filewct, std::ios_base::app);
383
384 // Create text file to which GPU time will be saved
385 std::string path_time_filegpu = path_output + "/time_GPU" + std::to_string(v_cl.rank()) + ".txt";
386 std::ofstream out_GPUtime_file(path_time_filegpu, std::ios_base::app);
387
388 // Iterative diffusion
389 const size_t iterations = 103;
390 const size_t number_write_outs = 10;
391 const size_t interval_write = iterations / number_write_outs; // Find interval for writing to reach
392 size_t iter = 0;
393 float t = 0;
394
395 // Copy from host to GPU for simulation
396 g_sparse.template hostToDevice<CONC_N, CONC_NPLUS1, DIFFUSION_COEFFICIENT, PHI_PHASE>();
397 g_sparse.template ghost_get<DIFFUSION_COEFFICIENT, PHI_PHASE>(RUN_ON_DEVICE | SKIP_LABELLING);
398 t_iterative_diffusion_total.start();
399 while(iter <= iterations)
400 {
401 if (iter % 2 == 0)
402 {
403 t_iteration_wct.start();
404 g_sparse.template ghost_get<CONC_N>(RUN_ON_DEVICE | SKIP_LABELLING);
405 t_GPU.startGPU();
406 g_sparse.template conv3_b<
407 CONC_N,
408 DIFFUSION_COEFFICIENT,
409 PHI_PHASE,
410 CONC_NPLUS1,
411 DIFFUSION_COEFFICIENT,
412 PHI_PHASE, 1>
413 ({0, 0, 0},
414 {(long int) sz[x]-1, (long int) sz[y]-1, (long int) sz[z]-1},
415 func_inhomogDiffusion);
416 t_GPU.stopGPU();
417 t_iteration_wct.stop();
418
419 // Write out time to text-file
420 out_wct_file << t_iteration_wct.getwct() << ",";
421 out_GPUtime_file << t_GPU.getwctGPU() << ",";
422 }
423 else
424 {
425 t_iteration_wct.start();
426 g_sparse.template ghost_get<CONC_NPLUS1>(RUN_ON_DEVICE | SKIP_LABELLING);
427 t_GPU.startGPU();
428 g_sparse.template conv3_b<
429 CONC_NPLUS1,
430 DIFFUSION_COEFFICIENT,
431 PHI_PHASE,
432 CONC_N,
433 DIFFUSION_COEFFICIENT,
434 PHI_PHASE, 1>
435 ({0, 0, 0},
436 {(long int) sz[x]-1, (long int) sz[y]-1, (long int) sz[z]-1},
437 func_inhomogDiffusion);
438 t_GPU.stopGPU();
439 t_iteration_wct.stop();
440 // Write out time to text-file
441 out_wct_file << t_iteration_wct.getwct() << ",";
442 out_GPUtime_file << t_GPU.getwctGPU() << ",";
443 }
444
445
446 if (iter % interval_write == 0)
447 {
448 // Copy from GPU to host for writing
449 g_sparse.template deviceToHost<CONC_N, CONC_NPLUS1, DIFFUSION_COEFFICIENT, PHI_PHASE>();
450
451 // Write g_sparse to vtk
452 g_sparse.write_frame(path_output + "g_sparse_diffuse", iter, FORMAT_BINARY);
453
454 // Write g_sparse to hdf5
455 g_sparse.save(path_output + "g_sparse_diffuse_" + std::to_string(iter) + ".hdf5");
456
457 if (iter % 2 == 0)
458 {
459 monitor_total_concentration<CONC_N>(g_sparse, t, iter, path_output,"total_conc.csv");
460 }
461 else
462 {
463 monitor_total_concentration<CONC_NPLUS1>(g_sparse, t, iter, path_output,"total_conc.csv");
464 }
465 }
466
467
468 t += dt;
469 iter += 1;
470 }
471 t_iterative_diffusion_total.stop();
472 out_GPUtime_file.close();
473
475 t_total.stop();
476 std::cout << "Rank " << v_cl.rank() << " total time for the whole programm : " << t_total.getwct() << std::endl;
477
479
490 openfpm_finalize(); // Finalize openFPM library
491 return 0;
492}
494
Header file containing help-functions that perform on OpenFPM-grids.
Header file containing functions for loading pixel values from 2D image or 3D image stack (volume) st...
Header file containing functions for creating files and folders.
static void create_directory_if_not_exist(std::string path, bool silent=0)
Creates a directory if not already existent.
static std::string get_cwd()
Gets the current working directory and returns path as string.
Class for reinitializing a level-set function into a signed distance function using Sussman redistanc...
This class represent an N-dimensional box.
Definition Box.hpp:61
This class decompose a space into sub-sub-domains and distribute them across processors.
Derivative second order on h (spacing)
This class implement the point shape in an N-dimensional space.
Definition Point.hpp:28
__device__ __host__ const T & get(unsigned int i) const
Get coordinate.
Definition Point.hpp:172
This is a distributed grid.
Implementation of 1-D std::vector like structure.
Class for cpu time benchmarking.
Definition timer.hpp:28
void stop()
Stop the timer.
Definition timer.hpp:119
void start()
Start the timer.
Definition timer.hpp:90
double getwct()
Return the elapsed real time.
Definition timer.hpp:130
KeyT const ValueT ValueT OffsetIteratorT OffsetIteratorT int
[in] The number of segments that comprise the sorting data
get the type of the block
aggregate of properties, from a list of object if create a struct that follow the OPENFPM native stru...