OpenFPM_pdata  4.1.0
Project that contain the implementation of distributed structures
 
Loading...
Searching...
No Matches
main.cpp
1//
2// Created by Abhinav Singh @absingh on 26/04/2021
3//
68
69// Include Vector Expression,Vector Expressions for Subset,DCPSE,Odeint header files
70#include "Operators/Vector/vector_dist_operators.hpp"
71#include "Vector/vector_dist_subset.hpp"
72#include "DCPSE/DCPSE_op/DCPSE_op.hpp"
73#include "OdeIntegrators/OdeIntegrators.hpp"
75
76
98
99constexpr int x = 0;
100constexpr int y = 1;
101
102double dt;
103
104void *PointerDistGlobal, *PointerDistSubset;
105
110
140template<typename DXX,typename DYY>
142{
143 //Intializing the operators
144 DXX &Dxx;
145 DYY &Dyy;
146 //Constructor
147 RHSFunctor(DXX &Dxx,DYY &Dyy):Dxx(Dxx),Dyy(Dyy)
148 {}
149
150 void operator()( const state_type_2d_ofp &X , state_type_2d_ofp &dxdt , const double t ) const
151 {
152 //Casting the pointers to OpenFPM vector distributions
153 dist_vector_type &Particles= *(dist_vector_type *) PointerDistGlobal;
154 dist_vector_subset_type &Particles_bulk= *(dist_vector_subset_type *) PointerDistSubset;
155
156 //Aliasing the properties.
157 auto C = getV<0>(Particles);
158 auto C_bulk = getV<0>(Particles_bulk);
159 auto dC = getV<2>(Particles);
160 auto dC_bulk = getV<2>(Particles_bulk);
161
162 //Since the state vector does not have enough information to compute derivates, we copy them back to the Particle Distribution.
163 //In this way we can compute all the things that are required to be computed in each stage of the ode integrator.
164 // 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).
165
166 //These expressions only update the bulk values of C.
167 C_bulk[x]=X.data.get<0>();
168 C_bulk[y]=X.data.get<1>();
169 Particles.ghost_get<0>();
170
171 // We do the RHS computations for the Laplacian (Updating bulk only).
172 dC_bulk[x] = 0.1*(Dxx(C[x])+Dyy(C[x]));
173 dC_bulk[y] = 0.1*(Dxx(C[y])+Dyy(C[y]));
174
175 //We copy back to the dxdt state_type for Odeint
176 dxdt.data.get<0>()=dC[x];
177 dxdt.data.get<1>()=dC[y];
178 }
179};
181
182
183
196int main(int argc, char* argv[]) {
197 // initialize library
198 openfpm_init(&argc, &argv);
200
202
214 const size_t sz[2] = {41, 41};
215 Box<2, double> box({-1, -1}, {1.0, 1.0});
216 size_t bc[2] = {NON_PERIODIC, NON_PERIODIC};
217 double spacing[2];
218 spacing[0] = 2.0 / (sz[0] - 1);
219 spacing[1] = 2.0 / (sz[1] - 1);
220 double rCut = 3.1 * spacing[0];
221 Ghost<2, double> ghost(rCut);
222
223 dist_vector_type Particles(0, box, bc, ghost);
224 Particles.setPropNames({"Concentration", "Velocity", "TempConcentration"});
225
226 auto it = Particles.getGridIterator(sz);
227 while (it.isNext()) {
228 Particles.add();
229 auto key = it.get();
230 double x = -1.0 + key.get(0) * spacing[0];
231 Particles.getLastPos()[0] = x;
232 double y = -1.0 + key.get(1) * spacing[1];
233 Particles.getLastPos()[1] = y;
234 //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.
235 if (x != -1.0 && x != 1.0 && y != -1.0 && y != 1) {
236 Particles.getLastSubset(0);
237 } else {
238 Particles.getLastSubset(1);
239 }
240 // Here fill the Initial value of the concentration.
241 if (x == 0.0 && y > -0.5 && y < 0) {
242 Particles.template getLastProp<0>()[0] = 1.0;
243 Particles.template getLastProp<0>()[1] = 0.0;
244 } else if (x == 0.0 && y > 0 && y < 0.5) {
245 Particles.template getLastProp<0>()[0] = 0.0;
246 Particles.template getLastProp<0>()[1] = 1.0;
247 } else {
248 Particles.template getLastProp<0>()[0] = 0.0;
249 Particles.template getLastProp<0>()[1] = 0.0;
250 }
251 ++it;
252 }
253 Particles.map();
254 Particles.ghost_get<0>();
255
256 //We write the particles to check if the initialization is correct.
257 Particles.write("Init");
259
272 // Now we initialize the grid with a filled circle. Outside the circle, the value of Phi_0 will be -1, inside +1.
273 //Now we construct the subsets based on the subset number.
274 dist_vector_subset_type Particles_bulk(Particles, 0);
275 dist_vector_subset_type Particles_boundary(Particles, 1);
276
277 //We cast the global pointers to Particles and Particles_bulk as expected by the RHS functor.
278 PointerDistGlobal = (void *) &Particles;
279 PointerDistSubset = (void *) &Particles_bulk;
281
295 //We create the DCPSE Based operators for discretization of the operators.
296 Derivative_xx Dxx(Particles, 2, rCut);
297 Derivative_yy Dyy(Particles, 2, rCut);
298 //We create aliases for referring to property and and positions.
299 auto Pos = getV<PROP_POS>(Particles);
300 auto C = getV<0>(Particles);
301 auto V = getV<1>(Particles);
302 auto dC = getV<2>(Particles);
303
304 //We create aliases for referring to the subset properties.
305 auto C_bulk = getV<0>(Particles_bulk);
306 auto V_bulk = getV<1>(Particles_bulk);
307 auto dC_bulk = getV<2>(Particles_bulk);
308
310
326 //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.
327 // 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)
328 boost::numeric::odeint::runge_kutta4<state_type_2d_ofp, double, state_type_2d_ofp, double, boost::numeric::odeint::vector_space_algebra_ofp> Odeint_rk4;
329
330 //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.
332
333 //Furhter, odeint needs data in a state type "state_type_2d_ofp", we create one and fill in the initial condition.
335 //Since we created a 2d state_type we initialize the two fields in the object data using the method get.
336 X.data.get<x>() = C[x];
337 X.data.get<y>() = C[y];
339
363 // Now we will create a time loop for controlling the step size ourselves but call odeint to do the 4 stages of RK4.
364 int ctr = 0;
365 double t = 0, tf = 1, dt = 1e-2;
366 while (t < tf) {
367 //computing the velocity at the current position and at time step t (only in the bulk, so boundary remains 0).
368 V_bulk[x] = -Pos[y] * exp(-10.0 * (Pos[x] * Pos[x] + Pos[y] * Pos[y]));
369 V_bulk[y] = Pos[x] * exp(-10.0 * (Pos[x] * Pos[x] + Pos[y] * Pos[y]));
370 //Observing the state
371 Particles.deleteGhost();
372 Particles.write_frame("PDE_Sol", ctr);
373 Particles.ghost_get<0>();
374 //calling rk4 with the function System, the state_type, current time t and the stepsize dt. It computes one step of the RK4.
375 Odeint_rk4.do_step(System, X, t, dt);
376 //Copying back the step, updating only the bulk values.
377 C_bulk[x] = X.data.get<0>();
378 C_bulk[y] = X.data.get<1>();
379 //We advect the particles with the velocity using an Euler step.
380 Pos = Pos + dt * V;
381 //We need to map and get the ghost, as we moved the particles.
382 Particles.map();
383 Particles.ghost_get<0>();
384 //We update the subset and operators as the particles moved.
385 Particles_bulk.update();
386 Particles_boundary.update();
387 Dxx.update(Particles);
388 Dyy.update(Particles);
389 //Reinitialzing as distributed size can change.
390 X.data.get<0>()=C[x];
391 X.data.get<1>()=C[y];
392 ctr++;
393 t += dt;
394 } //time loop end
395
396 //Deallocating the operators
397 Dxx.deallocate(Particles);
398 Dyy.deallocate(Particles);
399 openfpm_finalize(); // Finalize openFPM library
400 return 0;
401} //main end
403
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 ghost_get(size_t opt=WITH_POSITION)
It synchronize the properties and position of the ghost particles.
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.