API Reference#
OasisMove Solvers#
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.assemble_first_inner_iter(A, a_conv, dt, M, scalar_components, les_model, nn_model, a_scalar, K, nu, nut_, nunn_, u_components, LT, KT, NT, b_tmp, b0, x_1, x_2, u_ab, bcs, mesh, boundary, u, v, backflow_facets, backflow_beta, wx_, bw_tmp, bw0, use_supg, **NS_namespace)[source]#
Called on first inner iteration of velocity/pressure system.
Assemble convection matrix, compute rhs of tentative velocity and reset coefficient matrix for solve.
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.get_solvers(use_krylov_solvers, krylov_solvers, krylov_solvers_w, scalar_components, velocity_krylov_solver, pressure_krylov_solver, scalar_krylov_solver, mesh_velocity_krylov_solver, **NS_namespace)[source]#
Return linear solvers.
- We are solving for
mesh velocity
tentative velocity
pressure correction
and possibly: - scalars
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.mesh_velocity_assemble(A_mesh, ui, bw, bw_tmp, a_mesh, bc_mesh, A_cache, **NS_namespace)[source]#
Assemble rhs of mesh velocity equation.
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.mesh_velocity_solve(A_mesh, bw, wx_, w_, dof_map, dt, dx_, coordinates, w_sol, ui, bc_mesh, **NS_namespace)[source]#
Solve mesh equation.
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.pressure_assemble(b, x_, dt, Ap, divu, **NS_namespace)[source]#
Assemble rhs of pressure equation.
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.pressure_solve(dp_, x_, Ap, b, p_sol, bcs, **NS_namespace)[source]#
Solve pressure equation.
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.print_velocity_pressure_info(num_iter, print_velocity_pressure_convergence, norm, info_blue, inner_iter, udiff, dp_, **NS_namespace)[source]#
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.scalar_assemble(a_scalar, a_conv, Ta, dt, M, scalar_components, Schmidt_T, KT, nu, Schmidt, b, K, x_1, b0, les_model, nn_model, use_supg, u_ab, u, v, mesh, **NS_namespace)[source]#
Assemble scalar equation.
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.scalar_solve(ci, scalar_components, Ta, b, x_, bcs, c_sol, nu, Schmidt, K, use_supg, **NS_namespace)[source]#
Solve scalar equation.
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.setup(u_components, u, v, p, q, u_w, v_w, bcs, les_model, nut_, scalar_components, V, Q, x_, p_, u_, velocity_update_solver, assemble_matrix, homogenize, GradFunction, DivFunction, LESsource, nn_model, nunn_, NNsource, mesh, **NS_namespace)[source]#
Preassemble mass and diffusion matrices.
Set up and prepare all equations to be solved. Called once, before going into time loop.
- oasismove.solvers.NSfracStep.IPCS_ABCN_Move.velocity_tentative_assemble(ui, b, b_tmp, p_, gradp, **NS_namespace)[source]#
Add pressure gradient to rhs of tentative velocity system.
Common methods#
- oasismove.common.io.check_if_kill(folder, killtime, total_timer)[source]#
Check if user has put a file named killoasis in folder or if given killtime has been reached.
- oasismove.common.io.check_if_reset_statistics(folder)[source]#
Check if user has put a file named resetoasis in folder.
- oasismove.common.io.create_initial_folders(folder, restart_folder, sys_comp, tstep, info_red, scalar_components, output_timeseries_as_vector, **NS_namespace)[source]#
Create necessary folders.
- oasismove.common.io.init_from_restart(restart_folder, sys_comp, uc_comp, u_components, q_, q_1, q_2, w_, d_, tstep, velocity_degree, previous_velocity_degree, mesh, constrained_domain, V, Q, **NS_namespace)[source]#
Initialize solution from checkpoint files
- oasismove.common.io.save_checkpoint_solution_xdmf(q_, q_1, w_, d_, newfolder, u_components, mesh, NS_parameters)[source]#
Overwrite solution in Checkpoint folder
- oasismove.common.io.save_solution(tstep, t, q_, q_1, w_, d_, folder, newfolder, save_step, checkpoint, NS_parameters, tstepfiles, u_, u_components, output_timeseries_as_vector, mesh, AssignedVectorFunction, killtime, total_timer, **NS_namespace)[source]#
Called at end of timestep. Check for kill and save solution if required.
- oasismove.common.io.save_tstep_solution_xdmf(tstep, q_, u_, newfolder, tstepfiles, output_timeseries_as_vector, AssignedVectorFunction, NS_parameters)[source]#
Store solution on current timestep to XDMF file.
- class oasismove.common.utilities.AssignedVectorFunction(u, name='Assigned Vector Function')[source]#
Bases:
Function
Vector function used for postprocessing.
Assign data from ListTensor components using FunctionAssigner.
- class oasismove.common.utilities.CG1Function(form, mesh, bcs=[], name='CG1', method={}, bounded=False)[source]#
Bases:
OasisFunction
Function used for projecting a CG1 space, possibly using a weighted average.
Typically used for computing turbulent viscosity in LES.
- class oasismove.common.utilities.DivFunction(u_, Space, bcs=[], name='div', method={})[source]#
Bases:
OasisFunction
Function used for projecting divergence of vector.
Typically used for computing divergence of velocity on pressure function space.
- class oasismove.common.utilities.GradFunction(p_, Space, i=0, bcs=[], name='grad', method={})[source]#
Bases:
OasisFunction
Function used for projecting gradients.
Typically used for computing pressure gradient on velocity function space.
- class oasismove.common.utilities.LESsource(nut, u, Space, bcs=[], name='')[source]#
Bases:
Function
Function used for computing the transposed source to the LES equation.
- class oasismove.common.utilities.Mat_cache_dict[source]#
Bases:
dict
Items in dictionary are matrices stored for efficient reuse.
- class oasismove.common.utilities.NNsource(nunn, u, Space, bcs=[], name='')[source]#
Bases:
Function
Function used for computing the transposed source to the non-Newtonian equation.
- class oasismove.common.utilities.OasisFunction(form, Space, bcs=[], name='x', matvec=[None, None], method='default', solver_type='cg', preconditioner_type='default')[source]#
Bases:
Function
Function with more or less efficient projection methods of associated linear form.
The matvec option is provided for letting the right hand side be computed through a fast matrix vector product. Both the matrix and the Coefficient of the required vector must be provided.
- method = “default”
Solve projection with regular linear algebra using solver_type and preconditioner_type
- method = “lumping”
Solve through lumping of mass matrix
- class oasismove.common.utilities.Solver_cache_dict[source]#
Bases:
dict
Items in dictionary are Linear algebra solvers stored for efficient reuse.
OasisMove problems#
- oasismove.problems.__init__.add_function_to_tstepfiles(function, newfolder, tstepfiles, tstep)[source]#
- oasismove.problems.__init__.compute_flow_quantities(u, L, nu, mesh, t, tstep, dt, h, outlet_area=1, boundary=None, outlet_ids=[], inlet_ids=[], id_wall=0, period=1.0, newfolder=None, dynamic_mesh=False, write_to_file=False)[source]#
Compute max and mean Reynolds numbers, CFL numbers, and fluxes through boundaries.
- Parameters
u (Function) – Velocity field.
L (float) – Characteristic length.
nu (float) – Kinematic viscosity of fluid.
mesh (Mesh) – Computational mesh.
t (float) – Current time.
tstep (int) – Current time step.
dt (float) – Time step size.
h (float) – Mesh element size.
outlet_area (float) – Area of the outlet, default is 1.
boundary (MeshFunction) – Boundary mesh function, default is None.
outlet_ids (list of int) – List of outlet boundary IDs, default is empty list.
inlet_ids (list of int) – List of inlet boundary IDs, default is empty list.
id_wall (int) – Wall boundary ID, default is 0.
period (float) – Period of oscillation, default is 1.0.
newfolder (str) – Directory for output files, default is None.
dynamic_mesh (bool) – Flag for dynamic mesh, default is False.
write_to_file (bool) – Flag to write results to file, default is False.
- oasismove.problems.__init__.create_bcs(sys_comp, **NS_namespace)[source]#
Return dictionary of Dirichlet boundary conditions.
- oasismove.problems.__init__.post_import_problem(NS_parameters, mesh, commandline_kwargs, NS_expressions, **NS_namespace)[source]#
Called after importing from problem.
- oasismove.problems.__init__.pre_solve_hook(**NS_namespace)[source]#
Called just prior to entering time-loop. Must return a dictionary.
- oasismove.problems.__init__.problem_parameters(**NS_namespace)[source]#
Updates problem specific parameters, and handles restart
- oasismove.problems.__init__.recursive_update(dst, src)[source]#
Update dict dst with items from src deeply (“deep update”).
- oasismove.problems.__init__.scalar_source(scalar_components, **NS_namespace)[source]#
Return a dictionary of scalar sources.
- class oasismove.problems.NSfracStep.MovingCommon.Counter(**kwargs)[source]#
Bases:
UserExpression
A UserExpression designed for evaluating and counting mesh nodes, and associating each node with a unique integer, thus keeping a reference to the original mesh.
- map#
A dictionary that maps a unique integer to a node coordinate.
- Type
dict
- counter#
An integer that counts the number of evaluated nodes.
- Type
int
- oasismove.problems.NSfracStep.MovingCommon.get_coordinate_map(V, boundary, w_, facet_id)[source]#
Create a mapping of coordinates for a given function space and boundary.
- Parameters
V (dolfin.function.Functionspace) – A FunctionSpace object representing the field on which the map is created.
boundary (dolfin.cpp.mesh.MeshFunction) – A MeshFunction defining with marked boundaries for the mesh.
w (dict) – A dict containing the mesh velocity at the current time step.
facet_id (int) – An integer representing the id of the boundary facets.
- Returns
An integer count of the number of elements in the field V. x_hat_map (dict): A dict mapping coordinates to elements in the field V.
- Return type
counter_max (int)
- oasismove.problems.NSfracStep.MovingCommon.get_visualization_writers(newfolder, list_of_quantities)[source]#
Create a list of XDMFFile writers for different quantities.
- Parameters
newfolder (str) – The path of the folder where the XDMF files will be saved.
list_of_quantities (list) – A list of string names representing different quantities to be written,
'pressure'. (such as 'velocity' or) –
- Returns
A list of XDMFFile objects for writing each quantity to a separate file.
- Return type
writers (list)
- oasismove.problems.NSfracStep.MovingCylinder.create_bcs(V, Q, D, u_inf, St, F, A_ratio, sys_comp, boundary, NS_expressions, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingCylinder.pre_boundary_condition(D, Re, u_inf, mesh, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingCylinder.pre_solve_hook(V, p_, u_, velocity_degree, nu, mesh, newfolder, u_components, boundary, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingCylinder.problem_parameters(commandline_kwargs, NS_parameters, NS_expressions, **NS_namespace)[source]#
Problem file for running CFD simulation for the oscillating cylinder in a rectangular 2D domain, as described by Blackburn and Henderson [1].The cylinder is prescribed an oscillatory motion and is placed in a free stream, with a diameter of D cm. The kinetmatic viscosity is computed from the free stream velocity of 1 m/s for a Reynolds number of 500, which is well above the critical value for vortex shedding. Moreover, the oscillation is mainly controlled by the amplitude ratio A_ratio, the Strouhal number St, and the frequency ratio F.
[1] Blackburn, H. M., & Henderson, R. D. (1999). A study of two-dimensional flow past an oscillating cylinder. Journal of Fluid Mechanics, 385, 255-286.
- oasismove.problems.NSfracStep.MovingCylinder.temporal_hook(mesh, t, dt, nu, St, F, A_ratio, tstep, save_solution_frequency, forces, q_, u_inf, D, viz_u, viz_p, u_vec, newfolder, p_, u_, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingCylinder.update_boundary_conditions(t, NS_expressions, **NS_namespace)[source]#
- class oasismove.problems.NSfracStep.MovingTaylorGreen3D.PeriodicDomain[source]#
Bases:
SubDomain
- oasismove.problems.NSfracStep.MovingTaylorGreen3D.initialize(q_, q_1, q_2, VV, initial_fields, OasisFunction, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingTaylorGreen3D.pre_boundary_condition(mesh, Re, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingTaylorGreen3D.pre_solve_hook(V, mesh, A0, T_G, L, t, newfolder, initial_fields_w, velocity_degree, u_components, boundary, NS_expressions, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingTaylorGreen3D.problem_parameters(commandline_kwargs, NS_parameters, NS_expressions, **NS_namespace)[source]#
Problem file for running CFD simulation for the moving Taylor-Green problem in 3D as described by Fehn et al.[1]. The problem solves the N-S equations in the absence of body forces, and is commonly used to study transitional and turbulent flows. The problem initializes the solution at the two previous time steps, and applies periodic boundary condition on the domain walls in all coordinate directions. The domain boundaries are moving for the given mesh motion, while the fluid velocity is defined periodically in order to ensure consistency with the periodic boundary conditions. The main simulation parameters are the ampliture A, and period time duration T_G.
[1] Fehn, N., Heinz, J., Wall, W. A., & Kronbichler, M. (2021). High-order arbitrary Lagrangian–Eulerian discontinuous Galerkin methods for the incompressible Navier–Stokes equations. Journal of Computational Physics, 430, 110040.
- oasismove.problems.NSfracStep.MovingTaylorGreen3D.temporal_hook(nu, L, dt, mesh, u_, save_solution_frequency, u_vec, viz_u, tstep, t, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingTaylorGreen3D.update_boundary_conditions(t, NS_expressions, **NS_namespace)[source]#
- class oasismove.problems.NSfracStep.MovingWall.MovingSideWall(t, sigma, h0, eps, x_hat_map, counter_max, **kwargs)[source]#
Bases:
UserExpression
- class oasismove.problems.NSfracStep.MovingWall.MovingWall(t, sigma, h0, eps, **kwargs)[source]#
Bases:
UserExpression
- oasismove.problems.NSfracStep.MovingWall.create_bcs(V, Q, w_, sigma, h0, eps, sys_comp, boundary, NS_expressions, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingWall.mesh(eps, h0, scale, dt, Nx, Ny, u_max, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingWall.pre_boundary_condition(mesh, h0, eps, scale, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingWall.pre_solve_hook(V, mesh, newfolder, velocity_degree, u_components, boundary, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingWall.problem_parameters(NS_parameters, **NS_namespace)[source]#
Problem file for running CFD simulation for the MovingWall problem inspired by the Wall-induced channel flow presented by Chnafa [1]. The problem considers flow in a long, straight, time-dependent domain, where the flow is induced by a moving wall at y=h(t). The motion is mainly controlled by the pulsation (sigma) and the ampliture of oscillations (eps), with an initial height of h0. The flow has an analytic solution for small Reynolds numbers, and may be used as a validation case.
[1] Chnafa, C. (2014). Using image-based large-eddy simulations to investigate the intracardiac flow and its turbulent nature (Doctoral dissertation, Université Montpellier II-Sciences et Techniques du Languedoc).
- oasismove.problems.NSfracStep.MovingWall.temporal_hook(mesh, dt, h0, eps, nu, t, tstep, save_solution_frequency, p_, u_, viz_u, viz_p, u_vec, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingWall.update_boundary_conditions(t, NS_expressions, **NS_namespace)[source]#
- class oasismove.problems.NSfracStep.MovingVortex.Walls(t, coor, A, L, T_G, x_hat_map, counter_max, **kwargs)[source]#
Bases:
UserExpression
- oasismove.problems.NSfracStep.MovingVortex.create_bcs(V, t, dt, nu, sys_comp, boundary, initial_fields, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingVortex.initialize(q_, q_1, q_2, VV, t, nu, dt, initial_fields, **NS_namespace)[source]#
Initialize solution.
Use t=dt/2 for pressure since pressure is computed in between timesteps.
- oasismove.problems.NSfracStep.MovingVortex.pre_solve_hook(V, mesh, t, nu, L, w_, T_G, A0, newfolder, velocity_degree, u_components, boundary, exact_fields, **NS_namespace)[source]#
- oasismove.problems.NSfracStep.MovingVortex.problem_parameters(commandline_kwargs, NS_parameters, NS_expressions, **NS_namespace)[source]#
Problem file for running CFD simulation for the MovingVortex problem inspired by the problem by Fehn et al.[1], resembling the 2D Taylor-Green vortex. The problem solves the N-S equations in the absence of body forces, with a manufactured velocity solution. The mesh velocity is also described by an analytic displacement field, describing the oscillatory boundary movement. The movement is mainly controlled by the amplitude A0, and period length of the mesh motion T_G.
[1] Fehn, N., Heinz, J., Wall, W. A., & Kronbichler, M. (2021). High-order arbitrary Lagrangian–Eulerian discontinuous Galerkin methods for the incompressible Navier–Stokes equations. Journal of Computational Physics, 430, 110040.
- oasismove.problems.NSfracStep.MovingVortex.temporal_hook(u_, q_, t, nu, VV, dt, u_vec, ue_vec, p_, viz_u, viz_p, viz_ue, initial_fields, tstep, save_solution_frequency, sys_comp, compute_error, total_error, ue_x, ue_y, pe, testing, **NS_namespace)[source]#
Function called at end of timestep.
Plot solution and compute error by comparing to analytical solution. Remember pressure is computed in between timesteps.