aboutsummaryrefslogtreecommitdiff
path: root/src/InitialData.c
blob: 5cb3f4035b487e391f3befb94938d704b31f3a25 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
 /*@@
   @file    InitialData.c
   @date
   @author  Werner Benger
   @desc
            Initial data for the 3D Wave Equation
            Derived from Tom Goodale
   @enddesc
   @version $Id$
 @@*/

#include <math.h>

#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"

static const char *rcsid = "$Header$";

CCTK_FILEVERSION(CactusWave_IDScalarWaveC_InitialData_c)

#define SQR(val) ((val) * (val))

void IDScalarWaveC_InitialData(CCTK_ARGUMENTS);


 /*@@
   @routine    IDScalarWaveC_InitialData
   @date
   @author     Tom Goodale
   @desc
               Set up initial data for the wave equation
   @enddesc
   @history
   @hdate Mon Oct 11 11:48:03 1999 @hauthor Werner Benger
   @hdesc  Converted to C++
   @hdate Mon Oct 11 11:48:20 1999 @hauthor Tom Goodale
   @hdesc Added the rest of the initial data.
   @hdate Thu Feb 17 09:22:01 2000 @hauthor Tom Goodale
   @hdesc Converted to C
   @endhistory

@@*/

void IDScalarWaveC_InitialData(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  int i,j,k;

  CCTK_REAL dt;
  CCTK_REAL omega;
  int vindex;
  CCTK_REAL X, Y, Z, R;
  CCTK_REAL pi;

  dt = CCTK_DELTA_TIME;

  if(CCTK_Equals(initial_data, "plane"))
  {
    omega = sqrt(SQR(kx)+SQR(ky)+SQR(kz));

    for(k=0; k<cctk_lsh[2]; k++)
    {
      for(j=0; j<cctk_lsh[1]; j++)
      {
        for(i=0; i<cctk_lsh[0]; i++)
        {
          vindex =  CCTK_GFINDEX3D(cctkGH,i,j,k);

          phi[vindex]     = amplitude*cos(kx*x[vindex]+ky*y[vindex]+kz*z[vindex]+omega*cctk_time);
          phi_p[vindex] = amplitude*cos(kx*x[vindex]+ky*y[vindex]+kz*z[vindex]+omega*(cctk_time-dt));
        }
      }
    }
  }
  else if(CCTK_Equals(initial_data, "gaussian"))
  {
    for(k=0; k<cctk_lsh[2]; k++)
    {
      for(j=0; j<cctk_lsh[1]; j++)
      {
        for(i=0; i<cctk_lsh[0]; i++)
        {
          vindex =  CCTK_GFINDEX3D(cctkGH,i,j,k);

          X = x[vindex];
          Y = y[vindex];
          Z = z[vindex];

          R = sqrt(X*X + Y*Y + Z*Z);

          phi[vindex] = amplitude*exp( - SQR( (R - radius) / sigma ) );

          if (R == 0.0)
          {
            phi_p[vindex] = amplitude*(1.0 - 2.0*dt*dt/sigma/sigma)*exp(-dt*dt/sigma/sigma);
          }
          else
          {
            phi_p[vindex] = amplitude/2.0*(R-dt)/R*
              exp( - SQR( (R - radius - dt)/ sigma ) )
              + amplitude/2.0*(R+dt)/R*
              exp( - SQR( (R - radius + dt)/ sigma ) );
          }
        }
      }
    }
  }
  else if(CCTK_Equals(initial_data, "box"))
  {
    pi = 4.0*atan(1.0);
    omega = sqrt(SQR(kx)+SQR(ky)+SQR(kz));

    for(k=0; k<cctk_lsh[2]; k++)
    {
      for(j=0; j<cctk_lsh[1]; j++)
      {
        for(i=0; i<cctk_lsh[0]; i++)
        {
          vindex =  CCTK_GFINDEX3D(cctkGH,i,j,k);

          phi[vindex] = amplitude*sin(kx*(x[vindex]-0.5)*pi)*
                                 sin(ky*(y[vindex]-0.5)*pi)*
                                 sin(kz*(z[vindex]-0.5)*pi)*
                                 cos(omega*cctk_time*pi);

          phi_p[vindex] = amplitude*sin(kx*(x[vindex]-0.5)*pi)*
                                       sin(ky*(y[vindex]-0.5)*pi)*
                                       sin(kz*(z[vindex]-0.5)*pi)*
                                       cos(omega*(cctk_time-dt)*pi);
        }
      }
    }
  }
  else if (CCTK_Equals(initial_data, "none"))
  {
    for(k=0; k<cctk_lsh[2]; k++)
    {
      for(j=0; j<cctk_lsh[1]; j++)
      {
        for(i=0; i<cctk_lsh[0]; i++)
        {
          vindex =  CCTK_GFINDEX3D(cctkGH,i,j,k);

          phi[vindex] = phi_p[vindex] = 0.0;
        }
      }
    }
  }
}