Tutorial de modelamiento de ingresos a tajo con fallas de 1m usado MODFLOW6 y mf6Voronoi

Ejemplo aplicado de la simulación del flujo de agua subterránea hacia una mina a cielo abierto considerando fallas con un espesor de 1 metro. El modelo tiene 6 períodos de esfuerzo, uno de estado estacionario y 5 períodos transitorios que representan el desarrollo de la mina a lo largo de 4 años dentro de una vida total de la mina de 20 años. El tutorial muestra el proceso completo comenzando desde la creación de la malla, construcción del modelo, análisis del balance hídrico y representación tridimensional.

Tutorial

Código


from mf6Voronoi.utils import listTemplates, copyTemplate
#listTemplates()
copyTemplate('generateVoronoi','tajo')
copyTemplate('multilayeredTransient','tajo')
copyTemplate('vtkGeneration','tajo')

Part 1 : Voronoi mesh generation

#!pip install -U mf6Voronoi
import warnings ## Org
warnings.filterwarnings('ignore') ## Org

import os, sys ## Org
import geopandas as gpd ## Org
from mf6Voronoi.geoVoronoi import createVoronoi ## Org
from mf6Voronoi.meshProperties import meshShape ## Org
from mf6Voronoi.utils import initiateOutputFolder, getVoronoiAsShp ## Org
#Create mesh object specifying the coarse mesh and the multiplier
vorMesh = createVoronoi(meshName='openPit',maxRef = 200, multiplier=2.5, overlapping=False) ## Org

#Open limit layers and refinement definition layers
vorMesh.addLimit('basin','../shp/local/aoi.shp') ## Org
vorMesh.addLayer('faults','../shp/faults_v2.shp',1) #<======= Insert
vorMesh.addLayer('river','../shp/river_basin.shp',50) ## Org
vorMesh.addLayer('pit','../shp/minePlan/pitShell.shp',50) ## Org 
vorMesh.addLayer('regionalFlow','../shp/local/regionalFlow.shp',50)
#Generate point pair array
vorMesh.generateOrgDistVertices() ## Org

#Generate the point cloud and voronoi
vorMesh.createPointCloud() ## Org
vorMesh.generateVoronoi() ## Org
Hatarilabs

mf6Voronoi will have a web version in 2028

Follow us:

Hatarilabs Hatarilabs Hatarilabs Hatarilabs Hatarilabs Hatarilabs
/--------Layer faults discretization-------/
Progressive cell size list: [1, 3.5, 9.75, 25.375, 64.4375, 162.09375] m.

/--------Layer river discretization-------/
Progressive cell size list: [50, 175.0] m.

/--------Layer pit discretization-------/
Progressive cell size list: [50, 175.0] m.

/--------Layer regionalFlow discretization-------/
Progressive cell size list: [50, 175.0] m.

/----Sumary of points for voronoi meshing----/
Distributed points from layers: 4
Points from layer buffers: 26874
Points from max refinement areas: 384
Points from min refinement areas: 629
Total points inside the limit: 35446
/--------------------------------------------/

Time required for point generation: 15.48 seconds 


/----Generation of the voronoi mesh----/

Time required for voronoi generation: 5.12 seconds
#Uncomment the next two cells if you have strong differences on discretization or you have encounter an FORTRAN error while running MODFLOW6
vorMesh.checkVoronoiQuality(threshold=0.01)
/----Performing quality verification of voronoi mesh----/
Short side on polygon: 35445 with length = 0.00123
Short side on polygon: 35445 with length = 0.00123
Short side on polygon: 35445 with length = 0.00123
vorMesh.fixVoronoiShortSides()
vorMesh.generateVoronoi()
vorMesh.checkVoronoiQuality(threshold=0.01)
/----Generation of the voronoi mesh----/

Time required for voronoi generation: 5.67 seconds 


/----Performing quality verification of voronoi mesh----/
Your mesh has no edges shorter than your threshold
#Export generated voronoi mesh
initiateOutputFolder('../output') ## Org
getVoronoiAsShp(vorMesh.modelDis, shapePath='../output/'+vorMesh.modelDis['meshName']+'.shp') ## Org
The output folder ../output exists and has been cleared

