opencmp.models package

Submodules

opencmp.models.base_model module

class opencmp.models.base_model.Model(config, t_param)[source]

Bases: abc.ABC

Base class which other models will subclass from.

_add_multiple_parameters(config)[source]

Function to add multiple time-varying parameters (whose values are not known in closed form) to the model from the config file.

Parameters

config (ConfigParser) – The ConfigParser

Return type

None

abstract _construct_fes()[source]

Function to construct the finite element space for the model.

Return type

FESpace

Returns

The finite element space.

_construct_linearization_terms()[source]

Function to construct the list of linearization terms.

If a specific model does not need linearization terms, return None.

Return type

Optional[List[GridFunction], None]

Returns

None or a list of grid functions used to linearize a non-linear model.

_create_test_and_trial_functions()[source]

Function to create the trial and test (weighting) function(s) for the model.

Return type

Tuple[List[ProxyFunction], List[ProxyFunction]]

Returns

Returns two lists, one for the trial and test (weighting) function(s).

abstract _define_bc_types()[source]

Function to specify the types of boundary conditions (BCs) available for this model.

Return type

List[str]

Returns

List containing the names, as str, of all available/allowable BC types.

abstract _define_model_components()[source]

Function to specify the mapping, as a dictionary, between variable names and their position within the finite element space, and all other locations which derive their indexing order from the finite element space.

For a single variable model, e.g. Poisson, the index value MUST be None. This is due to how NGSolve defines gridfunctions for single vs mixed finite element spaces.

Return type

Dict[str, Optional[int, None]]

Returns

Returns a dictionary whose keys are the names of the trial functions and values are the index of a given trial function within the finite element space.

abstract _define_model_local_error_components()[source]

Function to specify which model components (trial functions) should be included in the local error calculation.

Return type

Dict[str, bool]

Returns

Returns a dictionary whose keys are the names of the trial functions and values is a bool which indicates if the local error should be calculated for that trial function.

abstract _define_num_weak_forms()[source]

Function to specify the number of seperate weak forms this model has.

This is mostly used if a linearization scheme is used in order to de-couple multiple equations. E.g. Solving INS by solving conservation of mass first, using the previous iteration/time-step’s conservation of momentum values, and then use that result to solve for conservation of momentum.

Return type

int

Returns

Int representing the number of weak forms.

abstract _define_time_derivative_components()[source]

Function to specify which trial functions have a time derivative associated with them. Additionally, it specifies in which of the model’s (potentially) multiple weak forms the time derivative should be present.

The weak form must be derived so that the time derivative has a coefficient of 1, i.e. not multiplied by anything.

Return type

List[Dict[str, bool]]

Returns

A list of dictionaries, one for each weak form of the model (in order). The keys for each are the name of the trial function and value is a bool indicating if this trial function has a time derivative in the weak form represented by the index of the dictionary within the list.

_ds(marker)[source]

Function to get the appropriate (DG vs non-DG) version of ds for the given marker.

Parameters

marker (str) – String representing the mesh marker to get ds on.

Return type

DifferentialSymbol

Returns

The appropriate ds.

abstract _post_init()[source]

Function to do extra model-specific work after the ENTIRE Model.__init__() has run.

This approach is used, instead of letting the user override __init__() and then calling super().__init__() from their custom model, so that we can make it explicit which parts are inherited from which parent.

Return type

None

abstract _pre_init()[source]

Function to do extra model-specific work after _define_model_components, _define_model_local_error_components, _define_time_derivative_components, _define_num_weak_forms, and _define_bc_types have run.

Return type

None

abstract _set_model_parameters()[source]

Function to initialize model parameters and functions from the configfile.

Return type

None

abstract accept_stationary_solution(gfu)[source]

nasser@nasser: TODO

abstract static allows_explicit_schemes()[source]

Function to specify whether a given model works with explicit time integration schemes.

Return type

bool

Returns

True if the model can be used with fully explicit time integration schemes, else False.

apply_dirichlet_bcs_to(gfu, time_step=0)[source]

Function to set the Dirichlet boundary conditions within the solution GridFunction.

