HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
ReadInputs.c
Go to the documentation of this file.
1 
6 #include <stdio.h>
7 #include <stdlib.h>
8 #include <string.h>
9 #include <basic.h>
10 #include <timeintegration.h>
11 #include <mpivars.h>
12 #include <simulation_object.h>
13 
93 int ReadInputs( void *s,
95  int nsims,
96  int rank
97  )
98 {
100  int n, ferr = 0;
101 
102  if (sim == NULL) {
103  printf("Error: simulation object array is NULL!\n");
104  printf("Please consider killing this run.\n");
105  return(1);
106  }
107 
108  if (!rank) {
109 
110  /* set some default values for optional inputs */
111  for (n = 0; n < nsims; n++) {
112  sim[n].solver.ndims = 1;
113  sim[n].solver.nvars = 1;
114  sim[n].solver.ghosts = 1;
115  sim[n].solver.dim_global = NULL;
116  sim[n].solver.dim_local = NULL;
117  sim[n].solver.dim_global_ex = NULL;
118  sim[n].mpi.iproc = NULL;
119  sim[n].mpi.N_IORanks = 1;
120  sim[n].solver.dt = 0.0;
121  sim[n].solver.n_iter = 0;
122  sim[n].solver.restart_iter = 0;
123  sim[n].solver.screen_op_iter = 1;
124  sim[n].solver.file_op_iter = 1000;
125  sim[n].solver.write_residual = 0;
126  sim[n].solver.flag_ib = 0;
127 #if defined(HAVE_CUDA)
128  sim[n].solver.use_gpu = 0;
129  sim[n].solver.gpu_device_no = -1;
130 #endif
131  strcpy(sim[n].solver.time_scheme ,"euler" );
132  strcpy(sim[n].solver.time_scheme_type ," " );
133  strcpy(sim[n].solver.spatial_scheme_hyp ,"1" );
134  strcpy(sim[n].solver.spatial_type_par ,_NC_1STAGE_ );
135  strcpy(sim[n].solver.spatial_scheme_par ,"2" );
136  strcpy(sim[n].solver.interp_type ,"characteristic");
137  strcpy(sim[n].solver.ip_file_type ,"ascii" );
138  strcpy(sim[n].solver.input_mode ,"serial" );
139  strcpy(sim[n].solver.output_mode ,"serial" );
140  strcpy(sim[n].solver.op_file_format ,"text" );
141  strcpy(sim[n].solver.op_overwrite ,"no" );
142  strcpy(sim[n].solver.plot_solution ,"no" );
143  strcpy(sim[n].solver.model ,"none" );
144  strcpy(sim[n].solver.ConservationCheck ,"no" );
145  strcpy(sim[n].solver.SplitHyperbolicFlux,"no" );
146  strcpy(sim[n].solver.ib_filename ,"none" );
147  }
148 
149  /* open the file */
150  FILE *in;
151  printf("Reading solver inputs from file \"solver.inp\".\n");
152  in = fopen("solver.inp","r");
153  if (!in) {
154  fprintf(stderr,"Error: File \"solver.inp\" not found.\n");
155  fprintf(stderr,"Please consider killing this run.\n");
156  return(1);
157  }
158 
159  /* reading solver inputs */
160  char word[_MAX_STRING_SIZE_];
161  ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
162 
163  if (!strcmp(word, "begin")){
164 
165  while (strcmp(word, "end")) {
166 
167  ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
168 
169  if (!strcmp(word, "ndims")) {
170 
171  ferr = fscanf(in,"%d",&(sim[0].solver.ndims)); if (ferr != 1) return(1);
172  sim[0].solver.dim_global = (int*) calloc (sim[0].solver.ndims,sizeof(int));
173  sim[0].mpi.iproc = (int*) calloc (sim[0].solver.ndims,sizeof(int));
174  sim[0].solver.dim_global_ex = (int*) calloc (sim[0].solver.ndims,sizeof(int));
175 
176  int n;
177  for (n = 1; n < nsims; n++) {
178  sim[n].solver.ndims = sim[0].solver.ndims;
179  sim[n].solver.dim_global = (int*) calloc (sim[n].solver.ndims,sizeof(int));
180  sim[n].mpi.iproc = (int*) calloc (sim[n].solver.ndims,sizeof(int));
181  sim[n].solver.dim_global_ex = (int*) calloc (sim[n].solver.ndims,sizeof(int));
182  }
183 
184  } else if (!strcmp(word, "nvars")) {
185 
186  ferr = fscanf(in,"%d",&(sim[0].solver.nvars));
187  for (int n = 1; n < nsims; n++) sim[n].solver.nvars = sim[0].solver.nvars;
188 
189  } else if (!strcmp(word, "size")) {
190 
191  for (int n = 0; n < nsims; n++) {
192  if (!sim[n].solver.dim_global) {
193  fprintf(stderr,"Error in ReadInputs(): dim_global not allocated for n=%d.\n", n);
194  fprintf(stderr,"Please specify ndims before dimensions.\n" );
195  return(1);
196  } else {
197  for (int i=0; i<sim[n].solver.ndims; i++) {
198  ferr = fscanf(in,"%d",&(sim[n].solver.dim_global[i]));
199  if (ferr != 1) {
200  fprintf(stderr,"Error in ReadInputs() while reading grid sizes for domain %d.\n", n);
201  return(1);
202  }
203  sim[n].solver.dim_global_ex[i] = sim[n].solver.dim_global[i];
204  }
205  }
206  }
207 
208  } else if (!strcmp(word, "size_exact")) {
209 
210  for (int n = 0; n < nsims; n++) {
211  if (!sim[n].solver.dim_global_ex) {
212  fprintf(stderr,"Error in ReadInputs(): dim_global_ex not allocated for n=%d.\n", n);
213  fprintf(stderr,"Please specify ndims before dimensions.\n" );
214  return(1);
215  } else {
216  for (int i=0; i<sim[n].solver.ndims; i++) {
217  ferr = fscanf(in,"%d",&(sim[n].solver.dim_global_ex[i]));
218  if (ferr != 1) {
219  fprintf(stderr,"Error in ReadInputs() while reading exact solution grid sizes for domain %d.\n", n);
220  return(1);
221  }
222  }
223  }
224  }
225 
226  } else if (!strcmp(word, "iproc")) {
227 
228  int n;
229  for (n = 0; n < nsims; n++) {
230  if (!sim[n].mpi.iproc) {
231  fprintf(stderr,"Error in ReadInputs(): iproc not allocated for n=%d.\n", n);
232  fprintf(stderr,"Please specify ndims before iproc.\n" );
233  return(1);
234  } else {
235  int i;
236  for (i=0; i<sim[n].solver.ndims; i++) {
237  ferr = fscanf(in,"%d",&(sim[n].mpi.iproc[i]));
238  if (ferr != 1) {
239  fprintf(stderr,"Error in ReadInputs() while reading iproc for domain %d.\n", n);
240  return(1);
241  }
242  }
243  }
244  }
245 
246  } else if (!strcmp(word, "ghost")) {
247 
248  ferr = fscanf(in,"%d",&(sim[0].solver.ghosts));
249 
250  int n;
251  for (n = 1; n < nsims; n++) sim[n].solver.ghosts = sim[0].solver.ghosts;
252 
253  } else if (!strcmp(word, "n_iter")) {
254 
255  ferr = fscanf(in,"%d",&(sim[0].solver.n_iter));
256 
257  int n;
258  for (n = 1; n < nsims; n++) sim[n].solver.n_iter = sim[0].solver.n_iter;
259 
260  } else if (!strcmp(word, "restart_iter")) {
261 
262  ferr = fscanf(in,"%d",&(sim[0].solver.restart_iter));
263 
264  int n;
265  for (n = 1; n < nsims; n++) sim[n].solver.restart_iter = sim[0].solver.restart_iter;
266 
267  } else if (!strcmp(word, "time_scheme")) {
268 
269  ferr = fscanf(in,"%s",sim[0].solver.time_scheme);
270 
271  int n;
272  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.time_scheme, sim[0].solver.time_scheme);
273 
274  } else if (!strcmp(word, "time_scheme_type" )) {
275 
276  ferr = fscanf(in,"%s",sim[0].solver.time_scheme_type);
277 
278  int n;
279  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.time_scheme_type, sim[0].solver.time_scheme_type);
280 
281  } else if (!strcmp(word, "hyp_space_scheme")) {
282 
283  ferr = fscanf(in,"%s",sim[0].solver.spatial_scheme_hyp);
284 
285  int n;
286  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.spatial_scheme_hyp, sim[0].solver.spatial_scheme_hyp);
287 
288  } else if (!strcmp(word, "hyp_flux_split")) {
289 
290  ferr = fscanf(in,"%s",sim[0].solver.SplitHyperbolicFlux);
291 
292  int n;
293  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.SplitHyperbolicFlux, sim[0].solver.SplitHyperbolicFlux);
294 
295  } else if (!strcmp(word, "hyp_interp_type")) {
296 
297  ferr = fscanf(in,"%s",sim[0].solver.interp_type);
298 
299  int n;
300  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.interp_type, sim[0].solver.interp_type);
301 
302  } else if (!strcmp(word, "par_space_type")) {
303 
304  ferr = fscanf(in,"%s",sim[0].solver.spatial_type_par);
305 
306  int n;
307  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.spatial_type_par, sim[0].solver.spatial_type_par);
308 
309  } else if (!strcmp(word, "par_space_scheme")) {
310 
311  ferr = fscanf(in,"%s",sim[0].solver.spatial_scheme_par);
312 
313  int n;
314  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.spatial_scheme_par, sim[0].solver.spatial_scheme_par);
315 
316  } else if (!strcmp(word, "dt")) {
317 
318  ferr = fscanf(in,"%lf",&(sim[0].solver.dt));
319 
320  int n;
321  for (n = 1; n < nsims; n++) sim[n].solver.dt = sim[0].solver.dt;
322 
323  } else if (!strcmp(word, "conservation_check" )) {
324 
325  ferr = fscanf(in,"%s",sim[0].solver.ConservationCheck);
326 
327  int n;
328  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.ConservationCheck, sim[0].solver.ConservationCheck);
329 
330  } else if (!strcmp(word, "screen_op_iter")) {
331 
332  ferr = fscanf(in,"%d",&(sim[0].solver.screen_op_iter));
333 
334  int n;
335  for (n = 1; n < nsims; n++) sim[n].solver.screen_op_iter = sim[0].solver.screen_op_iter;
336 
337  } else if (!strcmp(word, "file_op_iter")) {
338 
339  ferr = fscanf(in,"%d",&(sim[0].solver.file_op_iter));
340 
341  int n;
342  for (n = 1; n < nsims; n++) sim[n].solver.file_op_iter = sim[0].solver.file_op_iter;
343 
344  } else if (!strcmp(word, "op_file_format")) {
345 
346  ferr = fscanf(in,"%s",sim[0].solver.op_file_format);
347 
348  int n;
349  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.op_file_format, sim[0].solver.op_file_format);
350 
351  } else if (!strcmp(word, "ip_file_type")) {
352 
353  ferr = fscanf(in,"%s",sim[0].solver.ip_file_type);
354 
355  int n;
356  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.ip_file_type, sim[0].solver.ip_file_type);
357 
358  } else if (!strcmp(word, "input_mode")) {
359 
360  ferr = fscanf(in,"%s",sim[0].solver.input_mode);
361  if (strcmp(sim[0].solver.input_mode,"serial")) ferr = fscanf(in,"%d",&(sim[0].mpi.N_IORanks));
362 
363  int n;
364  for (n = 1; n < nsims; n++) {
365  strcpy(sim[n].solver.input_mode, sim[0].solver.input_mode);
366  if (strcmp(sim[n].solver.input_mode,"serial")) sim[n].mpi.N_IORanks = sim[0].mpi.N_IORanks;
367  }
368 
369  } else if (!strcmp(word, "output_mode")) {
370 
371  ferr = fscanf(in,"%s",sim[0].solver.output_mode);
372  if (strcmp(sim[0].solver.output_mode,"serial")) ferr = fscanf(in,"%d",&(sim[0].mpi.N_IORanks));
373 
374  int n;
375  for (n = 1; n < nsims; n++) {
376  strcpy(sim[n].solver.output_mode, sim[0].solver.output_mode);
377  if (strcmp(sim[n].solver.output_mode,"serial")) sim[n].mpi.N_IORanks = sim[0].mpi.N_IORanks;
378  }
379 
380  } else if (!strcmp(word, "op_overwrite")) {
381 
382  ferr = fscanf(in,"%s",sim[0].solver.op_overwrite);
383 
384  int n;
385  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.op_overwrite, sim[0].solver.op_overwrite);
386 
387  } else if (!strcmp(word, "plot_solution")) {
388 
389  ferr = fscanf(in,"%s",sim[0].solver.plot_solution);
390 
391  int n;
392  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.plot_solution, sim[0].solver.plot_solution);
393 
394  } else if (!strcmp(word, "model")) {
395 
396  ferr = fscanf(in,"%s",sim[0].solver.model);
397 
398  int n;
399  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.model, sim[0].solver.model);
400 
401  } else if (!strcmp(word, "immersed_body")) {
402 
403  ferr = fscanf(in,"%s",sim[0].solver.ib_filename);
404 
405  int n;
406  for (n = 1; n < nsims; n++) strcpy(sim[n].solver.ib_filename, sim[0].solver.ib_filename);
407 
408  }
409 #if defined(HAVE_CUDA)
410  else if (!strcmp(word, "use_gpu")) {
411  ferr = fscanf(in,"%s",word);
412  if (!strcmp(word, "yes") || !strcmp(word, "true")) sim[0].solver.use_gpu = 1;
413 
414  int n;
415  for (n = 1; n < nsims; n++) sim[n].solver.use_gpu = sim[0].solver.use_gpu;
416  } else if (!strcmp(word, "gpu_device_no")) {
417  ferr = fscanf(in,"%d", &sim[0].solver.gpu_device_no);
418 
419  int n;
420  for (n = 1; n < nsims; n++) sim[n].solver.gpu_device_no = sim[0].solver.gpu_device_no;
421  }
422 #endif
423  else if (strcmp(word, "end")) {
424 
425  char useless[_MAX_STRING_SIZE_];
426  ferr = fscanf(in,"%s",useless);
427  printf("Warning: keyword %s in file \"solver.inp\" with value %s not recognized or extraneous. Ignoring.\n",
428  word,useless);
429 
430  }
431  if (ferr != 1) return(1);
432 
433  }
434 
435  } else {
436 
437  fprintf(stderr,"Error: Illegal format in file \"solver.inp\".\n");
438  return(1);
439 
440  }
441 
442  /* close the file */
443  fclose(in);
444 
445  /* some checks */
446  for (n = 0; n < nsims; n++) {
447 
448  if (sim[n].solver.screen_op_iter <= 0) sim[n].solver.screen_op_iter = 1;
449  if (sim[n].solver.file_op_iter <= 0) sim[n].solver.file_op_iter = sim[n].solver.n_iter;
450 
451  if ((sim[n].solver.ndims != 3) && (strcmp(sim[n].solver.ib_filename,"none"))) {
452  printf("Warning: immersed boundaries not implemented for ndims = %d. ",sim[n].solver.ndims);
453  printf("Ignoring input for \"immersed_body\" (%s).\n",sim[n].solver.ib_filename);
454  strcpy(sim[n].solver.ib_filename,"none");
455  }
456  sim[n].solver.flag_ib = strcmp(sim[n].solver.ib_filename,"none");
457 
458  /* restart only supported for binary output files */
459  if ((sim[n].solver.restart_iter != 0) && strcmp(sim[n].solver.op_file_format,"binary")) {
460  if (!sim[n].mpi.rank) fprintf(stderr,"Error in ReadInputs(): Restart is supported only for binary output files.\n");
461  return(1);
462  }
463  }
464  }
465 
466 #ifndef serial
467  for (n = 0; n < nsims; n++) {
468 
469  /* Broadcast the input parameters */
470  MPIBroadcast_integer(&(sim[n].solver.ndims),1,0,&(sim[n].mpi.world));
471  if (sim[n].mpi.rank) {
472  sim[n].solver.dim_global = (int*) calloc (sim[n].solver.ndims,sizeof(int));
473  sim[n].mpi.iproc = (int*) calloc (sim[n].solver.ndims,sizeof(int));
474  sim[n].solver.dim_global_ex = (int*) calloc (sim[n].solver.ndims,sizeof(int));
475  }
476  MPIBroadcast_integer(&(sim[n].solver.nvars) ,1 ,0,&(sim[n].mpi.world));
477  MPIBroadcast_integer( sim[n].solver.dim_global ,sim[n].solver.ndims,0,&(sim[n].mpi.world));
478  MPIBroadcast_integer( sim[n].solver.dim_global_ex ,sim[n].solver.ndims,0,&(sim[n].mpi.world));
479  MPIBroadcast_integer( sim[n].mpi.iproc ,sim[n].solver.ndims,0,&(sim[n].mpi.world));
480  MPIBroadcast_integer(&(sim[n].mpi.N_IORanks) ,1 ,0,&(sim[n].mpi.world));
481  MPIBroadcast_integer(&(sim[n].solver.ghosts) ,1 ,0,&(sim[n].mpi.world));
482  MPIBroadcast_integer(&(sim[n].solver.n_iter) ,1 ,0,&(sim[n].mpi.world));
483  MPIBroadcast_integer(&(sim[n].solver.restart_iter) ,1 ,0,&(sim[n].mpi.world));
484  MPIBroadcast_integer(&(sim[n].solver.screen_op_iter),1 ,0,&(sim[n].mpi.world));
485  MPIBroadcast_integer(&(sim[n].solver.file_op_iter) ,1 ,0,&(sim[n].mpi.world));
486  MPIBroadcast_integer(&(sim[n].solver.flag_ib) ,1 ,0,&(sim[n].mpi.world));
487 #if defined(HAVE_CUDA)
488  MPIBroadcast_integer(&(sim[n].solver.use_gpu) ,1 ,0,&(sim[n].mpi.world));
489  MPIBroadcast_integer(&(sim[n].solver.gpu_device_no) ,1 ,0,&(sim[n].mpi.world));
490 #endif
491  MPIBroadcast_character(sim[n].solver.time_scheme ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
492  MPIBroadcast_character(sim[n].solver.time_scheme_type ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
493  MPIBroadcast_character(sim[n].solver.spatial_scheme_hyp ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
494  MPIBroadcast_character(sim[n].solver.interp_type ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
495  MPIBroadcast_character(sim[n].solver.spatial_type_par ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
496  MPIBroadcast_character(sim[n].solver.spatial_scheme_par ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
497  MPIBroadcast_character(sim[n].solver.ConservationCheck ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
498  MPIBroadcast_character(sim[n].solver.SplitHyperbolicFlux,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
499  MPIBroadcast_character(sim[n].solver.op_file_format ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
500  MPIBroadcast_character(sim[n].solver.ip_file_type ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
501  MPIBroadcast_character(sim[n].solver.input_mode ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
502  MPIBroadcast_character(sim[n].solver.output_mode ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
503  MPIBroadcast_character(sim[n].solver.op_overwrite ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
504  MPIBroadcast_character(sim[n].solver.plot_solution ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
505  MPIBroadcast_character(sim[n].solver.model ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
506  MPIBroadcast_character(sim[n].solver.ib_filename ,_MAX_STRING_SIZE_,0,&(sim[n].mpi.world));
507 
508  MPIBroadcast_double(&(sim[n].solver.dt),1,0,&(sim[n].mpi.world));
509  }
510 #endif
511 
512  return 0;
513 }
char plot_solution[_MAX_STRING_SIZE_]
Definition: hypar.h:194
char ib_filename[_MAX_STRING_SIZE_]
Definition: hypar.h:439
int nvars
Definition: hypar.h:29
MPI related function definitions.
Structure defining a simulation.
char SplitHyperbolicFlux[_MAX_STRING_SIZE_]
Definition: hypar.h:92
int ReadInputs(void *s, int nsims, int rank)
Definition: ReadInputs.c:93
int * dim_global_ex
Definition: hypar.h:75
int MPIBroadcast_integer(int *, int, int, void *)
Definition: MPIBroadcast.c:23
char ip_file_type[_MAX_STRING_SIZE_]
Definition: hypar.h:180
int n_iter
Definition: hypar.h:55
Some basic definitions and macros.
Contains function declarations for time integration.
Simulation object.
int flag_ib
Definition: hypar.h:441
char output_mode[_MAX_STRING_SIZE_]
Definition: hypar.h:183
char spatial_scheme_par[_MAX_STRING_SIZE_]
Definition: hypar.h:99
int restart_iter
Definition: hypar.h:58
char input_mode[_MAX_STRING_SIZE_]
Definition: hypar.h:177
int ndims
Definition: hypar.h:26
int MPIBroadcast_character(char *, int, int, void *)
Definition: MPIBroadcast.c:37
int MPIBroadcast_double(double *, int, int, void *)
Definition: MPIBroadcast.c:9
#define _MAX_STRING_SIZE_
Definition: basic.h:14
char ConservationCheck[_MAX_STRING_SIZE_]
Definition: hypar.h:376
char spatial_scheme_hyp[_MAX_STRING_SIZE_]
Definition: hypar.h:84
#define _NC_1STAGE_
Definition: hypar.h:475
char time_scheme[_MAX_STRING_SIZE_]
Definition: hypar.h:78
char op_overwrite[_MAX_STRING_SIZE_]
Definition: hypar.h:191
MPI_Comm world
char interp_type[_MAX_STRING_SIZE_]
Definition: hypar.h:88
char spatial_type_par[_MAX_STRING_SIZE_]
Definition: hypar.h:96
int * dim_local
Definition: hypar.h:37
double dt
Definition: hypar.h:67
char model[_MAX_STRING_SIZE_]
Definition: hypar.h:263
int gpu_device_no
Definition: hypar.h:450
int ghosts
Definition: hypar.h:52
char time_scheme_type[_MAX_STRING_SIZE_]
Definition: hypar.h:81
int use_gpu
Definition: hypar.h:449
char op_file_format[_MAX_STRING_SIZE_]
Definition: hypar.h:186
int * dim_global
Definition: hypar.h:33
int write_residual
Definition: hypar.h:174
int file_op_iter
Definition: hypar.h:171
int screen_op_iter
Definition: hypar.h:168