The FARGO_THORIN code developer's guide
 All Data Structures Files Functions Variables Typedefs Macros
LowTasks.c
Go to the documentation of this file.
1 /** \file LowTasks.c
2 
3 Contains many low level short functions. The name of these functions
4 should be self-explanatory in most cases. The prefix 'prs_' stands
5 for 'personal'. The prefix 'master' means that only the process 0
6 executes the function [note that the architecture is not of the kind
7 master/slaves, all processes perform similar tasks, but a minor
8 number of tasks (like output of information on the standard output) do
9 not need to be performed by several processes.] The function fopenp()
10 is an upper layer of fopen(), which should be used only in the case
11 of writing or appending a file (and not reading a file). It tries to
12 create the output directory if it does not exist, and it issues an
13 error message if it fails, so that the calling function does not
14 need to worry about these details.
15 */
16 
17 #include "fargo.h"
18 #include <stdarg.h>
19 
20 
22  real r;
23 {
24  int i=0;
25  real ifrac;
26  if (r < GlobalRmed[0]) return 0.0;
27  if (r > GlobalRmed[GLOBALNRAD-1]) return (real)GLOBALNRAD-1.0;
28  while (GlobalRmed[i] <= r) i++;
29  ifrac = (real)i+(r-GlobalRmed[i-1])/(GlobalRmed[i]-GlobalRmed[i-1])-1.0;
30  return ifrac;
31 }
32 
33 void prs_exit (numb)
34  int numb;
35 {
36  MPI_Finalize ();
37  exit (numb);
38 }
39 
40 void masterprint (const char *template, ...)
41 {
42  va_list ap;
43  if (!CPU_Master) return;
44  va_start (ap, template);
45  vfprintf (stdout, template, ap);
46  va_end (ap);
47 }
48 
49 void mastererr (const char *template, ...)
50 {
51  va_list ap;
52  if (!CPU_Master) return;
53  va_start (ap, template);
54  vfprintf (stderr, template, ap);
55  va_end (ap);
56 }
57 
58 void
59 prs_error(string)
60 char *string;
61 {
62  fprintf(stderr, "%s\n", string);
63  prs_exit(1);
64 }
65 
66 void message (msg)
67 char *msg;
68 {
69  fprintf (stdout, "%s", msg);
70 }
71 
72 PolarGrid *
73 CreatePolarGrid(Nr, Ns, name)
74 int Nr, Ns;
75 char *name;
76 {
77  PolarGrid *array;
78  real *field;
79  char *string;
80  int i, j, l;
81 
82  array = (PolarGrid *) malloc(sizeof(PolarGrid));
83  if (array == NULL)
84  prs_error("Insufficient memory for PolarGrid creation");
85  field = (real *) malloc(sizeof(real) * (Nr + 3) * (Ns + 1) + 5);
86  if (field == NULL)
87  prs_error("Insufficient memory for PolarGrid creation");
88  string = (char *) malloc(sizeof(char) * 80);
89  if (string == NULL)
90  prs_error("Insufficient memory for PolarGrid creation");
91  sprintf(string, "gas%s", name);
92  array->Field = field;
93  array->Name = string;
94  array->Nrad = Nr;
95  array->Nsec = Ns;
96  for (i = 0; i <= Nr; i++) {
97  for (j = 0; j < Ns; j++) {
98  l = j + i*Ns;
99  field[l] = 0.;
100  }
101  }
102  return array;
103 }
104 
105 
106 void MultiplyPolarGridbyConstant (arraysrc, constant)
107 PolarGrid *arraysrc;
108 real constant;
109 {
110  int i, nr, ns;
111  real *fieldsrc;
112 
113  nr = arraysrc->Nrad;
114  ns = arraysrc->Nsec;
115 
116  fieldsrc = arraysrc->Field;
117 #pragma omp parallel for
118  for (i = 0; i < (nr+1)*ns; i++) {
119  fieldsrc[i] *= constant;
120  }
121 }
122 
123 void DumpSources (argc, argv)
124  int argc;
125  char *argv[];
126 {
127  char CommandLine[1024];
128  char filecom[1024];
129  int i;
130  FILE *COM;
131  if (!CPU_Master) return;
132  sprintf (filecom, "%srun.commandline", OUTPUTDIR);
133  COM = fopenp (filecom, "w");
134  for (i = 0; i < argc; i++) {
135  fprintf (COM, "%s ",argv[i]);
136  }
137  fclose (COM);
138  sprintf (CommandLine, "cp .source.tar.bz2 %ssrc.tar.bz2", OUTPUTDIR);
139  system (CommandLine);
140 }
141 
142 void MakeDir (string)
143  char *string;
144 {
145  int foo=0;
146  char command[MAX1D];
147  DIR *dir;
148  /* Each processor tries to create the directory, sequentially */
149  /* Silent if directory exists */
150  if (CPU_Rank) MPI_Recv (&foo, 1, MPI_INT, CPU_Rank-1, 53, MPI_COMM_WORLD, &fargostat);
151  dir = opendir (string);
152  if (dir) {
153  closedir (dir);
154  } else {
155  fprintf (stdout, "Process %d creates the directory %s\n", CPU_Rank, string);
156  sprintf (command, "mkdir -p %s", string);
157  system (command);
158  }
159  if (CPU_Rank < CPU_Number-1) MPI_Send (&foo, 1, MPI_INT, CPU_Rank+1, 53, MPI_COMM_WORLD);
160 }
161 
162 
163 FILE *fopenp (string, mode)
164  char *string, *mode;
165 {
166  FILE *f;
167  f = fopen (string, mode);
168  if (f == NULL) {
169  /* This should be redundant with the call to MakeDir () at the
170  beginning, from main.c; this is not a problem however */
171  printf ("Process %d could not open %s\n", CPU_Rank, string);
172  printf ("Trying to create %s\n", OUTPUTDIR);
173  MakeDir (OUTPUTDIR);
174  f = fopen (string, "w"); /* "w" instead of mode: at this stage we know the file does not exist */
175  if (f == NULL) {
176  fprintf (stdout, "I still cannot open %s.\n", string);
177  fprintf (stdout, "You should check that the permissions are correctly set.\n");
178  fprintf (stdout, "Run aborted\n");
179  prs_exit (1);
180  }
181  }
182  return f;
183 }
double real
Definition of the type 'real' used throughout the code.
Definition: types.h:20
MPI_Status fargostat
Definition: global.h:32
void MPI_Send()
Definition: mpi_dummy.c:56
void MultiplyPolarGridbyConstant(PolarGrid *arraysrc, real constant)
Definition: LowTasks.c:106
char * Name
Pointer to the array of Nrad*Nsec reals (e.g., density, etc.)
Definition: types.h:41
int CPU_Number
Definition: global.h:2
void message(char *msg)
Definition: LowTasks.c:66
char OUTPUTDIR[512]
Definition: param_noex.h:11
PolarGrid * CreatePolarGrid(int Nr, int Ns, char *name)
Definition: LowTasks.c:73
A structure used to store any scalar fied on the computational domain.
Definition: types.h:37
real GlobalRmed[MAX1D]
Definition: global.h:13
void prs_error(char *string)
Definition: LowTasks.c:59
FILE * fopenp(char *string, char *mode)
Definition: LowTasks.c:163
void DumpSources(int argc, argv)
Definition: LowTasks.c:123
real * Field
Definition: types.h:40
int Nrad
Radial size of the grid, in number of zones.
Definition: types.h:38
void MPI_Finalize()
Definition: mpi_dummy.c:40
int GLOBALNRAD
Definition: global.h:10
void MakeDir(char *string)
Definition: LowTasks.c:142
#define MPI_INT
Definition: mpi_dummy.h:14
void mastererr(const char *template,...)
Definition: LowTasks.c:49
void MPI_Recv()
Definition: mpi_dummy.c:60
Contains all the include directives requested by the code.
int Nsec
Azimuthal size of the grid, in number of zones.
Definition: types.h:39
real GetGlobalIFrac(real r)
Definition: LowTasks.c:21
int CPU_Rank
Definition: global.h:1
boolean CPU_Master
Definition: global.h:3
#define MPI_COMM_WORLD
Definition: mpi_dummy.h:10
void prs_exit(int numb)
Definition: LowTasks.c:33
#define MAX1D
Definition: types.h:65
void masterprint(const char *template,...)
Definition: LowTasks.c:40