Stocks and Lifetime Models

pydantic model flodym.Stock

Stock objects are components of an MFASystem, where materials can accumulate over time. They consist of three flodym.FlodymArray objects: stock (the accumulation), inflow, outflow.

The base class only allows to compute the stock from known inflow and outflow. The subclasses allows computations using a lifetime distribution function, which is necessary if not both inflow and outflow are known.

Config:
  • protected_namespaces: tuple = ()

  • arbitrary_types_allowed: bool = True

Fields:
  • dims (flodym.dimensions.DimensionSet)

  • inflow (flodym.flodym_arrays.StockArray | None)

  • name (str | None)

  • outflow (flodym.flodym_arrays.StockArray | None)

  • process (flodym.processes.Process | None)

  • stock (flodym.flodym_arrays.StockArray | None)

  • time_letter (str)

Validators:
  • validate_stock_arrays » all fields

  • validate_time_first_dim » all fields

field dims: DimensionSet [Required]

Dimensions of the stock, inflow, and outflow arrays. Time must be the first dimension.

field inflow: StockArray | None = None

Inflow into the stock

field name: str | None = 'unnamed'

Name of the stock

field outflow: StockArray | None = None

Outflow from the stock

field process: Process | None = None

Process the stock is associated with, if any. Needed for example for the mass balance.

field stock: StockArray | None = None

Accumulation of the stock

field time_letter: str = 't'

Letter of the time dimension in the dimensions set, to make sure it’s the first one.

check_stock_balance()
abstract compute()
get_stock_balance()

Check whether inflow, outflow, and stock are balanced. If possible, the method returns the vector ‘Balance’, where Balance = inflow - outflow - stock_change

Return type:

ndarray

to_stock_type(desired_stock_type, **kwargs)

Return an object of a new stock type with values and dimensions the same as the original. **kwargs can be used to pass additional model attributes as required by the desired stock type, if these are not contained in the original stock type.

Parameters:

desired_stock_type (type)

validator validate_stock_arrays  »  all fields
validator validate_time_first_dim  »  all fields
property process_id: int

ID of the process the stock is associated with.

property shape: tuple

Shape of the stock, inflow, outflow arrays, defined by the dimensions.

pydantic model flodym.SimpleFlowDrivenStock

Given inflows and outflows, the stock can be calculated without a lifetime model or cohorts.

Config:
  • protected_namespaces: tuple = ()

  • arbitrary_types_allowed: bool = True

Fields:
  • dims ()

  • inflow ()

  • name ()

  • outflow ()

  • process ()

  • stock ()

  • time_letter ()

Validators:
  • validate_stock_arrays » all fields

  • validate_time_first_dim » all fields

field dims: DimensionSet [Required]

Dimensions of the stock, inflow, and outflow arrays. Time must be the first dimension.

field inflow: StockArray | None = None

Inflow into the stock

field name: str | None = 'unnamed'

Name of the stock

field outflow: StockArray | None = None

Outflow from the stock

field process: Process | None = None

Process the stock is associated with, if any. Needed for example for the mass balance.

field stock: StockArray | None = None

Accumulation of the stock

field time_letter: str = 't'

Letter of the time dimension in the dimensions set, to make sure it’s the first one.

check_stock_balance()
compute()
get_stock_balance()

Check whether inflow, outflow, and stock are balanced. If possible, the method returns the vector ‘Balance’, where Balance = inflow - outflow - stock_change

Return type:

ndarray

to_stock_type(desired_stock_type, **kwargs)

Return an object of a new stock type with values and dimensions the same as the original. **kwargs can be used to pass additional model attributes as required by the desired stock type, if these are not contained in the original stock type.

Parameters:

desired_stock_type (type)

validator validate_stock_arrays  »  all fields
validator validate_time_first_dim  »  all fields
property process_id: int

ID of the process the stock is associated with.

property shape: tuple

Shape of the stock, inflow, outflow arrays, defined by the dimensions.

pydantic model flodym.InflowDrivenDSM

Inflow driven model. Given inflow and lifetime distribution calculate stocks and outflows.

Config:
  • protected_namespaces: tuple = ()

  • arbitrary_types_allowed: bool = True

Fields:
  • dims ()

  • inflow ()

  • lifetime_model ()

  • name ()

  • outflow ()

  • process ()

  • stock ()

  • time_letter ()

