def getPotential(self):
        """
        returns a list containing 3 lists one for each the primary, secondary
        and total potential.
        """


        primCon=self.primaryConductivity
        coords=self.domain.getX()
        pde=LinearPDE(self.domain, numEquations=1)
        tol=1e-8
        pde.getSolverOptions().setTolerance(tol)
        pde.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)
        pde.setValue(A=A,q=q)

        delPhiSecondary = []
        delPhiPrimary = []
        delPhiTotal = []
        if(len(self.electrodes[0])==3):

            for i in range(self.numElectrodes-1):
                analyticRs=es.Data(0,(3,),es.ContinuousFunction(self.domain))
                analyticRs[0]=(coords[0]-self.electrodes[i][0])
                analyticRs[1]=(coords[1]-self.electrodes[i][1])
                analyticRs[2]=(coords[2])
                rsMag=(analyticRs[0]**2+analyticRs[1]**2+analyticRs[2]**2)**0.5
                analyticPrimaryPot=(self.current*(1./primCon))/(2*pi*(rsMag+(es.whereZero(rsMag)*0.0000001))) #the magic number 0.0000001 is to avoid devide by 0
                analyticRsPolePower=(analyticRs[0]**2+analyticRs[1]**2+analyticRs[2]**2)**1.5
                analyticRsPolePower = analyticRsPolePower+(es.whereZero(analyticRsPolePower)*0.0000001)
                gradUPrimary = es.Data(0,(3,),es.ContinuousFunction(self.domain))
                gradUPrimary[0] =(self.current/(2*pi*primCon)) * (analyticRs[0]/analyticRsPolePower)
                gradUPrimary[1] =(self.current/(2*pi*primCon)) * (analyticRs[1]/analyticRsPolePower)
                gradUPrimary[2] =(self.current/(2*pi*primCon)) * (analyticRs[2]/analyticRsPolePower)
                gradUPrimary=-gradUPrimary
                X=(primCon-self.secondaryConductivity) * gradUPrimary
                pde.setValue(X=X)
                u=pde.getSolution()
                loc=Locator(self.domain,self.electrodes[i+1])
                delPhiSecondary.append(loc.getValue(u))
                delPhiPrimary.append(loc.getValue(analyticPrimaryPot))
        else:
            raise NotImplementedError("2d forward model is not yet implemented")

        self.delPhiSecondary = delPhiSecondary
        self.delPhiPrimary = delPhiPrimary
        for i in range(len(delPhiPrimary)):
            delPhiTotal.append(delPhiPrimary[i] + delPhiSecondary[i])
        self.delPhiTotal=delPhiTotal
        return [delPhiPrimary, delPhiSecondary, delPhiTotal]
 def test_replaceNaNConstant(self):
     dom=self.domain
     dat = es.Data(10,es.ContinuousFunction(dom))
     dat=(dat*0)/0
     self.assertTrue(dat.hasNaN(),"dat should contain NaN but its doesn't")
     dat.replaceNaN(10)
     self.assertEqual(es.Lsup(dat), 10)
     dat = es.Data(10,es.ContinuousFunction(dom))
     dat.promote()
     dat=(dat*0)/0
     self.assertTrue(dat.hasNaN(),"dat should contain NaN but its doesn't")
     dat.replaceNaN(4+3j)
     self.assertEqual(es.Lsup(dat), 5)
 def test_replaceNaNExpanded(self):
     dom=self.domain
     scl=es.Scalar(0,es.ContinuousFunction(dom))
     scl.expand()
     sclNaN=scl/0
     self.assertTrue(sclNaN.hasNaN(),"sclNaN should contain NaN but its doesn't")
     sclNaN.replaceNaN(15.0)
     self.assertEqual(es.Lsup(sclNaN), 15.0)
     scl=es.Scalar(0,es.ContinuousFunction(dom))
     scl.expand()
     scl.promote()
     if not es.getEscriptParamInt('AUTOLAZY')==1:            
         sclNaN=scl/0
         self.assertTrue(sclNaN.hasNaN(),"sclNaN should contain NaN but its doesn't")
         sclNaN.replaceNaN(3+4j)
         self.assertEqual(es.Lsup(sclNaN), 5.0)
 def test_replaceInfExpanded(self):
     dom=self.domain
     scl=es.Scalar(1,es.ContinuousFunction(dom))
     scl.expand()
     sclNaN=scl/0
     self.assertTrue(sclNaN.hasInf(),"sclNaN should contain Inf but its doesn't")
     sclNaN.replaceNaN(17.0) # Make sure it is distinguishing between NaN and Inf
     sclNaN.replaceInf(15.0)
     self.assertEqual(es.Lsup(sclNaN), 15.0)
     scl=es.Scalar(10,es.ContinuousFunction(dom))
     scl.expand()
     scl.promote()
     if not es.getEscriptParamInt('AUTOLAZY')==1:            
         sclNaN=scl/0
         self.assertTrue(sclNaN.hasInf(),"sclNaN should contain Inf but its doesn't")
         sclNaN.replaceInf(3+4j)
         sclNaN.replaceNaN(3+19j)    # Make sure it is distinguishing between NaN and Inf
         self.assertEqual(es.Lsup(sclNaN), 5.0)
 def checkBounds(self):
     X = es.ContinuousFunction(self.domain).getX()
     xDim=[es.inf(X[0]),es.sup(X[0])]
     yDim=[es.inf(X[1]),es.sup(X[1])]
     zDim=[es.inf(X[2]),es.sup(X[2])]
     for i in range(self.numElectrodes):
         if (self.electrodes[i][0] < xDim[0] or self.electrodes[i][0] > xDim[1]
                 or self.electrodes[i][1] < yDim[0] or self.electrodes[i][1] > yDim[1]):
             print (self.electrodes[i])
             raise ValueError("Electrode setup extents past domain dimentions")
