aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--doc/documentation.tex13
-rw-r--r--src/metrics/Minkowski_gauge_wave.F7736
2 files changed, 22 insertions, 27 deletions
diff --git a/doc/documentation.tex b/doc/documentation.tex
index 3467456..a10b8db 100644
--- a/doc/documentation.tex
+++ b/doc/documentation.tex
@@ -533,15 +533,16 @@ The line element is
\begin{equation}
ds^2=-H dt^2 +Hdx^2+dy^2+dz^2,
\end{equation}
-where $H=H(x-t)$, for instance $H=1-a*\sin\big((x-t)/d\big)$.
+where $H=H(x-t)$, for instance $H=1-a\sin\left((x-t)/\Lambda\right)$.
This is flat space but the slice is a planar wave travelling along the x axis.
This thorn implements several possible choices for the $H$ function,
controlled by the \verb|Minkowski_gauge_wave__what_fn| parameter:
\begin{eqnarray}
-H(x-t) &=& 1- A \sin \left(\frac{x-t}{d}\right) \\
-H(x-t) &=& \exp \left(A*sin(x)*cos(x)\right) \\
-H(x-t) &=& 1- A exp(-x^2) %%%\\
+H(x-t) &=& 1 - A \sin \left(\frac{x-\omega t}{\Lambda} - \delta\right) \\
+H(x-t) &=& \exp \left(- A \sin \left(\frac{x-\omega t}{\Lambda}
+ - \delta\right)\right) \\
+H(x-t) &=& 1 - A \exp(-x^2) %%%\\
\end{eqnarray}
The parameters are
\begin{itemize}
@@ -550,10 +551,6 @@ The parameters are
\item $\lambda = \verb|Minkowski_gauge_wave__lambda|$, the wavelength
\item $\delta = \verb|Minkowski_gauge_wave__phase|$, the phase shift
\end{itemize}
-Unfortunately, you have to look at the source code to see exactly what
-these do. :(
-FIXME: someone reverse-engineer the code and change the equations above
-to match whatever the code actually does.
If the Boolean parameter \verb|Minkowski_gauge_wave__diagonal| is true,
then we make the gauge wave travel diagonally across the grid by the
diff --git a/src/metrics/Minkowski_gauge_wave.F77 b/src/metrics/Minkowski_gauge_wave.F77
index 6ef07be..f4b6360 100644
--- a/src/metrics/Minkowski_gauge_wave.F77
+++ b/src/metrics/Minkowski_gauge_wave.F77
@@ -36,11 +36,8 @@ c output arguments
CCTK_REAL psi
LOGICAL Tmunu_flag
-c local static variables
- logical firstcall
+c local parameter copies
CCTK_REAL a, o, d, fs
- data firstcall /.true./
- save firstcall, a, o, d, fs
c local variables
CCTK_REAL H
@@ -54,29 +51,30 @@ C This is a vacuum spacetime with no cosmological constant
Tmunu_flag = .false.
C Get parameters of the exact solution.
-
- if (firstcall) then
- a = Minkowski_gauge_wave__amplitude
- o = Minkowski_gauge_wave__omega
- d = Minkowski_gauge_wave__lambda
- fs = Minkowski_gauge_wave__phase
-
- firstcall = .false.
- end if
+ a = Minkowski_gauge_wave__amplitude
+ o = Minkowski_gauge_wave__omega
+ d = Minkowski_gauge_wave__lambda
+ fs = Minkowski_gauge_wave__phase
C How should the wave look like.
- if (CCTK_Equals(Minkowski_gauge_wave__what_fn,"sin").eq.1) then
+ if (CCTK_EQUALS(Minkowski_gauge_wave__what_fn,"sin")) then
d = Minkowski_gauge_wave__lambda * half / Pi
if (Minkowski_gauge_wave__diagonal.ne.0) then
- H = one - a * sin((x-y)/d - o*t/d)
+ H = one - a * sin((x-y)/d - o*t/d - fs)
else
H = one - a * sin((x-o*t)/d - fs)
end if
- elseif (CCTK_Equals(Minkowski_gauge_wave__what_fn,"expsin").eq.1) then
+ elseif (CCTK_EQUALS(Minkowski_gauge_wave__what_fn,"expsin")) then
+c$$$ d = Minkowski_gauge_wave__lambda * half / Pi
+c$$$ H = exp(a*sin(x/d)*cos(t/d))
d = Minkowski_gauge_wave__lambda * half / Pi
- H = exp(a*sin(x/d)*cos(t/d))
- elseif (CCTK_Equals(Minkowski_gauge_wave__what_fn,"Gaussian").eq.1) then
- H = one - a*dexp(-(x-t)**2/d**2)
+ if (Minkowski_gauge_wave__diagonal.ne.0) then
+ H = exp(- a * sin((x-y)/d - o*t/d - fs))
+ else
+ H = exp(- a * sin((x-o*t)/d - fs))
+ end if
+ elseif (CCTK_EQUALS(Minkowski_gauge_wave__what_fn,"Gaussian")) then
+ H = one - a*exp(-(x-t)**2/d**2)
else
write (warn_buffer, '(a,a,a)')
$ 'Unknown Minkowski_gauge_wave__what_fn = "',