geodezyx.marine package

Submodules

geodezyx.marine.dac module

Created on Mon Jun 3 15:45:36 2024

@author: psakic

geodezyx.marine.dac.dac_extract(dac_files_lis, mlon=45.46, mlat=-12.78, epoch_min=datetime.datetime(1980, 1, 1, 0, 0), epoch_max=datetime.datetime(2099, 1, 1, 0, 0))

Extracts DAC (Doppler-derived Atmospheric Correction) values from a list of DAC files for a given location and date range.

Parameters:
  • dac_files_lis (list) – The list of DAC files from which the DAC values are to be extracted.

  • mlon (float, optional) – The longitude of the location for which the DAC values are to be extracted. Default is 45.46.

  • mlat (float, optional) – The latitude of the location for which the DAC values are to be extracted. Default is -12.78.

  • epoch_min (datetime, optional) – The start date of the range for which the DAC values are to be extracted. Default is January 1, 1980.

  • epoch_max (datetime, optional) – The end date of the range for which the DAC values are to be extracted. Default is January 1, 2099.

Returns:

A DataFrame with the extracted DAC values. The index of the DataFrame is the date and the column contains the DAC values.

Return type:

DataFrame

geodezyx.marine.dac.dac_get_aviso(date_srt_inp, date_end_inp, out_dir, user, passwd)

Downloads DAC (Doppler-derived Atmospheric Correction) files from the AVISO FTP server for a given date range.

Parameters:
  • date_srt_inp (datetime) – The start date of the range for which DAC files are to be downloaded.

  • date_end_inp (datetime) – The end date of the range for which DAC files are to be downloaded.

  • out_dir (str) – The local directory where the downloaded DAC files should be saved.

  • user (str) – The username for the FTP server.

  • passwd (str) – The password for the FTP server.

Returns:

A list of paths to the downloaded DAC files.

Return type:

list

geodezyx.marine.marine module

@author: psakic

This sub-module of geodezyx.marine contains functions to interpolate GEBCO bathymetry.

it can be imported directly with: from geodezyx import marine

The GeodeZYX Toolbox is a software for simple but useful functions for Geodesy and Geophysics under the GNU LGPL v3 License

Copyright (C) 2019 Pierre Sakic et al. (IPGP, sakic@ipgp.fr) GitHub repository : https://github.com/GeodeZYX/geodezyx-toolbox

geodezyx.marine.marine.gebco_bathy_grid_extractor(dataset, latmin, latmax, lonmin, lonmax)

for safety reasons, lat and lon input MUST BE in the dataset, replaced by the closest elsewhere

return latnew , lonnew , Znew

geodezyx.marine.obp module

Created on Sat Aug 5 21:08:05 2023

@author: psakic

This module regroups the functions for the exploitation of the A0A pressure sensors in the context of the REVOSIMA network

It is based on the work of Yann-Terden Tranchant, LIENSs, La Rochelle, France

geodezyx.marine.obp.butter_bandpass(lowcut, highcut, fs, order=3)

Designs a bandpass Butterworth filter.

Parameters:
  • lowcut (float) – The lower cutoff frequency of the filter.

  • highcut (float) – The upper cutoff frequency of the filter.

  • fs (float) – The sampling frequency of the signal.

  • order (int, optional) – The order of the filter. Default is 3.

Returns:

  • b (ndarray) – The numerator coefficient vector of the filter.

  • a (ndarray) – The denominator coefficient vector of the filter.

geodezyx.marine.obp.butter_bandpass_filtfilt(data, lowcut, highcut, fs, order=4)

Applies a bandpass Butterworth filter to the data using forward and backward filtering.

Parameters:
  • data (array_like) – The input data to be filtered.

  • lowcut (float) – The lower cutoff frequency of the filter.

  • highcut (float) – The upper cutoff frequency of the filter.

  • fs (float) – The sampling frequency of the signal.

  • order (int, optional) – The order of the filter. Default is 4.

Returns:

y – The filtered data.

Return type:

ndarray

geodezyx.marine.obp.butter_highpass(cutoff, fs, order=5)