Ejemplo n.º 6
0
    def __init__(self,
                 domain,
                 mode,
                 freq_def,
                 tags,
                 rho,
                 rho_1d,
                 ifc_1d,
                 xstep=100,
                 zstep=100,
                 maps=None,
                 plot=False,
                 limits=None):
        """
    DESCRIPTION:
    ------------
    Constructor which initialises the 2D magnetotelluric class:
    (*) check for argument type
    (*) check for valid argument values
    (*) initialises required values

    ARGUMENTS:
    ----------
    param  domain       :: the 2d mesh domain
    type   domain       :: ``escript data object``

    param  mode         :: TE or TM mode
    type   mode         :: ``string``

    param  freq_def     :: highest/lowest frequency & points per decade
    type   freq_def     :: ``dictionary``

    param  tags         :: the tag names of the regions defined in the mesh
    type   tags         :: ``list``

    param  rho          :: the resistivity values of the regions in the mesh
    type   rho          :: ``list``

    param  rho_1d       :: the resistivity values at the left & right boundary
    type   rho_1d       :: ``dictionary``

    param  ifc_1d       :: the layer interface depths of the left & right boundary
    type   ifc_1d       :: ``dictionary``

    param  xstep        :: user-defined step size for horizontal sample list
    type   xstep        :: ``number``  (optional)

    param  zstep        :: user-defined step size for vertical sample list
    type   zstep        :: ``number``  (optional)

    param  maps         :: list with user-defined  functions which map the resistivity for each region
    type   maps         :: ``list``    (optional)

    param  plot         :: user-defined flag to show a plot of apparent resistivity and phase at each frequency
    type   plot         :: ``boolean`` (optional)



    DATA ATTRIBUTES:
    ----------------
    self.domain         :: escript data object of mesh
    self.X              :: escript data object with all mesh coordinates
    self.mode           :: string with TE or TM mode
    self.xmin           :: float with x-coordinate minimum
    self.xmax           :: float with x-coordinate maximum
    self.zmin           :: float with z-coordinate minimum
    self.zmax           :: float with z-coordinate maximum
    self.zstep          :: number with sample step in vertical direction
    self.xstep          :: number with sample step in horizontal direction
    self.rho            :: list with resistivity values of all regions
    self.rho_1d         :: dictionary with resistivity values at boundaries left/right
    self.ifc_1d         :: dictionary with interface depths at boundaries left/right
    self.plot           :: boolean flag to show plots of apparent resistivity and phase
    self.sigma          :: escript data object with the conductivity model (based on 'rho' and 'maps')
    self.frequencies    :: list of sounding frequencies
    self.boundary_mask  :: Dirichlet mask at boundaries
    """

        if not HAVE_FINLEY:
            raise ImportError("Finley module not available")
        #make python3 compatible, since long disappeared in python 3
        if sys.version_info[0] == 3:
            long_type = int
        else:
            long_type = long
        # ---
        # Checks
        # ---

        # Types:
        if not isinstance(domain, escript.Domain):
            raise ValueError("Input parameter DOMAIN must be an Escript mesh")
        if not isinstance(mode, str):
            raise ValueError("Input parameter MODE must be a string")
        if not isinstance(freq_def, dict) or len(freq_def) != 3:
            raise ValueError(
                "Input parameter FREQ_DEF must be a dictionary with length 3")
        if not isinstance(tags, list) or not all(
                isinstance(t, str) for t in tags):
            raise ValueError("Input parameter TAGS must be a list of strings")
        if not isinstance(rho, list) or not all(
                isinstance(d, (int, long_type, float)) for d in rho):
            raise ValueError("Input parameter RHO must be a list of numbers")
        if not isinstance(rho_1d, dict) or len(rho_1d) != 2:
            raise ValueError(
                "Input parameter RHO_1D must be a dictionary with length 2")
        if not isinstance(ifc_1d, dict) or len(ifc_1d) != 2:
            raise ValueError(
                "Input parameter IFC_1D must be a dictionary with length 2")
        if not isinstance(xstep, (int, long_type, float)):
            raise ValueError("Optional input parameter XSTEP must be a number")
        if not isinstance(zstep, (int, long_type, float)):
            raise ValueError("Optional input parameter ZSTEP must be a number")
        if maps is not None:
            if not isinstance(maps, list) or not all(
                    isinstance(m, (types.FunctionType, types.NoneType))
                    for m in maps):
                raise ValueError(
                    "Optional input parameter MAPS must be a list of Functions or Nones"
                )
        if plot is not None:
            if not isinstance(plot, bool):
                raise ValueError(
                    "Optional input parameter PLOT must be True or False")

        # Values:
        if mode.upper() != "TE" and mode.upper() != "TM":  # Check mode:
            raise ValueError(
                "Input parameter mode must be either 'TE' or 'TM'")
        if not 'high' in freq_def and not 'low' in freq_def and not 'step' in freq_def:
            raise ValueError(
                "Input dictionary FREQ_DEF must have keys 'high', 'low' and 'step' defined"
            )
        if freq_def['high'] < freq_def['low']:
            raise ValueError(
                "High frequency value is smaller than low frequency value in input parameter FREQ_DEF"
            )
        if freq_def['step'] < 1:
            raise ValueError(
                "Step frequency value is smaller than 1 in input parameter FREQ_DEF"
            )
        if not all(r > 0 for r in rho):  # Check resistivity values:
            raise ValueError("Input parameter RHO must be all positive")
        if len(rho) != len(tags):  # Check resistivity list-length:
            raise ValueError(
                "Input parameter RHO must have the same length as input parameter TAGS"
            )
        if not 'left' in rho_1d and not 'right' in rho_1d:
            raise ValueError(
                "Input dictionary RHO_1D must have keys 'left' and 'right' defined"
            )
        if not 'left' in ifc_1d and not 'right' in ifc_1d:
            raise ValueError(
                "Input dictionary IFC_1D must have keys 'left' and 'right' defined"
            )
        if len(ifc_1d['left']) - 1 != len(rho_1d['left']) and len(
                ifc_1d['right']) - 1 != len(rho_1d['right']):
            raise ValueError(
                "Lists with values in input dictionary RHO_1D must have length equal to IFC_1D"
            )
        if xstep < 0.5:  # Step size should be non-zero but should not be smaller than 0.5m:
            raise ValueError("Input parameter XSTEP must be at least 0.5")
        if zstep < 0.5:  # Step size should be non-zero but should not be smaller than 0.5m:
            raise ValueError("Input parameter ZSTEP must be at least 0.5")

        # ---
        # Domain coordinates & tags:
        # ---

        # Extract the model coordinates..
        X = escript.Solution(domain).getX()

        # Get the Min/Max coordinates:
        xmin = escript.inf(X[0])
        xmax = escript.sup(X[0])
        zmin = escript.inf(X[1])
        zmax = escript.sup(X[1])

        # Get the tag names from the mesh file
        mesh_tags = escript.getTagNames(domain)

        if xmin >= xmax or zmin >= zmax:
            raise ValueError("The mesh limits are not valid (min >= max)")
        if zmin >= 0:
            raise ValueError(
                "The mesh must be defined with a negative vertical axis")
        if not set(mesh_tags) == set(tags):
            print("user-tags:", tags)
            print("mesh-tags:", mesh_tags)
            raise ValueError(
                "Input parameter TAGS does not match the tags defined in the mesh"
            )

        # ---
        # Define the boundary mask:
        # ---

        boundary_mask = self.__setBoundaryMask(X)

        # ---
        # Calculate list of sounding frequencies:
        # ---

        frequencies = self.__getSoundingFrequencies(freq_def)

        # ---
        # Tag the domain with conductivity maps:
        # ---

        sigma_model = self.__tagDomain(domain, X, tags, rho, maps)

        # Check for valid values
        if escript.inf(sigma_model) < 0 or escript.sup(sigma_model) < 0:
            raise ValueError("Negative conductivity encountered")
        if cmath.isnan( escript.inf(sigma_model) ) or \
           cmath.isnan( escript.sup(sigma_model) ) or \
           cmath.isinf( escript.sup(sigma_model) ):
            raise ValueError("The conductivity model contains NaNs or INFs")

        # ---
        # Projector and Locator objects.
        # ---

        print("Setting up Escript Locator and Projector objects...")

        # Setup a list with sample points along the vertical mesh extent, bottom to top:
        xsample = self.__getSamplePoints(escript.inf(X[0]),
                                         escript.sup(X[0]),
                                         xstep,
                                         constant=0.0)

        # Get the locations of mesh points at the surface via the Locator object
        # operating on the continuous function-space (i.e. nodes) of the domain.
        loc = pdetools.Locator(escript.ContinuousFunction(domain), xsample)

        # Instantiate the Projector class with smoothing on (fast=False);
        # the Projector is used to calculate the gradient correctly.
        proj = pdetools.Projector(domain, reduce=False, fast=False)

        # ---
        # Print information:
        # ---

        print("")
        print("=" * 72)
        print("Escript MT2D, version", self.__version)
        print("=" * 72)
        print("Mesh XMin/XMax       : ", xmin, xmax)
        print("Mesh ZMin/ZMax       : ", zmin, zmax)
        print("Number of Tags       : ", len(tags))
        print("Mapping              : ", {
            True: 'Yes',
            False: 'No'
        }[maps is not None])
        print("Conductivity Model   : ", sigma_model)
        print("Nr of Frequencies    : ", len(frequencies))
        print("Start/End/Step (Hz)  : ", freq_def["high"], freq_def["low"],
              freq_def["step"])
        print("Mode                 : ", mode.upper())
        print("Solver               : ", MT_2D._solver)
        print("Show plots           : ", plot)
        print("=" * 72)
        print("")

        if self._debug:
            print("Mesh-Info     : ")
            print(domain.print_mesh_info(full=False))

        # ---
        # Set all required variables as data attributes
        # ---

        self.domain = domain
        self.X = X
        self.mode = mode
        self.xmin = xmin
        self.xmax = xmax
        self.zmin = zmin
        self.zmax = zmax
        self.zstep = zstep
        self.xstep = xstep
        self.rho = rho
        self.rho_1d = rho_1d
        self.ifc_1d = ifc_1d
        self.plot = plot
        self.limits = limits
        self.sigma = sigma_model
        self.frequencies = frequencies
        self.boundary_mask = boundary_mask
        self.proj = proj
        self.loc = loc
    def getPotential(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()
        pde=LinearPDE(self.domain, numEquations=1)
        tol=1e-8
        pde.getSolverOptions().setTolerance(tol)
        pde.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)
        pde.setValue(A=A,q=q)

        delPhiSecondaryList = []
        delPhiPrimaryList = []
        delPhiTotalList = []
        for i in range(1,self.n+1): # 1 to n
            maxR = self.numElectrodes - 1 - i #max amount of readings that will fit in the survey
            delPhiSecondary = []
            delPhiPrimary = []
            delPhiTotal = []
            for j in range(maxR):
                analyticRs=es.Data(0,(3,),es.ContinuousFunction(self.domain))
                analyticRs[0]=(coords[0]-self.electrodes[j][0])
                analyticRs[1]=(coords[1]-self.electrodes[j][1])
                analyticRs[2]=(coords[2])
                rsMag=(analyticRs[0]**2+analyticRs[1]**2+analyticRs[2]**2)**0.5
                analyticPrimaryPot=(self.current*(1./primCon))/(2*pi*(rsMag+(es.whereZero(rsMag)*0.0000001))) #the magic number 0.0000001 is to avoid devide by 0

                analyticRsPolePower=(analyticRs[0]**2+analyticRs[1]**2+analyticRs[2]**2)**1.5
                analyticRsPolePower = analyticRsPolePower+(es.whereZero(analyticRsPolePower)*0.0000001)
                gradUPrimary = es.Data(0,(3,),es.ContinuousFunction(self.domain))
                gradUPrimary[0] =(self.current/(2*pi*primCon)) * (analyticRs[0]/analyticRsPolePower)
                gradUPrimary[1] =(self.current/(2*pi*primCon)) * (analyticRs[1]/analyticRsPolePower)
                gradUPrimary[2] =(self.current/(2*pi*primCon)) * (analyticRs[2]/analyticRsPolePower)
                gradUPrimary=-gradUPrimary
                X=(primCon-self.secondaryConductivity) * gradUPrimary
                pde.setValue(X=X)
                u=pde.getSolution()
                loc=Locator(self.domain,[self.electrodes[i+j],self.electrodes[i+j+1]])
                valPrimary=loc.getValue(analyticPrimaryPot)
                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]
    def getPotential(self):
        """
        returns a list containing 3 lists one for each the primary, secondary
        and total potential.
        """

        primCon=self.primaryConductivity
        coords=self.domain.getX()
        pde=LinearPDE(self.domain, numEquations=1)
        tol=1e-8
        pde.getSolverOptions().setTolerance(tol)
        pde.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)
        pde.setValue(A=A,q=q)

        delPhiSecondary = []
        delPhiPrimary = []
        delPhiTotal = []
        if(len(self.electrodes[0])==3):

            for i in range(self.numElectrodes-3):
                analyticRsOne=es.Data(0,(3,),es.ContinuousFunction(self.domain))
                analyticRsOne[0]=(coords[0]-self.electrodes[i][0])
                analyticRsOne[1]=(coords[1]-self.electrodes[i][1])
                analyticRsOne[2]=(coords[2])
                rsMagOne=(analyticRsOne[0]**2+analyticRsOne[1]**2+analyticRsOne[2]**2)**0.5
                analyticRsTwo=es.Data(0,(3,),es.ContinuousFunction(self.domain))
                analyticRsTwo[0]=(coords[0]-self.electrodes[i+3][0])
                analyticRsTwo[1]=(coords[1]-self.electrodes[i+3][1])
                analyticRsTwo[2]=(coords[2])
                rsMagTwo=(analyticRsTwo[0]**2+analyticRsTwo[1]**2+analyticRsTwo[2]**2)**0.5
                rsMagOne+=(es.whereZero(rsMagOne)*0.0000001)
                rsMagTwo+=(es.whereZero(rsMagTwo)*0.0000001)
                analyticPrimaryPot=(self.current/(2*pi*primCon*rsMagOne))-(self.current/(2*pi*primCon*rsMagTwo))

                analyticRsOnePower=(analyticRsOne[0]**2+analyticRsOne[1]**2+analyticRsOne[2]**2)**1.5
                analyticRsOnePower = analyticRsOnePower+(es.whereZero(analyticRsOnePower)*0.0001)
                analyticRsTwoPower=(analyticRsTwo[0]**2+analyticRsTwo[1]**2+analyticRsTwo[2]**2)**1.5
                analyticRsTwoPower = analyticRsTwoPower+(es.whereZero(analyticRsTwoPower)*0.0001)

                gradAnalyticPrimaryPot = es.Data(0,(3,),es.ContinuousFunction(self.domain))
                gradAnalyticPrimaryPot[0] =(self.current/(2*pi*primCon)) \
                        * ((-analyticRsOne[0]/analyticRsOnePower) \
                            + (analyticRsTwo[0]/analyticRsTwoPower))
                gradAnalyticPrimaryPot[1] =(self.current/(2*pi*primCon)) \
                        * ((-analyticRsOne[1]/analyticRsOnePower) \
                            + (analyticRsTwo[1]/analyticRsTwoPower))
                gradAnalyticPrimaryPot[2] =(self.current/(2*pi*primCon)) \
                        * ((-analyticRsOne[2]/analyticRsOnePower)
                            + (analyticRsTwo[2]/analyticRsTwoPower))
                X=(primCon-self.secondaryConductivity) * (gradAnalyticPrimaryPot)
                pde.setValue(X=X)
                u=pde.getSolution()
                loc=Locator(self.domain,[self.electrodes[i+1],self.electrodes[i+2]])
                valPrimary=loc.getValue(analyticPrimaryPot)
                valSecondary=loc.getValue(u)
                delPhiPrimary.append(valPrimary[1]-valPrimary[0])
                delPhiSecondary.append(valSecondary[1]-valSecondary[0])
                delPhiTotal.append(delPhiPrimary[i]+delPhiSecondary[i])
        else:
            raise NotImplementedError("2d forward model is not yet implemented")

        self.delPhiSecondary = delPhiSecondary
        self.delPhiPrimary = delPhiPrimary
        self.delPhiTotal=delPhiTotal
        return [delPhiPrimary, delPhiSecondary, delPhiTotal]
    def getPotentialAnalytic(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.
        """
        coords=self.domain.getX()
        pde=LinearPDE(self.domain, numEquations=1)
        tol=1e-8
        pde.getSolverOptions().setTolerance(tol)
        pde.setSymmetryOn()
        primCon=self.primaryConductivity
        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)
        pde.setValue(A=A,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):
                analyticRsOne=es.Data(0,(3,),es.ContinuousFunction(self.domain))
                analyticRsOne[0]=(coords[0]-self.electrodes[j][0])
                analyticRsOne[1]=(coords[1]-self.electrodes[j][1])
                analyticRsOne[2]=(coords[2])
                rsMagOne=(analyticRsOne[0]**2+analyticRsOne[1]**2+analyticRsOne[2]**2)**0.5
                analyticRsTwo=es.Data(0,(3,),es.ContinuousFunction(self.domain))
                analyticRsTwo[0]=(coords[0]-self.electrodes[j + ((2*i) + 1)][0])
                analyticRsTwo[1]=(coords[1]-self.electrodes[j + ((2*i) + 1)][1])
                analyticRsTwo[2]=(coords[2])
                rsMagTwo=(analyticRsTwo[0]**2+analyticRsTwo[1]**2+analyticRsTwo[2]**2)**0.5
                self.sources.append([self.electrodeTags[j], self.electrodeTags[j + ((2*i) + 1)]])
                rsMagOne+=(es.whereZero(rsMagOne)*0.0000001)
                rsMagTwo+=(es.whereZero(rsMagTwo)*0.0000001)
                
                analyticPrimaryPot=(self.current/(2*pi*primCon*rsMagOne))-(self.current/(2*pi*primCon*rsMagTwo))
                analyticRsOnePower=(analyticRsOne[0]**2+analyticRsOne[1]**2+analyticRsOne[2]**2)**1.5
                analyticRsOnePower = analyticRsOnePower+(es.whereZero(analyticRsOnePower)*0.0001)
                analyticRsTwoPower=(analyticRsTwo[0]**2+analyticRsTwo[1]**2+analyticRsTwo[2]**2)**1.5
                analyticRsTwoPower = analyticRsTwoPower+(es.whereZero(analyticRsTwoPower)*0.0001)

                gradAnalyticPrimaryPot = es.Data(0,(3,),es.ContinuousFunction(self.domain))
                gradAnalyticPrimaryPot[0] =(self.current/(2*pi*primCon)) * ((-analyticRsOne[0]/analyticRsOnePower) + (analyticRsTwo[0]/analyticRsTwoPower))
                gradAnalyticPrimaryPot[1] =(self.current/(2*pi*primCon)) * ((-analyticRsOne[1]/analyticRsOnePower) + (analyticRsTwo[1]/analyticRsTwoPower))
                gradAnalyticPrimaryPot[2] =(self.current/(2*pi*primCon)) * ((-analyticRsOne[2]/analyticRsOnePower) + (analyticRsTwo[2]/analyticRsTwoPower))
                X=(primCon-self.secondaryConductivity) * (gradAnalyticPrimaryPot)
                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(analyticPrimaryPot)
                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]
Ejemplo n.º 10
0
domain = finley.ReadGmsh(mesh_file, 3)
mesh_tags = escript.getTagNames(domain)
directionVector = [1, 0]

# -------------------------------------------------------------------------------------------------
# Define the primary and secondary resistivity model.
# -------------------------------------------------------------------------------------------------

#<Note>: the mesh was created with domains, which were all 'tagged' with unique id-s;
# based on the tagged domain id-s, the conductivity values are assigned to the domains.
#<Note>: the primary conductivity is the conductivity in the vicinity of the electrodes;
# the secondary conductivity is defined as the difference between primary and the
# actual input conductivity.

# Setup primary and secondary conductivity objects:
sig_p = escript.Scalar(0, escript.ContinuousFunction(domain))
sig_s = escript.Scalar(0, escript.ContinuousFunction(domain))

# Cycle tag_values dictionary and assign conductivity values;
# this is done for both, primary and secondary, potentials as
# both dictionaries are setup with identical dictionary-keys:
for tag in tag_p:
    # All initially defined tags must be in the tag list.
    # Print an error if it doesn't and exit the program.
    if tag in mesh_tags:
        # Assign value:
        sig_p.setTaggedValue(tag, tag_p[tag])
        sig_s.setTaggedValue(tag, tag_s[tag])
    else:
        print("Error: the defined tag is not defined in the mesh: " & tag)
        sys.exit()