Validators:
  • init_cohort_arrays » all fields

  • init_lifetime_model » all fields

  • validate_stock_arrays » all fields

  • validate_time_first_dim » all fields

field dims: DimensionSet [Required]

Dimensions of the stock, inflow, and outflow arrays. Time must be the first dimension.

field inflow: StockArray | None = None

Inflow into the stock

field lifetime_model: LifetimeModel | type [Required]

Lifetime model, which contains the lifetime distribution function. Can be input either as a LifetimeModel subclass, or as an instance of a LifetimeModel subclass. For available subclasses, see flodym.lifetime_models.

field name: str | None = 'unnamed'

Name of the stock

field outflow: StockArray | None = None

Outflow from the stock

field process: Process | None = None

Process the stock is associated with, if any. Needed for example for the mass balance.

field stock: StockArray | None = None

Accumulation of the stock

field time_letter: str = 't'

Letter of the time dimension in the dimensions set, to make sure it’s the first one.

check_stock_balance()
compute()

Determine stocks and outflows and store values in the class instance.

get_outflow_by_cohort()

Outflow by cohort, i.e. the outflow of each production year at each time step.

Return type:

ndarray

get_stock_balance()

Check whether inflow, outflow, and stock are balanced. If possible, the method returns the vector ‘Balance’, where Balance = inflow - outflow - stock_change

Return type:

ndarray

get_stock_by_cohort()

Stock by cohort, i.e. the stock of each production year at each time step.

Return type:

ndarray

validator init_cohort_arrays  »  all fields
validator init_lifetime_model  »  all fields
to_stock_type(desired_stock_type, **kwargs)

Return an object of a new stock type with values and dimensions the same as the original. **kwargs can be used to pass additional model attributes as required by the desired stock type, if these are not contained in the original stock type.

Parameters:

desired_stock_type (type)

validator validate_stock_arrays  »  all fields
validator validate_time_first_dim  »  all fields
property process_id: int

ID of the process the stock is associated with.

property shape: tuple

Shape of the stock, inflow, outflow arrays, defined by the dimensions.

pydantic model flodym.StockDrivenDSM

Stock driven model. Given total stock and lifetime distribution, calculate inflows and outflows. This involves solving the lower triangular equation system A*x=b, where A is the survival function matrix, x is the inflow vector, and b is the stock vector.

Config:
  • protected_namespaces: tuple = ()

  • arbitrary_types_allowed: bool = True

Fields:
  • dims ()

  • inflow ()

  • lifetime_model ()

  • name ()

  • outflow ()

  • process ()

  • solver (str)

  • stock ()

  • time_letter ()

Validators:
  • init_cohort_arrays » all fields

  • init_lifetime_model » all fields

  • init_solver » all fields

  • validate_stock_arrays » all fields

  • validate_time_first_dim » all fields

field dims: DimensionSet [Required]

Dimensions of the stock, inflow, and outflow arrays. Time must be the first dimension.

field inflow: StockArray | None = None

Inflow into the stock

field lifetime_model: LifetimeModel | type [Required]

Lifetime model, which contains the lifetime distribution function. Can be input either as a LifetimeModel subclass, or as an instance of a LifetimeModel subclass. For available subclasses, see flodym.lifetime_models.

field name: str | None = 'unnamed'

Name of the stock

field outflow: StockArray | None = None

Outflow from the stock

field process: Process | None = None

Process the stock is associated with, if any. Needed for example for the mass balance.

field solver: str = 'manual'

Algorithm to use for solving the equation system. Options are: “manual” (default), which uses an own python implementation, and “lapack”, which calls the lapack trtrs routine via scipy. The lapack implementation may be more precise. Speed depends on the dimensionality, but the manual implementation is usually faster.

field stock: StockArray | None = None

Accumulation of the stock

field time_letter: str = 't'

Letter of the time dimension in the dimensions set, to make sure it’s the first one.

check_stock_balance()
compute()

Determine inflows and outflows and store values in the class instance.

get_outflow_by_cohort()

Outflow by cohort, i.e. the outflow of each production year at each time step.

Return type:

ndarray

get_stock_balance()

Check whether inflow, outflow, and stock are balanced. If possible, the method returns the vector ‘Balance’, where Balance = inflow - outflow - stock_change

Return type:

ndarray

get_stock_by_cohort()

Stock by cohort, i.e. the stock of each production year at each time step.

Return type:

ndarray

validator init_cohort_arrays  »  all fields
validator init_lifetime_model  »  all fields
validator init_solver  »  all fields
to_stock_type(desired_stock_type, **kwargs)