Parameters
  • gfu (GridFunction) – The GridFunction to add the Dirichlet boundary condition values to.

  • time_step (int) – Specifies which time step’s boundary condition values to use. This would usually be 0 (the t^n+1 values) except in the case of adaptive_three_step and Runge Kutta schemes.

Return type

None

construct_and_run_solver(a_assembled, L_assembled, precond, gfu)[source]

Function to construct the solver and run one solve on the provided gridfunction.

Parameters
  • a_assembled (BilinearForm) – The assembled bilinear form.

  • L_assembled (LinearForm) – The assembled linear form.

  • precond (Preconditioner) – The preconditioner.

  • gfu (GridFunction) – The gridfunction holding information about any Dirichlet BCs which will be updated to hold the solution.

abstract construct_bilinear_time_ODE(U, V, dt=<ngsolve.fem.Parameter object>, time_step=0)[source]

Function to construct a portion of the bilinear form for model variables WITH time derivatives (basically any bilinear form terms not already in construct_bilinear_time_coefficient).

A given model with multiple model variables may not include time derivatives of all of its model variables. It is necessary to split the bilinear form into two portions as some time discretization schemes use different coefficients for terms with and without time derivatives.

Parameters
  • U (Union[List[ProxyFunction], List[GridFunction]]) – A list of trial functions for the model’s finite element space.

  • V (List[ProxyFunction]) – A list of testing (weighting) functions for the model’s finite element space.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[BilinearForm]

Returns

List of the described portion of the bilinear form for the model. The list will contain multiple bilinear forms if the model involves iterating between different versions of a linearized weak form.

abstract construct_bilinear_time_coefficient(U, V, dt, time_step)[source]

Function to construct a portion of the bilinear form for model variables WITHOUT time derivatives.

A given model with multiple model variables may not include time derivatives of all of its model variables. It is necessary to split the bilinear form into two portions as some time discretization schemes use different coefficients for terms with and without time derivatives.

Parameters
  • U (List[ProxyFunction]) – A list of trial functions for the model’s finite element space.

  • V (List[ProxyFunction]) – A list of testing (weighting) functions for the model’s finite element space.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[BilinearForm]

Returns

List of the described portion of the bilinear form for the model. The list will contain multiple bilinear forms if the model involves iterating between different versions of a linearized weak form.

construct_gfu()[source]

Function to construct a solution GridFunction.

Return type

GridFunction

Returns

A GridFunction initialized on the finite element space of the model.

abstract construct_imex_explicit(V, gfu_0, dt, time_step)[source]

Function to construct the linear form containing any terms that arise solely from an IMEX scheme.

This function should only include the terms specific to the IMEX scheme. E. g. for INS only the explicit form of the convection term would be included in this function. The other linear form terms like the body force and boundary condition terms (which exist regardless of the linearization scheme) would be handled by construct_linear.

Parameters
  • V (List[ProxyFunction]) – A list of test (weighting) functions for the model’s finite element space.

  • gfu_0 (Optional[List[GridFunction], None]) – List of gridfunction components from previous time steps.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[LinearForm]

Returns

List of the explicit IMEX terms to add to the linear form for the model. The list will contain multiple items if the model involves iterating between different versions of a linearized weak form.

abstract construct_linear(V, gfu_0, dt, time_step)[source]

Function to construct the linear form.

DO NOT include any terms that arise from and an IMEX scheme, see construct_imex_explicit for those terms.

Parameters
  • V (List[ProxyFunction]) – A list of test (weighting) functions for the model’s finite element space.

  • gfu_0 (Optional[List[GridFunction], None]) – List of gridfunction components from previous time steps.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[LinearForm]

Returns

List of the linear form for the model. The list will contain multiple linear forms if the model involves iterating between different versions of a linearized weak form.

construct_preconditioners(a_assembled)[source]

Function to construct the preconditioners needed by the model.

Parameters

a_assembled (List[BilinearForm]) – A list of the assembled bilinear form.

Return type

List[Preconditioner]

Returns

A list of the constructed preconditioners.

get_trial_and_test_functions()[source]

Function return the trial and test (weighting) function(s) for the model.

Return type

