![]() |
NLMech
0.1.0
|
Processes simulation results and computes post-processing quantities. More...
#include <compute.h>

Public Member Functions | |
| Compute (const std::string &filename) | |
| Constructor. More... | |
| void | finalize () |
| Finalize. | |
| void | init () |
| Initialize the class based on input data. | |
| void | initWriter (rw::writer::Writer *writer, std::vector< util::Point3 > *u) |
| Create new writer object if it is not already done so. More... | |
| void | readComputeInstruction (const std::string &set, tools::pp::InstructionData *data) |
| Read compute instruction. More... | |
| void | readCrackTipData (const std::string &filename, int crack_id, std::vector< tools::pp::CrackTipData > *data, bool is_inclined=false) |
| Read crack tip data. More... | |
Post-processing methods | |
| void | transformU (rw::writer::Writer *writer) |
| Transform displacement. More... | |
| void | transformV (rw::writer::Writer *writer) |
| Transform velocity. More... | |
| void | computeStrain (rw::writer::Writer *writer) |
| Compute strain and stress. More... | |
| void | computeDamage (rw::writer::Writer *writer, std::vector< double > *Z, bool perf_out=false) |
| Compute damage at nodes. More... | |
| void | findCrackTip (std::vector< double > *Z, rw::writer::Writer *writer) |
| Find crack tip location and velocity. More... | |
| void | computeJIntegral () |
| Compute J integral. More... | |
| void | computeJIntegralAngledCrack () |
| Compute J integral where crack line is angled as opposed to either horizontal or verticle line. | |
Private Member Functions | |
Utility methods | |
| void | addUniqueToList (size_t i, std::vector< size_t > *list) |
| Inserts element to list if not present. More... | |
| size_t | findNode (const util::Point3 &x, const std::vector< util::Point3 > *nodes, const std::vector< util::Point3 > *u=nullptr) |
| Find node closest to given point. More... | |
| void | listElemsAndNodesInDomain (const std::pair< util::Point3, util::Point3 > &cd, const double &tol, const double &tol_elem, std::vector< size_t > *nodes, std::vector< size_t > *elements, bool calc_in_ref) |
| Find node within contour and elements intersecting contour. More... | |
| void | decomposeSearchNodes (const std::pair< util::Point3, util::Point3 > &cd, std::vector< size_t > *nodes, std::vector< size_t > *nodes_new) |
| Decompose node list into two lists. More... | |
| bool | triCheckAndInterpolateUV (const util::Point3 &p, util::Point3 &up, util::Point3 &vp, const std::vector< size_t > &ids, bool calc_in_ref, bool check_only=false) |
| Interpolates displacement and velocity at given point in triangle. More... | |
| void | interpolateUV (const util::Point3 &p, util::Point3 &up, util::Point3 &vp, const std::vector< size_t > *nodes, const std::vector< size_t > *elements, bool calc_in_ref) |
| Interpolates displacement and velocity at given point. More... | |
| size_t | interpolateUVNodes (const util::Point3 &p, util::Point3 &up, util::Point3 &vp, const std::vector< size_t > *nodes, bool calc_in_ref) |
| void | getContourContribJInt (const util::Point3 &p, const std::vector< size_t > *nodes, const std::vector< size_t > *elements, const util::Point3 &normal, double &pd_energy, double &kinetic_energy, double &elastic_energy, util::Point3 &dot_u, bool calc_in_ref, const tools::pp::CrackTipData &ctip) |
| Computes contribution to energy into crack from the quadrature point on contour. More... | |
| void | updateCrack (const double &time, const std::vector< double > *Z) |
| Updates crack tip location and crack velocity. More... | |
| void | crackOutput () |
| Performs output of crack tip data. | |
| void | getRectsAndNodesForCrackTip (inp::EdgeCrack &crack, std::vector< std::pair< util::Point3, util::Point3 >> &rects_t, std::vector< std::pair< util::Point3, util::Point3 >> &rects_b, std::vector< std::vector< size_t >> &nodes_t, std::vector< std::vector< size_t >> &nodes_b, std::vector< std::vector< double >> &Z_t, std::vector< std::vector< double >> &Z_b, const std::vector< double > *Z) |
| Creates sequence of rectangles for crack tip search, and finds the ids of nodes in each rectangle and damage associated to the node. More... | |
| void | addNewCrackTip (inp::EdgeCrack &crack, util::Point3 pnew, double time, bool is_top) |
| Adds new tip to crack, computes velocity and crack length. More... | |
| util::Point3 | findTipInRects (inp::EdgeCrack &crack, const double &max_Z, const std::vector< std::pair< util::Point3, util::Point3 >> &rects, const std::vector< std::vector< size_t >> &nodes, const std::vector< std::vector< double >> &Zs, const std::vector< double > *Z, bool is_top) |
| Finds new crack tip location from the rectangle data. More... | |
| void | safeExit (const std::string &err_message, rw::writer::Writer *writer=nullptr) |
| Exits code safely in the event there is error. More... | |
Private Attributes | |
| std::vector< tools::pp::InstructionData > | d_computeData |
| List of compute sets. | |
| tools::pp::InstructionData * | d_currentData |
| Current active compute data. | |
| data::DataManager * | d_dataManager_p |
| Pointer to the Data Mananger. | |
| size_t | d_dtEnd |
| Global output start step. Default is maximum output step for the simulation, i.e. d_dtN. | |
| size_t | d_dtN |
| Total number of output files to process. | |
| size_t | d_dtOutChange |
| Time step at which output interval changes. More... | |
| size_t | d_dtStart |
| Global output start step. Default is 1. | |
| std::string | d_fileZ |
| Filename (if any) to read damage at nodes. | |
| std::string | d_fnErrMsg |
| bool | d_fnSuccess |
| Save success/failur of function. | |
| inp::FractureDeck * | d_fractureDeck_p |
| Fracture deck. | |
| std::string | d_inpFilename |
| Input filename for compute instruction. | |
| inp::Input * | d_input_p |
| Pointer to Fracture object. More... | |
| inp::MaterialDeck * | d_matDeck_p |
| Material deck. | |
| material::pd::BaseMaterial * | d_material_p |
| Pointer to Material object. | |
| std::vector< double > | d_mX |
| Weighted volume. More... | |
| size_t | d_nC |
| Current active compute set. | |
| bool | d_needDamageZ |
| Need damage at nodes. | |
| bool | d_needNeighborList |
| Need neighbor list. | |
| int | d_nOut |
| Current active simulation output step. | |
| std::string | d_outFilename |
| Current active output filename to save postprocessed data. | |
| std::string | d_outPath |
| Output path where postprocessed files are created. | |
| std::string | d_outPreTag |
| Initial tag for each post-processing output. | |
| inp::OutputDeck * | d_outputDeck_p |
| Model deck. More... | |
| std::string | d_simInpFilename |
| Simulation input filename. More... | |
| std::string | d_simOutFilename |
| Mesh filename to get displacement and velocity (results of simulation) | |
| std::string | d_tagZ |
| Tag of damage function Z in vtu file. | |
| std::vector< double > | d_thetaX |
| Dilation. More... | |
| double | d_time |
| Current time. | |
| std::vector< util::Point3 > | d_u |
| Displacement of nodes. | |
| bool | d_uPlus |
| Specify if we consider u(n+1) or u(n) | |
| std::vector< util::Point3 > | d_v |
| Velocity of nodes. | |
| bool | d_writerReady |
| State of writer class. | |
| std::vector< double > | d_Z |
| Damage of nodes. | |
Processes simulation results and computes post-processing quantities.
Currently we have implemented
|
explicit |
Constructor.
| filename | Name of input (yaml) file |