Return an object of a new stock type with values and dimensions the same as the original. **kwargs can be used to pass additional model attributes as required by the desired stock type, if these are not contained in the original stock type.

Parameters:

desired_stock_type (type)

validator validate_stock_arrays  »  all fields
validator validate_time_first_dim  »  all fields
property process_id: int

ID of the process the stock is associated with.

property shape: tuple

Shape of the stock, inflow, outflow arrays, defined by the dimensions.

pydantic model flodym.LifetimeModel

Contains shared functionality across the various lifetime models.

Fields:
  • dims (flodym.dimensions.DimensionSet)

  • inflow_at (str)

  • n_pts_per_interval (int)

  • time_letter (str)

Validators:
  • cast_prms » all fields

  • check_inflow_at » all fields

  • init_t » all fields

field dims: DimensionSet [Required]
field inflow_at: str = 'middle'

If no quadrature is used, all inflow happens at one point in time, either at the beginning of the time period (start), in the middle (middle) or at the end (end).

field n_pts_per_interval: int = 1

Inflow into the stock is in reality quite uniform over time, while survival factors only take into account a single point in time. This can be alleviated here by using numerical integration over each inflow time period with n points. A value of 1 means that the inflow is only evaluated once, depending on the inflow_at parameter. This is the default and should be used for long life times, e.g. >> 1 year. For n_pts > 1, the inflow is evaluated at n points in time within the time period. inflow_at is ignored in this case. A value of 2 means evaluation at the beginning and end of the time period, while higher values additionally evaluate at more points within the time period. Higher point numbers are only needed for short life times, e.g. < 1 year. Default is 1, meaning that the inflow is evaluated only once per time period.

field time_letter: str = 't'
cast_any_to_np_array(prm_in)
validator cast_prms  »  all fields
validator check_inflow_at  »  all fields
compute_outflow_pdf()

Returns an array year-by-cohort of the probability that an item added to stock in year m (aka cohort m) leaves in in year n. This value equals pdf(n,m).

compute_survival_factor()

Survival table self.sf(m,n) denotes the share of an inflow in year n (age-cohort) still present at the end of year m (after m-n years). The computation is self.sf(m,n) = ProbDist.sf(m-n), where ProbDist is the appropriate scipy function for the lifetime model chosen. For lifetimes 0 the sf is also 0, meaning that the age-cohort leaves during the same year of the inflow. The method compute outflow_sf returns an array year-by-cohort of the surviving fraction of a flow added to stock in year m (aka cohort m) in in year n. This value equals sf(n,m). This is the only method for the inflow-driven model where the lifetime distribution directly enters the computation. All other stock variables are determined by mass balance. The shape of the output sf array is NoofYears * NoofYears, and the meaning is years by age-cohorts. The method does nothing if the sf alreay exists. For example, sf could be assigned to the dynamic stock model from an exogenous computation to save time.

get_quad_points_and_weights()

Returns the quadrature points and weights for the inflow time periods.

validator init_t  »  all fields
abstract set_prms()
property pdf
abstract property prms: dict[str, ndarray | None]
property sf
property shape
pydantic model flodym.FixedLifetime

Fixed lifetime, age-cohort leaves the stock in the model year when a certain age, specified as ‘Mean’, is reached.

Fields:
  • dims ()

  • inflow_at ()

  • mean (Any)

  • n_pts_per_interval ()

  • time_letter ()

Validators:
  • cast_prms » all fields

  • check_inflow_at » all fields

  • init_t » all fields

field dims: DimensionSet [Required]
field inflow_at: str = 'middle'

If no quadrature is used, all inflow happens at one point in time, either at the beginning of the time period (start), in the middle (middle) or at the end (end).

field mean: Any = None
field n_pts_per_interval: int = 1

Inflow into the stock is in reality quite uniform over time, while survival factors only take into account a single point in time. This can be alleviated here by using numerical integration over each inflow time period with n points. A value of 1 means that the inflow is only evaluated once, depending on the inflow_at parameter. This is the default and should be used for long life times, e.g. >> 1 year. For n_pts > 1, the inflow is evaluated at n points in time within the time period. inflow_at is ignored in this case. A value of 2 means evaluation at the beginning and end of the time period, while higher values additionally evaluate at more points within the time period. Higher point numbers are only needed for short life times, e.g. < 1 year. Default is 1, meaning that the inflow is evaluated only once per time period.

