The Brown Book
Documentation for the PooPyLab project
PooPyLab.unit_procs.bio.asm_reactor Class Reference

Bioreactor using ASM kinetics, derived as a "pipe" w/ active volume. More...

Inheritance diagram for PooPyLab.unit_procs.bio.asm_reactor:
PooPyLab.unit_procs.streams.pipe PooPyLab.unit_procs.streams.splitter PooPyLab.unit_procs.base.poopy_lab_obj

Public Member Functions

def __init__ (self, act_vol=38000, swd=3.5, ww_temp=20, DO=2, *args, **kw)
 Init w/ active volume, water depth, water temperature, & dissolved O2. More...
 
def is_converged (self, limit=0.01)
 Return the convergence status of the unit. More...
 
def discharge (self, method_name="BDF")
 Pass the total flow and blended components to the downstreams. More...
 
def assign_initial_guess (self, initial_guess)
 Assign the intial guess to the unit before simulation. More...
 
def set_active_vol (self, vol=380)
 Set the active process volume. More...
 
def get_active_vol (self)
 Return the active process volume. More...
 
def set_model_condition (self, ww_temp, DO)
 Set the wastewater temperature and dissolved O2 for the model. More...
 
def get_model_params (self)
 Return the kinetic parameters of the applied model. More...
 
def get_model_stoichs (self)
 Return the stoichiometrics of the applied model. More...
 
- Public Member Functions inherited from PooPyLab.unit_procs.streams.pipe
def __init__ (self)
 A few special steps are taken when initializing a "pipe": More...
 
def set_downstream_side (self, receiver)
 Define the downstream side outlet's connection. More...
 
def set_sidestream_flow (self, flow)
 Define the flow rate for the sidestream. More...
 
- Public Member Functions inherited from PooPyLab.unit_procs.streams.splitter
def set_flow_data_src (self, branch='Main', flow_ds=flow_data_src.TBD)
 Set the flow data source of the branch specified by the user. More...
 
def get_flow_data_src (self)
 Return the flow data source tags of all three branches.
 
def get_type (self)
 Return the type string of the process unit. More...
 
def has_sidestream (self)
 Return whether the unit has a sidestream. More...
 
def add_upstream (self, discharger, upst_branch='Main')
 Add the discharger's branch to inlet. More...
 
def has_discharger (self)
 Return whether the unit's inlet has been connected. More...
 
def get_upstream (self)
 Return the _inlet {} of the unit. More...
 
def totalize_inflow (self)
 Combine the individual flows specified in the self._inlet into one. More...
 
def blend_inlet_comps (self)
 Calculate the flow weighted average model component concentrations. More...
 
def update_combined_input (self)
 Update both total inflow and blended concentrations (model components).
 
def remove_upstream (self, discharger)
 Remove an existing discharger from inlet. More...
 
def set_downstream_main (self, rcvr)
 Define the main outlet by specifying the receiving process unit.
 
def main_outlet_connected (self)
 Return whether the mainstream outlet is connected. More...
 
def get_downstream_main (self)
 Return the process unit connected at the mainstream outlet. More...
 
def set_mainstream_flow_by_upstream (self, f=True)
 Set whether the mainstream flow = (total inflow - side outflow).
 
def set_mainstream_flow (self, flow=0)
 Define the mainstream outlet flow. More...
 
def get_main_outflow (self)
 Return the mainstream outflow. More...
 
def get_main_outlet_concs (self)
 Return a copy of the mainstream outlet concentrations. More...
 
def side_outlet_connected (self)
 Return True if the main outlet is connected, False if not. More...
 
def get_downstream_side (self)
 Return the process unit connected to the side outlet. More...
 
def sidestream_flow_defined (self)
 Return whether the sidestream flow rate has been defined. More...
 
def get_side_outflow (self)
 Return the sidestream outlet flow rate. More...
 
def get_side_outlet_concs (self)
 Return a copy of the sidestream outlet concentrations. More...
 
def set_flow (self, dschgr, flow)
 Specify the flow from the discharger. More...
 
def discharge (self)
 Pass the total flow and blended components to the downstreams. More...
 
def get_TSS (self, br='Main')
 Return the Total Suspended Solids of the specified branch.
 
