Source code for pkpd.biopharmaceutics

# **************************************************************************
# *
# * Authors:     Carlos Oscar Sorzano (info@kinestat.com)
# *
# * Kinestat Pharma
# *
# * This program is free software; you can redistribute it and/or modify
# * it under the terms of the GNU General Public License as published by
# * the Free Software Foundation; either version 2 of the License, or
# * (at your option) any later version.
# *
# * This program is distributed in the hope that it will be useful,
# * but WITHOUT ANY WARRANTY; without even the implied warranty of
# * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
# * GNU General Public License for more details.
# *
# * You should have received a copy of the GNU General Public License
# * along with this program; if not, write to the Free Software
# * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA
# * 02111-1307  USA
# *
# *  All comments concerning this program package may be sent to the
# *  e-mail address 'info@kinestat.com'
# *
# **************************************************************************
"""
Biopharmaceutics: Drug sources and how they dissolve
"""
import copy
import math
import numpy as np
from .pkpd_units import PKPDUnit, changeRateToWeight, divideUnits, inverseUnits
from pkpd.utils import uniqueFloatValues, excelWriteRow
from scipy.interpolate import InterpolatedUnivariateSpline, PchipInterpolator
from scipy.special import lambertw, wrightomega

[docs]class BiopharmaceuticsModel: def __init__(self): self.parameters = [] self.ptrExperiment = None
[docs] def getNumberOfParameters(self): return len(self.getParameterNames())
def getDescription(self): pass
[docs] def getParameterNames(self): pass
[docs] def calculateParameterUnits(self,sample): pass
[docs] def setExperiment(self, ptrExperiment): self.ptrExperiment = ptrExperiment
[docs] def getExperiment(self): return self.ptrExperiment
[docs] def getDoseUnits(self): if self.ptrExperiment is None: return PKPDUnit.UNIT_WEIGHT_mg else: return self.ptrExperiment.getDoseUnits()
[docs] def setParameters(self, parameters): self.parameters = parameters
[docs] def getAg(self,t): # Total amount of drug that is available at time t return 0.0
[docs] def getEquation(self): return ""
[docs] def getModelEquation(self): return ""
[docs] def getDescription(self): return ""
[docs] def areParametersSignificant(self, lowerBound, upperBound): retval=[] for i in range(len(self.parameters)): lower = lowerBound[i] upper = upperBound[i] if lower<0 and upper>0: retval.append("False") elif lower>0 or upper<0: retval.append("True") else: retval.append("NA") return retval
[docs] def areParametersValid(self, p): return np.sum(p<0)==0
[docs]class BiopharmaceuticsModelOrder0(BiopharmaceuticsModel): def getDescription(self): return ['Constant absorption rate']
[docs] def getParameterNames(self): return ['Rin']
[docs] def calculateParameterUnits(self,sample): if self.ptrExperiment is None: dunits=PKPDUnit.UNIT_WEIGHT_mg tunits=PKPDUnit.UNIT_TIME_MIN else: dunits=self.ptrExperiment.getDoseUnits() tunits = self.ptrExperiment.getTimeUnits().unit print(dunits,tunits) rateUnits = divideUnits(dunits,tunits) self.parameterUnits = [rateUnits] return self.parameterUnits
[docs] def getAg(self,t): if t<0: return 0.0 Rin = self.parameters[0] return max(self.Amax-Rin*t,0.0)
[docs] def getEquation(self): Rin = self.parameters[0] return "D(t)=(%f)*t"%(Rin)
[docs] def getModelEquation(self): return "D(t)=Rin*t"
[docs] def getDescription(self): return "Zero order absorption (%s)"%self.__class__.__name__
[docs]class BiopharmaceuticsModelOrder01(BiopharmaceuticsModel): def getDescription(self): return ['Constant absorption rate','Constant absorption time','Absorption rate']
[docs] def getParameterNames(self): return ['Rin','t0','Ka']
[docs] def calculateParameterUnits(self,sample): tunits = self.ptrExperiment.getTimeUnits().unit if tunits == PKPDUnit.UNIT_TIME_MIN: rateUnits = PKPDUnit.UNIT_WEIGHTINVTIME_mg_MIN else: rateUnits = PKPDUnit.UNIT_WEIGHTINVTIME_mg_H self.parameterUnits = [rateUnits, tunits, inverseUnits(tunits)] return self.parameterUnits
[docs] def getAg(self,t): if t<0: return 0.0 Rin = self.parameters[0] t0 = self.parameters[1] Ka = self.parameters[2] if t<t0: return max(self.Amax-Rin*t,0.0) else: A0=max(self.Amax-Rin*t0,0.0) return A0*math.exp(-Ka*(t-t0))
[docs] def getEquation(self): Rin = self.parameters[0] t0 = self.parameters[1] Ka = self.parameters[2] return "D(t)=(%f)*t if t<%f and (%f)*t+((%f)-(%f)*(%f))*(1-exp(-(%f)*t) if t>%f"%(Rin, t0, Rin, self.Amax, Rin, t0, Ka, t0)
[docs] def getModelEquation(self): return "D(t)=Rin*t if t<t0 and Rin*t0+(Amax-Rin*t0)*(1-exp(-Ka*(t-t0)) if t>t0"
[docs] def getDescription(self): return "Zero-First Mixed order absorption (%s)"%self.__class__.__name__
[docs]class BiopharmaceuticsModelOrder01Tlag1(BiopharmaceuticsModel): def getDescription(self): return ['Fraction absorbed at order 0','Constant absorption rate','Constant absorption time','Absorption rate']
[docs] def getParameterNames(self): return ['F0','Rin','tlag1','Ka']
[docs] def calculateParameterUnits(self,sample): dUnits = self.getDoseUnits() tUnits = self.ptrExperiment.getTimeUnits().unit self.parameterUnits = [PKPDUnit.UNIT_NONE, divideUnits(dUnits,tUnits), tUnits, inverseUnits(tUnits)] return self.parameterUnits
[docs] def getAg(self,t): if t<0: return 0.0 F0 = self.parameters[0] Rin = self.parameters[1] if Rin<0.0: return 0.0 tlag1 = self.parameters[2] Ka = self.parameters[3] A=self.Amax A-=min(self.Amax*F0,Rin*t) if t>tlag1: A-=self.Amax*(1-F0)*(1-math.exp(-Ka*(t-tlag1))) return A
[docs] def getEquation(self): F0 = self.parameters[0] Rin = self.parameters[1] tlag1 = self.parameters[2] Ka = self.parameters[3] return "D(t)=(%f)*t if 0<t<%f and (%f)*(%f)*(1-exp(-(%f)*(t-(%f))) if t>%f"%(Rin, F0*self.Amax/Rin, self.Amax, 1-F0, Ka, tlag1, tlag1)
[docs] def getModelEquation(self): return "D(t)=Rin*t if t<Amax*F0/Rin and Amax*(1-F0)*(1-exp(-Ka*(t-tlag1)) if t>tlag1"
[docs] def getDescription(self): return "Zero-First Mixed order absorption (%s)"%self.__class__.__name__
[docs]class BiopharmaceuticsModelOrder1(BiopharmaceuticsModel): def getDescription(self): return ['Absorption rate']
[docs] def getParameterNames(self): return ['Ka']
[docs] def calculateParameterUnits(self,sample): self.parameterUnits = [inverseUnits(self.ptrExperiment.getTimeUnits().unit)] return self.parameterUnits
[docs] def getAg(self,t): if t<0: return 0.0 Ka = self.parameters[0] return self.Amax*math.exp(-Ka*t)
[docs] def getEquation(self): Ka = self.parameters[0] return "D(t)=(%f)*(1-exp(-(%f)*t)"%(self.Amax,Ka)
[docs] def getModelEquation(self): return "D(t)=Amax*(1-exp(-Ka*t))"
[docs] def getDescription(self): return "First order absorption (%s)"%self.__class__.__name__
[docs]class BiopharmaceuticsModelOrderFractional(BiopharmaceuticsModel): def getDescription(self): return ['Initial amount','Constant absorption rate', 'alpha']
[docs] def getParameterNames(self): return ['Amax','K','alpha']
[docs] def calculateParameterUnits(self,sample): self.parameterUnits = [PKPDUnit.UNIT_WEIGHT_mg,PKPDUnit.UNIT_WEIGHTINVTIME_mg_MIN,PKPDUnit.UNIT_NONE] return self.paramterUnits
[docs] def getAg(self,t): if t<=0: return 0.0 Amax = self.parameters[0] K = self.parameters[1] alpha = self.parameters[2] aux = alpha*K*t if aux>Amax: return Amax return Amax-math.pow(math.pow(Amax,alpha)-aux,1.0/alpha)
[docs] def getEquation(self): Amax = self.parameters[0] K = self.parameters[1] alpha = self.parameters[2] return "D(t)=(%f)-((%f)^(%f)-(%f)*(%f)*t)^(1/(%f))"%(Amax,Amax,alpha,K,alpha,alpha)
[docs] def getModelEquation(self): return "D(t)=Amax-(Amax^alpha-alpha*K*t)^(1/alpha)"
[docs] def getDescription(self): return "Fractional order absorption (%s)"%self.__class__.__name__
[docs] def areParametersValid(self, p): return np.sum(p<0)==0 and p[2]>0 and p[2]<1
[docs]class BiopharmaceuticsModelImmediateAndOrder1(BiopharmaceuticsModel): def getDescription(self): return ['Absorption rate','Immediate fraction']
[docs] def getParameterNames(self): return ['Ka','F']
[docs] def calculateParameterUnits(self,sample): self.parameterUnits = [inverseUnits(self.ptrExperiment.getTimeUnits().unit),PKPDUnit.UNIT_NONE] return self.parameterUnits
[docs] def getAg(self,t): if t<0: return 0.0 Ka = self.parameters[0] F = self.parameters[1] return self.Amax*(1-F)*math.exp(-Ka*t)
[docs] def getEquation(self): Ka = self.parameters[0] F = self.parameters[1] return "D(t)=(%f)*delta(t)+(%f)*(1-exp(-(%f)*t)"%(self.Amax*F,self.Amax*(1-F),Ka)
[docs] def getModelEquation(self): return "D(t)=Amax*F*delta(t)+(1-F)*Amax*(1-exp(-Ka*t))"
[docs] def getDescription(self): return "Immediate and First order absorption (%s)"%self.__class__.__name__
[docs]class BiopharmaceuticsModelOrder1AndOrder1(BiopharmaceuticsModel): def getDescription(self): return ['Absorption rate1','Absorption rate2','tlag2','Fraction 1']
[docs] def getParameterNames(self): return ['Ka1','Ka2','tlag12','F1']
[docs] def calculateParameterUnits(self,sample): tunits = self.ptrExperiment.getTimeUnits().unit self.parameterUnits = [inverseUnits(tunits),inverseUnits(tunits),tunits,PKPDUnit.UNIT_NONE] return self.parameterUnits
[docs] def getAg(self,t): if t<0: return 0.0 Ka1 = self.parameters[0] Ka2 = self.parameters[1] tlag12 = self.parameters[2] F1 = self.parameters[3] A1=F1*math.exp(-Ka1*t) A2=1-F1 if t>tlag12: A2=(1-F1)*math.exp(-Ka2*(t-tlag12)) return self.Amax*(A1+A2)
[docs] def getEquation(self): Ka1 = self.parameters[0] Ka2 = self.parameters[1] tlag12 = self.parameters[2] F1 = self.parameters[3] return "D(t)=(%f)*(1-exp(-(%f)*t)+(%f)*(1-exp(-(%f)*(t-%f))"%(self.Amax*F1,Ka1,self.Amax*(1-F1),Ka2,tlag12)
[docs] def getModelEquation(self): return "D(t)=Amax*F1*(1-exp(-Ka1*t))+(1-F1)*Amax*(1-exp(-Ka2*(t-tlag12)))"
[docs] def getDescription(self): return "First and First order absorption (%s)"%self.__class__.__name__
[docs]class BiopharmaceuticsModelOrder1AndOrder1Saturable(BiopharmaceuticsModel): """ Srihari Gopal M.D., M.H.S., An Vermeulen Ph.D., Partha Nandy Ph.D., Paulien Ravenstijn Ph.D., Isaac Nuamah Ph.D., J.A. Buron Vidal M.D., Joris Berwaerts M.D., Adam Savitz M.D., Ph.D., David Hough M.D. & Mahesh N. Samtani Ph.D. Practical Guidance for Dosing and Switching from Paliperidone Palmitate 1-Monthly to 3-Monthly Formulation in Schizophrenia, Current Medical Research and Opinion, 2015, 31: 2043-2054. There are two absorptions of order 1: one is slow (fraction=1-F3), the other is rapid (fraction=F3). Both absorptions are saturable by the amount present in the absorption compartment so that the instantaneous absorption rate depends on the amount remaining. """ def getDescription(self): return ['Maximal absorption rate slow','Amount at 50% of slow absorption rate','Hill factor', 'Maximal absorption rate rapid','Amount at 50% of rapid absorption rate', 'Fraction 3']
[docs] def getParameterNames(self): return ['ka1max','kamt150','gamma','ka3max','kamt350','F3']
[docs] def calculateParameterUnits(self,sample): tunits = self.ptrExperiment.getTimeUnits().unit Dunits = self.ptrExperiment.getDoseUnits() kaUnits = divideUnits(Dunits,tunits) self.parameterUnits = [kaUnits, Dunits, PKPDUnit.UNIT_NONE, kaUnits, Dunits, PKPDUnit.UNIT_NONE] return self.parameterUnits
[docs] def getAg(self,t): if t<0: return 0.0 ka1max = self.parameters[0] kamt150 = self.parameters[1] gamma = self.parameters[2] ka3max = self.parameters[3] kamt350 = self.parameters[4] F3 = self.parameters[5] if gamma<0 or F3<0 or F3>1: return 0.0 try: A50 = kamt150 Ka = ka1max A0 = F3*self.Amax arg = math.exp(-(Ka * t - A50 * math.log(A0 * math.exp(A0 / A50))) / A50) / A50 Aslow = A50 * lambertw(arg) except: Aslow = 0.0 try: A50 = kamt350 Ka = ka3max A0 = (1-F3)*self.Amax A0gamma = math.pow(A0,gamma) A50gamma = math.pow(A50,gamma) A0log = math.log(A0) Arapid = math.exp(((A0gamma + A50gamma * gamma * A0log) / gamma - Ka * t) / A50gamma - wrightomega(math.log(1 / A50gamma) + (gamma * ((A0gamma + A50gamma * gamma * A0log) / gamma - Ka * t)) / A50gamma) / gamma) except: Arapid = 0.0 retval = Aslow+Arapid if np.iscomplexobj(retval): retval = np.real(retval) return retval
[docs] def getEquation(self): ka1max = self.parameters[0] kamt150 = self.parameters[1] gamma = self.parameters[2] ka3max = self.parameters[3] kamt350 = self.parameters[4] F3 = self.parameters[5] A50 = kamt150 Ka = ka1max A0 = F3*self.Amax C = A50 * math.log(A0 * math.exp(A0 / A50)) Aslow = "%f * lambertw(exp(-(%f * t - (%f)) / %f) / %f)"%(A50,Ka,C,A50,A50) A50 = kamt350 Ka = ka3max A0 = (1-F3)*self.Amax A0gamma = math.pow(A0,gamma) A50gamma = math.pow(A50,gamma) A0log = math.log(A0) C1 = (A0gamma + A50gamma * gamma * A0log) / gamma Arapid = "exp((%f - Ka * t) / %f - "%(C1,A50gamma) C2 = math.log(1 / A50gamma) C3 = (A0gamma + A50gamma * gamma * A0log) / gamma Arapid += "wrightOmega((%f) + (%f * (%f - (%f) * t)) / %f) / %f)"%(C2,gamma,C3,Ka,A50gamma,gamma) return "D(t)=(%s)+(%s)"%(Aslow,Arapid)
[docs] def getModelEquation(self): return "D(t)=kamt150 * lambertw(exp(-(ka1max * t - kamt150 * log(F3*Amax * exp(F3*Amax / kamt150))) / kamt150) / kamt150) "+\ "exp(((((1-F3)*Amax)^gamma + kamt350^gamma*gamma*log((1-F3)*Amax))/gamma - ka3max*t)/kamt350^gamma - wrightOmega(log(1/kamt350^gamma) + (gamma*((((1-F3)*Amax)^gamma + kamt350^gamma*gamma*log((1-F3)*Amax))/gamma - ka3max*t))/kamt350^gamma)/gamma)"
[docs] def getDescription(self): return "First and First order absorption with saturation (%s)"%self.__class__.__name__
[docs]class BiopharmaceuticsModelDoubleWeibull(BiopharmaceuticsModel): def getDescription(self): return ['Time1 63%','Exponential power1','Fraction1', 'Time2 63%','Exponential power2','tlag2']
[docs] def getParameterNames(self): return ['td1','b1','F1','td2','b2','tlag2']
[docs] def calculateParameterUnits(self,sample): tunits = self.ptrExperiment.getTimeUnits().unit self.parameterUnits = [tunits,PKPDUnit.UNIT_NONE,PKPDUnit.UNIT_NONE, tunits,PKPDUnit.UNIT_NONE,tunits] return self.parameterUnits
[docs] def getAg(self,t): if t<=0: return 0.0 td1 = self.parameters[0] b1 = self.parameters[1] F1 = self.parameters[2] td2 = self.parameters[3] b2 = self.parameters[4] tlag2 = self.parameters[5] t2 = t-tlag2 f1 = F1*math.exp(-math.pow(t/td1,b1)) if t2>0: f2=(1-F1)*math.exp(-math.pow(t2/td2,b2)) else: f2=0 return self.Amax*(f1+f2)
[docs] def getEquation(self): td1 = self.parameters[0] b1 = self.parameters[1] F1 = self.parameters[2] td2 = self.parameters[3] b2 = self.parameters[4] tlag2 = self.parameters[5] return "D(t)=%f*(1-(%f)*exp(-(t/%f)^(%f))-(%f)*exp(-((t-%f)/%f)^(%f))"%(self.Amax,F1,td1,b1,1-F1,tlag2,td2,b2)
[docs] def getModelEquation(self): return "D(t)=Amax*(1-F1*exp(-(t/td1)^b1)-(1-F1)*exp(-((t-tlag2)/td2)^b2)"
[docs] def getDescription(self): return "Double Weibull absorption (%s)"%self.__class__.__name__
[docs] def areParametersValid(self, p): return np.sum(p<0)==0 and p[2]>0 and p[2]<1
[docs]class BiopharmaceuticsModelTripleWeibull(BiopharmaceuticsModel): def getDescription(self): return ['Time1 63%','Exponential power1','Fraction1', 'Time2 63%','Exponential power2','tlag2', 'Fraction2', 'Time3 63%','Exponential power3','tlag3']
[docs] def getParameterNames(self): return ['td1','b1','F1','td2','b2','tlag2','F2','td3','b3','tlag3']
[docs] def calculateParameterUnits(self,sample): tunits = self.ptrExperiment.getTimeUnits().unit self.parameterUnits = [tunits,PKPDUnit.UNIT_NONE,PKPDUnit.UNIT_NONE, tunits,PKPDUnit.UNIT_NONE,tunits, PKPDUnit.UNIT_NONE, tunits,PKPDUnit.UNIT_NONE,tunits] return self.parameterUnits
[docs] def getAg(self,t): if t<=0: return 0.0 td1 = self.parameters[0] b1 = self.parameters[1] F1 = self.parameters[2] td2 = self.parameters[3] b2 = self.parameters[4] tlag2 = self.parameters[5] F2 = self.parameters[6] td3 = self.parameters[7] b3 = self.parameters[8] tlag3 = self.parameters[9] f1 = F1*math.exp(-math.pow(t/td1,b1)) t2 = t-tlag2 if t2>0: f2=F2*math.exp(-math.pow(t2/td2,b2)) else: f2=0 t3 = t-tlag3 if t3>0: f3=(1-F1-F2)*math.exp(-math.pow(t3/td3,b3)) else: f3=0 return self.Amax*(f1+f2+f3)
[docs] def getEquation(self): td1 = self.parameters[0] b1 = self.parameters[1] F1 = self.parameters[2] td2 = self.parameters[3] b2 = self.parameters[4] tlag2 = self.parameters[5] F2 = self.parameters[6] td3 = self.parameters[7] b3 = self.parameters[8] tlag3 = self.parameters[9] return "D(t)=%f*(1-(%f)*exp(-(t/%f)^(%f))-(%f)*exp(-((t-%f)/%f)^(%f)-(%f)*exp(-((t-%f)/%f)^(%f)))"%\ (self.Amax,F1,td1,b1,F2,tlag2,td2,b2,1-F1-F2,tlag3,td3,b3)
[docs] def getModelEquation(self): return "D(t)=Amax*(1-F1*exp(-(t/td1)^b1)-F2*exp(-((t-tlag2)/td2)^b2)-(1-F1-F2)*exp(-((t-tlag3)/td3)^b3))"
[docs] def getDescription(self): return "Triple Weibull absorption (%s)"%self.__class__.__name__
[docs] def areParametersValid(self, p): return np.sum(p<0)==0 and p[2]>0 and p[2]<1 and p[6]>0 and p[6]<1
[docs]class BiopharmaceuticsModelOrder1Multiple4(BiopharmaceuticsModel): """ Br J Clin Pharmacol. 2013 Dec;76(6):868-79. doi: 10.1111/bcp.12118. Determination of the pharmacokinetics of glycopyrronium in the lung using a population pharmacokinetic modelling approach. Christian Bartels 1, Michael Looby, Romain Sechaud, Guenther Kaiser There are 4 entry vias: - F1, Order1: F1*(1-exp(-Ka1*t)) - (1-F1), Remaining: - Slow: (1-F1)*Fslow*(1-exp(-Kaslow*t)) - Medium: (1-F1)*Fmed*(1-exp(-Kamed*t)) - Fast: (1-F1)*(1-Fmed-Fslow) """ def getDescription(self): return ['Fraction 1','Absorption rate 1','Fraction medium','Absorption rate medium','Fraction slow', 'Absorption rate slow']
[docs] def getParameterNames(self): return ['F1','Ka1','Fmed','Kamed','Fslow','Kaslow']
[docs] def calculateParameterUnits(self,sample): Kunits = inverseUnits(self.ptrExperiment.getTimeUnits().unit) self.parameterUnits = [PKPDUnit.UNIT_NONE, Kunits, PKPDUnit.UNIT_NONE, Kunits, PKPDUnit.UNIT_NONE, Kunits] return self.parameterUnits
[docs] def getAg(self,t): if t<0: return 0.0 F1 = self.parameters[0] Ka1 = self.parameters[1] Fmed = self.parameters[2] Kamed = self.parameters[3] Fslow = self.parameters[4] Kaslow = self.parameters[5] via1 = F1*math.exp(-Ka1*t) viafast = (1-F1)*(1-Fmed-Fslow) viamed = (1-F1)*Fmed*math.exp(-Kamed*t) viaslow = (1-F1)*Fslow*math.exp(-Kaslow*t) return self.Amax*(via1 + viafast + viamed + viaslow)
[docs] def getEquation(self): F1 = self.parameters[0] Ka1 = self.parameters[1] Fmed = self.parameters[2] Kamed = self.parameters[3] Fslow = self.parameters[4] Kaslow = self.parameters[5] return "D(t)=(%f)*(%f*(1-exp(-%f*t)) + (1-%f)*(1-%f*exp(-%f*t)-%f*exp(-%f*t)))"%\ (self.Amax,F1,Ka1,F1,Fmed,Kamed,Fslow,Kaslow)
[docs] def getModelEquation(self): return "D(t)=F1*(1-exp(-Ka1*t)) + (1-F1)*(1-Fmed*exp(-Kamed*t)-Fslow*exp(-Kaslow*t))"
[docs] def getDescription(self): return "Multiple first order absorption (%s)"%self.__class__.__name__
[docs]class BiopharmaceuticsModelSplineGeneric(BiopharmaceuticsModel): def __init__(self): self.nknots=0 self.parametersPrepared=None def getDescription(self): return ['B-spline model with %d knots'%self.nknots]
[docs] def getParameterNames(self): retval = ['tmax'] if self.nknots>=10: retval += ['spline%d_A%02d' % (self.nknots, i) for i in range(self.nknots)] else: retval+=['spline%d_A%d'%(self.nknots,i) for i in range(self.nknots)] return retval
[docs] def calculateParameterUnits(self,sample): self.parameterUnits = [self.ptrExperiment.getTimeUnits().unit] self.parameterUnits += [PKPDUnit.UNIT_NONE]*(self.nknots) return self.parameterUnits
[docs] def rearrange(self,parameters): retval = parameters retval[1:]=np.sort(retval[1:]) return retval
[docs] def getAg(self,t): if t<=0: return self.Amax self.tmax=self.parameters[0] if t>=self.tmax or self.tmax<=0: return 0.0 if self.parametersPrepared is None or not np.array_equal(self.parametersPrepared,self.parameters): self.knots = np.linspace(0, self.tmax, self.nknots+2) self.parameters[1:] = np.sort(self.parameters[1:]) self.knotsY = np.append(np.insert(self.parameters[1:],0,0),1) self.knotsY = np.sort(self.knotsY) knotsUnique, knotsYUnique=uniqueFloatValues(self.knots, self.knotsY) try: self.B=PchipInterpolator(knotsUnique, knotsYUnique) # self.B=InterpolatedUnivariateSpline(knotsUnique, knotsYUnique, k=1) except: print("self.tmax",self.tmax) print("self.nknots",self.nknots) print("self.parameters[1:]",self.parameters[1:]) print("Error en spline",self.knots, self.knotsY, knotsUnique, knotsYUnique) raise Exception("Bug in spline") self.parametersPrepared=copy.copy(self.parameters) fraction=self.B(t) fraction=np.clip(fraction,0.0,1.0) #print("getAg t= %f B(t)= %f fraction= %f Ag= %f"%(t,self.B(t),fraction,self.Amax*(1-fraction))) return self.Amax*(1-fraction)
[docs] def getEquation(self): self.knotsY=np.sort(self.knotsY) retval="D(t) interpolating spline at x=%s and y=%s"%(np.array2string(self.knots,max_line_width=10000),np.array2string(self.knotsY*self.Amax,max_line_width=10000)) return retval
[docs] def getModelEquation(self): # https://en.wikipedia.org/wiki/De_Boor%27s_algorithm return "D(t)=interpolating BSpline1 with %d knots distributed until tmax"%self.nknots
[docs] def getDescription(self): return "BSplines with %d knots (%s)"%(self.nknots,self.__class__.__name__)
[docs] def areParametersValid(self, p): return np.sum(p<0)==0 and np.sum(p[1:]>1)==0
[docs]class BiopharmaceuticsModelSpline2(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 2
[docs]class BiopharmaceuticsModelSpline3(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 3
[docs]class BiopharmaceuticsModelSpline4(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 4
[docs]class BiopharmaceuticsModelSpline5(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 5
[docs]class BiopharmaceuticsModelSpline6(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 6
[docs]class BiopharmaceuticsModelSpline7(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 7
[docs]class BiopharmaceuticsModelSpline8(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 8
[docs]class BiopharmaceuticsModelSpline9(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 9
[docs]class BiopharmaceuticsModelSpline10(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 10
[docs]class BiopharmaceuticsModelSpline11(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 11
[docs]class BiopharmaceuticsModelSpline12(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 12
[docs]class BiopharmaceuticsModelSpline13(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 13
[docs]class BiopharmaceuticsModelSpline14(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 14
[docs]class BiopharmaceuticsModelSpline15(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 15
[docs]class BiopharmaceuticsModelSpline16(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 16
[docs]class BiopharmaceuticsModelSpline17(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 17
[docs]class BiopharmaceuticsModelSpline18(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 18
[docs]class BiopharmaceuticsModelSpline19(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 19
[docs]class BiopharmaceuticsModelSpline20(BiopharmaceuticsModelSplineGeneric): def __init__(self): BiopharmaceuticsModelSplineGeneric.__init__(self) self.nknots = 20
[docs]class BiopharmaceuticsModelSplineXYGeneric(BiopharmaceuticsModel): def __init__(self): self.nknots=0 self.parametersPrepared=None def getDescription(self): return ['B-splineXY model with %d knots'%self.nknots]
[docs] def getParameterNames(self): retval = ['tmax'] for i in range(self.nknots): retval += ['spline%d_X%d' % (self.nknots, i)] retval += ['spline%d_A%d' % (self.nknots, i)] return retval
[docs] def calculateParameterUnits(self,sample): self.parameterUnits = [self.ptrExperiment.getTimeUnits().unit] self.parameterUnits += [self.ptrExperiment.getTimeUnits().unit,PKPDUnit.UNIT_NONE]*(self.nknots) return self.parameterUnits
[docs] def rearrange(self,parameters): retval = parameters retval[1::2]=np.sort(retval[1::2]) retval[2::2]=np.sort(retval[2::2]) return retval
[docs] def getAg(self,t): if t<=0: return self.Amax self.tmax=self.parameters[0] if t>=self.tmax or self.tmax<=0: return 0.0 if self.parametersPrepared is None or not np.array_equal(self.parametersPrepared,self.parameters): self.parameters[1::2]=np.sort(self.parameters[1::2]) self.parameters[2::2]=np.sort(self.parameters[2::2]) self.knots = np.append(np.insert(self.parameters[1::2],0,0),1)*self.tmax self.knotsY = np.append(np.insert(self.parameters[2::2],0,0),1) self.knots=np.sort(self.knots) self.knotsY=np.sort(self.knotsY) knotsUnique, knotsYUnique=uniqueFloatValues(self.knots, self.knotsY) try: self.B=PchipInterpolator(knotsUnique, knotsYUnique) # self.B=InterpolatedUnivariateSpline(knotsUnique, knotsYUnique, k=1) except: print("self.tmax",self.tmax) print("self.nknots",self.nknots) print("self.parameters[1:]",self.parameters[1:]) print("Error en splineXY",self.knots, self.knotsY, knotsUnique, knotsYUnique) raise Exception("Bug in spline") self.parametersPrepared=copy.copy(self.parameters) fraction=self.B(t) fraction=np.clip(fraction,0.0,1.0) #print("getAg t= %f B(t)= %f fraction= %f Ag= %f"%(t,self.B(t),fraction,self.Amax*(1-fraction))) return self.Amax*(1-fraction)
[docs] def getEquation(self): self.knotsY=np.sort(self.knotsY) retval="D(t) interpolating spline at x=%s and y=%s"%(np.array2string(self.knots,max_line_width=10000),np.array2string(self.knotsY*self.Amax,max_line_width=10000)) return retval
[docs] def getModelEquation(self): # https://en.wikipedia.org/wiki/De_Boor%27s_algorithm return "D(t)=interpolating BSpline1 with %d knots distributed until tmax"%self.nknots
[docs] def getDescription(self): return "BSplinesXY with %d knots (%s)"%(self.nknots,self.__class__.__name__)
[docs] def areParametersValid(self, p): return np.sum(p<0)==0 and np.sum(p[1:]>1)==0
[docs]class BiopharmaceuticsModelSplineXY2(BiopharmaceuticsModelSplineXYGeneric): def __init__(self): BiopharmaceuticsModelSplineXYGeneric.__init__(self) self.nknots = 2
[docs]class BiopharmaceuticsModelSplineXY3(BiopharmaceuticsModelSplineXYGeneric): def __init__(self): BiopharmaceuticsModelSplineXYGeneric.__init__(self) self.nknots = 3
[docs]class BiopharmaceuticsModelSplineXY4(BiopharmaceuticsModelSplineXYGeneric): def __init__(self): BiopharmaceuticsModelSplineXYGeneric.__init__(self) self.nknots = 4
[docs]class BiopharmaceuticsModelSplineXY5(BiopharmaceuticsModelSplineXYGeneric): def __init__(self): BiopharmaceuticsModelSplineXYGeneric.__init__(self) self.nknots = 5
[docs]class BiopharmaceuticsModelSplineXY6(BiopharmaceuticsModelSplineXYGeneric): def __init__(self): BiopharmaceuticsModelSplineXYGeneric.__init__(self) self.nknots = 6
[docs]class BiopharmaceuticsModelSplineXY7(BiopharmaceuticsModelSplineXYGeneric): def __init__(self): BiopharmaceuticsModelSplineXYGeneric.__init__(self) self.nknots = 7
[docs]class BiopharmaceuticsModelSplineXY8(BiopharmaceuticsModelSplineXYGeneric): def __init__(self): BiopharmaceuticsModelSplineXYGeneric.__init__(self) self.nknots = 8
[docs]class BiopharmaceuticsModelSplineXY9(BiopharmaceuticsModelSplineXYGeneric): def __init__(self): BiopharmaceuticsModelSplineXYGeneric.__init__(self) self.nknots = 9
[docs]class BiopharmaceuticsModelSplineXY10(BiopharmaceuticsModelSplineXYGeneric): def __init__(self): BiopharmaceuticsModelSplineXYGeneric.__init__(self) self.nknots = 10
[docs]class BiopharmaceuticsModelNumerical(BiopharmaceuticsModel):
[docs] def setXYValues(self,t,A): # A is the accumulated fraction released tUnique, Aunique = uniqueFloatValues(t,np.asarray(A,dtype=np.float64)/100.0) # self.B = InterpolatedUnivariateSpline(tUnique, Aunique,k=1) self.B = PchipInterpolator(tUnique, Aunique) self.tmin=np.min(t) self.tmax=np.max(t)
def getDescription(self): return ['Numerical source with t and A']
[docs] def getParameterNames(self): return []
[docs] def calculateParameterUnits(self,sample): self.parameterUnits = [] return self.parameterUnits
[docs] def getAg(self,t): if t<=0: return self.Amax tToUse=t if t>=self.tmax: tToUse=self.tmax fraction=self.B(tToUse) fraction=np.clip(fraction,0.0,1.0) return self.Amax*(1-fraction)
[docs] def getEquation(self): return 'Numerical source with t and A'
[docs] def getModelEquation(self): return 'Numerical source with t and A'
[docs] def getDescription(self): return "Numerical source (%s)"%self.__class__.__name__
[docs] def areParametersValid(self, p): return True
[docs]class PKPDVia: def __init__(self, ptrExperiment = None): self.viaName = None self.via = None self.viaProfile = None self.tlag = 0 self.tunits = None self.bioavailability = 1 self.paramsToOptimize = [] self.paramsUnitsToOptimize=[] self.ptrExperiment = ptrExperiment
[docs] def parseTokens(self,tokens): # Intravenous; iv; [tlag=0 min]; [bioavailability=1] # Oral; ev1; tlag=0 min; bioavailability=1 # Default values self.tlag = 0 self.bioavailability = 1 if self.ptrExperiment is not None: self.tunits = self.ptrExperiment.getTimeUnits() else: self.tunits = PKPDUnit("min") # Get name currentToken = 0 self.viaName = tokens[currentToken].strip() currentToken+=1 # Get via self.via = tokens[currentToken].strip() currentToken+=1 while currentToken<len(tokens): optionalTokens=tokens[currentToken].strip() if '=' in optionalTokens: optionalTokens=optionalTokens.split('=') optionalVar=optionalTokens[0].strip() if optionalVar=="tlag": optionalTokens=optionalTokens[1].split() self.tlag=float(optionalTokens[0].strip()) unitString = optionalTokens[1].strip() self.tunits = PKPDUnit(unitString) if not self.tunits.unit: raise Exception("Unrecognized unit: %s"%unitString) if not self.tunits.isTime(): raise Exception("Time unit is not valid") elif optionalVar=="bioavailability": self.bioavailability=float(optionalTokens[1].strip()) else: optionalTokens=optionalTokens.split()[0] self.paramsToOptimize.append(optionalTokens) if optionalTokens=="tlag": self.paramsUnitsToOptimize.append(self.tunits.unit) elif optionalTokens=="bioavailability": self.paramsUnitsToOptimize.append(PKPDUnit.UNIT_NONE) currentToken+=1
def _printToStrig(self): outStr="%s; %s; "%(self.viaName,self.via) if "tlag" in self.paramsToOptimize: outStr+=" tlag min; " else: outStr+="tlag=%f %s; "%(self.tlag,self.tunits._toString()) if "bioavailability" in self.paramsToOptimize: outStr+=" bioavailability" else: outStr+="bioavailability=%f"%self.bioavailability return outStr def _printToStream(self,fh): fh.write("%s\n"%self._printToStrig()) def _printToExcel(self, wb, row): excelWriteRow(self._printToStrig(), wb, row) return row+1
[docs] def prepare(self): if self.via=="iv": self.viaProfile=None else: if self.viaProfile==None: if self.via=="iv-ev1": self.viaProfile=BiopharmaceuticsModelImmediateAndOrder1() elif self.via=="ev0": self.viaProfile=BiopharmaceuticsModelOrder0() elif self.via=="ev01": self.viaProfile=BiopharmaceuticsModelOrder01() elif self.via == "ev0tlag1": self.viaProfile = BiopharmaceuticsModelOrder01Tlag1() elif self.via=="ev1": self.viaProfile=BiopharmaceuticsModelOrder1() elif self.via == "ev1x4": self.viaProfile = BiopharmaceuticsModelOrder1Multiple4() elif self.via=="evFractional": self.viaProfile=BiopharmaceuticsModelOrderFractional() elif self.via=="ev1-ev1": self.viaProfile=BiopharmaceuticsModelOrder1AndOrder1() elif self.via == "ev1-ev1Saturable": self.viaProfile = BiopharmaceuticsModelOrder1AndOrder1Saturable() elif self.via == "doubleWeibull": self.viaProfile = BiopharmaceuticsModelDoubleWeibull() elif self.via == "tripleWeibull": self.viaProfile = BiopharmaceuticsModelTripleWeibull() elif self.via=="spline2": self.viaProfile=BiopharmaceuticsModelSpline2() elif self.via=="spline3": self.viaProfile=BiopharmaceuticsModelSpline3() elif self.via=="spline4": self.viaProfile=BiopharmaceuticsModelSpline4() elif self.via=="spline5": self.viaProfile=BiopharmaceuticsModelSpline5() elif self.via=="spline6": self.viaProfile=BiopharmaceuticsModelSpline6() elif self.via=="spline7": self.viaProfile=BiopharmaceuticsModelSpline7() elif self.via=="spline8": self.viaProfile=BiopharmaceuticsModelSpline8() elif self.via=="spline9": self.viaProfile=BiopharmaceuticsModelSpline9() elif self.via=="spline10": self.viaProfile=BiopharmaceuticsModelSpline10() elif self.via == "spline11": self.viaProfile = BiopharmaceuticsModelSpline11() elif self.via == "spline12": self.viaProfile = BiopharmaceuticsModelSpline12() elif self.via == "spline13": self.viaProfile = BiopharmaceuticsModelSpline13() elif self.via == "spline14": self.viaProfile = BiopharmaceuticsModelSpline14() elif self.via == "spline15": self.viaProfile = BiopharmaceuticsModelSpline15() elif self.via == "spline16": self.viaProfile = BiopharmaceuticsModelSpline16() elif self.via == "spline17": self.viaProfile = BiopharmaceuticsModelSpline17() elif self.via == "spline18": self.viaProfile = BiopharmaceuticsModelSpline18() elif self.via == "spline19": self.viaProfile = BiopharmaceuticsModelSpline19() elif self.via == "spline20": self.viaProfile = BiopharmaceuticsModelSpline20() elif self.via == "splineXY2": self.viaProfile = BiopharmaceuticsModelSplineXY2() elif self.via == "splineXY3": self.viaProfile = BiopharmaceuticsModelSplineXY3() elif self.via == "splineXY4": self.viaProfile = BiopharmaceuticsModelSplineXY4() elif self.via == "splineXY5": self.viaProfile = BiopharmaceuticsModelSplineXY5() elif self.via == "splineXY6": self.viaProfile = BiopharmaceuticsModelSplineXY6() elif self.via == "splineXY7": self.viaProfile = BiopharmaceuticsModelSplineXY7() elif self.via == "splineXY8": self.viaProfile = BiopharmaceuticsModelSplineXY8() elif self.via == "splineXY9": self.viaProfile = BiopharmaceuticsModelSplineXY9() elif self.via == "splineXY10": self.viaProfile = BiopharmaceuticsModelSplineXY10() elif self.via=="numerical": self.viaProfile=BiopharmaceuticsModelNumerical() self.viaProfile.setExperiment(self.ptrExperiment)
[docs] def changeTimeUnitsToMinutes(self): if self.tunits.unit==PKPDUnit.UNIT_TIME_MIN: pass elif self.tunits.unit==PKPDUnit.UNIT_TIME_H: self.tlag *= 60 else: raise Exception("Time for doses must be hours or minutes")
[docs] def getEquation(self): if self.via == "iv": retval = "D=Div(t)" else: retval = self.viaProfile.getEquation() return "%s (tlag=%f, bioavailability=%f)"%(retval,self.tlag,self.bioavailability)
[docs] def getModelEquation(self): if self.via == "iv": return "D=Div(t)" else: return self.viaProfile.getModelEquation()
[docs] def getDescription(self): if self.via == "iv": return "Intravenous dose" else: return self.viaProfile.getDescription()
[docs] def getParameterNames(self): names=copy.copy(self.paramsToOptimize) if self.via != "iv": names+=self.viaProfile.getParameterNames() return [self.viaName+"_"+x for x in names]
[docs] def getNumberOfParameters(self): return len(self.getParameterNames())
[docs] def calculateParameterUnits(self,sample): if self.via == "iv": return self.paramsUnitsToOptimize else: return self.paramsUnitsToOptimize+self.viaProfile.calculateParameterUnits(sample)
[docs] def areParametersSignificant(self, lowerBound, upperBound): retval=[] currentIdx=0 if self.paramsToOptimize: for paramName in self.paramsToOptimize: if paramName=="tlag": retval+=[str(lowerBound[currentIdx]>0)] elif paramName=="bioavailability": retval+=[str(upperBound[currentIdx]<1)] currentIdx+=1 if self.via == "iv": return retval else: return np.concatenate((np.asarray(retval),self.viaProfile.areParametersSignificant(lowerBound[currentIdx:], upperBound[currentIdx:])))
[docs] def areParametersValid(self, p): retval=True currentIdx=0 if self.paramsToOptimize: for paramName in self.paramsToOptimize: if paramName=="tlag": retval=retval and p[currentIdx]>=0 elif paramName=="bioavailability": retval=retval and p[currentIdx]>=0 and p[currentIdx]<=1 currentIdx+=1 if self.via == "iv": return retval else: return retval and self.viaProfile.areParametersValid(p[currentIdx:])
[docs] def setParameters(self, p): currentIdx=0 if self.paramsToOptimize: for paramName in self.paramsToOptimize: if paramName=="tlag": self.tlag=p[currentIdx] elif paramName=="bioavailability": self.bioavailability=p[currentIdx] currentIdx+=1 if self.via != "iv": self.viaProfile.setParameters(p[currentIdx:])
[docs]def createVia(line, ptrExperiment=None): via = PKPDVia(ptrExperiment) via.parseTokens(line.split(';')) via.prepare() return via
[docs]class PKPDDose: TYPE_BOLUS = 1 TYPE_REPEATED_BOLUS = 2 TYPE_INFUSION = 3 def __init__(self): self.doseName = None self.via = None self.doseType = None self.doseAmount = None self.t0 = None self.tF = None self.every = None self.tunits = None self.dunits = None self.paramsToOptimize = [] self.paramsUnitsToOptimize=[]
[docs] def parseTokens(self,tokens,vias): # Dose1; via=Intravenous; bolus; t=0 min; d=60*$(weight)/1000 mg # Dose1; via=Oral; repeated_bolus; t=0:8:48 h; d=60*$(weight)/1000 mg # Dose1; via=Intravenous; infusion; t=0:59 min; d=1 mg/min # Get name currentToken = 0 self.doseName = tokens[currentToken].strip() currentToken+=1 # Get via viaName = tokens[currentToken].strip().split('=')[1] if viaName in vias: self.via=vias[viaName] else: raise Exception("Unrecognized via %s"%viaName) currentToken+=1 # Get type doseTypeString = tokens[currentToken].strip() if doseTypeString=="bolus": self.doseType = PKPDDose.TYPE_BOLUS elif doseTypeString=="repeated_bolus": self.doseType = PKPDDose.TYPE_REPEATED_BOLUS elif doseTypeString=="infusion": self.doseType = PKPDDose.TYPE_INFUSION else: raise Exception("Unrecognized dose type %s"%doseTypeString) currentToken+=1 # Get time description timeUnitsString = tokens[currentToken].strip().lower() timeTokens=timeUnitsString.split() if len(timeTokens)!=2: raise Exception("Time description is badly formed %s"%timeUnitsString) timeString = timeTokens[0].strip().split('=')[1] unitString = timeTokens[1].strip() if doseTypeString=="bolus": self.t0 = float(timeString) elif doseTypeString=="repeated_bolus": timeTokens = timeString.split(":") self.t0 = float(timeTokens[0].strip()) self.every = float(timeTokens[1].strip()) self.tF = float(timeTokens[2].strip()) elif doseTypeString=="infusion": timeTokens = timeString.split(":") self.t0 = float(timeTokens[0].strip()) self.tF = float(timeTokens[1].strip()) else: raise Exception("Unrecognized dose type %s"%doseTypeString) # Get time units self.tunits = PKPDUnit(unitString) if not self.tunits.unit: raise Exception("Unrecognized unit: %s"%unitString) if not self.tunits.isTime(): raise Exception("Time unit is not valid") currentToken+=1 # Get dose units doseUnitsString = tokens[currentToken].strip().lower() doseTokens=doseUnitsString.split() if len(doseTokens)!=2: raise Exception("Dose description is badly formed %s"%doseUnitsString) self.doseAmount = doseTokens[0].strip().lower().split("=")[1] unitString = doseTokens[1].strip() self.dunits = PKPDUnit(unitString) if not self.dunits.unit: raise Exception("Unrecognized unit: %s"%unitString) if doseTypeString=="infusion": if not self.dunits.isWeightInvTime(): raise Exception("After normalization, the dose must be a weight/time (=rate)") else: if not self.dunits.isWeight(): raise Exception("After normalization, the dose must be a weight") currentToken+=1
def _printToStream(self,fh): outStr=self.doseName+"; "+self.getDoseString2() fh.write("%s\n"%outStr) def _printToExcel(self, wb, row): outStr=self.doseName+"; "+self.getDoseString2() excelWriteRow(outStr,wb,row) return row+1 def __str__(self): return self.getDoseString2()
[docs] def getDoseString(self): if self.doseType == PKPDDose.TYPE_BOLUS: doseString = "bolus; t=%f" % self.t0 elif self.doseType == PKPDDose.TYPE_REPEATED_BOLUS: doseString = "repeated_bolus; t=%f:%f:%f" % (self.t0, self.every, self.tF) elif self.doseType == PKPDDose.TYPE_INFUSION: doseString = "infusion; t=%f:%f" % (self.t0, self.tF) else: doseString = "" return doseString+" "+self.tunits._toString()
[docs] def getDoseString2(self): outStr="via=%s; %s; d=%s %s" % (self.via.viaName, self.getDoseString(), self.doseAmount, self.dunits._toString()) return outStr
[docs] def changeTimeUnitsTo(self, targetUnit): if targetUnit==PKPDUnit.UNIT_TIME_MIN: if self.tunits.unit==PKPDUnit.UNIT_TIME_MIN: pass elif self.tunits.unit==PKPDUnit.UNIT_TIME_H: if self.doseType == PKPDDose.TYPE_BOLUS: self.t0 *= 60 elif self.doseType == PKPDDose.TYPE_REPEATED_BOLUS: self.t0 *= 60 self.tF *= 60 self.every *= 60 elif self.doseType == PKPDDose.TYPE_INFUSION: self.t0 *= 60 self.tF *= 60 else: raise Exception("Time for doses must be hours or minutes") elif targetUnit==PKPDUnit.UNIT_TIME_H: if self.tunits.unit==PKPDUnit.UNIT_TIME_H: pass elif self.tunits.unit==PKPDUnit.UNIT_TIME_MIN: if self.doseType == PKPDDose.TYPE_BOLUS: self.t0 /= 60 elif self.doseType == PKPDDose.TYPE_REPEATED_BOLUS: self.t0 /= 60 self.tF /= 60 self.every /= 60 elif self.doseType == PKPDDose.TYPE_INFUSION: self.t0 /= 60 self.tF /= 60 else: raise Exception("Time for doses must be hours or minutes") else: raise Exception("Target time unit must be hours or minutes")
[docs] def prepare(self): self.via.prepare()
[docs] def getDoseAt(self,t0,dt=0.5): """Dose between t0<=t<t0+dt, t0 is in the units of the dose""" t0-=self.via.tlag t1=t0+dt-self.via.tlag if self.doseType == PKPDDose.TYPE_BOLUS: if t0<=self.t0 and self.t0<t1: return self.doseAmount else: return 0.0 elif self.doseType == PKPDDose.TYPE_REPEATED_BOLUS: doseAmount=0 for t in np.arange(self.t0,self.tF,self.every): if t0<=t and t<t1: doseAmount+=self.doseAmount return doseAmount elif self.doseType == PKPDDose.TYPE_INFUSION: if t0>self.tF or t1<self.t0: return 0.0 else: tLeft=max(t0,self.t0) tRight=min(t1,self.tF) return self.doseAmount*(tRight-tLeft)
[docs] def getAmountReleasedAt(self,t0,dt=0.5): doseAmount = 0.0 if self.via.viaProfile == None: doseAmount += self.getDoseAt(t0,dt) else: if self.doseType!=PKPDDose.TYPE_INFUSION: self.via.viaProfile.Amax = self.via.bioavailability*self.doseAmount doseAmount += self.via.viaProfile.getAg(t0-self.t0-self.via.tlag)-\ self.via.viaProfile.getAg(t0-self.t0-self.via.tlag+dt) else: doseAmount += self.getDoseAt(t0,dt) if doseAmount<0: doseAmount=0 # print("t0=%f self.t0=%f self.tlag=%f self.dt=%f -> released=%f"%(t0,self.t0,self.via.tlag,dt,doseAmount)) return doseAmount
[docs] def getAmountReleasedUpTo(self, t0): doseAmount = 0.0 if self.via.viaProfile == None: doseAmount += self.getDoseAt(0.0,t0) else: if self.doseType!=PKPDDose.TYPE_INFUSION: self.via.viaProfile.Amax = self.via.bioavailability*self.doseAmount doseAmount += self.via.viaProfile.getAg(0.0)-\ self.via.viaProfile.getAg(t0-self.t0-self.via.tlag) else: doseAmount += self.getDoseAt(0,t0) if doseAmount<0: doseAmount=0 return doseAmount
[docs] def isDoseABolus(self): if self.doseType != PKPDDose.TYPE_BOLUS: return False if self.t0 != 0: return False return True
[docs] def getTUnitsString(self): return self.tunits._toString()
[docs] def getDUnitsString(self): return self.dunits._toString()
[docs] def getDoseUnits(self): return self.dunits
[docs]def createDeltaDose(doseAmount,via,t=0,dunits="mg"): dose = PKPDDose() dose.doseName = "Bolus" dose.doseType = PKPDDose.TYPE_BOLUS dose.doseAmount = doseAmount dose.via = via dose.t0 = t dose.tunits = PKPDUnit("min") dose.dunits = PKPDUnit(dunits) return dose
[docs]class DrugSource: def __init__(self): self.parsedDoseList = [] self.parameterNames = [] self.parameterUnits = [] self.vias = []
[docs] def setDoses(self, parsedDoseList, t0, tF): self.originalDoseList = parsedDoseList self.parsedDoseList = [] collectedVias = [] self.parameterNames = [] self.vias = [] for dose in parsedDoseList: if dose.doseType == PKPDDose.TYPE_BOLUS: self.parsedDoseList.append(dose) elif dose.doseType == PKPDDose.TYPE_INFUSION: newDose = copy.copy(dose) newDose.dunits = copy.copy(dose.dunits) newDose.dunits.unit = changeRateToWeight(dose.dunits.unit) self.parsedDoseList.append(newDose) else: for t in np.arange(dose.t0,dose.tF,dose.every): if t0<=t and t<=tF: newDose = copy.copy(dose) newDose.doseType = PKPDDose.TYPE_BOLUS newDose.t0 = t self.parsedDoseList.append(newDose) if not dose.via.viaName in collectedVias: collectedVias.append(dose.via.viaName) viaParameterNames = dose.via.getParameterNames() self.parameterNames += viaParameterNames self.vias.append((dose.via,len(viaParameterNames)))
[docs] def getDoseUnits(self): if len(self.parsedDoseList)==0: return PKPDUnit.UNIT_NONE else: return self.parsedDoseList[0].getDoseUnits()
[docs] def getAmountReleasedAt(self,t0,dt=0.5): doseAmount = 0.0 for dose in self.parsedDoseList: doseAmount+=dose.getAmountReleasedAt(t0,dt) return doseAmount
[docs] def getAmountReleasedUpTo(self,t0): doseAmount = 0.0 for dose in self.parsedDoseList: doseAmount+=dose.getAmountReleasedUpTo(t0) return doseAmount
[docs] def getEquation(self): retval = "" for via,_ in self.vias: retval+=via.getEquation()+" " return retval.strip()
[docs] def getModelEquation(self): retval = "" for via,_ in self.vias: retval+=via.getModelEquation()+" " return retval.strip()
[docs] def getDescription(self): retval = "" for via,_ in self.vias: retval+=via.getDescription()+" " return retval.strip()
[docs] def getParameterNames(self): return self.parameterNames
[docs] def getNumberOfParameters(self): return len(self.getParameterNames())
[docs] def calculateParameterUnits(self,sample): if len(self.parameterUnits)>0: return self.parameterUnits else: retval = [] for via,_ in self.vias: retval+=via.calculateParameterUnits(sample) self.parameterUnits=retval return retval
[docs] def areParametersSignificant(self, lowerBound, upperBound): retval = [] currentToken=0 for via, viaPrmNo in self.vias: retval+=[x for x in via.areParametersSignificant(lowerBound[currentToken:(currentToken+viaPrmNo)], upperBound[currentToken:(currentToken+viaPrmNo)])] currentToken+=viaPrmNo if retval: return retval else: return True
[docs] def areParametersValid(self, p): retval = True currentToken=0 for via,viaPrmNo in self.vias: retval=retval and via.areParametersValid(p[currentToken:currentToken+viaPrmNo]) currentToken+=viaPrmNo return retval
[docs] def setParameters(self, p): currentToken=0 for via,viaPrmNo in self.vias: via.setParameters(p[currentToken:currentToken+viaPrmNo]) currentToken+=viaPrmNo
[docs] def getVia(self): return self.vias[0][0] # Only the first one is accessible through this function
[docs] def getDprofile(self,t): D = np.zeros(t.shape) Nsamples = D.shape if type(Nsamples)==tuple: Nsamples=Nsamples[0] for i in range(0,Nsamples-1): D[i]=self.getAmountReleasedAt(t[i],t[i+1]-t[i]) return D