.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "tutorials/a_getting_started/get_started.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_a_getting_started_get_started.py: Getting Started =============== .. GENERATED FROM PYTHON SOURCE LINES 9-19 .. code-block:: Python # Importing GemPy and viewer import gempy as gp import gempy_viewer as gpv # Auxiliary libraries import numpy as np import matplotlib.pyplot as plt import matplotlib.image as mpimg .. GENERATED FROM PYTHON SOURCE LINES 20-28 Initializing the model: ~~~~~~~~~~~~~~~~~~~~~~~ Create a gempy.Model object. This object will contain all other data structures and necessary functionality. We'll also define a regular grid for this example. This grid will be used for interpolating the 3D geological model. GemPy offers different grids for various purposes. For visualization, a regular grid is most appropriate. .. GENERATED FROM PYTHON SOURCE LINES 30-39 .. code-block:: Python geo_model: gp.data.GeoModel = gp.create_geomodel( project_name='Model1', extent=[0, 791, -200, 200, -582, 0], resolution=[50, 50, 50], refinement=4, structural_frame=gp.data.StructuralFrame.initialize_default_structure() ) geo_model .. rst-class:: sphx-glr-script-out .. code-block:: none {'grid': , 'interpolation_options': 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)>, 'compute_condition_number': False, 'kernel_solver': }, number_octree_levels=4, current_octree_level=0, compute_scalar_gradient=False, mesh_extraction=True, mesh_extraction_masking_options=MeshExtractionMaskingOptions.INTERSECT, mesh_extraction_fancy=True, debug=True, debug_water_tight=False, sigmoid_slope=50000, _number_octree_levels_surface=4), 'meta': GeoModelMeta(name='Model1', creation_date=None, last_modification_date=None, owner=None), 'structural_frame': StructuralFrame( structural_groups=[ StructuralGroup( name=default_formations, structural_relation=StackRelationType.ERODE, elements=[ Element( name=surface1, color=#015482, is_active=True ) ] ) ], fault_relations= [[False]], , 'transform': {'_is_default_transform': True, 'position': array([0., 0., 0.]), 'rotation': array([0., 0., 0.]), 'scale': array([1., 1., 1.])}} .. GENERATED FROM PYTHON SOURCE LINES 40-47 Creating a figure: ~~~~~~~~~~~~~~~~~~ GemPy utilizes matplotlib for 2D and pyvista-vtk for 3D visualizations. One design goal of GemPy is real-time model construction. This means as input data is added, you can see the 3D surfaces update in real-time. Let's initialize the visualization windows. First, the 2D figure: .. GENERATED FROM PYTHON SOURCE LINES 49-51 .. code-block:: Python p2d = gpv.plot_2d(geo_model) .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_001.png :alt: Cell Number: mid Direction: y :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_001.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 52-57 Adding a model section: ^^^^^^^^^^^^^^^^^^^^^^^ In the 2D renderer, we can add several cross sections of the model. For simplicity, we'll add just one, perpendicular to y. .. GENERATED FROM PYTHON SOURCE LINES 59-64 Loading a cross-section image: ------------------------------ GemPy uses standard matplotlib axes, allowing for flexibility. Let's load an image showing the details of a couple of boreholes: .. GENERATED FROM PYTHON SOURCE LINES 66-71 .. code-block:: Python img = mpimg.imread('wells.png') p2d = gpv.plot_2d(geo_model, show=False) p2d.axes[0].imshow(img, origin='upper', alpha=.8, extent=(0, 791, -582, 0)) plt.show() .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_002.png :alt: Cell Number: mid Direction: y :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_002.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 72-73 Similarly, we can visualize in 3D using pyvista and vtk: .. GENERATED FROM PYTHON SOURCE LINES 75-77 .. code-block:: Python p3d = gpv.plot_3d(geo_model, image=True) .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_003.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_003.png :class: sphx-glr-single-img .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_004.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_004.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 78-88 Building the model: ------------------- With everything initialized, we can begin constructing the geological model. Surfaces: ~~~~~~~~~ GemPy is a surface-based interpolator. All input data must be referred to a surface, which marks the bottom of a unit. By default, GemPy surfaces are empty: .. GENERATED FROM PYTHON SOURCE LINES 90-92 .. code-block:: Python geo_model.structural_frame.structural_elements .. rst-class:: sphx-glr-script-out .. code-block:: none [Element( name=surface1, color=#015482, is_active=True ), Element( name=basement, color=#9f0052, is_active=True )] .. GENERATED FROM PYTHON SOURCE LINES 93-96 Let's begin by adding data. GemPy input data consists of surface points and orientations (perpendicular to the layers). The 2D plot provides X and Z coordinates on mouse hover (in qt5 backend). We can add a surface point like this: .. GENERATED FROM PYTHON SOURCE LINES 98-109 .. code-block:: Python gp.add_surface_points( geo_model=geo_model, x=[223], y=[0.01], z=[-94], elements_names=['surface1'] ) gpv.plot_2d(geo_model, cell_number=11) gpv.plot_3d(geo_model, image=True) .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_005.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_005.png :class: sphx-glr-single-img .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_006.png :alt: Cell Number: 11 Direction: y :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_006.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 110-111 We can now add other points for the layer: .. GENERATED FROM PYTHON SOURCE LINES 113-124 .. code-block:: Python gp.add_surface_points( geo_model=geo_model, x=[458, 612], y=[0, 0], z=[-107, -14], elements_names=['surface1', 'surface1'] ) gpv.plot_2d(geo_model, cell_number=11) gpv.plot_3d(geo_model, image=True) .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_007.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_007.png :class: sphx-glr-single-img .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_008.png :alt: Cell Number: 11 Direction: y :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_008.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 125-132 To interpolate in gempy, the minimum data needed is: a) 2 surface points per surface b) One orientation per series Let's add an orientation: .. GENERATED FROM PYTHON SOURCE LINES 134-148 .. code-block:: Python gp.add_orientations( geo_model=geo_model, x=[350], y=[1], z=[-300], elements_names=['surface1'], pole_vector=[[0, 0, 1]] ) gpv.plot_2d(geo_model, cell_number=5) gpv.plot_3d(geo_model, image=True) .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_009.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_009.png :class: sphx-glr-single-img .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_010.png :alt: Cell Number: 5 Direction: y :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_010.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 149-152 Update and Recompute Model Transform: ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ Removing auto anisotropy for this 2.5D model. .. GENERATED FROM PYTHON SOURCE LINES 152-154 .. code-block:: Python geo_model.update_transform(gp.data.GlobalAnisotropy.NONE) .. GENERATED FROM PYTHON SOURCE LINES 155-158 Interpolation: ~~~~~~~~~~~~~~ With the provided data, we can now interpolate the 3D surface. .. GENERATED FROM PYTHON SOURCE LINES 158-161 .. code-block:: Python gp.compute_model(geo_model, engine_config=gp.data.GemPyEngineConfig()) .. rst-class:: sphx-glr-script-out .. code-block:: none Setting Backend To: AvailableBackends.numpy /home/leguark/gempy/gempy/core/data/geo_model.py:164: UserWarning: You are using refinement and passing a regular grid. The resolution of the regular grid will be overwritten warnings.warn( .. raw:: html
Solutions: 4 Octree Levels, 1 DualContouringMeshes