/----Generation of the voronoi shapefile----/

Time required for voronoi shapefile: 10.76 seconds
# Show the resulting voronoi mesh

#open the mesh file
mesh=gpd.read_file('../output/'+vorMesh.modelDis['meshName']+'.shp') ## Org
#plot the mesh
mesh.plot(figsize=(35,25), fc='crimson', alpha=0.3, ec='teal') ## Org

Part 2 generate disv properties

# open the mesh file
mesh=meshShape('../output/'+vorMesh.modelDis['meshName']+'.shp') ## Org
# get the list of vertices and cell2d data
gridprops=mesh.get_gridprops_disv() ## Org
Creating a unique list of vertices [[x1,y1],[x2,y2],...]


100%|█████████████████████████████████████████████████████████████████████████| 35685/35685 [00:02<00:00, 12119.34it/s]



Extracting cell2d data and grid index


100%|██████████████████████████████████████████████████████████████████████████| 35685/35685 [00:17<00:00, 2097.42it/s]
#create folder
initiateOutputFolder('../json') ## Org

#export disv
mesh.save_properties('../json/disvDict.json') ## Org
The output folder ../json exists and has been cleared

Part 2a: generate disv properties

import sys, json, os ## Org
import rasterio, flopy ## Org
import numpy as np ## Org
import matplotlib.pyplot as plt ## Org
import geopandas as gpd ## Org
from mf6Voronoi.meshProperties import meshShape ## Org
from shapely.geometry import MultiLineString ## Org
from mf6Voronoi.tools.cellWork import getLayCellElevTupleFromRaster, getLayCellElevTupleFromElev
C:\Users\saulm\anaconda3\Lib\site-packages\pyvista\examples\downloads.py:93: DeprecationWarning: support for supplying keyword arguments to pathlib.PurePath is deprecated and scheduled for removal in Python 3.14
  Path(USER_DATA_PATH, exist_ok=True).mkdir()
