class TTIWave(WaveBase): """ Solving the 2D TTI wave equation with `sigma_xx= c11*e_xx + c13*e_zz + c15*e_xz` `sigma_zz= c13*e_xx + c33*e_zz + c35*e_xz` `sigma_xz= c15*e_xx + c35*e_zz + c55*e_xz` the coefficients `c11`, `c13`, etc are calculated from the tompsen parameters `eps`, `delta` and the tilt `theta` :note: currently only the 2D case is supported. """ def __init__(self, domain, v_p, v_s, wavelet, source_tag, source_vector = [0.,1.], eps=0., delta=0., theta=0., rho=1., dt=None, u0=None, v0=None, absorption_zone=300*U.m, absorption_cut=1e-2, lumping=True): """ initialize the TTI wave solver :param domain: domain of the problem :type domain: `Domain` :param v_p: vertical p-velocity field :type v_p: `Scalar` :param v_s: vertical s-velocity field :type v_s: `Scalar` :param wavelet: wavelet to describe the time evolution of source term :type wavelet: `Wavelet` :param source_tag: tag of the source location :type source_tag: 'str' or 'int' :param source_vector: source orientation vector :param eps: first Thompsen parameter :param delta: second Thompsen parameter :param theta: tilting (in Rad) :param rho: density :param dt: time step size. If not present a suitable time step size is calculated. :param u0: initial solution. If not present zero is used. :param v0: initial solution change rate. If not present zero is used. :param absorption_zone: thickness of absorption zone :param absorption_cut: boundary value of absorption decay factor :param lumping: if True mass matrix lumping is being used. This is accelerates the computing but introduces some diffusion. """ DIM=domain.getDim() if not DIM == 2: raise ValueError("Only 2D is supported.") f=createAbsorptionLayerFunction(Function(domain).getX(), absorption_zone, absorption_cut) v_p=v_p*f v_s=v_s*f if u0 == None: u0=Vector(0.,Solution(domain)) else: u0=interpolate(p0, Solution(domain )) if v0 == None: v0=Vector(0.,Solution(domain)) else: v0=interpolate(v0, Solution(domain )) if dt == None: dt=min((1./5.)*min(inf(domain.getSize()/v_p), inf(domain.getSize()/v_s)), wavelet.getTimeScale()) super(TTIWave, self).__init__( dt, u0=u0, v0=v0, t0=0.) self.__wavelet=wavelet self.__mypde=LinearPDESystem(domain) if lumping: self.__mypde.getSolverOptions().setSolverMethod(SolverOptions.HRZ_LUMPING) self.__mypde.setSymmetryOn() self.__mypde.setValue(D=rho*kronecker(DIM), X=self.__mypde.createCoefficient('X')) self.__source_tag=source_tag self.__r=Vector(0, DiracDeltaFunctions(self.__mypde.getDomain())) self.__r.setTaggedValue(self.__source_tag, source_vector) c0_33=v_p**2 * rho c0_66=v_s**2 * rho c0_11=(1+2*eps) * c0_33 c0_13=sqrt(2*c0_33*(c0_33-c0_66) * delta + (c0_33-c0_66)**2)-c0_66 self.c11= c0_11*cos(theta)**4 - 2*c0_13*cos(theta)**4 + 2*c0_13*cos(theta)**2 + c0_33*sin(theta)**4 - 4*c0_66*cos(theta)**4 + 4*c0_66*cos(theta)**2 self.c13= -c0_11*cos(theta)**4 + c0_11*cos(theta)**2 + c0_13*sin(theta)**4 + c0_13*cos(theta)**4 - c0_33*cos(theta)**4 + c0_33*cos(theta)**2 + 4*c0_66*cos(theta)**4 - 4*c0_66*cos(theta)**2 self.c16= (-2*c0_11*cos(theta)**2 - 4*c0_13*sin(theta)**2 + 2*c0_13 + 2*c0_33*sin(theta)**2 - 8*c0_66*sin(theta)**2 + 4*c0_66)*sin(theta)*cos(theta)/2 self.c33= c0_11*sin(theta)**4 - 2*c0_13*cos(theta)**4 + 2*c0_13*cos(theta)**2 + c0_33*cos(theta)**4 - 4*c0_66*cos(theta)**4 + 4*c0_66*cos(theta)**2 self.c36= (2*c0_11*cos(theta)**2 - 2*c0_11 + 4*c0_13*sin(theta)**2 - 2*c0_13 + 2*c0_33*cos(theta)**2 + 8*c0_66*sin(theta)**2 - 4*c0_66)*sin(theta)*cos(theta)/2 self.c66= -c0_11*cos(theta)**4 + c0_11*cos(theta)**2 + 2*c0_13*cos(theta)**4 - 2*c0_13*cos(theta)**2 - c0_33*cos(theta)**4 + c0_33*cos(theta)**2 + c0_66*sin(theta)**4 + 3*c0_66*cos(theta)**4 - 2*c0_66*cos(theta)**2 def _getAcceleration(self, t, u): """ returns the acceleraton for time `t` and solution `u` at time `t` """ du = grad(u) sigma=self.__mypde.getCoefficient('X') e_xx=du[0,0] e_zz=du[1,1] e_xz=du[0,1]+du[1,0] sigma[0,0]= self.c11 * e_xx + self.c13 * e_zz + self.c16 * e_xz sigma[1,1]= self.c13 * e_xx + self.c33 * e_zz + self.c36 * e_xz sigma_xz = self.c16 * e_xx + self.c36 * e_zz + self.c66 * e_xz sigma[0,1]=sigma_xz sigma[1,0]=sigma_xz self.__mypde.setValue(X=-sigma, y_dirac= self.__r * self.__wavelet.getValue(t)) return self.__mypde.getSolution()
class SonicHTIWave(WaveBase): """ Solving the HTI wave equation (along the x_0 axis) with azimuth (rotation around verticle axis) under the assumption of zero shear wave velocities The unknowns are the transversal (along x_0) and vertial stress (Q, P) :note: In case of a two dimensional domain the second spatial dimenion is depth. """ def __init__(self, domain, v_p, wavelet, source_tag, source_vector = [1.,0.], eps=0., delta=0., azimuth=0., dt=None, p0=None, v0=None, absorption_zone=300*U.m, absorption_cut=1e-2, lumping=True): """ initialize the HTI wave solver :param domain: domain of the problem :type domain: `Doamin` :param v_p: vertical p-velocity field :type v_p: `Scalar` :param v_s: vertical s-velocity field :type v_s: `Scalar` :param wavelet: wavelet to describe the time evolution of source term :type wavelet: `Wavelet` :param source_tag: tag of the source location :type source_tag: 'str' or 'int' :param source_vector: source orientation vector :param eps: first Thompsen parameter :param azimuth: azimuth (rotation around verticle axis) :param gamma: third Thompsen parameter :param rho: density :param dt: time step size. If not present a suitable time step size is calculated. :param p0: initial solution (Q(t=0), P(t=0)). If not present zero is used. :param v0: initial solution change rate. If not present zero is used. :param absorption_zone: thickness of absorption zone :param absorption_cut: boundary value of absorption decay factor :param lumping: if True mass matrix lumping is being used. This is accelerates the computing but introduces some diffusion. """ DIM=domain.getDim() f=createAbsorptionLayerFunction(v_p.getFunctionSpace().getX(), absorption_zone, absorption_cut) self.v2_p=v_p**2 self.v2_t=self.v2_p*sqrt(1+2*delta) self.v2_n=self.v2_p*(1+2*eps) if p0 == None: p0=Data(0.,(2,),Solution(domain)) else: p0=interpolate(p0, Solution(domain )) if v0 == None: v0=Data(0.,(2,),Solution(domain)) else: v0=interpolate(v0, Solution(domain )) if dt == None: dt=min(min(inf(domain.getSize()/sqrt(self.v2_p)), inf(domain.getSize()/sqrt(self.v2_t)), inf(domain.getSize()/sqrt(self.v2_n))) , wavelet.getTimeScale())*0.2 super(SonicHTIWave, self).__init__( dt, u0=p0, v0=v0, t0=0.) self.__wavelet=wavelet self.__mypde=LinearPDESystem(domain) if lumping: self.__mypde.getSolverOptions().setSolverMethod(SolverOptions.HRZ_LUMPING) self.__mypde.setSymmetryOn() self.__mypde.setValue(D=kronecker(2), X=self.__mypde.createCoefficient('X')) self.__source_tag=source_tag self.__r=Vector(0, DiracDeltaFunctions(self.__mypde.getDomain())) self.__r.setTaggedValue(self.__source_tag, source_vector) def _getAcceleration(self, t, u): """ returns the acceleraton for time `t` and solution `u` at time `t` """ dQ = grad(u[0])[0] dP = grad(u[1])[1:] sigma=self.__mypde.getCoefficient('X') sigma[0,0] = self.v2_n*dQ sigma[0,1:] = self.v2_t*dP sigma[1,0] = self.v2_t*dQ sigma[1,1:] = self.v2_p*dP self.__mypde.setValue(X=-sigma, y_dirac= self.__r * self.__wavelet.getValue(t)) return self.__mypde.getSolution()
class HTIWave(WaveBase): """ Solving the HTI wave equation (along the x_0 axis) :note: In case of a two dimensional domain a horizontal domain is considered, i.e. the depth component is dropped. """ def __init__(self, domain, v_p, v_s, wavelet, source_tag, source_vector = [1.,0.,0.], eps=0., gamma=0., delta=0., rho=1., dt=None, u0=None, v0=None, absorption_zone=None, absorption_cut=1e-2, lumping=True, disable_fast_assemblers=False): """ initialize the VTI wave solver :param domain: domain of the problem :type domain: `Domain` :param v_p: vertical p-velocity field :type v_p: `Scalar` :param v_s: vertical s-velocity field :type v_s: `Scalar` :param wavelet: wavelet to describe the time evolution of source term :type wavelet: `Wavelet` :param source_tag: tag of the source location :type source_tag: 'str' or 'int' :param source_vector: source orientation vector :param eps: first Thompsen parameter :param delta: second Thompsen parameter :param gamma: third Thompsen parameter :param rho: density :param dt: time step size. If not present a suitable time step size is calculated. :param u0: initial solution. If not present zero is used. :param v0: initial solution change rate. If not present zero is used. :param absorption_zone: thickness of absorption zone :param absorption_cut: boundary value of absorption decay factor :param lumping: if True mass matrix lumping is being used. This is accelerates the computing but introduces some diffusion. :param disable_fast_assemblers: if True, forces use of slower and more general PDE assemblers """ DIM=domain.getDim() self.fastAssembler = hasattr(domain, "createAssembler") and not disable_fast_assemblers f=createAbsorptionLayerFunction(v_p.getFunctionSpace().getX(), absorption_zone, absorption_cut) v_p=v_p*f v_s=v_s*f if u0 == None: u0=Vector(0.,Solution(domain)) else: u0=interpolate(p0, Solution(domain )) if v0 == None: v0=Vector(0.,Solution(domain)) else: v0=interpolate(v0, Solution(domain )) if dt == None: dt=min((1./5.)*min(inf(domain.getSize()/v_p), inf(domain.getSize()/v_s)), wavelet.getTimeScale()) super(HTIWave, self).__init__( dt, u0=u0, v0=v0, t0=0.) self.__wavelet=wavelet self.c33 = v_p**2 * rho self.c44 = v_s**2 * rho self.c11 = (1+2*eps) * self.c33 self.c66 = (1+2*gamma) * self.c44 self.c13 = sqrt(2*self.c33*(self.c33-self.c44) * delta + (self.c33-self.c44)**2)-self.c44 self.c23 = self.c33-2*self.c66 if self.fastAssembler: C = [("c11", self.c11), ("c23", self.c23), ("c13", self.c13), ("c33", self.c33), ("c44", self.c44), ("c66", self.c66)] if "speckley" in domain.getDescription().lower(): C = [(n, interpolate(d, ReducedFunction(domain))) for n,d in C] self.__mypde=WavePDE(domain, C) else: self.__mypde=LinearPDESystem(domain) self.__mypde.setValue(X=self.__mypde.createCoefficient('X')) if lumping: self.__mypde.getSolverOptions().setSolverMethod(SolverOptions.HRZ_LUMPING) self.__mypde.setSymmetryOn() self.__mypde.setValue(D=rho*kronecker(DIM)) self.__source_tag=source_tag if DIM == 2: source_vector= [source_vector[0],source_vector[2]] self.__r=Vector(0, DiracDeltaFunctions(self.__mypde.getDomain())) self.__r.setTaggedValue(self.__source_tag, source_vector) def setQ(self,q): """ sets the PDE q value :param q: the value to set """ self.__mypde.setValue(q=q) def _getAcceleration(self, t, u): """ returns the acceleraton for time `t` and solution `u` at time `t` """ du = grad(u) if self.fastAssembler: self.__mypde.setValue(du=du, y_dirac= self.__r * self.__wavelet.getValue(t)) else: sigma=self.__mypde.getCoefficient('X') if self.__mypde.getDim() == 3: e11=du[0,0] e22=du[1,1] e33=du[2,2] sigma[0,0]=self.c11*e11+self.c13*(e22+e33) sigma[1,1]=self.c13*e11+self.c33*e22+self.c23*e33 sigma[2,2]=self.c13*e11+self.c23*e22+self.c33*e33 s=self.c44*(du[2,1]+du[1,2]) sigma[1,2]=s sigma[2,1]=s s=self.c66*(du[2,0]+du[0,2]) sigma[0,2]=s sigma[2,0]=s s=self.c66*(du[0,1]+du[1,0]) sigma[0,1]=s sigma[1,0]=s else: e11=du[0,0] e22=du[1,1] sigma[0,0]=self.c11*e11+self.c13*e22 sigma[1,1]=self.c13*e11+self.c33*e22 s=self.c66*(du[1,0]+du[0,1]) sigma[0,1]=s sigma[1,0]=s self.__mypde.setValue(X=-sigma, y_dirac= self.__r * self.__wavelet.getValue(t)) return self.__mypde.getSolution()