aboutsummaryrefslogtreecommitdiff
path: root/src/ParamChecker.c
blob: 7eb4d107efc25e7d010203dbdfabf41685335be2 (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
 /*@@
   @file      ParamChecker.F
   @date      March 1999
   @author    Gabrielle Allen
   @desc 
      Check parameters for black hole initial data and give some
      information
   @enddesc 
 @@*/

#include <stdio.h>
#include <stdlib.h>

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

#include "CactusEinstein/Einstein/src/Einstein.h"

 /*@@
   @routine    ParamChecker
   @date       March 1999
   @author     Gabrielle Allen
   @desc 
      Check parameters for black hole initial data and give some
      information
   @enddesc 
   @calls     
   @history 
 
   @endhistory 

@@*/

void ParamChecker(CCTK_CARGUMENTS)
{
  DECLARE_CCTK_CARGUMENTS
  DECLARE_CCTK_PARAMETERS

  char *message;

  if (CCTK_Equals(initial_data,"schwarzschild") == 1)
  {
    CCTK_INFO("Schwarzschild black hole");
    message = (char *)malloc(200*sizeof(char));
    sprintf(message,"  throat at %f",mass/2.0);
    CCTK_INFO(message);
    free(message);
  }
  else if (CCTK_Equals(initial_data,"bl_bh") == 1)
  {
    CCTK_INFO("Brill Lindquist black holes");
    message = (char *)malloc(200*sizeof(char));
    sprintf(message,"  %d black holes",bl_nbh);
    CCTK_INFO(message);
    if (bl_nbh > 0)
    {
      sprintf(message,  "  mass %f at (%f,%f,%f)",
	      bl_M_1,bl_x0_1,bl_y0_1,bl_z0_1);
      CCTK_INFO(message);
    }
    if (bl_nbh > 1)
    {
      sprintf(message,  "  mass %f at (%f,%f,%f)",
	      bl_M_2,bl_x0_2,bl_y0_2,bl_z0_2);
      CCTK_INFO(message);
    }
    if (bl_nbh > 2)
    {
      sprintf(message,  "  mass %f at (%f,%f,%f)",
	      bl_M_3,bl_x0_3,bl_y0_3,bl_z0_3);
      CCTK_INFO(message);
    }
    if (bl_nbh > 3)
    {
      sprintf(message,  "  mass %f at (%f,%f,%f)",
	      bl_M_4,bl_x0_4,bl_y0_4,bl_z0_4);
      CCTK_INFO(message);
    }
    free(message);
  }
  else if (CCTK_Equals(initial_data,"multiple_misner") == 1)
  {
    CCTK_INFO("Setting up Misner solution for multiple holes");
  }
  else if (CCTK_Equals(initial_data,"misner_bh")==1)
  {
    CCTK_INFO("Two Misner black holes (on z-axis)");
    message = (char *)malloc(200*sizeof(char));
    sprintf(message,"  mu is %f",mu);
    CCTK_INFO(message);
    free(message);
  }

  /*     Remind users about the conformal metric
   *     ---------------------------------------
   */
  
  if (use_conformal == 1)
  {
    CCTK_INFO("Black hole initial data uses conformal metric");
    if (use_conformal_derivs == 1)
    {
      CCTK_INFO("  and conformal derivatives");
    }
    else
    {
      CCTK_INFO("  but no conformal derivatives");
    }
  }
  else
  { 
    CCTK_INFO("Implements non-conformal metric");
    CCTK_INFO("  (Not usually a good idea!)");
  }        

  
  

}