Skip to content

Commit

Permalink
Fixed the various damage model issues
Browse files Browse the repository at this point in the history
  • Loading branch information
reverendbedford committed Aug 30, 2022
1 parent cfa05c7 commit ac7d7b8
Show file tree
Hide file tree
Showing 5 changed files with 30 additions and 12 deletions.
12 changes: 6 additions & 6 deletions src/damage.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -113,7 +113,7 @@ void NEMLScalarDamagedModel_sd::update_sd_actual(
std::copy(s_n, s_n+6, s_prime_n);
for (int i=0; i<6; i++) s_prime_n[i] /= (1-h_n[0]);

base_->update_sd(e_np1, e_n, T_np1, T_n,
base_->update_sd_actual(e_np1, e_n, T_np1, T_n,
t_np1, t_n, s_prime_np1, s_prime_n,
&h_np1[1], &h_n[1],
A_prime_np1, u_np1, u_n, p_np1, p_n);
Expand Down Expand Up @@ -180,15 +180,15 @@ void NEMLScalarDamagedModel_sd::RJ(const double * const x, TrialState * ts,
double s_prime_np1[6];
double s_prime_n[6];
double A_prime_np1[36];
std::vector<double> h_np1_v(base_->nhist());
std::vector<double> h_np1_v(base_->nstate());
double * h_np1 = &h_np1_v[0];
double u_np1;
double p_np1;

std::copy(tss->s_n, tss->s_n+6, s_prime_n);
for (int i=0; i<6; i++) s_prime_n[i] /= (1-tss->w_n);

base_->update_sd(tss->e_np1, tss->e_n, tss->T_np1, tss->T_n,
base_->update_sd_actual(tss->e_np1, tss->e_n, tss->T_np1, tss->T_n,
tss->t_np1, tss->t_n, s_prime_np1, s_prime_n,
h_np1, &tss->h_n[0],
A_prime_np1, u_np1, tss->u_n, p_np1, tss->p_n);
Expand Down Expand Up @@ -238,8 +238,8 @@ void NEMLScalarDamagedModel_sd::make_trial_state(
tss.t_np1 = t_np1;
tss.t_n = t_n;
std::copy(s_n, s_n+6, tss.s_n);
tss.h_n.resize(base_->nhist());
std::copy(h_n+1, h_n+base_->nhist()+1, tss.h_n.begin());
tss.h_n.resize(base_->nstate());
std::copy(h_n+1, h_n+base_->nstate()+1, tss.h_n.begin());
tss.u_n = u_n;
tss.p_n = p_n;
tss.w_n = h_n[0];
Expand Down Expand Up @@ -292,7 +292,7 @@ void NEMLScalarDamagedModel_sd::ekill_update_(double T_np1,
double & u_np1, double u_n,
double & p_np1, double p_n)
{
std::copy(h_n, h_n + nhist(), h_np1);
std::copy(h_n, h_n + nstate(), h_np1);
h_np1[0] = 1.0;
elastic_->C(T_np1, A_np1);
for (int i=0; i<36; i++) A_np1[i] /= sfact_;
Expand Down
9 changes: 7 additions & 2 deletions src/history.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -418,8 +418,13 @@ std::vector<std::string> History::formatted_names() const
std::vector<std::string> names;

for (auto name : order_)
for (size_t i = 0; i < size_of_entry(name); i++)
names.push_back(name + "_" + std::to_string(i));
if (size_of_entry(name) == 1) {
names.push_back(name);
}
else {
for (size_t i = 0; i < size_of_entry(name); i++)
names.push_back(name + "_" + std::to_string(i));
}

return names;
}
Expand Down
1 change: 1 addition & 0 deletions src/models.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -33,6 +33,7 @@ void NEMLModel::populate_hist(History & h) const

void NEMLModel::init_hist(History & h) const
{
h.zero(); // This shuts up valgrind errors.
init_static(h);
init_state(h);
}
Expand Down
16 changes: 14 additions & 2 deletions test/test_damage.py
Original file line number Diff line number Diff line change
Expand Up @@ -674,8 +674,9 @@ def test_is_damage(self):
self.assertTrue(self.model.is_damage_model())

def test_current_damage(self):
print(self.model.report_internal_variable_names())
self.assertEqual(self.model.get_damage(self.hist_n), self.d_n)
hist = self.model.init_store()
hist[6] = self.d_n
self.assertEqual(self.model.get_damage(hist), self.d_n)

def test_kill(self):
self.assertFalse(self.model.should_del_element(self.hist_n))
Expand Down Expand Up @@ -746,6 +747,17 @@ def test_function(self):
f_calcd = self.A * self.effective(self.stress) ** self.a
self.assertTrue(np.isclose(f_model, f_calcd))

def test_is_damage(self):
self.assertTrue(self.model.is_damage_model())

def test_current_damage(self):
hist = self.model.init_store()
hist[6] = self.d_n
self.assertEqual(self.model.get_damage(hist), self.d_n)

def test_kill(self):
self.assertFalse(self.model.should_del_element(self.hist_n))

class TestExponentialDamage(unittest.TestCase, CommonStandardDamageModel,
CommonScalarDamageModel, CommonDamagedModel):
def setUp(self):
Expand Down
4 changes: 2 additions & 2 deletions test/test_neml.py
Original file line number Diff line number Diff line change
Expand Up @@ -285,7 +285,7 @@ def setUp(self):

def test_internal_variable_names(self):
should = self.model.report_internal_variable_names()
correct = ['small_stress_0', 'small_stress_1', 'small_stress_2', 'small_stress_3', 'small_stress_4', 'small_stress_5', 'alpha_0']
correct = ['small_stress_0', 'small_stress_1', 'small_stress_2', 'small_stress_3', 'small_stress_4', 'small_stress_5', 'alpha']
self.assertEqual(should, correct)

def gen_hist(self):
Expand Down Expand Up @@ -498,7 +498,7 @@ def setUp(self):

def test_internal_variable_names(self):
should = self.model.report_internal_variable_names()
correct = ['small_stress_0', 'small_stress_1', 'small_stress_2', 'small_stress_3', 'small_stress_4', 'small_stress_5', 'alpha_0', 'backstress_0_0', 'backstress_0_1', 'backstress_0_2', 'backstress_0_3', 'backstress_0_4', 'backstress_0_5', 'backstress_1_0', 'backstress_1_1', 'backstress_1_2', 'backstress_1_3', 'backstress_1_4', 'backstress_1_5', 'backstress_2_0', 'backstress_2_1', 'backstress_2_2', 'backstress_2_3', 'backstress_2_4', 'backstress_2_5']
correct = ['small_stress_0', 'small_stress_1', 'small_stress_2', 'small_stress_3', 'small_stress_4', 'small_stress_5', 'alpha', 'backstress_0_0', 'backstress_0_1', 'backstress_0_2', 'backstress_0_3', 'backstress_0_4', 'backstress_0_5', 'backstress_1_0', 'backstress_1_1', 'backstress_1_2', 'backstress_1_3', 'backstress_1_4', 'backstress_1_5', 'backstress_2_0', 'backstress_2_1', 'backstress_2_2', 'backstress_2_3', 'backstress_2_4', 'backstress_2_5']
self.assertEqual(should, correct)

def gen_hist(self):
Expand Down

0 comments on commit ac7d7b8

Please sign in to comment.