PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
modelData.h
Go to the documentation of this file.
1/*
2 * -------------------------------------------
3 * Copyright (c) 2021 - 2024 Prashant K. Jha
4 * -------------------------------------------
5 * PeriDEM https://github.com/prashjha/PeriDEM
6 *
7 * Distributed under the Boost Software License, Version 1.0. (See accompanying
8 * file LICENSE)
9 */
10
11#ifndef MODEL_MODELDATA_H
12#define MODEL_MODELDATA_H
13
14#include "util/point.h"
15#include "util/matrix.h"
16#include "util/methods.h"
18#include "inp/input.h"
21#include "nsearch/nsearch.h"
22#include "geometry/fracture.h"
23#include <cstdint> // uint8_t type
24#include <cstring> // string and size_t type
25#include <vector>
26#include <map>
27#include <fstream>
28#include <iostream>
29
31
32// forward declare particle and wall
33namespace particle {
34class BaseParticle;
35class RefParticle;
36class Particle;
37}
38
39namespace model {
40
47class ModelData {
48
49public:
55 : d_n(0),
56 d_time(0.),
57 d_currentDt(0.),
58 d_infoN(1),
59 d_input_p(deck),
60 d_modelDeck_p(deck->getModelDeck()),
61 d_restartDeck_p(deck->getRestartDeck()),
62 d_outputDeck_p(deck->getOutputDeck()),
63 d_pDeck_p(deck->getParticleDeck()), d_cDeck_p(deck->getContactDeck()),
64 d_stop(false), d_hMax(0.), d_hMin(0.), d_maxContactR(0.),
68 d_uLoading_p(nullptr), d_fLoading_p(nullptr),
69 d_fracture_p(nullptr), d_nsearch_p(nullptr) {}
70
76 const particle::BaseParticle* getParticleFromAllList(size_t i) const { return
78
82
90
94
102
106
112 double getDensity(size_t i);
113
119 double getHorizon(size_t i);
120
126 size_t &getPtId(size_t i) { return d_ptId[i]; };
127
129 const size_t &getPtId(size_t i) const { return d_ptId[i]; };
130
136 void setPtId(size_t i, const size_t &id) { d_ptId[i] = id; };
137
138
145 double getKeyData(std::string key, bool issue_err = false) {
146 return util::methods::getKeyData<double>(key, d_dbgData, issue_err);
147 };
148
154 void appendKeyData(std::string key, double data, bool issue_err = false) {
155 util::methods::appendKeyData<double>(key, data, d_dbgData, issue_err);
156 };
157
163 void setKeyData(std::string key, double data, bool issue_err = false) {
164 util::methods::setKeyData<double>(key, data, d_dbgData, issue_err);
165 };
166
177 util::Point &getXRef(size_t i) { return d_xRef[i]; };
178
180 const util::Point &getXRef(size_t i) const { return d_xRef[i]; };
181
187 void setXRef(size_t i, const util::Point &x) { d_xRef[i] = x; };
188
194 void addXRef(size_t i, const util::Point &x) { d_xRef[i] += x; };
195
202 void setXRef(size_t i, int dof, double x) { d_xRef[i][dof] = x; };
203
210 void addXRef(size_t i, int dof, double x) { d_xRef[i][dof] += x; };
211
224 util::Point &getX(size_t i) { return d_x[i]; };
225
227 const util::Point &getX(size_t i) const { return d_x[i]; };
228
234 void setX(size_t i, const util::Point &x) { d_x[i] = x; };
235
241 void addX(size_t i, const util::Point &x) { d_x[i] += x; };
242
249 void setX(size_t i, int dof, double x) { d_x[i][dof] = x; };
250
257 void addX(size_t i, int dof, double x) { d_x[i][dof] += x; };
258
271 util::Point &getU(size_t i) { return d_u[i]; };
272
274 const util::Point &getU(size_t i) const { return d_u[i]; };
275
281 void setU(size_t i, const util::Point &u) { d_u[i] = u; };
282
288 void addU(size_t i, const util::Point &u) { d_u[i] += u; };
289
296 void setU(size_t i, int dof, double u) { d_u[i][dof] = u; };
297
304 void addU(size_t i, int dof, double u) { d_u[i][dof] += u; };
305
318 util::Point &getV(size_t i) { return d_v[i]; };
319
321 const util::Point &getV(size_t i) const { return d_v[i]; };
322
328 void setV(size_t i, const util::Point &v) { d_v[i] = v; };
329
335 void addV(size_t i, const util::Point &v) { d_v[i] += v; };
336
343 void setV(size_t i, int dof, double v) { d_v[i][dof] = v; };
344
351 void addV(size_t i, int dof, double v) { d_v[i][dof] += v; };
352
365 util::Point &getF(size_t i) { return d_f[i]; };
366
368 const util::Point &getF(size_t i) const { return d_f[i]; };
369
375 void setF(size_t i, const util::Point &f) { d_f[i] = f; };
376
382 void addF(size_t i, const util::Point &f) { d_f[i] += f; };
383
390 void setF(size_t i, int dof, double f) { d_f[i][dof] = f; };
391
398 void addF(size_t i, int dof, double f) { d_f[i][dof] += f; };
399
412 double &getVol(size_t i) { return d_vol[i]; };
413
415 const double &getVol(size_t i) const { return d_vol[i]; };
416
422 void setVol(size_t i, const double &vol) { d_vol[i] = vol; };
423
429 void addVol(size_t i, const double &vol) { d_vol[i] += vol; };
430
443 uint8_t &getFix(size_t i) { return d_fix[i]; };
444
446 const uint8_t &getFix(size_t i) const { return d_fix[i]; };
447
454 void setFix(size_t i, const unsigned int &dof,
455 const bool &flag) {
456 // to set i^th bit as true of integer a,
457 // a |= 1UL << (i % 8)
458
459 // to set i^th bit as false of integer a,
460 // a &= ~(1UL << (i % 8))
461
462 flag ? (d_fix[i] |= 1UL << dof) : (d_fix[i] &= ~(1UL << dof));
463 };
464
477 double &getMx(size_t i) { return d_mX[i]; };
478
480 const double &getMx(size_t i) const { return d_mX[i]; };
481
487 void setMx(size_t i, const double &mx) { d_mX[i] = mx; };
488
494 void addMx(size_t i, const double &mx) { d_mX[i] += mx; };
495
508 double &getThetax(size_t i) { return d_thetaX[i]; };
509
511 const double &getThetax(size_t i) const { return d_thetaX[i]; };
512
518 void setThetax(size_t i, const double &thetax) { d_thetaX[i] = thetax; };
519
525 void addThetax(size_t i, const double &thetax) { d_thetaX[i] += thetax; };
526
529public:
531 size_t d_n;
532
534 double d_time;
535
538
540 size_t d_infoN;
541
543 std::map<std::string, double> d_dbgData;
544
546 std::ofstream d_ppFile;
547
550
552 std::shared_ptr<inp::ModelDeck> d_modelDeck_p;
553
555 std::shared_ptr<inp::RestartDeck> d_restartDeck_p;
556
558 std::shared_ptr<inp::OutputDeck> d_outputDeck_p;
559
561 std::shared_ptr<inp::ParticleDeck> d_pDeck_p;
562
564 std::shared_ptr<inp::ContactDeck> d_cDeck_p;
565
566
568 bool d_stop;
569
571 double d_hMax;
572
574 double d_hMin;
575
578
581
584
587
589 std::vector<std::shared_ptr<particle::RefParticle>> d_referenceParticles;
590
592 std::vector<particle::BaseParticle*> d_particlesListTypeAll;
593
595 std::vector<particle::BaseParticle*> d_particlesListTypeParticle;
596
598 std::vector<particle::BaseParticle*> d_particlesListTypeWall;
599
604 std::vector<inp::MatData> d_particlesMatDataList;
605
608
611
615 std::vector<std::vector<size_t>> d_zInfo;
616
618 std::unique_ptr<loading::ParticleULoading> d_uLoading_p;
619
621 std::unique_ptr<loading::ParticleFLoading> d_fLoading_p;
622
624 std::unique_ptr<geometry::Fracture> d_fracture_p;
625
627 std::unique_ptr<NSearch> d_nsearch_p;
628
630 std::vector<util::Point> d_xRef;
631
633 std::vector<util::Point> d_x;
634
636 std::vector<util::Point> d_u;
637
639 std::vector<util::Point> d_v;
640
642 std::vector<double> d_vMag;
643
645 std::vector<util::Point> d_f;
646
648 std::vector<double> d_vol;
649
652 std::vector<size_t> d_ptId;
653
655 std::vector<std::vector<size_t>> d_neighC;
656
658 std::vector<std::vector<size_t>> d_neighPd;
659
661 std::vector<std::vector<float>> d_neighPdSqdDist;
662
664 std::vector<std::vector<std::vector<size_t>>> d_neighWallNodes;
665
667 std::vector<std::vector<std::vector<double>>> d_neighWallNodesDistance;
668
670 std::vector<std::vector<size_t>> d_neighWallNodesCondensed;
671
680 std::vector<uint8_t> d_fix;
681
683 std::vector<uint8_t> d_forceFixity;
684
690 std::vector<double> d_thetaX;
691
696 std::vector<double> d_mX;
697
700 std::vector<size_t> d_fPdCompNodes;
701
703 std::vector<size_t> d_fContCompNodes;
704
706 std::vector<float> d_Z;
707
716 std::vector<float> d_e;
717
719 std::vector<float> d_w;
720
722 std::vector<float> d_phi;
723
725 std::vector<float> d_eF;
726
728 std::vector<float> d_eFB;
729
731 std::vector<util::Point> d_xQuadCur;
732
734 std::vector<util::SymMatrix3> d_strain;
735
737 std::vector<util::SymMatrix3> d_stress;
738
740 float d_te;
741
743 float d_tw;
744
746 float d_tk;
747
749 float d_teF;
750
752 float d_teFB;
753
755};
756
759} // namespace model
760
761#endif // MODEL_MODELDATA_H
A class to read input file.
Definition input.h:61
A class to store model data.
Definition modelData.h:47
std::vector< uint8_t > d_fix
Vector of fixity mask of each node.
Definition modelData.h:680
std::vector< size_t > d_fContCompNodes
List of global nodes on which force (contact) is to be computed.
Definition modelData.h:703
std::shared_ptr< inp::RestartDeck > d_restartDeck_p
Restart deck.
Definition modelData.h:555
std::shared_ptr< inp::ContactDeck > d_cDeck_p
Contact deck.
Definition modelData.h:564
double & getThetax(size_t i)
Get volumetric deformation (thetax) of the node.
Definition modelData.h:508
void addX(size_t i, int dof, double x)
Add to specific current coordinate of the node.
Definition modelData.h:257
void setFix(size_t i, const unsigned int &dof, const bool &flag)
Set fixity of the node.
Definition modelData.h:454
util::Point & getXRef(size_t i)
Get reference coordinate of the node.
Definition modelData.h:177
std::vector< std::vector< float > > d_neighPdSqdDist
Square distance neighbor data for peridynamic forces.
Definition modelData.h:661
void setF(size_t i, int dof, double f)
Set force of the node.
Definition modelData.h:390
void setVol(size_t i, const double &vol)
Set volume of the node.
Definition modelData.h:422
void addXRef(size_t i, const util::Point &x)
Add reference coordinate of the node.
Definition modelData.h:194
double getHorizon(size_t i)
Get horizon of particle.
Definition modelData.cpp:17
std::vector< util::Point > d_f
Total force on the nodes.
Definition modelData.h:645
const double & getThetax(size_t i) const
Get volumetric deformation (thetax) of the node.
Definition modelData.h:511
std::vector< double > d_maxVelocityParticlesListTypeAll
Maximum velocity among all nodes in the particle for each particle.
Definition modelData.h:607
const particle::BaseParticle * getParticleFromWallList(size_t i) const
Get pointer to wall.
Definition modelData.h:100
particle::BaseParticle *& getParticleFromWallList(size_t i)
Definition modelData.h:104
std::vector< float > d_eFB
Bond-based fracture energy of the nodes.
Definition modelData.h:728
std::unique_ptr< NSearch > d_nsearch_p
Pointer to nsearch.
Definition modelData.h:627
const util::Point & getV(size_t i) const
Get velocity of the node.
Definition modelData.h:321
void addU(size_t i, int dof, double u)
Add to displacement of the node.
Definition modelData.h:304
std::vector< std::vector< size_t > > d_neighC
Neighbor data for contact forces.
Definition modelData.h:655
double d_hMax
Minimum mesh over all particles and walls.
Definition modelData.h:571
double getKeyData(std::string key, bool issue_err=false)
Get data for a key.
Definition modelData.h:145
std::vector< particle::BaseParticle * > d_particlesListTypeAll
List of particles + walls.
Definition modelData.h:592
const util::Point & getU(size_t i) const
Get displacement of the node.
Definition modelData.h:274
const double & getVol(size_t i) const
Get volume of the node.
Definition modelData.h:415
double & getVol(size_t i)
Get volume of the node.
Definition modelData.h:412
std::vector< util::Point > d_v
Velocity of the nodes.
Definition modelData.h:639
std::vector< float > d_e
Energy of the nodes.
Definition modelData.h:716
particle::BaseParticle *& getParticleFromParticleList(size_t i)
Definition modelData.h:92
std::vector< util::Point > d_xRef
reference positions of the nodes
Definition modelData.h:630
std::vector< std::vector< std::vector< size_t > > > d_neighWallNodes
Neighbor data for contact between particle and walls.
Definition modelData.h:664
std::unique_ptr< loading::ParticleFLoading > d_fLoading_p
Pointer to force Loading object.
Definition modelData.h:621
void addU(size_t i, const util::Point &u)
Add to displacement of the node.
Definition modelData.h:288
std::vector< std::vector< size_t > > d_neighWallNodesCondensed
Neighbor data for contact between particle and walls condensed into single vector for each particle.
Definition modelData.h:670
uint8_t & getFix(size_t i)
Get fixity of the node.
Definition modelData.h:443
size_t d_infoN
Print log step interval.
Definition modelData.h:540
std::vector< double > d_vol
Nodal volumes.
Definition modelData.h:648
void addMx(size_t i, const double &mx)
Add to weighted-volume (mx) of the node.
Definition modelData.h:494
std::vector< float > d_Z
Damage at nodes.
Definition modelData.h:706
void setThetax(size_t i, const double &thetax)
Set volumetric deformation (thetax) of the node.
Definition modelData.h:518
void setU(size_t i, int dof, double u)
Set displacement of the node.
Definition modelData.h:296
std::vector< util::Point > d_u
Displacement of the nodes.
Definition modelData.h:636
void setU(size_t i, const util::Point &u)
Set displacement of the node.
Definition modelData.h:281
double d_contNeighSearchRadius
Neighborlist contact search radius (multiple of d_maxContactR). This variable will be updated during ...
Definition modelData.h:586
void setXRef(size_t i, int dof, double x)
Set specific reference coordinate of the node.
Definition modelData.h:202
std::vector< std::vector< std::vector< double > > > d_neighWallNodesDistance
Neighbor data (distance) for contact between particle and walls.
Definition modelData.h:667
util::Point & getV(size_t i)
Get velocity of the node.
Definition modelData.h:318
std::shared_ptr< inp::ModelDeck > d_modelDeck_p
Model deck.
Definition modelData.h:552
std::vector< float > d_phi
Damage function at the nodes.
Definition modelData.h:722
const particle::BaseParticle * getParticleFromParticleList(size_t i) const
Get pointer to particle (excluding wall)
Definition modelData.h:88
void addV(size_t i, int dof, double v)
Add to velocity of the node.
Definition modelData.h:351
std::vector< double > d_thetaX
Dilation.
Definition modelData.h:690
const util::Point & getX(size_t i) const
Get current coordinate of the node.
Definition modelData.h:227
std::vector< util::SymMatrix3 > d_stress
Stress in elements (values at quadrature points)
Definition modelData.h:737
ModelData(inp::Input *deck)
Constructor.
Definition modelData.h:54
const size_t & getPtId(size_t i) const
Get particle id given the location in particle list.
Definition modelData.h:129
size_t d_contNeighUpdateInterval
Neighborlist update interval.
Definition modelData.h:580
void setV(size_t i, const util::Point &v)
Set velocity of the node.
Definition modelData.h:328
const util::Point & getXRef(size_t i) const
Get reference coordinate of the node.
Definition modelData.h:180
particle::BaseParticle *& getParticleFromAllList(size_t i)
Get pointer to base particle.
Definition modelData.h:80
size_t d_n
Current time step.
Definition modelData.h:531
std::vector< util::Point > d_x
Current positions of the nodes.
Definition modelData.h:633
std::vector< float > d_w
Work done on each of the nodes.
Definition modelData.h:719
std::unique_ptr< loading::ParticleULoading > d_uLoading_p
Pointer to displacement Loading object.
Definition modelData.h:618
void setXRef(size_t i, const util::Point &x)
Set reference coordinate of the node.
Definition modelData.h:187
std::ofstream d_ppFile
File stream to output information.
Definition modelData.h:546
std::vector< util::Point > d_xQuadCur
Current position of quadrature points.
Definition modelData.h:731
std::vector< double > d_vMag
Magnitude of velocity of the nodes.
Definition modelData.h:642
void addXRef(size_t i, int dof, double x)
Add specific reference coordinate of the node.
Definition modelData.h:210
void setV(size_t i, int dof, double v)
Set velocity of the node.
Definition modelData.h:343
const double & getMx(size_t i) const
Get weighted-volume (mx) of the node.
Definition modelData.h:480
void setX(size_t i, int dof, double x)
Set specific current coordinate of the node.
Definition modelData.h:249
float d_tk
Total kinetic energy.
Definition modelData.h:746
size_t & getPtId(size_t i)
Get particle id given the location in particle list.
Definition modelData.h:126
inp::Input * d_input_p
Pointer to Input object.
Definition modelData.h:549
void setPtId(size_t i, const size_t &id)
Set particle id given the location in particle list.
Definition modelData.h:136
std::vector< size_t > d_ptId
Global node to particle id (walls are assigned id after last particle id)
Definition modelData.h:652
bool d_stop
flag to stop the simulation midway
Definition modelData.h:568
void setF(size_t i, const util::Point &f)
Set force of the node.
Definition modelData.h:375
std::shared_ptr< inp::OutputDeck > d_outputDeck_p
Output deck.
Definition modelData.h:558
void addThetax(size_t i, const double &thetax)
Add to volumetric deformation (thetax) of the node.
Definition modelData.h:525
double d_maxContactR
Maximum contact radius between over pairs of particles and walls.
Definition modelData.h:577
void addVol(size_t i, const double &vol)
Add to volume of the node.
Definition modelData.h:429
void setKeyData(std::string key, double data, bool issue_err=false)
Set value to data associated with key.
Definition modelData.h:163
double d_maxVelocity
Maximum velocity among all nodes.
Definition modelData.h:610
std::map< std::string, double > d_dbgData
Debug data.
Definition modelData.h:543
void addX(size_t i, const util::Point &x)
Add current coordinate of the node.
Definition modelData.h:241
double d_time
Current time.
Definition modelData.h:534
size_t d_contNeighTimestepCounter
Contact neighborlist time step counter.
Definition modelData.h:583
std::vector< uint8_t > d_forceFixity
Vector of fixity mask of each node for force.
Definition modelData.h:683
double & getMx(size_t i)
Get weighted-volume (mx) of the node.
Definition modelData.h:477
util::Point & getX(size_t i)
Get current coordinate of the node.
Definition modelData.h:224
util::Point & getU(size_t i)
Get displacement of the node.
Definition modelData.h:271
double getDensity(size_t i)
Get density of particle.
Definition modelData.cpp:14
void addF(size_t i, int dof, double f)
Add to force of the node.
Definition modelData.h:398
float d_teFB
Total bond-based fracture energy.
Definition modelData.h:752
const uint8_t & getFix(size_t i) const
Get fixity of the node.
Definition modelData.h:446
const particle::BaseParticle * getParticleFromAllList(size_t i) const
Get pointer to base particle.
Definition modelData.h:76
std::vector< std::vector< size_t > > d_zInfo
Zone information of particles. For zone 0, d_zInfo[0] = [n1, n2] where n1 is the index at which the p...
Definition modelData.h:615
util::Point & getF(size_t i)
Get force of the node.
Definition modelData.h:365
std::vector< float > d_eF
Fracture energy of the nodes.
Definition modelData.h:725
double d_currentDt
Current timestep.
Definition modelData.h:537
std::shared_ptr< inp::ParticleDeck > d_pDeck_p
Particle deck.
Definition modelData.h:561
float d_tw
Total work done.
Definition modelData.h:743
std::unique_ptr< geometry::Fracture > d_fracture_p
Fracture state of bonds.
Definition modelData.h:624
std::vector< util::SymMatrix3 > d_strain
Strain in elements (values at quadrature points)
Definition modelData.h:734
void appendKeyData(std::string key, double data, bool issue_err=false)
Append value to data associated with key.
Definition modelData.h:154
double d_hMin
Maximum mesh over all particles and walls.
Definition modelData.h:574
std::vector< particle::BaseParticle * > d_particlesListTypeParticle
List of particles.
Definition modelData.h:595
std::vector< std::shared_ptr< particle::RefParticle > > d_referenceParticles
Pointer to reference particle.
Definition modelData.h:589
void addF(size_t i, const util::Point &f)
Add to force of the node.
Definition modelData.h:382
const util::Point & getF(size_t i) const
Get force of the node.
Definition modelData.h:368
void setMx(size_t i, const double &mx)
Set weighted-volume (mx) of the node.
Definition modelData.h:487
std::vector< double > d_mX
Weighted volume.
Definition modelData.h:696
std::vector< inp::MatData > d_particlesMatDataList
List of particle material data. Only populated if needed for calculation of stress or other quantitie...
Definition modelData.h:604
std::vector< std::vector< size_t > > d_neighPd
Neighbor data for peridynamic forces.
Definition modelData.h:658
std::vector< particle::BaseParticle * > d_particlesListTypeWall
List of walls.
Definition modelData.h:598
void setX(size_t i, const util::Point &x)
Set current coordinate of the node.
Definition modelData.h:234
std::vector< size_t > d_fPdCompNodes
List of global nodes on which force (peridynamic/internal) is to be computed.
Definition modelData.h:700
float d_teF
Total fracture energy.
Definition modelData.h:749
void addV(size_t i, const util::Point &v)
Add to velocity of the node.
Definition modelData.h:335
float d_te
Total internal energy.
Definition modelData.h:740
A class for nearest neighbor search using nanoflann library.
Definition nsearch.h:178
A class to store particle geometry, nodal discretization, and methods.
nsearch::NFlannSearchKd< 3 > NSearch
Definition modelData.h:30
Collection of methods and data related to particle object.
A structure to represent 3d vectors.
Definition point.h:30