Designs a highpass Butterworth filter.

Parameters:
  • cutoff (float) – The cutoff frequency of the filter.

  • fs (float) – The sampling frequency of the signal.

  • order (int, optional) – The order of the filter. Default is 5.

Returns:

  • b (ndarray) – The numerator coefficient vector of the filter.

  • a (ndarray) – The denominator coefficient vector of the filter.

geodezyx.marine.obp.butter_highpass_filtfilt(data, cutoff, fs, order=4)

Applies a highpass Butterworth filter to the data using forward and backward filtering.

Parameters:
  • data (array_like) – The input data to be filtered.

  • cutoff (float) – The cutoff frequency of the filter.

  • fs (float) – The sampling frequency of the signal.

  • order (int, optional) – The order of the filter. Default is 4.

Returns:

y – The filtered data.

Return type:

ndarray

geodezyx.marine.obp.butter_lowpass(cutoff, fs, order=5)

Designs a lowpass Butterworth filter.

Parameters:
  • cutoff (float) – The cutoff frequency of the filter.

  • fs (float) – The sampling frequency of the signal.

  • order (int, optional) – The order of the filter. Default is 5.

Returns:

  • b (ndarray) – The numerator coefficient vector of the filter.

  • a (ndarray) – The denominator coefficient vector of the filter.

geodezyx.marine.obp.butter_lowpass_filtfilt(data, cutoff, fs, order=4)

Applies a lowpass Butterworth filter to the data using forward and backward filtering.

Parameters:
  • data (array_like) – The input data to be filtered.

  • cutoff (float) – The cutoff frequency of the filter.

  • fs (float) – The sampling frequency of the signal.

  • order (int, optional) – The order of the filter. Default is 4.

Returns:

y – The filtered data.

Return type:

ndarray

geodezyx.marine.obp.butterworth(df, t0=259200, t1=864000, kind='bandpass', order=4)

Applies a Butterworth filter to the data.

Parameters:
  • df (pandas.DataFrame) – The input data to be filtered.

  • t0 (int, optional) – The lower cutoff period in seconds. Default is 3 days (3*24*3600 seconds).

  • t1 (int, optional) – The upper cutoff period in seconds. Default is 10 days (10*24*3600 seconds).

  • kind (str, optional) – The type of filter to apply. Options are ‘highpass’, ‘lowpass’, and ‘bandpass’. Default is ‘bandpass’.

  • order (int, optional) – The order of the filter. Default is 4.

Returns:

The filtered data.

Return type:

pandas.Series

geodezyx.marine.obp.compute_dens_profile(profile)

Computes the density profile from the given dataset.

Parameters:

profile (xarray.Dataset) – The input dataset containing temperature (theta) and salinity (salt) variables.

Returns:

The computed density profile.

Return type:

numpy.ndarray

geodezyx.marine.obp.compute_phibot(profile, ssh=None, integration='forward', rho=10.35, remove_median=True)

Computes the ocean hydrostatic bottom pressure anomaly (phibot).

Parameters:
  • profile (xarray.Dataset) – The input dataset containing temperature (theta) and salinity (salt) variables.

  • ssh (xarray.DataArray, optional) – The sea surface height data. If None, it is taken from the profile. Default is None.

  • integration (str, optional) – The integration method to use for steric height computation. Options are ‘forward’, ‘backward’, and ‘averaged’. Default is ‘forward’.

  • rho (float, optional) – The reference density in kg/m^3. Default is 10.35.

  • remove_median (bool, optional) – Whether to remove the median from the computed phibot values. Default is True.

Returns:

The computed phibot values.

Return type:

pandas.DataFrame

Notes

phibot = Ocean hydrostatic bottom pressure anomaly https://cmr.earthdata.nasa.gov/search/concepts/V2028471168-POCLOUD.html https://cmr.earthdata.nasa.gov/search/concepts/V2146301108-POCLOUD.html

PHIBOT Bottom Pressure Pot. Anomaly (p/rhonil, m^2/s^2)

