Skip to content

Commit

Permalink
Need a way to get active orientations back out of models
Browse files Browse the repository at this point in the history
  • Loading branch information
reverendbedford committed Nov 10, 2022
1 parent a94df3d commit 607703d
Show file tree
Hide file tree
Showing 5 changed files with 31 additions and 1 deletion.
3 changes: 3 additions & 0 deletions include/cp/batch.h
Original file line number Diff line number Diff line change
Expand Up @@ -27,6 +27,9 @@ NEML_EXPORT void set_orientation_passive_batch(SingleCrystalModel & model, size_
NEML_EXPORT void get_orientation_passive_batch(SingleCrystalModel & model, size_t n,
double * const hist,
std::vector<Orientation> & orientations);
NEML_EXPORT void get_orientation_active_batch(SingleCrystalModel & model, size_t n,
double * const hist,
std::vector<Orientation> & orientations);

} // namespace neml

Expand Down
1 change: 1 addition & 0 deletions include/cp/polycrystal.h
Original file line number Diff line number Diff line change
Expand Up @@ -33,6 +33,7 @@ class NEML_EXPORT PolycrystalModel: public NEMLModel_ldi
const double * w(const double * const store, size_t i) const;

virtual std::vector<Orientation> orientations(double * const store) const;
virtual std::vector<Orientation> orientations_active(double * const store) const;

protected:
std::shared_ptr<SingleCrystalModel> model_;
Expand Down
12 changes: 12 additions & 0 deletions src/cp/batch.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -70,4 +70,16 @@ void get_orientation_passive_batch(SingleCrystalModel & model, size_t n,
}
}

void get_orientation_active_batch(SingleCrystalModel & model, size_t n,
double * const hist,
std::vector<Orientation> & orientations)
{
orientations.resize(n);
size_t nh = model.nstore();

for (size_t i=0; i<n; i++) {
orientations[i] = model.get_active_orientation(&hist[i*nh]);
}
}

} // namespace neml
9 changes: 9 additions & 0 deletions src/cp/polycrystal.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -111,6 +111,15 @@ std::vector<Orientation> PolycrystalModel::orientations(double * const store) co
return res;
}

std::vector<Orientation> PolycrystalModel::orientations_active(double * const store) const
{
std::vector<Orientation> res;

get_orientation_active_batch(*model_, n(),
history(store, 0), res);
return res;
}

TaylorModel::TaylorModel(ParameterSet & params) :
PolycrystalModel(params)
{
Expand Down
7 changes: 6 additions & 1 deletion src/cp/polycrystal_wrap.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -21,7 +21,12 @@ PYBIND11_MODULE(polycrystal, m) {
[](PolycrystalModel & m, py::array_t<double, py::array::c_style> h) -> std::vector<Orientation>
{
return m.orientations(arr2ptr<double>(h));
}, "Return the current vector of orientations")
}, "Return the current vector of orientations in the passive convention")
.def("orientations_active",
[](PolycrystalModel & m, py::array_t<double, py::array::c_style> h) -> std::vector<Orientation>
{
return m.orientations_active(arr2ptr<double>(h));
}, "Return the current vector of orientations in the active convention")
;

py::class_<TaylorModel, PolycrystalModel, std::shared_ptr<TaylorModel>>(m, "TaylorModel")
Expand Down

0 comments on commit 607703d

Please sign in to comment.