Esempio n. 1
0
    def __init__(self,
                 domain,
                 v_p,
                 wavelet,
                 source_tag,
                 dt=None,
                 p0=None,
                 p0_t=None,
                 absorption_zone=300 * U.m,
                 absorption_cut=1e-2,
                 lumping=True):
        """
           initialize the sonic wave solver

           :param domain: domain of the problem
           :type domain: `Domain`
           :param v_p: p-velocity field
           :type v_p: `escript.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 dt: time step size. If not present a suitable time step size is calculated.
           :param p0: initial solution. If not present zero is used.
           :param p0_t: 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.
           """
        f = createAbsorptionLayerFunction(
            escript.Function(domain).getX(), absorption_zone, absorption_cut)
        v_p = v_p * f

        if p0 == None:
            p0 = escript.Scalar(0., escript.Solution(domain))
        else:
            p0 = escript.interpolate(p0, escript.Solution(domain))

        if p0_t == None:
            p0_t = escript.Scalar(0., escript.Solution(domain))
        else:
            p0_t = escript.interpolate(p0_t, escript.Solution(domain))

        if dt == None:
            dt = min(escript.inf((1. / 5.) * domain.getSize() / v_p),
                     wavelet.getTimeScale())

        super(SonicWave, self).__init__(dt, u0=p0, v0=p0_t, t0=0.)

        self.__wavelet = wavelet
        self.__mypde = lpde.LinearSinglePDE(domain)
        if lumping:
            self.__mypde.getSolverOptions().setSolverMethod(
                lpde.SolverOptions.HRZ_LUMPING)
        self.__mypde.setSymmetryOn()
        self.__mypde.setValue(D=1. / v_p**2)
        self.__source_tag = source_tag
        self.__r = escript.Scalar(
            0., escript.DiracDeltaFunctions(self.__mypde.getDomain()))
    def getPotentialNumeric(self):
        """
        Returns 3 list each made up of a number of list containing primary, secondary and total
        potentials diferences. Each of the lists contain a list for each value of n.
        """
        primCon=self.primaryConductivity
        coords=self.domain.getX()
        tol=1e-8
        pde=LinearPDE(self.domain, numEquations=1)
        pde.getSolverOptions().setTolerance(tol)
        pde.setSymmetryOn()
        primaryPde=LinearPDE(self.domain, numEquations=1)
        primaryPde.getSolverOptions().setTolerance(tol)
        primaryPde.setSymmetryOn()
        DIM=self.domain.getDim()
        x=self.domain.getX()
        q=es.whereZero(x[DIM-1]-es.inf(x[DIM-1]))
        for i in xrange(DIM-1):
            xi=x[i]
            q+=es.whereZero(xi-es.inf(xi))+es.whereZero(xi-es.sup(xi))
        A = self.secondaryConductivity * es.kronecker(self.domain)
        APrimary = self.primaryConductivity * es.kronecker(self.domain)
        pde.setValue(A=A,q=q)
        primaryPde.setValue(A=APrimary,q=q)

        delPhiSecondaryList = []
        delPhiPrimaryList = []
        delPhiTotalList = []
        for i in range(1,self.n+1): # 1 to n
            maxR = self.numElectrodes - 1 - (2*i) #max amount of readings that will fit in the survey
            delPhiSecondary = []
            delPhiPrimary = []
            delPhiTotal = []
            for j in range(maxR):
                y_dirac=es.Scalar(0,es.DiracDeltaFunctions(self.domain))
                y_dirac.setTaggedValue(self.electrodeTags[j],self.current)
                y_dirac.setTaggedValue(self.electrodeTags[j + ((2*i) + 1)],-self.current)
                self.sources.append([self.electrodeTags[j], self.electrodeTags[j + ((2*i) + 1)]])
                primaryPde.setValue(y_dirac=y_dirac)
                numericPrimaryPot = primaryPde.getSolution()
                X=(primCon-self.secondaryConductivity) * es.grad(numericPrimaryPot)
                pde.setValue(X=X)
                u=pde.getSolution()
                loc=Locator(self.domain,[self.electrodes[j+i],self.electrodes[j+i+1]])
                self.samples.append([self.electrodeTags[j+i],self.electrodeTags[j+i+1]])
                valPrimary=loc.getValue(numericPrimaryPot)
                valSecondary=loc.getValue(u)
                delPhiPrimary.append(valPrimary[1]-valPrimary[0])
                delPhiSecondary.append(valSecondary[1]-valSecondary[0])
                delPhiTotal.append(delPhiPrimary[j]+delPhiSecondary[j])
            delPhiPrimaryList.append(delPhiPrimary)
            delPhiSecondaryList.append(delPhiSecondary)
            delPhiTotalList.append(delPhiTotal)

        self.delPhiPrimaryList=delPhiPrimaryList
        self.delPhiSecondaryList=delPhiSecondaryList
        self.delPhiTotalList = delPhiTotalList
        return [delPhiPrimaryList, delPhiSecondaryList, delPhiTotalList]