.. GENERATED FROM PYTHON SOURCE LINES 162-163 Display interpolation kernel options: .. GENERATED FROM PYTHON SOURCE LINES 163-165 .. code-block:: Python geo_model.interpolation_options.kernel_options .. raw:: html
KernelOptions
range5
c_o10
uni_degree1
i_res4
gi_res2
number_dimensions3
kernel_functionAvailableKernelFunctions.cubic
compute_condition_numberFalse
kernel_solverSolvers.DEFAULT


.. GENERATED FROM PYTHON SOURCE LINES 166-169 Visualization: ~~~~~~~~~~~~~~ Interpolated 3D surface can be visualized both in 2D and 3D. .. GENERATED FROM PYTHON SOURCE LINES 169-176 .. code-block:: Python # 2D visualization: gpv.plot_2d(geo_model, cell_number=[5]) # 3D visualization: gpv.plot_3d(geo_model, show_surfaces=True, image=True) .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_011.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_011.png :class: sphx-glr-single-img .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_012.png :alt: Cell Number: 5 Direction: y :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_012.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none .. GENERATED FROM PYTHON SOURCE LINES 177-180 Expanding the Model with More Layers: ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ Our cross-section image displays 4 layers, yet we only defined 2. Let's add two more. .. GENERATED FROM PYTHON SOURCE LINES 180-184 .. code-block:: Python # Display current structural frame: geo_model.structural_frame .. raw:: html
Structural Groups: StructuralGroup:
Name:default_formations
Structural Relation:StackRelationType.ERODE
Elements:
StructuralElement:
Name:surface1
Fault Relations:
default_fo...
default_formations
True
False


