Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add a general lattice class #143

Merged
merged 1 commit into from
Sep 14, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
27 changes: 27 additions & 0 deletions doc/sphinx/cp/crystallography/GeneralLattice.rst
Original file line number Diff line number Diff line change
@@ -0,0 +1,27 @@
GeneralLattice
==============

Overview
--------

A general :doc:`Lattice` class where the user can provide the lattice
directions and the symmetry class.

Parameters
----------

.. csv-table::
:header: "Parameter", "Object type", "Description", "Default"
:widths: 12, 30, 50, 8

``a1``, :code:`std::vector<double>`, First lattice vector, No
``a2``, :code:`std::vector<double>`, Second lattice vector, No
``a3``, :code:`std::vector<double>`, Third lattice vector, No
``symmetry_group``, :cpp:class:`neml::SymmetryGroup`, Symmetry group, No
``slip_systems``, :code:`list_systems`, Initial list of slip systems, ``{}``
``twin_systems``, :code:`list_systems`, Initial list of twin systems, ``{}``

Class description
-----------------

.. doxygenclass:: neml::GeneralLattice
3 changes: 2 additions & 1 deletion doc/sphinx/cp/crystallography/Lattice.rst
Original file line number Diff line number Diff line change
Expand Up @@ -24,7 +24,8 @@ the lattice vectors and symmetry group.

.. toctree::
:maxdepth: 1


GeneralLattice
CubicLattice
HCPLattice

Expand Down
18 changes: 18 additions & 0 deletions include/cp/crystallography.h
Original file line number Diff line number Diff line change
Expand Up @@ -182,6 +182,24 @@ class NEML_EXPORT Lattice {
std::vector<std::vector<size_t>> normal_map_;
};

/// General lattice with no specialization
class NEML_EXPORT GeneralLattice: public NEMLObject, public Lattice {
public:
GeneralLattice(ParameterSet & params);

/// String type for the object system
static std::string type();
/// Initialize from parameter set
static std::unique_ptr<NEMLObject> initialize(ParameterSet & params);
/// Default parameters
static ParameterSet parameters();

/// Override serialization to account for dynamic changes
virtual ParameterSet & current_parameters();
};

static Register<GeneralLattice> regGeneralLattice;

class NEML_EXPORT CubicLattice: public NEMLObject, public Lattice {
public:
/// Specialized Lattice for cubic systems, initialize with the lattice
Expand Down
46 changes: 46 additions & 0 deletions src/cp/crystallography.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -609,6 +609,52 @@ void Lattice::update_normals_(const std::vector<Vector> & new_planes)
normal_map_.push_back(new_indices);
}

GeneralLattice::GeneralLattice(ParameterSet & params) :
NEMLObject(params),
Lattice(
Vector(params.get_parameter<std::vector<double>>("a1")),
Vector(params.get_parameter<std::vector<double>>("a2")),
Vector(params.get_parameter<std::vector<double>>("a3")),
params.get_object_parameter<SymmetryGroup>("symmetry_group"),
params.get_parameter<list_systems>("slip_systems"),
params.get_parameter<twin_systems>("twin_systems"))
{

}

std::string GeneralLattice::type()
{
return "GeneralLattice";
}

ParameterSet GeneralLattice::parameters()
{
ParameterSet pset(GeneralLattice::type());

pset.add_parameter<std::vector<double>>("a1");
pset.add_parameter<std::vector<double>>("a2");
pset.add_parameter<std::vector<double>>("a3");
pset.add_parameter<NEMLObject>("symmetry_group");
pset.add_optional_parameter<list_systems>("slip_systems",
list_systems());
pset.add_optional_parameter<list_systems>("twin_systems",
twin_systems());

return pset;
}

ParameterSet & GeneralLattice::current_parameters()
{
current_params_.assign_parameter("slip_systems", current_slip_);
current_params_.assign_parameter("twin_systems", current_twin_);
return current_params_;
}

std::unique_ptr<NEMLObject> GeneralLattice::initialize(ParameterSet & params)
{
return neml::make_unique<GeneralLattice>(params);
}

CubicLattice::CubicLattice(ParameterSet & params) :
NEMLObject(params),
Lattice(
Expand Down
8 changes: 8 additions & 0 deletions src/cp/crystallography_wrap.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -80,6 +80,14 @@ PYBIND11_MODULE(crystallography, m) {
.def("plane_index", &Lattice::plane_index)
;

py::class_<GeneralLattice, Lattice, NEMLObject, std::shared_ptr<GeneralLattice>>(m, "GeneralLattice")
.def(py::init([](py::args args, py::kwargs kwargs)
{
return create_object_python<GeneralLattice>(
args, kwargs, {"a1", "a2", "a3", "symmetry_group"});
}))
;

py::class_<CubicLattice, Lattice, NEMLObject, std::shared_ptr<CubicLattice>>(m, "CubicLattice")
.def(py::init([](py::args args, py::kwargs kwargs)
{
Expand Down
25 changes: 25 additions & 0 deletions test/test_crystallography.py
Original file line number Diff line number Diff line change
Expand Up @@ -301,6 +301,31 @@ def test_slip(self):
self.assertAlmostEqual(n.norm(), 1.0)
self.assertAlmostEqual(d.dot(n), 0.0)

class TestGeneralLattice(unittest.TestCase, LTests, ShearTests, SlipTest):
def setUp(self):
self.a = 1.3
self.lattice = crystallography.GeneralLattice(
[self.a,0,0],
[0,self.a,0],
[0,0,self.a],
crystallography.SymmetryGroup("432"))

self.S = np.array([[20.0,-15.0,12.0],[-15.0,-40.0,5.0],[12.0,5.0,60.0]])
self.ST = tensors.Symmetric(self.S)
self.Q = rotations.Orientation(30.0,43.0,10.0, angle_type = "degrees")
self.QM = self.Q.to_matrix()

self.lattice.add_slip_system([1,1,0],[1,1,1])

self.correct_numbers = [12,]

def test_planes(self):
self.assertEqual(self.lattice.nplanes, 4)
for i in range(12):
from_normals = self.lattice.unique_planes[self.lattice.plane_index(0, i)]
from_direct = self.lattice.slip_planes[0][i]
self.assertAlmostEqual(np.abs(from_normals.dot(from_direct)), 1.0)

class TestCubicFCC(unittest.TestCase, LTests, ShearTests, SlipTest):
def setUp(self):
self.a = 1.3
Expand Down