Esempio n. 3
0
    def getGradient(self, sigma, u, uTar, uTai, uTu):
        """
        Returns the gradient of the defect with respect to density.

        :param sigma: a suggestion for complex 1/V**2
        :type sigma: ``escript.Data`` of shape (2,)
        :param u: a u vector
        :type u: ``escript.Data`` of shape (2,)
        :param uTar: equals `integrate( w  * (data[0]*u[0]+data[1]*u[1]))`
        :type uTar: `float`
        :param uTai: equals `integrate( w  * (data[1]*u[0]-data[0]*u[1]))`
        :type uTa: `float`
        :param uTu: equals `integrate( w  * (u,u))`
        :type uTu: `float`
        """
        pde = self.setUpPDE()

        if self.scaleF and abs(uTu) > 0:
            Z = ((uTar**2 + uTai**2) / uTu**2) * escript.interpolate(
                u, self.__data.getFunctionSpace())
            Z[0] += (-uTar / uTu) * self.__data[0] + (-uTai /
                                                      uTu) * self.__data[1]
            Z[1] += (-uTar /
                     uTu) * self.__data[1] + uTai / uTu * self.__data[0]

        else:
            Z = u - self.__data
        if Z.getFunctionSpace() == escript.DiracDeltaFunctions(
                self.getDomain()):
            pde.setValue(y_dirac=self.__weight * Z)
        else:
            pde.setValue(y=self.__weight * Z)
        D = pde.getCoefficient('D')
        D[0, 0] = -self.__omega**2 * sigma[0]
        D[0, 1] = -self.__omega**2 * sigma[1]
        D[1, 0] = self.__omega**2 * sigma[1]
        D[1, 1] = -self.__omega**2 * sigma[0]
        pde.setValue(D=D)
        ZTo2 = pde.getSolution() * self.__omega**2
        return escript.inner(
            ZTo2, u) * [1, 0] + (ZTo2[1] * u[0] - ZTo2[0] * u[1]) * [0, 1]