|
private |
Adds new tip to crack, computes velocity and crack length.
| crack | Crack data which will be updated |
| pnew | New crack tip |
| time | Time of update |
| is_top | True if this is top (right) tip otherwise false |
|
private |
Inserts element to list if not present.
| i | Element to be inserted |
| list | Pointer to list |
| void tools::pp::Compute::computeDamage | ( | rw::writer::Writer * | writer, |
| std::vector< double > * | Z, | ||
| bool | perf_out = false |
||
| ) |
Compute damage at nodes.
We compute the damage function
where
is the ball centered at
,
is the bond-strain between points
, and
is the critical strain.
| writer | Pointer to vtk writer |
| Z | Pointer to nodal damage |
| perf_out | Flag to perform output of damage data |
Referenced by Compute().

| void tools::pp::Compute::computeJIntegral | ( | ) |
Compute J integral.
We consider rectangle near the crack tip and compute J-integral using the contour of rectangle. Schematic for horizontal crack (similar for vertical crack)
D C
+ + + + + + + + + + +
+ +
+ +
------------------+ +
+ +
+ +
+ + + + + + + + + + +
A B
Note: Contour is formed by lines A-B, B-C, C-D, D-A
Let contour is denoted as
, where
indicates contour moves with crack tip. Let domain inside contour is defined as
. Let the outward normal to the domain
is
, crack velocity direction is
, and the crack velocity is
. Since
is the crack velocity direction, we have crack velocity
, where
is the magnitude of crack velocity.
We compute following seven kinds of energy associated to crack tip:
is the energy density at point
given by
is the pairwise energy density. For regularized bond based model (see material::pd::RNPBond), it is given by
is the magnitude of crack velocity,
is the crack velocity.
is velocity of material point
.
is the elasticity tensor.
is the internal pairwise strain energy density as described above.
is the derivative of displacement along the crack velocity direction.
is the critical energy release rate.
is the damage at node x and
is the peridynamic energy density at point x.Method:
To compute
we loop over nodes which are in complement of domain
, i.e.
. We compute the contribution of each node in
using pdForceWork.
as we require gradient of displacement along the normal to contour.
as we require gradient of displacement at nodes.Referenced by Compute().


