FrictionQPotSpringBlock v0.22.7
Loading...
Searching...
No Matches
Line1d.h
Go to the documentation of this file.
1
7#ifndef FRICTIONQPOTSPRINGBLOCK_LINE1D_H
8#define FRICTIONQPOTSPRINGBLOCK_LINE1D_H
9
10#include "config.h"
11#include "detail.h"
12
13#include <prrng.h>
14
15#include <GMatTensor/version.h>
16
18
22namespace Line1d {
23
34inline std::vector<std::string> version_dependencies()
35{
36 return GMatTensor::version_dependencies();
37}
38
43inline std::vector<std::string> version_compiler()
44{
45 return GMatTensor::version_compiler();
46}
47
51using Generator =
52 prrng::pcg32_tensor_cumsum<array_type::tensor<double, 2>, array_type::tensor<ptrdiff_t, 1>, 1>;
53
113 : public detail::System<1, detail::Cuspy<Generator>, Generator, detail::Laplace1d> {
114protected:
118
119public:
135 double m,
136 double eta,
137 double mu,
138 double k_interactions,
139 double k_frame,
140 double dt,
141 const std::array<size_t, 1>& shape,
142 uint64_t seed,
143 const std::string& distribution,
144 const std::vector<double>& parameters,
145 double offset = -100.0,
146 size_t nchunk = 5000
147 )
148 : m_gen(
149 std::array<size_t, 1>{nchunk},
150 xt::eval(seed + xt::arange<uint64_t>(shape[0])),
151 xt::eval(xt::zeros<uint64_t>(shape)),
152 detail::string_to_distribution(distribution),
153 parameters,
154 prrng::alignment(/*buffer*/ 2, /*margin*/ 30, /*min_margin*/ 6, /*strict*/ false)
155 )
156 {
157 m_gen += offset;
159 m_int = detail::Laplace1d(k_interactions, shape[0]);
160 this->initSystem(m, eta, k_frame, mu, dt, &m_pot, &m_gen, &m_int);
161 }
162};
163
174 1,
175 detail::Cuspy<Generator>,
176 Generator,
177 detail::Laplace1d,
178 void,
179 detail::Overdamped> {
180protected:
184
185public:
200 double mu,
201 double k_interactions,
202 double k_frame,
203 const std::array<size_t, 1>& shape,
204 uint64_t seed,
205 const std::string& distribution,
206 const std::vector<double>& parameters,
207 double offset = -100.0,
208 size_t nchunk = 5000,
209 double eta = 0.0,
210 double dt = 0.0
211 )
212 : m_gen(
213 std::array<size_t, 1>{nchunk},
214 xt::eval(seed + xt::arange<uint64_t>(shape[0])),
215 xt::eval(xt::zeros<uint64_t>(shape)),
216 detail::string_to_distribution(distribution),
217 parameters,
218 prrng::alignment(/*buffer*/ 2, /*margin*/ 30, /*min_margin*/ 6, /*strict*/ false)
219 )
220 {
221 m_gen += offset;
223 m_int = detail::Laplace1d(k_interactions, shape[0]);
224 this->initSystem(1.0, eta, k_frame, mu, dt, &m_pot, &m_gen, &m_int);
225 }
226
227protected:
231 void timeStep();
232 void timeSteps(size_t);
233 size_t timeStepsUntilEvent(double, size_t, size_t);
234 void flowSteps(size_t, double);
238};
239
262 1,
263 detail::Cuspy<Generator>,
264 Generator,
265 detail::Laplace1d,
266 detail::RandomNormalForcing<1>,
267 detail::None> {
268protected:
273
274public:
284 double m,
285 double eta,
286 double mu,
287 double k_interactions,
288 double k_frame,
289 double dt,
290 double mean,
291 double stddev,
292 uint64_t seed_forcing,
293 const array_type::tensor<ptrdiff_t, 1>& dinc_init,
295 const std::array<size_t, 1>& shape,
296 uint64_t seed,
297 const std::string& distribution,
298 const std::vector<double>& parameters,
299 double offset = -100.0,
300 size_t nchunk = 5000
301 )
302 : m_gen(
303 std::array<size_t, 1>{nchunk},
304 xt::eval(seed + xt::arange<uint64_t>(shape[0])),
305 xt::eval(xt::zeros<uint64_t>(shape)),
306 detail::string_to_distribution(distribution),
307 parameters,
308 prrng::alignment(/*buffer*/ 2, /*margin*/ 30, /*min_margin*/ 6, /*strict*/ false)
309 )
310 {
311 m_gen += offset;
313 m_int = detail::Laplace1d(k_interactions, shape[0]);
315 m_gen.generators().shape(), mean, stddev, seed_forcing, dinc_init, dinc
316 );
317 this->initSystem(m, eta, k_frame, mu, dt, &m_pot, &m_gen, &m_int, &m_ext);
318 }
319
320protected:
324 size_t quasistaticActivityFirst() const;
325 size_t quasistaticActivityLast() const;
326 double eventDrivenStep(double, bool, int);
330};
331
337 : public detail::System<1, detail::SemiSmooth<Generator>, Generator, detail::Laplace1d> {
338protected:
342
343public:
349 double m,
350 double eta,
351 double mu,
352 double kappa,
353 double k_interactions,
354 double k_frame,
355 double dt,
356 const std::array<size_t, 1>& shape,
357 uint64_t seed,
358 const std::string& distribution,
359 const std::vector<double>& parameters,
360 double offset = -100.0,
361 size_t nchunk = 5000
362 )
363 : m_gen(
364 std::array<size_t, 1>{nchunk},
365 xt::eval(seed + xt::arange<uint64_t>(shape[0])),
366 xt::eval(xt::zeros<uint64_t>(shape)),
367 detail::string_to_distribution(distribution),
368 parameters,
369 prrng::alignment(/*buffer*/ 2, /*margin*/ 30, /*min_margin*/ 6, /*strict*/ false)
370 )
371 {
372 m_gen += offset;
374 m_int = detail::Laplace1d(k_interactions, shape[0]);
375 this->initSystem(m, eta, k_frame, mu, dt, &m_pot, &m_gen, &m_int);
376 }
377};
378
384 : public detail::System<1, detail::Smooth<Generator>, Generator, detail::Laplace1d> {
385protected:
389
390public:
395 double m,
396 double eta,
397 double mu,
398 double k_interactions,
399 double k_frame,
400 double dt,
401 const std::array<size_t, 1>& shape,
402 uint64_t seed,
403 const std::string& distribution,
404 const std::vector<double>& parameters,
405 double offset = -100.0,
406 size_t nchunk = 5000
407 )
408 : m_gen(
409 std::array<size_t, 1>{nchunk},
410 xt::eval(seed + xt::arange<uint64_t>(shape[0])),
411 xt::eval(xt::zeros<uint64_t>(shape)),
412 detail::string_to_distribution(distribution),
413 parameters,
414 prrng::alignment(/*buffer*/ 2, /*margin*/ 30, /*min_margin*/ 6, /*strict*/ false)
415 )
416 {
417 m_gen += offset;
419 m_int = detail::Laplace1d(k_interactions, shape[0]);
420 this->initSystem(m, eta, k_frame, mu, dt, &m_pot, &m_gen, &m_int);
421 }
422};
423
429 : public detail::System<1, detail::Cuspy<Generator>, Generator, detail::Quartic1d> {
430protected:
434
435public:
452 double m,
453 double eta,
454 double mu,
455 double a1,
456 double a2,
457 double k_frame,
458 double dt,
459 const std::array<size_t, 1>& shape,
460 uint64_t seed,
461 const std::string& distribution,
462 const std::vector<double>& parameters,
463 double offset = -100.0,
464 size_t nchunk = 5000
465 )
466 : m_gen(
467 std::array<size_t, 1>{nchunk},
468 xt::eval(seed + xt::arange<uint64_t>(shape[0])),
469 xt::eval(xt::zeros<uint64_t>(shape)),
470 detail::string_to_distribution(distribution),
471 parameters,
472 prrng::alignment(/*buffer*/ 2, /*margin*/ 30, /*min_margin*/ 6, /*strict*/ false)
473 )
474 {
475 m_gen += offset;
477 m_int = detail::Quartic1d(a1, a2, shape[0]);
478 this->initSystem(m, eta, k_frame, mu, dt, &m_pot, &m_gen, &m_int);
479 }
480};
481
487 1,
488 detail::Cuspy<Generator>,
489 Generator,
490 detail::Quartic1d,
491 detail::RandomNormalForcing<1>,
492 detail::None> {
493protected:
498
499public:
509 double m,
510 double eta,
511 double mu,
512 double a1,
513 double a2,
514 double k_frame,
515 double dt,
516 double mean,
517 double stddev,
518 uint64_t seed_forcing,
519 const array_type::tensor<ptrdiff_t, 1>& dinc_init,
521 const std::array<size_t, 1>& shape,
522 uint64_t seed,
523 const std::string& distribution,
524 const std::vector<double>& parameters,
525 double offset = -100.0,
526 size_t nchunk = 5000
527 )
528 : m_gen(
529 std::array<size_t, 1>{nchunk},
530 xt::eval(seed + xt::arange<uint64_t>(shape[0])),
531 xt::eval(xt::zeros<uint64_t>(shape)),
532 detail::string_to_distribution(distribution),
533 parameters,
534 prrng::alignment(/*buffer*/ 2, /*margin*/ 30, /*min_margin*/ 6, /*strict*/ false)
535 )
536 {
537 m_gen += offset;
539 m_int = detail::Quartic1d(a1, a2, shape[0]);
541 m_gen.generators().shape(), mean, stddev, seed_forcing, dinc_init, dinc
542 );
543 this->initSystem(m, eta, k_frame, mu, dt, &m_pot, &m_gen, &m_int, &m_ext);
544 }
545
546protected:
550 size_t quasistaticActivityFirst() const;
551 size_t quasistaticActivityLast() const;
552 double eventDrivenStep(double, bool, int);
556};
557
563 : public detail::System<1, detail::Cuspy<Generator>, Generator, detail::QuarticGradient1d> {
564protected:
568
569public:
586 double m,
587 double eta,
588 double mu,
589 double k2,
590 double k4,
591 double k_frame,
592 double dt,
593 const std::array<size_t, 1>& shape,
594 uint64_t seed,
595 const std::string& distribution,
596 const std::vector<double>& parameters,
597 double offset = -100.0,
598 size_t nchunk = 5000
599 )
600 : m_gen(
601 std::array<size_t, 1>{nchunk},
602 xt::eval(seed + xt::arange<uint64_t>(shape[0])),
603 xt::eval(xt::zeros<uint64_t>(shape)),
604 detail::string_to_distribution(distribution),
605 parameters,
606 prrng::alignment(/*buffer*/ 2, /*margin*/ 30, /*min_margin*/ 6, /*strict*/ false)
607 )
608 {
609 m_gen += offset;
612 this->initSystem(m, eta, k_frame, mu, dt, &m_pot, &m_gen, &m_int);
613 }
614};
615
621 : public detail::System<1, detail::Cuspy<Generator>, Generator, detail::LongRange1d> {
622protected:
626
627public:
644 double m,
645 double eta,
646 double mu,
647 double k_interactions,
648 double alpha,
649 double k_frame,
650 double dt,
651 const std::array<size_t, 1>& shape,
652 uint64_t seed,
653 const std::string& distribution,
654 const std::vector<double>& parameters,
655 double offset = -100.0,
656 size_t nchunk = 5000
657 )
658 : m_gen(
659 std::array<size_t, 1>{nchunk},
660 xt::eval(seed + xt::arange<uint64_t>(shape[0])),
661 xt::eval(xt::zeros<uint64_t>(shape)),
662 detail::string_to_distribution(distribution),
663 parameters,
664 prrng::alignment(/*buffer*/ 2, /*margin*/ 30, /*min_margin*/ 6, /*strict*/ false)
665 )
666 {
667 m_gen += offset;
669 m_int = detail::LongRange1d(k_interactions, alpha, shape[0]);
670 this->initSystem(m, eta, k_frame, mu, dt, &m_pot, &m_gen, &m_int);
671 }
672};
673
674} // namespace Line1d
675} // namespace FrictionQPotSpringBlock
676
677#endif
System_Cuspy_Laplace() assuming overdamped dynamics.
Definition Line1d.h:179
Generator m_gen
Pointer to chunk of yield 'positions' (automatically updated if needed)
Definition Line1d.h:181
System_Cuspy_Laplace_Nopassing(double mu, double k_interactions, double k_frame, const std::array< size_t, 1 > &shape, uint64_t seed, const std::string &distribution, const std::vector< double > &parameters, double offset=-100.0, size_t nchunk=5000, double eta=0.0, double dt=0.0)
Definition Line1d.h:199
detail::Cuspy< Generator > m_pot
Class to get the forces from the local potential energy landscape.
Definition Line1d.h:182
detail::Laplace1d m_int
Class to get the forces from particle interaction.
Definition Line1d.h:183
System in which the effect of temperature in mimicked by random forcing.
Definition Line1d.h:267
Generator m_gen
Pointer to chunk of yield 'positions' (automatically updated if needed)
Definition Line1d.h:269
detail::Laplace1d m_int
Class to get the forces from particle interaction.
Definition Line1d.h:271
System_Cuspy_Laplace_RandomForcing(double m, double eta, double mu, double k_interactions, double k_frame, double dt, double mean, double stddev, uint64_t seed_forcing, const array_type::tensor< ptrdiff_t, 1 > &dinc_init, const array_type::tensor< ptrdiff_t, 1 > &dinc, const std::array< size_t, 1 > &shape, uint64_t seed, const std::string &distribution, const std::vector< double > &parameters, double offset=-100.0, size_t nchunk=5000)
Definition Line1d.h:283
detail::RandomNormalForcing< 1 > m_ext
Add extra random force to the residual.
Definition Line1d.h:272
detail::Cuspy< Generator > m_pot
Class to get the forces from the local potential energy landscape.
Definition Line1d.h:270
Standard system with a cuspy potential energy landscape and short range interactions.
Definition Line1d.h:113
Generator m_gen
Pointer to chunk of yield 'positions' (automatically updated if needed)
Definition Line1d.h:115
detail::Cuspy< Generator > m_pot
Class to get the forces from the local potential energy landscape.
Definition Line1d.h:116
detail::Laplace1d m_int
Class to get the forces from particle interaction.
Definition Line1d.h:117
System_Cuspy_Laplace(double m, double eta, double mu, double k_interactions, double k_frame, double dt, const std::array< size_t, 1 > &shape, uint64_t seed, const std::string &distribution, const std::vector< double > &parameters, double offset=-100.0, size_t nchunk=5000)
Definition Line1d.h:134
Same as System_Cuspy_Laplace() but with a quartic interactions.
Definition Line1d.h:621
detail::LongRange1d m_int
Class to get the forces from particle interaction.
Definition Line1d.h:625
Generator m_gen
Pointer to chunk of yield 'positions' (automatically updated if needed)
Definition Line1d.h:623
detail::Cuspy< Generator > m_pot
Class to get the forces from the local potential energy landscape.
Definition Line1d.h:624
System_Cuspy_LongRange(double m, double eta, double mu, double k_interactions, double alpha, double k_frame, double dt, const std::array< size_t, 1 > &shape, uint64_t seed, const std::string &distribution, const std::vector< double > &parameters, double offset=-100.0, size_t nchunk=5000)
Definition Line1d.h:643
Same as System_Cuspy_Laplace() but with a quartic interactions.
Definition Line1d.h:563
detail::QuarticGradient1d m_int
Class to get the forces from particle interaction.
Definition Line1d.h:567
Generator m_gen
Pointer to chunk of yield 'positions' (automatically updated if needed)
Definition Line1d.h:565
System_Cuspy_QuarticGradient(double m, double eta, double mu, double k2, double k4, double k_frame, double dt, const std::array< size_t, 1 > &shape, uint64_t seed, const std::string &distribution, const std::vector< double > &parameters, double offset=-100.0, size_t nchunk=5000)
Definition Line1d.h:585
detail::Cuspy< Generator > m_pot
Class to get the forces from the local potential energy landscape.
Definition Line1d.h:566
See System_Cuspy_Quartic() and System_Cuspy_Laplace_RandomForcing()
Definition Line1d.h:492
detail::Cuspy< Generator > m_pot
Class to get the forces from the local potential energy landscape.
Definition Line1d.h:495
detail::Quartic1d m_int
Class to get the forces from particle interaction.
Definition Line1d.h:496
detail::RandomNormalForcing< 1 > m_ext
Add extra random force to the residual.
Definition Line1d.h:497
Generator m_gen
Pointer to chunk of yield 'positions' (automatically updated if needed)
Definition Line1d.h:494
System_Cuspy_Quartic_RandomForcing(double m, double eta, double mu, double a1, double a2, double k_frame, double dt, double mean, double stddev, uint64_t seed_forcing, const array_type::tensor< ptrdiff_t, 1 > &dinc_init, const array_type::tensor< ptrdiff_t, 1 > &dinc, const std::array< size_t, 1 > &shape, uint64_t seed, const std::string &distribution, const std::vector< double > &parameters, double offset=-100.0, size_t nchunk=5000)
Definition Line1d.h:508
Same as System_Cuspy_Laplace() but with a quartic interactions.
Definition Line1d.h:429
Generator m_gen
Pointer to chunk of yield 'positions' (automatically updated if needed)
Definition Line1d.h:431
System_Cuspy_Quartic(double m, double eta, double mu, double a1, double a2, double k_frame, double dt, const std::array< size_t, 1 > &shape, uint64_t seed, const std::string &distribution, const std::vector< double > &parameters, double offset=-100.0, size_t nchunk=5000)
Definition Line1d.h:451
detail::Cuspy< Generator > m_pot
Class to get the forces from the local potential energy landscape.
Definition Line1d.h:432
detail::Quartic1d m_int
Class to get the forces from particle interaction.
Definition Line1d.h:433
Same as System_Cuspy_Laplace() but with a semi-smooth potential.
Definition Line1d.h:337
detail::SemiSmooth< Generator > m_pot
Class to get the forces from the local potential energy landscape.
Definition Line1d.h:340
detail::Laplace1d m_int
Class to get the forces from particle interaction.
Definition Line1d.h:341
Generator m_gen
Pointer to chunk of yield 'positions' (automatically updated if needed)
Definition Line1d.h:339
System_SemiSmooth_Laplace(double m, double eta, double mu, double kappa, double k_interactions, double k_frame, double dt, const std::array< size_t, 1 > &shape, uint64_t seed, const std::string &distribution, const std::vector< double > &parameters, double offset=-100.0, size_t nchunk=5000)
Definition Line1d.h:348
Same as System_Cuspy_Laplace() but with a smooth potential.
Definition Line1d.h:384
detail::Smooth< Generator > m_pot
Class to get the forces from the local potential energy landscape.
Definition Line1d.h:387
detail::Laplace1d m_int
Class to get the forces from particle interaction.
Definition Line1d.h:388
System_Smooth_Laplace(double m, double eta, double mu, double k_interactions, double k_frame, double dt, const std::array< size_t, 1 > &shape, uint64_t seed, const std::string &distribution, const std::vector< double > &parameters, double offset=-100.0, size_t nchunk=5000)
Definition Line1d.h:394
Generator m_gen
Pointer to chunk of yield 'positions' (automatically updated if needed)
Definition Line1d.h:386
A piece-wise quadratic local potential energy.
Definition detail.h:114
Short range elastic interactions with other particles.
Definition detail.h:446
Short range interaction based on a quartic potential.
Definition detail.h:757
Short range interaction based on a quartic potential.
Definition detail.h:616
Each particle experiences a random force representing the effect of temperature.
Definition detail.h:882
A potential energy landscape of each particle that is piecewise smooth.
Definition detail.h:213
A potential energy landscape of each particle that is smooth.
Definition detail.h:357
System in generic number of dimensions.
Definition detail.h:1053
void flowSteps(size_t n, double v_frame)
Make a number of steps with the frame moving at a constant velocity.
Definition detail.h:1637
void initSystem(double m, double eta, double k_frame, double mu, double dt, detail::Cuspy< Generator > *potential, Generator *chunk, detail::Laplace1d *interactions=nullptr, void *external=nullptr)
Initialise the system.
Definition detail.h:1096
size_t timeStepsUntilEvent(double tol=1e-5, size_t niter_tol=10, size_t max_iter=1e9)
Perform a series of time-steps until the next plastic event, or equilibrium.
Definition detail.h:1595
std::vector< std::string > version_compiler()
Return information on the compiler, platform, C++ standard, and the compilation data.
Definition Line1d.h:43
prrng::pcg32_tensor_cumsum< array_type::tensor< double, 2 >, array_type::tensor< ptrdiff_t, 1 >, 1 > Generator
Chunked storage of the cumulative sum of random numbers, used in all classes.
Definition Line1d.h:51
std::vector< std::string > version_dependencies()
Return versions of this library and of all of its dependencies.
Definition Line1d.h:34
xt::xarray< T > array
Arbitrary rank array.
Definition config.h:178
Tensor products / operations.
Definition config.h:145