Esempio n. 4
0
    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: `escript.Scalar`
           :param v_s: vertical s-velocity field
           :type v_s: `escript.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 * escript.sqrt(1 + 2 * delta)
        self.v2_n = self.v2_p * (1 + 2 * eps)

        if p0 == None:
            p0 = escript.Data(0., (2, ), escript.Solution(domain))
        else:
            p0 = escript.interpolate(p0, escript.Solution(domain))

        if v0 == None:
            v0 = escript.Data(0., (2, ), escript.Solution(domain))
        else:
            v0 = escript.interpolate(v0, escript.Solution(domain))

        if dt == None:
            dt = min(
                min(escript.inf(domain.getSize() / escript.sqrt(self.v2_p)),
                    escript.inf(domain.getSize() / escript.sqrt(self.v2_t)),
                    escript.inf(domain.getSize() / escript.sqrt(self.v2_n))),
                wavelet.getTimeScale()) * 0.2

        super(SonicHTIWave, self).__init__(dt, u0=p0, v0=v0, t0=0.)

        self.__wavelet = wavelet

        self.__mypde = lpde.LinearPDESystem(domain)
        if lumping:
            self.__mypde.getSolverOptions().setSolverMethod(
                lpde.SolverOptions.HRZ_LUMPING)
        self.__mypde.setSymmetryOn()
        self.__mypde.setValue(D=escript.kronecker(2),
                              X=self.__mypde.createCoefficient('X'))
        self.__source_tag = source_tag

        self.__r = escript.Vector(
            0, escript.DiracDeltaFunctions(self.__mypde.getDomain()))
        self.__r.setTaggedValue(self.__source_tag, source_vector)
Esempio n. 5
0
    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: `escript.Scalar`
           :param v_s: vertical s-velocity field
           :type v_s: `escript.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.
           """
        cos = escript.cos
        sin = escript.sin
        DIM = domain.getDim()
        if not DIM == 2:
            raise ValueError("Only 2D is supported.")
        f = createAbsorptionLayerFunction(
            escript.Function(domain).getX(), absorption_zone, absorption_cut)

        v_p = v_p * f
        v_s = v_s * f

        if u0 == None:
            u0 = escript.Vector(0., escript.Solution(domain))
        else:
            u0 = escript.interpolate(p0, escript.Solution(domain))

        if v0 == None:
            v0 = escript.Vector(0., escript.Solution(domain))
        else:
            v0 = escript.interpolate(v0, escript.Solution(domain))

        if dt == None:
            dt = min((1. / 5.) * min(escript.inf(domain.getSize() / v_p),
                                     escript.inf(domain.getSize() / v_s)),
                     wavelet.getTimeScale())

        super(TTIWave, self).__init__(dt, u0=u0, v0=v0, t0=0.)

        self.__wavelet = wavelet

        self.__mypde = lpde.LinearPDESystem(domain)
        if lumping:
            self.__mypde.getSolverOptions().setSolverMethod(
                lpde.SolverOptions.HRZ_LUMPING)
        self.__mypde.setSymmetryOn()
        self.__mypde.setValue(D=rho * escript.kronecker(DIM),
                              X=self.__mypde.createCoefficient('X'))
        self.__source_tag = source_tag

        self.__r = escript.Vector(
            0, escript.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 = escript.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
Esempio n. 6
0
    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: `escript.Scalar`
       :param v_s: vertical s-velocity field
       :type v_s: `escript.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 = escript.Vector(0., escript.Solution(domain))
        else:
            u0 = escript.interpolate(p0, escript.Solution(domain))

        if v0 == None:
            v0 = escript.Vector(0., escript.Solution(domain))
        else:
            v0 = escript.interpolate(v0, escript.Solution(domain))

        if dt == None:
            dt = min((1. / 5.) * min(escript.inf(domain.getSize() / v_p),
                                     escript.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 = escript.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, escript.interpolate(d,
                                             escript.ReducedFunction(domain)))
                     for n, d in C]
            self.__mypde = lpde.WavePDE(domain, C)
        else:
            self.__mypde = lpde.LinearPDESystem(domain)
            self.__mypde.setValue(X=self.__mypde.createCoefficient('X'))

        if lumping:
            self.__mypde.getSolverOptions().setSolverMethod(
                lpde.SolverOptions.HRZ_LUMPING)
        self.__mypde.setSymmetryOn()
        self.__mypde.setValue(D=rho * escript.kronecker(DIM))
        self.__source_tag = source_tag

        if DIM == 2:
            source_vector = [source_vector[0], source_vector[2]]

        self.__r = escript.Vector(
            0, escript.DiracDeltaFunctions(self.__mypde.getDomain()))
        self.__r.setTaggedValue(self.__source_tag, source_vector)
Esempio n. 7
0
    def __init__(self,
                 domain,
                 omega,
                 w,
                 data,
                 F,
                 coordinates=None,
                 fixAtBottom=False,
                 tol=1e-10,
                 saveMemory=True,
                 scaleF=True):
        """
        initializes a new forward model with acoustic wave form inversion.

        :param domain: domain of the model
        :type domain: `Domain`
        :param w: weighting factors
        :type w: ``Scalar``
        :param data: real and imaginary part of data
        :type data: ``escript.Data`` of shape (2,)
        :param F: real and imaginary part of source given at Dirac points,
                  on surface or at volume.
        :type F: ``escript.Data`` of shape (2,)
        :param coordinates: defines coordinate system to be used (not supported yet)
        :type coordinates: `ReferenceSystem` or `SpatialCoordinateTransformation`
        :param tol: tolerance of underlying PDE
        :type tol: positive ``float``
        :param saveMemory: if true stiffness matrix is deleted after solution
                           of PDE to minimize memory requests. This will
                           require more compute time as the matrix needs to be
                           reallocated.
        :type saveMemory: ``bool``
        :param scaleF: if true source F is scaled to minimize defect.
        :type scaleF: ``bool``
        :param fixAtBottom: if true pressure is fixed to zero at the bottom of
                            the domain
        :type fixAtBottom: ``bool``
        """
        super(AcousticWaveForm, self).__init__()
        self.__trafo = edc.makeTransformation(domain, coordinates)
        if not self.getCoordinateTransformation().isCartesian():
            raise ValueError(
                "Non-Cartesian Coordinates are not supported yet.")
        if not isinstance(data, escript.Data):
            raise ValueError("data must be an escript.Data object.")
        if not data.getFunctionSpace() == escript.FunctionOnBoundary(domain):
            raise ValueError("data must be defined on boundary")
        if not data.getShape() == (2, ):
            raise ValueError(
                "data must have shape (2,) (real and imaginary part).")
        if w is None:
            w = 1.
        if not isinstance(w, escript.Data):
            w = escript.Data(w, escript.FunctionOnBoundary(domain))
        else:
            if not w.getFunctionSpace() == escript.FunctionOnBoundary(domain):
                raise ValueError("Weights must be defined on boundary.")
            if not w.getShape() == ():
                raise ValueError("Weights must be scalar.")

        self.__domain = domain
        self.__omega = omega
        self.__weight = w
        self.__data = data
        self.scaleF = scaleF
        if scaleF:
            A = escript.integrate(self.__weight *
                                  escript.length(self.__data)**2)
            if A > 0:
                self.__data *= 1. / escript.sqrt(A)

        self.__BX = escript.boundingBox(domain)
        self.edge_lengths = np.asarray(escript.boundingBoxEdgeLengths(domain))

        if not isinstance(F, escript.Data):
            F = escript.interpolate(F, escript.DiracDeltaFunctions(domain))
        if not F.getShape() == (2, ):
            raise ValueError(
                "Source must have shape (2,) (real and imaginary part).")

        self.__F = escript.Data()
        self.__f = escript.Data()
        self.__f_dirac = escript.Data()

        if F.getFunctionSpace() == escript.DiracDeltaFunctions(domain):
            self.__f_dirac = F
        elif F.getFunctionSpace() == escript.FunctionOnBoundary(domain):
            self.__f = F
        else:
            self.__F = F
        self.__tol = tol
        self.__fixAtBottom = fixAtBottom
        self.__pde = None
        if not saveMemory:
            self.__pde = self.setUpPDE()