| void tools::pp::Compute::computeStrain | ( | rw::writer::Writer * | writer | ) |
Compute strain and stress.
We perform following post-processing calculation
and its shape functions
, the displacement at quadrature point is given by
is displacement of node
. Linear strain is given by
. Otherwise,
| writer | Pointer to vtk writer |
Referenced by Compute().

|
private |
Decompose node list into two lists.
We remove those nodes from the list which are outside the specified rectangle and add them to the new list.
| cd | Rectangle defined by two corner points |
| nodes | Pointer to ids of nodes |
| nodes_new | Pointer to new list of ids of nodes |

| void tools::pp::Compute::findCrackTip | ( | std::vector< double > * | Z, |
| rw::writer::Writer * | writer | ||
| ) |
Find crack tip location and velocity.
To compute the crack tip location at current output time, we first compute the damage at nodes. We then search for node which has smallest damage among nodes with damage above 1.
To compute crack tip velocity, we compare the tip at two different times . If not specified, then the two different times are the previous output time and current output time.
Since we compute damage in this method, we store the damage in global variable which can then be used in method computeDamage().
| Z | Pointer to nodal damage |
| writer | Pointer to vtk writer |
Referenced by Compute().

|
private |
Find node closest to given point.
If u is not null then given point is assumed to be in current configuration and therefore closest node in current configuration is searched.
| x | Point in reference/current configuration |
| nodes | Pointer to list of nodes |
| u | Pointer to nodal displacements |

|
private |
Finds new crack tip location from the rectangle data.
In this method we implement Step 2 to 5 of Compute::updateCrack.
| crack | Crack data |
| max_Z | Maximum value of damage |
| rects | Sequence of rectangles for search of new tip |
| nodes | Vector of nodal ids in rectangles |
| Zs | Vector of nodal damages in rectangles |
| Z | Pointer to nodal damage data |
| is_top | True if searching for new tip for top side of crack |

|
private |
Computes contribution to energy into crack from the quadrature point on contour.
Computes
at given point, see computeJIntegral() for more details.
| p | Point |
| nodes | Pointer to ids of nodes to perform search |
| elements | Pointer to ids of elements to perform search |
| normal | Normal to the edge of contour |
| pd_energy | Peridynamic energy density |
| kinetic_energy | Kinetic energy density |
| elastic_energy | Elastic energy density |
| dot_u | Material velocity at the quadrature point |
| calc_in_ref | Calculate in reference configuration |
| ctip | Current crack tip data |

