StressStrainControlBoundary Class Reference

A cuboid box consists of periodic boundaries that can be strain/stress controlled and achieve different deformation modes. User needs to define target stress/strainrate matrix, gain_factor and a boolean parameter isStrainRateControlled to True/False to activate/deactivate strainrate control. More...

#include <StressStrainControlBoundary.h>

+ Inheritance diagram for StressStrainControlBoundary:

Public Member Functions

 StressStrainControlBoundary ()
 default constructor. More...
 
 StressStrainControlBoundary (const StressStrainControlBoundary &b)=default
 copy constructor. More...
 
virtual ~StressStrainControlBoundary ()=default
 destructor. More...
 
void read (std::istream &is) override
 Reads the object's id_ from given istream. More...
 
void write (std::ostream &os) const override
 Adds object's id_ to given ostream. More...
 
std::string getName () const override
 Sets the name of the boundary. More...
 
StressStrainControlBoundarycopy () const override
 Used to create a copy of the object. More...
 
void checkBoundaryAfterParticlesMove (ParticleHandler &particleHandler) override
 Virtual function that does things to particles, each timestep after particles have moved. More...
 
void set (const Matrix3D &stressGoal, const Matrix3D &strainRate, const Matrix3D &pGain, bool isStrainRateControlled, const Matrix3D &iGain={0, 0, 0, 0, 0, 0, 0, 0, 0})
 Sets all boundary inputs at once and determines which deformation mode it is, then combine the right boundaries to assemble the cuboid box/domain. More...
 
void setStrainRate (const Matrix3D &strainRate)
 
void createPeriodicParticles (ParticleHandler &particleHandler) override
 Create the periodic particles after read in from a restart file to attain right information. More...
 
void checkPeriodicLeesEdwardsBoundariesAfterParticlesMove (ParticleHandler &particleHandler)
 Call the boundary and update them based on the new domain size after the stress-control movement. More...
 
void determineLengthAndCentre ()
 Determine the length in x,y,z and center location of domain. More...
 
void activateStrainRateControl (const ParticleHandler &particleHandler)
 Activate the strainrate control for particle movement based on user's boolean input. More...
 
void computeStrainRate ()
 Compute the change of strainrate tensor based on the stress difference and then update the tensor. More...
 
void determineStressControlledShearBoundaries ()
 Determines stress-controlled shear Lees-Edwards boundary in x-y direction and normal periodic in z direction. More...
 
void updateDomainSize ()
 Update the domain to new sizes. More...
 
Matrix3D getStrainRate () const
 Accesses the strainrate tensor. More...
 
Matrix3D getStressGoal () const
 Accesses the target stress tensor. More...
 
Matrix3D getGainFactor () const
 Accesses the gainFactor. More...
 
Mdouble computeStressError ()
 
double getIntegratedShift () const
 
const PIControllergetXX () const
 getter for xx More...
 
const PIControllergetXY () const
 getter for xy More...
 
const PIControllergetYY () const
 getter for yy More...
 
const PIControllergetZZ () const
 getter for zz More...
 
- Public Member Functions inherited from BaseBoundary
 BaseBoundary ()
 default constructor. More...
 
 BaseBoundary (const BaseBoundary &b)
 copy constructor More...
 
 ~BaseBoundary () override
 destructor More...
 
virtual void createPeriodicParticle (BaseParticle *p UNUSED, ParticleHandler &pH UNUSED)
 Creates a periodic particle in case of periodic boundaries in serial build. More...
 
virtual void createPeriodicParticles (ParticleHandler &pH UNUSED)
 Creates periodic copies of given particle in case of periodic boundaries. More...
 
virtual void checkBoundaryBeforeTimeStep (DPMBase *md)
 Virtual function that does things before each time step. More...
 
virtual void actionsBeforeTimeLoop ()
 Virtual function that does something after DPMBase::setupInitialConditions but before the first time step. More...
 
virtual void modifyGhostAfterCreation (BaseParticle *particle, int i)
 
virtual void writeVTK (std::fstream &file)
 
void setHandler (BoundaryHandler *handler)
 Sets the boundary's BoundaryHandler. More...
 
BoundaryHandlergetHandler () const
 Returns the boundary's BoundaryHandler. More...
 
- Public Member Functions inherited from BaseObject
 BaseObject ()=default
 Default constructor. More...
 
 BaseObject (const BaseObject &p)=default
 Copy constructor, copies all the objects BaseObject contains. More...
 
virtual ~BaseObject ()=default
 virtual destructor More...
 
virtual void moveInHandler (unsigned int index)
 Except that it is virtual, it does the same thing as setIndex() does. More...
 
void setIndex (unsigned int index)
 Allows one to assign an index to an object in the handler/container. More...
 
void setId (unsigned long id)
 Assigns a unique identifier to each object in the handler (container) which remains constant even after the object is deleted from the container/handler. More...
 
unsigned int getIndex () const
 Returns the index of the object in the handler. More...
 
unsigned int getId () const
 Returns the unique identifier of any particular object. More...
 
void setGroupId (unsigned groupId)
 
unsigned getGroupId () const
 

Private Attributes