field time_letter: str = 't'
cast_any_to_np_array(prm_in)
validator cast_prms  »  all fields
validator check_inflow_at  »  all fields
compute_outflow_pdf()

Returns an array year-by-cohort of the probability that an item added to stock in year m (aka cohort m) leaves in in year n. This value equals pdf(n,m).

compute_survival_factor()

Survival table self.sf(m,n) denotes the share of an inflow in year n (age-cohort) still present at the end of year m (after m-n years). The computation is self.sf(m,n) = ProbDist.sf(m-n), where ProbDist is the appropriate scipy function for the lifetime model chosen. For lifetimes 0 the sf is also 0, meaning that the age-cohort leaves during the same year of the inflow. The method compute outflow_sf returns an array year-by-cohort of the surviving fraction of a flow added to stock in year m (aka cohort m) in in year n. This value equals sf(n,m). This is the only method for the inflow-driven model where the lifetime distribution directly enters the computation. All other stock variables are determined by mass balance. The shape of the output sf array is NoofYears * NoofYears, and the meaning is years by age-cohorts. The method does nothing if the sf alreay exists. For example, sf could be assigned to the dynamic stock model from an exogenous computation to save time.

get_quad_points_and_weights()

Returns the quadrature points and weights for the inflow time periods.

validator init_t  »  all fields
set_prms(mean)
Parameters:

mean (FlodymArray)

property pdf
property prms
property sf
property shape
pydantic model flodym.NormalLifetime

Normally distributed lifetime with mean and standard deviation. Watch out for nonzero values, for negative ages, no correction or truncation done here. NOTE: As normal distributions have nonzero pdf for negative ages, which are physically impossible, these outflow contributions can either be ignored ( violates the mass balance) or allocated to the zeroth year of residence, the latter being implemented in the method compute compute_o_c_from_s_c. As alternative, use lognormal or folded normal distribution options.

Fields:
  • dims ()

  • inflow_at ()

  • mean ()

  • n_pts_per_interval ()

  • std ()

  • time_letter ()

Validators:
  • cast_prms » all fields

  • check_inflow_at » all fields

  • init_t » all fields

field dims: DimensionSet [Required]
field inflow_at: str = 'middle'

If no quadrature is used, all inflow happens at one point in time, either at the beginning of the time period (start), in the middle (middle) or at the end (end).

field mean: Any = None
field n_pts_per_interval: int = 1

Inflow into the stock is in reality quite uniform over time, while survival factors only take into account a single point in time. This can be alleviated here by using numerical integration over each inflow time period with n points. A value of 1 means that the inflow is only evaluated once, depending on the inflow_at parameter. This is the default and should be used for long life times, e.g. >> 1 year. For n_pts > 1, the inflow is evaluated at n points in time within the time period. inflow_at is ignored in this case. A value of 2 means evaluation at the beginning and end of the time period, while higher values additionally evaluate at more points within the time period. Higher point numbers are only needed for short life times, e.g. < 1 year. Default is 1, meaning that the inflow is evaluated only once per time period.

field std: Any = None
field time_letter: str = 't'
cast_any_to_np_array(prm_in)
validator cast_prms  »  all fields
validator check_inflow_at  »  all fields
compute_outflow_pdf()

Returns an array year-by-cohort of the probability that an item added to stock in year m (aka cohort m) leaves in in year n. This value equals pdf(n,m).

compute_survival_factor()

Survival table self.sf(m,n) denotes the share of an inflow in year n (age-cohort) still present at the end of year m (after m-n years). The computation is self.sf(m,n) = ProbDist.sf(m-n), where ProbDist is the appropriate scipy function for the lifetime model chosen. For lifetimes 0 the sf is also 0, meaning that the age-cohort leaves during the same year of the inflow. The method compute outflow_sf returns an array year-by-cohort of the surviving fraction of a flow added to stock in year m (aka cohort m) in in year n. This value equals sf(n,m). This is the only method for the inflow-driven model where the lifetime distribution directly enters the computation. All other stock variables are determined by mass balance. The shape of the output sf array is NoofYears * NoofYears, and the meaning is years by age-cohorts. The method does nothing if the sf alreay exists. For example, sf could be assigned to the dynamic stock model from an exogenous computation to save time.

get_quad_points_and_weights()

Returns the quadrature points and weights for the inflow time periods.