def get_VSS (self, br='Main')
 Return the Volatile Suspended Solids of the specified branch.
 
def get_COD (self, br='Main')
 Return the Chemical Oxygen Demand (total) of the specified branch.
 
def get_sCOD (self, br='Main')
 Return the soluble COD of the specified branch.
 
def get_pCOD (self, br='Main')
 Return the particultate COD of the specified branch.
 
def get_TN (self, br='Main')
 Return the total nitrogen of the specified branch. More...
 
def get_orgN (self, br='Main')
 Return the organic nitrogen of the specified branch.
 
def get_inorgN (self, br='Main')
 Return the inorganic nitrogen of the specified branch.
 
def get_pN (self, br='Main')
 Return the particulate nitrogen of the specified branch.
 
def get_sN (self, br='Main')
 Return the soluble nitrogen of the specified branch.
 
def set_as_SRT_controller (self, setting=False)
 Set the current splitter as an Solids Retention Time controller. More...
 
def is_SRT_controller (self)
 Return whether a splitter is an SRT controller.
 

Private Member Functions

def _RKF45_ks (self)
 Calculate k1...k6 used in RKF45 method. More...
 
def _RKF45_err (self, k1, k3, k4, k5, k6)
 Calculate the norm of the error vector in RKF45 method. More...
 
def _runge_kutta_fehlberg_45 (self, tol=1e-4)
 Integration by using the Runge-Kutta-Fehlberg (RKF45) method. More...
 

Private Attributes

 __name__
 
 _type
 
 _active_vol
 active volume, m3
 
 _swd
 side water depth, m
 
 _area
 plan section area, m2
 
 _sludge
 sludge mixed liquor contained in the reactor
 
 _del_C_del_t
 
 _in_comps
 
 _mo_comps
 
 _prev_mo_comps
 
 _prev_so_comps
 
 _upstream_set_mo_flow
 
 _step
 
 _prev_local_err
 
 _atol
 
 _rtol
 
 _solultion
 
 _so_comps
 

Static Private Attributes

int __id
 

Detailed Description

Bioreactor using ASM kinetics, derived as a "pipe" w/ active volume.

Current design only uses ASM 1 model. Will need to add flexibility for using ASM 2d, ASM 3, and user revised versions of them.

The "asm_reactor" contain sludge mixed liquor of a certain kinetics described by the chosen model.

The integration of the model is also done by the "asm_reactor".

Constructor & Destructor Documentation

◆ __init__()

def PooPyLab.unit_procs.bio.asm_reactor.__init__ (   self,
  act_vol = 38000,
  swd = 3.5,
  ww_temp = 20,
  DO = 2,
args,
**  kw 
)

Init w/ active volume, water depth, water temperature, & dissolved O2.

Parameters
act_volactive process volume, m3
swdside water depth, m
ww_tempwastewater temperature, degC
DOdissolved oxygen, mg/L
args(provision for other parameters for different models)
kw(provision for other parameters)
Return
None

Member Function Documentation

◆ is_converged()

def PooPyLab.unit_procs.bio.asm_reactor.is_converged (   self,
  limit = 0.01 
)

Return the convergence status of the unit.

   This function checks the absolute changes of flows and concentrations
   between last round of simulation and current. If all changes are within
   the tolerance limit, it flags the self._converged True.

   The same absolute changes will be used in evaluting the convergence of
   splitter, clarifier, pipe, effluent, WAS, and other unit processes that
   do not have splicit definitions of change rate terms (dy/dt = f(...)).

   This function is redefined for unit processes like asm_reactor whose
   models have specific dy/dt terms.
Parameters
limitLimit within which the simulation is considered converged
Return
True/False
See
_check_conc_cnvg()

Reimplemented from PooPyLab.unit_procs.streams.splitter.

◆ discharge()

def PooPyLab.unit_procs.bio.asm_reactor.discharge (   self,
  method_name = "BDF" 
)

Pass the total flow and blended components to the downstreams.

   This function is re-implemented for "asm_reactor". Because of the
   biological reactions "happening" in the "asm_reactor", integration of
   the model (Note 1) is carried out here before sending the results to
   the down stream.
