tests.test_tensors

PyDRex: Tests for tensor operations.

  1"""> PyDRex: Tests for tensor operations."""
  2
  3import numpy as np
  4from pydrex import tensors as _tensors
  5from pydrex.minerals import StiffnessTensors
  6
  7
  8def test_voigt_decompose():
  9    """Test decomposition of Voigt 6x6 matrix into distinct contractions."""
 10    olivine_tensor = _tensors.voigt_to_elastic_tensor(StiffnessTensors().olivine)
 11    dilat, voigt = _tensors.voigt_decompose(StiffnessTensors().olivine)
 12    np.testing.assert_array_equal(np.einsum("ijkk", olivine_tensor), dilat)
 13    np.testing.assert_array_equal(np.einsum("ijkj", olivine_tensor), voigt)
 14
 15
 16def test_voigt_tensor():
 17    """Test elasticity tensor <-> 6x6 Voigt matrix conversions."""
 18    olivine_tensor = np.array(
 19        [
 20            [
 21                [[320.71, 0.0, 0.0], [0.0, 69.84, 0.0], [0.0, 0.0, 71.22]],
 22                [[0.0, 78.36, 0.0], [78.36, 0.0, 0.0], [0.0, 0.0, 0.0]],
 23                [[0.0, 0.0, 77.67], [0.0, 0.0, 0.0], [77.67, 0.0, 0.0]],
 24            ],
 25            [
 26                [[0.0, 78.36, 0.0], [78.36, 0.0, 0.0], [0.0, 0.0, 0.0]],
 27                [[69.84, 0.0, 0.0], [0.0, 197.25, 0.0], [0.0, 0.0, 74.8]],
 28                [[0.0, 0.0, 0.0], [0.0, 0.0, 63.77], [0.0, 63.77, 0.0]],
 29            ],
 30            [
 31                [[0.0, 0.0, 77.67], [0.0, 0.0, 0.0], [77.67, 0.0, 0.0]],
 32                [[0.0, 0.0, 0.0], [0.0, 0.0, 63.77], [0.0, 63.77, 0.0]],
 33                [[71.22, 0.0, 0.0], [0.0, 74.8, 0.0], [0.0, 0.0, 234.32]],
 34            ],
 35        ]
 36    )
 37    np.testing.assert_array_equal(
 38        _tensors.voigt_to_elastic_tensor(StiffnessTensors().olivine),
 39        olivine_tensor,
 40    )
 41    np.testing.assert_array_equal(
 42        _tensors.elastic_tensor_to_voigt(
 43            _tensors.voigt_to_elastic_tensor(StiffnessTensors().olivine)
 44        ),
 45        StiffnessTensors().olivine,
 46    )
 47    np.testing.assert_array_equal(
 48        _tensors.voigt_to_elastic_tensor(
 49            _tensors.elastic_tensor_to_voigt(olivine_tensor)
 50        ),
 51        olivine_tensor,
 52    )
 53    np.testing.assert_array_equal(
 54        _tensors.elastic_tensor_to_voigt(
 55            _tensors.voigt_to_elastic_tensor(StiffnessTensors().enstatite)
 56        ),
 57        StiffnessTensors().enstatite,
 58    )
 59
 60
 61def test_voigt_to_vector():
 62    """Test Voigt vector construction."""
 63    np.testing.assert_allclose(
 64        _tensors.voigt_matrix_to_vector(StiffnessTensors().enstatite),
 65        np.array(
 66            [
 67                236.9,
 68                180.5,
 69                230.4,
 70                80.32733034,
 71                89.37829714,
 72                112.57139956,
 73                168.6,
 74                158.8,
 75                160.2,
 76                0.0,
 77                0.0,
 78                0.0,
 79                0.0,
 80                0.0,
 81                0.0,
 82                0.0,
 83                0.0,
 84                0.0,
 85                0.0,
 86                0.0,
 87                0.0,
 88            ]
 89        ),
 90        atol=1e-9,
 91    )
 92    np.testing.assert_array_equal(
 93        StiffnessTensors().olivine,
 94        _tensors.voigt_vector_to_matrix(
 95            _tensors.voigt_matrix_to_vector(StiffnessTensors().olivine),
 96        ),
 97    )
 98    r = np.random.rand(6, 6)
 99    np.testing.assert_array_equal(
100        _tensors.voigt_matrix_to_vector(r),
101        np.array(
102            [
103                r[0, 0],
104                r[1, 1],
105                r[2, 2],
106                np.sqrt(2) * r[1, 2],
107                np.sqrt(2) * r[2, 0],
108                np.sqrt(2) * r[0, 1],
109                2 * r[3, 3],
110                2 * r[4, 4],
111                2 * r[5, 5],
112                2 * r[0, 3],
113                2 * r[1, 4],
114                2 * r[2, 5],
115                2 * r[2, 3],
116                2 * r[0, 4],
117                2 * r[1, 5],
118                2 * r[1, 3],
119                2 * r[2, 4],
120                2 * r[0, 5],
121                2 * np.sqrt(2) * r[4, 5],
122                2 * np.sqrt(2) * r[5, 3],
123                2 * np.sqrt(2) * r[3, 4],
124            ]
125        ),
126    )
def test_voigt_decompose():
 9def test_voigt_decompose():