Tuple[List[ProxyFunction], List[ProxyFunction]]

Returns

Returns two lists, one for the trial and one for the test (weighting) function(s).

load_mesh_fes(mesh=True, fes=True)[source]

Function to load the model’s mesh.

Parameters
  • mesh (bool) – If True, reload the original mesh from file.

  • fes (bool) – If True, reconstruct the finite element space.

abstract solve_single_step(a_lst, L_lst, precond_lst, gfu, time_step=0)[source]

Function to solve the model.

This function produces the model solution after one complete time step (or just the steady state solution for a stationary solve). It includes an inner iterations needed for ex: Picard iterations with Oseen linearization or iteration between weak forms for a model with multiple weak forms.

Parameters
  • a_lst (List[BilinearForm]) – A list of the bilinear forms.

  • L_lst (List[LinearForm]) – A list of the linear forms.

  • precond_lst (List[Preconditioner]) – A list of preconditioners to use.

  • gfu (GridFunction) – The gridfunction to store the solution from the current iteration.

  • time_step (int) – What time step values to use if _apply_dirichlet_bcs_to must be called.

Return type

None

abstract solve_stationary_iteration(a_lst, L_lst, precond_lst, gfu)[source]

nasser@nasser: TODO

Return type

None

time_derivative_terms(gfu_lst, scheme, step=1)[source]

Function to produce the time derivative terms for the linear and bilinear forms.

Parameters
  • gfu_lst (List[List[GridFunction]]) – List of the solutions of previous time steps in reverse chronological order.

  • scheme (str) – The name of the time integration scheme being used.

  • step (int) – Which intermediate step the time derivatives are needed for. This is specific to Runge Kutta schemes.

Returns

  • a: The time derivative terms for the bilinear forms.

  • L: The time derivative terms for the linear forms.

Return type

Tuple[CoefficientFunction, CoefficientFunction]

update_bcs(bc_dict_patch)[source]

Function to update BCs to arbitrary values.

This function is used to implement controllers. It lets them manipulate the manipulated variable.

Parameters

bc_dict_patch (Dict[str, Dict[str, Dict[str, List[Optional[CoefficientFunction, None]]]]]) – Dictionary containing new values for the BCs being updated.

Return type

None

update_linearization_terms(gfu)[source]

Function to update the values in the linearization terms.

Required for Runge Kutta-type solvers if the model’s weak form involves linearization terms. Assumes that linearization terms (self.W) have already been created for the model.

Parameters

gfu (GridFunction) – The gridfunction to use to update the linearization terms.

Return type

None

update_model_variables(updated_gfu, ic_update=False, ref_sol_update=False, time_step=None)[source]

Function to update the values of the model variables and then re-parse any expressions that contain those model variables in ex: the model functions, the model boundary conditions…

Parameters
  • updated_gfu (Union[GridFunction, List[ProxyFunction]]) – Gridfunction containing updated values for all model variables.

  • ic_update (bool) – If True update the initial conditions.

  • ref_sol_update (bool) – If True update the reference solutions.

  • time_step (Optional[int, None]) – If specified, only update the model component values for one particular time step.

Return type

None

update_timestep(gfu, gfu_0)[source]

Function to update the previous time-step gridfunction with result of the current time-step.

Parameters
  • gfu (GridFunction) – The gridfunction containing the result of the current time-step.

  • gfu_0 (GridFunction) – The gridfunction containing the result of the previous time-step.

Return type

None

opencmp.models.ins module

class opencmp.models.ins.INS(config, t_param)[source]

Bases: opencmp.models.base_model.Model

A single phase incompressible Navier-Stokes model.

_contruct_fes_helper()[source]

Helper function for creating the FESpace.

This function exists in order to allow multicomponent_ins (and any further additions) to simplify their FES creation by being able to use this helper function to create the finite element space for velocity and pressure.

Return type

List[FESpace]

Returns

A list containing individual finite element spaces

_get_wind(gfu, time_step)[source]

Function to obtain the wind term used to linearize the convection term.

Parameters
  • gfu (Union[List[ProxyFunction], List[GridFunction]]) – List of

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

Optional[ProxyFunction, None]

Returns