To convert to m, divide by g (g=9.81 m/s^2) PHIBOT is the anomaly relative to Depth * rhonil * g The absolute bottom pressure in Pa is: Depth * rhonil * g + PHIBOT * rhonil (rhonil=1027.5 kg/m^3)

http://apdrc.soest.hawaii.edu/doc/Readme_ecco2_cube92

geodezyx.marine.obp.compute_spectrogram(df, max_period=45, nchunks=21600)

Computes the spectrogram of the given data.

Parameters:
  • df (pandas.DataFrame) – The input data to compute the spectrogram for.

  • max_period (int, optional) – The maximum period to consider in the spectrogram. Default is 45.

  • nchunks (int, optional) – The number of chunks to divide the data into. Default is 3600*6.

Returns:

  • p (numpy.ndarray) – The periods corresponding to the spectrogram.

  • t (pandas.DatetimeIndex) – The time points corresponding to the spectrogram.

  • sxx (numpy.ndarray) – The spectrogram values.

geodezyx.marine.obp.compute_steric(profile, integration='forward')

Computes the steric height anomaly from the given profile.

Parameters:
  • profile (xarray.Dataset) – The input dataset containing temperature (theta) and salinity (salt) variables.

  • integration (str, optional) – The integration method to use. Options are ‘forward’, ‘backward’, and ‘averaged’. Default is ‘forward’.

Returns:

The steric height anomaly.

Return type:

pandas.Series

geodezyx.marine.obp.exp(x, a, b, c, d, **kwargs)

Exponential model function.

Parameters:
  • x (array-like) – The input data.

  • a (float) – Coefficient for the exponential term.

  • b (float) – Coefficient for the exponential term.

  • c (float) – Constant term.

  • **kwargs (dict) – Additional keyword arguments.

Returns:

The computed exponential value.

Return type:

float

geodezyx.marine.obp.exp_linear(x, a, b, c, d, **kwargs)

Exponential-linear model function.

Parameters:
  • x (array-like) – The input data.

  • a (float) – Coefficient for the exponential term.

  • b (float) – Coefficient for the exponential term.

  • c (float) – Linear term coefficient.

  • d (float) – Constant term.

  • **kwargs (dict) – Additional keyword arguments.

Returns:

The computed exponential-linear value.

Return type:

float

geodezyx.marine.obp.exp_polynomial_deg2(x, a, b, c, **kwargs)

Computes an exponential function combined with a second-degree polynomial.

This function calculates the value of an exponential function combined with a second-degree polynomial given the input data and coefficients.

Parameters:
  • x (array-like) – The input data.

  • a (float) – Coefficient for the polynomial term.

  • b (float) – Coefficient for the exponential term.

  • c (float) – Constant term.

  • **kwargs (dict) – Additional keyword arguments.

Returns:

The computed value from the exponential function combined with a second-degree polynomial.

Return type:

float

geodezyx.marine.obp.exp_power_combined(x, a, b, c, d, **kwargs)

Generalized nonlinear model function combining exponential decay and power law.

This function models data using a combination of an exponential decay term and a power law term.

Parameters:
  • x (array-like) – The input data.

  • a (float) – Coefficient for the exponential term.

  • b (float) – Exponent for the exponential term.

  • c (float) – Coefficient for the power law term.

  • d (float) – Exponent for the power law term.

  • **kwargs (dict) – Additional keyword arguments.

Returns:

The computed value from the combined exponential and power law model.

Return type:

float

geodezyx.marine.obp.exp_with_offset(x, a, b, c, **kwargs)

Exponential model function with an offset.

This function calculates an exponential decay with an added constant offset.

Parameters:
  • x (array-like) – The input data.

  • a (float) – Coefficient for the exponential term.

  • b (float) – Exponent for the exponential term.

  • c (float) – Constant offset term.

  • **kwargs (dict) – Additional keyword arguments.

Returns:

The computed value from the exponential model with offset.

Return type:

float

geodezyx.marine.obp.extract_profile(ds, x, y, method='linear', zbottom=None)

