Transient simulation of the trench excavation impact on a coastal aquifer with mf6Voronoi - Tutorial
/Local scale modeling of a 4 stage trench excavation on a coastal aquifer with mf6Voronoi. The main objective of this modeling work is assessing the impact of the excavation on the nearby groundwater flow regime.
Tutorial
Code
#!pip install -U mf6Voronoi
from mf6Voronoi.utils import listTemplates, copyTemplate
#listTemplates()
copyTemplate('generateVoronoi','trench')
copyTemplate('multilayeredTransient','trench')
copyTemplate('vtkGeneration','trench')
Mesh generation
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='trenchExcavation',maxRef = 50, multiplier=2.5) ## Org
#Open limit layers and refinement definition layers
vorMesh.addLimit('basin','../shp/modelAoi.shp') ## <=== update
vorMesh.addLayer('trench','../shp/trenchExcavationDissolved.shp',1) ## <=== update
vorMesh.addLayer('ghb','../shp/compoundGhb.shp',10) ## <=== update
vorMesh.addLayer('wells','../shp/pumpingWells.shp',2) ## <=== update
#Generate point pair array
vorMesh.generateOrgDistVertices() ## Org
#Generate the point cloud and voronoi
vorMesh.createPointCloud() ## Org
vorMesh.generateVoronoi() ## Org
Follow us: |
|
|
|
|
|
|
/--------Layer trench discretization-------/
Progressive cell size list: [1, 3.5, 9.75, 25.375] m.
/--------Layer ghb discretization-------/
Progressive cell size list: [10, 35.0] m.
/--------Layer wells discretization-------/
Progressive cell size list: [2, 7.0, 19.5] m.
/----Sumary of points for voronoi meshing----/
Distributed points from layers: 3
Points from layer buffers: 2427
Points from max refinement areas: 1121
Points from min refinement areas: 1011
Total points inside the limit: 5234
/--------------------------------------------/
Time required for point generation: 0.65 seconds
/----Generation of the voronoi mesh----/
Time required for voronoi generation: 0.46 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)
#vorMesh.fixVoronoiShortSides()
#vorMesh.generateVoronoi()
#vorMesh.checkVoronoiQuality(threshold=0.01)
#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: 1.33 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%|███████████████████████████████████████████████████████████████████████████| 5234/5234 [00:00<00:00, 19900.33it/s]
Extracting cell2d data and grid index
100%|████████████████████████████████████████████████████████████████████████████| 5234/5234 [00:01<00:00, 3461.15it/s]
#create folder
initiateOutputFolder('../json') ## Org
#export disv
mesh.save_properties('../json/disvDict.json') ## Org
The output folder ../json exists and has been cleared
Multilayer and transient model
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
from mf6Voronoi.tools.graphs2d import generateRasterFromArray, generateContoursFromRaster
# 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/elevDem.tif') ## Org
elevation=[x for x in src.sample(centroids)] ## Org
nlay = 5 ## Org
mtop=np.array([elev[0] for i,elev in enumerate(elevation)]) ## Org
zbot=np.zeros((nlay,ncpl)) ## Org
AcuifInf_Bottom = -20 ## Org
zbot[0,] = AcuifInf_Bottom + (0.8 * (mtop - AcuifInf_Bottom)) ## Org
zbot[1,] = AcuifInf_Bottom + (0.6 * (mtop - AcuifInf_Bottom)) ## Org
zbot[2,] = AcuifInf_Bottom + (0.4 * (mtop - AcuifInf_Bottom)) ## Org
zbot[3,] = AcuifInf_Bottom + (0.2 * (mtop - AcuifInf_Bottom)) ## Org
zbot[4,] = 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 = [(1.0, 1, 1.0)] + [(86400, 1, 1.0) for level in range(4)] ## Org
print(tdis_rc[:3]) ## Org
tdis = flopy.mf6.ModflowTdis(sim, pname='tdis', time_units='SECONDS', ## Org
perioddata=tdis_rc, ## Org
nper=5) ## Org
[(1.0, 1, 1.0), (86400, 1, 1.0), (86400, 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)
Kx =[7E-4 for x in range(5)] ## <=== updated
icelltype = [1, 1, 1, 0, 0] ## Org
# node property flow
npf = flopy.mf6.ModflowGwfnpf(gwf, ## Org
save_specific_discharge=True, ## Org
icelltype=icelltype, ## Org
k=Kx, ## Org
k33=np.array(Kx)/2) ## Org
# 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
},
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
# regional flow as ghb ## <=== updated
layCellTupleList = getLayCellElevTupleFromElev(gwf,interIx,-1,'../shp/compoundGhb.shp') #elev of -1 to get all cells
ghbSpd = {} ## <=== updated
ghbSpd[0] = [] ## <=== updated
for index, layCellTuple in enumerate(layCellTupleList): ## <=== updated
ghbSpd[0].append([layCellTuple,0,0.01]) ## <=== updated
ghbSpd[0][:5] ## <=== updated
You have inserted a fixed elevation
[[(0, 4), 0, 0.01],
[(0, 9), 0, 0.01],
[(0, 10), 0, 0.01],
[(0, 21), 0, 0.01],
[(0, 26), 0, 0.01]]
ghb = flopy.mf6.ModflowGwfghb(gwf, stress_period_data=ghbSpd) ## <===== modified
#regional flow plot
ghb.plot(mflay=0, kper=0) ## <===== modified
# trench as drain package from stress period 1 flow as ghb ## <=== updated
drainSpd = {} ## Org
drainSpd[0] = [] ## Org
drainDf = gpd.read_file('../shp/trenchExcavationDissolved.shp')
tempTrench = drainDf.iloc[1].geometry
tempTrench
i = 1
for index, row in drainDf.iterrows():
tempTrench = row.geometry
tempLayCellTupleList = getLayCellElevTupleFromElev(gwf,interIx,-4.5,tempTrench) #elev of -1 to get all cells
drainSpd[i] = [] # start a new list for a given stress period
for layCellTuple in tempLayCellTupleList:
drainSpd[i].append([layCellTuple,-4.5,0.01])
i+=1
drn = flopy.mf6.ModflowGwfdrn(gwf, stress_period_data=drainSpd) ## Org
You have inserted a fixed elevation
You have inserted a fixed elevation
You have inserted a fixed elevation
You have inserted a fixed elevation
#trench plot
drn.plot(mflay=1, kper=2) ## 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
xsect = flopy.plot.PlotCrossSection(model=gwf, line={'Line': sectionLine}) ## Org
lc = xsect.plot_grid(lw=0.5, alpha=0.3) ## Org
xsect.plot_bc('DRN',kper=3) ## <== updated
ax.grid(alpha=0.2) ## Org
#river package
#layCellTupleList, cellElevList = getLayCellElevTupleFromRaster(gwf,interIx,'rst/elevWgs18S.tif','shp/riversSpart.shp') ## Org
#riverSpd = {} ## Org
#riverSpd[0] = [] ## Org
#for index, layCellTuple in enumerate(layCellTupleList): ## Org
# riverSpd[0].append([layCellTuple,cellElevList[index],0.01]) ## Org
#riv = flopy.mf6.ModflowGwfdrn(gwf, stress_period_data=riverSpd) ## Org
#river plot
#riv.plot(mflay=0, kper=1) ## 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
#xsect = flopy.plot.PlotCrossSection(model=gwf, line={'Line': sectionLine}) ## Org
#lc = xsect.plot_grid(lw=0.5) ## Org
#xsect.plot_bc('DRN',kper=4) ## Org
#ax.grid() ## 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_0...
writing model mf6Model...
writing model name file...
writing package disv...
writing package ic...
writing package npf...
writing package sto...
writing package ghb_0...
INFORMATION: maxbound in ('gwf6', 'ghb', 'dimensions') changed to 1486 based on size of stress_period_data
writing package drn_0...
INFORMATION: maxbound in ('gwf6', 'drn', 'dimensions') changed to 102 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/12 9:50:17
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
Run end date and time (yyyy/mm/dd hh:mm:ss): 2025/09/12 9:50:25
Elapsed run time: 7.248 Seconds
Normal termination of simulation.
Model output visualization
headObj = gwf.output.head() ## Org
headObj.get_kstpkper() ## Org
[(np.int32(0), np.int32(0)),
(np.int32(0), np.int32(1)),
(np.int32(0), np.int32(2)),
(np.int32(0), np.int32(3)),
(np.int32(0), np.int32(4))]
kper = 3 ## 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)
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()
generateRasterFromArray(gwf,
waterTable,
meshLayer=None,
rasterRes=2,
epsg=crossSection.crs.to_epsg(),
outputPath='../output/waterTable.tif',
limitLayer=None)
Raster X Dim: 2220.00, Raster Y Dim: 1350.00
Number of cols: 1111, Number of rows: 676
generateContoursFromRaster('../output/waterTable.tif',
interval = 1,
outputPath = '../output/waterTable_0_5_m.tif')
3d geometry generation on Vtk format
#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 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 exists and has been cleared
Follow us: |
|
|
|
|
|
|
/---------------------------------------/
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', 'GHB_0', 'DRN_0', 'OC']
#show output data
headObj = mf6Vtk.gwf.output.head() ## Org
headObj.get_kstpkper() ## Org
[(np.int32(0), np.int32(0)),
(np.int32(0), np.int32(1)),
(np.int32(0), np.int32(2)),
(np.int32(0), np.int32(3)),
(np.int32(0), np.int32(4))]
#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=1) ## Org
/--------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.0549 seconds to be generated.
/--------GHB_0 vtk generated-------/
/--------DRN_0 vtk generation-------/
Working for DRN_0 package, creating the datasets: ('elev', 'cond')
Vtk file took 0.2655 seconds to be generated.
/--------DRN_0 vtk generated-------/
mf6Vtk.generateHeadVtk(nper=1, crop=True) ## Org
mf6Vtk.generateWaterTableVtk(nper=1) ## Org
Input data
You can download the input data from this link:
owncloud.hatarilabs.com/s/KIJP2EC3SKuFMF3
Password: Hatarilabs