MDStressLab++
Loading...
Searching...
No Matches
Mls Class Reference

Computes deformation gradients and stress tensors using Moving Least Squares (MLS) interpolation on a spatial grid. More...

#include <Mls.h>

Collaboration diagram for Mls:

Public Member Functions

 Mls (const BoxConfiguration &body, const Grid< Reference > *pgrid, double radiusMls, const std::string name)
 Constructs the deformation gradient field from a body configuration and grid.
 
 ~Mls ()
 
void pushToCauchy (const std::vector< Matrix3d > &piolaStress, std::vector< Matrix3d > &cauchyStress)
 Converts Piola-Kirchhoff stress to Cauchy stress at grid points.
 
void writeDeformationGradient ()
 Writes the deformation gradient at each grid point to a file.
 
void writeGridPushed ()
 Writes the pushed (deformed) grid point coordinates to a file.
 
void writePushedCauchyStress (std::vector< Matrix3d > &cauchyStress)
 Writes Cauchy stresses at grid points to a file.
 

Public Attributes

std::string name
 Identifier name used for output file prefixes.
 
double radiusMls
 Radius of influence for the MLS weighting function.
 
std::vector< Matrix3ddeformationGradient
 Deformation gradient at each grid point (one per grid coordinate).
 
std::vector< Vector3dgridPushed
 Positions of grid points after applying MLS displacements.
 

Detailed Description

Computes deformation gradients and stress tensors using Moving Least Squares (MLS) interpolation on a spatial grid.

This class implements a 3D MLS approach to compute the local deformation gradient and map Piola-Kirchhoff stress tensors to Cauchy stresses on a grid. It accounts for periodic boundary conditions and constructs local neighbor lists for accuracy.

Typical use: compute deformations and push stresses from atomistic data onto grid points.

Examples
testMls.cpp.

Definition at line 26 of file Mls.h.

Constructor & Destructor Documentation

◆ Mls()

Mls::Mls ( const BoxConfiguration body,
const Grid< Reference > *  pgrid,
double  radiusMls,
const std::string  name 
)

Constructs the deformation gradient field from a body configuration and grid.

Performs MLS interpolation to compute deformation gradients at each grid point using displacement fields derived from a reference and current configuration.

Handles periodic boundary conditions, builds neighbor lists, checks sanity of interpolation points, and stores the final deformation gradient and displaced grid points.

Parameters
bodyThe configuration containing reference and current atom positions.
pgridPointer to the spatial grid (reference coordinates).
radiusMlsThe cutoff radius for MLS interpolation.
nameAn identifier used in output filenames.

Definition at line 21 of file Mls.cpp.

◆ ~Mls()

Mls::~Mls ( )

Definition at line 791 of file Mls.cpp.

Member Function Documentation

◆ pushToCauchy()

void Mls::pushToCauchy ( const std::vector< Matrix3d > &  piolaStress,
std::vector< Matrix3d > &  cauchyStress 
)

Converts Piola-Kirchhoff stress to Cauchy stress at grid points.

Uses the computed deformation gradients and their determinants to push stress tensors from the reference configuration to the current spatial grid.

Parameters
piolaStressInput Piola-Kirchhoff stresses at reference points.
cauchyStressOutput Cauchy stresses on the deformed grid.
Examples
testMls.cpp.

Definition at line 796 of file Mls.cpp.

◆ writeDeformationGradient()

void Mls::writeDeformationGradient ( )

Writes the deformation gradient at each grid point to a file.

File name is constructed as <name>.DeformationGradient. Includes both grid positions and corresponding 3x3 deformation gradient matrices.

Examples
testMls.cpp.

Definition at line 832 of file Mls.cpp.

◆ writeGridPushed()

void Mls::writeGridPushed ( )

Writes the pushed (deformed) grid point coordinates to a file.

File name is <name>.pushedGrids. One line per grid point with x/y/z coordinates.

Examples
testMls.cpp.

Definition at line 872 of file Mls.cpp.

◆ writePushedCauchyStress()

void Mls::writePushedCauchyStress ( std::vector< Matrix3d > &  cauchyStress)

Writes Cauchy stresses at grid points to a file.

File name is <name>.stress. Outputs 6 unique components per symmetric 3x3 stress tensor.

Parameters
cauchyStressVector of Cauchy stress tensors (must match grid size).
Examples
testMls.cpp.

Definition at line 887 of file Mls.cpp.

Member Data Documentation

◆ deformationGradient

std::vector<Matrix3d> Mls::deformationGradient

Deformation gradient at each grid point (one per grid coordinate).

Definition at line 42 of file Mls.h.

◆ gridPushed

std::vector<Vector3d> Mls::gridPushed

Positions of grid points after applying MLS displacements.

Definition at line 47 of file Mls.h.

◆ name

std::string Mls::name

Identifier name used for output file prefixes.

Definition at line 32 of file Mls.h.

◆ radiusMls

double Mls::radiusMls

Radius of influence for the MLS weighting function.

Definition at line 37 of file Mls.h.


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