10    """Test decomposition of Voigt 6x6 matrix into distinct contractions."""
11    olivine_tensor = _tensors.voigt_to_elastic_tensor(StiffnessTensors().olivine)
12    dilat, voigt = _tensors.voigt_decompose(StiffnessTensors().olivine)
13    np.testing.assert_array_equal(np.einsum("ijkk", olivine_tensor), dilat)
14    np.testing.assert_array_equal(np.einsum("ijkj", olivine_tensor), voigt)

Test decomposition of Voigt 6x6 matrix into distinct contractions.

def test_voigt_tensor():
17def test_voigt_tensor():
18    """Test elasticity tensor <-> 6x6 Voigt matrix conversions."""
19    olivine_tensor = np.array(
20        [
21            [
22                [[320.71, 0.0, 0.0], [0.0, 69.84, 0.0], [0.0, 0.0, 71.22]],
23                [[0.0, 78.36, 0.0], [78.36, 0.0, 0.0], [0.0, 0.0, 0.0]],
24                [[0.0, 0.0, 77.67], [0.0, 0.0, 0.0], [77.67, 0.0, 0.0]],
25            ],
26            [
27                [[0.0, 78.36, 0.0], [78.36, 0.0, 0.0], [0.0, 0.0, 0.0]],
28                [[69.84, 0.0, 0.0], [0.0, 197.25, 0.0], [0.0, 0.0, 74.8]],
29                [[0.0, 0.0, 0.0], [0.0, 0.0, 63.77], [0.0, 63.77, 0.0]],
30            ],
31            [
32                [[0.0, 0.0, 77.67], [0.0, 0.0, 0.0], [77.67, 0.0, 0.0]],
33                [[0.0, 0.0, 0.0], [0.0, 0.0, 63.77], [0.0, 63.77, 0.0]],
34                [[71.22, 0.0, 0.0], [0.0, 74.8, 0.0], [0.0, 0.0, 234.32]],
35            ],
36        ]
37    )
38    np.testing.assert_array_equal(
39        _tensors.voigt_to_elastic_tensor(StiffnessTensors().olivine),
40        olivine_tensor,
41    )
42    np.testing.assert_array_equal(
43        _tensors.elastic_tensor_to_voigt(
44            _tensors.voigt_to_elastic_tensor(StiffnessTensors().olivine)
45        ),
46        StiffnessTensors().olivine,
47    )
48    np.testing.assert_array_equal(
49        _tensors.voigt_to_elastic_tensor(
50            _tensors.elastic_tensor_to_voigt(olivine_tensor)
51        ),
52        olivine_tensor,
53    )
54    np.testing.assert_array_equal(
55        _tensors.elastic_tensor_to_voigt(
56            _tensors.voigt_to_elastic_tensor(StiffnessTensors().enstatite)
57        ),
58        StiffnessTensors().enstatite,
59    )

Test elasticity tensor <-> 6x6 Voigt matrix conversions.

def test_voigt_to_vector():
 62def test_voigt_to_vector():
 63    """Test Voigt vector construction."""
 64    np.testing.assert_allclose(
 65        _tensors.voigt_matrix_to_vector(StiffnessTensors().enstatite),
 66        np.array(
 67            [
 68                236.9,
 69                180.5,
 70                230.4,
 71                80.32733034,
 72                89.37829714,
 73                112.57139956,
 74                168.6,
 75                158.8,
 76                160.2,
 77                0.0,
 78                0.0,
 79                0.0,
 80                0.0,
 81                0.0,
 82                0.0,
 83                0.0,
 84                0.0,
 85                0.0,
 86                0.0,
 87                0.0,
 88                0.0,
 89            ]
 90        ),
 91        atol=1e-9,
 92    )
 93    np.testing.assert_array_equal(
 94        StiffnessTensors().olivine,
 95        _tensors.voigt_vector_to_matrix(
 96            _tensors.voigt_matrix_to_vector(StiffnessTensors().olivine),
 97        ),
 98    )
 99    r = np.random.rand(6, 6)
100    np.testing.assert_array_equal(
101        _tensors.voigt_matrix_to_vector(r),
102        np.array(
103            [
104                r[0, 0],
105                r[1, 1],
106                r[2, 2],
107                np.sqrt(2) * r[1, 2],
108                np.sqrt(2) * r[2, 0],
109                np.sqrt(2) * r[0, 1],
110                2 * r[3, 3],
111                2 * r[4, 4],
112                2 * r[5, 5],
113                2 * r[0, 3],
114                2 * r[1, 4],
115                2 * r[2, 5],
116                2 * r[2, 3],
117                2 * r[0, 4],
118                2 * r[1, 5],
119                2 * r[1, 3],
120                2 * r[2, 4],
121                2 * r[0, 5],
122                2 * np.sqrt(2) * r[4, 5],
123                2 * np.sqrt(2) * r[5, 3],
124                2 * np.sqrt(2) * r[3, 4],
125            ]
126        ),
127    )

Test Voigt vector construction.