Stokes integrators
stokes_integrators
This module provides a fine-tuned solver class for the Stokes system of conservation
equations and a function to automatically set the associated null spaces. Users
instantiate the StokesSolver
class by providing relevant parameters and call the
solve
method to request a solver update.
iterative_stokes_solver_parameters = {'mat_type': 'matfree', 'ksp_type': 'preonly', 'pc_type': 'fieldsplit', 'pc_fieldsplit_type': 'schur', 'pc_fieldsplit_schur_type': 'full', 'fieldsplit_0': {'ksp_type': 'cg', 'ksp_rtol': 1e-05, 'ksp_max_it': 1000, 'pc_type': 'python', 'pc_python_type': 'gadopt.SPDAssembledPC', 'assembled_pc_type': 'gamg', 'assembled_mg_levels_pc_type': 'sor', 'assembled_pc_gamg_threshold': 0.01, 'assembled_pc_gamg_square_graph': 100, 'assembled_pc_gamg_coarse_eq_limit': 1000, 'assembled_pc_gamg_mis_k_minimum_degree_ordering': True}, 'fieldsplit_1': {'ksp_type': 'fgmres', 'ksp_rtol': 0.0001, 'ksp_max_it': 200, 'pc_type': 'python', 'pc_python_type': 'firedrake.MassInvPC', 'Mp_pc_type': 'ksp', 'Mp_ksp_ksp_rtol': 1e-05, 'Mp_ksp_ksp_type': 'cg', 'Mp_ksp_pc_type': 'sor'}}
module-attribute
Default iterative solver parameters for solution of stokes system.
We configure the Schur complement approach as described in Section of 4.3 of Davies et al. (2022), using PETSc's fieldsplit preconditioner type, which provides a class of preconditioners for mixed problems that allows a user to apply different preconditioners to different blocks of the system.
The fieldsplit_0 entries configure solver options for the velocity block. The linear systems associated with this matrix are solved using a combination of the Conjugate Gradient (cg) method and an algebraic multigrid preconditioner (GAMG).
The fieldsplit_1 entries contain solver options for the Schur complement solve itself. For preconditioning, we approximate the Schur complement matrix with a mass matrix scaled by viscosity, with the viscosity provided through the optional mu keyword argument to Stokes solver. Since this preconditioner step involves an iterative solve, the Krylov method used for the Schur complement needs to be of flexible type, and we use FGMRES by default.
We note that our default solver parameters can be augmented or adjusted by accessing the solver_parameter dictionary, for example:
stokes_solver.solver_parameters['fieldsplit_0']['ksp_converged_reason'] = None
stokes_solver.solver_parameters['fieldsplit_0']['ksp_rtol'] = 1e-3
stokes_solver.solver_parameters['fieldsplit_0']['assembled_pc_gamg_threshold'] = -1
stokes_solver.solver_parameters['fieldsplit_1']['ksp_converged_reason'] = None
stokes_solver.solver_parameters['fieldsplit_1']['ksp_rtol'] = 1e-2
Note
G-ADOPT defaults to iterative solvers in 3-D.
direct_stokes_solver_parameters = {'mat_type': 'aij', 'ksp_type': 'preonly', 'pc_type': 'lu', 'pc_factor_mat_solver_type': 'mumps'}
module-attribute
Default direct solver parameters for solution of Stokes system.
Configured to use LU factorisation, using the MUMPS library.
Note
G-ADOPT defaults to direct solvers in 2-D.
newton_stokes_solver_parameters = {'snes_type': 'newtonls', 'snes_linesearch_type': 'l2', 'snes_max_it': 100, 'snes_atol': 1e-10, 'snes_rtol': 1e-05}
module-attribute
Default solver parameters for non-linear systems.
We use a setup based on Newton's method (newtonls) with a secant line search over the L2-norm of the function.
StokesSolver(z, T, approximation, bcs={}, mu=1, quad_degree=6, solver_parameters=None, J=None, constant_jacobian=False, **kwargs)
Solves the Stokes system in place.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
z |
Function
|
Firedrake function representing mixed Stokes system |
required |
T |
Function
|
Firedrake function representing temperature |
required |
approximation |
BaseApproximation
|
Approximation describing system of equations |
required |
bcs |
dict[int, dict[str, Number]]
|
Dictionary of identifier-value pairs specifying boundary conditions |
{}
|
mu |
Function | Number
|
Firedrake function representing dynamic viscosity |
1
|
quad_degree |
int
|
Quadrature degree. Default value is |
6
|
solver_parameters |
Optional[dict[str, str | Number] | str]
|
Either a dictionary of PETSc solver parameters or a string specifying a default set of parameters defined in G-ADOPT |
None
|
J |
Optional[Function]
|
Firedrake function representing the Jacobian of the system |
None
|
constant_jacobian |
bool
|
Whether the Jacobian of the system is constant |
False
|
Source code in g-adopt/gadopt/stokes_integrators.py
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 |
|
setup_solver()
Sets up the solver.
Source code in g-adopt/gadopt/stokes_integrators.py
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 |
|
solve()
Solves the system.
Source code in g-adopt/gadopt/stokes_integrators.py
325 326 327 328 329 330 |
|
create_stokes_nullspace(Z, closed=True, rotational=False, translations=None, ala_approximation=None, top_subdomain_id=None)
Create a null space for the mixed Stokes system.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
Z |
WithGeometry
|
Firedrake mixed function space associated with the Stokes system |
required |
closed |
bool
|
Whether to include a constant pressure null space |
True
|
rotational |
bool
|
Whether to include all rotational modes |
False
|
translations |
Optional[list[int]]
|
List of translations to include |
None
|
ala_approximation |
Optional[AnelasticLiquidApproximation]
|
AnelasticLiquidApproximation for calculating (non-constant) right nullspace |
None
|
top_subdomain_id |
Optional[str | int]
|
Boundary id of top surface. Required when providing ala_approximation. |
None
|
Returns:
Type | Description |
---|---|
MixedVectorSpaceBasis
|
A Firedrake mixed vector space basis incorporating the null space components |
Source code in g-adopt/gadopt/stokes_integrators.py
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 |
|
ala_right_nullspace(W, approximation, top_subdomain_id)
Compute pressure nullspace for Anelastic Liquid Approximation.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
W |
WithGeometry
|
pressure function space |
required |
approximation |
AnelasticLiquidApproximation
|
AnelasticLiquidApproximation with equation parameters |
required |
top_subdomain_id |
str | int
|
boundary id of top surface |
required |
Returns:
Type | Description |
---|---|
pressure nullspace solution |
To obtain the pressure nullspace solution for the Stokes equation in Anelastic Liquid Approximation, which includes a pressure-dependent buoyancy term, we try to solve the equation:
Taking the divergence:
then testing it with q:
followed by integration by parts:
This elliptic equation can be solved with natural boundary conditions by imposing our original equation above, which eliminates all boundary terms:
However, if we do so on all boundaries we end up with a system that has the same nullspace, as the one we are after (note that we ended up merely testing the original equation with \(nabla q\)). Instead we use the fact that the gradient of the null mode is always vertical, and thus the null mode is constant at any horizontal level (geoid), specifically the top surface. Choosing any nonzero constant for this surface fixes the arbitrary scalar multiplier of the null mode. We choose the value of one and apply it as a Dirichlet boundary condition.
Note that this procedure does not necessarily compute the exact nullspace of the discretised Stokes system. In particular, since not every test function \(v in V\), the velocity test space, can be written as \(v=nabla q\) with \(q in W\), the pressure test space, the two terms do not necessarily exactly cancel when tested with \(v\) instead of \(nabla q\) as in our final equation. However, in practice the discrete error appears to be small enough, and providing this nullspace gives an improved convergence of the iterative Stokes solver.
Source code in g-adopt/gadopt/stokes_integrators.py
333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 |
|