Modelamiento de ingresos en galerías filtrantes con MODFLOW6 y mf6Voronoi - Tutorial

Infiltration galleries are a low cost and low maintenance option for domestic water supply. The amount of inflow water and the interaction with the water bodies are main concerns on the evaluation and design of infiltration galleries. We have done an applied case of inflow simulation to infiltration galleries with MODFLOW6 based on Voronoi meshes. The example covers all steps from mesh creation, steady state model construction, simulation of infiltration galleries and the inflow calculation per gallery group. Finally a 3D representation of the model geometry, boundary conditions and head distribution is performed on Paraview.

Tutorial

Código


from mf6Voronoi.utils import listTemplates, copyTemplate
#listTemplates()
copyTemplate('generateVoronoi','galeria')
copyTemplate('modelCreation','galeria')
copyTemplate('vtkGeneration','galeria')

Part 1 : Voronoi mesh generation

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='infiltrationGalleries',maxRef = 25, multiplier=1.5) ## Org

#Open limit layers and refinement definition layers
vorMesh.addLimit('basin','../shp/modelLimit.shp') ## Org
vorMesh.addLayer('galleries','../shp/infiltrationGalleries.shp',1) ## Org
vorMesh.addLayer('river','../shp/river.shp',5) ## Org
vorMesh.addLayer('ghb','../shp/regionalFlow.shp',10) ## Org
#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 galleries discretization-------/
Progressive cell size list: [1, 2.5, 4.75, 8.125, 13.1875, 20.78125] m.

/--------Layer river discretization-------/
Progressive cell size list: [5, 12.5, 23.75] m.

/--------Layer ghb discretization-------/
Progressive cell size list: [10, 25.0] m.

/----Sumary of points for voronoi meshing----/
Distributed points from layers: 3
Points from layer buffers: 2291
Points from max refinement areas: 312
Points from min refinement areas: 1042
Total points inside the limit: 4140
/--------------------------------------------/

Time required for point generation: 0.77 seconds 


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

Time required for voronoi generation: 0.57 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: 4139 with length = 0.00080
Short side on polygon: 4139 with length = 0.00080
Short side on polygon: 4139 with length = 0.00316
Short side on polygon: 4139 with length = 0.00316
Short side on polygon: 4139 with length = 0.00985
Short side on polygon: 4139 with length = 0.00985
Short side on polygon: 4139 with length = 0.00685
Short side on polygon: 4139 with length = 0.00685
Short side on polygon: 4139 with length = 0.00610
Short side on polygon: 4139 with length = 0.00610
Short side on polygon: 4139 with length = 0.00492
Short side on polygon: 4139 with length = 0.00492
Short side on polygon: 4139 with length = 0.00440
Short side on polygon: 4139 with length = 0.00440
Short side on polygon: 4139 with length = 0.00440
Short side on polygon: 4139 with length = 0.00440
vorMesh.fixVoronoiShortSides()
vorMesh.generateVoronoi()
vorMesh.checkVoronoiQuality(threshold=0.01)
/----Generation of the voronoi mesh----/

Time required for voronoi generation: 0.49 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: 0.93 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')
# 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%|███████████████████████████████████████████████████████████████████████████| 4148/4148 [00:00<00:00, 18028.78it/s]



Extracting cell2d data and grid index


100%|████████████████████████████████████████████████████████████████████████████| 4148/4148 [00:01<00:00, 2972.32it/s]
#create folder
initiateOutputFolder('../json') 

#export disv
mesh.save_properties('../json/disvDict.json')
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 shapely.geometry import MultiLineString, MultiPolygon ## <==== updated
from mf6Voronoi.tools.cellWork import getLayCellElevTupleFromRaster, getLayCellElevTupleFromElev ## <==== updated
from mf6Voronoi.tools.graphs2d import generateRasterFromArray ## <==== updated
# 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/dem44N_v2.tif')  ## Org
elevation=[x for x in src.sample(centroids)] ## Org
nlay = 10   ## Org

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