Extracts a profile from the dataset at the specified coordinates and interpolates to the given depth.

Parameters:
  • ds (xarray.Dataset) – The input dataset containing the data.

  • x (float) – The longitude coordinate.

  • y (float) – The latitude coordinate.

  • method (str, optional) – The interpolation method to use. Default is ‘linear’.

  • zbottom (int, float, list, or numpy.ndarray, optional) – The depth(s) to interpolate to. If None, no depth interpolation is performed. Default is None.

Returns:

The extracted and interpolated profile.

Return type:

xarray.Dataset

geodezyx.marine.obp.fit_model(data, t_fit, model='log_linear', offset_last=True, maxfev=1000, pn=2, normalize_t=False)

Fits a model to the given data.

Parameters:
  • data (pandas.Series) – The input data to fit the model to. Is a Series with time as Index

  • t_fit (pandas.DatetimeIndex) – Time points to get fitted values once the fit model has been determeined.

  • model (str, optional) – The model to use for fitting. Options are ‘log_linear’, ‘log’, ‘exp’, ‘linear’, ‘exp_linear’, and ‘poly’. Default is ‘log_linear’.

  • offset_last (bool, optional) – Whether to offset the last value. Default is True.

  • maxfev (int, optional) – The maximum number of function evaluations. Default is 1000.

  • pn (int, optional) – The degree of the polynomial if the model is ‘poly’. Default is 2.

  • normalize_t (bool, optional) – Whether to normalize the time values. Default is False.

Returns:

The fitted model values.

Return type:

pandas.Series

geodezyx.marine.obp.get_bottom_depth(ds)

Computes the bottom depth of the dataset.

Parameters:

ds (xarray.Dataset) – The input dataset containing depth and theta variables.

Returns:

The bottom depth values.

Return type:

numpy.ndarray

geodezyx.marine.obp.interp_time(ds, time, method='linear')

Interpolates the dataset along the time dimension.

Parameters:
  • ds (xarray.Dataset) – The input dataset to be interpolated.

  • time (array-like or scalar) – The new time points to interpolate to.

  • method (str or Litteral, optional) – The interpolation method to use. Default is ‘linear’.

Returns:

The interpolated dataset.

Return type:

xarray.Dataset

geodezyx.marine.obp.interp_xy(ds, x=None, y=None, method='linear')

frontend for xarray spatial interpolation

Parameters:
  • ds (xarray.Dataset) – input xarray grid.

  • x (float, optional) – longitude in deg (but depends on ds units). The default is None.

  • y (float, optional) – latitude in deg (but depends on ds units). The default is None.

  • method (str or Litteral, optional) – interpolation method. The default is ‘linear’.

Returns:

Interpolated values.

Return type:

ds

geodezyx.marine.obp.interp_z(ds, depth, method='linear')

Interpolates the dataset along the depth dimension.

Parameters:
  • ds (xarray.Dataset) – The input dataset to be interpolated.

  • depth (array-like or scalar) – The new depth points to interpolate to.

  • method (str or Litteral, optional) – The interpolation method to use. Default is ‘linear’.

Returns:

The interpolated dataset.

Return type:

xarray.Dataset

geodezyx.marine.obp.linear(x, a1, a2, a3, a4, **kwargs)

Linear model function.

Parameters:
  • x (array-like) – The input data.

  • a1 (float) – Slope of the linear function.

  • a2 (float) – Intercept of the linear function.

  • a3 (float) – Additional coefficient.

  • a4 (float) – Additional coefficient.

  • **kwargs (dict) – Additional keyword arguments.

Returns:

The computed linear value.

Return type:

float

geodezyx.marine.obp.log(x, a1, a2, a3, a4, **kwargs)

Logarithmic model function.

Parameters:
  • x (array-like) – The input data.

  • a1 (float) – Coefficient for the logarithmic term.

  • a2 (float) – Coefficient for the logarithmic term.

  • a3 (float) – Constant term.

  • **kwargs (dict) – Additional keyword arguments.

Returns:

The computed logarithmic value.

Return type:

float

