src.Libs.Simulation package
Subpackages
- src.Libs.Simulation.Plots package
- src.Libs.Simulation.Telescope package
- Subpackages
- src.Libs.Simulation.Telescope.Functions package
- src.Libs.Simulation.Telescope.Specifications package
- Submodules
- src.Libs.Simulation.Telescope.Specifications.APEX module
- src.Libs.Simulation.Telescope.Specifications.Effelsberg module
- src.Libs.Simulation.Telescope.Specifications.EffelsbergSec module
- src.Libs.Simulation.Telescope.Specifications.Lovell module
- src.Libs.Simulation.Telescope.Specifications.MRO module
- src.Libs.Simulation.Telescope.Specifications.SRT_F1 module
- src.Libs.Simulation.Telescope.Specifications.SRT_F2 module
- Module contents
- Submodules
- src.Libs.Simulation.Telescope.TelescopeDataClass module
Telescope
Telescope.GregorianOrCass()
Telescope.Parabolic()
Telescope.addOptics()
Telescope.attemptTelescopeDataImport()
Telescope.createFromSpecs()
Telescope.createWithFunction()
Telescope.extractET()
Telescope.getAngularSizes()
Telescope.importAtmosphere()
Telescope.importTelescopeDict()
Telescope.makePrimary()
Telescope.makeSecondary()
Telescope.makeSky()
Telescope.makeSpill()
Telescope.plotTelescope()
Telescope.printOut()
Telescope.setTelescope()
Telescope.simulate()
Telescope.simulateAndAnalyze()
- Module contents
- Subpackages
- src.Libs.Simulation.UFO package
- Submodules
- src.Libs.Simulation.UFO.OpticsDataClass module
- src.Libs.Simulation.UFO.UFODataClass module
UFOData
UFOData.applyLosses()
UFOData.correlate()
UFOData.fitBeam()
UFOData.fromPoints()
UFOData.getAeff()
UFOData.getCenter()
UFOData.getCutData()
UFOData.getCutDataInter()
UFOData.getData()
UFOData.getPolarization()
UFOData.getPower()
UFOData.getRot()
UFOData.makeAperture()
UFOData.makeAsphLens()
UFOData.makeConicSection()
UFOData.makeFarField()
UFOData.makeFlat()
UFOData.makeLens()
UFOData.makeParabola()
UFOData.makeSphere()
UFOData.makeWaist()
UFOData.move()
UFOData.moveZero()
UFOData.plotCut()
UFOData.plotIntensity()
UFOData.plotMirror()
UFOData.plotPhase()
UFOData.printOut()
UFOData.readFile()
UFOData.resetPosition()
UFOData.rotX()
UFOData.rotY()
UFOData.rotZ()
UFOData.setI0()
UFOData.setLog()
UFOData.simulate()
UFOData.surfaceError()
UFOData.transform()
UFOData.writeFile()
- src.Libs.Simulation.UFO.UFOGaussCalc module
- src.Libs.Simulation.UFO.UFOTrace module
GeometricTraceData
GeometricTraceData.Trace()
GeometricTraceData.getPattern()
GeometricTraceData.getPower()
GeometricTraceData.getTraceLine()
GeometricTraceData.plotGeometricTrace()
GeometricTraceData.plotGeometricTrace3D()
GeometricTraceData.plotIntensity()
GeometricTraceData.plotPattern()
GeometricTraceData.printOut()
UFORay
- src.Libs.Simulation.UFO.UFO_CCode module
- src.Libs.Simulation.UFO.constants module
- Module contents
Submodules
src.Libs.Simulation.ArrayDataClass module
src.Libs.Simulation.AtmosphereDataClass module
This file is part of the PAFFrontendSim.
The PAFFrontendSim is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License v3 as published by the Free Software Foundation.
The PAFFrontendSim is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with the PAFFrontendSim. If not, see https://www.gnu.org/licenses/.
- class src.Libs.Simulation.AtmosphereDataClass.Atmosphere(log=None)[source]
Bases:
object
This module is used to read atmospheric transmission data from a .csv-file. The data is an optional component of the telescopes environment and meant to improve accuracy.
Initialize class.
- Parameters:
log (LogClass) – Documentation Class used by the calling application
- SkyTransmission(frequency)[source]
Return the interpolated sky temperature and sky transmittance at a given wavelength.
- Parameters:
frequency (float) – frequency of interpolated data
- plotData(freq_list=[], figname='', show=False)[source]
Plot the atmospheric data.
Transmittance and temperature over frequency.
- Parameters:
freq_list (list) – List of frequencies [GHz] at which the atmospheric data is calculated
figname (string) – file-name of the figure (if given figure is only shown on screen if show == True
show (boolean) – boolean flag to show plot on screen even when being written to disk (True)
src.Libs.Simulation.BeamWeightDataClass module
This file is part of the PAFFrontendSim.
The PAFFrontendSim is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License v3 as published by the Free Software Foundation.
The PAFFrontendSim is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with the PAFFrontendSim. If not, see https://www.gnu.org/licenses/.
- class src.Libs.Simulation.BeamWeightDataClass.BeamWeights(log=None, wl=214.0, t=0.0, p=0.0)[source]
Bases:
object
Class to handle beam-weight data
Initialize class.
- plot(figname='', show=False, freq=-1)[source]
Plotting the Beam-weights.
- Parameters:
figname (string) – filename of the created plot
show (boolean) – flag if plot is shown (True) or safed to disk
- printOut(line, message_type, noSpace=True, lineFeed=True)[source]
Wrapper for the printOut method of the actual log class
src.Libs.Simulation.CSTDataClass module
This file is part of the PAFFrontendSim.
The PAFFrontendSim is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License v3 as published by the Free Software Foundation.
The PAFFrontendSim is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with the PAFFrontendSim. If not, see https://www.gnu.org/licenses/.
- class src.Libs.Simulation.CSTDataClass.CST(log=None, project='')[source]
Bases:
object
This module provides two methods to read farfield data generated by CST (https://www.3ds.com/products/simulia/cst-studio-suite).
- Method 1 (readFarFieldSources): Reading farfield files (.ffs) exported by CST
This method requires no additional libraries and no local installation of the CST suite. It does require .ffs-files. CST has an export feature that allows farfield data to be stored in said format.
- Method 2: Using the official CST Python library to access project files (.cst)
An installation of CST contains Python libraries which are able to read the proprietary project files. This requires a local installation of CST and the project files containing your desired farfield data. Be aware that the Python libraries provided by CST only support specific Python versions. Make sure that you run this module with a supported version of Python when using the official CST Python library.
The PAF Frontend simulator requires some form of farfield data, though it doesn’t have to be generated by CST. For this reason, this module is optional and internal methods to generate farfield data may be used instead (see FarFieldDataClass.py).
Initialize class.
- Parameters :
project (string): Name of the CST project (folder name within the _CSTProjectPath folder). If empty, an empty object is created. log (LogClass): Logging class to be used
- getItemData(itemName, minX=None, maxX=None)[source]
Get the data of the item with name ‘itemName’.
- Parameters:
itemName (string) – Full or part name of the item (first occurance of the string in the list will be returned).
- Returns:
List (list) – of x-data, y-data, plot-title, x-axis label, y-axis label
or None
- plotItem(itemName, unit='Mag', figname='', show=False)[source]
Plot the data of the item with name ‘itemName’.
Breakdown of the units: “Mag” : using a linear scale of the magnutude of the value (if complex, amplitude is used) “dB” : using dB scale (power) of the magnutude of the value (if complex, amplitude is used) “Phase” : using the phase (only for complex y-data)
- Parameters:
itemName (string) – Full or part name of the item (first occurance of the string in the list will be returned).
unit (string) – Identifier to switch between amplitude, phase and liner / dB scale
figname (string) – File-name of the figure (if given figure is only shown on screen if show == True
show (boolean) – Flag to show plot on screen even when being written to disk (True)
- printOut(line, message_type, noSpace=True, lineFeed=True)[source]
Wrapper for the printOut method of the actual log class
- readFarFieldSources(cstImpedance=-1)[source]
Method to read a CST farfield source file (.ffs).
Currently, each farfield source file must contain only one frequency.
- Parameters:
cstImpedance (int) – Override CST simulation reference impedance (negative, use simulation setting)
- Returns:
dataValid (boolean) – error flag (true on success, false on error)
FarFields (list) – list of far-field data classes
src.Libs.Simulation.DistributionDataClass module
This file is part of the PAFFrontendSim.
The PAFFrontendSim is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License v3 as published by the Free Software Foundation.
The PAFFrontendSim is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with the PAFFrontendSim. If not, see https://www.gnu.org/licenses/.
- class src.Libs.Simulation.DistributionDataClass.DistributionData(log=None)[source]
Bases:
object
This module is used to create a distributions of points on a plane with boundary conditions.
It can be used to generate a list of antenna positions on the focal plane of a PAF. Input parameters are the size of the plane (XY), the type of distribution (square, hexagonal, random, etc), spacing between antennas and more.
Initialize class.
- createDistribution()[source]
Creates the distribution as set.
- Returns:
step (int?) – Spacing between the elements NOTE?
- getDisorder()[source]
Calculates the number variance var^2_N(R) of the point-counts in various size circles. var ~ R^(2-alpha) Algorithm calculates number variance in a sample of 500 circles of one size. Radius of the circle is then changed. The curve variance versus circle radius is fitted to f(x) = a*x^b and 2-b is returned This only works reliable for very large sets of distributions!
- Returns:
alpha (float) – ?
- plotDist(figname='', show=False)[source]
Plot the element usage of the beam-forming over the FoV.
- Parameters:
figname (string) – file-name of the figure (if given figure is only shown on screen if show == True
show (boolean) – flag to show plot on screen even when being written to disk (True)
src.Libs.Simulation.FarFieldDataClass module
This file is part of the PAFFrontendSim.
The PAFFrontendSim is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License v3 as published by the Free Software Foundation.
The PAFFrontendSim is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with the PAFFrontendSim. If not, see https://www.gnu.org/licenses/.
- class src.Libs.Simulation.FarFieldDataClass.FarFieldData(wavelength, raster=5, log=None)[source]
Bases:
object
Class to handle far field data
Initialize class. A wavelength can be given.
- adjustFarfield()[source]
Re-calculate internal rotation matrices, used after importing data including coordinate system (implemented into import function above).
- calculateCurrent(ufo, add=True)[source]
Calculate current distribution induced by the farfield on a UFO file.
- Parameters:
ufo (UFODataClass) – If Type is _CURRENT the farfield inducted currents will be added
add (boolean) – If the additional field and current distribution should be added to the existing
- Returns:
ufoNew (UFODataClass) – with the new data
- calculateDipoleFarfieldPoint(Phi, Theta, antennaLength=100)[source]
Calculates one point (Phi,Theta) of a farfield for a dipole with groundplane. The dipole is aligned with the X-axis and has a distance of wavelength/4 from the groundplane.
- Parameters:
Phi (float) – Phi angle of data point
Thete (float) – Theta angle of data point
antennaLength (float) – length of the dipole section (oriented along the y-axis)
- Returns:
list (list) – phi angle, theta angle, Phi component of farfield, Theta component of farfield
- calculateGaussianFarfieldPoint(Phi, Theta, intensity=1, width=90)[source]
Calculates one point (Phi,Theta) of a farfield for a Gaussian distribution.
With linear polarization vector in cross-Theta direction.
- Parameters:
Phi (float) – Phi angle of data point
Thete (float) – Theta angle of data point
intensity (float) – intensity of the distribution
width (float) – 1/e width of the resulting field (in degrees)
- Returns:
list (list) – phi angle, theta angle, Phi component of farfield, Theta component of farfield
- combineFarField(FFData, dPhi=5, dTheta=5)[source]
Combines far-field data of this instance with the once of an other instance and creates a new class instance. Returns FarFieldDataClass with the resulting field.
- Parameters:
FFData (FarFieldDataClass) – with the far-field to combine. For phase variations use setPhase in before hand
dPhi (float) – stepping in phi-direction of the new far-field [deg]
dTheta (float) – stepping in theta-direction of the new far-field [deg]
- Returns:
newFF (FarFieldDataClass) – of the combined port
- createFarfield(Type='Dipole', intensity=1, width=80, antennaLength=100, thetaMax=180)[source]
Creates a farfield of either a dipole or with a Gaussian distribution.
- Parameters:
intensity (float) – intensity modifier for Gaussian farfield.
width (float) – maximum angle of the Gaussian farfield (in degrees).
antennaLength (float) – length of the dipole section (oriented along the x-axis).
thetaMax (float) – maximum value for theta (degree)
- createWaist(waist=10, intensity=1, stepsPerWaist=50)[source]
creates a Gaussian far field distribution corresponding to the far field of the given waist size.
- estimatePhaseCenter()[source]
Try to estimate the phase center of the given far field distribution The result is stored to PhaseCenter
- getField(phi, th, interpolate=False)[source]
Method to interpolate a point within the loaded FarField data array.
- Parameters:
phi (float) – Angle Phi for which the data is interpolated
th (float) – Angle Theta for which the data is interpolated
interpolate (bool) – Whether or not the point should be interpolated or nearest
- Returns:
farfieldPoint (list) – The point on the farfield at phi,th in format [phi, th, Ephi, Eth]
- getFieldPoint(xp, yp, zp)[source]
Calculate E-field at position (xp, yp, zp).
- Parameters:
xp (float) – x coordinate
yp (float) – y coordinate
zp (float) – z coordinate
- Returns:
complex 3d E-vector (?)
normalized pointing vector (from FarField center-coordinate to surface point)
- getPower()[source]
Calculates the actual power in the far-field distribution and compares it with the excitation power via the excitation current.
- Returns:
p (float) – Far-Field power [W]
- importFarfieldFile(filename, setPhaseCenter=False, single=True, fileformat='')[source]
Method to import a Farfield data file in either CST or UFO format
- Parameters:
filename (string) – Name of the input file
setPhaseCenter (boolean) – Use phase center (or fit it in case of CST-files) – True; or set it to zero
single (boolean) – Handling a folder setup (False) or read a single file (True)
fileformat (string) – File-Format indicator, empty string “”: try out, ‘“UFO”: internal fime format, “CST”: CST fileformat
- Returns:
idOK (boolean) – error flag (true on success, false on error)
- makeMovie(movName='farfield.mov', folder='')[source]
Creates a movie on all far-field data in the folder set (if a folder was set).
- Parameters:
movName (string) – name of the movie
folder (string) – folder to store the individual plots and the movie
- Returns:
boolReturn (boolean) – flag to indicate whether a movie was created.
- moveCenter(x, y, z)[source]
Moves the FarField coordinate system to point (x, y, z)
- Parameters:
x (float) – x coordinate of new point
y (float) – y coordinate of new point
z (float) – z coordinate of new point
- movePhaseCenter(x, y, z)[source]
Moves phase center (rotation-point of the far-field distribution) to point (x, y, z)
- Parameters:
x (float) – x coordinate of new point
y (float) – y coordinate of new point
z (float) – z coordinate of new point
- plotFarField(figname='', show=False, logScale=True)[source]
Plotting the farfield distribution.
- Parameters:
figname (string) – filename of the created plot
show (boolean) – flag if plot is shown (True) or safed to disk
- plotFarField2d(figname='', show=False, rotate=False, upside=False, vmin=-30)[source]
Plotting the farfield distribution as 2d-false color plot including rotation.
- Parameters:
figname (string) – Filename of the created plot.
show (boolean) – Flag if plot is shown (True) or saved to disk.
rotate (boolean) – Flag indicating if rotation angels are included or not
upside (boolean) – Flip z-axis (has only effect if rotate is set to True)
vmin (float) – Indicating lower bound for dB scale
- plotFarFieldCut(phi=0, polAngle=0, dB=True, figname='', show=False)[source]
Plotting the farfield distribution as cut on a given phi-angle.
- Parameters:
phi (float or list of floats) – phi angle of the cut
polAngle (float) – Phi-angle of the polarization to look at
dB (boolean) – Indicating if a dB scale (True) or a lineare scale is used
figname (string) – Filename of the created plot
show (boolean) – Flag if plot is shown (True) or safed to disk
- plotOpeningAngle(figname='', show=False)[source]
Plots the average opening angle over frequency (if a folder is set).
- Parameters:
figname (string) – Filename of the created plot.
show (boolean) – Flag if plot is shown (True) or saved to disk.
- plotSParam(figname='', show=False)[source]
Plots the s-parameter data of the folder as set.
- Parameters:
figname (string) – Filename of the created plot.
show (boolean) – Flag if plot is shown (True) or saved to disk.
- printOut(line, message_type, noSpace=True, lineFeed=True)[source]
Wrapper for the printOut method of the actual log class.
- readCSTFarfieldFile(filename, cstImpedance=50.0)[source]
Method to read a CST ascii file with Far-Field data.
- Parameters:
filename (string) – Name of the input file
single (boolean?) – Handling a folder setup (False) or read a single file (True)
- Returns:
idOK (boolean) – error flag (true on success, false on error)
- readUFOFarfieldFile(filename)[source]
A method to read a farfield file in UFO (?) format.
This fileformat seems to correspond to the one found in the writeFile-Method of this class.
- Parameters:
filename (string) – the path of the file to be read
- Returns:
idOK (boolean) – whether or not reading the file was successful
- rotate(alpha_x, alpha_y, alpha_z)[source]
Rotates the actual coordinate system relative to x-, y- and z-axis
- Parameters:
alpha_x (float) – angle of rotation relative to x axis [Deg]
alpha_y (float) – angle of rotation relative to y axis [Deg]
alpha_z (float) – angle of rotation relative to z axis [Deg]
- scaleI0(newI0=200000.0)[source]
Scales the field-amplitudes such that they fullfill the new excitation current without changing the input impedance.
- Parameters:
newI0 (float) – new excitation current
- setCST(CSTproject, folder='', fileNameBase='ffx_', PolAngle=0.0)[source]
Read CST-files from the project, store them in an internal format to folder ‘folder’ and select this folder for data access.
- Parameters:
CSTproject (string) – name of the CST project to access
folder (string) – name of the folder to save the far-field sources (if empty, no access to far-field source data)
fileNameBase (string) – base-name for the new farfield files. Available farfields will be numbered in sequence.
PolAngle (float) – gives th polarization angle [deg] of the created source (0 deg = along x-axis)
- setFolder(folder, wildcard='*.*')[source]
Set a folder with far-field pattern for various frequencies.
- Parameters:
folder (string?) – name of the folder to be used as far-field source
wildcard (string?) – wildcard of the filename
- Returns:
list of frequencies available
- setFrequency(frequency)[source]
Loads the far-field file belonging to the given frequency (closest is chosen).
- Parameters:
frequency (float) – frequency [GHz] of the far-field pattern to be loaded (it is the nearest available frequency chosen!)
- setGain(gain, relative=True)[source]
Sets the gain of the system, which means the interpolated field is multiplied with the scalar factor gain.
- Parameters:
gain (float) – Factor to multiply the field with (set to zero gain, means correct back)
relative (boolean) – Correct the I0 value (True) or introduce a loss (False)
- setPhase(phase, absolute=False)[source]
Sets the phase of the system, which means that the field points are shifted in phase by this value.
- Parameters:
phase (float) – new phase [deg]
absolute (boolean) – phase is added to the actual value (False) or taken as new absolute phase value (True)
src.Libs.Simulation.FoVDataClass module
src.Libs.Simulation.ImpedanceMatrixDataClass module
This file is part of the PAFFrontendSim.
The PAFFrontendSim is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License v3 as published by the Free Software Foundation.
The PAFFrontendSim is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with the PAFFrontendSim. If not, see https://www.gnu.org/licenses/.
- class src.Libs.Simulation.ImpedanceMatrixDataClass.ImpedanceMatrix(wl=214.0, z=50.0, i0=200000.0, bw=1000000.0, dg=1.0, aa=0.0, log=None, name='')[source]
Bases:
object
Module to handle mutual impedance matrix (MIM) data.
In a PAF, the radiation pattern of each antenna is affected by every other antenna in the array. As a result there is a distinction between isolated element patterns (radiation pattern of an antenna in isolation) and embedded element patterns. The latter includes the effects of the other antennas.
The embedded element patterns can be determined through the MIM. The MIM (usually denoted as Z in literature) is an N x N matrix, where N is the amount of antennas and the entries are complex values in Ohm. It describes the relationship between voltage and current at each port: V = Z * I V and I are column vectors with length N.
The MIM can also be transformed into S-Parameters.
This class has a method (correlateFields) to determine the MIM, as well as methods to read and write MIM data to a file.
Initialize class.
- correlate(baseName, detectorGain=10000.0)[source]
Calculate mutual impedance matrix out of the individual element far-field data or spill-over data.
Method will be removed in future!!!! (not necessary anymore to work on files)
- Returns:
boolean (boolean) – True on success.
- correlateFields(e_fields, df_array)[source]
Correlate the fields given in e_fields with the corresponding surface array in df_array and fill the matrix.
- Parameters:
e_fields (numpy.array) – [element-no., sky position point no., 3 (E_vector complex, 3d)] of the e-fields
df_array (numpy.array) – [sky position point no.] of surface element sizes
- Returns:
boolean (boolean) – True if successful.
- getMutualPower(element)[source]
Calculates the power coupled to the other antennas.
- Parameters:
element (int) – Number of the element to calculate.
- Returns:
p (float) – Denoting the power in the mutual coupling to other elements.
float (float?) – Power in the element self overlap.
- plotMatrix(name='', figname='', show=False)[source]
Plot the far-field response of the full array on sky.
- Parameters:
name (string) – Name of the matrix (shown in the title)
figname (string) – Filename of the figure (if given figure is only shown on screen if show == True
show (boolean) – Flag to show plot on screen even when being written to disk (True)
- printOut(line, message_type, noSpace=True, lineFeed=True)[source]
Wrapper for the printOut method of the actual log class.
src.Libs.Simulation.NoiseDataClass module
src.Libs.Simulation.SignalChainDataClass module
This file is part of the PAFFrontendSim.
The PAFFrontendSim is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License v3 as published by the Free Software Foundation.
The PAFFrontendSim is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with the PAFFrontendSim. If not, see https://www.gnu.org/licenses/.
- class src.Libs.Simulation.SignalChainDataClass.SignalChain(log=None)[source]
Bases:
object
Class to handle signal chain data.
It is also able to read (readFile) and write (writeFile) signal chain data in its own format (PAF-Simulator signal chain noise data file). Alternatively, LNA data can be read if formatted as a Touchstone file. In this case, the module SnPDataClass is used. If no input file is provided, the LNA data must be passed as lists to the internal method “makeDataset”.
Initialize class.
- getData(f)[source]
Linear interpolation of the data set to get value for frequency f.
- Parameters:
f (float) – Frequency [GHz]
- Returns:
tuple (tuple) – tmin, rn, zopt, s11, s12, s21, s22
- makeDataset(Name, F, Tmin=[], Zopt=[], Rn=[], S11=[], S12=[], S21=[], S22=[])[source]
Make a data-set from pre-defined data points.
- Parameters:
Name (string) – Name of the signal-chain
F (list) – All other data will be interpolated to these frequencies [GHz] and are assumed to be equally distributed in the range Fmin to Fmax
Tmin (list) – List of minimum tempertures [K]
Rn (list) – List of noise resistances [Ohm]
Zopt (list) – List of input impedances [Ohm]
S11 (list) – List of S11 data [dB]
S12 (list) – List of S12 data [dB]
S21 (list) – List of S21 data [dB]
S22 (list) – List of S22 data [dB]
- plot(figname='', show=False, onlyNoise=False, Zmatch=-1.0)[source]
Plot the actual data set.
- Parameters:
figname (string) – Name of the graphics file written
show (boolean) – If True, plot will be shown on screen without being written to file
onlyNoise (boolean) – Plot only noise data
Zmatch (float) – Signal impedance. If > 0, the noise temperature belonging to this impedance is plotted as well
- plotSParam(figname='', show=False)[source]
Plot the actual S-Parameter data.
- Parameters:
figname (string) – Name of the graphics file written
show (boolean) – If True, plot will be shown on screen without being written to file
- printOut(line, message_type, noSpace=True, lineFeed=True)[source]
Wrapper for the printOut method of the actual log class.
- readFile(filename)[source]
Method to read an signal-chain noise data file.
- Parameters:
filename (string) – Name of the input file
- Returns:
idOK (boolean) – error flag (true on success, false on error)
src.Libs.Simulation.SnPDataClass module
This file is part of the PAFFrontendSim.
The PAFFrontendSim is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License v3 as published by the Free Software Foundation.
The PAFFrontendSim is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with the PAFFrontendSim. If not, see https://www.gnu.org/licenses/.
- class src.Libs.Simulation.SnPDataClass.SnPData(log=None)[source]
Bases:
object
Class to read LNA data input files.
The expected format is a Touchstone file for a 2 port network (.s2p). The data must be blocks of S-Parameters (9 columns) and/or Noise parameters (5 columns). The header of the file must specify the units as well as the impedance (”# ghz S db R 50” as an example).
Initialize class.
src.Libs.Simulation.constants module
This file is part of the PAFFrontendSim.
The PAFFrontendSim is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License v3 as published by the Free Software Foundation.
The PAFFrontendSim is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with the PAFFrontendSim. If not, see https://www.gnu.org/licenses/.
Module contents
This file is part of the PAFFrontendSim.
The PAFFrontendSim is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License v3 as published by the Free Software Foundation.
The PAFFrontendSim is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with the PAFFrontendSim. If not, see https://www.gnu.org/licenses/.