AcuifInf_Bottom = 1262 ## Org
zbot[0,] = AcuifInf_Bottom + (0.97 * (mtop - AcuifInf_Bottom)) ## Org
zbot[1,] = AcuifInf_Bottom + (0.94 * (mtop - AcuifInf_Bottom)) ## Org
zbot[2,] = AcuifInf_Bottom + (0.91 * (mtop - AcuifInf_Bottom)) ## Org
zbot[3,] = AcuifInf_Bottom + (0.88 * (mtop - AcuifInf_Bottom)) ## Org
zbot[4,] = AcuifInf_Bottom + (0.82 * (mtop - AcuifInf_Bottom)) ## Org
zbot[5,] = AcuifInf_Bottom + (0.75 * (mtop - AcuifInf_Bottom)) ## Org
zbot[6,] = AcuifInf_Bottom + (0.65 * (mtop - AcuifInf_Bottom)) ## Org
zbot[7,] = AcuifInf_Bottom + (0.45 * (mtop - AcuifInf_Bottom)) ## Org
zbot[8,] = AcuifInf_Bottom + (0.15 * (mtop - AcuifInf_Bottom)) ## Org
zbot[9,] = 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
                             sim_ws=modelWs) ## Org
# create tdis package
tdis_rc = [(1000.0, 1, 1.0)] ## Org
tdis = flopy.mf6.ModflowTdis(sim, pname='tdis', time_units='SECONDS', ## Org
                             perioddata=tdis_rc) ## Org
# 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=30, ## 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
# initial conditions
ic = flopy.mf6.ModflowGwfic(gwf, strt=np.stack([mtop for i in range(nlay)])) ## Org
Kx =[4E-4 for i in range(5)] + [5E-6 for i in range(5)] ## Org
icelltype = [1 for i in range(2)] + [0 for i in range(8)] ## Org

# node property flow
npf = flopy.mf6.ModflowGwfnpf(gwf, ## Org
                              save_specific_discharge=True, ## Org
                              icelltype=icelltype, ## Org
                              k=Kx,
                              k33=np.array(Kx)/10) ## Org
# define storage and transient stress periods
sto = flopy.mf6.ModflowGwfsto(gwf, ## Org
                              iconvert=1, ## Org
                              steady_state={ ## Org
                                0:True, ## Org
                              } ## Org
                              ) ## Org

Working with rechage, evapotranspiration

rchr = 0.15/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
#open the river shapefile
rivers =gpd.read_file('../shp/river.shp') ## Org
list_rivers=[] ## Org
for i in range(rivers.shape[0]): ## Org
    list_rivers.append(rivers['geometry'].loc[i]) ## Org
    
riverMls = MultiPolygon(polygons=list_rivers) ## Org

#intersec rivers with our grid
riverCells=interIx.intersect(riverMls).cellids ## Org
riverCells[:10] ## Org
array([4, 90, 132, 133, 134, 136, 141, 142, 144, 186], dtype=object)
#river package
riverSpd = {} ## Org
riverSpd[0] = [] ## Org
for cell in riverCells: ## Org
    riverSpd[0].append([(0,cell),mtop[cell]-0.3,0.01,mtop[cell]-0.6]) ## Org