The wind/linearization term used for linearizing the convection term.

accept_stationary_solution(gfu)[source]

nasser@nasser: TODO

static allows_explicit_schemes()[source]

Function to specify whether a given model works with explicit time integration schemes.

Return type

bool

Returns

True if the model can be used with fully explicit time integration schemes, else False.

construct_bilinear_time_ODE(U, V, dt=<ngsolve.fem.Parameter object>, time_step=0)[source]

Function to construct a portion of the bilinear form for model variables WITH time derivatives (basically any bilinear form terms not already in construct_bilinear_time_coefficient).

A given model with multiple model variables may not include time derivatives of all of its model variables. It is necessary to split the bilinear form into two portions as some time discretization schemes use different coefficients for terms with and without time derivatives.

Parameters
  • U (Union[List[ProxyFunction], List[GridFunction]]) – A list of trial functions for the model’s finite element space.

  • V (List[ProxyFunction]) – A list of testing (weighting) functions for the model’s finite element space.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[BilinearForm]

Returns

List of the described portion of the bilinear form for the model. The list will contain multiple bilinear forms if the model involves iterating between different versions of a linearized weak form.

construct_bilinear_time_coefficient(U, V, dt, time_step)[source]

Function to construct a portion of the bilinear form for model variables WITHOUT time derivatives.

A given model with multiple model variables may not include time derivatives of all of its model variables. It is necessary to split the bilinear form into two portions as some time discretization schemes use different coefficients for terms with and without time derivatives.

Parameters
  • U (List[ProxyFunction]) – A list of trial functions for the model’s finite element space.

  • V (List[ProxyFunction]) – A list of testing (weighting) functions for the model’s finite element space.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[BilinearForm]

Returns

List of the described portion of the bilinear form for the model. The list will contain multiple bilinear forms if the model involves iterating between different versions of a linearized weak form.

construct_imex_explicit(V, gfu_0, dt, time_step)[source]

Function to construct the linear form containing any terms that arise solely from an IMEX scheme.

This function should only include the terms specific to the IMEX scheme. E. g. for INS only the explicit form of the convection term would be included in this function. The other linear form terms like the body force and boundary condition terms (which exist regardless of the linearization scheme) would be handled by construct_linear.

Parameters
  • V (List[ProxyFunction]) – A list of test (weighting) functions for the model’s finite element space.

  • gfu_0 (Optional[List[GridFunction], None]) – List of gridfunction components from previous time steps.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[LinearForm]

Returns

List of the explicit IMEX terms to add to the linear form for the model. The list will contain multiple items if the model involves iterating between different versions of a linearized weak form.

construct_linear(V, gfu_0, dt, time_step)[source]

Function to construct the linear form.

DO NOT include any terms that arise from and an IMEX scheme, see construct_imex_explicit for those terms.

Parameters
  • V (List[ProxyFunction]) – A list of test (weighting) functions for the model’s finite element space.

  • gfu_0 (Optional[List[GridFunction], None]) – List of gridfunction components from previous time steps.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[LinearForm]

Returns

List of the linear form for the model. The list will contain multiple linear forms if the model involves iterating between different versions of a linearized weak form.

solve_single_step(a_lst, L_lst, precond_lst, gfu, time_step=0)[source]

Function to solve the model.

This function produces the model solution after one complete time step (or just the steady state solution for a stationary solve). It includes an inner iterations needed for ex: Picard iterations with Oseen linearization or iteration between weak forms for a model with multiple weak forms.

Parameters
  • a_lst (List[BilinearForm]) – A list of the bilinear forms.

  • L_lst (List[LinearForm]) – A list of the linear forms.

  • precond_lst (List[Preconditioner]) – A list of preconditioners to use.

  • gfu (GridFunction) – The gridfunction to store the solution from the current iteration.

  • time_step (int) – What time step values to use if _apply_dirichlet_bcs_to must be called.

Return type

None

solve_stationary_iteration(a_lst, L_lst, precond_lst, gfu)[source]

nasser@nasser: TODO

Return type

None

update_linearization_terms(gfu)[source]

Function to update the values in the linearization terms.