geodezyx.marine.obp.log_linear(x, a1, a2, a3, a4, **kwargs)

Log-linear model function.

Parameters:
  • x (array-like) – The input data.

  • a1 (float) – Coefficient for the logarithmic term.

  • a2 (float) – Coefficient for the logarithmic term.

  • a3 (float) – Coefficient for the linear term.

  • a4 (float) – Constant term.

  • **kwargs (dict) – Additional keyword arguments.

Returns:

The computed log-linear value.

Return type:

float

geodezyx.marine.obp.polynomial_3(x, a, b, c, d, **kwargs)

Computes a 3rd degree polynomial.

This function calculates the value of a 3rd degree polynomial given the input data and coefficients.

Parameters:
  • x (array-like) – The input data.

  • a (float) – The constant term of the polynomial.

  • b (float) – The coefficient of the linear term.

  • c (float) – The coefficient of the quadratic term.

  • d (float) – The coefficient of the cubic term.

  • **kwargs (dict) – Additional keyword arguments.

Returns:

The computed value of the 3rd degree polynomial.

Return type:

float

geodezyx.marine.obp.read_duacs(file)

Reads a DUACS dataset and renames its variables.

Parameters:

file (str) – The path to the DUACS dataset file.

Returns:

The renamed dataset.

Return type:

xarray.Dataset

geodezyx.marine.obp.read_ecco2(file)

Reads an ECCO2 dataset and renames its variables.

Parameters:

file (str) – The path to the ECCO2 dataset file.

Returns:

The renamed dataset.

Return type:

xarray.Dataset

geodezyx.marine.obp.read_glorys(file)

Reads a GLORYS dataset and renames its variables.

Parameters:

file (str) – The path to the GLORYS dataset file.

Returns:

The renamed dataset.

Return type:

xarray.Dataset

geodezyx.marine.obp.read_hycom(file)

Reads a HYCOM dataset and renames its variables.

Parameters:

file (str) – The path to the HYCOM dataset file.

Returns:

The renamed dataset.

Return type:

xarray.Dataset

geodezyx.marine.obscom module

Created on Thu Nov 17 17:41:32 2022

@author: psakic

This module regroups the functions for the resolution estimation of the REVOSIMA’s OBSCOM-embeded pressure sensor

geodezyx.marine.obscom.counter2freq(n_counted_by_clk, count_sensor, freq_clk, integ_timespan=1, out='freq')

Convert count number to sensor frequency

Parameters:
  • n_counted_by_clk (int or float) – number of counts counted by the sensor.

  • count_sensor (int or float) – number of clock periods counted for one sensor count

  • freq_clk (int or float) – primary frequency of the clock that count period.

  • integ_timespan (int or float, optional) – integration timespan in seconds. The default is 1.

  • out (str, optional) – output ‘freq’ for freqency or ‘tau’ for period (=1/frequency). The default is ‘freq’.

Returns:

freqency (‘freq’) or period (‘tau’) of the sensor

Return type:

float

Note

usually, count_sensor=30e3 for pressure, count_sensor=168e3 for temperature and freq_clk=4.096e6

geodezyx.marine.obscom.freq2U(ft_in, U0)

Convert temperature sensor frequency to coefficient U

Parameters:
  • ft_in (float) – frequency of the temperature sensor in Hertz.

  • U0 (float) – input coefficient U0.

Returns:

coefficient U.

Return type:

float

geodezyx.marine.obscom.freq2counter(freq_or_tau_sensor, count_sensor, freq_or_tau_clk, integ_timespan=1, inp='freq', round_fct=<ufunc 'floor'>)

Convert sensor frequency to count number

Parameters:
  • freq_or_tau_sensor – freqency (‘freq’) or period (‘tau’) of the sensor.

  • count_sensor (int or float) – number of clock periods counted for one sensor count

  • freq_or_tau_clk – primary frequency/period of the clock that count period..

  • integ_timespan (int or float, optional) – integration timespan in seconds. The default is 1.

  • inp (str, optional) – input ‘freq’ for freqency or ‘tau’ for period (=1/frequency). The default is ‘freq’.

  • round_fct (function, optional) – the function that round the count values. The default is np.floor.