validator init_t  »  all fields
set_prms(mean, std)
Parameters:
property pdf
property prms
property sf
property shape
pydantic model flodym.FoldedNormalLifetime

Folded normal distribution, cf. https://en.wikipedia.org/wiki/Folded_normal_distribution NOTE: call this with the parameters of the normal distribution mu and sigma of curve BEFORE folding, curve after folding will have different mu and sigma.

Fields:
  • dims ()

  • inflow_at ()

  • mean ()

  • n_pts_per_interval ()

  • std ()

  • time_letter ()

Validators:
  • cast_prms » all fields

  • check_inflow_at » all fields

  • init_t » all fields

field dims: DimensionSet [Required]
field inflow_at: str = 'middle'

If no quadrature is used, all inflow happens at one point in time, either at the beginning of the time period (start), in the middle (middle) or at the end (end).

field mean: Any = None
field n_pts_per_interval: int = 1

Inflow into the stock is in reality quite uniform over time, while survival factors only take into account a single point in time. This can be alleviated here by using numerical integration over each inflow time period with n points. A value of 1 means that the inflow is only evaluated once, depending on the inflow_at parameter. This is the default and should be used for long life times, e.g. >> 1 year. For n_pts > 1, the inflow is evaluated at n points in time within the time period. inflow_at is ignored in this case. A value of 2 means evaluation at the beginning and end of the time period, while higher values additionally evaluate at more points within the time period. Higher point numbers are only needed for short life times, e.g. < 1 year. Default is 1, meaning that the inflow is evaluated only once per time period.

field std: Any = None
field time_letter: str = 't'
cast_any_to_np_array(prm_in)
validator cast_prms  »  all fields
validator check_inflow_at  »  all fields
compute_outflow_pdf()

Returns an array year-by-cohort of the probability that an item added to stock in year m (aka cohort m) leaves in in year n. This value equals pdf(n,m).

compute_survival_factor()

Survival table self.sf(m,n) denotes the share of an inflow in year n (age-cohort) still present at the end of year m (after m-n years). The computation is self.sf(m,n) = ProbDist.sf(m-n), where ProbDist is the appropriate scipy function for the lifetime model chosen. For lifetimes 0 the sf is also 0, meaning that the age-cohort leaves during the same year of the inflow. The method compute outflow_sf returns an array year-by-cohort of the surviving fraction of a flow added to stock in year m (aka cohort m) in in year n. This value equals sf(n,m). This is the only method for the inflow-driven model where the lifetime distribution directly enters the computation. All other stock variables are determined by mass balance. The shape of the output sf array is NoofYears * NoofYears, and the meaning is years by age-cohorts. The method does nothing if the sf alreay exists. For example, sf could be assigned to the dynamic stock model from an exogenous computation to save time.

get_quad_points_and_weights()

Returns the quadrature points and weights for the inflow time periods.

validator init_t  »  all fields
set_prms(mean, std)
Parameters:
property pdf
property prms
property sf
property shape
pydantic model flodym.LogNormalLifetime

Lognormal distribution Here, the mean and stddev of the lognormal curve, not those of the underlying normal distribution, need to be specified! Values chosen according to description on https://docs.scipy.org/doc/scipy-0.13.0/reference/generated/scipy.stats.lognorm.html Same result as EXCEL function “=LOGNORM.VERT(x;LT_LN;SG_LN;TRUE)”

Fields:
  • dims ()

  • inflow_at ()

  • mean ()

  • n_pts_per_interval ()

  • std ()

  • time_letter ()

Validators:
  • cast_prms » all fields

  • check_inflow_at » all fields

  • init_t » all fields

field dims: DimensionSet [Required]
field inflow_at: str = 'middle'

If no quadrature is used, all inflow happens at one point in time, either at the beginning of the time period (start), in the middle (middle) or at the end (end).

field mean: Any = None
field n_pts_per_interval: int = 1

Inflow into the stock is in reality quite uniform over time, while survival factors only take into account a single point in time. This can be alleviated here by using numerical integration over each inflow time period with n points. A value of 1 means that the inflow is only evaluated once, depending on the inflow_at parameter. This is the default and should be used for long life times, e.g. >> 1 year. For n_pts > 1, the inflow is evaluated at n points in time within the time period. inflow_at is ignored in this case. A value of 2 means evaluation at the beginning and end of the time period, while higher values additionally evaluate at more points within the time period. Higher point numbers are only needed for short life times, e.g. < 1 year. Default is 1, meaning that the inflow is evaluated only once per time period.