Required for Runge Kutta-type solvers if the model’s weak form involves linearization terms. Assumes that linearization terms (self.W) have already been created for the model.

Parameters

gfu (GridFunction) – The gridfunction to use to update the linearization terms.

Return type

None

opencmp.models.misc module

opencmp.models.misc.get_model_class(model_name, dim_used)[source]

Function to find the correct model to to find, initialize, and return an instance of the desired model class(es).

Will raise a ValueError if a model could not be found.

Parameters
  • model_name (str) – String representing the model to use, must be the name of the model class.

  • dim_used (bool) – bool indicating if the Diffuse Interface Method implemetation of a class should be used.

Return type

Type[Model]

Returns

An initialized instance of the desired model.

opencmp.models.multi_component_ins module

class opencmp.models.multi_component_ins.MultiComponentINS(config, t_param)[source]

Bases: opencmp.models.ins.INS

A single phase multicomponent incompressible Navier-Stokes model.

_add_multiple_components(config)[source]

Function to add multiple components to the model from the config file.

Parameters

config (ConfigParser) – The ConfigParser

Return type

None

construct_bilinear_time_ODE(U, V, dt=<ngsolve.fem.Parameter object>, time_step=0)[source]

Function to construct a portion of the bilinear form for model variables WITH time derivatives (basically any bilinear form terms not already in construct_bilinear_time_coefficient).

A given model with multiple model variables may not include time derivatives of all of its model variables. It is necessary to split the bilinear form into two portions as some time discretization schemes use different coefficients for terms with and without time derivatives.

Parameters
  • U (Union[List[ProxyFunction], List[GridFunction]]) – A list of trial functions for the model’s finite element space.

  • V (List[ProxyFunction]) – A list of testing (weighting) functions for the model’s finite element space.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[BilinearForm]

Returns

List of the described portion of the bilinear form for the model. The list will contain multiple bilinear forms if the model involves iterating between different versions of a linearized weak form.

construct_bilinear_time_coefficient(U, V, dt, time_step)[source]

Function to construct a portion of the bilinear form for model variables WITHOUT time derivatives.

A given model with multiple model variables may not include time derivatives of all of its model variables. It is necessary to split the bilinear form into two portions as some time discretization schemes use different coefficients for terms with and without time derivatives.

Parameters
  • U (List[ProxyFunction]) – A list of trial functions for the model’s finite element space.

  • V (List[ProxyFunction]) – A list of testing (weighting) functions for the model’s finite element space.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[BilinearForm]

Returns

List of the described portion of the bilinear form for the model. The list will contain multiple bilinear forms if the model involves iterating between different versions of a linearized weak form.

construct_imex_explicit(V, gfu_0, dt, time_step)[source]

Function to construct the linear form containing any terms that arise solely from an IMEX scheme.

This function should only include the terms specific to the IMEX scheme. E. g. for INS only the explicit form of the convection term would be included in this function. The other linear form terms like the body force and boundary condition terms (which exist regardless of the linearization scheme) would be handled by construct_linear.

Parameters
  • V (List[ProxyFunction]) – A list of test (weighting) functions for the model’s finite element space.

  • gfu_0 (Optional[List[GridFunction], None]) – List of gridfunction components from previous time steps.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[LinearForm]

Returns

List of the explicit IMEX terms to add to the linear form for the model. The list will contain multiple items if the model involves iterating between different versions of a linearized weak form.

construct_linear(V, gfu_0, dt, time_step)[source]

Function to construct the linear form.

DO NOT include any terms that arise from and an IMEX scheme, see construct_imex_explicit for those terms.

Parameters
  • V (List[ProxyFunction]) – A list of test (weighting) functions for the model’s finite element space.

  • gfu_0 (Optional[List[GridFunction], None]) – List of gridfunction components from previous time steps.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[LinearForm]

Returns

List of the linear form for the model. The list will contain multiple linear forms if the model involves iterating between different versions of a linearized weak form.

opencmp.models.poisson module

class opencmp.models.poisson.Poisson(config, t_param)[source]

Bases: opencmp.models.base_model.Model

A poisson model.

accept_stationary_solution(gfu)[source]

