.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "tutorials/ch1_fundamentals/ch1_1_basics.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code .. rst-class:: sphx-glr-example-title .. _sphx_glr_tutorials_ch1_fundamentals_ch1_1_basics.py: 1.1 -Basics of geological modeling with GemPy ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ .. GENERATED FROM PYTHON SOURCE LINES 6-8 .. code-block:: Python import numpy as np .. GENERATED FROM PYTHON SOURCE LINES 9-11 Importing Necessary Libraries """""""""""""""""""""""""""""" .. GENERATED FROM PYTHON SOURCE LINES 11-14 .. code-block:: Python import gempy as gp import gempy_viewer as gpv .. GENERATED FROM PYTHON SOURCE LINES 15-48 Importing and Defining Input Data """"""""""""""""""""""""""""""""" :obj:`gempy.core.data.GeoModel` GemPy uses Python objects to store the data that builds the geological model. The main data classes include: - :obj:`gempy.core.data.GeoModel` - :obj:`gempy.core.data.StructuralFrame` - :obj:`gempy.core.data.StructuralGroup` - :obj:`gempy.core.data.StructuralElement` - :obj:`gempy.core.data.SurfacePointsTable` - :obj:`gempy.core.data.OrientationsTable` - :obj:`gempy.core.data.Grid` You can also create data from raw CSV files (comma-separated values). This could be useful if you are exporting model data from a different program or creating it in a spreadsheet software like Microsoft Excel or LibreOffice Calc. In this tutorial, we'll use CSV files to generate input data. You can find these example files in the `gempy data` repository on GitHub. The data consists of x, y, and z positional values for all surface points and orientation measurements. Additional data includes poles, azimuth and polarity (or the gradient components). Surface points are assigned a formation, which can be a lithological unit (like "Sandstone") or a structural feature (like "Main Fault"). It's important to note that, in GemPy, interface position points mark the **bottom** of a layer. If you need points to represent the top of a formation (for example, when modeling an intrusion), you can define an inverted orientation measurement. While generating data from CSV files, we also need to define the model's real extent in x, y, and z. This extent defines the area used for interpolation and many of the plotting functions. We also set a resolution to establish a regular grid right away. This resolution will dictate the number of voxels used during modeling. We're using a medium resolution of 50x50x50 here, which results in 125,000 voxels. The model extent should enclose all relevant data in a representative space. As our model voxels are prisms rather than cubes, the resolution can differ from the extent. However, it is recommended to avoid going beyond 100 cells in each direction (1,000,000 voxels) to prevent excessive computational costs. .. GENERATED FROM PYTHON SOURCE LINES 50-64 .. code-block:: Python data_path = 'https://raw.githubusercontent.com/cgre-aachen/gempy_data/master/' geo_model: gp.data.GeoModel = gp.create_geomodel( project_name='Tutorial_ch1_1_Basics', extent=[0, 2000, 0, 2000, 0, 750], refinement=4, # * Here we define the number of octree levels. If octree levels are defined, the resolution is ignored. importer_helper=gp.data.ImporterHelper( path_to_orientations=data_path + "/data/input_data/getting_started/simple_fault_model_orientations.csv", path_to_surface_points=data_path + "/data/input_data/getting_started/simple_fault_model_points.csv", hash_surface_points="4cdd54cd510cf345a583610585f2206a2936a05faaae05595b61febfc0191563", hash_orientations="7ba1de060fc8df668d411d0207a326bc94a6cdca9f5fe2ed511fd4db6b3f3526" ) ) .. rst-class:: sphx-glr-script-out .. code-block:: none Surface points hash: 4cdd54cd510cf345a583610585f2206a2936a05faaae05595b61febfc0191563 Orientations hash: 7ba1de060fc8df668d411d0207a326bc94a6cdca9f5fe2ed511fd4db6b3f3526 .. GENERATED FROM PYTHON SOURCE LINES 65-72 .. admonition:: New in GemPy 3! GemPy 3 has introduced the ``ImporterHelper`` class to streamline importing data from various sources. This class simplifies the process of passing multiple arguments needed for importing data and will likely see further extensions in the future. Currently, one of its uses is to handle `pooch` arguments for downloading data from the internet. .. GENERATED FROM PYTHON SOURCE LINES 72-76 .. code-block:: Python # The input data can be reviewed using the properties `surface_points` and `orientations`. However, note that at this point, # the sequence of formations and their assignment to series are still arbitrary. We will rectify this in the subsequent steps. .. GENERATED FROM PYTHON SOURCE LINES 77-79 .. code-block:: Python geo_model.surface_points_copy .. raw:: html
XYZidnugget
700.001000.00300.00269027820.00
600.001000.00200.00269027820.00
500.001000.00100.00269027820.00
1000.001000.00600.00269027820.00
1100.001000.00700.00269027820.00
1000.0050.00350.001171800310.00
1000.00150.00333.331171800310.00
1000.00300.00333.331171800310.00
1000.00500.00366.671171800310.00
1000.001000.00433.331171800310.00
...............
1100.001700.00473.334123465060.00
1100.001950.00490.004123465060.00
0.001000.00433.334123465060.00
300.001000.00400.004123465060.00
600.001000.00366.674123465060.00
1300.001000.00566.674123465060.00
1600.001000.00550.004123465060.00
1900.001000.00566.674123465060.00
1700.00500.00533.334123465060.00
1700.001500.00516.674123465060.00


.. GENERATED FROM PYTHON SOURCE LINES 80-82 .. code-block:: Python geo_model.orientations_copy .. raw:: html
XYZG_xG_yG_zidnugget
500.001000.00300.00-0.95-0.000.32269027820.01
400.001000.00420.000.320.000.952733860110.01
1000.001000.00300.000.320.000.953897245790.01


.. GENERATED FROM PYTHON SOURCE LINES 83-111 Declaring the Sequential Order of Geological Formations ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ In our model, we want the geological units to appear in the correct chronological order. Such order could be determined by a sequence of stratigraphic deposition, unconformities due to erosion, or other lithological genesis events like igneous intrusions. A similar age-related order is declared for faults in our model. In GemPy, we use the function `gempy.map_stack_to_surfaces` to assign formations or faults to different sequential series by declaring them in a Python dictionary. The correct ordering of series is crucial for model construction! It's possible to assign several surfaces to one series. The order of units within a series only affects the color code, so we recommend maintaining consistency. The order can be defined by simply changing the order of the lists within `gempy.core.data.StructuralFrame.structural_groups` and `gempy.core.data.StructuralGroups.elements` attributes. Faults are treated as independent groups and must be younger than the groups they affect. The relative order between different faults defines their tectonic relationship (the first entry is the youngest). For a model with simple sequential stratigraphy, all layer formations can be assigned to one series without an issue. All unit boundaries and their order would then be determined by interface points. However, to model more complex lithostratigraphical relations and interactions, separate series definition becomes important. For example, modeling an unconformity or an intrusion that disrupts older stratigraphy would require declaring a "newer" series. By default, we create a simple sequence inferred from the data: .. GENERATED FROM PYTHON SOURCE LINES 113-115 .. code-block:: Python geo_model.structural_frame .. raw:: html
Structural Groups: StructuralGroup:
Name:default_formation
Structural Relation:StackRelationType.ERODE
Elements:
StructuralElement:
Name:Main_Fault

StructuralElement:
Name:Sandstone_1

StructuralElement:
Name:Sandstone_2

StructuralElement:
Name:Shale

StructuralElement:
Name:Siltstone
Fault Relations:
default_fo...
default_formation
True
False


.. GENERATED FROM PYTHON SOURCE LINES 116-126 Our example model comprises four main layers (plus an underlying basement that is automatically generated by GemPy) and one main normal fault displacing those layers. Assuming a simple stratigraphy where each younger unit was deposited onto the underlying older one, we can assign these layer formations to one series called "Strat\_Series". For the fault, we declare a respective "Fault\_Series" as the first key entry in the mapping dictionary. We could give any other names to these series, the formations however have to be referred to as named in the input data. .. GENERATED FROM PYTHON SOURCE LINES 129-140 .. code-block:: Python gp.map_stack_to_surfaces( gempy_model=geo_model, mapping_object= # TODO: This mapping I do not like it too much. We should be able to do it passing the data objects directly { "Fault_Series": 'Main_Fault', "Strat_Series": ('Sandstone_2', 'Siltstone', 'Shale', 'Sandstone_1') } ) geo_model.structural_frame # Display the resulting structural frame .. raw:: html
Structural Groups: StructuralGroup:
Name:Fault_Series
Structural Relation:StackRelationType.ERODE
Elements:
StructuralElement:
Name:Main_Fault

StructuralGroup:
Name:Strat_Series
Structural Relation:StackRelationType.ERODE
Elements:
StructuralElement:
Name:Sandstone_2

StructuralElement:
Name:Siltstone

StructuralElement:
Name:Shale

StructuralElement:
Name:Sandstone_1
Fault Relations:
Fault_Seri...Strat_Seri...
Fault_Series
Strat_Series
True
False


.. GENERATED FROM PYTHON SOURCE LINES 141-146 .. code-block:: Python gp.set_is_fault( frame=geo_model.structural_frame, fault_groups=['Fault_Series'] ) .. raw:: html
Structural Groups: StructuralGroup:
Name:Fault_Series
Structural Relation:StackRelationType.FAULT
Elements:
StructuralElement:
Name:Main_Fault

StructuralGroup:
Name:Strat_Series
Structural Relation:StackRelationType.ERODE
Elements:
StructuralElement:
Name:Sandstone_2

StructuralElement:
Name:Siltstone

StructuralElement:
Name:Shale

StructuralElement:
Name:Sandstone_1
Fault Relations:
Fault_Seri...Strat_Seri...
Fault_Series
Strat_Series
True
False


.. GENERATED FROM PYTHON SOURCE LINES 147-155 Now, all surfaces have been assigned to a series and are displayed in the correct order (from young to old). Returning Information from Our Input Data ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ Our model input data, named "geo_model", contains all essential information for constructing our model. You can access different types of information by accessing the attributes. For instance, you can retrieve the coordinates of our modeling grid as follows: .. GENERATED FROM PYTHON SOURCE LINES 158-160 .. code-block:: Python geo_model.grid .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 161-170 Visualizing input data ~~~~~~~~~~~~~~~~~~~~~~ We can also visualize our input data. This might for example be useful to check if all points and measurements are defined the way we want them to. Using the function :obj:`gempy_viewer.plot2d`, we attain a 2D projection of our data points onto a plane of chosen *direction* (we can choose this attribute to be either :math:`x`, :math:`y` or :math:`z`). .. GENERATED FROM PYTHON SOURCE LINES 173-175 .. code-block:: Python plot = gpv.plot_2d(geo_model, show_lith=False, show_boundaries=False) .. image-sg:: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_001.png :alt: Cell Number: mid Direction: y :srcset: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_001.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 176-180 Using :obj:`gempy_viewer.plot_3d`, # we can also visualize this data in 3D. Note that direct 3D visualization in GemPy requires `the Visualization Toolkit `__ (VTK) to be installed. .. GENERATED FROM PYTHON SOURCE LINES 182-184 .. code-block:: Python gpv.plot_3d(geo_model, image=False, plotter_type='basic') .. image-sg:: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_002.png :alt: ch1 1 basics :srcset: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_002.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 185-198 Model Generation ~~~~~~~~~~~~~~~~ Once we've correctly defined all our primary information in our `gempy.core.data.GeoModel` object (referred to as `geo_model` in these tutorials), we can proceed to the next step: preparing the input data for interpolation. .. admonition:: New in GemPy 3! Numpy and TensorFlow backend Unlike previous versions, GemPy 3 doesn't rely on `theano` or `asera`. Instead, it utilizes `numpy` or `tensorflow`. Consequently, we no longer need to recompile all theano functions (TensorFlow uses eager execution; we found no notable speed difference after profiling the XLA compiler). .. GENERATED FROM PYTHON SOURCE LINES 200-204 The parameters used for the interpolation are stored in `gempy.core.data.GeoModel.interpolation_options`. These parameters have sensible default values that you can modify if necessary. However, we advise caution when changing these parameters unless you fully understand their implications. .. GENERATED FROM PYTHON SOURCE LINES 206-207 Display the current interpolation options .. GENERATED FROM PYTHON SOURCE LINES 207-209 .. code-block:: Python geo_model.interpolation_options .. raw:: html
InterpolationOptions
kernel_options{'range': 5, 'c_o': 10, 'uni_degree': 1, 'i_res': 4, 'gi_res': 2, 'number_dimensions': 3, 'kernel_function': , derivative_div_r=, second_derivative=, consume_sq_distance=False)>, 'kernel_solver': , 'compute_condition_number': False, 'optimizing_condition_number': False, 'condition_number': None}
_number_octree_levels4
current_octree_level0
block_solutions_typeBlockSolutionType.OCTREE
compute_scalar_gradientFalse
mesh_extractionTrue
mesh_extraction_masking_optionsMeshExtractionMaskingOptions.INTERSECT
mesh_extraction_fancyTrue
debugTrue
debug_water_tightFalse
sigmoid_slope50000
evaluation_chunk_size50000
_number_octree_levels_surface4
_model_nameTutorial_ch1_1_Basics


.. GENERATED FROM PYTHON SOURCE LINES 210-214 With all our prerequisites in place, we can now compute our complete geological model using :func:`gempy.compute_model`. This function returns a :obj:`gempy.core.data.Solutions` object. The following sections illustrate these different model solutions and how to utilize them. .. GENERATED FROM PYTHON SOURCE LINES 216-217 Compute the geological model and get the solutions .. GENERATED FROM PYTHON SOURCE LINES 217-220 .. code-block:: Python sol = gp.compute_model(geo_model) sol .. rst-class:: sphx-glr-script-out .. code-block:: none Setting Backend To: AvailableBackends.numpy .. raw:: html
Solutions: 4 Octree Levels, 5 DualContouringMeshes


.. GENERATED FROM PYTHON SOURCE LINES 221-222 Solutions are also stored within the :obj:`gempy.core.data.GeoModel` object, for future reference. .. GENERATED FROM PYTHON SOURCE LINES 224-226 .. code-block:: Python geo_model.solutions .. raw:: html
Solutions: 4 Octree Levels, 5 DualContouringMeshes


.. GENERATED FROM PYTHON SOURCE LINES 227-233 Direct model visualization in GemPy ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ Model solutions can be easily visualized in 2D sections in GemPy directly. Let's take a look at our lithology block: .. GENERATED FROM PYTHON SOURCE LINES 235-237 .. code-block:: Python gpv.plot_2d(geo_model, show_data=True, cell_number="mid", direction='y') .. image-sg:: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_003.png :alt: Cell Number: mid Direction: y :srcset: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_003.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 238-248 With ``cell_number=mid``, we have chosen a section going through the middle of our block. We have moved in ``direction='y'``, the plot thus depicts a plane parallel to the :math:`x`- and :math:`y`-axes. Setting ``show_data=True``, we could plot original data together with the results. Changing the values for ``cell_number`` and ``direction``, we can move through our 3D block model and explore it by looking at different 2D planes. We can do the same with the underlying scalar-field solution: .. GENERATED FROM PYTHON SOURCE LINES 251-259 .. code-block:: Python gpv.plot_2d( model=geo_model, series_n=0, # This will plot the scalar field used for the fault show_data=False, show_scalar=True, show_lith=False ) .. image-sg:: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_004.png :alt: Cell Number: mid Direction: y :srcset: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_004.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 260-268 .. code-block:: Python gpv.plot_2d( model=geo_model, series_n=1, # This will plot the scalar field used for the stratigraphy show_data=False, show_scalar=True, show_lith=False ) .. image-sg:: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_005.png :alt: Cell Number: mid Direction: y :srcset: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_005.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 269-284 Dual Contouring and vtk visualization ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ In addition to 2D sections we can extract surfaces to visualize in 3D renderers. Surfaces can be visualized as 3D triangle complexes in VTK (see function plot\_surfaces\_3D below). To create these triangles, we need to extract respective vertices and simplices from the potential fields of lithologies and faults. This process is automatized in GemPy using dual contouring in the :obj:`gempy_engine`. .. admonition:: New in GemPy 3! Dual Contouring GemPy 3 uses dual contouring to extract surfaces from the scalar fields. The method is completely coded in :obj:`gempy_engine` what also enables further improvements in the midterm. This method is more efficient to use together with octrees and suited better the new capabilities of gempy3. .. GENERATED FROM PYTHON SOURCE LINES 286-288 .. code-block:: Python gpv.plot_3d(geo_model, show_data=False, image=False, plotter_type='basic') .. image-sg:: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_006.png :alt: ch1 1 basics :srcset: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_006.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 289-294 Adding topography ~~~~~~~~~~~~~~~~~~~~~~~~~~~ In gempy we can add more grid types for different purposes. We will explore this concept in more detail in the next tutorials (:doc:`ch1_3a_grids`). For now, we will just add a topography grid to our model. This grid allows us to intersect the surfaces as well as compute a high resolution geological mal. .. GENERATED FROM PYTHON SOURCE LINES 297-316 .. code-block:: Python gp.set_topography_from_random( grid=geo_model.grid, fractal_dimension=1.2, d_z=np.array([350, 750]), topography_resolution=np.array([50, 50]), ) gp.compute_model(geo_model) gpv.plot_2d(geo_model, show_topography=True) gpv.plot_3d( model=geo_model, plotter_type='basic', show_topography=True, show_surfaces=True, show_lith=True, image=False ) .. image-sg:: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_007.png :alt: ch1 1 basics :srcset: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_007.png :class: sphx-glr-single-img .. image-sg:: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_008.png :alt: Cell Number: mid Direction: y :srcset: /tutorials/ch1_fundamentals/images/sphx_glr_ch1_1_basics_008.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none Active grids: ['topography'] Setting Backend To: AvailableBackends.numpy .. GENERATED FROM PYTHON SOURCE LINES 317-322 Compute at a given location ~~~~~~~~~~~~~~~~~~~~~~~~~~~ This is done by modifying the grid to a custom grid and recomputing. .. GENERATED FROM PYTHON SOURCE LINES 324-331 .. code-block:: Python x_i = np.array([[1000, 1000, 1000]]) lith_values_at_coords: np.ndarray = gp.compute_model_at( gempy_model=geo_model, at=x_i ) lith_values_at_coords .. rst-class:: sphx-glr-script-out .. code-block:: none Active grids: ['custom' 'topography'] Setting Backend To: AvailableBackends.numpy array([5.]) .. GENERATED FROM PYTHON SOURCE LINES 332-333 Therefore if we just want the value at **x\_i**: .. GENERATED FROM PYTHON SOURCE LINES 335-337 .. code-block:: Python geo_model.solutions.raw_arrays.custom .. rst-class:: sphx-glr-script-out .. code-block:: none array([5.]) .. GENERATED FROM PYTHON SOURCE LINES 338-345 .. admonition:: Work in progress GemPy3 model serialization is currently being redisigned. Therefore, at the current version, there is not a build in method to save the model. However, since now the data model should be completely robust, you should be able to save the :obj:`gempy.core.data.GeoModel` and all its attributes using the standard python library [pickle](https://docs.python.org/3/library/pickle.html) sphinx_gallery_thumbnail_number = -2 .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 6.424 seconds) .. _sphx_glr_download_tutorials_ch1_fundamentals_ch1_1_basics.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: ch1_1_basics.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: ch1_1_basics.py ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_