1#include "Grid/grid_dist_id.hpp"
2#include "data_type/aggregate.hpp"
54extern "C" void update_new(
const int* lo,
const int* hi,
55 double* u,
const int* ulo,
const int* uhi,
56 double* v,
const int* vlo,
const int* vhi,
57 double* flu,
const int* fulo,
const int* fuhi,
58 double* flv,
const int* fvlo,
const int* fvhi,
59 const double * dt,
const double * uFactor,
const double * vFactor,
const double *
F,
71 auto it = OldU.getDomainIterator();
79 OldU.get<0>(key) = 1.0;
80 OldV.get<0>(key) = 0.0;
83 NewU.get<0>(key) = 0.0;
84 NewV.get<0>(key) = 0.0;
89 long int x_start = OldU.size(0)*1.55f/domain.getHigh(0);
90 long int y_start = OldU.size(1)*1.55f/domain.getHigh(1);
91 long int z_start = OldU.size(1)*1.55f/domain.getHigh(2);
93 long int x_stop = OldU.size(0)*1.85f/domain.getHigh(0);
94 long int y_stop = OldU.size(1)*1.85f/domain.getHigh(1);
95 long int z_stop = OldU.size(1)*1.85f/domain.getHigh(2);
99 auto it_init = OldU.getSubDomainIterator(start,stop);
101 while (it_init.isNext())
103 auto key = it_init.
get();
105 OldU.get<0>(key) = 0.5 + (((
double)std::rand())/RAND_MAX -0.5)/10.0;
106 OldV.get<0>(key) = 0.25 + (((
double)std::rand())/RAND_MAX -0.5)/20.0;
120 double uFactor_s,
double vFactor_s,
double deltaT,
double F,
double K)
123#ifndef FORTRAN_UPDATE
127 Vc::double_v uFactor = uFactor_s;
128 Vc::double_v vFactor = vFactor_s;
130 WHILE_M(OldU,star_stencil_3D)
131 auto & U_old = GET_GRID_M(OldU);
132 auto & V_old = GET_GRID_M(OldV);
134 auto & U_new = GET_GRID_M(NewU);
135 auto & V_new = GET_GRID_M(NewV);
136 ITERATE_3D_M(Vc::double_v::Size)
139 auto Cp = it.getStencil<0>();
142 auto mx = it.getStencil<1>();
143 auto px = it.getStencil<2>();
144 auto my = it.getStencil<3>();
145 auto py = it.getStencil<4>();
146 auto mz = it.getStencil<5>();
147 auto pz = it.getStencil<6>();
150 Vc::double_v u_c(&U_old.get<0>(Cp),Vc::Unaligned);
151 Vc::double_v u_mz(&U_old.get<0>(mz),Vc::Unaligned);
152 Vc::double_v u_pz(&U_old.get<0>(pz),Vc::Unaligned);
153 Vc::double_v u_my(&U_old.get<0>(my),Vc::Unaligned);
154 Vc::double_v u_py(&U_old.get<0>(py),Vc::Unaligned);
155 Vc::double_v u_mx(&U_old.get<0>(mx),Vc::Unaligned);
156 Vc::double_v u_px(&U_old.get<0>(px),Vc::Unaligned);
159 Vc::double_v v_c(&V_old.get<0>(Cp),Vc::Unaligned);
160 Vc::double_v v_mz(&V_old.get<0>(mz),Vc::Unaligned);
161 Vc::double_v v_pz(&V_old.get<0>(pz),Vc::Unaligned);
162 Vc::double_v v_my(&V_old.get<0>(my),Vc::Unaligned);
163 Vc::double_v v_py(&V_old.get<0>(py),Vc::Unaligned);
164 Vc::double_v v_mx(&V_old.get<0>(mx),Vc::Unaligned);
165 Vc::double_v v_px(&V_old.get<0>(px),Vc::Unaligned);
167 Vc::double_v out1 = u_c + uFactor * (u_mz + u_pz +
171 - deltaT * u_c * v_c * v_c
172 - deltaT *
F * (u_c - 1.0);
174 Vc::double_v out2 = v_c + vFactor * (v_mz + v_pz +
178 deltaT * u_c * v_c * v_c +
179 - deltaT * (
F+K) * v_c;
181 out1.store(&U_new.get<0>(Cp),Vc::Unaligned);
182 out2.store(&V_new.get<0>(Cp),Vc::Unaligned);
183 END_LOOP_M(Vc::double_v::Size)
191 auto & ginfo = OldU.getLocalGridsInfo();
193 for (
size_t i = 0 ; i < OldU.getN_loc_grid() ; i++)
195 auto & U_old = OldU.get_loc_grid(i);
196 auto & V_old = OldV.get_loc_grid(i);
198 auto & U_new = NewU.get_loc_grid(i);
199 auto & V_new = NewV.get_loc_grid(i);
201 int lo[3] = {(
int)ginfo.get(i).Dbox.getLow(0),(
int)ginfo.get(i).Dbox.getLow(1),(
int)ginfo.get(i).Dbox.getLow(2)};
202 int hi[3] = {(
int)ginfo.get(i).Dbox.getHigh(0),(
int)ginfo.get(i).Dbox.getHigh(1),(
int)ginfo.get(i).Dbox.getHigh(2)};
204 int ulo[3] = {0,0,0};
205 int uhi[3] = {(
int)ginfo.get(i).GDbox.getHigh(0),(
int)ginfo.get(i).GDbox.getHigh(1),(
int)ginfo.get(i).GDbox.getHigh(2)};
206 int nulo[3] = {0,0,0};
207 int nuhi[3] = {(
int)ginfo.get(i).GDbox.getHigh(0),(
int)ginfo.get(i).GDbox.getHigh(1),(
int)ginfo.get(i).GDbox.getHigh(2)};
209 int vlo[3] = {0,0,0};
210 int vhi[3] = {(
int)ginfo.get(i).GDbox.getHigh(0),(
int)ginfo.get(i).GDbox.getHigh(1),(
int)ginfo.get(i).GDbox.getHigh(2)};
211 int nvlo[3] = {0,0,0};
212 int nvhi[3] = {(
int)ginfo.get(i).GDbox.getHigh(0),(
int)ginfo.get(i).GDbox.getHigh(1),(
int)ginfo.get(i).GDbox.getHigh(2)};
215 (
double *)U_old.getPointer(),ulo,uhi,
216 (
double *)V_old.getPointer(),vlo,vhi,
217 (
double *)U_new.getPointer(),nulo,nuhi,
218 (
double *)V_new.getPointer(),nulo,nvhi,
219 &deltaT, &uFactor_s, &vFactor_s,&
F,&K);
229int main(
int argc,
char* argv[])
231 openfpm_init(&argc,&argv);
237 size_t sz[3] = {128,128,128};
255 size_t timeSteps = 5000;
285 double spacing[3] = {OldU.spacing(0),OldU.spacing(1),OldU.spacing(2)};
287 init(OldU,OldV,NewU,NewV,domain);
293 OldU.template ghost_get<0>();
294 OldV.template ghost_get<0>();
298 double uFactor = deltaT * du/(spacing[x]*spacing[x]);
299 double vFactor = deltaT * dv/(spacing[x]*spacing[x]);
304 auto & v_cl = create_vcluster();
314 for (
size_t i = 0; i < timeSteps; ++i)
316 if (i % 100 == 0 && v_cl.rank() == 0)
318 std::cout <<
"STEP: " << i << std::endl;
356 step(OldU,OldV,NewU,NewV,star_stencil_3D,uFactor,vFactor,deltaT,
F,K);
363 step(NewU,NewV,OldU,OldV,star_stencil_3D,uFactor,vFactor,deltaT,
F,K);
395 if (create_vcluster().rank() == 0)
396 {std::cout <<
"Total simulation: " << tot_sim.
getwct() << std::endl;}
This class represent an N-dimensional box.
This is a distributed grid.
grid_key_dx is the key to access any element in the grid
__device__ __host__ index_type get(index_type i) const
Get the i index.
Class for cpu time benchmarking.
void stop()
Stop 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
OutputIteratorT OffsetT ReductionOpT OuputT init
< [in] The initial value of the reduction
[v_transform metafunction]
aggregate of properties, from a list of object if create a struct that follow the OPENFPM native stru...