nasser@nasser: TODO

static allows_explicit_schemes()[source]

Function to specify whether a given model works with explicit time integration schemes.

Return type

bool

Returns

True if the model can be used with fully explicit time integration schemes, else False.

construct_bilinear_time_ODE(U, V, dt=<ngsolve.fem.Parameter object>, time_step=0)[source]

Function to construct a portion of the bilinear form for model variables WITH time derivatives (basically any bilinear form terms not already in construct_bilinear_time_coefficient).

A given model with multiple model variables may not include time derivatives of all of its model variables. It is necessary to split the bilinear form into two portions as some time discretization schemes use different coefficients for terms with and without time derivatives.

Parameters
  • U (Union[List[ProxyFunction], List[GridFunction]]) – A list of trial functions for the model’s finite element space.

  • V (List[ProxyFunction]) – A list of testing (weighting) functions for the model’s finite element space.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[BilinearForm]

Returns

List of the described portion of the bilinear form for the model. The list will contain multiple bilinear forms if the model involves iterating between different versions of a linearized weak form.

construct_bilinear_time_coefficient(U, V, dt, time_step)[source]

Function to construct a portion of the bilinear form for model variables WITHOUT time derivatives.

A given model with multiple model variables may not include time derivatives of all of its model variables. It is necessary to split the bilinear form into two portions as some time discretization schemes use different coefficients for terms with and without time derivatives.

Parameters
  • U (List[ProxyFunction]) – A list of trial functions for the model’s finite element space.

  • V (List[ProxyFunction]) – A list of testing (weighting) functions for the model’s finite element space.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[BilinearForm]

Returns

List of the described portion of the bilinear form for the model. The list will contain multiple bilinear forms if the model involves iterating between different versions of a linearized weak form.

construct_imex_explicit(V, gfu_0, dt, time_step)[source]

Function to construct the linear form containing any terms that arise solely from an IMEX scheme.

This function should only include the terms specific to the IMEX scheme. E. g. for INS only the explicit form of the convection term would be included in this function. The other linear form terms like the body force and boundary condition terms (which exist regardless of the linearization scheme) would be handled by construct_linear.

Parameters
  • V (List[ProxyFunction]) – A list of test (weighting) functions for the model’s finite element space.

  • gfu_0 (Optional[List[GridFunction], None]) – List of gridfunction components from previous time steps.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[LinearForm]

Returns

List of the explicit IMEX terms to add to the linear form for the model. The list will contain multiple items if the model involves iterating between different versions of a linearized weak form.

construct_linear(V, gfu_0, dt, time_step)[source]

Function to construct the linear form.

DO NOT include any terms that arise from and an IMEX scheme, see construct_imex_explicit for those terms.

Parameters
  • V (List[ProxyFunction]) – A list of test (weighting) functions for the model’s finite element space.

  • gfu_0 (Optional[List[GridFunction], None]) – List of gridfunction components from previous time steps.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[LinearForm]

Returns

List of the linear form for the model. The list will contain multiple linear forms if the model involves iterating between different versions of a linearized weak form.

solve_single_step(a_lst, L_lst, precond_lst, gfu, time_step=0)[source]

Function to solve the model.

This function produces the model solution after one complete time step (or just the steady state solution for a stationary solve). It includes an inner iterations needed for ex: Picard iterations with Oseen linearization or iteration between weak forms for a model with multiple weak forms.

Parameters
  • a_lst (List[BilinearForm]) – A list of the bilinear forms.

  • L_lst (List[LinearForm]) – A list of the linear forms.

  • precond_lst (List[Preconditioner]) – A list of preconditioners to use.

  • gfu (GridFunction) – The gridfunction to store the solution from the current iteration.

  • time_step (int) – What time step values to use if _apply_dirichlet_bcs_to must be called.

Return type

None

solve_stationary_iteration(a_lst, L_lst, precond_lst, gfu)[source]

nasser@nasser: TODO

Return type

None

opencmp.models.stokes module

class opencmp.models.stokes.Stokes(config, t_param)[source]

Bases: opencmp.models.ins.INS

A single phase Stokes model.

accept_stationary_solution(gfu)[source]