Matrix3D stressGoal_
 Stores the stress value the boundary should attain. More...
 
Matrix3D strainRate_
 
Matrix3D pGainFactor_
 
Matrix3D iGainFactor_
 
PIController xx
 PI-Controllers for the strain rate in xx, xy, yy, zz directions. More...
 
PIController xy
 
PIController yy
 
PIController zz
 
bool isStrainRateControlled_
 The boolean input, true means switch on the strain rate control for particles affine movements. More...
 
Mdouble integratedShift_
 Shift integrated for all the time when using Lees-Edwards Boundary. More...
 
std::vector< LeesEdwardsBoundaryleesEdwardsBoundaries_
 
std::vector< PeriodicBoundaryperiodicBoundaries_
 

Detailed Description

A cuboid box consists of periodic boundaries that can be strain/stress controlled and achieve different deformation modes. User needs to define target stress/strainrate matrix, gain_factor and a boolean parameter isStrainRateControlled to True/False to activate/deactivate strainrate control.

Inherits from BaseObject, combining the PeriodicBoundary and LeesEdwardsBoundary. This boundary takes matrices inputs from user and determine which boundaries should be combined together and therefore achieve different deformation mode: e.g. triaxial, bi-aixial, uni-axial and simple shear. The Lees-Edwards boundary is only activated when there is a target shear stress/strainrate set in the matrix. Note that the same component can only be set in either target stress matrix or strainrate matrix, i.e. if the target shear stress XY is set to a constant, the strainrate_XY has to be set to 0, or just leave it free, it is by default set to 0. Note also there are 9 components in the stress/strainrate tensor but we only support the following components due to the boundary inheritance: XX, YY, ZZ, XY

Constructor & Destructor Documentation

◆ StressStrainControlBoundary() [1/2]

StressStrainControlBoundary::StressStrainControlBoundary ( )

default constructor.

constructor, set all the parameters that boundary needs as inputs to zero.

