aboutsummaryrefslogtreecommitdiff
path: root/CarpetExtra/IDScalarWaveMoL
diff options
context:
space:
mode:
authorschnetter <>2003-08-14 15:22:00 +0000
committerschnetter <>2003-08-14 15:22:00 +0000
commit9424a56b281ab58546e77ee8177661a9f0b03ffc (patch)
treeb73645afc2cfee81a16217c1a74074be0bc3705a /CarpetExtra/IDScalarWaveMoL
parentb8ed9935ab73be68264d06f930a2a9f633e43ad8 (diff)
Add new kind of initial data "Gaussian pulse".
darcs-hash:20030814152209-07bb3-442fd3095423f18d9ad4c961b6745b65e5328a05.gz
Diffstat (limited to 'CarpetExtra/IDScalarWaveMoL')
-rw-r--r--CarpetExtra/IDScalarWaveMoL/param.ccl32
-rw-r--r--CarpetExtra/IDScalarWaveMoL/src/initialdata.F7747
2 files changed, 61 insertions, 18 deletions
diff --git a/CarpetExtra/IDScalarWaveMoL/param.ccl b/CarpetExtra/IDScalarWaveMoL/param.ccl
index 67d341ef5..b7d933e69 100644
--- a/CarpetExtra/IDScalarWaveMoL/param.ccl
+++ b/CarpetExtra/IDScalarWaveMoL/param.ccl
@@ -1,5 +1,15 @@
# Parameter definitions for thorn IDScalarWaveMoL
-# $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDScalarWaveMoL/param.ccl,v 1.1 2003/06/18 18:24:29 schnetter Exp $
+# $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDScalarWaveMoL/param.ccl,v 1.2 2003/08/14 17:22:09 schnetter Exp $
+
+KEYWORD initial_data "Initial data"
+{
+ "plane wave" :: "Plane wave"
+ "Gaussian pulse" :: "Gaussian pulse"
+} "plane wave"
+
+
+
+# Plane wave specification
CCTK_REAL wave_number[3] "Wave number"
{
@@ -11,12 +21,30 @@ CCTK_REAL phase_offset[3] "Phase offset"
*:* :: ""
} 0.0
+
+
+# Gaussian pulse specification
+
+CCTK_REAL pulse_direction[3] "Pulse width"
+{
+ *:* :: ""
+} 0.0
+
+CCTK_REAL pulse_offset[3] "Pulse offset"
+{
+ *:* :: ""
+} 0.0
+
+
+
+# Generic specifications
+
CCTK_REAL time_offset "Time offset"
{
*:* :: ""
} 0.0
-CCTK_REAL amplitude "Wave amplitude"
+CCTK_REAL amplitude "Amplitude"
{
*:* :: ""
} 1.0
diff --git a/CarpetExtra/IDScalarWaveMoL/src/initialdata.F77 b/CarpetExtra/IDScalarWaveMoL/src/initialdata.F77
index 1e7f5923c..638030eb0 100644
--- a/CarpetExtra/IDScalarWaveMoL/src/initialdata.F77
+++ b/CarpetExtra/IDScalarWaveMoL/src/initialdata.F77
@@ -1,4 +1,4 @@
-c $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDScalarWaveMoL/src/initialdata.F77,v 1.1 2003/06/18 18:24:29 schnetter Exp $
+c $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDScalarWaveMoL/src/initialdata.F77,v 1.2 2003/08/14 17:22:09 schnetter Exp $
#include "cctk.h"
#include "cctk_Arguments.h"
@@ -12,22 +12,37 @@ c $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDScalarWaveMoL/src/in
CCTK_REAL pi
parameter (pi = 3.141592653589793238462643383279502884197169399375105820974944592307816406286208998628034825342117068d0)
CCTK_REAL omega
+ CCTK_REAL tmp
integer i, j, k
- omega = sqrt(wave_number(1)**2 + wave_number(2)**2 + wave_number(3)**2)
- do k=1,cctk_lsh(3)
- do j=1,cctk_lsh(2)
- do i=1,cctk_lsh(1)
- phi(i,j,k) = amplitude * cos (2*pi *
- $ ( wave_number(1)*(x(i,j,k)-phase_offset(1))
- $ + wave_number(2)*(y(i,j,k)-phase_offset(2))
- $ + wave_number(3)*(z(i,j,k)-phase_offset(3))
- $ + omega*(cctk_time-time_offset)))
- psi(i,j,k) = - amplitude * 2*pi * omega * sin (2*pi *
- $ ( wave_number(1)*(x(i,j,k)-phase_offset(1))
- $ + wave_number(2)*(y(i,j,k)-phase_offset(2))
- $ + wave_number(3)*(z(i,j,k)-phase_offset(3))
- $ + omega*(cctk_time-time_offset)))
+ if (CCTK_EQUALS(initial_data, "plane wave")) then
+ omega = sqrt(wave_number(1)**2 + wave_number(2)**2
+ $ + wave_number(3)**2)
+ do k=1,cctk_lsh(3)
+ do j=1,cctk_lsh(2)
+ do i=1,cctk_lsh(1)
+ tmp = wave_number(1)*(x(i,j,k)-phase_offset(1))
+ $ + wave_number(2)*(y(i,j,k)-phase_offset(2))
+ $ + wave_number(3)*(z(i,j,k)-phase_offset(3))
+ $ + omega*(cctk_time-time_offset)
+ phi(i,j,k) = amplitude * cos (2*pi * tmp)
+ psi(i,j,k) = - amplitude * 2*pi * omega * sin (2*pi * tmp)
+ end do
end do
end do
- end do
+ else if (CCTK_EQUALS(initial_data, "Gaussian pulse")) then
+ omega = sqrt(pulse_direction(1)**2 + pulse_direction(2)**2
+ $ + pulse_direction(3)**2)
+ do k=1,cctk_lsh(3)
+ do j=1,cctk_lsh(2)
+ do i=1,cctk_lsh(1)
+ tmp = pulse_direction(1)*(x(i,j,k)-pulse_offset(1))
+ $ + pulse_direction(2)*(y(i,j,k)-pulse_offset(2))
+ $ + pulse_direction(3)*(z(i,j,k)-pulse_offset(3))
+ $ + omega*(cctk_time-time_offset)
+ phi(i,j,k) = amplitude * exp (-0.5d0 * tmp**2)
+ psi(i,j,k) = - tmp * omega * phi(i,j,k)
+ end do
+ end do
+ end do
+ end if
end