Parameters
method_name"BDF", "RK45", "Radau", etc.(see Note 2 below)
Retrun
self._sludge._comps
Notes
1) It is highly recommended the model components are arranged such that all the soluble ones are ahead of the particulate ones in the array. Generally, soluble components requires smaller time steps than particulate ones. This kind of arrangement will enable quick identification of soluble/particulate components that may have very different suitable time step during integration. Using appropriate but different time steps for the soluble and particulate components is required for fast integrations with correct results. This is how the ODE partitioning method suggested in the IWA ASM1 report works. Although PooPyLab doesn't apply this relaxation scheme as of now, arranging the model components in such partitioned way will allow future exploration of optimization approaches.
2) There are a few integration methods attempted for PooPyLab: Euler, Runge-Kutta 4th order, Runge-Kutta-Felhberg 4/5, RK-Dormand-Prince-4/5, and the ODE system partitioning scheme suggested in the IWA ASM1 report. After much study, it is decided to settle with scipy.integrate.solve_ivp routine for now so that the rest of the PooPyLab development can progress, while KZ continues in his study of BDF methods and attempts for a home brew version. Euler, RK4, RKF45, RKDP45, and Partitioned ODE methods have been coded and tested in the past but no longer in use as of now, except for RKF45. The unused code is moved to bio_py_funcs_not_used.txt for archiving.
See
_runge_kutta_fehlberg_45()

◆ assign_initial_guess()

def PooPyLab.unit_procs.bio.asm_reactor.assign_initial_guess (   self,
  initial_guess 
)

Assign the intial guess to the unit before simulation.

   This function is re-implemented for "asm_reactor" which contains the
   "sludge" whose kinetics are described by the model.

   When passing the initial guess into an "asm_reactor", the reactor's
   inlet, mainstream outlet, and the "sludge" in it all get the same list
   of model component concentrations.
Parameters
initial_guesslist of model components
Return
None

Reimplemented from PooPyLab.unit_procs.streams.splitter.

◆ set_active_vol()

def PooPyLab.unit_procs.bio.asm_reactor.set_active_vol (   self,
  vol = 380 
)

Set the active process volume.

Parameters
volactive volume to be used. (m3)
Return
None

◆ get_active_vol()

def PooPyLab.unit_procs.bio.asm_reactor.get_active_vol (   self)

Return the active process volume.

(m3)

◆ set_model_condition()

def PooPyLab.unit_procs.bio.asm_reactor.set_model_condition (   self,
  ww_temp,
  DO 
)

Set the wastewater temperature and dissolved O2 for the model.

   This function updates the model conditions for the "sludge" the
   "asm_reactor" contains.
Parameters
ww_tempwastewtaer temperature in degC;
DOdissolved O2 concentration in mg/L.
Return
None
See
ASMModel.ASM_1.update().

◆ get_model_params()

def PooPyLab.unit_procs.bio.asm_reactor.get_model_params (   self)

Return the kinetic parameters of the applied model.

Return
{param_name, param_val_adjusted}
See
ASMModel.ASM_1.get_params().

◆ get_model_stoichs()

def PooPyLab.unit_procs.bio.asm_reactor.get_model_stoichs (   self)

Return the stoichiometrics of the applied model.

Return
{id_of_stoich, val}
See
ASMModel.ASM_1.get_stoichs().

◆ _RKF45_ks()

def PooPyLab.unit_procs.bio.asm_reactor._RKF45_ks (   self)
private

Calculate k1...k6 used in RKF45 method.

See
_runge_kutta_fehlberg_45();
_RKF45_err().

◆ _RKF45_err()

def PooPyLab.unit_procs.bio.asm_reactor._RKF45_err (   self,
  k1,
  k3,
  k4,
  k5,
  k6 
)
private

Calculate the norm of the error vector in RKF45 method.

        k1, k3, ... ,k6: intermediate step vectors of RKF45
Return
Norm of the error vector
See
_runge_kutta_fehlberg_45();
_RKF45_ks().

◆ _runge_kutta_fehlberg_45()

def PooPyLab.unit_procs.bio.asm_reactor._runge_kutta_fehlberg_45 (   self,
  tol = 1e-4 
)
private

Integration by using the Runge-Kutta-Fehlberg (RKF45) method.

Parameters
toluser defined tolerance of error
Return
step size used
See
_RKF45_ks();
_RKF45_err().

The documentation for this class was generated from the following file: