OpenFPM_pdata  4.1.0
Project that contain the implementation of distributed structures
 
Loading...
Searching...
No Matches
main2.cpp
1//
2// Created by Abhinav Singh @absingh on 26/04/2021
3//
41
42
43// Include Vector Expression,Vector Expressions for Subset,DCPSE,Odeint header files
44#include "Operators/Vector/vector_dist_operators.hpp"
45#include "Vector/vector_dist_subset.hpp"
46#include "DCPSE/DCPSE_op/DCPSE_op.hpp"
47#include "OdeIntegrators/OdeIntegrators.hpp"
49
50
73
74constexpr int x = 0;
75constexpr int y = 1;
76
77double dt=1e-2,tf=1.0,vf=1.0;
78
79void *PointerDistGlobal, *PointerDistSubset,*PointerDistSubset2;
80
85
115template<typename DXX,typename DYY>
116struct RHSFunctor
117{
118 //Intializing the operators
119 DXX &Dxx;
120 DYY &Dyy;
121 //Constructor
122 RHSFunctor(DXX &Dxx,DYY &Dyy):Dxx(Dxx),Dyy(Dyy)
123 {}
124
125 void operator()( const state_type_2d_ofp &X , state_type_2d_ofp &dxdt , const double t ) const
126 {
127 //Casting the pointers to OpenFPM vector distributions
128 dist_vector_type &Particles= *(dist_vector_type *) PointerDistGlobal;
129 dist_vector_subset_type &Particles_bulk= *(dist_vector_subset_type *) PointerDistSubset;
130
131 //Aliasing the properties.
132 auto C = getV<0>(Particles);
133 auto C_bulk = getV<0>(Particles_bulk);
134 auto dC = getV<2>(Particles);
135 auto dC_bulk = getV<2>(Particles_bulk);
136
137 //Since the state vector does not have enough information to compute derivates, we copy them back to the Particle Distribution.
138 //In this way we can compute all the things that are required to be computed in each stage of the ode integrator.
139 // Note that we are going to use bulk expressions as needed. Whenever we don't want to update the boundaries (Due to the boundary conditions).
140
141 //These expressions only update the bulk values of C.
142 C_bulk[x]=X.data.get<0>();
143 C_bulk[y]=X.data.get<1>();
144 Particles.ghost_get<0>();
145
146 // We do the RHS computations for the Laplacian (Updating bulk only).
147 dC_bulk[x] = 0.1*(Dxx(C[x])+Dyy(C[x]));
148 dC_bulk[y] = 0.1*(Dxx(C[y])+Dyy(C[y]));
149
150 //We copy back to the dxdt state_type for Odeint
151 dxdt.data.get<0>()=dC[x];
152 dxdt.data.get<1>()=dC[y];
153 }
154};
156
157
186template<typename DXX,typename DYY>
188
189 DXX &Dxx;
190 DYY &Dyy;
191 int ctr;
192 double t_old;
193
194 //Constructor
195 ObserverFunctor(DXX &Dxx, DYY &Dyy) : Dxx(Dxx), Dyy(Dyy) {
196 //a counter for counting the np. of steps
197 ctr = 0;
198 //Starting with t=0, we compute the step size take by t-t_old. So for the first observed step size is what we provide. Which will be 0-(-dt)=dt.
199 t_old = -dt;
200 }
201
202 void operator()(state_type_2d_ofp &X, double t) {
203 //Casting the pointers to OpenFPM vector distributions
204 dist_vector_type &Particles = *(dist_vector_type *) PointerDistGlobal;
205 dist_vector_subset_type &Particles_bulk = *(dist_vector_subset_type *) PointerDistSubset;
206 dist_vector_subset_type &Particles_boundary = *(dist_vector_subset_type *) PointerDistSubset2;
207
208 //Aliasing the position and properties.
209 auto Pos = getV<PROP_POS>(Particles);
210 auto Concentration = getV<0>(Particles);
211 auto Velocity = getV<1>(Particles);
212 auto Concentration_bulk = getV<0>(Particles_bulk);
213 auto Velocity_bulk = getV<1>(Particles_bulk);
214
215 //We would like to move the particles after t=0. (Remember, odeint calls the observer before taking the step.)
216 if (t != 0) {
217 //copy back the state after the time step into data structure. This is required as we are going to move the particles and the distributed state can be resized correctly (by copying back after map). Also this expression preserves the boundary condition on concentration.
218 Concentration_bulk[x] = X.data.get<0>();
219 Concentration_bulk[y] = X.data.get<1>();
220 //computing the velocity and move the particles
221 Velocity_bulk[x] = -vf*Pos[y] * exp(-10.0 * (Pos[x] * Pos[x] + Pos[y] * Pos[y]));
222 Velocity_bulk[y] = vf*Pos[x] * exp(-10.0 * (Pos[x] * Pos[x] + Pos[y] * Pos[y]));
223 Pos = Pos + dt * Velocity;
224 //Map and ghost_get is required after moving particles.
225 Particles.map();
226 Particles.ghost_get<0>();
227 //Updating the subset and operators based on new positions
228 Particles_bulk.update();
229 Particles_boundary.update();
230 Dxx.update(Particles);
231 Dyy.update(Particles);
232
233 //Since we did Map, we assign the Odeint state again.
234 X.data.get<0>() = Concentration[x];
235 X.data.get<1>() = Concentration[y];
236 //counting the step number
237
238 }
239 ctr++;
240 std::cout<<"Taking a step at t="<<t<<" with dt="<<t-t_old<<std::endl;
241 t_old=t;
242 //Writing the data after every step
243 Particles.deleteGhost();
244 Particles.write_frame("PDE_sol", ctr);
245 Particles.ghost_get<0>();
246
247 }
248};
250
262int main(int argc, char *argv[])
263{
264 // initialize library
265 openfpm_init(&argc, &argv);
266
267 dt=std::atof(argv[1]);
268 tf=std::atof(argv[2]);
269 vf=std::atof(argv[3]);
271
285 const size_t sz[2] = {41, 41};
286 Box<2, double> box({-1, -1}, {1.0, 1.0});
287 size_t bc[2] = {NON_PERIODIC, NON_PERIODIC};
288 double spacing[2];
289 spacing[0] = 2.0 / (sz[0] - 1);
290 spacing[1] = 2.0 / (sz[1] - 1);
291 double rCut = 3.1 * spacing[0];
292 Ghost<2, double> ghost(rCut);
293
294 dist_vector_type Particles(0, box, bc, ghost);
295 Particles.setPropNames({"Concentration", "Velocity", "TempConcentration"});
296
297 auto it = Particles.getGridIterator(sz);
298 while (it.isNext()) {
299 Particles.add();
300 auto key = it.get();
301 double x = -1.0 + key.get(0) * spacing[0];
302 Particles.getLastPos()[0] = x;
303 double y = -1.0 + key.get(1) * spacing[1];
304 Particles.getLastPos()[1] = y;
305 //We detect Boundaries. By labelling the subset to be 0 for bulk and 1 for the boundary. Subsets will be constructed later based on the label.
306 if (x != -1.0 && x != 1.0 && y != -1.0 && y != 1) {
307 Particles.getLastSubset(0);
308 } else {
309 Particles.getLastSubset(1);
310 }
311 // Here fill the Initial value of the concentration.
312 if (x == 0.0 && y > -0.5 && y < 0) {
313 Particles.template getLastProp<0>()[0] = 1.0;
314 Particles.template getLastProp<0>()[1] = 0.0;
315 } else if (x == 0.0 && y > 0 && y < 0.5) {
316 Particles.template getLastProp<0>()[0] = 0.0;
317 Particles.template getLastProp<0>()[1] = 1.0;
318 } else {
319 Particles.template getLastProp<0>()[0] = 0.0;
320 Particles.template getLastProp<0>()[1] = 0.0;
321 }
322 ++it;
323 }
324 Particles.map();
325 Particles.ghost_get<0>();
326
327 //We write the particles to check if the initialization is correct.
328 Particles.write("Init");
330
343 // Now we initialize the grid with a filled circle. Outside the circle, the value of Phi_0 will be -1, inside +1.
344 //Now we construct the subsets based on the subset number.
345 dist_vector_subset_type Particles_bulk(Particles, 0);
346 dist_vector_subset_type Particles_boundary(Particles, 1);
347
348 //We cast the global pointers to Particles and Particles_bulk as expected by the RHS functor.
349 PointerDistGlobal = (void *) &Particles;
350 PointerDistSubset = (void *) &Particles_bulk;
351 PointerDistSubset2 = (void *) &Particles_boundary;
353
365 //We create the DCPSE Based operators for discretization of the operators.
366 Derivative_xx Dxx(Particles, 2, rCut);
367 Derivative_yy Dyy(Particles, 2, rCut);
368 //We create aliases for referring to property and and positions.
369 auto Pos = getV<PROP_POS>(Particles);
370 auto C = getV<0>(Particles);
371 auto V = getV<1>(Particles);
372 auto dC = getV<2>(Particles);
373
374 //We create aliases for referring to the subset properties.
375 auto C_bulk = getV<0>(Particles_bulk);
376 auto dC_bulk = getV<2>(Particles_bulk);
378
394 //Now we create a odeint stepper object (RK4). Since we are in 2d, we are going to use "state_type_2d_ofp". Which is a structure or state_type compatible with odeint. We further pass all the parameters including "boost::numeric::odeint::vector_space_algebra_ofp",which tell odeint to use openfpm algebra.
395 // The template parameters are: state_type_2d_ofp (state type of X), double (type of the value inside the state), state_type_2d_ofp (state type of DxDt), double (type of the time), boost::numeric::odeint::vector_space_algebra_ofp (our algebra)
396 boost::numeric::odeint::runge_kutta4<state_type_2d_ofp, double, state_type_2d_ofp, double, boost::numeric::odeint::vector_space_algebra_ofp> Odeint_rk4;
397
398 //The method Odeint_rk4 from Odeint, requires system (a function which computes RHS of the PDE), an instance of the Compute RHS functor. We create the System with the correct types and parameteres for the operators as declared before.
400
401 //Since we are using Odeint to control the time steps, we also create a observer instance. Which also updates the position via an euler step for moving thr particles.
402 ObserverFunctor<Derivative_xx, Derivative_yy> ObserveAndUpdate(Dxx, Dyy);
403
404
405 //Furhter, odeint needs data in a state type "state_type_2d_ofp", we create one and fill in the initial condition.
407 //Since we created a 2d state_type we initialize the two fields in the object data using the method get.
408 X.data.get<x>() = C[0];
409 X.data.get<y>() = C[1];
411
433 std::vector<double> inter_times; // vector to store intermediate time steps taken by odeint.
434
435 size_t steps = boost::numeric::odeint::integrate_const(Odeint_rk4, System, X, 0.0, tf, dt, ObserveAndUpdate);
436 /* From Odeint: (Please refer to Odeint on more details about these steppers.)
437 runge_kutta_dopri5 is maybe the best default stepper. It has step size control as well as dense-output functionality. Simple create a dense-output stepper by make_dense_output( 1.0e-6 , 1.0e-5 , runge_kutta_dopri5< state_type >() ).
438 runge_kutta4 is a good stepper for constant step sizes. It is widely used and very well known. If you need to create artificial time series this stepper should be the first choice.
439 'runge_kutta_fehlberg78' is similar to the 'runge_kutta4' with the advantage that it has higher precision. It can also be used with step size control.
440 adams_bashforth_moulton is very well suited for ODEs where the r.h.s. is expensive (in terms of computation time). It will calculate the system function only once during each step.*/
441 /*Below are listed examples of using different steppers and their usage.*/
442 //size_t steps = integrate_adaptive(controlled_stepper_type() , System , X , t , tf , dt, ObserveAndUpdate);
443 //size_t steps = boost::numeric::odeint::integrate_adaptive( boost::numeric::odeint::make_controlled( 1.0e-7 , 1.0e-7 ,boost::numeric::odeint::runge_kutta_cash_karp54< state_type_2d_ofp,double,state_type_2d_ofp,double,boost::numeric::odeint::vector_space_algebra_ofp>() ) , System , X , t , tf , dt, ObserveAndUpdate );
444 //size_t steps = boost::numeric::odeint::integrate_const( boost::numeric::odeint::make_dense_output(1.0e-7, 1.0e-7, boost::numeric::odeint::runge_kutta_dopri5< state_type_2d_ofp,double,state_type_2d_ofp,double,boost::numeric::odeint::vector_space_algebra_ofp >() ) , System , X , t , tf , dt, ObserveAndUpdate );
445 //size_t steps = boost::numeric::odeint::integrate_adaptive( boost::numeric::odeint::make_dense_output( 1.0e-8 , 1.0e-8 , boost::numeric::odeint::runge_kutta_dopri5< state_type_2d_ofp,double,state_type_2d_ofp,double,boost::numeric::odeint::vector_space_algebra_ofp >() ) , System , X , t , tf , dt, ObserveAndUpdate );
446 //size_t steps = boost::numeric::odeint::integrate_adaptive( boost::numeric::odeint::make_controlled( 1.0e-7 , 1.0e-7 , boost::numeric::odeint::runge_kutta_dopri5< state_type_2d_ofp,double,state_type_2d_ofp,double,boost::numeric::odeint::vector_space_algebra_ofp >() ) , System , X , t , tf , dt, ObserveAndUpdate );
447
448
449 std::cout << "No. of Time steps taken: " << steps << std::endl;
450 //Copying back the final solution and outputting the number of steps taken by Odeint.
451 C_bulk[x] = X.data.get<0>();
452 C_bulk[y] = X.data.get<1>();
453 //Deallocating the operators
454 Dxx.deallocate(Particles);
455 Dyy.deallocate(Particles);
456 openfpm_finalize(); // Finalize openFPM library
457 return 0;
458} //main end
460
This class represent an N-dimensional box.
Definition Box.hpp:61
This class implement the point shape in an N-dimensional space.
Definition Point.hpp:28
System of equations.
Definition System.hpp:24
void update()
Update the subset indexes.
void ghost_get(size_t opt=WITH_POSITION)
It synchronize the properties and position of the ghost particles.
void map(size_t opt=NONE)
It move all the particles that does not belong to the local processor to the respective processor.
void deleteGhost()
Delete the particles on the ghost.
aggregate of properties, from a list of object if create a struct that follow the OPENFPM native stru...
A 2d Odeint and Openfpm compatible structure.