#include "smd_material_models.h"
#include <cmath>
#include <cstdlib>
#include <utility>
#include <iostream>
#include <cstdio>
#include "math_special.h"
#include <Eigen/Eigen>
using namespace LAMMPS_NS::MathSpecial;
using namespace std;
using namespace Eigen;
#define MIN(A,B) ((A) < (B) ? (A) : (B))
#define MAX(A,B) ((A) > (B) ? (A) : (B))
void LinearEOS(double lambda, double pInitial, double d, double dt, double &pFinal, double &p_rate) {
p_rate = lambda * d;
pFinal = pInitial + dt * p_rate;
}
void ShockEOS(double rho, double rho0, double e, double e0, double c0, double S, double Gamma, double pInitial, double dt,
double &pFinal, double &p_rate) {
double mu = rho / rho0 - 1.0;
double pH = rho0 * square(c0) * mu * (1.0 + mu) / square(1.0 - (S - 1.0) * mu);
pFinal = (-pH + rho * Gamma * (e - e0));
p_rate = (pFinal - pInitial) / dt;
}
void polynomialEOS(double rho, double rho0, double , double C0, double C1, double C2, double C3, double , double , double ,
double pInitial, double dt, double &pFinal, double &p_rate) {
double mu = rho / rho0 - 1.0;
if (mu > 0.0) {
pFinal = C0 + C1 * mu + C2 * mu * mu + C3 * mu * mu * mu; } else {
pFinal = C0 + C1 * mu + C3 * mu * mu * mu; }
pFinal = -pFinal;
p_rate = (pFinal - pInitial) / dt;
}
void TaitEOS_density(const double exponent, const double c0_reference, const double rho_reference, const double rho_current,
double &pressure, double &sound_speed) {
double B = rho_reference * c0_reference * c0_reference / exponent;
double tmp = pow(rho_current / rho_reference, exponent);
pressure = B * (tmp - 1.0);
double bulk_modulus = B * tmp * exponent; sound_speed = sqrt(bulk_modulus / rho_current);
}
void PerfectGasEOS(const double gamma, const double vol, const double mass, const double energy, double &pFinal, double &c0) {
if (energy > 0.0) {
pFinal = (1.0 - gamma) * energy / vol;
c0 = sqrt((gamma - 1.0) * energy / mass);
} else {
pFinal = c0 = 0.0;
}
}
void LinearStrength(const double mu, const Matrix3d sigmaInitial_dev, const Matrix3d d_dev, const double dt,
Matrix3d &sigmaFinal_dev__, Matrix3d &sigma_dev_rate__) {
sigma_dev_rate__ = 2.0 * mu * d_dev;
sigmaFinal_dev__ = sigmaInitial_dev + dt * sigma_dev_rate__;
}
void LinearPlasticStrength(const double G, const double yieldStress, const Matrix3d sigmaInitial_dev, const Matrix3d d_dev,
const double dt, Matrix3d &sigmaFinal_dev__, Matrix3d &sigma_dev_rate__, double &plastic_strain_increment) {
Matrix3d sigmaTrial_dev, dev_rate;
double J2;
dev_rate = 2.0 * G * d_dev;
sigmaTrial_dev = sigmaInitial_dev + dt * dev_rate;
J2 = sqrt(3. / 2.) * sigmaTrial_dev.norm();
if (J2 < yieldStress) {
sigma_dev_rate__ = dev_rate;
sigmaFinal_dev__ = sigmaTrial_dev;
plastic_strain_increment = 0.0;
} else {
plastic_strain_increment = (J2 - yieldStress) / (3.0 * G);
sigmaFinal_dev__ = (yieldStress / J2) * sigmaTrial_dev;
sigma_dev_rate__ = sigmaFinal_dev__ - sigmaInitial_dev;
}
}
void JohnsonCookStrength(const double G, const double cp, const double espec, const double A, const double B, const double a,
const double C, const double epdot0, const double T0, const double Tmelt, const double , const double dt, const double ep,
const double epdot, const Matrix3d sigmaInitial_dev, const Matrix3d d_dev, Matrix3d &sigmaFinal_dev__,
Matrix3d &sigma_dev_rate__, double &plastic_strain_increment) {
Matrix3d sigmaTrial_dev, dev_rate;
double J2, yieldStress;
double deltaT = espec / cp;
double TH = deltaT / (Tmelt - T0);
TH = MAX(TH, 0.0);
double epdot_ratio = epdot / epdot0;
epdot_ratio = MAX(epdot_ratio, 1.0);
yieldStress = (A + B * pow(ep, a)) * (1.0 + C * log(epdot_ratio));
dev_rate = 2.0 * G * d_dev;
sigmaTrial_dev = sigmaInitial_dev + dt * dev_rate;
J2 = sqrt(3. / 2.) * sigmaTrial_dev.norm();
if (J2 < yieldStress) {
sigma_dev_rate__ = dev_rate;
sigmaFinal_dev__ = sigmaTrial_dev;
plastic_strain_increment = 0.0;
} else {
plastic_strain_increment = (J2 - yieldStress) / (3.0 * G);
sigmaFinal_dev__ = (yieldStress / J2) * sigmaTrial_dev;
sigma_dev_rate__ = sigmaFinal_dev__ - sigmaInitial_dev;
}
}
bool IsotropicMaxStrainDamage(const Matrix3d E, const double maxStrain) {
SelfAdjointEigenSolver < Matrix3d > es;
es.compute(E);
double max_eigenvalue = es.eigenvalues().maxCoeff();
if (max_eigenvalue > maxStrain) {
return true;
} else {
return false;
}
}
bool IsotropicMaxStressDamage(const Matrix3d S, const double maxStress) {
SelfAdjointEigenSolver < Matrix3d > es;
es.compute(S);
double max_eigenvalue = es.eigenvalues().maxCoeff();
if (max_eigenvalue > maxStress) {
return true;
} else {
return false;
}
}
double JohnsonCookFailureStrain(const double p, const Matrix3d Sdev, const double d1, const double d2, const double d3,
const double d4, const double epdot0, const double epdot) {
double vm = sqrt(3. / 2.) * Sdev.norm(); if (vm < 0.0) {
cout << "this is sdev " << endl << Sdev << endl;
printf("vm=%f < 0.0, surely must be an error\n", vm);
exit(1);
}
double triax = p / (vm + 0.01 * fabs(p)); if (triax < 0.0) {
triax = 0.0;
} else if (triax > 3.0) {
triax = 3.0;
}
double jc_failure_strain = d1 + d2 * exp(d3 * triax);
if (d4 > 0.0) { if (epdot > epdot0) {
double epdot_ratio = epdot / epdot0;
jc_failure_strain *= (1.0 + d4 * log(epdot_ratio));
}
}
return jc_failure_strain;
}