field std: Any = None
field time_letter: str = 't'
cast_any_to_np_array(prm_in)
validator cast_prms  »  all fields
validator check_inflow_at  »  all fields
compute_outflow_pdf()

Returns an array year-by-cohort of the probability that an item added to stock in year m (aka cohort m) leaves in in year n. This value equals pdf(n,m).

compute_survival_factor()

Survival table self.sf(m,n) denotes the share of an inflow in year n (age-cohort) still present at the end of year m (after m-n years). The computation is self.sf(m,n) = ProbDist.sf(m-n), where ProbDist is the appropriate scipy function for the lifetime model chosen. For lifetimes 0 the sf is also 0, meaning that the age-cohort leaves during the same year of the inflow. The method compute outflow_sf returns an array year-by-cohort of the surviving fraction of a flow added to stock in year m (aka cohort m) in in year n. This value equals sf(n,m). This is the only method for the inflow-driven model where the lifetime distribution directly enters the computation. All other stock variables are determined by mass balance. The shape of the output sf array is NoofYears * NoofYears, and the meaning is years by age-cohorts. The method does nothing if the sf alreay exists. For example, sf could be assigned to the dynamic stock model from an exogenous computation to save time.

get_quad_points_and_weights()

Returns the quadrature points and weights for the inflow time periods.

validator init_t  »  all fields
set_prms(mean, std)
Parameters:
property pdf
property prms
property sf
property shape
pydantic model flodym.WeibullLifetime

Weibull distribution with standard definition of scale and shape parameters.

Fields:
  • dims ()

  • inflow_at ()

  • n_pts_per_interval ()

  • time_letter ()

  • weibull_scale (Any)

  • weibull_shape (Any)

Validators:
  • cast_prms » all fields

  • check_inflow_at » all fields

  • init_t » all fields

field dims: DimensionSet [Required]
field inflow_at: str = 'middle'

If no quadrature is used, all inflow happens at one point in time, either at the beginning of the time period (start), in the middle (middle) or at the end (end).

field n_pts_per_interval: int = 1

Inflow into the stock is in reality quite uniform over time, while survival factors only take into account a single point in time. This can be alleviated here by using numerical integration over each inflow time period with n points. A value of 1 means that the inflow is only evaluated once, depending on the inflow_at parameter. This is the default and should be used for long life times, e.g. >> 1 year. For n_pts > 1, the inflow is evaluated at n points in time within the time period. inflow_at is ignored in this case. A value of 2 means evaluation at the beginning and end of the time period, while higher values additionally evaluate at more points within the time period. Higher point numbers are only needed for short life times, e.g. < 1 year. Default is 1, meaning that the inflow is evaluated only once per time period.

field time_letter: str = 't'
field weibull_scale: Any = None
field weibull_shape: Any = None
cast_any_to_np_array(prm_in)
validator cast_prms  »  all fields
validator check_inflow_at  »  all fields
compute_outflow_pdf()

Returns an array year-by-cohort of the probability that an item added to stock in year m (aka cohort m) leaves in in year n. This value equals pdf(n,m).

compute_survival_factor()

Survival table self.sf(m,n) denotes the share of an inflow in year n (age-cohort) still present at the end of year m (after m-n years). The computation is self.sf(m,n) = ProbDist.sf(m-n), where ProbDist is the appropriate scipy function for the lifetime model chosen. For lifetimes 0 the sf is also 0, meaning that the age-cohort leaves during the same year of the inflow. The method compute outflow_sf returns an array year-by-cohort of the surviving fraction of a flow added to stock in year m (aka cohort m) in in year n. This value equals sf(n,m). This is the only method for the inflow-driven model where the lifetime distribution directly enters the computation. All other stock variables are determined by mass balance. The shape of the output sf array is NoofYears * NoofYears, and the meaning is years by age-cohorts. The method does nothing if the sf alreay exists. For example, sf could be assigned to the dynamic stock model from an exogenous computation to save time.

get_quad_points_and_weights()

Returns the quadrature points and weights for the inflow time periods.

validator init_t  »  all fields
set_prms(weibull_shape, weibull_scale)
Parameters:
property pdf
property prms
property sf
property shape
flodym.make_empty_stocks(stock_definitions, processes, dims)

Initialise empty Stock objects for each of the stocks listed in stock definitions.

Parameters:
Return type:

dict[str, Stock]