Code documentation
This part of the documentation is generated automatically from the python source using SPHINX.
Folders
Distributed port-Hamiltonian system
file: dphs.py
authors: Giuseppe Ferraro, Ghislain Haine, Florian Monteghetti
date: 22 nov. 2022
brief: class for distributed port-Hamiltonian system
- class scrimp.dphs.DPHS(basis_field='real')
Bases:
object
A generic class handling distributed pHs using the GetFEM tools
This is a wrapper in order to simplify the coding process
Access to fine tunings is preserved as much as possible
- F
A PETSc Vec for residual computation
- IFunction(TS, t, z, zd, F)
IFunction for the time-resolution of the dphs with PETSc TS fully implicit integrator
- Args:
TS (PETSc TS): the PETSc TS object handling time-resolution t (float): time parameter z (PETSc Vec): the state zd (PETSc Vec): the time-derivative of the state F (PETSc Vec): the rhs vector
- IJacobian(TS, t, z, zd, sig, A, P)
IJacobian for the time-resolution of the dphs with PETSc TS fully implicit integrator
- Args:
TS (PETSc TS): the PETSc TS object handling time-resolution t (float): time parameter z (PETSc Vec): the state zd (PETSc Vec): the time-derivative of the state sig (float): a shift-parameter, depends on dt A (PETSc Mat): the jacobian matrix P:(PETSc Mat) the jacobian matrix to use for pre-conditionning
- J
A PETSc Mat for Jacobian computation
- add_FEM(fem: FEM)
This function adds a FEM (Finite Element Method) for the variables associated to a port of the dphs
TODO: handle tensor-field
- Args:
fem (FEM): the FEM to use
- add_brick(brick: Brick)
This function adds a brick in the getfem Model thanks to a form in GWFL getfem language
The form may be non-linear
- Args:
brick (Brick): the brick
- add_control_port(control_port: Control_Port)
This function adds a control port to the dphs
- Args:
control_port (Control_Port): the Control Port.
- add_costate(costate: CoState)
This function adds a costate to the costate dict of the dphs, then defines and adds the dynamical port gathering the couple (State, CoState).
- Args:
costate (CoState): the costate
- add_parameter(parameter: Parameter)
This function adds a time-independent (possibly space-varying parameter: x, y and z are the space variables) associated to a port of the dphs
- Args:
parameter (Parameter): the parameter
- add_port(port: Port)
This function adds a port to the port dict of the dphs
- Args:
port (Port): the port
- add_state(state: State)
This functions adds a state to state dict of the dphs.
- Args:
state (State): the state
- allocate_memory()
Pre-allocate memory for matrices and vectors
- assemble_mass()
This function performs the assembly of the bricks dt=True and linear=True and set the PETSc.Mat attribute mass
- assemble_nl_mass()
This function performs the assembly of the bricks dt=True and linear=False and set the PETSc.Mat attribute nl_mass
- assemble_nl_stiffness()
This function performs the assembly of the bricks dt=False and linear=False and set the PETSc.Mat attribute nl_stiffness
- assemble_rhs()
This function performs the assembly of the rhs and set the PETSc.Vec attribute rhs
- assemble_stiffness()
This function performs the assembly of the bricks dt=False and linear=True and set the PETSc.Mat attribute stiffness
- bricks
The dict of bricks, associating getfem bricks to petsc matrices obtained by the PFEM, store many infos for display()
- buffer
A PETSc Vec buffering computation
- compute_Hamiltonian()
Compute each term constituting the Hamiltonian
Wrapper to the function compute of Hamiltonian class
- compute_powers()
Compute each power associated to each algebraic port if it is not substituted (because of the parameter-dependency if it is)
Wrapper to the function compute of Port class (loop on self.ports)
- controls
The dict of controls, collecting information about control ports. Also appear in ports entries
- costates
The dict of costates, store many infos for display()
- disable_all_bricks()
This function disables all bricks in the Model
- display(verbose=2)
A method giving infos about the dphs
- Args:
verbose (int): the level of verbosity (defaults to 2)
- domain
The domain of a dphs is an object that handle mesh(es) and dict of regions with getfem indices (for each mesh), useful to define bricks (i.e. forms) in the getfem syntax
- enable_all_bricks()
This function enables all bricks in the Model
- event(TS, t, z, fvalue)
Check if the time step is not too small
- exclude_algebraic_var_from_lte(TS)
Exclude the algebraic variable from the local error troncature in the time-resolution
- Args:
TS (PETSc TS): the PETSc TS object handling time-resolution
- export_matrices(t=None, state=None, path=None, to='matlab')
TODO:
- export_to_pv(name_variable, path=None, t='All')
Export the solution to .vtu file(s) (for ParaView), with associated .pvd if t=’All’
name_variable (str): the variable to export path (str): the path for the output file (default: in the outputs/pv folder next to your .py file) t (Numpy array): the time values of extraction (default: All the stored times), All, Init, Final
- get_Hamiltonian()
This function returns the Hamiltonian in function of time
- Returns:
numpy array: time array of the values of the Hamiltonian.
- get_cleared_TS_options()
To ensure a safe database for the PETSc TS environment
- get_quantity(expression, region=-1, order=0, CN=False, mesh_id=0) list
This functions computes the integral over region of the expression at each time step
- Args:
expression (str): the GFWL expression to compute
region (int, optional): the id of the region (defaults to -1)
order (int, optional): the order of the quantity to be computed (0: scalar, 1: vector, 2: tensor) (defaults to 0)
CN (bool, optional): take half the sum 0.5*(z^{n+1}+z^n) instead of z^n in the computation (defaults to False)
mesh_id (int, optional): the id of the mesh (defaults to 0)
- Returns:
list(float): a list of float at each time step (according to self.solution[“t”])
- get_solution(name_variable) list
This functions is useful (especially in 1D) to handle post-processing by extracting the variable of interest from the PETSc Vec self.solution[“z”]
- Args:
name_variable (str): the name of the variable to extract
- Returns:
list(numpy array): a list of numpy array (dofs of name_variable) at each time step (according to self.solution[“t”])
- gf_model
A getfem Model object that is use as core for the dphs
- hamiltonian
The Hamiltonian of a dphs is a list of dict containing several useful information for each term
- init_matrix(A)
Set all diagonal entries of A without modifying existing entries.
Parameters
A : PETSc.Mat
- init_parameter(name, name_port)
This function initializes the parameter name in the FEM of the port name_port of the dphs and adds it to the getfem Model
- Args:
name (str): the name of the parameter as defined with add_parameter() name_port (str): the name of the port where the parameter belongs
- init_step()
Perform a first initial step with a pseudo bdf
It needs set_time_scheme(init_step=True)
- initial_value_set
To check if the initial values have been set before time-resolution
- linear_mass
To speed-up IFunction calls for linear systems
- mass
Linear mass matrix of the system in PETSc CSR format
- monitor(TS, i, t, z, dt_save=1.0, t_0=0.0, initial_step=False)
Monitor to use during iterations of time-integration at each successful time step
- Args:
TS (PETSc TS): the PETSc TS object handling time-resolution i (int): the iteration in the time-resolution t (float): time parameter z (PETSc Vec): the state dt_save (float): save the solution each dt_save s (different from the time-step dt used for resolution) t_0 (float): the initial time initial_step (bool): True if this is the initial consistency step (default=`False`)
- nl_mass
Non-linear mass matrix of the system in PETSc CSR format
- nl_stiffness
Non-linear stiffness matrix of the system in PETSc CSR format
- plot_Hamiltonian(with_powers=True, save_figure=False, filename='Hamiltonian.png')
Plot each term constituting the Hamiltonian and the Hamiltonian
May include the power terms, i.e. the sum over [t_0, t_f] of the flow/effort product of algebraic ports
- Args:
with_powers (bool): if True (default), the plot will also contains the power of each algebraic ports
save_figure (bool): if ‘True’ (defaults: False), save the plot
filename (str): the name of the file where the plot is saved (defaults: Hamiltonian.png)
- plot_powers(ax=None, HamTot=None)
Plot each power associated to each algebraic port
The time integration for visual comparison with the Hamiltonian is done using a midpoint method
If HamTot is provided, a Balance showing structure-preserving property is shown: must be constant on the plot
ax (Matplotlib axis): the gca of matplotlib when this function is called from plot_Hamiltonian() HamTot (Numpy array): the values of the Hamiltonian over time
- ports
The dict of ports, store many infos for display()
- postevent(TS, event, t, z, forward)
If the time step is too small, ask for the end of the simulation
- powers_computed
To check if the powers have been computed
- rhs
rhs of the system in PETSc Vec
- set_control(name, expression)
This function applies a source term expression to the control port name
- Args:
name (str): the name of the port expression (str): the expression of the source term
- set_domain(domain: Domain)
This function sets a domain for the dphs.
TODO: If not built_in, given from a script ‘name.py’ or a .geo file with args in the dict ‘parameters’ should be able to handle several meshes e.g. for interconnections, hence the list type
- Args:
name (str): id of the domain, either for built in, or user-defined auxiliary script parameters (dict): parameters for the construction, either for built in, or user-defined auxiliary script
- set_from_vector(name_variable, x)
This function sets the value of the variable name_variable in the getfem Model from a numpy vector
- Args:
name_variable (str): the name of the variable to set x (numpy array): the vector of values
- set_initial_value(name_variable, expression)
This function sets the initial value of the variable name_variable of the dphs from an expression
- Args:
name_variable (str): the name of the variable to set expression (str): the expression of the function to use
- set_linear_flags()
This function set the “linear flags” to speed-up IFunction calls for linear systems.
- set_time_scheme(**kwargs)
Allows an easy setting of the PETSc TS environment
- Args:
**kwargs: PETSc TS options and more (see examples)
- solution
Will contain both time t and solution z
- solve()
Perform the time-resolution of the dphs thanks to PETSc TS
The options database is set in the time_scheme attribute
- solve_done
To check if the system has been solved
- spy_Dirac(t=None, state=None)
!TO DO: improve a lot with position of bricks and no constitutive relations!!!
- states
The dict of states, store many infos for display()
- stiffness
Linear stiffness matrix of the system in PETSc CSR format
- stop_TS
To stop TS integration by keeping already computed timesteps if one step fails
- tangent_mass
Tangent (non-linear + linear) mass matrix of the system in PETSc CSR format
- tangent_stiffness
Tangent (non-linear + linear) stiffness matrix of the system in PETSc CSR format
- ts_start
For monitoring time in TS resolution
Domain
file: domain.py
authors: Giuseppe Ferraro, Ghislain Haine
date: 31 may 2023
brief: class for domain object
- class scrimp.domain.Domain(name: str, parameters: dict, refine=0, terminal=1)
Bases:
object
A class handling meshes and indices for regions
Lists are used to handle interconnection of pHs, allowing for several meshes in the dpHs.
- display()
A method giving infos about the domain
- get_boundaries() list
This function gets the list of the bounderies for the domain.
- Returns:
list: list of the boundaries for the domain
- get_dim() list
This function gets the list of dimensions for the domain.
- Returns:
list: list of dimensions for the domain
- get_isSet() bool
This function gets the boolean vale indicating wether if a Mesh has been set for the domain or not.
- Returns:
bool: boolean indicating if a Mesh has been set
- get_mesh() list
This function gets the list of mesh for the domain.
- Returns:
list: list of mesh for the domain
- get_name() str
This function get the name of the domain.
- Returns:
str: name of the domain.
- get_subdomains() list
This function gets the list of subdomains for the domain.
- Returns:
list: list of the subdomains for the domain
- set_mim_auto()
Define the integration method to a default choice
Hamiltonian / Term
file: hamiltonian.py
authors: Giuseppe Ferraro, Ghislain Haine
date: 31 may 2023
brief: class for hamiltonian and term objects
- class scrimp.hamiltonian.Hamiltonian(name: str)
Bases:
object
This class defines the Hamiltoninan.
- add_term(term: Term)
This function adds a term to the term list of the Hamiltonian
- Args:
term (Term): term for the Hamiltonian
- compute(solution: dict, gf_model: getfem.Model, domain: Domain)
Compute each term constituting the Hamiltonian
- Args:
solutions (dict): The solution of the dphs
gf_model (GetFEM Model): The model getfem of the dphs
domain (Domain): The domain of the dphs
- get_is_computed() bool
This function returns True if the Hamiltonian is computed, False otherwise.
- Returns:
bool: flag that indicates if the hamiltonian terms have been computed
- get_name() str
This function returns the name of the Hamiltonion.
- Args:
name (str): name of the Hamiltonian
- get_terms() list
This function returns a copy of the list of terms of the Hamiltonian.
- Returns:
list: list of terms of the Hamiltonian.
- set_is_computed()
This function sets the Hamiltonian as computed.
- set_name(name: str)
This function set the name for the Hamiltonion. For plotting purposes.
- Args:
name (str): name of the Hamiltonian
- class scrimp.hamiltonian.Term(description: str, expression: str, regions: str, mesh_id: int = 0)
Bases:
object
This class defines a term for the Hamiltoninan.
- get_description() str
This function gets the description of the term.
- Returns:
str: description of the term
- get_expression() str
This function gets the matematical expression of the term.
- Returns:
str: matematical expression of the term.
- get_mesh_id() int
This function gets the mesh id of the mesh where the regions belong to..
- Returns:
int: the mesh id of the mesh where the regions belong to.
- get_regions() str
This function gets the regions of the term.
- Returns:
str: regions of the termthe region IDs of the mesh where the expression has to be evaluated
- get_values() list
This function gets the valeus of the term.
- Returns:
list: list of values of the term
- set_value(value)
This function sets a value for the term.
- Args:
value : a value for the term
Port / Parameter
file: port.py
authors: Giuseppe Ferraro, Ghislain Haine
date: 31 may 2023
brief: class for port and parameter objects
- class scrimp.port.Parameter(name: str, description: str, kind: str, expression: str, name_port: str)
Bases:
object
This class describes the Parameter for a Port.
- get_description() str
This function gets the description of the parameter.
- Returns:
str: description of the parameter
- get_expression() str
This function gets the matematical expression of the parameter.
- Returns:
str: matematical expression of the parameter.
- get_kind() str
This function gets the kind of the parameter.
- Returns:
str: kind of the parameter
- get_name() str
This function gets the name of the parameter.
- Returns:
str: name of the parameter
- get_name_port() str
This function gets the name of the port whose the parameter is bounded.
- Returns:
str: name of the port whose the parameter is bounded.
- class scrimp.port.Port(name: str, flow: str, effort: str, kind: str, mesh_id: int = 0, algebraic: bool = True, substituted: bool = False, dissipative: bool = True, region: int = None)
Bases:
object
A class to handle a port of a dpHs
It is mainly constituted of a flow variable, an effort variable, and a fem.
- add_parameter(parameter: Parameter) bool
This function adds a Parameter object that is acting on the variables of the port.
- Args:
parameter (Parameter): parameter for the port.:
- Returns:
bool: True if the insertion has been complete correctly, False otherwise
- compute(solution: dict, gf_model: getfem.Model, domain: Domain)
Compute the power flowing through the algebraic port if it is not substituted (because of the parameter-dependency if it is)
- Args:
solutions (dict): The solution of the dphs
gf_model (GetFEM Model): The model getfem of the dphs
domain (Domain): The domain of the dphs
- get_algebraic() bool
- This function gets boolean value of the algebraic parameter of the port.
If True, the equation associated to this port is algebraic, otherwise dynamic and the flow is derivated in time at resoltuion
- Returns:
bool: value of the algebraic parameter
- get_dissipative() bool
This function gets boolean value for the dissipativeness flag of the port. If True, the power associated to the port gets a negative sign
- Returns:
bool: value of the dissipativeness flag
- get_effort() str
This function gets the name of the effort variable.
- Returns:
str: name of the effort variable
- get_fem()
This function returns the fetfem Meshfem object to discretize the port.
- Returns:
_type_: the getfem Meshfem object to discretize the port
- get_flow() str
This function gets the name of the flow variable.
- Returns:
str: name of the flow variable
- get_isSet() bool
This funcion gets the boolean value that indicates wether the port is set or not.
- Returns:
bool: value that indicates if the port is set.
- get_is_computed() bool
This function returns True if the power of the Port is computed, False otherwise.
- Returns:
bool: flag that indicates if the power flowing through the Port has been computed
- get_kind() str
This function gets the type of the variables (e.g. scalar-field)
- Returns:
str: type of the variables (e.g. scalar-field)
- get_mesh_id() int
This function gets he id of the mesh where the variables belong
- Returns:
int: The id of the mesh where the variables belong
- get_name() str
This function gets the name of the port.
- Returns:
str: name of the port
- get_parameter(name) Parameter
This function return the parameter with a specific name.
- Args:
name (str): the name of the parameter of interest
- Returns:
Parameter: the desired parameter, None otherwise
- get_parameters() list
This function returns the list of all the parameters inserted for the port.
- Returns:
list(Parameter): list of all the parameters inserted for the port
- get_power()
Gives access to the computed power of the port.
- get_region() int
This function gets the region of the mesh. If any, the int of the region of mesh_id where the flow/effort variables belong
- Returns:
int: region of the mesh
- get_substituted() bool
This function gets boolean value of the subtituted parameter. If True, the getfem `Model will only have an unknown variable for the effort: the constitutive relation is substituted into the mass matrix on the flow side
- Returns:
bool: value of the subtituted parameter
- init_parameter(name: str, expression: str)
This function sets the chosen parameter object for the current port by initialization in the FE basis.
- Args:
name (str): the name of the parameter object expression (str):
- Returns:
out (numpy array): the evaluation of the parameter in the fem of the port.
- set_fem(fem: FEM)
This function sets the Meshfem getfem object defining the finite element method to use to discretize the port.
- Args:
fem (FEM): the FEM object to use
- set_isSet() bool
This funcion sets the boolean value that indicates the port is set.
- Returns:
bool: value that indicates if the port is set.
- set_is_computed()
This function sets the power of the Port as computed.
- set_power(power)
This function sets the power along time of the Port.
State
file: state.py
authors: Giuseppe Ferraro, Ghislain Haine
date: 31 may 2023
brief: class for state object
- class scrimp.state.State(name: str, description: str, kind: str, region: int = None, mesh_id: int = 0)
Bases:
object
This class defines a State.
- get_costate() object
This function gets the Co-state of the state
- Returns:
object: Costate
- get_description() str
This function gets the description of the State.
- Returns:
str: description of the state
- get_kind() str
This function gets the kind of the State.
- Returns:
str: kind of the state
- get_mesh_id() int
This funtion gets the integer number of the mesh of the state.
- Returns:
int: id of the mesh
- get_name() str
This function gets the name of the State.
- Returns:
str: name of the state
- get_port() object
This function gets the port of the state
- Returns:
object: Port
- get_region() int
This function gets the integer number of the region of the state.
- Returns:
int: region of the State.
- set_costate(costate)
This function sets a Co-state to the State.
- Args:
costate (Costate): Co-state
- set_port(port)
This function sets a port to the state.
- Args:
port (Port): Port
Co-state
file: costate.py
authors: Giuseppe Ferraro, Ghislain Haine
date: 31 may 2023
brief: class for co-state object
- class scrimp.costate.CoState(name: str, description: str, state: State, substituted=False)
Bases:
State
This class defines a Co-State.
- get_state() object
This function gets the State of the Costate.
- Returns:
object: State
- get_substituted() bool
This function gets the boolean indicating wether to substitute the variable or not.
- Returns:
bool: boolean indicating wether to substitute the variable
Control
file: control.py
authors: Giuseppe Ferraro, Ghislain Haine
date: 31 may 2023
brief: class for control port object
- class scrimp.control.Control_Port(name: str, name_control: str, description_control: str, name_observation: str, description_observation: str, kind: str, region: int = None, position: str = 'effort', mesh_id: int = 0)
Bases:
Port
This class defines a Control Port.
- get_description_control() str
This function gets the description of the control.
- Returns:
str: the description of the control
- get_description_observation() str
This function gets the description of the observation.
- Returns:
str: the description of the observation
- get_name_control() str
This function gets the name of the control.
- Returns:
str: the name of the control
- get_name_obervation() str
This function gets the name of the obervation.
- Returns:
str: the name of the observation
- get_position() str
This function gets the position of the control in the Dirac structure.
- Returns:
str: the position of the control
FEM
file: fem.py
authors: Giuseppe Ferraro, Ghislain Haine
date: 31 may 2023
brief: class for fem object
- class scrimp.fem.FEM(name, order, FEM='CG')
Bases:
object
This class defines what is a FEM object in SCRIMP.
An negative order allows to access to GetFEM syntax for the FEM, e.g., by setting FEM=”FEM_HERMITE(2)” for Hermite FE in dimension 2.
- get_dim() int
This function gets the dimension of the FEM.
- Returns:
int: the dimension of FEM
- get_fem()
This function returns the the FEM of getfem.
- Returns:
gf.MeshFem: the FEM
- get_isSet() bool
This function gets the flag to know if the FEM are set in getfem.
- Returns:
bool: the flag to assert setting in getfem
- get_mesh()
This function gets the mesh of the FEM
- Returns:
mesh: the mesh of FEM
- get_name() str
This function gets the name of the FEM.
- Returns:
str: name of the FEM
- get_order() int
This function gets the order for the FEM.
- Returns:
int: dim of the flow FEM
- get_type() str
This function gets the tyoe of the FEM.
- Returns:
str: type of the FEM
- set_dim(dim: int)
This function sets the dimension for the FEM.
- Args:
dim (int): the dimension fro the FEM
- set_fem()
This function sets the Meshfem getfem object defining the finite element method to use to discretize the port.
- set_mesh(mesh)
This function sets the Meshfem getfem object FEM object of scrimp.
- Args:
mesh (Mesh): the mesh where the FE are define
Brick
file: brick.py
authors: Giuseppe Ferraro, Ghislain Haine
date: 31 may 2023
brief: class for brick object
- class scrimp.brick.Brick(name: str, form: str, regions: list, linear: bool = True, dt: bool = False, position: str = 'constitutive', explicit: bool = False, mesh_id: int = 0)
Bases:
object
This class defines a Brick.
- add_id_brick_to_list(id_brick: int)
This function adds a brick ID to the brick ID list.
- Args:
id_brick (int): the id of the brick
- disable_id_bricks(gf_model)
This function disable the brick in the getfem model.
- enable_id_bricks(gf_model)
This function enable the brick in the getfem model.
- get_dt() bool
This function returns the boolean that defines wether the matrices applied to time-derivative of a variable or not.
- Returns:
bool: parameter to help easy identification of matrices applied to time-derivative of a variable (e.g. mass matrices).
- get_explicit() bool
This function returns the boolean that defines wether the brick is explicit or not.
- Returns:
bool: parameter to help easy identification of explicit bricks.
- get_form() str
This function returns the form of the brick.
- Returns:
str: the form in GWFL getfem language.
- get_id_bricks() list
This function returns the list of integers related to the ids of the bricks.
- Returns:
list: the list of integers related to the ids of the bricks.
- get_linear() bool
This function returns the boolean that defines wether the brick is linear or not.
- Returns:
bool: parameter to help easy identification of linear bricks.
- get_mesh_id() int
This function returns the ID of the brick.
- Returns:
int: the id of the mesh where the form applies.
- get_name() str
This function returns the name of the brick.
- Returns:
str: the name of the brick, will be used mainly for plotting purpose
- get_position() int
This function returns the id of the position where the form of the brick applies.
- Returns:
int: the id of the mesh where the form applies. Defaults to 0.
- get_regions() list
This function returns the regions of the brick.
- Returns:
list: the regions of mesh where the form applies.