nasser@nasser: TODO

construct_bilinear_time_ODE(U, V, dt=<ngsolve.fem.Parameter object>, time_step=0)[source]

Function to construct a portion of the bilinear form for model variables WITH time derivatives (basically any bilinear form terms not already in construct_bilinear_time_coefficient).

A given model with multiple model variables may not include time derivatives of all of its model variables. It is necessary to split the bilinear form into two portions as some time discretization schemes use different coefficients for terms with and without time derivatives.

Parameters
  • U (Union[List[ProxyFunction], List[GridFunction]]) – A list of trial functions for the model’s finite element space.

  • V (List[ProxyFunction]) – A list of testing (weighting) functions for the model’s finite element space.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[BilinearForm]

Returns

List of the described portion of the bilinear form for the model. The list will contain multiple bilinear forms if the model involves iterating between different versions of a linearized weak form.

construct_bilinear_time_coefficient(U, V, dt=<ngsolve.fem.Parameter object>, time_step=0)[source]

Function to construct a portion of the bilinear form for model variables WITHOUT time derivatives.

A given model with multiple model variables may not include time derivatives of all of its model variables. It is necessary to split the bilinear form into two portions as some time discretization schemes use different coefficients for terms with and without time derivatives.

Parameters
  • U (List[ProxyFunction]) – A list of trial functions for the model’s finite element space.

  • V (List[ProxyFunction]) – A list of testing (weighting) functions for the model’s finite element space.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[BilinearForm]

Returns

List of the described portion of the bilinear form for the model. The list will contain multiple bilinear forms if the model involves iterating between different versions of a linearized weak form.

construct_imex_explicit(V, gfu_0=None, dt=<ngsolve.fem.Parameter object>, time_step=0)[source]

Function to construct the linear form containing any terms that arise solely from an IMEX scheme.

This function should only include the terms specific to the IMEX scheme. E. g. for INS only the explicit form of the convection term would be included in this function. The other linear form terms like the body force and boundary condition terms (which exist regardless of the linearization scheme) would be handled by construct_linear.

Parameters
  • V (List[ProxyFunction]) – A list of test (weighting) functions for the model’s finite element space.

  • gfu_0 (Optional[List[GridFunction], None]) – List of gridfunction components from previous time steps.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[LinearForm]

Returns

List of the explicit IMEX terms to add to the linear form for the model. The list will contain multiple items if the model involves iterating between different versions of a linearized weak form.

construct_linear(V, gfu_0=None, dt=<ngsolve.fem.Parameter object>, time_step=0)[source]

Function to construct the linear form.

DO NOT include any terms that arise from and an IMEX scheme, see construct_imex_explicit for those terms.

Parameters
  • V (List[ProxyFunction]) – A list of test (weighting) functions for the model’s finite element space.

  • gfu_0 (Optional[List[GridFunction], None]) – List of gridfunction components from previous time steps.

  • dt (Parameter) – Time step parameter (in the case of a stationary solve is just one).

  • time_step (int) – What time step values to use for ex: boundary conditions. The value corresponds to the index of the time step in t_param and dt_param.

Return type

List[LinearForm]

Returns

List of the linear form for the model. The list will contain multiple linear forms if the model involves iterating between different versions of a linearized weak form.

solve_single_step(a_lst, L_lst, precond_lst, gfu, time_step=0)[source]

Function to solve the model.

This function produces the model solution after one complete time step (or just the steady state solution for a stationary solve). It includes an inner iterations needed for ex: Picard iterations with Oseen linearization or iteration between weak forms for a model with multiple weak forms.

Parameters
  • a_lst (List[BilinearForm]) – A list of the bilinear forms.

  • L_lst (List[LinearForm]) – A list of the linear forms.

  • precond_lst (List[Preconditioner]) – A list of preconditioners to use.

  • gfu (GridFunction) – The gridfunction to store the solution from the current iteration.

  • time_step (int) – What time step values to use if _apply_dirichlet_bcs_to must be called.

Return type

None

solve_stationary_iteration(a_lst, L_lst, precond_lst, gfu)[source]

nasser@nasser: TODO

Return type

None

Module contents