|
private |
Creates sequence of rectangles for crack tip search, and finds the ids of nodes in each rectangle and damage associated to the node.
We filter out the nodes which have damage above max threshold and below min threshold defined in tools::pp::FindCrackTip.
| crack | Crack data |
| rects_t | Sequence of rectangles for search of new top (right) tip |
| rects_b | Sequence of rectangles for search of new bottom (left) tip |
| nodes_t | Vector of nodal ids in rectangles for top (right) tip |
| nodes_b | Vector of nodal ids in rectangles for bottom (left) tip |
| Z_t | Vector of nodal damages in rectangles for top (right) tip |
| Z_b | Vector of nodal damages in rectangles for bottom (left) tip |
| Z | Pointer to nodal damage data |

| void tools::pp::Compute::initWriter | ( | rw::writer::Writer * | writer, |
| std::vector< util::Point3 > * | u | ||
| ) |
Create new writer object if it is not already done so.
| writer | Pointer to vtk writer |
| u | Pointer to nodal displacement vector |
Referenced by Compute().


|
private |
Interpolates displacement and velocity at given point.
are vertices of quadrangle, then we consider triangle
and
and find if point belongs to
or
and proceed similar to the case of triangle element (see triCheckAndInterpolateUV).| p | Point at which we want to interpolate |
| up | Displacement at the point p |
| vp | Velocity at the point p |
| nodes | Pointer to ids of nodes to perform search |
| elements | Pointer to ids of elements to perform search |
| calc_in_ref | Calculate in reference configuration |

|
private |
Find node within contour and elements intersecting contour.
For nodes list, we consider domain which envelopes contour. Thickness of domain is typically taken as size of horizon + 2 * mesh size.
For element list, we look for element which has at least 1 node in the domain enveloping contour. The thickness of domain is typically taken as mesh size.
| cd | Rectangle defining the contour |
| tol | Thickness for node list search |
| tol_elem | Thickness for element list search |
| nodes | Pointer to ids of nodes |
| elements | Pointer to ids of elements |
| calc_in_ref | Calculate in reference configuration |

| void tools::pp::Compute::readComputeInstruction | ( | const std::string & | set, |
| tools::pp::InstructionData * | data | ||
| ) |
Read compute instruction.
| set | Tag for given compute set |
| data | Pointer to instruction data |
| void tools::pp::Compute::readCrackTipData | ( | const std::string & | filename, |
| int | crack_id, | ||
| std::vector< tools::pp::CrackTipData > * | data, | ||
| bool | is_inclined = false |
||
| ) |
Read crack tip data.
| filename | Filename to read crack tip data |
| crack_id | Id of crack to read |
| data | Pointer to data where crack tip data will be stored |
|
private |
Exits code safely in the event there is error.
| err_message | Error message |
| writer | Writer object |

| void tools::pp::Compute::transformU | ( | rw::writer::Writer * | writer | ) |
Transform displacement.
Currently following transformation is implemented
| writer | Pointer to vtk writer |
Referenced by Compute().


| void tools::pp::Compute::transformV | ( | rw::writer::Writer * | writer | ) |
Transform velocity.
There are two transformation implemented
| writer | Pointer to vtk writer |
Referenced by Compute().


|
private |
Interpolates displacement and velocity at given point in triangle.
Triangle is specified by global ids of nodes. We first check if point is within the triangle
. If yes then we use TriElem::getShapes to compute the shape functions
at the given point. Using shape function we can interpolate displacement and velocity as follows:
where
are displacement and velocity of vertex
of triangle.
| p | Point at which we want to interpolate |
| up | Displacement at the point p |
| vp | Velocity at the point p |
| ids | Global ids of vertices of triangle |
| calc_in_ref | Calculate in reference configuration |
| check_only | True if only interested in whether the point belongs to triangle |

|
private |
Updates crack tip location and crack velocity.
Method:
Step 2 - 5 is implemented in Compute::findTipInRects.
| time | Current time |
| Z | Vector of damage at nodes |

|
private |
Time step at which output interval changes.
Since we now support criteria based output intervals, we need to know when the transition from one interval to another happens. See model::FDModel::checkOutputCriteria() for details about the criteria based output intervals.
Default value is total number of simulation steps. For this value, we always read output files which correspond to time steps in the interval Dt, where Dt is maximum of two intervals in the simulation input file.
Referenced by Compute().
|
private |
|
private |
Weighted volume.
In case of Rob's state based model, this data is not required
|
private |
|
private |
Simulation input filename.
This filename is same input file used in running simulation
|
private |
Dilation.
In case of Rob's state based model, this will give the spherical (hydrostatic) strain