-
Notifications
You must be signed in to change notification settings - Fork 0
/
solvers_mod.py
352 lines (283 loc) · 12.6 KB
/
solvers_mod.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
import logging, time
import scipy.sparse as sps
import scipy.io as sps_io
import numpy as np
import porepy as pp
# ------------------------------------------------------------------------------#
logger = logging.getLogger(__name__)
def run_flow(gb, node_discretization, source_discretization, folder, is_FV, discr_3d=None, discr_2d=None):
if discr_3d is None:
discr_3d = node_discretization
if discr_2d is None:
discr_2d = node_discretization
kw_t = "transport" # For storing fluxes
grid_variable = "pressure"
mortar_variable = "mortar_flux"
# Identifier of the discretization operator on each grid
diffusion_term = "diffusion"
# Identifier of the discretization operator between grids
coupling_operator_keyword = "coupling_operator"
edge_discretization_3d = pp.RobinCoupling(
"flow", discr_3d, discr_2d
)
edge_discretization_2d = pp.RobinCoupling(
"flow", discr_2d, node_discretization
)
edge_discretization = pp.RobinCoupling(
"flow", node_discretization, node_discretization
)
# Loop over the nodes in the GridBucket, define primary variables and discretization schemes
for g, d in gb:
# Assign primary variables on this grid. It has one degree of freedom per cell.
if g.dim == 3 and isinstance(discr_3d, pp.MVEM):
d[pp.PRIMARY_VARIABLES] = {grid_variable: {"cells": 1, "faces": 1}}
elif g.dim == 2 and isinstance(discr_2d, pp.MVEM):
d[pp.PRIMARY_VARIABLES] = {grid_variable: {"cells": 1, "faces": 1}}
else:
d[pp.PRIMARY_VARIABLES] = {grid_variable: {"cells": 1, "faces": int(not is_FV)}}
# Assign discretization operator for the variable.
if g.dim == 3:
d[pp.DISCRETIZATION] = {
grid_variable: {
diffusion_term: discr_3d,
}
}
elif g.dim == 2:
d[pp.DISCRETIZATION] = {
grid_variable: {
diffusion_term: discr_2d,
}
}
else:
d[pp.DISCRETIZATION] = {
grid_variable: {
diffusion_term: node_discretization,
}
}
# Loop over the edges in the GridBucket, define primary variables and discretizations
for e, d in gb.edges():
g1, g2 = gb.nodes_of_edge(e)
# The mortar variable has one degree of freedom per cell in the mortar grid
d[pp.PRIMARY_VARIABLES] = {mortar_variable: {"cells": 1}}
# The coupling discretization links an edge discretization with variables
# and discretization operators on each neighboring grid
if g2.dim == 3:
d[pp.COUPLING_DISCRETIZATION] = {
coupling_operator_keyword: {
g1: (grid_variable, diffusion_term),
g2: (grid_variable, diffusion_term),
e: (mortar_variable, edge_discretization_3d),
}
}
elif g2.dim == 2:
d[pp.COUPLING_DISCRETIZATION] = {
coupling_operator_keyword: {
g1: (grid_variable, diffusion_term),
g2: (grid_variable, diffusion_term),
e: (mortar_variable, edge_discretization_2d),
}
}
else:
d[pp.COUPLING_DISCRETIZATION] = {
coupling_operator_keyword: {
g1: (grid_variable, diffusion_term),
g2: (grid_variable, diffusion_term),
e: (mortar_variable, edge_discretization),
}
}
logger.info(
"Assembly and discretization using the discretizer " + str(node_discretization)
)
tic = time.time()
assembler = pp.Assembler()
# Assemble the linear system, using the information stored in the GridBucket
A, b, block_dof, full_dof = assembler.assemble_matrix_rhs(
gb, active_variables=[grid_variable, mortar_variable]
)
logger.info("Done. Elapsed time: " + str(time.time() - tic))
logger.info("Linear solver")
tic = time.time()
x = sps.linalg.spsolve(A, b)
logger.info("Done. Elapsed time " + str(time.time() - tic))
# if is_FV:
assembler.distribute_variable(gb, x, block_dof, full_dof)
if is_FV:
if isinstance(discr_3d, pp.MVEM):
g = gb.grids_of_dimension(3)[0]
d = gb.node_props(g)
discr_scheme = d[pp.DISCRETIZATION][grid_variable][diffusion_term]
d["pressure"] = discr_scheme.extract_pressure(g, d[grid_variable], d)
if isinstance(discr_2d, pp.MVEM):
for g in gb.grids_of_dimension(2):
d = gb.node_props(g)
discr_scheme = d[pp.DISCRETIZATION][grid_variable][diffusion_term]
d["pressure"] = discr_scheme.extract_pressure(g, d[grid_variable], d)
# pp.fvutils.compute_darcy_flux(gb, lam_name=mortar_variable, keyword_store=kw_t)
else:
for g, d in gb:
discr_scheme = d[pp.DISCRETIZATION][grid_variable][diffusion_term]
#d[pp.PARAMETERS][kw_t]["darcy_flux"] = discr_scheme.extract_flux(
# g, d[grid_variable], d
#)
# Note the order: we overwrite d["pressure"] so this has to be done after
# extracting the flux
d["pressure"] = discr_scheme.extract_pressure(g, d[grid_variable], d)
for e, d in gb.edges():
# g1, g2 = gb.nodes_of_edge(e)
d[pp.PARAMETERS][kw_t]["darcy_flux"] = d[mortar_variable].copy()
export_flow(gb, folder)
# sps_io.mmwrite(folder + "/matrix.mtx", A)
return A, b, block_dof, full_dof
def export_flow(gb, folder):
gb.add_node_props(["cell_volumes", "cell_centers"])
for g, d in gb:
d["cell_volumes"] = g.cell_volumes
d["cell_centers"] = g.cell_centers
save = pp.Exporter(gb, "sol", folder=folder)
props = ["pressure", "cell_volumes", "cell_centers"]
# extra properties, problem specific
if all(gb.has_nodes_prop(gb.get_grids(), "low_zones")):
gb.add_node_props("bottom_domain")
for g, d in gb:
d["bottom_domain"] = 1 - d["low_zones"]
props.append("bottom_domain")
if all(gb.has_nodes_prop(gb.get_grids(), "color")):
props.append("color")
if all(gb.has_nodes_prop(gb.get_grids(), "aperture")):
props.append("aperture")
save.write_vtk(props)
# ------------------------------------------------------------------------------#
def solve_rt0(gb, folder, discr_3d=None):
# Choose and define the solvers and coupler
flow_discretization = pp.RT0("flow")
source_discretization = pp.DualScalarSource("flow")
run_flow(gb, flow_discretization, source_discretization, folder, is_FV=False, discr_3d=discr_3d)
# ------------------------------------------------------------------------------#
def solve_tpfa(gb, folder, discr_3d=None):
# Choose and define the solvers and coupler
flow_discretization = pp.Tpfa("flow")
source_discretization = pp.ScalarSource("flow")
run_flow(gb, flow_discretization, source_discretization, folder, is_FV=True, discr_3d=discr_3d)
# ------------------------------------------------------------------------------#
def solve_mpfa(gb, folder, discr_3d=None, discr_2d=None):
# Choose and define the solvers and coupler
flow_discretization = pp.Mpfa("flow")
source_discretization = pp.ScalarSource("flow")
run_flow(gb, flow_discretization, source_discretization, folder, is_FV=True,
discr_3d=discr_3d, discr_2d=discr_2d)
# ------------------------------------------------------------------------------#
def solve_vem(gb, folder, discr_3d=None):
# Choose and define the solvers and coupler
flow_discretization = pp.MVEM("flow")
source_discretization = pp.DualScalarSource("flow")
run_flow(gb, flow_discretization, source_discretization, folder, is_FV=False, discr_3d=discr_3d)
# ------------------------------------------------------------------------------#
def transport(gb, data, solver_name, folder, callback=None, save_every=1):
grid_variable = "tracer"
mortar_variable = "mortar_tracer"
kw = "transport"
# Identifier of the discretization operator on each grid
advection_term = "advection"
source_term = "source"
mass_term = "mass"
# Identifier of the discretization operator between grids
advection_coupling_term = "advection_coupling"
# Discretization objects
node_discretization = pp.Upwind(kw)
source_discretization = pp.ScalarSource(kw)
mass_discretization = pp.MassMatrix(kw)
edge_discretization = pp.UpwindCoupling(kw)
# Loop over the nodes in the GridBucket, define primary variables and discretization schemes
for g, d in gb:
# Assign primary variables on this grid. It has one degree of freedom per cell.
d[pp.PRIMARY_VARIABLES] = {grid_variable: {"cells": 1, "faces": 0}}
# Assign discretization operator for the variable.
d[pp.DISCRETIZATION] = {
grid_variable: {
advection_term: node_discretization,
source_term: source_discretization,
mass_term: mass_discretization,
}
}
# Loop over the edges in the GridBucket, define primary variables and discretizations
for e, d in gb.edges():
g1, g2 = gb.nodes_of_edge(e)
# The mortar variable has one degree of freedom per cell in the mortar grid
d[pp.PRIMARY_VARIABLES] = {mortar_variable: {"cells": 1}}
d[pp.COUPLING_DISCRETIZATION] = {
advection_coupling_term: {
g1: (grid_variable, advection_term),
g2: (grid_variable, advection_term),
e: (mortar_variable, edge_discretization),
}
}
d[pp.DISCRETIZATION_MATRICES] = {kw: {}}
assembler = pp.Assembler()
# Assemble the linear system, using the information stored in the GridBucket. By
# not adding the matrices, we can arrange them at will to obtain the efficient
# solver defined below, which LU factorizes the system only once for all time steps.
A, b, block_dof, full_dof = assembler.assemble_matrix_rhs(
gb, active_variables=[grid_variable, mortar_variable], add_matrices=False
)
advection_coupling_term += (
"_" + mortar_variable + "_" + grid_variable + "_" + grid_variable
)
mass_term += "_" + grid_variable
advection_term += "_" + grid_variable
source_term += "_" + grid_variable
lhs = A[mass_term] + data["time_step"] * (
A[advection_term] + A[advection_coupling_term]
)
rhs_source_adv = b[source_term] + data["time_step"] * (
b[advection_term] + b[advection_coupling_term]
)
IEsolver = sps.linalg.factorized(lhs)
n_steps = int(np.round(data["t_max"] / data["time_step"]))
# Initial condition
tracer = np.zeros(rhs_source_adv.size)
assembler.distribute_variable(
gb, tracer, block_dof, full_dof, variable_names=[grid_variable, mortar_variable]
)
# Exporter
exporter = pp.Exporter(gb, name="tracer", folder=folder)
export_fields = ["tracer"]
# Keep track of the outflow for each time step
outflow = np.empty(0)
# Time loop
for i in range(n_steps):
# Export existing solution (final export is taken care of below)
assembler.distribute_variable(
gb,
tracer,
block_dof,
full_dof,
variable_names=[grid_variable, mortar_variable],
)
if np.isclose(i % save_every, 0):
exporter.write_vtk(export_fields, time_step=int(i // save_every))
tracer = IEsolver(A[mass_term] * tracer + rhs_source_adv)
outflow = compute_flow_rate(gb, grid_variable, outflow)
if callback is not None:
callback(gb)
exporter.write_vtk(export_fields, time_step=(n_steps // save_every))
time_steps = np.arange(
0, data["t_max"] + data["time_step"], save_every * data["time_step"]
)
exporter.write_pvd(time_steps)
return tracer, outflow, A, b, block_dof, full_dof
def compute_flow_rate(gb, grid_variable, outflow):
# this function is only for the first benchmark case
for g, d in gb:
if g.dim < 3:
continue
faces, cells, sign = sps.find(g.cell_faces)
index = np.argsort(cells)
faces, sign = faces[index], sign[index]
discharge = d[pp.PARAMETERS]["transport"]["darcy_flux"].copy()
tracer = d[grid_variable].copy()
discharge[faces] *= sign
discharge[g.get_internal_faces()] = 0
discharge[discharge < 0] = 0
val = np.dot(discharge, np.abs(g.cell_faces) * tracer)
outflow = np.r_[outflow, val]
return outflow