42 #include "CSVReader/CSVReader.hpp"
49 #include "include/RemoveLines.hpp"
51 #include "FiniteDifference/FD_simple.hpp"
53 #include "include/DiffusionSpace_sparseGrid.hpp"
54 #include "include/HelpFunctions_diffusion.hpp"
56 #include "Decomposition/Distribution/BoxDistribution.hpp"
58 const size_t dims = 3;
61 const size_t PHI_FULL = 0;
69 const size_t PHI_PHASE = 0;
71 const size_t U_NPLUS1 = 2;
72 const size_t K_SINK = 3;
83 constexpr
size_t x = 0, y = 1, z = 2;
86 const std::string path_to_redistancing_result =
87 "output_sussman_sparse_grid_porousCeramics_1216x1016x941/";
89 const std::string redistancing_filename =
"sparseGrid_initial.hdf5";
91 const std::string path_to_size = path_to_redistancing_result;
93 const std::string output_name =
"output_heat_conduction_porous_ceramics_1216x1016x941_sparse_grid";
96 int main(
int argc,
char* argv[])
99 openfpm_init(&argc, &argv);
100 auto & v_cl = create_vcluster();
103 timer t_iterative_diffusion_total;
104 timer t_iteration_wct;
113 const std::string path_output = cwd +
"/" + output_name +
"/";
116 if(v_cl.rank()==0) std::cout <<
"Redistancing result will be loaded from " << path_to_redistancing_result << redistancing_filename << std::endl;
118 if(v_cl.rank()==0) std::cout <<
"Outputs will be saved to " << path_output << std::endl;
130 read_csv_to_vector(path_to_size +
"/N.csv", v_sz, m, n);
131 read_csv_to_vector(path_to_size +
"/L.csv", v_L, m, n);
133 const float Lx_low = 0.0;
134 const float Lx_up = v_L.get(x);
135 const float Ly_low = 0.0;
136 const float Ly_up = v_L.get(y);
137 const float Lz_low = 0.0;
138 const float Lz_up = v_L.get(z);
141 size_t sz[dims] = {v_sz.get(x), v_sz.get(y), v_sz.get(z)};
148 grid_in_type g_dist(sz, box, ghost);
150 std::cout <<
"Rank " << v_cl.rank() <<
" starts loading redistancing result..." << std::endl;
151 g_dist.load(path_to_redistancing_result +
"/" + redistancing_filename);
152 std::cout <<
"Rank " << v_cl.rank() <<
" finished loading redistancing result." << std::endl;
159 std::cout <<
"Rank " << v_cl.rank() <<
" starts creating sparse grid." << std::endl;
162 typedef sgrid_dist_id_gpu<dims, float, props_sparse, CudaMemory, Dec> sparse_grid_type;
163 sparse_grid_type g_sparse(sz, box, ghost);
164 g_sparse.setPropNames(prop_names_sparse);
167 const float d_low = 0.0, d_up = Lx_up;
168 get_diffusion_domain_sparse_grid<PHI_FULL, PHI_PHASE>(g_dist, g_sparse, d_low, d_up);
173 std::cout <<
"Rank " << v_cl.rank() <<
" finished creating sparse grid." << std::endl;
178 const float u0 = 1.0;
180 const float sink_rate = 0.1;
181 std::cout <<
"Diffusion coefficient in s/(mm^2): " <<
D << std::endl;
182 std::cout <<
"Sink rate in 1/s: " << sink_rate << std::endl;
185 const float center[dims] = {(Lx_up+Lx_low)/(
float)2,
186 (Ly_up+Ly_low)/(
float)2,
187 (Lz_up+Lz_low)/(
float)2};
188 const float R = Lz_up / 4.0;
190 auto dom = g_sparse.getDomainIterator();
193 auto key = dom.get();
198 if(((coords.
get(x) - center[x]) * (coords.
get(x) - center[x])
199 + (coords.
get(y) - center[y]) * (coords.
get(y) - center[y])
200 + (coords.
get(z) - center[z]) * (coords.
get(z) - center[z])) <= R * R)
202 g_sparse.template insertFlush<U_N>(key) = u0;
204 else g_sparse.template insertFlush<U_N>(key) = 0.0;
217 typedef typename GetCpBlockType<decltype(g_sparse),0,1>::type CpBlockType;
219 const float dx = g_sparse.spacing(x), dy = g_sparse.spacing(y), dz = g_sparse.spacing(z);
221 const float dt = diffusion_time_step(g_sparse,
D);
222 std::cout <<
"dx, dy, dz = " << g_sparse.spacing(x) <<
"," << g_sparse.spacing(y) <<
"," << g_sparse.spacing(z) << std::endl;
223 std::cout <<
"dt = 1/(4D) * 1/(1/dx^2 + 1/dy^2 + 1/dz^2) = " << dt << std::endl;
225 auto func_heatDiffusion_withSink = [dx, dy, dz, dt, d_low,
D, sink_rate] __device__ (
240 float u_c = u(i, j, k);
241 float u_px = u(i+1, j, k);
242 float u_mx = u(i-1, j, k);
243 float u_py = u(i, j+1, k);
244 float u_my = u(i, j-1, k);
245 float u_pz = u(i, j, k+1);
246 float u_mz = u(i, j, k-1);
249 float phi_c = phi(i, j, k);
250 float phi_px = phi(i+1, j, k);
251 float phi_mx = phi(i-1, j, k);
252 float phi_py = phi(i, j+1, k);
253 float phi_my = phi(i, j-1, k);
254 float phi_pz = phi(i, j, k+1);
255 float phi_mz = phi(i, j, k-1);
261 if (phi_px <= d_low + std::numeric_limits<float>::epsilon()) {u_px = u_c; k_sink = sink_rate;}
262 if (phi_mx <= d_low + std::numeric_limits<float>::epsilon()) {u_mx = u_c; k_sink = sink_rate;}
263 if (phi_py <= d_low + std::numeric_limits<float>::epsilon()) {u_py = u_c; k_sink = sink_rate;}
264 if (phi_my <= d_low + std::numeric_limits<float>::epsilon()) {u_my = u_c; k_sink = sink_rate;}
265 if (phi_pz <= d_low + std::numeric_limits<float>::epsilon()) {u_pz = u_c; k_sink = sink_rate;}
266 if (phi_mz <= d_low + std::numeric_limits<float>::epsilon()) {u_mz = u_c; k_sink = sink_rate;}
268 block.template get<K_SINK>()[offset] = k_sink;
271 u_out = u_c + dt * (
D * ((u_px + u_mx - 2 * u_c)/(dx * dx)
272 + (u_py + u_my - 2 * u_c)/(dy * dy)
273 + (u_pz + u_mz - 2 * u_c)/(dz * dz))
286 std::string path_time_filewct = path_output +
"/time_wct_incl_ghostCommunication" + std::to_string(v_cl.rank()) +
".txt";
287 std::ofstream out_wct_file(path_time_filewct, std::ios_base::app);
290 std::string path_time_filegpu = path_output +
"/time_GPU" + std::to_string(v_cl.rank()) +
".txt";
291 std::ofstream out_GPUtime_file(path_time_filegpu, std::ios_base::app);
296 const size_t iterations = 103;
297 const size_t number_write_outs = 10;
298 const size_t interval_write = iterations / number_write_outs;
303 g_sparse.template hostToDevice<U_N, U_NPLUS1, K_SINK, PHI_PHASE>();
304 g_sparse.template ghost_get<PHI_PHASE, K_SINK>(RUN_ON_DEVICE | SKIP_LABELLING);
305 t_iterative_diffusion_total.
start();
306 while(iter <= iterations)
310 t_iteration_wct.
start();
311 g_sparse.template ghost_get<U_N>(RUN_ON_DEVICE | SKIP_LABELLING);
313 g_sparse.template conv2_b<
319 {(
long int) sz[x]-1, (
long int) sz[y]-1, (
long int) sz[z]-1},
320 func_heatDiffusion_withSink);
322 t_iteration_wct.
stop();
325 out_wct_file << t_iteration_wct.
getwct() <<
",";
326 out_GPUtime_file << t_GPU.getwctGPU() <<
",";
328 t_iteration_wct.
reset();
332 t_iteration_wct.
start();
333 g_sparse.template ghost_get<U_NPLUS1>(RUN_ON_DEVICE | SKIP_LABELLING);
335 g_sparse.template conv2_b<
341 {(
long int) sz[x]-1, (
long int) sz[y]-1, (
long int) sz[z]-1},
342 func_heatDiffusion_withSink);
344 t_iteration_wct.
stop();
347 out_wct_file << t_iteration_wct.
getwct() <<
",";
348 out_GPUtime_file << t_GPU.getwctGPU() <<
",";
350 t_iteration_wct.
reset();
354 if (iter % interval_write == 0)
357 g_sparse.template deviceToHost<U_N, U_NPLUS1, K_SINK, PHI_PHASE>();
360 g_sparse.write_frame(path_output +
"g_sparse_diffuse", iter, FORMAT_BINARY);
364 monitor_total_concentration<U_N>(g_sparse, t, iter, path_output,
"total_conc.csv");
368 monitor_total_concentration<U_NPLUS1>(g_sparse, t, iter, path_output,
"total_conc.csv");
375 t_iterative_diffusion_total.
stop();
376 std::cout <<
"Rank " << v_cl.rank() <<
" total time for iterative diffusion incl. write outs: " << t_iterative_diffusion_total.
getwct() << std::endl;
377 std::cout <<
"Rank " << v_cl.rank() <<
" finished diffusion." << std::endl;
384 std::cout <<
"Rank " << v_cl.rank() <<
" total time for the whole programm : " << t_total.
getwct() << std::endl;
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.
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.
__device__ __host__ const T & get(unsigned int i) const
Get coordinate.
This is a distributed grid.
Class for cpu time benchmarking.
void stop()
Stop the timer.
void reset()
Reset the timer.
void start()
Start the timer.
double getwct()
Return the elapsed real time.
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...