Returns:

number of counts counted by the sensor..

Return type:

n_count

Note

usually, count_sensor=30e3 for pressure, count_sensor=168e3 for temperature and freq_clk=4.096e6

geodezyx.marine.obscom.freq2pres(fp_in, t_or_ft_in, C1, C2, C3, D1, D2, U0, Y1, Y2, Y3, T1, T2, T3, T4, T5, out='psi', inp_temp='temp')

Convert pressure sensor frequency to pressure or equivalent distance

Parameters:
  • fp_in – frequency of the pressure sensor in Hertz.

  • t_or_ft_in (float) – temperature in Celsius or temperature sensor frequency in Hertz. depends on inp_temp

  • D2 (C1 C2 C3 D1) – Calibration coefficients provided by PAROS.

  • U0 (float) – Calibration coefficients provided by PAROS for the temperature sensor.

  • Y1 (float) – Calibration coefficients provided by PAROS for the temperature sensor.

  • Y2 (float) – Calibration coefficients provided by PAROS for the temperature sensor.

  • Y3 (float) – Calibration coefficients provided by PAROS for the temperature sensor.

  • T1 (float) – Calibration coefficients provided by PAROS for the temperature sensor.

  • T2 (float) – Calibration coefficients provided by PAROS for the temperature sensor.

  • T3 (float) – Calibration coefficients provided by PAROS for the temperature sensor.

  • T4 (float) – Calibration coefficients provided by PAROS for the temperature sensor.

  • T5 (float) – Calibration coefficients provided by PAROS for the temperature sensor.

  • out (str, optional) – output as pressure (‘psi’, ‘bar’, ‘mbar’, ‘pa’) or distance (‘meter’). The default is “psi”.

  • inp_temp (str) – ‘freq’ temperature sensor frequency in Hertz, ‘temp’ temperature in Celsius, default is ‘temp’

Returns:

pressure or equivalent distance.

Return type:

float

Note

Using temperature as input is unstable because of the root finding of U frequency as input is recommended

geodezyx.marine.obscom.freq2pres_old(f, U, F0, C1=-22682.65, C2=-1143.743, C3=70903.62, D1=0.040903, D2=0.0)

freqence capteur pression > pression OLD version

geodezyx.marine.obscom.freq2temp(ft_in, U0, Y1, Y2, Y3)

Convert temperature sensor frequency to temperature

Parameters:
  • ft_in – frequency of the temperature sensor in Hertz.

  • U0 (float) – Calibration coefficients provided by PAROS.

  • Y1-3 (float) – Calibration coefficients provided by PAROS.

Returns:

temperature in Celsius.

Return type:

float

geodezyx.marine.obscom.get_coeffs(sens_type='all', sens_id=158073)
geodezyx.marine.obscom.obscom_import_pres_df(csv_pres_inp, dic_coeffs, freq_clk=4096000.0, count_sensor_temp=168000.0, count_sensor_pres=30000.0, integ_timespan_temp=10, integ_timespan_pres=100, epoc_idx=True, out_pres='pa')

Imports pressure data from a CSV file and computes the frequency, temperature, pressure, and depth.

Parameters:
  • csv_pres_inp (str) – The path to the CSV file containing the pressure data.

  • dic_coeffs (dict) – The dictionary of coefficients used for the computations.

  • freq_clk (float, optional) – The frequency of the clock in Hertz. Default is 4.096e6.

  • count_sensor_temp (float, optional) – The count of the temperature sensor. Default is 168e3.

  • count_sensor_pres (float, optional) – The count of the pressure sensor. Default is 30e3.

  • integ_timespan_temp (int, optional) – The integration timespan for the temperature sensor. Default is 10.

  • integ_timespan_pres (int, optional) – The integration timespan for the pressure sensor. Default is 100.

  • epoc_idx (bool, optional) – If True, sets the epoch as the index of the DataFrame. Default is True.

  • out_pres (str, optional) – The unit of the output pressure. Default is “pa”.

