你已经派生过 gprMax
镜像自地址
https://gitee.com/sunhf/gprMax.git
已同步 2025-08-07 15:10:13 +08:00
Re-structuring package layout
这个提交包含在:
700
toolboxes/DebyeFit/Debye_Fit.py
普通文件
700
toolboxes/DebyeFit/Debye_Fit.py
普通文件
@@ -0,0 +1,700 @@
|
||||
# Copyright (C) 2015-2022, Iraklis Giannakis and Sylwia Majchrowska
|
||||
#
|
||||
# This module is licensed under the Creative Commons Attribution-ShareAlike 4.0 International License.
|
||||
# To view a copy of this license, visit http://creativecommons.org/licenses/by-sa/4.0/.
|
||||
#
|
||||
# Please use the attribution at http://dx.doi.org/10.1109/TAP.2014.2308549
|
||||
|
||||
import os
|
||||
import sys
|
||||
import warnings
|
||||
from pathlib import Path
|
||||
|
||||
import matplotlib.gridspec as gridspec
|
||||
import numpy as np
|
||||
import scipy.interpolate
|
||||
from matplotlib import pylab as plt
|
||||
|
||||
from optimization import DA_DLS, DE_DLS, PSO_DLS
|
||||
|
||||
|
||||
class Relaxation(object):
|
||||
""" Create Relaxation function object for complex material.
|
||||
|
||||
:param sigma: The conductivity (Siemens/metre).
|
||||
:type sigma: float, non-optional
|
||||
:param mu: The relative permeability.
|
||||
:type mu: float, non-optional
|
||||
:param mu_sigma: The magnetic loss.
|
||||
:type mu_sigma: float, non-optional
|
||||
:param material_name: A string containing the given name of
|
||||
the material (e.g. "Clay").
|
||||
:type material_name: str, non-optional
|
||||
:param: number_of_debye_poles: Number of Debye functions used to
|
||||
approximate the given electric
|
||||
permittivity.
|
||||
:type number_of_debye_poles: int, optional
|
||||
:param: fn: Number of frequency points in frequency grid.
|
||||
:type fn: int, optional (Default: 50)
|
||||
:param plot: if True will plot the actual and the approximated
|
||||
permittivity at the end (neglected as default: False).
|
||||
:type plot: bool, optional, default:False
|
||||
:param save: if True will save approximated material parameters
|
||||
(not neglected as default: True).
|
||||
:type save: bool, optional, default:True
|
||||
:param optimizer: chosen optimization method:
|
||||
Hybrid Particle Swarm-Damped Least-Squares (PSO_DLS),
|
||||
Dual Annealing (DA) or Differential Evolution (DE)
|
||||
(Default: PSO_DLS).
|
||||
:type optimizer: Optimizer class, optional
|
||||
:param optimizer_options: Additional keyword arguments passed to
|
||||
optimizer class (Default: empty dict).
|
||||
:type optimizer_options: dict, optional, default: empty dict
|
||||
"""
|
||||
|
||||
def __init__(self, sigma, mu, mu_sigma,
|
||||
material_name, f_n=50,
|
||||
number_of_debye_poles=-1,
|
||||
plot=True, save=False,
|
||||
optimizer=PSO_DLS,
|
||||
optimizer_options={}):
|
||||
self.name = 'Relaxation function'
|
||||
self.params = {}
|
||||
self.number_of_debye_poles = number_of_debye_poles
|
||||
self.f_n = f_n
|
||||
self.sigma = sigma
|
||||
self.mu = mu
|
||||
self.mu_sigma = mu_sigma
|
||||
self.material_name = material_name
|
||||
self.plot = plot
|
||||
self.save = save
|
||||
self.optimizer = optimizer(**optimizer_options)
|
||||
|
||||
def set_freq(self, f_min, f_max, f_n=50):
|
||||
""" Interpolate frequency vector using n equally logarithmicaly spaced frequencies.
|
||||
|
||||
Args:
|
||||
f_min (float): First bound of the frequency range
|
||||
used to approximate the given function (Hz).
|
||||
f_max (float): Second bound of the frequency range
|
||||
used to approximate the given function (Hz).
|
||||
f_n (int): Number of frequency points in frequency grid
|
||||
(Default: 50).
|
||||
Note:
|
||||
f_min and f_max must satisfied f_min < f_max
|
||||
"""
|
||||
if abs(f_min - f_max) > 1e12:
|
||||
warnings.warn(f'The chosen range is realy big. '
|
||||
f'Consider setting greater number of points '
|
||||
f'on the frequency grid!')
|
||||
self.freq = np.logspace(np.log10(f_min),
|
||||
np.log10(f_max),
|
||||
int(f_n))
|
||||
|
||||
def check_inputs(self):
|
||||
""" Check the validity of the inputs. """
|
||||
try:
|
||||
d = [float(i) for i in
|
||||
[self.number_of_debye_poles,
|
||||
self.sigma, self.mu, self.mu_sigma]]
|
||||
except ValueError:
|
||||
sys.exit("The inputs should be numeric.")
|
||||
if not isinstance(self.number_of_debye_poles, int):
|
||||
sys.exit("The number of Debye poles must be integer.")
|
||||
if (np.array(d[1:]) < 0).sum() != 0:
|
||||
sys.exit("The inputs should be positive.")
|
||||
|
||||
def calculation(self):
|
||||
""" Approximate the given relaxation function
|
||||
(Havriliak-Negami function, Crim, Jonscher) or based on raw data.
|
||||
"""
|
||||
raise NotImplementedError()
|
||||
|
||||
def print_info(self):
|
||||
"""Readable string of parameters for given approximation settings.
|
||||
|
||||
Returns:
|
||||
s (str): Info about chosen function and its parameters.
|
||||
"""
|
||||
print(f"Approximating {self.name}"
|
||||
f" using {self.number_of_debye_poles} Debye poles")
|
||||
print(f"{self.name} parameters: ")
|
||||
s = ''
|
||||
for k, v in self.params.items():
|
||||
s += f"{k:10s} = {v}\n"
|
||||
print(s)
|
||||
return f'{self.name}:\n{s}'
|
||||
|
||||
def optimize(self):
|
||||
""" Calling the main optimisation module with defined lower and upper boundaries of search.
|
||||
|
||||
Returns:
|
||||
tau (ndarray): The optimised relaxation times.
|
||||
weights (ndarray): Resulting optimised weights for the given relaxation times.
|
||||
ee (float): Average error between the actual and the approximated real part.
|
||||
rl (ndarray): Real parts of chosen relaxation function
|
||||
for given frequency points.
|
||||
im (ndarray): Imaginary parts of chosen relaxation function
|
||||
for given frequency points.
|
||||
"""
|
||||
# Define the lower and upper boundaries of search
|
||||
lb = np.full(self.number_of_debye_poles,
|
||||
-np.log10(np.max(self.freq)) - 3)
|
||||
ub = np.full(self.number_of_debye_poles,
|
||||
-np.log10(np.min(self.freq)) + 3)
|
||||
# Call optimizer to minimize the cost function
|
||||
tau, weights, ee, rl, im = self.optimizer.fit(func=self.optimizer.cost_function,
|
||||
lb=lb, ub=ub,
|
||||
funckwargs={'rl': self.rl,
|
||||
'im': self.im,
|
||||
'freq': self.freq}
|
||||
)
|
||||
return tau, weights, ee, rl, im
|
||||
|
||||
def run(self):
|
||||
""" Solve the problem described by the given relaxation function
|
||||
(Havriliak-Negami function, Crim, Jonscher)
|
||||
or data given from a text file.
|
||||
|
||||
Returns:
|
||||
avg_err (float): average fractional error
|
||||
for relative permittivity (sum)
|
||||
properties (list(str)): Given material nad Debye expnasion parameters
|
||||
in a gprMax format.
|
||||
"""
|
||||
# Check the validity of the inputs
|
||||
self.check_inputs()
|
||||
# Print information about chosen approximation settings
|
||||
self.print_info()
|
||||
# Calculate both real and imaginary parts
|
||||
# for the frequencies included in the vector freq
|
||||
q = self.calculation()
|
||||
# Set the real and the imaginary part of the relaxation function
|
||||
self.rl, self.im = q.real, q.imag
|
||||
|
||||
if self.number_of_debye_poles == -1:
|
||||
print("\n#########",
|
||||
"Try to automaticaly fit number of Debye poles, up to 20!",
|
||||
"##########\n", sep="")
|
||||
error = np.infty # artificial best error starting value
|
||||
self.number_of_debye_poles = 1
|
||||
iteration = 1
|
||||
# stop increasing number of Debye poles if error is smaller then 5%
|
||||
# or 20 debye poles is reached
|
||||
while error > 5 and iteration < 21:
|
||||
# Calling the main optimisation module
|
||||
tau, weights, ee, rl, im = self.optimize()
|
||||
err_real, err_imag = self.error(rl + ee, im)
|
||||
error = err_real + err_imag
|
||||
self.number_of_debye_poles += 1
|
||||
iteration += 1
|
||||
else:
|
||||
# Calling the main optimisation module
|
||||
# for choosen number of debye poles
|
||||
# if one of the weights is negative increase the stabiliser
|
||||
# and repeat the optimisation
|
||||
tau, weights, ee, rl, im = self.optimize()
|
||||
err_real, err_imag = self.error(rl + ee, im)
|
||||
|
||||
# Print the results in gprMax format style
|
||||
properties = self.print_output(tau, weights, ee)
|
||||
print(f'The average fractional error for:\n'
|
||||
f'- real part: {err_real}\n'
|
||||
f'- imaginary part: {err_imag}\n')
|
||||
if self.save:
|
||||
self.save_result(properties)
|
||||
# Plot the actual and the approximate dielectric properties
|
||||
if self.plot:
|
||||
self.plot_result(rl + ee, im)
|
||||
return err_real + err_imag, properties
|
||||
|
||||
def print_output(self, tau, weights, ee):
|
||||
""" Print out the resulting Debye parameters in a gprMax format.
|
||||
|
||||
Args:
|
||||
tau (ndarray): The best known position form optimization module
|
||||
(optimal design).
|
||||
weights (ndarray): Resulting optimised weights for the given relaxation times.
|
||||
ee (float): Average error between the actual and the approximated real part.
|
||||
|
||||
Returns:
|
||||
material_prop (list(str)): Given material nad Debye expnasion parameters
|
||||
in a gprMax format.
|
||||
"""
|
||||
print("Debye expansion parameters: ")
|
||||
print(f" |{'e_inf':^14s}|{'De':^14s}|{'log(tau_0)':^25s}|")
|
||||
print("_" * 65)
|
||||
for i in range(0, len(tau)):
|
||||
print("Debye {0:}|{1:^14.5f}|{2:^14.5f}|{3:^25.5f}|"
|
||||
.format(i + 1, ee/len(tau), weights[i],
|
||||
tau[i]))
|
||||
print("_" * 65)
|
||||
|
||||
# Print the Debye expnasion in a gprMax format
|
||||
material_prop = []
|
||||
material_prop.append("#material: {} {} {} {} {}\n".format(ee, self.sigma,
|
||||
self.mu,
|
||||
self.mu_sigma,
|
||||
self.material_name))
|
||||
print(material_prop[0], end="")
|
||||
dispersion_prop = "#add_dispersion_debye: {}".format(len(tau))
|
||||
for i in range(len(tau)):
|
||||
dispersion_prop += " {} {}".format(weights[i], 10**tau[i])
|
||||
dispersion_prop += " {}".format(self.material_name)
|
||||
print(dispersion_prop)
|
||||
material_prop.append(dispersion_prop + '\n')
|
||||
return material_prop
|
||||
|
||||
def plot_result(self, rl_exp, im_exp):
|
||||
""" Plot the actual and the approximated electric permittivity,
|
||||
along with relative error for real and imaginary parts
|
||||
using a semilogarithm X axes.
|
||||
|
||||
Args:
|
||||
rl_exp (ndarray): Real parts of optimised Debye expansion
|
||||
for given frequency points (plus average error).
|
||||
im_exp (ndarray): Imaginary parts of optimised Debye expansion
|
||||
for given frequency points.
|
||||
"""
|
||||
plt.close("all")
|
||||
fig = plt.figure(figsize=(16, 8), tight_layout=True)
|
||||
gs = gridspec.GridSpec(2, 1)
|
||||
ax = fig.add_subplot(gs[0])
|
||||
ax.grid(b=True, which="major", linewidth=0.2, linestyle="--")
|
||||
ax.semilogx(self.freq * 1e-6, rl_exp, "b-", linewidth=2.0,
|
||||
label="Debye Expansion: Real part")
|
||||
ax.semilogx(self.freq * 1e-6, -im_exp, "k-", linewidth=2.0,
|
||||
label="Debye Expansion: Imaginary part")
|
||||
ax.semilogx(self.freq * 1e-6, self.rl, "r.",
|
||||
linewidth=2.0, label=f"{self.name}: Real part")
|
||||
ax.semilogx(self.freq * 1e-6, -self.im, "g.", linewidth=2.0,
|
||||
label=f"{self.name}: Imaginary part")
|
||||
ax.set_ylim([-1, np.max(np.concatenate([self.rl, -self.im])) + 1])
|
||||
ax.legend()
|
||||
ax.set_xlabel("Frequency (MHz)")
|
||||
ax.set_ylabel("Relative permittivity")
|
||||
|
||||
ax = fig.add_subplot(gs[1])
|
||||
ax.grid(b=True, which="major", linewidth=0.2, linestyle="--")
|
||||
ax.semilogx(self.freq * 1e-6, (rl_exp - self.rl)/(self.rl + 1), "b-", linewidth=2.0,
|
||||
label="Real part")
|
||||
ax.semilogx(self.freq * 1e-6, (-im_exp + self.im)/(self.im + 1), "k-", linewidth=2.0,
|
||||
label="Imaginary part")
|
||||
ax.legend()
|
||||
ax.set_xlabel("Frequency (MHz)")
|
||||
ax.set_ylabel("Relative approximation error")
|
||||
plt.show()
|
||||
|
||||
def error(self, rl_exp, im_exp):
|
||||
""" Calculate the average fractional error separately for
|
||||
relative permittivity (real part) and conductivity (imaginary part)
|
||||
|
||||
Args:
|
||||
rl_exp (ndarray): Real parts of optimised Debye expansion
|
||||
for given frequency points (plus average error).
|
||||
im_exp (ndarray): Imaginary parts of optimised Debye expansion
|
||||
for given frequency points.
|
||||
Returns:
|
||||
avg_err_real (float): average fractional error
|
||||
for relative permittivity (real part)
|
||||
avg_err_imag (float): average fractional error
|
||||
for conductivity (imaginary part)
|
||||
"""
|
||||
avg_err_real = np.sum(np.abs((rl_exp - self.rl)/(self.rl + 1)) * 100)/len(rl_exp)
|
||||
avg_err_imag = np.sum(np.abs((-im_exp + self.im)/(self.im + 1)) * 100)/len(im_exp)
|
||||
return avg_err_real, avg_err_imag
|
||||
|
||||
@staticmethod
|
||||
def save_result(output, fdir="../materials"):
|
||||
""" Save the resulting Debye parameters in a gprMax format.
|
||||
|
||||
Args:
|
||||
output (list(str)): Material and resulting Debye parameters
|
||||
in a gprMax format.
|
||||
fdir (str): Path to saving directory.
|
||||
"""
|
||||
if fdir != "../materials" and os.path.isdir(fdir):
|
||||
file_path = os.path.join(fdir, "my_materials.txt")
|
||||
elif os.path.isdir("../materials"):
|
||||
file_path = os.path.join("../materials",
|
||||
"my_materials.txt")
|
||||
elif os.path.isdir("materials"):
|
||||
file_path = os.path.join("materials",
|
||||
"my_materials.txt")
|
||||
elif os.path.isdir("user_libs/materials"):
|
||||
file_path = os.path.join("user_libs", "materials",
|
||||
"my_materials.txt")
|
||||
else:
|
||||
sys.exit("Cannot save material properties "
|
||||
f"in {os.path.join(fdir, 'my_materials.txt')}!")
|
||||
fileH = open(file_path, "a")
|
||||
fileH.write(f"## {output[0].split(' ')[-1]}")
|
||||
fileH.writelines(output)
|
||||
fileH.write("\n")
|
||||
fileH.close()
|
||||
print(f"Material properties save at: {file_path}")
|
||||
|
||||
|
||||
class HavriliakNegami(Relaxation):
|
||||
""" Approximate a given Havriliak-Negami function
|
||||
Havriliak-Negami function = ε_∞ + Δε / (1 + (2πfjτ)**α)**β,
|
||||
where f is the frequency in Hz.
|
||||
|
||||
:param f_min: First bound of the frequency range
|
||||
used to approximate the given function (Hz).
|
||||
:type f_min: float
|
||||
:param f_max: Second bound of the frequency range
|
||||
used to approximate the given function (Hz).
|
||||
:type f_max: float
|
||||
:param e_inf: The real relative permittivity at infinity frequency
|
||||
:type e_inf: float
|
||||
:param alpha: Real positive float number which varies 0 < alpha < 1.
|
||||
For alpha = 1 and beta !=0 & beta !=1 Havriliak-Negami
|
||||
transforms to Cole-Davidson function.
|
||||
:type alpha: float
|
||||
:param beta: Real positive float number which varies 0 < beta < 1.
|
||||
For beta = 1 and alpha !=0 & alpha !=1 Havriliak-Negami
|
||||
transforms to Cole-Cole function.
|
||||
:type beta: float
|
||||
:param de: The difference of relative permittivity at infinite frequency
|
||||
and the relative permittivity at zero frequency.
|
||||
:type de: float
|
||||
:param tau_0: Real positive float number, tau_0 is the relaxation time.
|
||||
:type tau_0: float
|
||||
"""
|
||||
def __init__(self, f_min, f_max,
|
||||
alpha, beta, e_inf, de, tau_0,
|
||||
sigma, mu, mu_sigma, material_name,
|
||||
number_of_debye_poles=-1, f_n=50,
|
||||
plot=False, save=False,
|
||||
optimizer=PSO_DLS,
|
||||
optimizer_options={}):
|
||||
super(HavriliakNegami, self).__init__(sigma=sigma, mu=mu, mu_sigma=mu_sigma,
|
||||
material_name=material_name, f_n=f_n,
|
||||
number_of_debye_poles=number_of_debye_poles,
|
||||
plot=plot, save=save,
|
||||
optimizer=optimizer,
|
||||
optimizer_options=optimizer_options)
|
||||
self.name = 'Havriliak-Negami function'
|
||||
# Place the lower frequency bound at f_min and the upper frequency bound at f_max
|
||||
if f_min > f_max:
|
||||
self.f_min, self.f_max = f_max, f_min
|
||||
else:
|
||||
self.f_min, self.f_max = f_min, f_max
|
||||
# Choosing n frequencies logarithmicaly equally spaced between the bounds given
|
||||
self.set_freq(self.f_min, self.f_max, self.f_n)
|
||||
self.e_inf, self.alpha, self.beta, self.de, self.tau_0 = e_inf, alpha, beta, de, tau_0
|
||||
self.params = {'f_min': self.f_min, 'f_max': self.f_max,
|
||||
'eps_inf': self.e_inf, 'Delta_eps': self.de, 'tau_0': self.tau_0,
|
||||
'alpha': self.alpha, 'beta': self.beta}
|
||||
|
||||
def check_inputs(self):
|
||||
""" Check the validity of the Havriliak Negami model's inputs. """
|
||||
super(HavriliakNegami, self).check_inputs()
|
||||
try:
|
||||
d = [float(i) for i in self.params.values()]
|
||||
except ValueError:
|
||||
sys.exit("The inputs should be numeric.")
|
||||
if (np.array(d) < 0).sum() != 0:
|
||||
sys.exit("The inputs should be positive.")
|
||||
if self.alpha > 1:
|
||||
sys.exit("Alpha value must range between 0-1 (0 <= alpha <= 1)")
|
||||
if self.beta > 1:
|
||||
sys.exit("Beta value must range between 0-1 (0 <= beta <= 1)")
|
||||
if self.f_min == self.f_max:
|
||||
sys.exit("Null frequency range")
|
||||
|
||||
def calculation(self):
|
||||
"""Calculates the Havriliak-Negami function for
|
||||
the given parameters."""
|
||||
return self.e_inf + self.de / (
|
||||
1 + (1j * 2 * np.pi *
|
||||
self.freq * self.tau_0)**self.alpha
|
||||
)**self.beta
|
||||
|
||||
|
||||
class Jonscher(Relaxation):
|
||||
""" Approximate a given Jonsher function
|
||||
Jonscher function = ε_∞ - ap * (-1j * 2πf / omegap)**n_p,
|
||||
where f is the frequency in Hz
|
||||
|
||||
:param f_min: First bound of the frequency range
|
||||
used to approximate the given function (Hz).
|
||||
:type f_min: float
|
||||
:param f_max: Second bound of the frequency range
|
||||
used to approximate the given function (Hz).
|
||||
:type f_max: float
|
||||
:params e_inf: The real relative permittivity at infinity frequency.
|
||||
:type e_inf: float, non-optional
|
||||
:params a_p: Jonscher parameter. Real positive float number.
|
||||
:type a_p: float, non-optional
|
||||
:params omega_p: Jonscher parameter. Real positive float number.
|
||||
:type omega_p: float, non-optional
|
||||
:params n_p: Jonscher parameter, 0 < n_p < 1.
|
||||
:type n_p: float, non-optional
|
||||
"""
|
||||
def __init__(self, f_min, f_max,
|
||||
e_inf, a_p, omega_p, n_p,
|
||||
sigma, mu, mu_sigma,
|
||||
material_name, number_of_debye_poles=-1,
|
||||
f_n=50, plot=False, save=False,
|
||||
optimizer=PSO_DLS,
|
||||
optimizer_options={}):
|
||||
super(Jonscher, self).__init__(sigma=sigma, mu=mu, mu_sigma=mu_sigma,
|
||||
material_name=material_name, f_n=f_n,
|
||||
number_of_debye_poles=number_of_debye_poles,
|
||||
plot=plot, save=save,
|
||||
optimizer=optimizer,
|
||||
optimizer_options=optimizer_options)
|
||||
self.name = 'Jonsher function'
|
||||
# Place the lower frequency bound at f_min and the upper frequency bound at f_max
|
||||
if f_min > f_max:
|
||||
self.f_min, self.f_max = f_max, f_min
|
||||
else:
|
||||
self.f_min, self.f_max = f_min, f_max
|
||||
# Choosing n frequencies logarithmicaly equally spaced between the bounds given
|
||||
self.set_freq(self.f_min, self.f_max, self.f_n)
|
||||
self.e_inf, self.a_p, self.omega_p, self.n_p = e_inf, a_p, omega_p, n_p
|
||||
self.params = {'f_min': self.f_min, 'f_max': self.f_max,
|
||||
'eps_inf': self.e_inf, 'n_p': self.n_p,
|
||||
'omega_p': self.omega_p, 'a_p': self.a_p}
|
||||
|
||||
def check_inputs(self):
|
||||
""" Check the validity of the inputs. """
|
||||
super(Jonscher, self).check_inputs()
|
||||
try:
|
||||
d = [float(i) for i in self.params.values()]
|
||||
except ValueError:
|
||||
sys.exit("The inputs should be numeric.")
|
||||
if (np.array(d) < 0).sum() != 0:
|
||||
sys.exit("The inputs should be positive.")
|
||||
if self.n_p > 1:
|
||||
sys.exit("n_p value must range between 0-1 (0 <= n_p <= 1)")
|
||||
if self.f_min == self.f_max:
|
||||
sys.exit("Error: Null frequency range!")
|
||||
|
||||
def calculation(self):
|
||||
"""Calculates the Q function for the given parameters"""
|
||||
return self.e_inf + (self.a_p * (2 * np.pi *
|
||||
self.freq / self.omega_p)**(self.n_p-1)) * (
|
||||
1 - 1j / np.tan(self.n_p * np.pi/2))
|
||||
|
||||
|
||||
class Crim(Relaxation):
|
||||
""" Approximate a given CRIM function
|
||||
CRIM = (Σ frac_i * (ε_∞_i + Δε_i/(1 + 2πfj*τ_i))^a)^(1/a)
|
||||
|
||||
:param f_min: First bound of the frequency range
|
||||
used to approximate the given function (Hz).
|
||||
:type f_min: float
|
||||
:param f_max: Second bound of the frequency range
|
||||
used to approximate the given function (Hz).
|
||||
:type f_max: float
|
||||
:param a: Shape factor.
|
||||
:type a: float, non-optional
|
||||
:param: volumetric_fractions: Volumetric fraction for each material.
|
||||
:type volumetric_fractions: ndarray, non-optional
|
||||
:param materials: Arrays of materials properties, for each material [e_inf, de, tau_0].
|
||||
:type materials: ndarray, non-optional
|
||||
"""
|
||||
|
||||
def __init__(self, f_min, f_max, a, volumetric_fractions,
|
||||
materials, sigma, mu, mu_sigma, material_name,
|
||||
number_of_debye_poles=-1, f_n=50,
|
||||
plot=False, save=False,
|
||||
optimizer=PSO_DLS,
|
||||
optimizer_options={}):
|
||||
|
||||
super(Crim, self).__init__(sigma=sigma, mu=mu, mu_sigma=mu_sigma,
|
||||
material_name=material_name, f_n=f_n,
|
||||
number_of_debye_poles=number_of_debye_poles,
|
||||
plot=plot, save=save,
|
||||
optimizer=optimizer,
|
||||
optimizer_options=optimizer_options)
|
||||
self.name = 'CRIM function'
|
||||
# Place the lower frequency bound at f_min and the upper frequency bound at f_max
|
||||
if f_min > f_max:
|
||||
self.f_min, self.f_max = f_max, f_min
|
||||
else:
|
||||
self.f_min, self.f_max = f_min, f_max
|
||||
# Choosing n frequencies logarithmicaly equally spaced between the bounds given
|
||||
self.set_freq(self.f_min, self.f_max, self.f_n)
|
||||
self.a = a
|
||||
self.volumetric_fractions = np.array(volumetric_fractions)
|
||||
self.materials = np.array(materials)
|
||||
self.params = {'f_min': self.f_min, 'f_max': self.f_max,
|
||||
'a': self.a, 'volumetric_fractions': self.volumetric_fractions,
|
||||
'materials': self.materials}
|
||||
|
||||
def check_inputs(self):
|
||||
""" Check the validity of the inputs. """
|
||||
super(Crim, self).check_inputs()
|
||||
try:
|
||||
d = [float(i) for i in
|
||||
[self.f_min, self.f_max, self.a]]
|
||||
except ValueError:
|
||||
sys.exit("The inputs should be numeric.")
|
||||
if (np.array(d) < 0).sum() != 0:
|
||||
sys.exit("The inputs should be positive.")
|
||||
if len(self.volumetric_fractions) != len(self.materials):
|
||||
sys.exit("Number of volumetric volumes does not match the dielectric properties")
|
||||
# Check if the materials are at least two
|
||||
if len(self.volumetric_fractions) < 2:
|
||||
sys.exit("The materials should be at least 2")
|
||||
# Check if the frequency range is null
|
||||
if self.f_min == self.f_max:
|
||||
sys.exit("Null frequency range")
|
||||
# Check if the inputs are positive
|
||||
f = [i for i in self.volumetric_fractions if i < 0]
|
||||
if len(f) != 0:
|
||||
sys.exit("Error: The inputs should be positive")
|
||||
for i in range(len(self.volumetric_fractions)):
|
||||
f = [i for i in self.materials[i][:] if i < 0]
|
||||
if len(f) != 0:
|
||||
sys.exit("Error: The inputs should be positive")
|
||||
# Check if the summation of the volumetric fractions equal to one
|
||||
if np.sum(self.volumetric_fractions) != 1:
|
||||
sys.exit("Error: The summation of volumetric volumes should be equal to 1")
|
||||
|
||||
def print_info(self):
|
||||
""" Print information about chosen approximation settings """
|
||||
print(f"Approximating Complex Refractive Index Model (CRIM)"
|
||||
f" using {self.number_of_debye_poles} Debye poles")
|
||||
print("CRIM parameters: ")
|
||||
for i in range(len(self.volumetric_fractions)):
|
||||
print("Material {}.:".format(i+1))
|
||||
print("---------------------------------")
|
||||
print(f"{'Vol. fraction':>27s} = {self.volumetric_fractions[i]}")
|
||||
print(f"{'e_inf':>27s} = {self.materials[i][0]}")
|
||||
print(f"{'De':>27s} = {self.materials[i][1]}")
|
||||
print(f"{'log(tau_0)':>27s} = {np.log10(self.materials[i][2])}")
|
||||
|
||||
def calculation(self):
|
||||
"""Calculates the Crim function for the given parameters"""
|
||||
return np.sum(np.repeat(self.volumetric_fractions, len(self.freq)
|
||||
).reshape((-1, len(self.materials)))*(
|
||||
self.materials[:, 0] + self.materials[:, 1] / (
|
||||
1 + 1j * 2 * np.pi * np.repeat(self.freq, len(self.materials)
|
||||
).reshape((-1, len(self.materials))) * self.materials[:, 2]))**self.a,
|
||||
axis=1)**(1 / self.a)
|
||||
|
||||
|
||||
class Rawdata(Relaxation):
|
||||
""" Interpolate data given from a text file.
|
||||
|
||||
:param filename: text file which contains three columns:
|
||||
frequency (Hz),Real,Imaginary (separated by comma).
|
||||
:type filename: str, non-optional
|
||||
:param delimiter: separator for three data columns
|
||||
:type delimiter: str, optional (Deafult: ',')
|
||||
"""
|
||||
def __init__(self, filename,
|
||||
sigma, mu, mu_sigma,
|
||||
material_name, number_of_debye_poles=-1,
|
||||
f_n=50, delimiter=',',
|
||||
plot=False, save=False,
|
||||
optimizer=PSO_DLS,
|
||||
optimizer_options={}):
|
||||
|
||||
super(Rawdata, self).__init__(sigma=sigma, mu=mu, mu_sigma=mu_sigma,
|
||||
material_name=material_name, f_n=f_n,
|
||||
number_of_debye_poles=number_of_debye_poles,
|
||||
plot=plot, save=save,
|
||||
optimizer=optimizer,
|
||||
optimizer_options=optimizer_options)
|
||||
self.delimiter = delimiter
|
||||
self.filename = Path(filename).absolute()
|
||||
self.params = {'filename': self.filename}
|
||||
|
||||
def check_inputs(self):
|
||||
""" Check the validity of the inputs. """
|
||||
super(Rawdata, self).check_inputs()
|
||||
|
||||
if not os.path.isfile(self.filename):
|
||||
sys.exit("File doesn't exists!")
|
||||
|
||||
def calculation(self):
|
||||
""" Interpolate real and imaginary part from data.
|
||||
Column framework of the input file three columns comma-separated
|
||||
Frequency(Hz),Real,Imaginary
|
||||
"""
|
||||
# Read the file
|
||||
with open(self.filename) as f:
|
||||
try:
|
||||
array = np.array(
|
||||
[[float(x) for x in line.split(self.delimiter)] for line in f]
|
||||
)
|
||||
except ValueError:
|
||||
sys.exit("Error: The inputs should be numeric")
|
||||
|
||||
self.set_freq(min(array[:, 0]), max(array[:, 0]), self.f_n)
|
||||
rl_interp = scipy.interpolate.interp1d(array[:, 0], array[:, 1],
|
||||
fill_value="extrapolate")
|
||||
im_interp = scipy.interpolate.interp1d(array[:, 0], array[:, 2],
|
||||
fill_value="extrapolate")
|
||||
return rl_interp(self.freq) - 1j * im_interp(self.freq)
|
||||
|
||||
|
||||
if __name__ == "__main__":
|
||||
# Kelley et al. parameters
|
||||
setup = HavriliakNegami(f_min=1e7, f_max=1e11,
|
||||
alpha=0.91, beta=0.45,
|
||||
e_inf=2.7, de=8.6-2.7, tau_0=9.4e-10,
|
||||
sigma=0, mu=0, mu_sigma=0,
|
||||
material_name="Kelley", f_n=100,
|
||||
number_of_debye_poles=6,
|
||||
plot=True, save=False,
|
||||
optimizer_options={'swarmsize': 30,
|
||||
'maxiter': 100,
|
||||
'omega': 0.5,
|
||||
'phip': 1.4,
|
||||
'phig': 1.4,
|
||||
'minstep': 1e-8,
|
||||
'minfun': 1e-8,
|
||||
'seed': 111,
|
||||
'pflag': True})
|
||||
setup.run()
|
||||
setup = HavriliakNegami(f_min=1e7, f_max=1e11,
|
||||
alpha=0.91, beta=0.45,
|
||||
e_inf=2.7, de=8.6-2.7, tau_0=9.4e-10,
|
||||
sigma=0, mu=0, mu_sigma=0,
|
||||
material_name="Kelley", f_n=100,
|
||||
number_of_debye_poles=6,
|
||||
plot=True, save=False,
|
||||
optimizer=DA_DLS,
|
||||
optimizer_options={'seed': 111})
|
||||
setup.run()
|
||||
setup = HavriliakNegami(f_min=1e7, f_max=1e11,
|
||||
alpha=0.91, beta=0.45,
|
||||
e_inf=2.7, de=8.6-2.7, tau_0=9.4e-10,
|
||||
sigma=0, mu=0, mu_sigma=0,
|
||||
material_name="Kelley", f_n=100,
|
||||
number_of_debye_poles=6,
|
||||
plot=True, save=False,
|
||||
optimizer=DE_DLS,
|
||||
optimizer_options={'seed': 111})
|
||||
setup.run()
|
||||
# Testing setup
|
||||
setup = Rawdata("examples/Test.txt", 0.1, 1, 0.1, "M1",
|
||||
number_of_debye_poles=3, plot=True,
|
||||
optimizer_options={'seed': 111})
|
||||
setup.run()
|
||||
np.random.seed(111)
|
||||
setup = HavriliakNegami(1e12, 1e-3, 0.5, 1, 10, 5,
|
||||
1e-6, 0.1, 1, 0, "M2",
|
||||
number_of_debye_poles=6,
|
||||
plot=True)
|
||||
setup.run()
|
||||
setup = Jonscher(1e6, 1e-5, 50, 1, 1e5, 0.7,
|
||||
0.1, 1, 0.1, "M3",
|
||||
number_of_debye_poles=4,
|
||||
plot=True)
|
||||
setup.run()
|
||||
f = np.array([0.5, 0.5])
|
||||
material1 = [3, 25, 1e6]
|
||||
material2 = [3, 0, 1e3]
|
||||
materials = np.array([material1, material2])
|
||||
setup = Crim(1*1e-1, 1e-9, 0.5, f, materials, 0.1,
|
||||
1, 0, "M4", number_of_debye_poles=2,
|
||||
plot=True)
|
||||
setup.run()
|
248
toolboxes/DebyeFit/README.rst
普通文件
248
toolboxes/DebyeFit/README.rst
普通文件
@@ -0,0 +1,248 @@
|
||||
User libraries is a sub-package where useful Python modules contributed by users are stored.
|
||||
|
||||
********
|
||||
DebyeFit
|
||||
********
|
||||
|
||||
Information
|
||||
===========
|
||||
|
||||
**Author/Contact**: Iraklis Giannakis (iraklis.giannakis@abdn.ac.uk), University of Aberdeen, UK and Sylwia Majchrowska (Sylwia.Majchrowska1993@gmail.com)
|
||||
|
||||
This module was created as part of the `Google Summer of Code <https://summerofcode.withgoogle.com/>`_ programme 2021 which gprMax participated.
|
||||
|
||||
**License**: `Creative Commons Attribution-ShareAlike 4.0 International License <http://creativecommons.org/licenses/by-sa/4.0/>`_
|
||||
|
||||
**Attribution/cite**: Giannakis, I., & Giannopoulos, A. (2014). A novel piecewise linear recursive convolution approach for dispersive media using the finite-difference time-domain method. *IEEE Transactions on Antennas and Propagation*, 62(5), 2669-2678. (http://dx.doi.org/10.1109/TAP.2014.2308549)
|
||||
|
||||
Electric permittivity is a complex function with both real and imaginary parts.
|
||||
In general, as a hard and fast rule, the real part dictates the velocity of the medium while the imaginary part is related to the electromagnetic losses.
|
||||
The generic form of dispersive media is
|
||||
|
||||
.. math::
|
||||
|
||||
\epsilon(\omega) = \epsilon^{'}(\omega) - j\epsilon^{''}(\omega),
|
||||
|
||||
where :math:`\omega` is the angular frequency, :math:`\epsilon^{'}` and :math:`\epsilon^{''}` are the real and imaginary parts of the permittivity respectively.
|
||||
|
||||
This package provides scripts and tools which can be used to fit a multi-Debye expansion to dielectric data, defined as
|
||||
|
||||
.. math::
|
||||
|
||||
\epsilon(\omega) = \epsilon_{\infty} + \sum_{i=1}^{N}\frac{\Delta\epsilon_{i}}{1+j\omega t_{0,i}},
|
||||
|
||||
where :math:`\epsilon(\omega)` is frequency dependent dielectric permittivity, :math:`\Delta\epsilon` - difference between the real permittivity at zero and infinite frequency.
|
||||
:math:`\tau_{0}` is relaxation time (seconds), :math:`\epsilon_{\infty}` - real part of relative permittivity at infinite frequency, and :math:`N` is number of the Debye poles.
|
||||
|
||||
To fit the data to a multi-Debye expansion, you can choose between Havriliak-Negami, Jonscher, or Complex Refractive Index Mixing (CRIM) models, as well as arbitrary dielectric data derived experimentally or calculated using a different function.
|
||||
|
||||
.. figure:: ../../images_shared/epsilon.png
|
||||
:width: 600 px
|
||||
|
||||
Real and imaginary parts of frequency-dependent permittivity
|
||||
|
||||
|
||||
Package contents
|
||||
================
|
||||
|
||||
There are two main scripts:
|
||||
|
||||
* ```Debye_fit.py``` contains definitions of all Relaxation functions classes
|
||||
* ```optimization.py``` contains definitions of three choosen global optimization methods
|
||||
|
||||
|
||||
Relaxation Class
|
||||
----------------
|
||||
|
||||
This class is designed for modelling different relaxation functions, like Havriliak-Negami (```HavriliakNegami```), Jonscher (```Jonscher```), Complex Refractive Index Mixing (```CRIM```) models, and arbitrary dielectric data derived experimentally or calculated using some other function (```Rawdata```).
|
||||
|
||||
More about the ``Relaxation`` class structure can be found in the :doc:`Relaxation doc <relaxation.rst>`.
|
||||
|
||||
Havriliak-Negami Function
|
||||
^^^^^^^^^^^^^^^^^^^^^^^^^
|
||||
|
||||
The Havriliak-Negami relaxation is an empirical modification of the Debye relaxation model in electromagnetism, which in addition to the Debye equation has two exponential parameters
|
||||
|
||||
.. math::
|
||||
|
||||
\epsilon(\omega) = \epsilon_{\infty} + \frac{\Delta\epsilon}{\left(1+\left(j\omega t_{0}\right)^{a}\right)^{b}}
|
||||
|
||||
|
||||
The ``HavriliakNegami`` class has the following structure:
|
||||
|
||||
.. code-block:: none
|
||||
|
||||
HavriliakNegami(f_min, f_max,
|
||||
alpha, beta, e_inf, de, tau_0,
|
||||
sigma, mu, mu_sigma, material_name,
|
||||
number_of_debye_poles=-1, f_n=50,
|
||||
plot=False, save=True,
|
||||
optimizer=PSO_DLS,
|
||||
optimizer_options={})
|
||||
|
||||
|
||||
* ``f_min`` is first bound of the frequency range used to approximate the given function (Hz),
|
||||
* ``f_max`` is second bound of the frequency range used to approximate the given function (Hz),
|
||||
* ``alpha`` is real positive float number which varies :math:`0 < \alpha < 1`,
|
||||
* ``beta`` is real positive float number which varies :math:`0 < \beta < 1`,
|
||||
* ``e_inf`` is a real part of relative permittivity at infinite frequency,
|
||||
* ``de`` is a difference between the real permittivity at zero and infinite frequency,
|
||||
* ``tau_0`` is a relaxation time (seconds),
|
||||
* ``sigma`` is a conductivity (Siemens/metre),
|
||||
* ``mu`` is a relative permeability,
|
||||
* ``mu_sigma`` is a magnetic loss (Ohms/metre),
|
||||
* ``material_name`` is the material name,
|
||||
* ``number_of_debye_poles`` is the chosen number of Debye poles,
|
||||
* ``f_n`` is the chosen number of frequences,
|
||||
* ``plot`` is a switch to turn on the plotting,
|
||||
* ``save`` is a switch to turn on saving final material properties,
|
||||
* ``optimizer`` is a chosen optimizer to fit model to dielectric data,
|
||||
* ``optimizer_options`` is a dict for options of chosen optimizer.
|
||||
|
||||
Jonscher Function
|
||||
^^^^^^^^^^^^^^^^^
|
||||
|
||||
Jonscher function is mainly used to describe the dielectric properties of concrete and soils. The frequency domain expression of Jonscher function is given by
|
||||
|
||||
.. math::
|
||||
|
||||
\epsilon(\omega) = \epsilon_{\infty} + a_{p}*\left( -j*\frac{\omega}{\omega_{p}} \right)^{n}
|
||||
|
||||
|
||||
The ``Jonscher`` class has the following structure:
|
||||
|
||||
.. code-block:: none
|
||||
|
||||
Jonscher(f_min, f_max,
|
||||
e_inf, a_p, omega_p, n_p,
|
||||
sigma, mu, mu_sigma,
|
||||
material_name, number_of_debye_poles=-1,
|
||||
f_n=50, plot=False, save=True,
|
||||
optimizer=PSO_DLS,
|
||||
optimizer_options={})
|
||||
|
||||
|
||||
* ``f_min`` is first bound of the frequency range used to approximate the given function (Hz),
|
||||
* ``f_max`` is second bound of the frequency range used to approximate the given function (Hz),
|
||||
* ``e_inf`` is a real part of relative permittivity at infinite frequency,
|
||||
* ``a_p``` is a Jonscher parameter. Real positive float number,
|
||||
* ``omega_p`` is a Jonscher parameter. Real positive float number,
|
||||
* ``n_p`` Jonscher parameter, 0 < n_p < 1.
|
||||
|
||||
Complex Refractive Index Mixing (CRIM) Function
|
||||
^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
|
||||
|
||||
CRIM is the most mainstream approach for estimating the bulk permittivity of heterogeneous materials and has been widely applied for GPR applications. The function takes form of
|
||||
|
||||
.. math::
|
||||
|
||||
\epsilon(\omega)^{d} = \sum_{i=1}^{m}f_{i}\epsilon_{m,i}(\omega)^{d}
|
||||
|
||||
|
||||
The ``CRIM`` class has the following structure:
|
||||
|
||||
.. code-block:: none
|
||||
|
||||
CRIM(f_min, f_max, a, volumetric_fractions,
|
||||
materials, sigma, mu, mu_sigma, material_name,
|
||||
number_of_debye_poles=-1, f_n=50,
|
||||
plot=False, save=True,
|
||||
optimizer=PSO_DLS,
|
||||
optimizer_options={})
|
||||
|
||||
|
||||
* ``f_min`` is first bound of the frequency range used to approximate the given function (Hz),
|
||||
* ``f_max`` is second bound of the frequency range used to approximate the given function (Hz),
|
||||
* ``a`` is a shape factor,
|
||||
* ``volumetric_fractions`` is a volumetric fraction for each material,
|
||||
* ``materials`` are arrays of materials properties, for each material [e_inf, de, tau_0].
|
||||
|
||||
Rawdata Class
|
||||
^^^^^^^^^^^^^
|
||||
|
||||
This package also has the ability to model dielectric properties obtained experimentally by fitting multi-Debye functions to data given from a file.
|
||||
The format of the file should be three columns: the first column contains the frequencies (Hz) associated with the electric permittivity; the second column contains the real part of the relative permittivity; the third column contains the imaginary part of the relative permittivity. The columns should separated by a coma by default, but it is also possible to define a different separator.
|
||||
|
||||
The ``Rawdata`` class has the following structure:
|
||||
|
||||
.. code-block:: none
|
||||
|
||||
Rawdata(self, filename,
|
||||
sigma, mu, mu_sigma,
|
||||
material_name, number_of_debye_poles=-1,
|
||||
f_n=50, delimiter =',',
|
||||
plot=False, save=True,
|
||||
optimizer=PSO_DLS,
|
||||
optimizer_options={})
|
||||
|
||||
|
||||
* ``filename`` is a path to text file which contains three columns,
|
||||
* ``delimiter`` is a separator for three data columns.
|
||||
|
||||
Class Optimizer
|
||||
---------------
|
||||
|
||||
This class supports global optimization algorithms (particle swarm, dual annealing, evolutionary algorithms) for finding an optimal set of relaxation times that minimise the error between the actual and the approximated electric permittivity, and calculates optimised weights for the given relaxation times.
|
||||
Code written here is mainly based on external libraries, like ```scipy``` and ```pyswarm```.
|
||||
|
||||
More about the ``Optimizer`` class structure can be found in the :doc:`Optimisation doc <optimisation.rst>`.
|
||||
|
||||
PSO_DLS Class
|
||||
^^^^^^^^^^^^^
|
||||
|
||||
Creation of hybrid Particle Swarm-Damped Least Squares optimisation object with predefined parameters.
|
||||
The code is a modified version of the pyswarm package which can be found at https://pythonhosted.org/pyswarm/.
|
||||
|
||||
DA_DLS Class
|
||||
^^^^^^^^^^^^
|
||||
|
||||
Creation of Dual Annealing-Damped Least Squares optimisation object with predefined parameters. The class is a modified version of the scipy.optimize package which can be found at:
|
||||
https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.dual_annealing.html#scipy.optimize.dual_annealing.
|
||||
|
||||
DE_DLS Class
|
||||
^^^^^^^^^^^^
|
||||
|
||||
Creation of Differential Evolution-Damped Least Squares object with predefined parameters. The class is a modified version of the scipy.optimize package which can be found at:
|
||||
https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.differential_evolution.html#scipy.optimize.differential_evolution.
|
||||
|
||||
DLS function
|
||||
^^^^^^^^^^^^
|
||||
|
||||
Finding the weights using a non-linear least squares (LS) method, the Levenberg-Marquardt algorithm (LMA or just LM), also known as the damped least-squares (DLS) method.
|
||||
|
||||
How to use the package
|
||||
======================
|
||||
|
||||
Examples
|
||||
--------
|
||||
|
||||
In the examples directory you will find Jupyter notebooks, scripts, and data that demonstrate different cases of how to use the main script ```DebyeFit.py```:
|
||||
|
||||
* ```example_DebyeFitting.ipynb```: simple cases of using all available implemented relaxation functions,
|
||||
* ```example_BiologicalTissues.ipynb```: simple cases of using Cole-Cole function for biological tissues,
|
||||
* ```example_ColeCole.py```: simple cases of using Cole-Cole function in case of 3, 5 and automatically chosen number of Debye poles,
|
||||
* ```Test.txt```: raw data for testing ```Rawdata``` class, file contains 3 columns: the first column contains the frequencies (Hz) associated with the value of the permittivity; the second column contains the real part of the relative permittivity; and the third column contains the imaginary part of the relative permittivity.
|
||||
|
||||
The following code shows a basic example of how to use the Havriliak-Negami function:
|
||||
|
||||
.. code-block:: python
|
||||
|
||||
# set Havrilak-Negami function with initial parameters
|
||||
setup = HavriliakNegami(f_min=1e4, f_max=1e11,
|
||||
alpha=0.3, beta=1,
|
||||
e_inf=3.4, de=2.7, tau_0=.8e-10,
|
||||
sigma=0.45e-3, mu=1, mu_sigma=0,
|
||||
material_name="dry_sand", f_n=100,
|
||||
plot=True, save=False,
|
||||
number_of_debye_poles=3,
|
||||
optimizer_options={'swarmsize':30,
|
||||
'maxiter':100,
|
||||
'omega':0.5,
|
||||
'phip':1.4,
|
||||
'phig':1.4,
|
||||
'minstep':1e-8,
|
||||
'minfun':1e-8,
|
||||
'seed': 111,
|
||||
'pflag': True})
|
||||
# run optimization
|
||||
setup.run()
|
@@ -0,0 +1,41 @@
|
||||
10000000,29.9213536481435,1.25169556541143
|
||||
11220184.5430196,29.9010908161542,1.40299702432563
|
||||
12589254.1179417,29.8756399586587,1.57217537410568
|
||||
14125375.4462276,29.8436916400072,1.76117429474009
|
||||
15848931.9246111,29.8036168365262,1.972079339829
|
||||
17782794.1003892,29.7533952203391,2.20709802135445
|
||||
19952623.1496888,29.6905309907513,2.46852367160776
|
||||
22387211.3856834,29.611956285414,2.75867656083315
|
||||
25118864.3150958,29.5139238169428,3.0798138913592
|
||||
28183829.3126445,29.3918930811742,3.43399830398795
|
||||
31622776.6016838,29.2404187150832,3.82291275991739
|
||||
35481338.9233575,29.0530558254938,4.24760876499865
|
||||
39810717.0553497,28.8223057349563,4.70817602888384
|
||||
44668359.2150963,28.5396365567361,5.20332657814755
|
||||
50118723.3627272,28.1956253280391,5.729897510637
|
||||
56234132.5190349,27.7802793753847,6.28229678071347
|
||||
63095734.4480193,27.283598801889,6.85194778826586
|
||||
70794578.4384137,26.6964309897057,7.42683054876774
|
||||
79432823.4724282,26.0116311645087,7.99126340085894
|
||||
89125093.8133746,25.2254715223637,8.52610387978976
|
||||
100000000,24.339136006498,9.00954486736777
|
||||
112201845.430197,23.3600192468154,9.41861299029905
|
||||
125892541.179417,22.302462791172,9.73132391277098
|
||||
141253754.462276,21.1875720188749,9.92923323827099
|
||||
158489319.246111,20.0419098344834,9.99991217790304
|
||||
177827941.003892,18.8951464050085,9.93877751706085
|
||||
199526231.496888,17.7770654023544,9.74979803763085
|
||||
223872113.856834,16.714546101342,9.44488182455412
|
||||
251188643.150958,15.7291491674702,9.04211441899954
|
||||
281838293.126445,14.8357276987863,8.56331078490773
|
||||
316227766.016838,14.0421664572876,8.03145189099272
|
||||
354813389.233576,13.3500622522389,7.46848900052093
|
||||
398107170.553497,12.7559948799797,6.89379357256362
|
||||
446683592.150963,12.2530185094144,6.32331224790652
|
||||
501187233.627271,11.8320839792422,5.76932819987279
|
||||
562341325.190349,11.4832204965472,5.24065519659195
|
||||
630957344.480194,11.1964125953287,4.74308431385235
|
||||
707945784.384137,10.9621824341848,4.27993617324392
|
||||
794328234.724282,10.7719266477577,3.85261757843109
|
||||
891250938.133744,10.6180694756136,3.46112404221322
|
||||
1000000000,10.4940904606372,3.10446192269295
|
文件差异因一行或多行过长而隐藏
@@ -0,0 +1,44 @@
|
||||
# I. Giannakis, A. Giannopoulos and N. Davidson,
|
||||
# "Incorporating dispersive electrical properties in FDTD GPR models
|
||||
# using a general Cole-Cole dispersion function,"
|
||||
# 2012 14th International Conference on Ground Penetrating Radar (GPR), 2012, pp. 232-236
|
||||
import os, sys
|
||||
sys.path.append(os.path.join(os.path.dirname(__file__), '..'))
|
||||
from Debye_Fit import HavriliakNegami
|
||||
|
||||
if __name__ == "__main__":
|
||||
# set Havrilak-Negami function with initial parameters
|
||||
setup = HavriliakNegami(f_min=1e4, f_max=1e11,
|
||||
alpha=0.3, beta=1,
|
||||
e_inf=3.4, de=2.7, tau_0=.8e-10,
|
||||
sigma=0.45e-3, mu=1, mu_sigma=0,
|
||||
material_name="dry_sand", f_n=100,
|
||||
plot=True, save=False,
|
||||
optimizer_options={'swarmsize':30,
|
||||
'maxiter':100,
|
||||
'omega':0.5,
|
||||
'phip':1.4,
|
||||
'phig':1.4,
|
||||
'minstep':1e-8,
|
||||
'minfun':1e-8,
|
||||
'seed': 111,
|
||||
'pflag': True})
|
||||
### Dry Sand in case of 3, 5
|
||||
# and automatically set number of Debye poles (-1)
|
||||
for number_of_debye_poles in [3, 5, -1]:
|
||||
setup.number_of_debye_poles = number_of_debye_poles
|
||||
setup.run()
|
||||
|
||||
### Moist sand
|
||||
# set Havrilak-Negami function parameters
|
||||
setup.material_name="moist_sand"
|
||||
setup.alpha = 0.25
|
||||
setup.beta = 1
|
||||
setup.e_inf = 5.6
|
||||
setup.de = 3.3
|
||||
setup.tau_0 = 1.1e-10,
|
||||
setup.sigma = 2e-3
|
||||
# calculate for different number of Debye poles
|
||||
for number_of_debye_poles in [3, 5, -1]:
|
||||
setup.number_of_debye_poles = number_of_debye_poles
|
||||
setup.run()
|
文件差异因一行或多行过长而隐藏
462
toolboxes/DebyeFit/optimization.py
普通文件
462
toolboxes/DebyeFit/optimization.py
普通文件
@@ -0,0 +1,462 @@
|
||||
# Copyright (C) 2015-2022, Iraklis Giannakis and Sylwia Majchrowska
|
||||
#
|
||||
# This module is licensed under the Creative Commons Attribution-ShareAlike 4.0 International License.
|
||||
# To view a copy of this license, visit http://creativecommons.org/licenses/by-sa/4.0/.
|
||||
#
|
||||
# Please use the attribution at http://dx.doi.org/10.1109/TAP.2014.2308549
|
||||
|
||||
import numpy as np
|
||||
import scipy.optimize
|
||||
from matplotlib import pylab as plt
|
||||
from tqdm import tqdm
|
||||
|
||||
|
||||
class Optimizer(object):
|
||||
"""
|
||||
Create choosen optimisation object.
|
||||
|
||||
:param maxiter: The maximum number of iterations for the
|
||||
optimizer (Default: 1000).
|
||||
:type maxiter: int, optional
|
||||
:param seed: Seed for RandomState. Must be convertible to 32 bit
|
||||
unsigned integers (Default: None).
|
||||
:type seed: int, NoneType, optional
|
||||
"""
|
||||
def __init__(self, maxiter=1000, seed=None):
|
||||
self.maxiter = maxiter
|
||||
self.seed = seed
|
||||
self.calc_weights = None
|
||||
|
||||
def fit(self, func, lb, ub, funckwargs={}):
|
||||
"""Call the optimization function that tries to find an optimal set
|
||||
of relaxation times that minimise the error
|
||||
between the actual and the approximated electric permittivity
|
||||
and calculate optimised weights for the given relaxation times.
|
||||
|
||||
Args:
|
||||
func (function): The function to be minimized.
|
||||
lb (ndarray): The lower bounds of the design variable(s).
|
||||
ub (ndarray): The upper bounds of the design variable(s).
|
||||
funckwargs (dict): Additional keyword arguments passed to
|
||||
objective and constraint function
|
||||
(Default: empty dict).
|
||||
|
||||
Returns:
|
||||
tau (ndarray): The the best relaxation times.
|
||||
weights (ndarray): Resulting optimised weights for the given
|
||||
relaxation times.
|
||||
ee (float): Average error between the actual and the approximated
|
||||
real part.
|
||||
rl (ndarray): Real parts of chosen relaxation function
|
||||
for given frequency points.
|
||||
im (ndarray): Imaginary parts of chosen relaxation function
|
||||
for given frequency points.
|
||||
"""
|
||||
np.random.seed(self.seed)
|
||||
# find the relaxation frequencies using choosen optimization alghoritm
|
||||
tau, _ = self.calc_relaxation_times(func, lb, ub, funckwargs)
|
||||
# find the weights using a calc_weights method
|
||||
if self.calc_weights is None:
|
||||
raise NotImplementedError()
|
||||
_, _, weights, ee, rl_exp, im_exp = \
|
||||
self.calc_weights(tau, **funckwargs)
|
||||
return tau, weights, ee, rl_exp, im_exp
|
||||
|
||||
def calc_relaxation_times(self):
|
||||
""" Optimisation method that tries to find an optimal set
|
||||
of relaxation times that minimise the error
|
||||
between the actual and the approximated electric permittivity.
|
||||
"""
|
||||
raise NotImplementedError()
|
||||
|
||||
@staticmethod
|
||||
def cost_function(x, rl, im, freq):
|
||||
"""
|
||||
The cost function is the average error between
|
||||
the actual and the approximated electric permittivity.
|
||||
|
||||
Args:
|
||||
x (ndarray): The logarithm with base 10 of relaxation times
|
||||
of the Debyes poles.
|
||||
rl (ndarray): Real parts of chosen relaxation function
|
||||
for given frequency points.
|
||||
im (ndarray): Imaginary parts of chosen relaxation function
|
||||
for given frequency points.
|
||||
freq (ndarray): The frequencies vector for defined grid.
|
||||
|
||||
Returns:
|
||||
cost (float): Sum of mean absolute errors for real and
|
||||
imaginary parts.
|
||||
"""
|
||||
cost_i, cost_r, _, _, _, _ = DLS(x, rl, im, freq)
|
||||
return cost_i + cost_r
|
||||
|
||||
|
||||
class PSO_DLS(Optimizer):
|
||||
""" Create hybrid Particle Swarm-Damped Least Squares optimisation
|
||||
object with predefined parameters.
|
||||
|
||||
:param swarmsize: The number of particles in the swarm (Default: 40).
|
||||
:type swarmsize: int, optional
|
||||
:param maxiter: The maximum number of iterations for the swarm
|
||||
to search (Default: 50).
|
||||
:type maxiter: int, optional
|
||||
:param omega: Particle velocity scaling factor (Default: 0.9).
|
||||
:type omega: float, optional
|
||||
:param phip: Scaling factor to search away from the particle's
|
||||
best known position (Default: 0.9).
|
||||
:type phip: float, optional
|
||||
:param phig: Scaling factor to search away from the swarm's
|
||||
best known position (Default: 0.9).
|
||||
:type phig: float, optional
|
||||
:param minstep: The minimum stepsize of swarm's best position
|
||||
before the search terminates (Default: 1e-8).
|
||||
:type minstep: float, optional
|
||||
:param minfun: The minimum change of swarm's best objective value
|
||||
before the search terminates (Default: 1e-8)
|
||||
:type minfun: float, optional
|
||||
:param pflag: if True will plot the actual and the approximated
|
||||
value during optimization process (Default: False).
|
||||
:type pflag: bool, optional
|
||||
"""
|
||||
def __init__(self, swarmsize=40, maxiter=50,
|
||||
omega=0.9, phip=0.9, phig=0.9,
|
||||
minstep=1e-8, minfun=1e-8,
|
||||
pflag=False, seed=None):
|
||||
|
||||
super(PSO_DLS, self).__init__(maxiter, seed)
|
||||
self.swarmsize = swarmsize
|
||||
self.omega = omega
|
||||
self.phip = phip
|
||||
self.phig = phig
|
||||
self.minstep = minstep
|
||||
self.minfun = minfun
|
||||
self.pflag = pflag
|
||||
self.calc_weights = DLS
|
||||
|
||||
def calc_relaxation_times(self, func, lb, ub, funckwargs={}):
|
||||
"""
|
||||
A particle swarm optimisation that tries to find an optimal set
|
||||
of relaxation times that minimise the error
|
||||
between the actual and the approximated electric permittivity.
|
||||
The current code is a modified edition of the pyswarm package
|
||||
which can be found at https://pythonhosted.org/pyswarm/
|
||||
|
||||
Args:
|
||||
func (function): The function to be minimized.
|
||||
lb (ndarray): The lower bounds of the design variable(s).
|
||||
ub (ndarray): The upper bounds of the design variable(s).
|
||||
funckwargs (dict): Additional keyword arguments passed to
|
||||
objective and constraint function
|
||||
(Default: empty dict).
|
||||
|
||||
Returns:
|
||||
g (ndarray): The swarm's best known position (relaxation times).
|
||||
fg (float): The objective value at ``g``.
|
||||
"""
|
||||
np.random.seed(self.seed)
|
||||
# check input parameters
|
||||
assert len(lb) == len(ub), \
|
||||
'Lower- and upper-bounds must be the same length'
|
||||
assert hasattr(func, '__call__'), 'Invalid function handle'
|
||||
lb = np.array(lb)
|
||||
ub = np.array(ub)
|
||||
assert np.all(ub > lb), \
|
||||
'All upper-bound values must be greater than lower-bound values'
|
||||
|
||||
vhigh = np.abs(ub - lb)
|
||||
vlow = -vhigh
|
||||
|
||||
# Initialize objective function
|
||||
obj = lambda x: func(x=x, **funckwargs)
|
||||
|
||||
# Initialize the particle swarm
|
||||
d = len(lb) # the number of dimensions each particle has
|
||||
x = np.random.rand(self.swarmsize, d) # particle positions
|
||||
v = np.zeros_like(x) # particle velocities
|
||||
p = np.zeros_like(x) # best particle positions
|
||||
fp = np.zeros(self.swarmsize) # best particle function values
|
||||
g = [] # best swarm position
|
||||
fg = np.inf # artificial best swarm position starting value
|
||||
|
||||
for i in range(self.swarmsize):
|
||||
# Initialize the particle's position
|
||||
x[i, :] = lb + x[i, :] * (ub - lb)
|
||||
# Initialize the particle's best known position
|
||||
p[i, :] = x[i, :]
|
||||
# Calculate the objective's value at the current particle's
|
||||
fp[i] = obj(p[i, :])
|
||||
# At the start, there may not be any feasible starting point,
|
||||
# so just
|
||||
# give it a temporary "best" point since it's likely to change
|
||||
if i == 0:
|
||||
g = p[0, :].copy()
|
||||
# If the current particle's position is better than the swarm's,
|
||||
# update the best swarm position
|
||||
if fp[i] < fg:
|
||||
fg = fp[i]
|
||||
g = p[i, :].copy()
|
||||
# Initialize the particle's velocity
|
||||
v[i, :] = vlow + np.random.rand(d) * (vhigh - vlow)
|
||||
|
||||
# Iterate until termination criterion met
|
||||
for it in tqdm(range(self.maxiter), desc='Debye fitting'):
|
||||
rp = np.random.uniform(size=(self.swarmsize, d))
|
||||
rg = np.random.uniform(size=(self.swarmsize, d))
|
||||
for i in range(self.swarmsize):
|
||||
# Update the particle's velocity
|
||||
v[i, :] = self.omega * v[i, :] + self.phip * rp[i, :] * \
|
||||
(p[i, :] - x[i, :]) + \
|
||||
self.phig * rg[i, :] * (g - x[i, :])
|
||||
# Update the particle's position,
|
||||
# correcting lower and upper bound
|
||||
# violations, then update the objective function value
|
||||
x[i, :] = x[i, :] + v[i, :]
|
||||
mark1 = x[i, :] < lb
|
||||
mark2 = x[i, :] > ub
|
||||
x[i, mark1] = lb[mark1]
|
||||
x[i, mark2] = ub[mark2]
|
||||
fx = obj(x[i, :])
|
||||
# Compare particle's best position
|
||||
if fx < fp[i]:
|
||||
p[i, :] = x[i, :].copy()
|
||||
fp[i] = fx
|
||||
# Compare swarm's best position to current
|
||||
# particle's position
|
||||
if fx < fg:
|
||||
tmp = x[i, :].copy()
|
||||
stepsize = np.sqrt(np.sum((g - tmp) ** 2))
|
||||
if np.abs(fg - fx) <= self.minfun:
|
||||
print(f'Stopping search: Swarm best objective '
|
||||
f'change less than {self.minfun}')
|
||||
return tmp, fx
|
||||
elif stepsize <= self.minstep:
|
||||
print(f'Stopping search: Swarm best position '
|
||||
f'change less than {self.minstep}')
|
||||
return tmp, fx
|
||||
else:
|
||||
g = tmp.copy()
|
||||
fg = fx
|
||||
|
||||
# Dynamically plot the error as the optimisation takes place
|
||||
if self.pflag:
|
||||
if it == 0:
|
||||
xpp = [it]
|
||||
ypp = [fg]
|
||||
else:
|
||||
xpp.append(it)
|
||||
ypp.append(fg)
|
||||
PSO_DLS.plot(xpp, ypp)
|
||||
return g, fg
|
||||
|
||||
@staticmethod
|
||||
def plot(x, y):
|
||||
"""
|
||||
Dynamically plot the error as the optimisation takes place.
|
||||
|
||||
Args:
|
||||
x (array): The number of current iterations.
|
||||
y (array): The objective value at for all x points.
|
||||
"""
|
||||
# it clears an axes
|
||||
plt.cla()
|
||||
plt.plot(x, y, "b-", linewidth=1.0)
|
||||
plt.ylim(min(y) - 0.1 * min(y),
|
||||
max(y) + 0.1 * max(y))
|
||||
plt.xlim(min(x) - 0.1, max(x) + 0.1)
|
||||
plt.grid(b=True, which="major", color="k",
|
||||
linewidth=0.2, linestyle="--")
|
||||
plt.suptitle("Debye fitting process")
|
||||
plt.xlabel("Iteration")
|
||||
plt.ylabel("Average Error")
|
||||
plt.pause(0.0001)
|
||||
|
||||
|
||||
class DA_DLS(Optimizer):
|
||||
""" Create Dual Annealing object with predefined parameters.
|
||||
The current class is a modified edition of the scipy.optimize
|
||||
package which can be found at:
|
||||
https://docs.scipy.org/doc/scipy/reference/generated/
|
||||
scipy.optimize.dual_annealing.html#scipy.optimize.dual_annealing
|
||||
"""
|
||||
def __init__(self, maxiter=1000,
|
||||
local_search_options={}, initial_temp=5230.0,
|
||||
restart_temp_ratio=2e-05, visit=2.62, accept=-5.0,
|
||||
maxfun=1e7, no_local_search=False,
|
||||
callback=None, x0=None, seed=None):
|
||||
super(DA_DLS, self).__init__(maxiter, seed)
|
||||
self.local_search_options = local_search_options
|
||||
self.initial_temp = initial_temp
|
||||
self.restart_temp_ratio = restart_temp_ratio
|
||||
self.visit = visit
|
||||
self.accept = accept
|
||||
self.maxfun = maxfun
|
||||
self.no_local_search = no_local_search
|
||||
self.callback = callback
|
||||
self.x0 = x0
|
||||
self.calc_weights = DLS
|
||||
|
||||
def calc_relaxation_times(self, func, lb, ub, funckwargs={}):
|
||||
"""
|
||||
Find the global minimum of a function using Dual Annealing.
|
||||
The current class is a modified edition of the scipy.optimize
|
||||
package which can be found at:
|
||||
https://docs.scipy.org/doc/scipy/reference/generated/
|
||||
scipy.optimize.dual_annealing.html#scipy.optimize.dual_annealing
|
||||
|
||||
Args:
|
||||
func (function): The function to be minimized
|
||||
lb (ndarray): The lower bounds of the design variable(s)
|
||||
ub (ndarray): The upper bounds of the design variable(s)
|
||||
funckwargs (dict): Additional keyword arguments passed to
|
||||
objective and constraint function
|
||||
(Default: empty dict)
|
||||
|
||||
Returns:
|
||||
x (ndarray): The solution array (relaxation times).
|
||||
fun (float): The objective value at the best solution.
|
||||
"""
|
||||
np.random.seed(self.seed)
|
||||
result = scipy.optimize.dual_annealing(
|
||||
func,
|
||||
bounds=list(zip(lb, ub)),
|
||||
args=funckwargs.values(),
|
||||
maxiter=self.maxiter,
|
||||
local_search_options=self.local_search_options,
|
||||
initial_temp=self.initial_temp,
|
||||
restart_temp_ratio=self.restart_temp_ratio,
|
||||
visit=self.visit,
|
||||
accept=self.accept,
|
||||
maxfun=self.maxfun,
|
||||
no_local_search=self.no_local_search,
|
||||
callback=self.callback,
|
||||
x0=self.x0)
|
||||
print(result.message)
|
||||
return result.x, result.fun
|
||||
|
||||
|
||||
class DE_DLS(Optimizer):
|
||||
"""
|
||||
Create Differential Evolution-Damped Least Squares object
|
||||
with predefined parameters.
|
||||
The current class is a modified edition of the scipy.optimize
|
||||
package which can be found at:
|
||||
https://docs.scipy.org/doc/scipy/reference/generated/
|
||||
scipy.optimize.differential_evolution.html#scipy.optimize.differential_evolution
|
||||
"""
|
||||
def __init__(self, maxiter=1000,
|
||||
strategy='best1bin', popsize=15, tol=0.01, mutation=(0.5, 1),
|
||||
recombination=0.7, callback=None, disp=False, polish=True,
|
||||
init='latinhypercube', atol=0,
|
||||
updating='immediate', workers=1,
|
||||
constraints=(), seed=None):
|
||||
super(DE_DLS, self).__init__(maxiter, seed)
|
||||
self.strategy = strategy
|
||||
self.popsize = popsize
|
||||
self.tol = tol
|
||||
self.mutation = mutation
|
||||
self.recombination = recombination
|
||||
self.callback = callback
|
||||
self.disp = disp
|
||||
self.polish = polish
|
||||
self.init = init
|
||||
self.atol = atol
|
||||
self.updating = updating
|
||||
self.workers = workers
|
||||
self.constraints = constraints
|
||||
self.calc_weights = DLS
|
||||
|
||||
def calc_relaxation_times(self, func, lb, ub, funckwargs={}):
|
||||
"""
|
||||
Find the global minimum of a function using Differential Evolution.
|
||||
The current class is a modified edition of the scipy.optimize
|
||||
package which can be found at:
|
||||
https://docs.scipy.org/doc/scipy/reference/generated/
|
||||
scipy.optimize.differential_evolution.html#scipy.optimize.differential_evolution
|
||||
|
||||
Args:
|
||||
func (function): The function to be minimized
|
||||
lb (ndarray): The lower bounds of the design variable(s)
|
||||
ub (ndarray): The upper bounds of the design variable(s)
|
||||
funckwargs (dict): Additional keyword arguments passed to
|
||||
objective and constraint function
|
||||
(Default: empty dict)
|
||||
|
||||
Returns:
|
||||
x (ndarray): The solution array (relaxation times).
|
||||
fun (float): The objective value at the best solution.
|
||||
"""
|
||||
np.random.seed(self.seed)
|
||||
result = scipy.optimize.differential_evolution(
|
||||
func,
|
||||
bounds=list(zip(lb, ub)),
|
||||
args=funckwargs.values(),
|
||||
strategy=self.strategy,
|
||||
popsize=self.popsize,
|
||||
tol=self.tol,
|
||||
mutation=self.mutation,
|
||||
recombination=self.recombination,
|
||||
callback=self.callback,
|
||||
disp=self.disp,
|
||||
polish=self.polish,
|
||||
init=self.init,
|
||||
atol=self.atol,
|
||||
updating=self.updating,
|
||||
workers=self.workers,
|
||||
constraints=self.constraints)
|
||||
print(result.message)
|
||||
return result.x, result.fun
|
||||
|
||||
|
||||
def DLS(logt, rl, im, freq):
|
||||
"""
|
||||
Find the weights using a non-linear least squares (LS) method,
|
||||
the Levenberg–Marquardt algorithm (LMA or just LM),
|
||||
also known as the damped least-squares (DLS) method.
|
||||
|
||||
Args:
|
||||
logt (ndarray): The best known position form optimization module
|
||||
(optimal design),
|
||||
the logarithm with base 10 of relaxation times
|
||||
of the Debyes poles.
|
||||
rl (ndarray): Real parts of chosen relaxation function
|
||||
for given frequency points.
|
||||
im (ndarray): Imaginary parts of chosen relaxation function
|
||||
for given frequency points.
|
||||
freq (ndarray): The frequencies vector for defined grid.
|
||||
|
||||
Returns:
|
||||
cost_i (float): Mean absolute error between the actual and
|
||||
the approximated imaginary part.
|
||||
cost_r (float): Mean absolute error between the actual and
|
||||
the approximated real part (plus average error).
|
||||
x (ndarray): Resulting optimised weights for the given
|
||||
relaxation times.
|
||||
ee (float): Average error between the actual and the approximated
|
||||
real part.
|
||||
rp (ndarray): The real part of the permittivity for the optimised
|
||||
relaxation times and weights for the frequnecies included
|
||||
in freq.
|
||||
ip (ndarray): The imaginary part of the permittivity for the optimised
|
||||
relaxation times and weights for the frequnecies included
|
||||
in freq.
|
||||
"""
|
||||
# The relaxation time of the Debyes are given at as logarithms
|
||||
# logt=log10(t0) for efficiency during the optimisation
|
||||
# Here they are transformed back t0=10**logt
|
||||
tt = 10**logt
|
||||
# y = Ax, here the A matrix for the real and the imaginary part is builded
|
||||
d = 1 / (1 + 1j * 2 * np.pi * np.repeat(
|
||||
freq, len(tt)).reshape((-1, len(tt))) * tt)
|
||||
# Adding dumping (Levenberg–Marquardt algorithm)
|
||||
# Solving the overdetermined system y=Ax
|
||||
x = np.abs(np.linalg.lstsq(d.imag, im, rcond=None)[0])
|
||||
# x - absolute damped least-squares solution
|
||||
rp, ip = np.matmul(d.real, x[np.newaxis].T).T[0], np.matmul(
|
||||
d.imag, x[np.newaxis].T).T[0]
|
||||
cost_i = np.sum(np.abs(ip-im))/len(im)
|
||||
ee = np.mean(rl - rp)
|
||||
if ee < 1:
|
||||
ee = 1
|
||||
cost_r = np.sum(np.abs(rp + ee - rl))/len(im)
|
||||
return cost_i, cost_r, x, ee, rp, ip
|
@@ -0,0 +1,50 @@
|
||||
Optimization methods of multi-Debye fitting
|
||||
-------------------------------------------
|
||||
|
||||
The ``Optimizer`` class supports global optimization algorithms (particle swarm, dual annealing, evolutionary algorithms) for finding an optimal set of relaxation times that minimise the error between the actual and the approximated electric permittivity, and calculates optimised weights for the given relaxation times.
|
||||
Code written here is mainly based on external libraries, like ```scipy``` and ```pyswarm```.
|
||||
|
||||
Supported methods:
|
||||
- [x] hybrid Particle Swarm-Damped Least Squares
|
||||
- [x] hybrid Dual Annealing-Damped Least Squares
|
||||
- [x] hybrid Differential Evolution-Damped Least Squares
|
||||
|
||||
Methods
|
||||
^^^^^^^
|
||||
|
||||
1. __constructor__ - is called in all child classes.
|
||||
|
||||
It takes the following arguments:
|
||||
- `maxiter`: maximum number of iterations for the optimizer,
|
||||
- `seed`: Seed for RandomState.
|
||||
|
||||
In constructor the attributes:
|
||||
- `maxiter`,
|
||||
- `seed`,
|
||||
- `calc_weights` (used to fit weight, non-linear least squares (LS) method is used as a default)
|
||||
are set.
|
||||
|
||||
2. __fit__ - is inherited by all children classes. It calls the optimization function that tries to find an optimal set of relaxation times that minimise the error between the actual and the approximated electric permittivity and calculate optimised weights for the given relaxation times.
|
||||
It takes the following arguments:
|
||||
- `func`: objective function to be optimized,
|
||||
- `lb`: the lower bounds of the design variable(s),
|
||||
- `ub`: the upper bounds of the design variable(s),
|
||||
- `funckwargs`: optional arguments takien by objective function.
|
||||
|
||||
3. __cost_function__ - is inherited by all child classes. It calculates the cost function as the average error between the actual and the approximated electric permittivity (sum of real and imaginary part).
|
||||
It takes the following arguments:
|
||||
- `x`: the logarithm with base 10 of relaxation times of the Debyes poles,
|
||||
- `rl`: real parts of chosen relaxation function for given frequency points,
|
||||
- `im`: imaginary parts of chosen relaxation function for given frequency points,
|
||||
- `freq`: the frequencies vector for defined grid.
|
||||
|
||||
4. __calc_relaxation_times__ - it finds an optimal set of relaxation times that minimise an objective function using appropriate optimization procedure.
|
||||
It takes the following arguments:
|
||||
- `func`: objective function to be optimized,
|
||||
- `lb`: the lower bounds of the design variable(s),
|
||||
- `ub`: the upper bounds of the design variable(s),
|
||||
- `funckwargs`: optional arguments takien by objective function.
|
||||
|
||||
Each new class of optimizer should:
|
||||
- define constructor with appropriate arguments,
|
||||
- overload __calc_relaxation_times__ method (and optional define __calc_weights__ function in case of hybrid optimization procedure).
|
@@ -0,0 +1,71 @@
|
||||
Relaxation classes for multi-Debye fitting
|
||||
------------------------------------------
|
||||
|
||||
This class is designed for modelling different relaxation functions, like Havriliak-Negami (```HavriliakNegami```), Jonscher (```Jonscher```), Complex Refractive Index Mixing (```CRIM```) models, and arbitrary dielectric data derived experimentally or calculated using some other function (```Rawdata```).
|
||||
|
||||
Supported relaxation classes:
|
||||
|
||||
- [x] Havriliak-Negami,
|
||||
- [x] Jonscher,
|
||||
- [x] Complex Refractive Index Mixing,
|
||||
- [x] Experimental data,
|
||||
|
||||
Methods
|
||||
^^^^^^^
|
||||
|
||||
1. __constructor__ - is called in all child classes, creates Relaxation function object for complex material.
|
||||
|
||||
It takes the following arguments:
|
||||
- ``sigma`` is a conductivity (Siemens/metre),
|
||||
- ``mu`` is a relative permeability,
|
||||
- ``mu_sigma`` is a magnetic loss (Ohms/metre),
|
||||
- ``material_name`` is definition of material name,
|
||||
- ``number_of_debye_poles`` is choosen number of Debye poles,
|
||||
- ``f_n`` is the chosen number of frequences,
|
||||
- ``plot`` is a switch to turn on the plotting,
|
||||
- ``save`` is a switch to turn on the saving final material properties,
|
||||
- ``optimizer`` is a chosen optimizer to fit model to dielectric data,
|
||||
- ``optimizer_options`` is a dict for options of chosen optimizer.
|
||||
|
||||
Additional parameters:
|
||||
- ``rl`` calculated real part of chosen relaxation function for given frequency points,
|
||||
- ``im`` calculated imaginary part of chosen relaxation function for given frequency points.
|
||||
|
||||
2. __set_freq__ - is inherited by all child classes, interpolates frequency vector using n equally logarithmicaly spaced frequencies.
|
||||
It takes the following arguments:
|
||||
- `f_min_`: first bound of the frequency range used to approximate the given function (Hz),
|
||||
- `f_max`: second bound of the frequency range used to approximate the given function (Hz),
|
||||
- `f_n`: the number of frequency points in frequency grid (Default: 50).
|
||||
|
||||
3. __run__ - is inherited by all child classes, solves the problem described by the given relaxation function (main operational method).
|
||||
It consists of following steps:
|
||||
1) Check the validity of the inputs using ```check_inputs``` method.
|
||||
2) Print information about chosen approximation settings using ```print_info``` method.
|
||||
3) Calculate both real and imaginary parts using ```calculation``` method, and then set ```self.rl``` and ```self.im``` properties.
|
||||
4) Calling the main optimisation module using ```optimize``` method and calculate error based on ```error``` method.
|
||||
a) [OPTIONAL] If number of debye poles is set to -1 optimization procedure is repeated until the percentage error is les than 5% or 20 Debye poles is reached.
|
||||
5) Print the results in gprMax format style using ```print_output``` method.
|
||||
6) [OPTIONAL] Save results in gprMax style using ```save_result``` method.
|
||||
7) [OPTIONAL] Plot the actual and the approximate dielectric properties using ```plot_result``` method.
|
||||
|
||||
4. __check_inputs__ - is called in all child classes, finds an optimal set of relaxation times that minimise an objective function using appropriate optimization procedure.
|
||||
|
||||
5. __calculation__ - is inherited by all child classes, should be definied in all new chil classes, approximates the given relaxation function.
|
||||
|
||||
6. __print_info__ - is inherited by all child classes, prints readable string of parameters for given approximation settings.
|
||||
|
||||
7. __optimize__ - is inherited by all child classes, calls the main optimisation module with defined lower and upper boundaries of search.
|
||||
|
||||
8. __print_output__ - is inherited by all child classes, prints out the resulting Debye parameters in a gprMax format.
|
||||
|
||||
9. __plot_result__ - is inherited by all child classes, plots the actual and the approximated electric permittivity, along with relative error for real and imaginary parts using a semilogarithm X axes.
|
||||
|
||||
10. __save_result__ - is inherited by all child classes, saves the resulting Debye parameters in a gprMax format.
|
||||
|
||||
11. __error__ -is inherited by all child classes, calculates the average fractional error separately for relative permittivity (real part) and conductivity (imaginary part).
|
||||
|
||||
Each new class of relaxation object should:
|
||||
|
||||
- define constructor with appropriate arguments,
|
||||
- define __check_inputs__ method to check relaxation class specific parameters,
|
||||
- overload __calculation__ method.
|
在新工单中引用
屏蔽一个用户