C:\Users\saulm\anaconda3\Lib\site-packages\pyvista\examples\downloads.py:98: UserWarning: Unable to access C:\Users\saulm\AppData\Local\pyvista_3\pyvista_3\Cache. Manually specify the PyVistaexamples cache with the PYVISTA_USERDATA_PATH environment variable.
  warnings.warn(
# open the json file
with open('../json/disvDict.json') as file: ## Org
    gridProps = json.load(file) ## Org
cell2d = gridProps['cell2d']           #cellid, cell centroid xy, vertex number and vertex id list
vertices = gridProps['vertices']       #vertex id and xy coordinates
ncpl = gridProps['ncpl']               #number of cells per layer
nvert = gridProps['nvert']             #number of verts
centroids=gridProps['centroids']       #cell centroids xy

Part 2b: Model construction and simulation

#Extract dem values for each centroid of the voronois
src = rasterio.open('../rst/topoWgs12N.tif')  ## Org
elevation=[x for x in src.sample(centroids)] ## Org
nlay = 15 ## Org

mtop=np.array([elev[0] for i,elev in enumerate(elevation)]) ## Org
zbot=np.zeros((nlay,ncpl)) ## Org


AcuifInf_Bottom = 1300 ## Org
zbot[0,] = AcuifInf_Bottom + (0.95 * (mtop - AcuifInf_Bottom)) ## Org
zbot[1,] = AcuifInf_Bottom + (0.90 * (mtop - AcuifInf_Bottom)) ## Org
zbot[2,] = AcuifInf_Bottom + (0.85 * (mtop - AcuifInf_Bottom)) ## Org
zbot[3,] = AcuifInf_Bottom + (0.80 * (mtop - AcuifInf_Bottom)) ## Org
zbot[4,] = AcuifInf_Bottom + (0.75 * (mtop - AcuifInf_Bottom)) ## Org
zbot[5,] = AcuifInf_Bottom + (0.70 * (mtop - AcuifInf_Bottom)) ## Org
zbot[6,] = AcuifInf_Bottom + (0.65 * (mtop - AcuifInf_Bottom)) ## Org
zbot[7,] = AcuifInf_Bottom + (0.60 * (mtop - AcuifInf_Bottom)) ## Org
zbot[8,] = AcuifInf_Bottom + (0.55 * (mtop - AcuifInf_Bottom)) ## Org
zbot[9,] = AcuifInf_Bottom + (0.50 * (mtop - AcuifInf_Bottom)) ## Org
zbot[10,] = AcuifInf_Bottom + (0.45 * (mtop - AcuifInf_Bottom)) ## Org
zbot[11,] = AcuifInf_Bottom + (0.40 * (mtop - AcuifInf_Bottom)) ## Org
zbot[12,] = AcuifInf_Bottom + (0.35 * (mtop - AcuifInf_Bottom)) ## Org
zbot[13,] = AcuifInf_Bottom + (0.2 * (mtop - AcuifInf_Bottom)) ## Org
zbot[14,] = AcuifInf_Bottom ## Org

Create simulation and model

# create simulation
simName = 'mf6Sim' ## Org
modelName = 'mf6Model' ## Org
modelWs = '../modelFiles' ## Org
sim = flopy.mf6.MFSimulation(sim_name=modelName, version='mf6', ## Org
                             exe_name='../bin/mf6.exe', ## Org
                             continue_=True,
                             sim_ws=modelWs) ## Org
# create tdis package
tdis_rc = [(1.0, 1, 1.0)] + [(86400*365*4, 1, 1.0) for level in range(5)] ## Org
print(tdis_rc[:3]) ## Org

tdis = flopy.mf6.ModflowTdis(sim, pname='tdis', time_units='SECONDS', ## Org
                             perioddata=tdis_rc, ## Org
                            nper=6) ## Org
[(1.0, 1, 1.0), (126144000, 1, 1.0), (126144000, 1, 1.0)]
# create gwf model
gwf = flopy.mf6.ModflowGwf(sim, ## Org
                           modelname=modelName, ## Org
                           save_flows=True, ## Org
                           newtonoptions="NEWTON UNDER_RELAXATION") ## Org
# create iterative model solution and register the gwf model with it
ims = flopy.mf6.ModflowIms(sim, ## Org
                           complexity='COMPLEX', ## Org
                           outer_maximum=150, ## Org
                           inner_maximum=50, ## Org
                           outer_dvclose=0.1, ## Org
                           inner_dvclose=0.0001, ## Org
                           backtracking_number=20, ## Org
                           linear_acceleration='BICGSTAB') ## Org
sim.register_ims_package(ims,[modelName]) ## Org
# disv
disv = flopy.mf6.ModflowGwfdisv(gwf, nlay=nlay, ncpl=ncpl, ## Org
                                top=mtop, botm=zbot, ## Org
                                nvert=nvert, vertices=vertices, ## Org
                                cell2d=cell2d) ## Org
disv.top.plot(figsize=(12,8), alpha=0.8) ## Org
crossSection = gpd.read_file('../shp/crossSection.shp') ## Org
sectionLine =list(crossSection.iloc[0].geometry.coords) ## Org

fig, ax = plt.subplots(figsize=(12,8)) ## Org
modelxsect = flopy.plot.PlotCrossSection(model=gwf, line={'Line': sectionLine}) ## Org
linecollection = modelxsect.plot_grid(lw=0.5) ## Org
ax.grid() ## Org
# initial conditions

ic = flopy.mf6.ModflowGwfic(gwf, strt=np.stack([mtop for i in range(nlay)])) ## Org
#headsInitial = np.load('../npy/headCalibInitial.npy')
#ic = flopy.mf6.ModflowGwfic(gwf, strt=headsInitial)
ncplOnes = np.ones(ncpl) # <===== Inserted
Kx =[4E-4] + [7E-7 for x in range(6)] + [3E-7 for x in range(6)] + [1E-7 for x in range(nlay - 13)] ## Org
Kx = [kValue*ncplOnes for kValue in Kx] # <===== Inserted
Kx = np.hstack(Kx) # <===== Inserted

# Kx =[4E-4] + [7E-7 for x in range(4)] + [3E-7 for x in range(4)] + [1E-7 for x in range(nlay - 9)] ## Org
icelltype = [1 for x in range(10)] + [0 for x in range(nlay - 10)] ## Org


# Define intersection object # <===== Inserted
interIx = flopy.utils.gridintersect.GridIntersect(gwf.modelgrid) # <===== Inserted

faultDf = gpd.read_file('../shp/faults_v2.shp') # <===== Inserted

# for faults <==== Inserted
for index, row in faultDf.iterrows(): # <===== Inserted
    cellids = interIx.intersect(row.geometry).cellids # <===== Inserted
    cellidsListLy1_end = [list(cellids+ncpl*lay) for lay in range(1,gwf.modelgrid.nlay)] # <===== Inserted
    Kx[cellidsListLy1_end] = 9e-7 # <===== Inserted

# node property flow
npf = flopy.mf6.ModflowGwfnpf(gwf, ## Org
                              save_specific_discharge=True, ## Org
                              icelltype=icelltype, ## Org
                              k=Kx, ## Org
                            #   k33=np.array(Kx)/10) ## Org
                              k33 = Kx) ## Org
#plot cross section
import matplotlib.colors as mcolors # <===== Inserted
crossSection = gpd.read_file('../shp/crossSection.shp') # <===== Inserted
sectionLine =list(crossSection.iloc[0].geometry.coords) # <===== Inserted

fig, ax = plt.subplots(figsize=(12,8)) # <===== Inserted
modelxsect = flopy.plot.PlotCrossSection(model=gwf, line={'Line': sectionLine}) # <===== Inserted
linecollection = modelxsect.plot_array(npf.k.array, 
                                       alpha=0.5, 
                                       norm=mcolors.LogNorm(vmin=Kx.min(), 
                                                            vmax=Kx.max()))  # <===== Inserted
modelxsect.plot_grid(lw=0.5) # <===== Inserted 

plt.colorbar(linecollection, shrink=0.75) # <===== Inserted
# define storage and transient stress periods
sto = flopy.mf6.ModflowGwfsto(gwf, ## Org
                              iconvert=1, ## Org
                              steady_state={ ## Org
                                0:True, ## Org
                              },
                              transient={
                                  1:True, ## Org
                                  2:True, ## Org
                                  3:True, ## Org
                                  4:True, ## Org
                                  5:True, ## Org
                              },
                              ss=1e-06,
                              sy=0.001,
                              ) ## Org

Working with rechage, evapotranspiration

rchr = 0.2/365/86400 ## Org
rch = flopy.mf6.ModflowGwfrcha(gwf, recharge=rchr) ## Org
evtr = 1.2/365/86400 ## Org
evt = flopy.mf6.ModflowGwfevta(gwf,ievt=1,surface=mtop,rate=evtr,depth=1.0) ## Org

Definition of the intersect object

For the manipulation of spatial data to determine hydraulic parameters or boundary conditions

# Define intersection object
interIx = flopy.utils.gridintersect.GridIntersect(gwf.modelgrid) ## Org
#river package
layCellTupleList, cellElevList = getLayCellElevTupleFromRaster(gwf,interIx,'../rst/topoWgs12N.tif',
                                                               '../shp/river_basin.shp') ## Org
riverSpd = {} ## Org
riverSpd[0] = [] ## Org
for index, layCellTuple in enumerate(layCellTupleList): ## Org
    riverSpd[0].append([layCellTuple,cellElevList[index],0.01,'Drain']) ## Org
import copy ## Org
#for yr 4 
riverSpd[1] = copy.copy(riverSpd[0]) ## Org
layCellTupleList, cellElevList = getLayCellElevTupleFromRaster(gwf, ## Org
                                                               interIx, ## Org
                                                               '../rst/minePlan/pitElevYr04.tif', ## Org
                                                               '../shp/minePlan/pitShell.shp') ## Org
for index, layCellTuple in enumerate(layCellTupleList): ## Org
    riverSpd[1].append([layCellTuple,cellElevList[index],0.011, 'Pit']) # <===== Inserted

#for yr 8 
riverSpd[2] = copy.copy(riverSpd[0]) ## Org
layCellTupleList, cellElevList = getLayCellElevTupleFromRaster(gwf, ## Org
                                                               interIx, ## Org
                                                               '../rst/minePlan/pitElevYr08.tif', ## Org
                                                               '../shp/minePlan/pitShell.shp') ## Org
for index, layCellTuple in enumerate(layCellTupleList): ## Org
    riverSpd[2].append([layCellTuple,cellElevList[index],0.012, 'Pit']) # <===== Inserted

#for yr 12 
riverSpd[3] = copy.copy(riverSpd[0]) ## Org
layCellTupleList, cellElevList = getLayCellElevTupleFromRaster(gwf, ## Org
                                                               interIx, ## Org
                                                               '../rst/minePlan/pitElevYr12.tif', ## Org
                                                               '../shp/minePlan/pitShell.shp') ## Org
for index, layCellTuple in enumerate(layCellTupleList): ## Org
    riverSpd[3].append([layCellTuple,cellElevList[index],0.013, 'Pit']) # <===== Inserted

#for yr 16 
riverSpd[4] = copy.copy(riverSpd[0]) ## Org
layCellTupleList, cellElevList = getLayCellElevTupleFromRaster(gwf, ## Org
                                                               interIx, ## Org
                                                               '../rst/minePlan/pitElevYr16.tif', ## Org
                                                               '../shp/minePlan/pitShell.shp') ## Org
for index, layCellTuple in enumerate(layCellTupleList): ## Org
    riverSpd[4].append([layCellTuple,cellElevList[index],0.014, 'Pit']) # <===== Inserted 

#for yr 20 
riverSpd[5] = copy.copy(riverSpd[0]) ## Org
layCellTupleList, cellElevList = getLayCellElevTupleFromRaster(gwf, ## Org
                                                               interIx, ## Org
                                                               '../rst/minePlan/pitElevYr20.tif', ## Org
                                                               '../shp/minePlan/pitShell.shp') ## Org
for index, layCellTuple in enumerate(layCellTupleList): ## Org
    riverSpd[5].append([layCellTuple,cellElevList[index],0.015, 'Pit']) # <===== Inserted
The cell 2385 has a elevation of 1897.48 outside the model vertical domain
The cell 2395 has a elevation of 1898.64 outside the model vertical domain
The cell 2402 has a elevation of 1898.52 outside the model vertical domain
The cell 2403 has a elevation of 1898.80 outside the model vertical domain
The cell 2406 has a elevation of 1898.12 outside the model vertical domain
The cell 2407 has a elevation of 1898.32 outside the model vertical domain
drn = flopy.mf6.ModflowGwfdrn(gwf, stress_period_data=riverSpd, boundnames=True) # <===== Inserted 

# Observation package for Drain
obsDict = { # <===== Inserted 
    "{}.drn.obs.csv".format(modelName): [ # <===== Inserted 
        ("drain_flow", "drn", "Drain"), # <===== Inserted 
        ("pit_flow", "drn", "Pit") # <===== Inserted 
    ] # <===== Inserted 
} # <===== Inserted 

# Attach observation package to DRN package
drn.obs.initialize( # <===== Inserted 
    filename=gwf.name+".drn.obs", # <===== Inserted 
    digits=10, # <===== Inserted 
    print_input=True, # <===== Inserted 
    continuous=obsDict # <===== Inserted 
) # <===== Inserted
#river plot
drn.plot(mflay=0, kper=1) # <===== Inserted
crossSection = gpd.read_file('../shp/crossSection.shp') ## Org
sectionLine =list(crossSection.iloc[0].geometry.coords) ## Org

fig, ax = plt.subplots(figsize=(12,8)) ## Org
xsect = flopy.plot.PlotCrossSection(model=gwf, line={'Line': sectionLine}) ## Org
lc = xsect.plot_grid(lw=0.5) ## Org
xsect.plot_bc('DRN',kper=5) ## Org
ax.grid() ## Org
#regional flow package
layCellTupleList, elevList = getLayCellElevTupleFromRaster(gwf,interIx,'../rst/waterTable_v2.tif','../shp/local/regionalFlow.shp') ## <=== updated
ghbSpd = {} ## Org
ghbSpd[0] = [] ## Org
for index, layCellTuple in enumerate(layCellTupleList): ## <=== updated
    ghbSpd[0].append([layCellTuple,elevList[index],0.01]) ## <=== updated
ghbSpd[0][:5] ## <=== updated
The cell 184 has a elevation of 1769.59 outside the model vertical domain





[[(0, 0), 1778.0935269188267, 0.01],
 [(1, 3), 1783.0812529535838, 0.01],
 [(0, 6), 1769.4535360685973, 0.01],
 [(1, 17), 1808.8903920767075, 0.01],
 [(0, 18), 1779.0706182680847, 0.01]]
ghb = flopy.mf6.ModflowGwfghb(gwf, stress_period_data=ghbSpd) ## <===== modified
#regional flow plot
ghb.plot(mflay=0, kper=0) ## <===== modified

Set the Output Control and run simulation

#oc
head_filerecord = f"{gwf.name}.hds" ## Org
budget_filerecord = f"{gwf.name}.cbc" ## Org
oc = flopy.mf6.ModflowGwfoc(gwf, ## Org
                            head_filerecord=head_filerecord, ## Org
                            budget_filerecord = budget_filerecord, ## Org
                            saverecord=[("HEAD", "LAST"),("BUDGET","LAST")]) ## Org
# Run the simulation
sim.write_simulation() ## Org
success, buff = sim.run_simulation() ## Org
writing simulation...
  writing simulation name file...
  writing simulation tdis package...
  writing solution package ims_-1...
  writing model mf6Model...
    writing model name file...
    writing package disv...
    writing package ic...
    writing package npf...
    writing package sto...
    writing package rcha_0...
    writing package evta_0...
    writing package drn_0...
INFORMATION: maxbound in ('gwf6', 'drn', 'dimensions') changed to 29349 based on size of stress_period_data
    writing package obs_0...
    writing package ghb_0...
INFORMATION: maxbound in ('gwf6', 'ghb', 'dimensions') changed to 462 based on size of stress_period_data
    writing package oc...
FloPy is using the following executable to run the model: ..\bin\mf6.exe
                                   MODFLOW 6
                U.S. GEOLOGICAL SURVEY MODULAR HYDROLOGIC MODEL
                            VERSION 6.6.2 05/12/2025

   MODFLOW 6 compiled May 12 2025 12:42:18 with Intel(R) Fortran Intel(R) 64
   Compiler Classic for applications running on Intel(R) 64, Version 2021.7.0
                             Build 20220726_000000

This software has been approved for release by the U.S. Geological 
Survey (USGS). Although the software has been subjected to rigorous 
review, the USGS reserves the right to update the software as needed 
pursuant to further analysis and review. No warranty, expressed or 
implied, is made by the USGS or the U.S. Government as to the 
functionality of the software and related material nor shall the 
fact of release constitute any such warranty. Furthermore, the 
software is released on condition that neither the USGS nor the U.S. 
Government shall be held liable for any damages resulting from its 
authorized or unauthorized use. Also refer to the USGS Water 
Resources Software User Rights Notice for complete use, copyright, 
and distribution information.

 
 MODFLOW runs in SEQUENTIAL mode
 
 Run start date and time (yyyy/mm/dd hh:mm:ss): 2025/09/03 10:28:36
 
 Writing simulation list file: mfsim.lst
 Using Simulation name file: mfsim.nam
 
    Solving:  Stress period:     1    Time step:     1
    Solving:  Stress period:     2    Time step:     1
    Solving:  Stress period:     3    Time step:     1
    Solving:  Stress period:     4    Time step:     1
    Solving:  Stress period:     5    Time step:     1
    Solving:  Stress period:     6    Time step:     1
 
 Run end date and time (yyyy/mm/dd hh:mm:ss): 2025/09/03 10:49:00
 Elapsed run time: 20 Minutes, 23.857 Seconds
 

WARNING REPORT:

  1.  Simulation convergence failure occurred 1 time(s).
 Normal termination of simulation.

Model output visualization

headObj = gwf.output.head() ## Org
headObj.get_kstpkper() ## Org
[(0, 0), (0, 1), (0, 2), (0, 3), (0, 4), (0, 5)]
kper = 5 ## Org
lay = 0 ## Org
heads = headObj.get_data(kstpkper=(0,kper)) 
#heads[lay,0,:5] 
#heads = headObj.get_data(kstpkper=(0,0)) 
#np.save('../npy/headCalibInitial', heads)
### Plot the heads for a defined layer and boundary conditions
fig = plt.figure(figsize=(12,8)) ## Org
ax = fig.add_subplot(1, 1, 1, aspect='equal') ## Org
modelmap = flopy.plot.PlotMapView(model=gwf) ## Org

####
levels = np.linspace(heads[heads>-1e+30].min(),heads[heads>-1e+30].max(),num=50) ## Org
contour = modelmap.contour_array(heads[lay],ax=ax,levels=levels,cmap='PuBu') 
ax.clabel(contour) ## Org


quadmesh = modelmap.plot_bc('DRN') ## Org
cellhead = modelmap.plot_array(heads[lay],ax=ax, cmap='Blues', alpha=0.8) 

linecollection = modelmap.plot_grid(linewidth=0.3, alpha=0.5, color='cyan', ax=ax) ## Org

plt.colorbar(cellhead, shrink=0.75) ## Org

plt.show() ## Org
crossSection = gpd.read_file('../shp/crossSection.shp')
sectionLine =list(crossSection.iloc[0].geometry.coords)

waterTable = flopy.utils.postprocessing.get_water_table(heads)
#waterTable2 = flopy.utils.postprocessing.get_water_table(heads2)

fig, ax = plt.subplots(figsize=(12,8))
xsect = flopy.plot.PlotCrossSection(model=gwf, line={'Line': sectionLine})
lc = modelxsect.plot_grid(lw=0.5)
xsect.plot_array(heads, alpha=0.5)
xsect.plot_surface(waterTable)
xsect.plot_bc('drn', kper=kper, facecolor='none', edgecolor='teal')
plt.show()

#Vtk generation
import flopy ## Org
from mf6Voronoi.tools.vtkGen import Mf6VtkGenerator ## Org
from mf6Voronoi.utils import initiateOutputFolder ## Org
C:\Users\saulm\anaconda3\Lib\site-packages\pyvista\examples\downloads.py:93: DeprecationWarning: support for supplying keyword arguments to pathlib.PurePath is deprecated and scheduled for removal in Python 3.14
  Path(USER_DATA_PATH, exist_ok=True).mkdir()
C:\Users\saulm\anaconda3\Lib\site-packages\pyvista\examples\downloads.py:98: UserWarning: Unable to access C:\Users\saulm\AppData\Local\pyvista_3\pyvista_3\Cache. Manually specify the PyVistaexamples cache with the PYVISTA_USERDATA_PATH environment variable.
  warnings.warn(
# load simulation
simName = 'mf6Sim' ## Org
modelName = 'mf6Model' ## Org
modelWs = '../modelFiles' ## Org
sim = flopy.mf6.MFSimulation.load(sim_name=modelName, version='mf6', ## Org
                             exe_name='bin/mf6.exe', ## Org
                             sim_ws=modelWs) ## Org
loading simulation...
  loading simulation name file...
  loading tdis package...
  loading model gwf6...
    loading package disv...
    loading package ic...
    loading package npf...
    loading package sto...
    loading package rch...
    loading package evt...
    loading package drn...
    loading package ghb...
    loading package oc...
  loading solution package mf6model...
vtkDir = '../vtk' ## Org
initiateOutputFolder(vtkDir) ## Org

mf6Vtk = Mf6VtkGenerator(sim, vtkDir) ## Org
The output folder ../vtk has been generated.
Hatarilabs

mf6Voronoi will have a web version in 2028

Follow us:

Hatarilabs Hatarilabs Hatarilabs Hatarilabs Hatarilabs Hatarilabs
/---------------------------------------/

The Vtk generator engine has been started

/---------------------------------------/
#list models on the simulation
mf6Vtk.listModels() ## Org
Models in simulation: ['mf6model']
mf6Vtk.loadModel(modelName) ## Org
Package list: ['DISV', 'IC', 'NPF', 'STO', 'RCHA_0', 'EVTA_0', 'DRN_OBS', 'DRN_0', 'GHB_0', 'OC']
#show output data
headObj = mf6Vtk.gwf.output.head() ## Org
headObj.get_kstpkper() ## Org
[(0, 0), (0, 1), (0, 2), (0, 3), (0, 4), (0, 5)]
#generate model geometry as vtk and parameter array
mf6Vtk.generateGeometryArrays() ## Org
#generate parameter vtk
mf6Vtk.generateParamVtk() ## Org
Parameter Vtk Generated
#generate bc and obs vtk
mf6Vtk.generateBcObsVtk(nper=5) ## Org
/--------RCHA_0 vtk generation-------/
Working for RCHA_0 package, creating the datasets: dict_keys(['irch', 'recharge', 'aux'])
[WARNING] There is no data for the required stress period
Vtk file took 0.1722 seconds to be generated.
/--------RCHA_0 vtk generated-------/


/--------EVTA_0 vtk generation-------/
Working for EVTA_0 package, creating the datasets: dict_keys(['ievt', 'surface', 'rate', 'depth', 'aux'])
[WARNING] There is no data for the required stress period
Vtk file took 0.1562 seconds to be generated.
/--------EVTA_0 vtk generated-------/


/--------DRN_0 vtk generation-------/
Working for DRN_0 package, creating the datasets: ('elev', 'cond', 'boundname')
Vtk file took 743.5004 seconds to be generated.
/--------DRN_0 vtk generated-------/


/--------GHB_0 vtk generation-------/
Working for GHB_0 package, creating the datasets: ('bhead', 'cond')
[WARNING] There is no data for the required stress period
Vtk file took 0.4015 seconds to be generated.
/--------GHB_0 vtk generated-------/
mf6Vtk.generateHeadVtk(nper=5, crop=True) ## Org
mf6Vtk.generateWaterTableVtk(nper=5) ## Org

import pandas as pd
import matplotlib.pyplot as plt
drainDf = pd.read_csv('../modelFiles/mf6Model.drn.obs.csv', na_values=3e30)
drainDf.head()

time DRAIN_FLOW PIT_FLOW
0 1.0 -0.731761 NaN
1 126144001.0 -0.323163 -0.010561
2 252288001.0 -0.294534 -0.046803
3 378432001.0 -0.269111 -0.080022
4 504576001.0 -0.249139 -0.105717
inflowsDf = pd.DataFrame()
inflowsDf.index = drainDf['time']/86400/365
inflowsDf.index = inflowsDf.index.astype(int)
inflowsDf['Pit(l/s)'] =  -1*drainDf['PIT_FLOW'].values*1000
inflowsDf['Rivers(l/s)'] =  -1*drainDf['DRAIN_FLOW'].values*1000
fig, ax = plt.subplots(figsize=(12,8))
ax.plot(inflowsDf.index,inflowsDf['Pit(l/s)'],label='Pit Inflows')
ax.plot(inflowsDf.index,inflowsDf['Rivers(l/s)'],label='Rivers')
ax.legend()
ax.grid()
ax.set_xlabel('Year')
ax.set_ylabel('Flow (l/s)')
plt.show()

Datos de ingreso

Descargue los datos de ingreso desde este link:

owncloud.hatarilabs.com/s/ALDt8WUwRcAdX8h

Password: Hatarilabs

 

Suscríbete a nuestro boletín electrónico

Suscríbase a nuestro boletín gratuito para recibir noticias, datos interesantes y fechas de nuestros cursos en recursos hídricos.

 

Posted on September 3, 2025 and filed under TutorialPython, TutorialModflow, Modelamiento, Minería.