.. GENERATED FROM PYTHON SOURCE LINES 185-188 Defining Layer 2: ~~~~~~~~~~~~~~~~~ Adding points and properties for the next layer. .. GENERATED FROM PYTHON SOURCE LINES 188-208 .. code-block:: Python element2 = gp.data.StructuralElement( name='surface2', color=next(geo_model.structural_frame.color_generator), surface_points=gp.data.SurfacePointsTable.from_arrays( x=np.array([225, 459]), y=np.array([0, 0]), z=np.array([-269, -279]), names='surface2' ), orientations=gp.data.OrientationsTable.initialize_empty() ) geo_model.structural_frame.structural_groups[0].append_element(element2) # Compute and visualize the updated model: gp.compute_model(geo_model) gpv.plot_2d(geo_model, cell_number=5, legend='force') gpv.plot_3d(geo_model, image=True) .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_013.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_013.png :class: sphx-glr-single-img .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_014.png :alt: Cell Number: 5 Direction: y :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_014.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none Setting Backend To: AvailableBackends.numpy .. GENERATED FROM PYTHON SOURCE LINES 209-212 Defining Layer 3: ~~~~~~~~~~~~~~~~~ Adding points and properties for another layer. .. GENERATED FROM PYTHON SOURCE LINES 212-233 .. code-block:: Python element3 = gp.data.StructuralElement( name='surface3', color=next(geo_model.structural_frame.color_generator), surface_points=gp.data.SurfacePointsTable.from_arrays( x=np.array([225, 464, 619]), y=np.array([0, 0, 0]), z=np.array([-439, -456, -433]), names='surface3' ), orientations=gp.data.OrientationsTable.initialize_empty() ) geo_model.structural_frame.structural_groups[0].append_element(element3) # Compute and visualize with adjusted parameters: gp.compute_model(geo_model) gpv.plot_2d(geo_model, cell_number=5, legend='force') gpv.plot_3d(geo_model, kwargs_plot_structured_grid={'opacity': .2}, image=True) .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_015.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_015.png :class: sphx-glr-single-img .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_016.png :alt: Cell Number: 5 Direction: y :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_016.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none Setting Backend To: AvailableBackends.numpy .. GENERATED FROM PYTHON SOURCE LINES 234-241 Adding a Fault: ~~~~~~~~~~~~~~~ To date, our model represents a simple depositional unit. With GemPy, we can incorporate unconformities and faults for more intricate models. Relationships are depicted as: input data (surface points/ orientations) .. GENERATED FROM PYTHON SOURCE LINES 278-279 Compute and visualize the updated model: .. GENERATED FROM PYTHON SOURCE LINES 279-283 .. code-block:: Python gp.compute_model(geo_model) gpv.plot_2d(geo_model, cell_number=5, legend='force') gpv.plot_3d(geo_model, kwargs_plot_structured_grid={'opacity': .2}) .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_018.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_018.png :class: sphx-glr-single-img .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_019.png :alt: Cell Number: 5 Direction: y :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_019.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none Setting Backend To: AvailableBackends.numpy .. GENERATED FROM PYTHON SOURCE LINES 284-288 Advanced Features: ~~~~~~~~~~~~~~~~~~ Over time, numerous capabilities have been integrated with GemPy. Here, we'll showcase a few of them. .. GENERATED FROM PYTHON SOURCE LINES 288-308 .. code-block:: Python # Topography: # GemPy offers built-in tools to manage topographic data through gdal. # For demonstration, we'll create a random topography: gp.set_topography_from_random( grid=geo_model.grid, fractal_dimension=1.9, d_z=np.array([-150, 0]), topography_resolution=np.array([200, 200]) ) # Visualize the topography: gpv.plot_2d(geo_model, cell_number=5, legend='force') gpv.plot_3d(geo_model, kwargs_plot_structured_grid={'opacity': .2}) # Calculate and visualize the area's geological map: gp.compute_model(geo_model) gpv.plot_3d(geo_model, show_topography=True) .. rst-class:: sphx-glr-horizontal * .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_020.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_020.png :class: sphx-glr-multi-img * .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_021.png :alt: get started :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_021.png :class: sphx-glr-multi-img .. image-sg:: /tutorials/a_getting_started/images/sphx_glr_get_started_022.png :alt: Cell Number: 5 Direction: y :srcset: /tutorials/a_getting_started/images/sphx_glr_get_started_022.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 309-315 Gravity Inversion: ------------------ .. admonition:: Coming soon: Gravity inversion This feature is not yet available in the current version of GemPy. .. GENERATED FROM PYTHON SOURCE LINES 317-331 Assign density values to model units: geo_model.add_surface_values([0, 2.6, 2.4, 3.2, 3.6], ['density']) Generate a centered grid around a device for improved accuracy: geo_model.set_centered_grid(centers=[[400, 0, 0]], resolution=[10, 10, 100], radius=800) Adjust the compile code for gravity computation: gp.set_interpolator(geo_model, output=['gravity'], aesara_optimizer='fast_run') Besides the interpolation, compute the model's forward gravity: gp.compute_model(geo_model) geo_model.solutions.fw_gravity sphinx_gallery_thumbnail_number = -2 .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 5.821 seconds) .. _sphx_glr_download_tutorials_a_getting_started_get_started.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: get_started.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: get_started.py ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_