riv = flopy.mf6.ModflowGwfriv(gwf, stress_period_data=riverSpd) ## Org
#river plot
riv.plot(mflay=0) ## Org
#regional flow package
layCellTupleList, elevList = getLayCellElevTupleFromRaster(gwf,interIx,'../rst/waterTable.tif','../shp/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 854 has a elevation of 1294.16 outside the model vertical domain
The cell 889 has a elevation of 1294.20 outside the model vertical domain
The cell 919 has a elevation of 1294.22 outside the model vertical domain
The cell 949 has a elevation of 1294.24 outside the model vertical domain
The cell 1010 has a elevation of 1294.27 outside the model vertical domain
The cell 1059 has a elevation of 1294.31 outside the model vertical domain
The cell 1101 has a elevation of 1294.30 outside the model vertical domain
The cell 1130 has a elevation of 1294.17 outside the model vertical domain
The cell 1217 has a elevation of 1294.12 outside the model vertical domain
The cell 1265 has a elevation of 1294.16 outside the model vertical domain
The cell 1308 has a elevation of 1294.11 outside the model vertical domain
The cell 2682 has a elevation of 1294.09 outside the model vertical domain
The cell 2991 has a elevation of 1294.09 outside the model vertical domain
The cell 3103 has a elevation of 1294.15 outside the model vertical domain
The cell 3141 has a elevation of 1294.16 outside the model vertical domain
The cell 3165 has a elevation of 1294.12 outside the model vertical domain
The cell 3232 has a elevation of 1294.13 outside the model vertical domain
The cell 3260 has a elevation of 1294.18 outside the model vertical domain
The cell 3293 has a elevation of 1294.19 outside the model vertical domain
The cell 3321 has a elevation of 1294.15 outside the model vertical domain
The cell 3377 has a elevation of 1294.16 outside the model vertical domain





[[(2, 83), np.float64(1293.4989), 0.01],
 [(2, 135), np.float64(1293.5421), 0.01],
 [(2, 188), np.float64(1293.6545), 0.01],
 [(2, 217), np.float64(1293.7581), 0.01],
 [(2, 247), np.float64(1293.7969), 0.01]]
ghb = flopy.mf6.ModflowGwfghb(gwf, stress_period_data=ghbSpd) ## <===== modified
#regional flow plot
ghb.plot(mflay=0, kper=0) ## <===== modified
#drain package
drnSpd = {} ## Org
drnSpd[0] = [] ## Org
#gallery1
layCellTupleList = getLayCellElevTupleFromElev(gwf,interIx,1290,'../shp/infiltrationGallery1.shp') ## <=== updated
for index, layCellTuple in enumerate(layCellTupleList): ## <=== updated
    drnSpd[0].append([layCellTuple,1290,0.01, 'gallery1']) ## <=== updated
#gallery2
layCellTupleList = getLayCellElevTupleFromElev(gwf,interIx,1290,'../shp/infiltrationGallery2.shp') ## <=== updated
for index, layCellTuple in enumerate(layCellTupleList): ## <=== updated
    drnSpd[0].append([layCellTuple,1290,0.01, 'gallery2']) ## <=== updated
drnSpd[0][:5] ## <=== updated
You have inserted a fixed elevation
You have inserted a fixed elevation





[[(4, 1946), 1290, 0.01, 'gallery1'],
 [(4, 1957), 1290, 0.01, 'gallery1'],
 [(4, 1970), 1290, 0.01, 'gallery1'],
 [(4, 1971), 1290, 0.01, 'gallery1'],
 [(4, 1973), 1290, 0.01, 'gallery1']]
drn = flopy.mf6.ModflowGwfdrn(gwf, stress_period_data=drnSpd, boundnames=True) ## Org

drnDict = { # <===== Inserted 
    "{}.drn.obs.csv".format(modelName): [ # <===== Inserted 
        ("gallery1", "drn", "gallery1"), # <===== Inserted 
        ("gallery2", "drn", "gallery2") 
    ] # <===== 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=drnDict # <===== Inserted 
) # <===== Inserted
#gallery plot
drn.plot(mflay=4, kper=0) ## Org
crossSection = gpd.read_file('../shp/crossSection.shp') ## <=== updated
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
modelxsect.plot_bc('drn', color='crimson')
modelxsect.plot_bc('riv')
modelxsect.plot_grid(lw=0.5) ## Org

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 riv_0...
INFORMATION: maxbound in ('gwf6', 'riv', 'dimensions') changed to 1958 based on size of stress_period_data
    writing package ghb_0...
INFORMATION: maxbound in ('gwf6', 'ghb', 'dimensions') changed to 194 based on size of stress_period_data
    writing package drn_0...
INFORMATION: maxbound in ('gwf6', 'drn', 'dimensions') changed to 189 based on size of stress_period_data
    writing package obs_0...
    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.0 12/19/2024

   MODFLOW 6 compiled Dec 19 2024 21:59:25 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/08/26 17:16:10
 
 Writing simulation list file: mfsim.lst
 Using Simulation name file: mfsim.nam
 
    Solving:  Stress period:     1    Time step:     1
 
 Run end date and time (yyyy/mm/dd hh:mm:ss): 2025/08/26 17:16:19
 Elapsed run time:  9.566 Seconds
 
 Normal termination of simulation.

Model output visualization

headObj = gwf.output.head() ## Org
headObj.get_kstpkper() ## Org
[(np.int32(0), np.int32(0))]
heads = headObj.get_data() ## Org
heads[2,0,:5] ## Org
array([1292.2563048 , 1291.78482183, 1292.32213598, 1292.4317531 ,
       1292.61689701])
# 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=10) ## Org
contour = modelmap.contour_array(heads[3],ax=ax,levels=levels,cmap='PuBu') ## Org
ax.clabel(contour) ## Org


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

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
#plot cross section
crossSection = gpd.read_file('../shp/crossSection.shp') ## <==== inserted
sectionLine =list(crossSection.iloc[0].geometry.coords) ## <==== inserted

waterTable = flopy.utils.postprocessing.get_water_table(heads)

fig = plt.figure(figsize=(18, 5))
ax = fig.add_subplot(1, 1, 1)
ax.set_title("plot_array()")

xsect = flopy.plot.PlotCrossSection(model=gwf, line={"line": sectionLine})
patch_collection = xsect.plot_array(heads, alpha=0.5)
xsect.plot_surface(waterTable)
line_collection = xsect.plot_grid()
cb = plt.colorbar(patch_collection, shrink=0.75)
waterTable = flopy.utils.postprocessing.get_water_table(heads)
generateRasterFromArray(gwf, 
                        waterTable, 
                        rasterRes=2, 
                        epsg=32644, 
                        outputPath='../output/waterTable.tif', 
                        limitLayer='../shp/modelLimit.shp')
Raster X Dim: 810.00, Raster Y Dim: 430.00
Number of cols:  406, Number of rows: 216

#Vtk generation
import flopy ## Org
from mf6Voronoi.tools.vtkGen import Mf6VtkGenerator ## Org
from mf6Voronoi.utils import initiateOutputFolder ## Org
# 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 riv...
    loading package ghb...
    loading package drn...
    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', 'RIV_0', 'GHB_0', 'DRN_OBS', 'DRN_0', 'OC']
#show output data
headObj = mf6Vtk.gwf.output.head() ## Org
headObj.get_kstpkper() ## Org
[(np.int32(0), np.int32(0))]
#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=0) ## Org
/--------RCHA_0 vtk generation-------/
Working for RCHA_0 package, creating the datasets: dict_keys(['irch', 'recharge', 'aux'])
Vtk file took 0.0535 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'])
Vtk file took 0.0621 seconds to be generated.
/--------EVTA_0 vtk generated-------/


/--------RIV_0 vtk generation-------/
Working for RIV_0 package, creating the datasets: ('stage', 'cond', 'rbot')
Vtk file took 6.4831 seconds to be generated.
/--------RIV_0 vtk generated-------/


/--------GHB_0 vtk generation-------/
Working for GHB_0 package, creating the datasets: ('bhead', 'cond')
Vtk file took 0.6133 seconds to be generated.
/--------GHB_0 vtk generated-------/


/--------DRN_0 vtk generation-------/
Working for DRN_0 package, creating the datasets: ('elev', 'cond', 'boundname')
Vtk file took 0.5896 seconds to be generated.
/--------DRN_0 vtk generated-------/
mf6Vtk.generateHeadVtk(nper=0, crop=True) ## Org
mf6Vtk.generateWaterTableVtk(nper=0) ## Org

Datos de ingreso

Puedes descargar los datos de ingreso desde este link:

owncloud.hatarilabs.com/s/T6tSBaXrCZvloUQ

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 August 27, 2025 and filed under TutorialModflow, Modelamiento, Hidrogeología.