Skip to content
GitLab
Explore
Sign in
Register
Primary navigation
Search or go to…
Project
LISA Instrument
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Admin message
Gitlab has been updated. More info
here
and
there
.
Show more breadcrumbs
LISA Simulation
LISA Instrument
Commits
ce40e3a9
Commit
ce40e3a9
authored
4 years ago
by
Jean-Baptiste Bayle
Browse files
Options
Downloads
Patches
Plain Diff
Add instrumental noises
parent
ed31dce9
No related branches found
No related tags found
1 merge request
!1
Resolve "Add some basic noise generators"
Pipeline
#99608
passed
4 years ago
Stage: test
Changes
1
Pipelines
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
lisainstrument/noises.py
+111
-0
111 additions, 0 deletions
lisainstrument/noises.py
with
111 additions
and
0 deletions
lisainstrument/noises.py
+
111
−
0
View file @
ce40e3a9
...
...
@@ -65,3 +65,114 @@ def infrared(fs, size, psd):
return
0
red_noise
=
red
(
fs
,
size
,
psd
)
return
numpy
.
cumsum
(
red_noise
)
*
(
2
*
pi
/
fs
)
def
laser
(
fs
,
size
,
asd
=
28.2
):
"""
Generate laser noise [Hz].
This is a white noise,
S_p(f) = asd^2.
Args:
asd: amplitude spectral density [Hz/sqrt(Hz)]
"""
logging
.
debug
(
"
Generating laser noise (fs=%s Hz, size=%s, asd=%s Hz/sqrt(Hz))
"
,
fs
,
size
,
asd
)
return
white
(
fs
,
size
,
asd
**
2
)
def
clock
(
fs
,
size
,
asd
=
6.32E-14
):
"""
Generate clock noise fluctuations [ffd].
The power spectral density in fractional frequency deviations is a pink noise,
S_q(f) [ffd] = (asd)^2 f^(-1)
Args:
asd: amplitude spectral density [/sqrt(Hz)]
"""
logging
.
debug
(
"
Generating clock noise fluctuations (fs=%s Hz, size=%s, asd=%s /sqrt(Hz))
"
,
fs
,
size
,
asd
)
return
pink
(
fs
,
size
,
asd
**
2
)
def
modulation
(
fs
,
size
,
asd
=
1E-14
,
fknee
=
1.5E-2
):
"""
Generate modulation noise [ffd].
The power spectral density in timing jitter reads
S_M(f) [s] = (asd)^2 [ 1 + (fknee / f)^2 ].
Multiplying by (2π f)^2 to express it as fractional frequency deviations,
S_M(f) [ffd] = (2π asd)^2 [ f^2 + fknee^2 ]
= (2π asd fknee)^2 + (2π asd)^2 f^2.
Args:
asd: amplitude spectral density [s/sqrt(Hz)]
fknee: cutoff frequency [Hz]
"""
logging
.
debug
(
"
Generating modulation noise (fs=%s Hz, size=%s, asd=%s s/sqrt(Hz), fknee=%s Hz)
"
,
fs
,
size
,
asd
,
fknee
)
return
white
(
fs
,
size
,
(
2
*
pi
*
asd
*
fknee
)
**
2
)
\
+
violet
(
fs
,
size
,
(
2
*
pi
*
asd
)
**
2
)
def
backlink
(
fs
,
size
,
asd
=
3E-12
,
fknee
=
2E-3
):
"""
Generate backlink noise as fractional frequency deviation [ffd].
The power spectral density in displacement is given by
S_bl(f) [m] = asd^2 [ 1 + (fknee / f)^4 ].
Multiplying by (2πf / c)^2 to express it as fractional frequency deviations,
S_bl(f) [ffd] = (2π asd / c)^2 [ f^2 + (fknee^4 / f^2) ]
= (2π asd / c)^2 (f)^2 + (2π asd fknee^2 / c)^2 f^(-2)
Because this is a optical pathlength noise expressed as fractional frequency deviation, it should
be multiplied by the beam frequency to obtain the beam frequency fluctuations.
Args:
asd: amplitude spectral density [m/sqrt(Hz)]
fknee: cutoff frequency [Hz]
"""
logging
.
debug
(
"
Generating modulation noise (fs=%s Hz, size=%s, asd=%s m/sqrt(Hz), fknee=%s Hz)
"
,
fs
,
size
,
asd
,
fknee
)
return
violet
(
fs
,
size
,
(
2
*
pi
*
asd
/
c
)
**
2
)
\
+
red
(
fs
,
size
,
(
2
*
pi
*
asd
*
fknee
**
2
/
c
)
**
2
)
def
ranging
(
fs
,
size
,
asd
=
3E-9
):
"""
Generate stochastic ranging noise [s].
This is a white noise as a timing jitter,
S_R(f) [s] = asd.
Args:
asd: amplitude spectral density [s/sqrt(Hz)]
"""
logging
.
debug
(
"
Generating ranging noise (fs=%s Hz, size=%s, asd=%s s/sqrt(Hz))
"
,
fs
,
size
,
asd
)
return
white
(
fs
,
size
,
asd
**
2
)
def
testmass
(
fs
,
size
,
asd
=
2.4E-15
,
fknee
=
0.4E-3
):
"""
Generate test-mass acceleration noise [ffd].
Expressed in acceleration, the noise power spectrum reads
S_delta(f) [ms^(-2)] = (asd)^2 [ 1 + (fknee / f)^2 ].
Multiplying by 1 / (2π f c)^2 yields the noise as a fractional frequency deviations,
S_delta(f) [ms^(-2)] = (asd / (2π c))^2 [ f^(-2) + (fknee^2 / f^4) ]
= (asd / (2π c))^2 f^(-2) + (asd fknee / (2π c))^2 f^(-4) ].
Args:
asd: amplitude spectral density [ms^(-2)/sqrt(Hz)]
fknee: cutoff frequency [Hz]
"""
logging
.
debug
(
"
Generating test-mass noise (fs=%s Hz, size=%s, asd=%s ms^(-2)/sqrt(Hz), fknee=%s Hz)
"
,
fs
,
size
,
asd
,
fknee
)
return
red
(
fs
,
size
,
(
asd
/
(
2
*
pi
*
c
))
**
2
)
\
+
infrared
(
fs
,
size
,
(
asd
*
fknee
/
(
2
*
pi
*
c
))
**
2
)
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment