Suspended Growth Bioreactor¶
QSDsan: Quantitative Sustainable Design for sanitation and resource recovery systems
This module is developed by:
Joy Zhang <joycheung1994@gmail.com>
This module is under the University of Illinois/NCSA Open Source License. Please refer to https://github.com/QSD-Group/QSDsan/blob/main/LICENSE.txt for license details.
- class qsdsan.unit_operations.dynamic._bioreactor.AerobicDigester(ID='', ins: Sequence[AbstractStream] | None = None, outs: Sequence[AbstractStream] | None = (), thermo=None, init_with='WasteStream', V_max=1000, activated_sludge_model=None, organic_particulate_inert_degradation_process=None, aeration=1.0, DO_ID='S_O2', isdynamic=True, **kwargs)¶
Models an aerobic digester with activated sludge models, is a subclass of CSTR. An additional degradation process of particulate inert organic materials is considered in addition to typical activated sludge processes.
- Parameters:
activated_sludge_model (
CompiledProcesses, optional) – The activated sludge model used for the biochemical reactions. The default is None.organic_particulate_inert_degradation_process (
Process, optional) – The degradation process of the particulate inert organic materials. The default is None.
See also
qsdsan.process_models.ASM1,qsdsan.process_models.ASM2d,qsdsan.process_models.mASM2d,qsdsan.process_models.ASM_AeDigAddOn,qsdsan.sanunites.CSTRExamples
>>> from qsdsan import WasteStream, process_models as pc, unit_operations as su >>> cmps = pc.create_asm1_cmps() >>> twas = WasteStream('thickened_WAS') >>> twas.set_flow_by_concentration( ... flow_tot=50, ... concentrations=dict( ... S_I=30, S_S=1, X_I=17000, X_S=800, X_BH=38000, X_BA=2300, ... X_P=6500, S_O=0.5, S_NH=2, S_ND=1, X_ND=65, S_NO=10, ... S_N2=25, S_ALK=84), ... units=('m3/d', 'mg/L')) >>> asm = pc.ASM1() >>> AED = su.AerobicDigester('AED', ins=twas, outs=('digested_WAS',), ... V_max=3000, activated_sludge_model=asm, ... DO_ID='S_O', aeration=1.0) >>> AED.simulate(t_span=(0, 400), method='BDF') >>> AED.show() AerobicDigester: AED ins... [0] thickened_WAS phase: 'l', T: 298.15 K, P: 101325 Pa flow (g/hr): S_I 62.5 S_S 2.08 X_I 3.54e+04 X_S 1.67e+03 X_BH 7.92e+04 X_BA 4.79e+03 X_P 1.35e+04 S_O 1.04 S_NO 20.8 S_NH 4.17 S_ND 2.08 X_ND 135 S_ALK 175 S_N2 52.1 H2O 1.99e+06 WasteStream-specific properties: pH : 7.0 Alkalinity : 2.5 mmol/L COD : 64631.0 mg/L BOD : 23300.3 mg/L TC : 22923.3 mg/L TOC : 22839.3 mg/L TN : 4712.0 mg/L TP : 1009.0 mg/L TK : 223.2 mg/L TSS : 48450.0 mg/L outs... [0] digested_WAS phase: 'l', T: 298.15 K, P: 101325 Pa flow (g/hr): S_I 62.5 S_S 3.58e+04 X_I 1.04e+04 X_S 123 X_BH 9.6e+03 X_BA 1.59e+03 X_P 2.77e+04 S_O 2.08 S_NO 4.17e+03 S_NH 0.101 S_ND 0.975 X_ND 8.74 S_N2 2.51e+03 H2O 2e+06 WasteStream-specific properties: pH : 7.0 Alkalinity : 2.5 mmol/L COD : 40987.6 mg/L BOD : 15416.1 mg/L TC : 13985.1 mg/L TOC : 13985.1 mg/L TN : 3534.1 mg/L TP : 458.2 mg/L TK : 89.3 mg/L TSS : 17813.2 mg/L
- line: str = 'Aerobic digester'¶
class-attribute Name denoting the type of Unit class. Defaults to the class name of the first child class
- property organic_particulate_inert_degradation_process¶
[
Processor NoneType] Process object for degradation of particulate inert organic materials in the aerobic digester. If none specified, will attempt to create a Process model according to components by default.
- run()¶
Run mass and energy balance. This method also runs specifications user defined specifications unless it is being run within a specification (to avoid infinite loops).
See also
_run,specifications,add_specification,add_bounded_numerical_specification
- class qsdsan.unit_operations.dynamic._bioreactor.BatchExperiment(ID='', ins: Sequence[AbstractStream] | None = None, outs: Sequence[AbstractStream] | None = (), thermo=None, init_with='WasteStream', model=None, isdynamic=True, exogenous_vars=(), **kwargs)¶
A batch reactor in experimental settings.
- Parameters:
model (
CompiledProcesses, optional) – Process model that describes the dynamics of state variables. The state of the batch reactor is entirely determined by the stoichiometry and rate function in this model.
Examples
>>> import qsdsan.unit_operations as su, qsdsan.process_models as pc >>> cmps = pc.create_asm1_cmps() >>> asm1 = pc.ASM1() >>> BE = su.BatchExperiment('BE', model=asm1) >>> BE.set_init_conc(S_S=20, X_BH=500, S_O=8, S_ND=3, S_ALK=84) >>> BE.simulate(t_span=(0,10), method='BDF') >>> for k,v in BE.state.items(): ... if v != 0: ... print(f'{k}{" "*(7-len(k))}{v:.2f}') S_S 0.93 X_S 446.16 X_BH 25.66 X_P 39.25 S_O -0.00 S_NO -0.00 S_NH 2.12 S_ND 0.00 X_ND 36.47 S_ALK 85.82 S_N2 0.00
- line: str = 'Batch experiment'¶
class-attribute Name denoting the type of Unit class. Defaults to the class name of the first child class
- property model¶
[
CompiledProcessesor NoneType] Suspended growth model.
- run()¶
Run mass and energy balance. This method also runs specifications user defined specifications unless it is being run within a specification (to avoid infinite loops).
See also
_run,specifications,add_specification,add_bounded_numerical_specification
- set_init_conc(**kwargs)¶
set the initial concentrations of the BatchExperiment.
- property state¶
The state of the BatchExperiment, i.e., component concentrations.
- class qsdsan.unit_operations.dynamic._bioreactor.CSTR(ID='', ins: Sequence[AbstractStream] | None = None, outs: Sequence[AbstractStream] | None = (), split=None, thermo=None, init_with='WasteStream', V_max=1000, W_tank=6.4, D_tank=3.65, freeboard=0.61, t_wall=None, t_slab=None, aeration=2.0, DO_ID='S_O2', suspended_growth_model=None, gas_stripping=False, gas_IDs=None, stripping_kLa_min=None, K_Henry=None, D_gas=None, p_gas_atm=None, isdynamic=True, exogenous_vars=(), **kwargs)¶
A single continuous stirred tank reactor.
- Parameters:
ID (str) – ID for the reactor.
ins (
WasteStream) – Influents to the reactor. Can be an array of up to 3 WasteStream objects by default, typically wastewater to be treated, recycled effluent, recycled activated sludge.outs (
WasteStream) – Treated effluent.split (iterable of float) – Volumetric splits of effluent flows if there are more than one effluent. The default is None.
V_max (float) – Designed volume, in [m^3]. The default is 1000.
W_tank (#)
tank (# Freeboard added to the depth of the reactor)
[1 ([m]. The default is 0.61 m (2 ft).)
code] (Yalin's adaptation of)
D_tank (#)
[1
code]
freeboard (#)
tank
[1
code]
aeration (float or
Process, optional) – Aeration setting. Either specify a targeted dissolved oxygen concentration in [mg O2/L] or provide aProcessobject to represent aeration, or None for no aeration. The default is 2.0.DO_ID (str, optional) – The
ComponentID for dissolved oxygen, only relevant when the reactor is aerated. The default is ‘S_O2’.suspended_growth_model (
Processes, optional) – The suspended growth biokinetic model. The default is None.exogenous_var (iterable[
ExogenousDynamicVariable], optional) – Any exogenous dynamic variables that affect the process mass balance, e.g., temperature, sunlight irradiance. Must be independent of state variables of the suspended_growth_model (if has one).
References
- [1] Shoener, B. D.; Zhong, C.; Greiner, A. D.; Khunjar, W. O.; Hong, P.-Y.; Guest, J. S.
Design of Anaerobic Membrane Bioreactors for the Valorization of Dilute Organic Carbon Waste Streams. Energy Environ. Sci. 2016, 9 (3), 1102–1112. https://doi.org/10.1039/C5EE03715H.
- property DO_ID¶
[str] The Component ID for dissolved oxygen used in the suspended growth model and the aeration model.
- property V_max¶
[float] The designed maximum liquid volume, not accounting for increased volume due to aeration, in m^3.
- property aeration¶
[
Processor float or NoneType] Aeration model.
- line: str = 'CSTR'¶
class-attribute Name denoting the type of Unit class. Defaults to the class name of the first child class
- run()¶
Run mass and energy balance. This method also runs specifications user defined specifications unless it is being run within a specification (to avoid infinite loops).
See also
_run,specifications,add_specification,add_bounded_numerical_specification
- set_init_conc(**kwargs)¶
set the initial concentrations [mg/L] of the CSTR.
- property split¶
[numpy.1darray or NoneType] The volumetric split of outs.
- property state¶
The state of the CSTR, including component concentrations [mg/L] and flow rate [m^3/d].
- property suspended_growth_model¶
[
CompiledProcessesor NoneType] Suspended growth model.
- class qsdsan.unit_operations.dynamic._bioreactor.PFR(ID='', ins: Sequence[AbstractStream] | None = None, outs: Sequence[AbstractStream] | None = (), thermo=None, init_with='WasteStream', N_tanks_in_series=5, V_tanks=[1500, 1500, 3000, 3000, 3000], influent_fractions=[[1.0, 0, 0, 0, 0]], internal_recycles=[(4, 0, 35000)], DO_setpoints=[], kLa=[0, 0, 120, 120, 60], DO_sat=8.0, DO_ID='S_O2', suspended_growth_model=None, gas_stripping=False, gas_IDs=None, stripping_kLa_min=None, K_Henry=None, D_gas=None, p_gas_atm=None, isdynamic=True, **kwargs)¶
A plug flow reactor discretized into CSTRs in series with internal recycles and multiple influents.
- Parameters:
N_tanks_in_series (int, optional) – The number of CSTRs or zones in which the PFR is discretized. The default is 5.
V_tanks (iterable[float], optional) – The volume [m3] for each zone, length must match the number of CSTRs. The default is [1500, 1500, 3000, 3000, 3000].
influent_fractions (iterable[float], optional) – The volume fractions fed to each zone for each influent. Number of rows must match the number of influents. Number of columns must match the number of zones. The default is [[1.0, 0,0,0,0]].
internal_recycles (list[3-tuple[int, int, float]], optional) – A list of three-element tuples (i, j, Q) indicating internal recycling streams from zone i to zone j with a flowrate of Q [m3/d], respectively. Indices i,j start from 0 not 1. The default is [(4,0,35000),].
DO_setpoints (iterable[float], optional) – Dissolve oxygen setpoints of each zone [mg-O2/L]. Length must match the number of zones. 0 is treated as no active aeration. DO setpoints take priority over kLa values. The default is [].
kLa (iterable[float], optional) – Oxygen transfer rate constant values of each zone [d^(-1)]. If DO setpoints are specified, kLa values would be ignored in process simulation. The default is [0, 0, 120, 120, 60].
DO_sat (float, optional) – Saturation dissolved oxygen concentration [mg/L]. The default is 8.0.
gas_stripping (bool, optional) – Whether to model gas stripping. The default is False.
gas_IDs (iterable[str], optional) – Component IDs of stripped gases. The default is None.
stripping_kLa_min (iterable[float], optional) – Minimum gas transfer rate constants [d^(-1)] of each stripped gas component. The default is None.
K_Henry (iterable[float], optional) – Henry’s law constants [(conc)/atm], where “conc” indicate the concentration unit for state variables in the suspended growth model. The default is None.
D_gas (iterable[float], optional) – Gas diffusion coefficients in water [m2/s]. The default is None.
p_gas_atm (iterable[float], optional) – Partial pressure of the stripped gases in the air [atm]. The default is None.
Examples
>>> import qsdsan.unit_operations as su, qsdsan.process_models as pc >>> from qsdsan import WasteStream >>> cmps = pc.create_asm1_cmps(adjust_MW_to_measured_as=False) >>> asm1 = pc.ASM1() >>> inf = WasteStream('inf', H2O=1.53e6, S_I=46, S_S=54, X_I=1770, X_S=230, ... X_BH=3870, X_BA=225, X_P=680, S_O=0.377, S_NO=7.98, ... S_NH=25.6, S_ND=5.87, X_ND=13.4, S_ALK=103) >>> AS = su.PFR('AS', ins=(inf,), outs=('eff',), ... N_tanks_in_series=5, V_tanks=[1000]*2+[1333]*3, ... influent_fractions=[[1.0, 0,0,0,0]], DO_setpoints=[0]*2+[1.7, 2.4, 0.5], ... internal_recycles=[(4,0,55338)], DO_ID='S_O', ... suspended_growth_model=asm1) >>> AS.set_init_conc(S_I=30, S_S=5, X_I=1000, X_S=100, X_BH=500, X_BA=100, ... X_P=100, S_O=2, S_NO=20, S_NH=2, S_ND=1, X_ND=1, S_ALK=84) >>> AS.simulate(t_span=(0,100), method='BDF') >>> eff, = AS.outs >>> eff.show() WasteStream: eff from <PFR: AS> phase: 'l', T: 298.15 K, P: 101325 Pa flow (g/hr): S_I 4.6e+04 S_S 2.24e+03 X_I 1.77e+06 X_S 7.59e+04 X_BH 3.94e+06 X_BA 2.3e+05 X_P 6.95e+05 S_O 770 S_NO 1.6e+04 S_NH 2.68e+03 S_ND 1.06e+03 X_ND 5.42e+03 S_ALK 7.65e+04 S_N2 2.11e+04 H2O 1.53e+09 WasteStream-specific properties: pH : 7.0 Alkalinity : 2.5 mmol/L COD : 4389.0 mg/L BOD : 1563.3 mg/L TC : 1599.8 mg/L TOC : 1550.1 mg/L TN : 329.0 mg/L TP : 68.2 mg/L TK : 15.1 mg/L TSS : 3268.3 mg/L Component concentrations (mg/L): S_I 29.9 S_S 1.5 X_I 1150.0 X_S 49.3 X_BH 2557.0 X_BA 149.5 X_P 451.9 S_O 0.5 S_NO 10.4 S_NH 1.7 S_ND 0.7 X_ND 3.5 S_ALK 49.7 S_N2 13.7 H2O 994138.1
See also
qsdsan.unit_operations.CSTR- property DO_ID¶
[str] The Component ID for dissolved oxygen used in the suspended growth model and the aeration model.
- property DO_setpoints¶
[iterable[float]] Dissolve oxygen setpoints of CSTRs in series [mg-O2/L]. 0 is treated as no active aeration. DO setpoints take priority over kLa values.
- property V_tanks¶
[iterable[float]] Volumes of CSTRs in series [m3]
- property influent_fractions¶
[iterable[float]] Fractions of influents fed into different zones [unitless]
- property internal_recycles¶
[list[3-tuple[int, int, float]]] A list of three-element tuples (i, j, Q) indicating internal recycling streams from zone i to zone j with a flowrate of Q [m3/d], respectively. Indices i,j start from 0 not 1.
- property kLa¶
[iterable[float]] Aeration kLa values of CSTRs in series [d^(-1)]. If DO setpoints are specified, kLa values would be ignored in process simulation.
- line: str = 'PFR'¶
class-attribute Name denoting the type of Unit class. Defaults to the class name of the first child class
- run()¶
Run mass and energy balance. This method also runs specifications user defined specifications unless it is being run within a specification (to avoid infinite loops).
See also
_run,specifications,add_specification,add_bounded_numerical_specification
- set_init_conc(concentrations=None, i_zone=None, **kwargs)¶
set the initial concentrations [mg/L] of specific zones.
- property state¶
The state of the PFR, including component concentrations [mg/L] and flow rate [m^3/d] for each zone.
- property suspended_growth_model¶
[
CompiledProcessesor NoneType] Suspended growth model.
Example¶
import qsdsan.process_models as pc, qsdsan.unit_operations as su
from qsdsan import System
cmps = pc.create_asm1_cmps()
asm1 = pc.ASM1()
BE = su.BatchExperiment('E1', model=asm1)
BE.set_init_conc(S_S=5, X_S=100, X_BH=500, X_BA=100, X_P=100, S_O=8,
S_NO=20, S_NH=2, S_ND=1, X_ND=1, S_ALK=7*12)
sys = System('sys', path=(BE,))
sys.set_dynamic_tracker(BE)
sys.simulate(t_span=(0, 3))
BE.scope.plot_time_series(('X_BH', 'X_BA', 'X_S'))