3.1. anisotropy.effective_modelling
The anisotropy.effective_modelling module provides a suite of functions
to model the effective anisotropic properties of composite materials. This
includes:
Hudson’s model for a material containing penny-shaped cracks
Tandon and Weng’s model for a material containing ellipsoidal inclusions
Silver and Savage’s model for N-layer modelling
Backus’ model for a material comprised of alternating layers of elastic material
- copyright
2021–2022, AnisotroPy developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
3.1.1. anisotropy.effective_modelling.backus
Modelling the effective elastic properties of a material consisting of repeated alternating layers of isotropic materials.
For more information, please read:
Backus, G.E., 1962. Long‐wave elastic anisotropy produced by horizontal layering. Journal of Geophysical Research, 67(11), pp.4427-4440.
If you use this code, please cite this article.
- copyright
2021–2022, AnisotroPy developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- anisotropy.effective_modelling.backus.model(material1, material2, relative_fraction)[source]
Calculate the averaged elastic parameters as described in Backus (1962). The variable naming system used is consistent with the elastic coefficients in the paper.
L = C44 = 1/⟨1/μ⟩ M = C66 = ⟨μ⟩ R = 1 / C33 = ⟨θ/μ⟩ S = ⟨θ * μ⟩ T = ⟨θ⟩
These map into the entries of the effective elastic matrix, C’ij, as:
C’11 = 4M - 4S + (1-2T)^2/R C’12 = 2M - 4S + (1-2T)^2/R C’13 = (1-2T)/R C’33 = 1/R C’44 = L C’66 = M = (C’11 - C’12)/2
- Parameters
material1 (anisotropy.Material object) – An isotropic material that comprises layer 1, whose elastic properties are fully described by its stiffness tensor and density.
material2 (anisotropy.Material object) – An isotropic material that comprises layer 2, whose elastic properties are fully described by its stiffness tensor and density.
relative_fraction – Volumetric fraction of layer 1.
- Returns
effective_material – The composite material with effective elastic properties.
- Return type
anisotropy.Material object
3.1.2. anisotropy.effective_modelling.hudson
Modelling the effective elastic properties of a material consisting of an isotropic host matrix with unidirectionally aligned cracks.
For more information, please read:
Hudson, J.A., 1980, September. Overall properties of a cracked solid. In Mathematical Proceedings of the Cambridge Philosophical Society (Vol. 88, No. 2, pp. 371-384). Cambridge University Press.
Hudson, J.A., 1981. Wave speeds and attenuation of elastic waves in material containing cracks. Geophysical Journal International, 64(1), pp.133-150.
Hudson, J.A., 1986. A higher order approximation to the wave propagation constants for a cracked solid. Geophysical Journal International, 87(1), pp.265-274.
If you use this code, please cite these articles.
- copyright
2021–2022, AnisotroPy developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- anisotropy.effective_modelling.hudson.model(matrix, inclusions, aspect_ratio, crack_density)[source]
Calculate the effective elastic properties of a material comprised of an isotropic matrix hosting penny-shaped cracks (ellipsoidal inclusions) that may contain a different isotropic material. The effective anisotropy is controlled by the crack density and the aspect ratio of the inclusions.
The theory is valid when the crack density (= Na^2 / nu) << 1, where N is the number of cracks of radius a in a volume nu.
The expressions used describe a system that is transversely isotropic about x1, as written in Crampin (1984). The original expressions in Hudson (1982) correspond to a system that is transversely isotropic about x3. They give the most general expressions, but can be tuned to model various different crack types:
Dry cracks: Set vpi, vsi, and rhoi to 0. Water-filled: Set vsi to 0 and vpi*rhoi = 2.25x10E9
- Parameters
matrix (anisotropy.Material object) – An isotropic material that forms the host matrix for the inclusions, whose elastic properties are fully described by its stiffness tensor and density.
inclusions (anisotropy.Material object) – An isotropic material that occupies the inclusions, whose elastic properties are fully described by its stiffness tensor and density.
aspect_ratio (float) – Aspect ratio of inclusions.
crack_density (float) – Volumetric fraction of the inclusions within the composite material.
- Returns
effective_material – The composite material with effective elastic properties.
- Return type
anisotropy.Material object
3.1.3. anisotropy.effective_modelling.n_layers
Modelling of the effective seismic anisotropy of a system of N arbitrarily anisotropic layers, based either on the analytical system of equations outlined by Silver and Savage (1994) or by directly applying a sequence of splitting operators to the first-derivative of a Gaussian wavelet.
For more information, please read:
Silver, P.G. and Savage, M.K., 1994. The interpretation of shear-wave splitting parameters in the presence of two anisotropic layers. Geophysical Journal International, 119(3), pp.949-963.
If you use this code, please cite this article.
- copyright
2021–2022, AnisotroPy developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- class anisotropy.effective_modelling.n_layers.ElasticLayer(material, thickness, ψ_a, dip=0.0, fractional_alignment=0.3)[source]
Bases:
objectSmall class to encapsulate information that defines a layer composed of an (an)isotropic material and provide a simple interface to a suite of useful functions to query its properties.
- material
(An)isotropic material, whose elastic properties are fully described by its stiffness tensor and density.
- Type
anisotropy.Material object
- thickness
Thickness of the elastic layer, in km.
- Type
float
- ψ_a
Azimuthal alignment of crystalline lattice a-axis of the material.
- Type
float
- dip
Angle of dip of the layer, relative to horizontal, in degrees. Defaults to 0.
- Type
float, optional
- fractional_alignment
Degree of mixing between completely aligned material and its isotropic component. Must be between 0 and 1.
- Type
float
- effective_fast()[source]
Calculate the effective fast orientation direction for a ray passing through the layer.
- effective_dt(angle_of_incidence, azimuth=array([0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 28, 29, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39, 40, 41, 42, 43, 44, 45, 46, 47, 48, 49, 50, 51, 52, 53, 54, 55, 56, 57, 58, 59, 60, 61, 62, 63, 64, 65, 66, 67, 68, 69, 70, 71, 72, 73, 74, 75, 76, 77, 78, 79, 80, 81, 82, 83, 84, 85, 86, 87, 88, 89, 90, 91, 92, 93, 94, 95, 96, 97, 98, 99, 100, 101, 102, 103, 104, 105, 106, 107, 108, 109, 110, 111, 112, 113, 114, 115, 116, 117, 118, 119, 120, 121, 122, 123, 124, 125, 126, 127, 128, 129, 130, 131, 132, 133, 134, 135, 136, 137, 138, 139, 140, 141, 142, 143, 144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155, 156, 157, 158, 159, 160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171, 172, 173, 174, 175, 176, 177, 178, 179, 180, 181, 182, 183, 184, 185, 186, 187, 188, 189, 190, 191, 192, 193, 194, 195, 196, 197, 198, 199, 200, 201, 202, 203, 204, 205, 206, 207, 208, 209, 210, 211, 212, 213, 214, 215, 216, 217, 218, 219, 220, 221, 222, 223, 224, 225, 226, 227, 228, 229, 230, 231, 232, 233, 234, 235, 236, 237, 238, 239, 240, 241, 242, 243, 244, 245, 246, 247, 248, 249, 250, 251, 252, 253, 254, 255, 256, 257, 258, 259, 260, 261, 262, 263, 264, 265, 266, 267, 268, 269, 270, 271, 272, 273, 274, 275, 276, 277, 278, 279, 280, 281, 282, 283, 284, 285, 286, 287, 288, 289, 290, 291, 292, 293, 294, 295, 296, 297, 298, 299, 300, 301, 302, 303, 304, 305, 306, 307, 308, 309, 310, 311, 312, 313, 314, 315, 316, 317, 318, 319, 320, 321, 322, 323, 324, 325, 326, 327, 328, 329, 330, 331, 332, 333, 334, 335, 336, 337, 338, 339, 340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351, 352, 353, 354, 355, 356, 357, 358, 359, 360]))[source]
Calculates the effective delay times of a ray traversing the layer with a given angle of incidence with the surface, as measured from the vertical.
- Parameters
angle_of_incidence (float) – Angle of incidence of a ray with the surface, relative to vertical, in degrees.
azimuth (float or list of floats, optional) – Azimuth(s) of ray(s) at which to calculate traversal distances. Defaults to the full range of azimuths at 1 degree increments.
- Returns
effective_dt – Effective delay times.
- Return type
numpy.ndarray of float
- effective_fast(angle_of_incidence, azimuth=array([0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 28, 29, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39, 40, 41, 42, 43, 44, 45, 46, 47, 48, 49, 50, 51, 52, 53, 54, 55, 56, 57, 58, 59, 60, 61, 62, 63, 64, 65, 66, 67, 68, 69, 70, 71, 72, 73, 74, 75, 76, 77, 78, 79, 80, 81, 82, 83, 84, 85, 86, 87, 88, 89, 90, 91, 92, 93, 94, 95, 96, 97, 98, 99, 100, 101, 102, 103, 104, 105, 106, 107, 108, 109, 110, 111, 112, 113, 114, 115, 116, 117, 118, 119, 120, 121, 122, 123, 124, 125, 126, 127, 128, 129, 130, 131, 132, 133, 134, 135, 136, 137, 138, 139, 140, 141, 142, 143, 144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155, 156, 157, 158, 159, 160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171, 172, 173, 174, 175, 176, 177, 178, 179, 180, 181, 182, 183, 184, 185, 186, 187, 188, 189, 190, 191, 192, 193, 194, 195, 196, 197, 198, 199, 200, 201, 202, 203, 204, 205, 206, 207, 208, 209, 210, 211, 212, 213, 214, 215, 216, 217, 218, 219, 220, 221, 222, 223, 224, 225, 226, 227, 228, 229, 230, 231, 232, 233, 234, 235, 236, 237, 238, 239, 240, 241, 242, 243, 244, 245, 246, 247, 248, 249, 250, 251, 252, 253, 254, 255, 256, 257, 258, 259, 260, 261, 262, 263, 264, 265, 266, 267, 268, 269, 270, 271, 272, 273, 274, 275, 276, 277, 278, 279, 280, 281, 282, 283, 284, 285, 286, 287, 288, 289, 290, 291, 292, 293, 294, 295, 296, 297, 298, 299, 300, 301, 302, 303, 304, 305, 306, 307, 308, 309, 310, 311, 312, 313, 314, 315, 316, 317, 318, 319, 320, 321, 322, 323, 324, 325, 326, 327, 328, 329, 330, 331, 332, 333, 334, 335, 336, 337, 338, 339, 340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351, 352, 353, 354, 355, 356, 357, 358, 359, 360]))[source]
Calculates the effective orientation of polarisation of the fast shear wave traversing the layer with a given angle of incidence with the surface, as measured from the vertical.
- Parameters
angle_of_incidence (float) – Angle of incidence of a ray with the surface, relative to vertical, in degrees.
azimuth (float or list of floats, optional) – Azimuth(s) of ray(s) at which to calculate traversal distances. Defaults to the full range of azimuths at 1 degree increments.
- Returns
effective_fast – Effective orientation of polarisation of the fast shear wave.
- Return type
numpy.ndarray of float
- property material
Get and set the constitutive material of the layer.
- traversal_distance(angle_of_incidence, azimuth=array([0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 28, 29, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39, 40, 41, 42, 43, 44, 45, 46, 47, 48, 49, 50, 51, 52, 53, 54, 55, 56, 57, 58, 59, 60, 61, 62, 63, 64, 65, 66, 67, 68, 69, 70, 71, 72, 73, 74, 75, 76, 77, 78, 79, 80, 81, 82, 83, 84, 85, 86, 87, 88, 89, 90, 91, 92, 93, 94, 95, 96, 97, 98, 99, 100, 101, 102, 103, 104, 105, 106, 107, 108, 109, 110, 111, 112, 113, 114, 115, 116, 117, 118, 119, 120, 121, 122, 123, 124, 125, 126, 127, 128, 129, 130, 131, 132, 133, 134, 135, 136, 137, 138, 139, 140, 141, 142, 143, 144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155, 156, 157, 158, 159, 160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171, 172, 173, 174, 175, 176, 177, 178, 179, 180, 181, 182, 183, 184, 185, 186, 187, 188, 189, 190, 191, 192, 193, 194, 195, 196, 197, 198, 199, 200, 201, 202, 203, 204, 205, 206, 207, 208, 209, 210, 211, 212, 213, 214, 215, 216, 217, 218, 219, 220, 221, 222, 223, 224, 225, 226, 227, 228, 229, 230, 231, 232, 233, 234, 235, 236, 237, 238, 239, 240, 241, 242, 243, 244, 245, 246, 247, 248, 249, 250, 251, 252, 253, 254, 255, 256, 257, 258, 259, 260, 261, 262, 263, 264, 265, 266, 267, 268, 269, 270, 271, 272, 273, 274, 275, 276, 277, 278, 279, 280, 281, 282, 283, 284, 285, 286, 287, 288, 289, 290, 291, 292, 293, 294, 295, 296, 297, 298, 299, 300, 301, 302, 303, 304, 305, 306, 307, 308, 309, 310, 311, 312, 313, 314, 315, 316, 317, 318, 319, 320, 321, 322, 323, 324, 325, 326, 327, 328, 329, 330, 331, 332, 333, 334, 335, 336, 337, 338, 339, 340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351, 352, 353, 354, 355, 356, 357, 358, 359, 360]))[source]
Calculates effective distance that a ray traverses when passing through the layer at some given angle of incidence
- Parameters
angle_of_incidence (float) – Angle of incidence of a ray with the surface, relative to vertical, in degrees.
azimuth (float or list of floats, optional) – Azimuth(s) of ray(s) at which to calculate traversal distances. Defaults to the full range of azimuths at 1 degree increments.
- Returns
distance – Effective distance travelled within the layer, in km.
- Return type
float
- anisotropy.effective_modelling.n_layers.fit_2layer_model(observations, material, angle_of_incidence)[source]
Perform a gridsearch over the 4 free parameters that define a 2 layer model, namely (T_upper, φ_upper) and (T_lower, φ_lower)
- TODO
Write tests
Testing the free parameters in increments of 5/10 km (thicknesses) and 2.5/5 degrees for azimuth? Already doing some additional model sweeping in the _calculate_misfit function.
Allow the code to smartly account for the phase? e.g. calculate the angle of incidence for different phases at different distances - we are searching over how the phase interacts with the model, so makes sense to include this information.
- Parameters
observations (pandas.DataFrame object) – This should either be a DataFrame or perhaps a custom data format develop as part of the overall package. This would mean that results generated with other splitting codes could be used by simply writing a parser that converts the outputs of the external codes into the custom data file format.
material (anisotropy.Material object) – (An)isotropic material, whose elastic properties are fully described by its stiffness tensor and density.
angle_of_incidence (float) – Angle of incidence of a ray with the surface, relative to vertical, in degrees.
- anisotropy.effective_modelling.n_layers.model(frequency, angle_of_incidence, layers, azimuth=array([0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 28, 29, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39, 40, 41, 42, 43, 44, 45, 46, 47, 48, 49, 50, 51, 52, 53, 54, 55, 56, 57, 58, 59, 60, 61, 62, 63, 64, 65, 66, 67, 68, 69, 70, 71, 72, 73, 74, 75, 76, 77, 78, 79, 80, 81, 82, 83, 84, 85, 86, 87, 88, 89, 90, 91, 92, 93, 94, 95, 96, 97, 98, 99, 100, 101, 102, 103, 104, 105, 106, 107, 108, 109, 110, 111, 112, 113, 114, 115, 116, 117, 118, 119, 120, 121, 122, 123, 124, 125, 126, 127, 128, 129, 130, 131, 132, 133, 134, 135, 136, 137, 138, 139, 140, 141, 142, 143, 144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155, 156, 157, 158, 159, 160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171, 172, 173, 174, 175, 176, 177, 178, 179, 180, 181, 182, 183, 184, 185, 186, 187, 188, 189, 190, 191, 192, 193, 194, 195, 196, 197, 198, 199, 200, 201, 202, 203, 204, 205, 206, 207, 208, 209, 210, 211, 212, 213, 214, 215, 216, 217, 218, 219, 220, 221, 222, 223, 224, 225, 226, 227, 228, 229, 230, 231, 232, 233, 234, 235, 236, 237, 238, 239, 240, 241, 242, 243, 244, 245, 246, 247, 248, 249, 250, 251, 252, 253, 254, 255, 256, 257, 258, 259, 260, 261, 262, 263, 264, 265, 266, 267, 268, 269, 270, 271, 272, 273, 274, 275, 276, 277, 278, 279, 280, 281, 282, 283, 284, 285, 286, 287, 288, 289, 290, 291, 292, 293, 294, 295, 296, 297, 298, 299, 300, 301, 302, 303, 304, 305, 306, 307, 308, 309, 310, 311, 312, 313, 314, 315, 316, 317, 318, 319, 320, 321, 322, 323, 324, 325, 326, 327, 328, 329, 330, 331, 332, 333, 334, 335, 336, 337, 338, 339, 340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351, 352, 353, 354, 355, 356, 357, 358, 359, 360]), method='silver_savage')[source]
Calculate the effective splitting of a ray passing through N arbitrarily anisotropic layers.
- Parameters
frequency (int or float) – Frequency of effective wavelet for purposes of modelling.
angle_of_incidence (float) – Angle of incidence of a ray with the surface, relative to vertical, in degrees.
layers (list of anisotropy.ElasticLayer object) – Layers in the model.
azimuth (float or list of floats, optional) – Azimuth(s) of ray(s) at which to calculate traversal distances. Defaults to the full range of azimuths at 1 degree increments.
method (str, {"silver_savage", "gaussian_wavelet"}) – Choose the method of modelling - either analytical or forward modelling.
- Returns
effective_results – Effective fast direction and delay times for the N-layer system.
- Return type
list of list of float
3.1.4. anisotropy.effective_modelling.tandon_weng
Modelling the effective elastic properties of a material consisting of an isotropic host matrix with unidirectionally aligned isotropic spheroid inclusions.
For more information, please read:
Tandon, G.P. and Weng, G.J., 1984. The effect of aspect ratio of inclusions on the elastic properties of unidirectionally aligned composites. Polymer composites, 5(4), pp.327-333.
If you use this code, please cite this article.
- copyright
2021–2022, AnisotroPy developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- anisotropy.effective_modelling.tandon_weng.model(matrix, inclusions, alpha, fraction_inc)[source]
Calculate the averaged elastic parameters as described in Tandon and Weng (1984), hereafter TW84.
Their theory describes the effective elastic properties of a medium consisting of an isotropic host matrix with unidirectionally aligned isotropic spheroid inclusions.
- Parameters
matrix (anisotropy.Material object) – An isotropic material that forms the host matrix for the inclusions, whose elastic properties are fully described by its stiffness tensor and density.
inclusions (anisotropy.Material object) – An isotropic material that occupies the inclusions, whose elastic properties are fully described by its stiffness tensor and density.
alpha (float) – Aspect ratio of inclusions.
fraction_inc (float) – Volumetric fraction of the inclusions within the composite material.
- Returns
effective_material – The composite material with effective elastic properties.
- Return type
anisotropy.Material object