diff options
author | schnetter <> | 2004-01-11 11:20:00 +0000 |
---|---|---|
committer | schnetter <> | 2004-01-11 11:20:00 +0000 |
commit | 29c34c256c04683b6958a6ab7d10fbb6df1a1b5a (patch) | |
tree | fba312c5a89e039b9e496c3b5a9beaaf5bc5bf9b /CarpetExtra | |
parent | 660d24114d2d704e48051ef4d7771c9c2e39f34b (diff) |
Allow different kinds of initial data.
darcs-hash:20040111112039-07bb3-7ffb487a26951fb48a0b090f5aaa6d160a4620f0.gz
Diffstat (limited to 'CarpetExtra')
-rw-r--r-- | CarpetExtra/IDScalarWaveFO/param.ccl | 11 | ||||
-rw-r--r-- | CarpetExtra/IDScalarWaveFO/src/initialdata.F77 | 54 |
2 files changed, 38 insertions, 27 deletions
diff --git a/CarpetExtra/IDScalarWaveFO/param.ccl b/CarpetExtra/IDScalarWaveFO/param.ccl index 2283a5aca..58c4c07ca 100644 --- a/CarpetExtra/IDScalarWaveFO/param.ccl +++ b/CarpetExtra/IDScalarWaveFO/param.ccl @@ -1,5 +1,14 @@ # Parameter definitions for thorn IDScalarWaveFO -# $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDScalarWaveFO/param.ccl,v 1.1 2003/06/18 18:24:29 schnetter Exp $ +# $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDScalarWaveFO/param.ccl,v 1.2 2004/01/11 12:20:39 schnetter Exp $ + +RESTRICTED: + +KEYWORD initial_data "Type of initial data" +{ + "plane" :: "Plane wave" +} "plane" + +PRIVATE: CCTK_REAL wave_number[3] "Wave number" { diff --git a/CarpetExtra/IDScalarWaveFO/src/initialdata.F77 b/CarpetExtra/IDScalarWaveFO/src/initialdata.F77 index f65fc7084..9e193be6b 100644 --- a/CarpetExtra/IDScalarWaveFO/src/initialdata.F77 +++ b/CarpetExtra/IDScalarWaveFO/src/initialdata.F77 @@ -1,4 +1,4 @@ -c $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDScalarWaveFO/src/initialdata.F77,v 1.3 2003/11/05 16:18:40 schnetter Exp $ +c $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDScalarWaveFO/src/initialdata.F77,v 1.4 2004/01/11 12:20:39 schnetter Exp $ #include "cctk.h" #include "cctk_Arguments.h" @@ -14,31 +14,33 @@ c $Header: /home/eschnett/C/carpet/Carpet/CarpetExtra/IDScalarWaveFO/src/ini parameter (pi = 3.141592653589793238462643383279502884197169399375105820974944592307816406286208998628034825342117068d0) CCTK_REAL omega 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))) - psix(i,j,k) = amplitude * wave_number(1) / omega * 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))) - psiy(i,j,k) = amplitude * wave_number(2) / omega * 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))) - psiz(i,j,k) = amplitude * wave_number(3) / omega * 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))) + if (CCTK_EQUALS(initial_data, "plane")) 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) + 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))) + psix(i,j,k) = amplitude * wave_number(1) / omega * 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))) + psiy(i,j,k) = amplitude * wave_number(2) / omega * 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))) + psiz(i,j,k) = amplitude * wave_number(3) / omega * 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))) + end do end do end do - end do + end if end |