Returns:

A DataFrame with the epoch, frequency, temperature, pressure, and depth.

Return type:

DataFrame

Notes

This function first reads the CSV file and creates a DataFrame. It then computes the frequency from the counter for both the temperature and pressure sensors. The temperature is computed from the frequency using the coefficients from the dictionary. The pressure is computed from the frequency of both the pressure and temperature sensors. The depth is also computed from the frequency of both the pressure and temperature sensors. If epoc_idx is True, the epoch is set as the index of the DataFrame.

geodezyx.marine.obscom.obscom_plot_pres_df(pres_df_in, plot_count=False)

Plots the pressure and temperature data from a DataFrame.

Parameters:
  • pres_df_in (DataFrame) – The DataFrame containing the pressure and temperature data. The index of the DataFrame should be the time.

  • plot_count (bool, optional) – If True, plots the count of the temperature and pressure sensors. If False, plots the pressure in pa and temperature in °C. Default is False.

Returns:

The Figure object with the plot of the pressure and temperature data.

Return type:

Figure

Notes

This function creates a plot with two y-axes. The left y-axis corresponds to the pressure data and the right y-axis corresponds to the temperature data. The x-axis corresponds to the time. The function uses different colors for the pressure and temperature data for better distinction.

geodezyx.marine.obscom.pres2freq(p_in, t_in, C1, C2, C3, D1, D2, U0, Y1, Y2, Y3, T1, T2, T3, T4, T5, inp='psi', return_optimize_object=False)

Convert pressure to pressure sensor frequency

There is no analytic solution, thus an root find is required

Parameters:
  • p_in – pressure in psi.

  • t_in – temperature in Celsius.

  • inp (optional) – input for p_in as pressure (‘psi’, ‘bar’, ‘mbar’, ‘pa’) or distance (‘meter’). The default is ‘psi’.

  • return_optimize_object (bool, optional) – if True, return the object of scipy’s optimize.root function. if False, return the root value The default is False.

Returns:

pressure sensor frequency (Hertz).

Return type:

float

geodezyx.marine.obscom.pres_resolution(val_presin, val_temp_in, val_clkin, count_presin, count_temp_in, input_val='freq', output_val='meter', relative_delta=True)
geodezyx.marine.obscom.resolution_grid_compute(Tau_pres, Tau_temp, fe, count_pres, count_temp, output_val='meter', relative_delta=False)
geodezyx.marine.obscom.resolution_plot_as_gradient_grid(PresVals, TempVals, ResVals, temp_ref, freq_input=False, coef_res=1000, relative_delta=False)
geodezyx.marine.obscom.temp2freq(t_in, U0, Y1, Y2, Y3, out='freq', force_root=False)

Convert temperature to temperature sensor frequency

Parameters:
  • t_in – temperature in Celsius.

  • U0 (float, optional) – Calibration coefficients provided by PAROS.

  • Y1-3 (float, optional) – Calibration coefficients provided by PAROS.

  • out (float, optional) – freqency (‘freq’) or period (‘tau’) or coefficient U (‘U’) of the sensor. The default is ‘freq’.

  • force_root (bool or int) – force the finded root for U if integer (1 or 2), use the first or second root default is False

Returns:

temperature sensor frequency or period.

Return type:

float

Note

This function is unstable because of the root finding

You can force it with force_root

geodezyx.marine.obscom.temp2freq_legacy(T, U0, Y, output='freq')
geodezyx.marine.obscom.temp2freq_old(T, out='freq')

temperature > frequence du capteur de temperature (output == “freq”) OU coef U (output == “U”) OU periode du capteur U (output == “tau”)

geodezyx.marine.obscom.unitary_tests(sens_id=158073)

return the correct values for pressure/temperature and intermediates parameters @ ftemp=172600.0 Hz and fpres=36300.0 Hz with the cofficients of the sensor #158073 (1st OBSCOM)

Return type:

ftemp (Hz), fpres (Hz), temp (°C), f0 (Hz), p (psi)