17  : BaseBoundary()
18 {
24  integratedShift_ = 0.0;
25  //
26  logger(DEBUG, "StressStrainControlBoundary::StressStrainControlBoundary() finished");
27 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ DEBUG
BaseBoundary()
default constructor.
Definition: BaseBoundary.cc:11
void setZero()
Sets all elements to zero.
Definition: Matrix.cc:54
bool isStrainRateControlled_
The boolean input, true means switch on the strain rate control for particles affine movements.
Definition: StressStrainControlBoundary.h:136
Mdouble integratedShift_
Shift integrated for all the time when using Lees-Edwards Boundary.
Definition: StressStrainControlBoundary.h:139
Matrix3D pGainFactor_
Definition: StressStrainControlBoundary.h:124
Matrix3D strainRate_
Definition: StressStrainControlBoundary.h:124
Matrix3D stressGoal_
Stores the stress value the boundary should attain.
Definition: StressStrainControlBoundary.h:124
Matrix3D iGainFactor_
Definition: StressStrainControlBoundary.h:130

References DEBUG, iGainFactor_, integratedShift_, isStrainRateControlled_, logger, pGainFactor_, Matrix3D::setZero(), strainRate_, and stressGoal_.

Referenced by copy().

◆ StressStrainControlBoundary() [2/2]

StressStrainControlBoundary::StressStrainControlBoundary ( const StressStrainControlBoundary b)
default

copy constructor.

◆ ~StressStrainControlBoundary()

virtual StressStrainControlBoundary::~StressStrainControlBoundary ( )
virtualdefault

destructor.

Member Function Documentation

◆ activateStrainRateControl()

void StressStrainControlBoundary::activateStrainRateControl ( const ParticleHandler particleHandler)

Activate the strainrate control for particle movement based on user's boolean input.

This function activate the strain rate control based on user inputs, it takes only the strainRate_ tensor and move particles based on this tensor, also move boundaries based on the new domain size.

177 {
178  const DPMBase* const dpm = getHandler()->getDPMBase();
179  // update the domain size based on the new strain rate tensor
181 
182  // Move the boundaries to next time step
183  if (strainRate_.XY != 0 || stressGoal_.XY != 0)
184  {
186  }
187  else
188  {
189  // Move the boundary in x direction to next time step
190  periodicBoundaries_[0].set(Vec3D(1.0, 0.0, 0.0), dpm->getXMin(), dpm->getXMax());
191  // Move the boundary in y direction to next time step
192  periodicBoundaries_[1].set(Vec3D(0.0, 1.0, 0.0), dpm->getYMin(), dpm->getYMax());
193  // Move the boundary in z direction to next time step
194  periodicBoundaries_[2].set(Vec3D(0.0, 0.0, 1.0), dpm->getZMin(), dpm->getZMax());
195  }
196 
197  // Give the strain-rate for all particles and move them to next time step before integration
199  {
200  for (auto& p : particleHandler)
201  {
202  // Box length Lx, Ly and Lz, and center point C of the box
203  Vec3D centerBox;
204  centerBox.X = (dpm->getXMax() + dpm->getXMin()) / 2.0;
205  centerBox.Y = (dpm->getYMax() + dpm->getYMin()) / 2.0;
206  centerBox.Z = (dpm->getZMax() + dpm->getZMin()) / 2.0;
207  Vec3D relativeToCenter;
208  relativeToCenter.X = p->getPosition().X - centerBox.X;
209  relativeToCenter.Y = p->getPosition().Y - centerBox.Y;
210  relativeToCenter.Z = p->getPosition().Z - centerBox.Z;
211  p->move(Vec3D(strainRate_.XX * dpm->getTimeStep() * relativeToCenter.X +
212  strainRate_.XY * dpm->getTimeStep() * relativeToCenter.Y,
213  strainRate_.YY * dpm->getTimeStep() * relativeToCenter.Y,
214  strainRate_.ZZ * dpm->getTimeStep() * relativeToCenter.Z));
215  }
216  }
217 }
float * p
Definition: Tutorial_Map_using.cpp:9
BoundaryHandler * getHandler() const
Returns the boundary's BoundaryHandler.
Definition: BaseBoundary.cc:122
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:733
The DPMBase header includes quite a few header files, defining all the handlers, which are essential....
Definition: DPMBase.h:56
Mdouble getXMin() const
If the length of the problem domain in x-direction is XMax - XMin, then getXMin() returns XMin.
Definition: DPMBase.h:603
Mdouble getXMax() const
If the length of the problem domain in x-direction is XMax - XMin, then getXMax() returns XMax.
Definition: DPMBase.h:610
Mdouble getYMin() const
If the length of the problem domain in y-direction is YMax - YMin, then getYMin() returns YMin.
Definition: DPMBase.h:616
Mdouble getTimeStep() const
Returns the simulation time step.
Definition: DPMBase.cc:1241
Mdouble getYMax() const
If the length of the problem domain in y-direction is YMax - YMin, then getYMax() returns XMax.
Definition: DPMBase.h:622
Mdouble getZMax() const
If the length of the problem domain in z-direction is ZMax - ZMin, then getZMax() returns ZMax.
Definition: DPMBase.h:634
Mdouble getZMin() const
If the length of the problem domain in z-direction is ZMax - ZMin, then getZMin() returns ZMin.
Definition: DPMBase.h:628
Mdouble XY
Definition: Kernel/Math/Matrix.h:22
Mdouble YY
Definition: Kernel/Math/Matrix.h:22
Mdouble ZZ
Definition: Kernel/Math/Matrix.h:22
Mdouble XX
all nine matrix elements
Definition: Kernel/Math/Matrix.h:22
std::vector< PeriodicBoundary > periodicBoundaries_
Definition: StressStrainControlBoundary.h:147
void determineStressControlledShearBoundaries()
Determines stress-controlled shear Lees-Edwards boundary in x-y direction and normal periodic in z di...
Definition: StressStrainControlBoundary.cc:245
void updateDomainSize()
Update the domain to new sizes.
Definition: StressStrainControlBoundary.cc:223
Definition: Kernel/Math/Vector.h:30
Mdouble Y
Definition: Kernel/Math/Vector.h:45
Mdouble Z
Definition: Kernel/Math/Vector.h:45
Mdouble X
the vector components
Definition: Kernel/Math/Vector.h:45

References determineStressControlledShearBoundaries(), BaseHandler< T >::getDPMBase(), BaseBoundary::getHandler(), DPMBase::getTimeStep(), DPMBase::getXMax(), DPMBase::getXMin(), DPMBase::getYMax(), DPMBase::getYMin(), DPMBase::getZMax(), DPMBase::getZMin(), isStrainRateControlled_, p, periodicBoundaries_, strainRate_, stressGoal_, updateDomainSize(), Vec3D::X, Matrix3D::XX, Matrix3D::XY, Vec3D::Y, Matrix3D::YY, Vec3D::Z, and Matrix3D::ZZ.

Referenced by checkBoundaryAfterParticlesMove().

◆ checkBoundaryAfterParticlesMove()

void StressStrainControlBoundary::checkBoundaryAfterParticlesMove ( ParticleHandler particleHandler)
overridevirtual

Virtual function that does things to particles, each timestep after particles have moved.

This is where the stress-strain control is implemented and the boundary will be checked at each time step to make sure the target stress/strain rate are achieved.

Reimplemented from BaseBoundary.

98 {
100 
101  //Real Stress and StrainRate Control every time step
102  logger.assert_debug(getHandler() != nullptr,
103  "you need to set the handler of this boundary before the parameters can be set");
104 
106 
107  //this checks if stressGoal matrix is non-zero and then activate only the stress control
108  if (stressGoal_.XX != 0 || stressGoal_.YY != 0 || stressGoal_.ZZ != 0 ||
109  stressGoal_.XY != 0)
110  {
111  //the strain rate in the corresponding direction will be calculated
113  }
114  activateStrainRateControl(particleHandler);
115 }
void checkPeriodicLeesEdwardsBoundariesAfterParticlesMove(ParticleHandler &particleHandler)
Call the boundary and update them based on the new domain size after the stress-control movement.
Definition: StressStrainControlBoundary.cc:121
void activateStrainRateControl(const ParticleHandler &particleHandler)
Activate the strainrate control for particle movement based on user's boolean input.
Definition: StressStrainControlBoundary.cc:176
void determineLengthAndCentre()
Determine the length in x,y,z and center location of domain.
Definition: StressStrainControlBoundary.cc:138
void computeStrainRate()
Compute the change of strainrate tensor based on the stress difference and then update the tensor.
Definition: StressStrainControlBoundary.cc:149

References activateStrainRateControl(), checkPeriodicLeesEdwardsBoundariesAfterParticlesMove(), computeStrainRate(), determineLengthAndCentre(), BaseBoundary::getHandler(), logger, stressGoal_, Matrix3D::XX, Matrix3D::XY, Matrix3D::YY, and Matrix3D::ZZ.

◆ checkPeriodicLeesEdwardsBoundariesAfterParticlesMove()

void StressStrainControlBoundary::checkPeriodicLeesEdwardsBoundariesAfterParticlesMove ( ParticleHandler particleHandler)

Call the boundary and update them based on the new domain size after the stress-control movement.

This function checks the boundaries after the particles are moved, to serve the stress-control for the following steps.

122 {
123  // Call Boundaries
125  {
126  b.checkBoundaryAfterParticlesMove(particleHandler);
127  }
129  {
130  b.checkBoundaryAfterParticlesMove(particleHandler);
131  }
132 }
Scalar * b
Definition: benchVecAdd.cpp:17
Class which creates a boundary with Lees-Edwards type periodic boundary conditions.
Definition: LeesEdwardsBoundary.h:26
Defines a pair of periodic walls. Inherits from BaseBoundary.
Definition: PeriodicBoundary.h:20
std::vector< LeesEdwardsBoundary > leesEdwardsBoundaries_
Definition: StressStrainControlBoundary.h:146

References b, leesEdwardsBoundaries_, and periodicBoundaries_.

Referenced by checkBoundaryAfterParticlesMove().

◆ computeStrainRate()

void StressStrainControlBoundary::computeStrainRate ( )

Compute the change of strainrate tensor based on the stress difference and then update the tensor.

This function is used to compute the new strain rate tensor based on the stress differences between the actual and user set target values. This is based on a PI-controller, strainRate = pGain * (stress-stressGoal) + iGain * integral(stress-stressGoal)

150 {
151  // calculate the stress total and average over the volume
152  Matrix3D stress = getHandler()->getDPMBase()->getTotalStress();
153 
154  // controller, sets the amount by which the strain rate tensor has to be changed
155  double timeStep = getHandler()->getDPMBase()->getTimeStep();
156  // amount by which the strain rate tensor has to be changed
157  if (stressGoal_.XX != 0) {
158  strainRate_.XX = xx.apply(stress.XX - stressGoal_.XX,timeStep);
159  }
160  if (stressGoal_.YY != 0) {
161  strainRate_.YY = yy.apply(stress.YY - stressGoal_.YY,timeStep);
162  }
163  if (stressGoal_.ZZ != 0) {
164  strainRate_.ZZ = zz.apply(stress.ZZ - stressGoal_.ZZ,timeStep);
165  }
166  if (stressGoal_.XY != 0) {
167  strainRate_.XY = xy.apply(stress.XY - stressGoal_.XY,timeStep);
168  }
169 }
Matrix3D getTotalStress() const
Calculate the total stress tensor in the system averaged over the whole volume.
Definition: DPMBase.cc:5510
Implementation of a 3D matrix.
Definition: Kernel/Math/Matrix.h:17
Mdouble apply(Mdouble error, Mdouble timeStep)
Definition: PIController.cc:40
PIController xy
Definition: StressStrainControlBoundary.h:133
PIController xx
PI-Controllers for the strain rate in xx, xy, yy, zz directions.
Definition: StressStrainControlBoundary.h:133
PIController zz
Definition: StressStrainControlBoundary.h:133
PIController yy
Definition: StressStrainControlBoundary.h:133

References PIController::apply(), BaseHandler< T >::getDPMBase(), BaseBoundary::getHandler(), DPMBase::getTimeStep(), DPMBase::getTotalStress(), strainRate_, stressGoal_, xx, Matrix3D::XX, xy, Matrix3D::XY, yy, Matrix3D::YY, zz, and Matrix3D::ZZ.

Referenced by checkBoundaryAfterParticlesMove().

◆ computeStressError()

Mdouble StressStrainControlBoundary::computeStressError ( )

◆ copy()

StressStrainControlBoundary * StressStrainControlBoundary::copy ( ) const
overridevirtual

Used to create a copy of the object.

Copy method; creates a copy on the boundary and returns its pointer.

Implements BaseBoundary.

36 {
37  return new StressStrainControlBoundary(*this);
38 }
StressStrainControlBoundary()
default constructor.
Definition: StressStrainControlBoundary.cc:16

References StressStrainControlBoundary().

◆ createPeriodicParticles()

void StressStrainControlBoundary::createPeriodicParticles ( ParticleHandler particleHandler)
override

Create the periodic particles after read in from a restart file to attain right information.

This function is used to create ghost particles such that the stress are calculated correctly after restart from another configuration.

375 {
376  // Call Boundaries
378  {
379  b.createPeriodicParticles(particleHandler);
380  }
382  {
383  b.createPeriodicParticles(particleHandler);
384  }
385 }

References b, leesEdwardsBoundaries_, and periodicBoundaries_.

◆ determineLengthAndCentre()

void StressStrainControlBoundary::determineLengthAndCentre ( )

Determine the length in x,y,z and center location of domain.

This function determines both the length in x,y,z direction and center location of the domain.

139 {
140 // const DPMBase* const dpm = getHandler()->getDPMBase();
141 }

Referenced by checkBoundaryAfterParticlesMove().

◆ determineStressControlledShearBoundaries()

void StressStrainControlBoundary::determineStressControlledShearBoundaries ( )

Determines stress-controlled shear Lees-Edwards boundary in x-y direction and normal periodic in z direction.

This function determines the boundary for stress-controlled shear situation. In this case, the Lees-Edwards boundary in x-y directions and normal periodic boundary in z direction are combined together as a cuboid shear box.

246 {
247  const DPMBase* const dpm = getHandler()->getDPMBase();
248  //Determine the current shear velocity and shift of the Lees-Edwards boundary
249  const Mdouble velocity_xy = strainRate_.XY * (dpm->getXMax() - dpm->getXMin());
250  integratedShift_ += velocity_xy * dpm->getTimeStep();
251 
252  // Determine how to move boundaries based on strain rate control or boundary control
254  {
255  // Move the Lees-Edwards boundary in z direction to next time step
256  leesEdwardsBoundaries_[0].set(
257  [this](Mdouble time)
258  { return integratedShift_; },
259  [](Mdouble time UNUSED)
260  { return 0; },
261  dpm->getXMin(), dpm->getXMax(), dpm->getYMin(), dpm->getYMax());
262  }
263  else
264  {
265  // Update the velocity of Lees-Edwards boundary in x-y direction to next time step
266  leesEdwardsBoundaries_[0].set(
267  [this](Mdouble time)
268  { return integratedShift_; },
269  [velocity_xy](Mdouble time UNUSED)
270  { return velocity_xy; },
271  dpm->getXMin(), dpm->getXMax(), dpm->getYMin(), dpm->getYMax());
272  }
273  // Move the boundary in z direction to next time step
274  periodicBoundaries_[0].set(Vec3D(0.0, 0.0, 1.0), dpm->getZMin(), dpm->getZMax());
275 }
#define UNUSED
Definition: GeneralDefine.h:18

References BaseHandler< T >::getDPMBase(), BaseBoundary::getHandler(), DPMBase::getTimeStep(), DPMBase::getXMax(), DPMBase::getXMin(), DPMBase::getYMax(), DPMBase::getYMin(), DPMBase::getZMax(), DPMBase::getZMin(), integratedShift_, isStrainRateControlled_, leesEdwardsBoundaries_, periodicBoundaries_, strainRate_, UNUSED, and Matrix3D::XY.

Referenced by activateStrainRateControl().

◆ getGainFactor()

Matrix3D StressStrainControlBoundary::getGainFactor ( ) const
inline

Accesses the gainFactor.

99 {return pGainFactor_;}

References pGainFactor_.

◆ getIntegratedShift()

double StressStrainControlBoundary::getIntegratedShift ( ) const
inline
106 {return integratedShift_;}

References integratedShift_.

◆ getName()

std::string StressStrainControlBoundary::getName ( ) const
overridevirtual

Sets the name of the boundary.

Returns the name of the object class.

Returns
the object class' name.

Implements BaseObject.

88 {
89  return "StressStrainControlBoundary";
90 }

◆ getStrainRate()

Matrix3D StressStrainControlBoundary::getStrainRate ( ) const
inline

Accesses the strainrate tensor.

93 {return strainRate_;}

References strainRate_.

Referenced by ShearStage::actionsAfterTimeStep(), ShearStage::printTime(), and ShearStage::ShearStage().

◆ getStressGoal()

Matrix3D StressStrainControlBoundary::getStressGoal ( ) const
inline

Accesses the target stress tensor.

96 {return stressGoal_;}

References stressGoal_.

◆ getXX()

const PIController& StressStrainControlBoundary::getXX ( ) const
inline

getter for xx

109 { return xx; }

References xx.

◆ getXY()

const PIController& StressStrainControlBoundary::getXY ( ) const
inline

getter for xy

111 { return xy; }

References xy.

◆ getYY()

const PIController& StressStrainControlBoundary::getYY ( ) const
inline

getter for yy

113 { return yy; }

References yy.

◆ getZZ()

const PIController& StressStrainControlBoundary::getZZ ( ) const
inline

getter for zz

115 { return zz; }

References zz.

◆ read()

void StressStrainControlBoundary::read ( std::istream &  is)
overridevirtual

Reads the object's id_ from given istream.

Reads the boundary properties from an istream.

Parameters
[in]isthe istream.

Implements BaseBoundary.

71 {
73  std::string dummy;
74  is >> dummy >> stressGoal_;
75  is >> dummy >> strainRate_;
76  is >> dummy >> pGainFactor_;
77  is >> dummy >> iGainFactor_;
78  is >> dummy >> isStrainRateControlled_;
80  is >> dummy >> integratedShift_;
81 }
void read(std::istream &is) override=0
Reads the object's id_ from given istream NB: purely virtual function, overriding the version of Base...
Definition: BaseBoundary.cc:40
void set(const Matrix3D &stressGoal, const Matrix3D &strainRate, const Matrix3D &pGain, bool isStrainRateControlled, const Matrix3D &iGain={0, 0, 0, 0, 0, 0, 0, 0, 0})
Sets all boundary inputs at once and determines which deformation mode it is, then combine the right ...
Definition: StressStrainControlBoundary.cc:288
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286

References iGainFactor_, integratedShift_, isStrainRateControlled_, pGainFactor_, BaseBoundary::read(), set(), strainRate_, stressGoal_, and oomph::Global_string_for_annotation::string().

◆ set()

void StressStrainControlBoundary::set ( const Matrix3D stressGoal,
const Matrix3D strainRate,
const Matrix3D pGain,
bool  isStrainRateControlled,
const Matrix3D iGain = {0, 0, 0, 0, 0, 0, 0, 0, 0} 
)

Sets all boundary inputs at once and determines which deformation mode it is, then combine the right boundaries to assemble the cuboid box/domain.

This function sets the inputs for the whole StressStrainControlBoundary based on the user inputs.

Parameters
[in]stressGoalThe target stress tensor that needs to achieve at the end of the deformation.
[in]strainRateThe target strain rate tensor, cannot be set non-zero with target Stress, i.e. stressGoal.XX != 0, then strainRate.XX = 0.
[in]pGainThe incremental factor for the stress control, usually a very small value ~ 0.0001.
[in]isStrainRateControlledThe boolean key to determine whether particles are moved by strain rate affine movement (true) or dragged by boundary itself (false).
[in]iGainThe incremental factor for the integral stress control
289 {
290  periodicBoundaries_.clear();
291  leesEdwardsBoundaries_.clear();
292 
293  isStrainRateControlled_ = isStrainRateControlled;
294  stressGoal_ = stressGoal;
295  strainRate_ = strainRate;
296  pGainFactor_ = pGain;
297  iGainFactor_ = iGain;
298 
299  logger.assert_always(stressGoal.XZ == 0, "Shear stress in XZ cannot be controlled; use shear stress in XY instead");
300  logger.assert_always(stressGoal.ZX == 0, "Shear stress in ZX cannot be controlled; use shear stress in XY instead");
301  logger.assert_always(stressGoal.YZ == 0, "Shear stress in YZ cannot be controlled; use shear stress in XY instead");
302  logger.assert_always(stressGoal.ZY == 0, "Shear stress in ZY cannot be controlled; use shear stress in XY instead");
303  logger.assert_always(strainRate.XZ == 0, "Strain rate in XZ cannot be controlled; use strain rate in XY instead");
304  logger.assert_always(strainRate.ZX == 0, "Strain rate in ZX cannot be controlled; use strain rate in XY instead");
305  logger.assert_always(strainRate.YZ == 0, "Strain rate in YZ cannot be controlled; use strain rate in XY instead");
306  logger.assert_always(strainRate.ZY == 0, "Strain rate in ZY cannot be controlled; use strain rate in XY instead");
307  logger.assert_always(stressGoal.XZ != 0 ? (pGain.XZ != 0 || iGain.XZ != 0) : true,
308  "You need to set a gain factor in XZ in order to control stress");
309  logger.assert_always(stressGoal.XX != 0 ? (pGain.XX != 0 || iGain.XX != 0) : true,
310  "You need to set a gain factor in XX in order to control stress");
311  logger.assert_always(stressGoal.YY != 0 ? (pGain.YY != 0 || iGain.YY != 0) : true,
312  "You need to set a gain factor in YY in order to control stress");
313  logger.assert_always(stressGoal.ZZ != 0 ? (pGain.ZZ != 0 || iGain.ZZ != 0) : true,
314  "You need to set a gain factor in ZZ in order to control stress");
315 
316  logger.assert_always(getHandler() != nullptr, "You need to set the handler of this boundary");
317  const DPMBase *dpm = getHandler()->getDPMBase();
318 
319  // set gains for the strain-rate controllers
324  // reset the integrated error to zero
325  xx.reset();
326  xy.reset();
327  yy.reset();
328  zz.reset();
329 
330  // if there is no shear (compression only)
331  if (stressGoal_.XY == 0 && strainRate_.XY == 0)
332  {
333  logger(INFO, "Shear rate is zero, setting up three periodic boundaries");
334  // Set up new box periodic boundaries, no simple shear activated
335  PeriodicBoundary boundary;
336  boundary.setHandler(getHandler());
337  boundary.set(Vec3D(1.0, 0.0, 0.0), dpm->getXMin(), dpm->getXMax());
338  periodicBoundaries_.push_back(boundary);
339  boundary.set(Vec3D(0.0, 1.0, 0.0), dpm->getYMin(), dpm->getYMax());
340  periodicBoundaries_.push_back(boundary);
341  boundary.set(Vec3D(0.0, 0.0, 1.0), dpm->getZMin(), dpm->getZMax());
342  periodicBoundaries_.push_back(boundary);
343  }
344  else
345  {
346  logger(INFO, "Shear rate is not zero, setting up Lees-Edwards in xy and a periodic boundary in z");
347 
348  PeriodicBoundary boundary;
349  boundary.setHandler(getHandler());
350  boundary.set(Vec3D(0.0, 0.0, 1.0), dpm->getZMin(), dpm->getZMax());
351  periodicBoundaries_.push_back(boundary);
352 
353  // Lees Edwards bc in y direction & periodic boundary in x direction, simple shear boundary activated
354  LeesEdwardsBoundary leesEdwardsBoundary;
355  leesEdwardsBoundary.setHandler(getHandler());
356  double integratedShift = integratedShift_;
357  leesEdwardsBoundary.set(
358  [&integratedShift](Mdouble time UNUSED) { return integratedShift; },
359  [](Mdouble time UNUSED) { return 0; },
360  dpm->getXMin(), dpm->getXMax(), dpm->getYMin(), dpm->getYMax());
361  leesEdwardsBoundaries_.push_back(leesEdwardsBoundary);
362  }
363 }
@ INFO
void setHandler(BoundaryHandler *handler)
Sets the boundary's BoundaryHandler.
Definition: BaseBoundary.cc:113
void set(std::function< Mdouble(Mdouble)> shift, std::function< Mdouble(Mdouble)> velocity, Mdouble left, Mdouble right, Mdouble down, Mdouble up)
Sets all boundary properties.
Definition: LeesEdwardsBoundary.cc:29
Mdouble ZX
Definition: Kernel/Math/Matrix.h:22
Mdouble ZY
Definition: Kernel/Math/Matrix.h:22
Mdouble YZ
Definition: Kernel/Math/Matrix.h:22
Mdouble XZ
Definition: Kernel/Math/Matrix.h:22
void set(Mdouble pGain, Mdouble iGain)
Definition: PIController.cc:11
void reset()
Definition: PIController.cc:46
void set(Vec3D normal, Mdouble distanceLeft, Mdouble distanceRight)
Defines a PeriodicBoundary by its normal and positions.
Definition: PeriodicBoundary.cc:63

References BaseHandler< T >::getDPMBase(), BaseBoundary::getHandler(), DPMBase::getXMax(), DPMBase::getXMin(), DPMBase::getYMax(), DPMBase::getYMin(), DPMBase::getZMax(), DPMBase::getZMin(), iGainFactor_, INFO, integratedShift_, isStrainRateControlled_, leesEdwardsBoundaries_, logger, periodicBoundaries_, pGainFactor_, PIController::reset(), PIControllerBasic::set(), LeesEdwardsBoundary::set(), PeriodicBoundary::set(), BaseBoundary::setHandler(), strainRate_, stressGoal_, UNUSED, xx, Matrix3D::XX, xy, Matrix3D::XY, Matrix3D::XZ, yy, Matrix3D::YY, Matrix3D::YZ, Matrix3D::ZX, Matrix3D::ZY, zz, and Matrix3D::ZZ.

Referenced by ShearStage::actionsAfterTimeStep(), read(), StressStrainControl::setupInitialConditions(), and ShearStage::ShearStage().

◆ setStrainRate()

void StressStrainControlBoundary::setStrainRate ( const Matrix3D strainRate)
366 {
367  strainRate_ = strainRate;
368 }

References strainRate_.

◆ updateDomainSize()

void StressStrainControlBoundary::updateDomainSize ( )

Update the domain to new sizes.

This function is used to update the domain size based on the strainRate tensor, note that the system is symmetric and therefore we have to update boundary in both Min and Max.

224 {
225  DPMBase* const dpm = getHandler()->getDPMBase();
226  // Box length Lx, Ly and Lz for next time step
227  Vec3D lengthBox;
228  lengthBox.X = (dpm->getXMax() - dpm->getXMin());
229  lengthBox.Y = (dpm->getYMax() - dpm->getYMin());
230  lengthBox.Z = (dpm->getZMax() - dpm->getZMin());
231  // Change the system size according to next time step
232  dpm->setXMax(dpm->getXMax() + 0.5 * lengthBox.X * strainRate_.XX * dpm->getTimeStep());
233  dpm->setXMin(dpm->getXMin() - 0.5 * lengthBox.X * strainRate_.XX * dpm->getTimeStep());
234  dpm->setYMax(dpm->getYMax() + 0.5 * lengthBox.Y * strainRate_.YY * dpm->getTimeStep());
235  dpm->setYMin(dpm->getYMin() - 0.5 * lengthBox.Y * strainRate_.YY * dpm->getTimeStep());
236  dpm->setZMax(dpm->getZMax() + 0.5 * lengthBox.Z * strainRate_.ZZ * dpm->getTimeStep());
237  dpm->setZMin(dpm->getZMin() - 0.5 * lengthBox.Z * strainRate_.ZZ * dpm->getTimeStep());
238 }
void setYMin(Mdouble newYMin)
Sets the value of YMin, the lower bound of the problem domain in the y-direction.
Definition: DPMBase.cc:1025
void setYMax(Mdouble newYMax)
Sets the value of YMax, the upper bound of the problem domain in the y-direction.
Definition: DPMBase.cc:1182
void setZMin(Mdouble newZMin)
Sets the value of ZMin, the lower bound of the problem domain in the z-direction.
Definition: DPMBase.cc:1049
void setXMax(Mdouble newXMax)
Sets the value of XMax, the upper bound of the problem domain in the x-direction.
Definition: DPMBase.cc:1156
void setZMax(Mdouble newZMax)
Sets the value of ZMax, the upper bound of the problem domain in the z-direction.
Definition: DPMBase.cc:1208
void setXMin(Mdouble newXMin)
Sets the value of XMin, the lower bound of the problem domain in the x-direction.
Definition: DPMBase.cc:1001

References BaseHandler< T >::getDPMBase(), BaseBoundary::getHandler(), DPMBase::getTimeStep(), DPMBase::getXMax(), DPMBase::getXMin(), DPMBase::getYMax(), DPMBase::getYMin(), DPMBase::getZMax(), DPMBase::getZMin(), DPMBase::setXMax(), DPMBase::setXMin(), DPMBase::setYMax(), DPMBase::setYMin(), DPMBase::setZMax(), DPMBase::setZMin(), strainRate_, Vec3D::X, Matrix3D::XX, Vec3D::Y, Matrix3D::YY, Vec3D::Z, and Matrix3D::ZZ.

Referenced by activateStrainRateControl().

◆ write()

void StressStrainControlBoundary::write ( std::ostream &  os) const
overridevirtual

Adds object's id_ to given ostream.

Writes boundary's properties to an ostream.

Parameters
[in]osthe ostream.

Implements BaseBoundary.

45 {
47  os << " stressGoal " << stressGoal_;
48  os << " strainRate " << strainRate_;
49  os << " pGainFactor " << pGainFactor_;
50  os << " iGainFactor " << iGainFactor_;
51  os << " isStrainRateControlled " << isStrainRateControlled_;
52  os << " integratedShift " << integratedShift_;
53 // os << leesEdwardsBoundaries_.size() << ' ';
54 // for (const LeesEdwardsBoundary& b : leesEdwardsBoundaries_)
55 // {
56 // os << b << ' ';
57 // }
58 //
59 // os << periodicBoundaries_.size() << ' ';
60 // for (const PeriodicBoundary& b : periodicBoundaries_)
61 // {
62 // os << b << ' ';
63 // }
64 }
void write(std::ostream &os) const override=0
Adds object's id_ to given ostream NB: purely virtual function, overriding the version of BaseObject.
Definition: BaseBoundary.cc:49

References iGainFactor_, integratedShift_, isStrainRateControlled_, pGainFactor_, strainRate_, stressGoal_, and BaseBoundary::write().

Member Data Documentation

◆ iGainFactor_

Matrix3D StressStrainControlBoundary::iGainFactor_
private

◆ integratedShift_

Mdouble StressStrainControlBoundary::integratedShift_
private

Shift integrated for all the time when using Lees-Edwards Boundary.

Referenced by determineStressControlledShearBoundaries(), getIntegratedShift(), read(), set(), StressStrainControlBoundary(), and write().

◆ isStrainRateControlled_

bool StressStrainControlBoundary::isStrainRateControlled_
private

The boolean input, true means switch on the strain rate control for particles affine movements.

Referenced by activateStrainRateControl(), determineStressControlledShearBoundaries(), read(), set(), StressStrainControlBoundary(), and write().

◆ leesEdwardsBoundaries_

std::vector<LeesEdwardsBoundary> StressStrainControlBoundary::leesEdwardsBoundaries_
private

Store boundaries into a vector for the pushback. Note, there is always either 0 LeesEdwardsBoundary (XY) and 3 PeriodicBoundary (XYZ), or 1 LeesEdwardsBoundary (XY) and 1 PeriodicBoundary (Z).

Referenced by checkPeriodicLeesEdwardsBoundariesAfterParticlesMove(), createPeriodicParticles(), determineStressControlledShearBoundaries(), and set().

◆ periodicBoundaries_

◆ pGainFactor_

Matrix3D StressStrainControlBoundary::pGainFactor_
private

◆ strainRate_

◆ stressGoal_

Matrix3D StressStrainControlBoundary::stressGoal_
private

Stores the stress value the boundary should attain.

Unused if the all stressGoal values are set to zero.

Referenced by activateStrainRateControl(), checkBoundaryAfterParticlesMove(), computeStrainRate(), computeStressError(), getStressGoal(), read(), set(), StressStrainControlBoundary(), and write().

◆ xx

PIController StressStrainControlBoundary::xx
private

PI-Controllers for the strain rate in xx, xy, yy, zz directions.

Referenced by computeStrainRate(), getXX(), and set().

◆ xy

PIController StressStrainControlBoundary::xy
private

Referenced by computeStrainRate(), getXY(), and set().

◆ yy

PIController StressStrainControlBoundary::yy
private

Referenced by computeStrainRate(), getYY(), and set().

◆ zz

PIController StressStrainControlBoundary::zz
private

Referenced by computeStrainRate(), getZZ(), and set().


The documentation for this class was generated from the following files: