libMesh
femparameters.C
Go to the documentation of this file.
1 
2 #include "femparameters.h"
3 
4 #include "libmesh/parsed_function.h"
5 #include "libmesh/zero_function.h"
6 #include "libmesh/auto_ptr.h" // libmesh_make_unique
7 
8 #include <unordered_set>
9 
10 using namespace libMesh;
11 
12 #define GETPOT_INPUT(A) { A = input(#A, A); \
13  variable_names.insert(#A); \
14  const std::string stringval = input(#A, std::string()); \
15  variable_assignments.push_back(std::string(#A "=") + stringval); }
16 #define GETPOT_INT_INPUT(A) { A = input(#A, (int)A); \
17  variable_names.insert(#A); \
18  const std::string stringval = input(#A, std::string()); \
19  variable_assignments.push_back(std::string(#A "=") + stringval); }
20 
21 #define GETPOT_REGISTER(A) { \
22  variable_names.insert(#A); \
23  std::string stringval = input(#A, std::string()); \
24  variable_assignments.push_back(std::string(#A "=") + stringval); }
25 
26 FEMParameters::FEMParameters(const Parallel::Communicator & comm_in) :
27  ParallelObject(comm_in),
28  initial_timestep(0), n_timesteps(100),
29  transient(true),
30  deltat_reductions(0),
31  timesolver_core("euler"),
32  end_time(std::numeric_limits<Real>::max()),
33  deltat(0.0001), timesolver_theta(0.5),
34  timesolver_maxgrowth(0.), timesolver_tolerance(0.),
35  timesolver_upper_tolerance(0.),
36  steadystate_tolerance(0.),
37  timesolver_norm(0, L2),
38 
39  dimension(2),
40  domaintype("square"), domainfile("mesh.xda"), elementtype("quad"),
41  elementorder(2),
42  domain_xmin(0.0), domain_ymin(0.0), domain_zmin(0.0),
43  domain_edge_width(1.0), domain_edge_length(1.0), domain_edge_height(1.0),
44  coarsegridx(1), coarsegridy(1), coarsegridz(1),
45  coarserefinements(0), extrarefinements(0),
46  mesh_redistribute_func("0"),
47 
48  mesh_partitioner_type("Default"),
49 
50  nelem_target(8000), global_tolerance(0.0),
51  refine_fraction(0.3), coarsen_fraction(0.3), coarsen_threshold(10),
52  max_adaptivesteps(1),
53  initial_adaptivesteps(0),
54 
55  write_interval(10),
56  write_gmv_error(false), write_tecplot_error(false),
57  write_exodus_error(false),
58  output_xda(false), output_xdr(false),
59  output_bz2(true), output_gz(true),
60  output_gmv(false), output_tecplot(false),
61  output_exodus(false), output_nemesis(false),
62 
63  system_types(0),
64 
65 #ifdef LIBMESH_ENABLE_PERIODIC
66  periodic_boundaries(0),
67 #endif
68 
69  run_simulation(true), run_postprocess(false),
70 
71  fe_family(1, "LAGRANGE"), fe_order(1, 1),
72  extra_quadrature_order(0),
73 
74  analytic_jacobians(true), verify_analytic_jacobians(0.0),
75  numerical_jacobian_h(TOLERANCE),
76  print_solution_norms(false), print_solutions(false),
77  print_residual_norms(false), print_residuals(false),
78  print_jacobian_norms(false), print_jacobians(false),
79  print_element_solutions(false),
80  print_element_residuals(false),
81  print_element_jacobians(false),
82 
83  use_petsc_snes(false),
84  time_solver_quiet(true), solver_quiet(true), solver_verbose(false),
85  reuse_preconditioner(true),
86  require_residual_reduction(true),
87  min_step_length(1e-5),
88  max_linear_iterations(200000), max_nonlinear_iterations(20),
89  relative_step_tolerance(1.e-7), relative_residual_tolerance(1.e-10),
90  absolute_residual_tolerance(1.e-10),
91  initial_linear_tolerance(1.e-3), minimum_linear_tolerance(TOLERANCE*TOLERANCE),
92  linear_tolerance_multiplier(1.e-3),
93 
94  initial_sobolev_order(1),
95  initial_extra_quadrature(0),
96  refine_uniformly(false),
97  indicator_type("kelly"), patch_reuse(true), sobolev_order(1)
98 {
99 }
100 
102 {
103  for (std::map<subdomain_id_type, FunctionBase<Number> *>::iterator
104  i = initial_conditions.begin(); i != initial_conditions.end();
105  ++i)
106  delete i->second;
107 
108  for (std::map<boundary_id_type, FunctionBase<Number> *>::iterator
109  i = dirichlet_conditions.begin(); i != dirichlet_conditions.end();
110  ++i)
111  delete i->second;
112 
113  for (std::map<boundary_id_type, FunctionBase<Number> *>::iterator
114  i = neumann_conditions.begin(); i != neumann_conditions.end();
115  ++i)
116  delete i->second;
117 
118  for (std::map<int,
120  >::iterator
121  i = other_boundary_functions.begin(); i != other_boundary_functions.end();
122  ++i)
123  for (std::map<boundary_id_type, FunctionBase<Number> *>::iterator
124  j = i->second.begin(); j != i->second.end();
125  ++j)
126  delete j->second;
127 
128  for (std::map<int,
130  >::iterator
131  i = other_interior_functions.begin(); i != other_interior_functions.end();
132  ++i)
133  for (std::map<subdomain_id_type, FunctionBase<Number> *>::iterator
134  j = i->second.begin(); j != i->second.end();
135  ++j)
136  delete j->second;
137 }
138 
139 
140 std::unique_ptr<FunctionBase<Number>> new_function_base(const std::string & func_type,
141  const std::string & func_value)
142 {
143  if (func_type == "parsed")
144  return libmesh_make_unique<ParsedFunction<Number>>(func_value);
145  else if (func_type == "zero")
146  return libmesh_make_unique<ZeroFunction<Number>>();
147  else
148  libmesh_not_implemented();
149 
150  return std::unique_ptr<FunctionBase<Number>>();
151 }
152 
153 
154 void FEMParameters::read(GetPot & input,
155  const std::vector<std::string> * other_variable_names)
156 {
157  std::vector<std::string> variable_assignments;
158  std::unordered_set<std::string> variable_names;
159  if (other_variable_names)
160  for (std::size_t i=0; i != other_variable_names->size(); ++i)
161  {
162  const std::string & name = (*other_variable_names)[i];
163  const std::string stringval = input(name, std::string());
164  variable_assignments.push_back(name + "=" + stringval);
165  }
166 
167  GETPOT_INT_INPUT(initial_timestep);
168  GETPOT_INT_INPUT(n_timesteps);
169  GETPOT_INPUT(transient);
170  GETPOT_INT_INPUT(deltat_reductions);
171  GETPOT_INPUT(timesolver_core);
172  GETPOT_INPUT(end_time);
173  GETPOT_INPUT(deltat);
174  GETPOT_INPUT(timesolver_theta);
175  GETPOT_INPUT(timesolver_maxgrowth);
176  GETPOT_INPUT(timesolver_tolerance);
177  GETPOT_INPUT(timesolver_upper_tolerance);
178  GETPOT_INPUT(steadystate_tolerance);
179 
180  GETPOT_REGISTER(timesolver_norm);
181  const unsigned int n_timesolver_norm = input.vector_variable_size("timesolver_norm");
182  timesolver_norm.resize(n_timesolver_norm, L2);
183  for (unsigned int i=0; i != n_timesolver_norm; ++i)
184  {
185  int current_norm = 0; // L2
186  if (timesolver_norm[i] == H1)
187  current_norm = 1;
188  if (timesolver_norm[i] == H2)
189  current_norm = 2;
190  current_norm = input("timesolver_norm", current_norm, i);
191  if (current_norm == 0)
192  timesolver_norm[i] = L2;
193  else if (current_norm == 1)
194  timesolver_norm[i] = H1;
195  else if (current_norm == 2)
196  timesolver_norm[i] = H2;
197  else
199  }
200 
201 
202  GETPOT_INT_INPUT(dimension);
203  GETPOT_INPUT(domaintype);
204  GETPOT_INPUT(domainfile);
205  GETPOT_INPUT(elementtype);
206  GETPOT_INPUT(elementorder);
207  GETPOT_INPUT(domain_xmin);
208  GETPOT_INPUT(domain_ymin);
209  GETPOT_INPUT(domain_zmin);
210  GETPOT_INPUT(domain_edge_width);
211  GETPOT_INPUT(domain_edge_length);
212  GETPOT_INPUT(domain_edge_height);
213  GETPOT_INT_INPUT(coarsegridx);
214  GETPOT_INT_INPUT(coarsegridy);
215  GETPOT_INT_INPUT(coarsegridz);
216  GETPOT_INT_INPUT(coarserefinements);
217  GETPOT_INT_INPUT(extrarefinements);
218  GETPOT_INPUT(mesh_redistribute_func);
219 
220 
221  GETPOT_INPUT(mesh_partitioner_type);
222 
223 
224  GETPOT_INT_INPUT(nelem_target);
225  GETPOT_INPUT(global_tolerance);
226  GETPOT_INPUT(refine_fraction);
227  GETPOT_INPUT(coarsen_fraction);
228  GETPOT_INPUT(coarsen_threshold);
229  GETPOT_INT_INPUT(max_adaptivesteps);
230  GETPOT_INT_INPUT(initial_adaptivesteps);
231 
232 
233  GETPOT_INT_INPUT(write_interval);
234  GETPOT_INPUT(output_xda);
235  GETPOT_INPUT(output_xdr);
236  GETPOT_INPUT(output_gz);
237 #ifndef LIBMESH_HAVE_GZSTREAM
238  output_gz = false;
239 #endif
240  GETPOT_INPUT(output_bz2);
241 #ifndef LIBMESH_HAVE_BZ2
242  output_bz2 = false;
243 #endif
244  GETPOT_INPUT(output_gmv);
245  GETPOT_INPUT(write_gmv_error);
246 #ifndef LIBMESH_HAVE_GMV
247  output_gmv = false;
248  write_gmv_error = false;
249 #endif
250  GETPOT_INPUT(output_tecplot);
251  GETPOT_INPUT(write_tecplot_error);
252 #ifndef LIBMESH_HAVE_TECPLOT_API
253  output_tecplot = false;
254  write_tecplot_error = false;
255 #endif
256  GETPOT_INPUT(output_exodus);
257  GETPOT_INPUT(write_exodus_error);
258 #ifndef LIBMESH_HAVE_EXODUS_API
259  output_exodus = false;
260  write_exodus_error = false;
261 #endif
262  GETPOT_INPUT(output_nemesis);
263 #ifndef LIBMESH_HAVE_NEMESIS_API
264  output_nemesis = false;
265 #endif
266 
267 
268  GETPOT_REGISTER(system_types);
269  const unsigned int n_system_types =
270  input.vector_variable_size("system_types");
271  if (n_system_types)
272  {
273  system_types.resize(n_system_types, "");
274  for (unsigned int i=0; i != n_system_types; ++i)
275  {
276  system_types[i] = input("system_types", system_types[i], i);
277  }
278  }
279 
280 
281 #ifdef LIBMESH_ENABLE_PERIODIC
282  GETPOT_REGISTER(periodic_boundaries);
283  const unsigned int n_periodic_bcs =
284  input.vector_variable_size("periodic_boundaries");
285 
286  if (n_periodic_bcs)
287  {
288  if (domaintype != "square" &&
289  domaintype != "cylinder" &&
290  domaintype != "file" &&
291  domaintype != "od2")
292  {
293  libMesh::out << "Periodic boundaries need rectilinear domains" << std::endl;;
294  libmesh_error();
295  }
296  for (unsigned int i=0; i != n_periodic_bcs; ++i)
297  {
298  const unsigned int myboundary =
299  input("periodic_boundaries", -1, i);
300  boundary_id_type pairedboundary = 0;
301  RealVectorValue translation_vector;
302  if (dimension == 2)
303  switch (myboundary)
304  {
305  case 0:
306  pairedboundary = 2;
307  translation_vector = RealVectorValue(0., domain_edge_length);
308  break;
309  case 1:
310  pairedboundary = 3;
311  translation_vector = RealVectorValue(-domain_edge_width, 0);
312  break;
313  default:
314  libMesh::out << "Unrecognized periodic boundary id " <<
315  myboundary << std::endl;;
316  libmesh_error();
317  }
318  else if (dimension == 3)
319  switch (myboundary)
320  {
321  case 0:
322  pairedboundary = 5;
323  translation_vector = RealVectorValue(Real(0), Real(0), domain_edge_height);
324  break;
325  case 1:
326  pairedboundary = 3;
327  translation_vector = RealVectorValue(Real(0), domain_edge_length, Real(0));
328  break;
329  case 2:
330  pairedboundary = 4;
331  translation_vector = RealVectorValue(-domain_edge_width, Real(0), Real(0));
332  break;
333  default:
334  libMesh::out << "Unrecognized periodic boundary id " <<
335  myboundary << std::endl;;
336  libmesh_error();
337  }
338  periodic_boundaries.push_back(PeriodicBoundary(translation_vector));
339  periodic_boundaries[i].myboundary = cast_int<boundary_id_type>(myboundary);
340  periodic_boundaries[i].pairedboundary = pairedboundary;
341  }
342  }
343 #endif // LIBMESH_ENABLE_PERIODIC
344 
345  // Use std::string inputs so GetPot doesn't have to make a bunch
346  // of internal C string copies
347  std::string zero_string = "zero";
348  std::string empty_string = "";
349 
350  GETPOT_REGISTER(dirichlet_condition_types);
351  GETPOT_REGISTER(dirichlet_condition_values);
352  GETPOT_REGISTER(dirichlet_condition_boundaries);
353  GETPOT_REGISTER(dirichlet_condition_variables);
354 
355  const unsigned int n_dirichlet_conditions=
356  input.vector_variable_size("dirichlet_condition_types");
357 
358  if (n_dirichlet_conditions !=
359  input.vector_variable_size("dirichlet_condition_values"))
360  {
361  libMesh::out << "Error: " << n_dirichlet_conditions
362  << " Dirichlet condition types does not match "
363  << input.vector_variable_size("dirichlet_condition_values")
364  << " Dirichlet condition values." << std::endl;
365 
366  libmesh_error();
367  }
368 
369  if (n_dirichlet_conditions !=
370  input.vector_variable_size("dirichlet_condition_boundaries"))
371  {
372  libMesh::out << "Error: " << n_dirichlet_conditions
373  << " Dirichlet condition types does not match "
374  << input.vector_variable_size("dirichlet_condition_boundaries")
375  << " Dirichlet condition boundaries." << std::endl;
376 
377  libmesh_error();
378  }
379 
380  if (n_dirichlet_conditions !=
381  input.vector_variable_size("dirichlet_condition_variables"))
382  {
383  libMesh::out << "Error: " << n_dirichlet_conditions
384  << " Dirichlet condition types does not match "
385  << input.vector_variable_size("dirichlet_condition_variables")
386  << " Dirichlet condition variables sets." << std::endl;
387 
388  libmesh_error();
389  }
390 
391  for (unsigned int dc=0; dc != n_dirichlet_conditions; ++dc)
392  {
393  const std::string func_type =
394  input("dirichlet_condition_types", zero_string, dc);
395 
396  const std::string func_value =
397  input("dirichlet_condition_values", empty_string, dc);
398 
399  const boundary_id_type func_boundary =
400  input("dirichlet_condition_boundaries", boundary_id_type(0), dc);
401 
402  dirichlet_conditions[func_boundary] =
403  (new_function_base(func_type, func_value).release());
404 
405  const std::string variable_set =
406  input("dirichlet_condition_variables", empty_string, dc);
407 
408  for (unsigned int i=0; i != variable_set.size(); ++i)
409  {
410  if (variable_set[i] == '1')
411  dirichlet_condition_variables[func_boundary].push_back(i);
412  else if (variable_set[i] != '0')
413  {
414  libMesh::out << "Unable to understand Dirichlet variable set"
415  << variable_set << std::endl;
416  libmesh_error();
417  }
418  }
419  }
420 
421  GETPOT_REGISTER(neumann_condition_types);
422  GETPOT_REGISTER(neumann_condition_values);
423  GETPOT_REGISTER(neumann_condition_boundaries);
424  GETPOT_REGISTER(neumann_condition_variables);
425 
426  const unsigned int n_neumann_conditions=
427  input.vector_variable_size("neumann_condition_types");
428 
429  if (n_neumann_conditions !=
430  input.vector_variable_size("neumann_condition_values"))
431  {
432  libMesh::out << "Error: " << n_neumann_conditions
433  << " Neumann condition types does not match "
434  << input.vector_variable_size("neumann_condition_values")
435  << " Neumann condition values." << std::endl;
436 
437  libmesh_error();
438  }
439 
440  if (n_neumann_conditions !=
441  input.vector_variable_size("neumann_condition_boundaries"))
442  {
443  libMesh::out << "Error: " << n_neumann_conditions
444  << " Neumann condition types does not match "
445  << input.vector_variable_size("neumann_condition_boundaries")
446  << " Neumann condition boundaries." << std::endl;
447 
448  libmesh_error();
449  }
450 
451  if (n_neumann_conditions !=
452  input.vector_variable_size("neumann_condition_variables"))
453  {
454  libMesh::out << "Error: " << n_neumann_conditions
455  << " Neumann condition types does not match "
456  << input.vector_variable_size("neumann_condition_variables")
457  << " Neumann condition variables sets." << std::endl;
458 
459  libmesh_error();
460  }
461 
462  for (unsigned int nc=0; nc != n_neumann_conditions; ++nc)
463  {
464  const std::string func_type =
465  input("neumann_condition_types", zero_string, nc);
466 
467  const std::string func_value =
468  input("neumann_condition_values", empty_string, nc);
469 
470  const boundary_id_type func_boundary =
471  input("neumann_condition_boundaries", boundary_id_type(0), nc);
472 
473  neumann_conditions[func_boundary] =
474  (new_function_base(func_type, func_value).release());
475 
476  const std::string variable_set =
477  input("neumann_condition_variables", empty_string, nc);
478 
479  for (unsigned int i=0; i != variable_set.size(); ++i)
480  {
481  if (variable_set[i] == '1')
482  neumann_condition_variables[func_boundary].push_back(i);
483  else if (variable_set[i] != '0')
484  {
485  libMesh::out << "Unable to understand Neumann variable set"
486  << variable_set << std::endl;
487  libmesh_error();
488  }
489  }
490  }
491 
492  GETPOT_REGISTER(initial_condition_types);
493  GETPOT_REGISTER(initial_condition_values);
494  GETPOT_REGISTER(initial_condition_subdomains);
495 
496  const unsigned int n_initial_conditions=
497  input.vector_variable_size("initial_condition_types");
498 
499  if (n_initial_conditions !=
500  input.vector_variable_size("initial_condition_values"))
501  {
502  libMesh::out << "Error: " << n_initial_conditions
503  << " initial condition types does not match "
504  << input.vector_variable_size("initial_condition_values")
505  << " initial condition values." << std::endl;
506 
507  libmesh_error();
508  }
509 
510  if (n_initial_conditions !=
511  input.vector_variable_size("initial_condition_subdomains"))
512  {
513  libMesh::out << "Error: " << n_initial_conditions
514  << " initial condition types does not match "
515  << input.vector_variable_size("initial_condition_subdomains")
516  << " initial condition subdomains." << std::endl;
517 
518  libmesh_error();
519  }
520 
521  for (unsigned int i=0; i != n_initial_conditions; ++i)
522  {
523  const std::string func_type =
524  input("initial_condition_types", zero_string, i);
525 
526  const std::string func_value =
527  input("initial_condition_values", empty_string, i);
528 
529  const subdomain_id_type func_subdomain =
530  input("initial_condition_subdomains", subdomain_id_type(0), i);
531 
532  initial_conditions[func_subdomain] =
533  (new_function_base(func_type, func_value).release());
534  }
535 
536  GETPOT_REGISTER(other_interior_function_types);
537  GETPOT_REGISTER(other_interior_function_values);
538  GETPOT_REGISTER(other_interior_function_subdomains);
539  GETPOT_REGISTER(other_interior_function_ids);
540 
541  const unsigned int n_other_interior_functions =
542  input.vector_variable_size("other_interior_function_types");
543 
544  if (n_other_interior_functions !=
545  input.vector_variable_size("other_interior_function_values"))
546  {
547  libMesh::out << "Error: " << n_other_interior_functions
548  << " other interior function types does not match "
549  << input.vector_variable_size("other_interior_function_values")
550  << " other interior function values." << std::endl;
551 
552  libmesh_error();
553  }
554 
555  if (n_other_interior_functions !=
556  input.vector_variable_size("other_interior_function_subdomains"))
557  {
558  libMesh::out << "Error: " << n_other_interior_functions
559  << " other interior function types does not match "
560  << input.vector_variable_size("other_interior_function_subdomains")
561  << " other interior function subdomains." << std::endl;
562 
563  libmesh_error();
564  }
565 
566  if (n_other_interior_functions !=
567  input.vector_variable_size("other_interior_function_ids"))
568  {
569  libMesh::out << "Error: " << n_other_interior_functions
570  << " other interior function types does not match "
571  << input.vector_variable_size("other_interior_function_ids")
572  << " other interior function ids." << std::endl;
573 
574  libmesh_error();
575  }
576 
577  for (unsigned int i=0; i != n_other_interior_functions; ++i)
578  {
579  const std::string func_type =
580  input("other_interior_function_types", zero_string, i);
581 
582  const std::string func_value =
583  input("other_interior_function_values", empty_string, i);
584 
585  const subdomain_id_type func_subdomain =
586  input("other_interior_condition_subdomains", subdomain_id_type(0), i);
587 
588  const int func_id =
589  input("other_interior_condition_ids", int(0), i);
590 
591  other_interior_functions[func_id][func_subdomain] =
592  (new_function_base(func_type, func_value).release());
593  }
594 
595  GETPOT_REGISTER(other_boundary_function_types);
596  GETPOT_REGISTER(other_boundary_function_values);
597  GETPOT_REGISTER(other_boundary_function_boundaries);
598  GETPOT_REGISTER(other_boundary_function_ids);
599 
600  const unsigned int n_other_boundary_functions =
601  input.vector_variable_size("other_boundary_function_types");
602 
603  if (n_other_boundary_functions !=
604  input.vector_variable_size("other_boundary_function_values"))
605  {
606  libMesh::out << "Error: " << n_other_boundary_functions
607  << " other boundary function types does not match "
608  << input.vector_variable_size("other_boundary_function_values")
609  << " other boundary function values." << std::endl;
610 
611  libmesh_error();
612  }
613 
614  if (n_other_boundary_functions !=
615  input.vector_variable_size("other_boundary_function_boundaries"))
616  {
617  libMesh::out << "Error: " << n_other_boundary_functions
618  << " other boundary function types does not match "
619  << input.vector_variable_size("other_boundary_function_boundaries")
620  << " other boundary function boundaries." << std::endl;
621 
622  libmesh_error();
623  }
624 
625  if (n_other_boundary_functions !=
626  input.vector_variable_size("other_boundary_function_ids"))
627  {
628  libMesh::out << "Error: " << n_other_boundary_functions
629  << " other boundary function types does not match "
630  << input.vector_variable_size("other_boundary_function_ids")
631  << " other boundary function ids." << std::endl;
632 
633  libmesh_error();
634  }
635 
636  for (unsigned int i=0; i != n_other_boundary_functions; ++i)
637  {
638  const std::string func_type =
639  input("other_boundary_function_types", zero_string, i);
640 
641  const std::string func_value =
642  input("other_boundary_function_values", empty_string, i);
643 
644  const boundary_id_type func_boundary =
645  input("other_boundary_function_boundaries", boundary_id_type(0), i);
646 
647  const int func_id =
648  input("other_boundary_function_ids", int(0), i);
649 
650  other_boundary_functions[func_id][func_boundary] =
651  (new_function_base(func_type, func_value).release());
652  }
653 
654  GETPOT_INPUT(run_simulation);
655  GETPOT_INPUT(run_postprocess);
656 
657 
658  GETPOT_REGISTER(fe_family);
659  const unsigned int n_fe_family =
660  std::max(1u, input.vector_variable_size("fe_family"));
661  fe_family.resize(n_fe_family, "LAGRANGE");
662  for (unsigned int i=0; i != n_fe_family; ++i)
663  fe_family[i] = input("fe_family", fe_family[i].c_str(), i);
664  GETPOT_REGISTER(fe_order);
665  const unsigned int n_fe_order =
666  input.vector_variable_size("fe_order");
667  fe_order.resize(n_fe_order, 1);
668  for (unsigned int i=0; i != n_fe_order; ++i)
669  fe_order[i] = input("fe_order", (int)fe_order[i], i);
670  GETPOT_INPUT(extra_quadrature_order);
671 
672 
673  GETPOT_INPUT(analytic_jacobians);
674  GETPOT_INPUT(verify_analytic_jacobians);
675  GETPOT_INPUT(numerical_jacobian_h);
676  GETPOT_INPUT(print_solution_norms);
677  GETPOT_INPUT(print_solutions);
678  GETPOT_INPUT(print_residual_norms);
679  GETPOT_INPUT(print_residuals);
680  GETPOT_INPUT(print_jacobian_norms);
681  GETPOT_INPUT(print_jacobians);
682  GETPOT_INPUT(print_element_solutions);
683  GETPOT_INPUT(print_element_residuals);
684  GETPOT_INPUT(print_element_jacobians);
685 
686 
687  GETPOT_INPUT(use_petsc_snes);
688  GETPOT_INPUT(time_solver_quiet);
689  GETPOT_INPUT(solver_quiet);
690  GETPOT_INPUT(solver_verbose);
691  GETPOT_INPUT(reuse_preconditioner);
692  GETPOT_INPUT(require_residual_reduction);
693  GETPOT_INPUT(min_step_length);
694  GETPOT_INT_INPUT(max_linear_iterations);
695  GETPOT_INT_INPUT(max_nonlinear_iterations);
696  GETPOT_INPUT(relative_step_tolerance);
697  GETPOT_INPUT(relative_residual_tolerance);
698  GETPOT_INPUT(absolute_residual_tolerance);
699  GETPOT_INPUT(initial_linear_tolerance);
700  GETPOT_INPUT(minimum_linear_tolerance);
701  GETPOT_INPUT(linear_tolerance_multiplier);
702 
703 
704  GETPOT_INT_INPUT(initial_sobolev_order);
705  GETPOT_INT_INPUT(initial_extra_quadrature);
706  GETPOT_INPUT(refine_uniformly);
707  GETPOT_INPUT(indicator_type);
708  GETPOT_INPUT(patch_reuse);
709  GETPOT_INT_INPUT(sobolev_order);
710 
711  GETPOT_INPUT(system_config_file);
712 
713  std::vector<std::string> bad_variables =
714  input.unidentified_arguments(variable_assignments);
715 
716  // The way unidentified_arguments() works can give us false
717  // positives from repeated (overridden) variable assignments or from
718  // other (e.g. PETSc) command line arguments.
719  std::vector<std::string> actually_bad_variables;
720  for (std::size_t i = 0; i < bad_variables.size(); ++i)
721  {
722  // If any of our ufo arguments start with a -, that's a false
723  // positive from an unrelated command line argument.
724  if (bad_variables[i].empty() || bad_variables[i][0] != '-')
725  {
726  std::string bad_variable_name =
727  bad_variables[i].substr(0, bad_variables[i].find('='));
728  if (!variable_names.count(bad_variable_name))
729  actually_bad_variables.push_back(bad_variables[i]);
730  }
731  // Skip any option variable and (to be safe from false
732  // positives, though it can create false negatives) any
733  // subsequent potential argument
734  else
735  if (bad_variables.size() > (i+1) &&
736  !bad_variables[i+1].empty() &&
737  bad_variables[i+1][0] != '-')
738  ++i;
739  }
740 
741  if (this->comm().rank() == 0 && !actually_bad_variables.empty())
742  {
743  libMesh::err << "ERROR: Unrecognized variables:" << std::endl;
744  for (auto var : actually_bad_variables)
745  libMesh::err << var << std::endl;
746  libMesh::err << "Not found among recognized variables:" << std::endl;
747  for (std::size_t i = 0; i != variable_names.size(); ++i)
748  libMesh::err << variable_assignments[i] << std::endl;
749  libmesh_error();
750  }
751 }
FEMParameters::reuse_preconditioner
bool reuse_preconditioner
Definition: femparameters.h:128
FEMParameters::use_petsc_snes
bool use_petsc_snes
Definition: femparameters.h:127
FEMParameters::print_jacobians
bool print_jacobians
Definition: femparameters.h:118
FEMParameters::timesolver_norm
std::vector< libMesh::FEMNormType > timesolver_norm
Definition: femparameters.h:41
FEMParameters::dirichlet_conditions
std::map< libMesh::boundary_id_type, libMesh::FunctionBase< libMesh::Number > * > dirichlet_conditions
Definition: femparameters.h:90
libMesh::subdomain_id_type
TestClass subdomain_id_type
Based on the 4-byte comment warning above, this probably doesn't work with exodusII at all....
Definition: id_types.h:43
libMesh::FunctionBase< Number >
FEMParameters::domain_edge_length
libMesh::Real domain_edge_length
Definition: femparameters.h:49
FEMParameters::n_timesteps
unsigned int n_timesteps
Definition: femparameters.h:34
FEMParameters::domain_xmin
libMesh::Real domain_xmin
Definition: femparameters.h:48
FEMParameters::system_config_file
std::string system_config_file
Definition: femparameters.h:151
FEMParameters::coarsen_fraction
libMesh::Real coarsen_fraction
Definition: femparameters.h:61
FEMParameters::steadystate_tolerance
libMesh::Real steadystate_tolerance
Definition: femparameters.h:38
FEMParameters::domain_ymin
libMesh::Real domain_ymin
Definition: femparameters.h:48
FEMParameters::max_nonlinear_iterations
unsigned int max_nonlinear_iterations
Definition: femparameters.h:131
FEMParameters::output_xda
bool output_xda
Definition: femparameters.h:68
FEMParameters::print_residuals
bool print_residuals
Definition: femparameters.h:118
FEMParameters::end_time
libMesh::Real end_time
Definition: femparameters.h:38
FEMParameters::print_solution_norms
bool print_solution_norms
Definition: femparameters.h:118
FEMParameters::print_residual_norms
bool print_residual_norms
Definition: femparameters.h:118
FEMParameters::initial_adaptivesteps
unsigned int initial_adaptivesteps
Definition: femparameters.h:63
FEMParameters::timesolver_maxgrowth
libMesh::Real timesolver_maxgrowth
Definition: femparameters.h:38
FEMParameters::dimension
unsigned int dimension
Definition: femparameters.h:45
libMesh::DISCRETE_L2
Definition: enum_norm_type.h:53
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
FEMParameters::coarserefinements
unsigned int coarserefinements
Definition: femparameters.h:51
FEMParameters::timesolver_tolerance
libMesh::Real timesolver_tolerance
Definition: femparameters.h:38
FEMParameters::mesh_partitioner_type
std::string mesh_partitioner_type
Definition: femparameters.h:55
FEMParameters::sobolev_order
unsigned int sobolev_order
Definition: femparameters.h:147
FEMParameters::output_gmv
bool output_gmv
Definition: femparameters.h:68
FEMParameters::refine_fraction
libMesh::Real refine_fraction
Definition: femparameters.h:61
FEMParameters::elementorder
libMesh::Real elementorder
Definition: femparameters.h:47
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
FEMParameters::run_postprocess
bool run_postprocess
Definition: femparameters.h:104
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
FEMParameters::output_exodus
bool output_exodus
Definition: femparameters.h:68
FEMParameters::initial_linear_tolerance
double initial_linear_tolerance
Definition: femparameters.h:134
libMesh::H2
Definition: enum_norm_type.h:38
FEMParameters::fe_order
std::vector< unsigned int > fe_order
Definition: femparameters.h:109
libMesh::RealVectorValue
VectorValue< Real > RealVectorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
Definition: hp_coarsentest.h:46
libMesh::boundary_id_type
int8_t boundary_id_type
Definition: id_types.h:51
FEMParameters::timesolver_core
std::string timesolver_core
Definition: femparameters.h:37
FEMParameters::solver_quiet
bool solver_quiet
Definition: femparameters.h:128
FEMParameters::write_gmv_error
bool write_gmv_error
Definition: femparameters.h:68
FEMParameters::require_residual_reduction
bool require_residual_reduction
Definition: femparameters.h:128
new_function_base
std::unique_ptr< FunctionBase< Number > > new_function_base(const std::string &func_type, const std::string &func_value)
Definition: femparameters.C:140
FEMParameters::coarsen_threshold
libMesh::Real coarsen_threshold
Definition: femparameters.h:61
FEMParameters::run_simulation
bool run_simulation
Definition: femparameters.h:104
FEMParameters::max_linear_iterations
unsigned int max_linear_iterations
Definition: femparameters.h:131
FEMParameters::domain_zmin
libMesh::Real domain_zmin
Definition: femparameters.h:48
libMesh::VectorValue< Real >
FEMParameters::neumann_conditions
std::map< libMesh::boundary_id_type, libMesh::FunctionBase< libMesh::Number > * > neumann_conditions
Definition: femparameters.h:90
FEMParameters::indicator_type
std::string indicator_type
Definition: femparameters.h:145
FEMParameters::refine_uniformly
bool refine_uniformly
Definition: femparameters.h:144
FEMParameters::periodic_boundaries
std::vector< libMesh::PeriodicBoundary > periodic_boundaries
Definition: femparameters.h:84
FEMParameters::print_solutions
bool print_solutions
Definition: femparameters.h:118
FEMParameters::fe_family
std::vector< std::string > fe_family
Definition: femparameters.h:108
FEMParameters::solver_verbose
bool solver_verbose
Definition: femparameters.h:128
FEMParameters::domainfile
std::string domainfile
Definition: femparameters.h:46
FEMParameters::max_adaptivesteps
unsigned int max_adaptivesteps
Definition: femparameters.h:62
FEMParameters::output_bz2
bool output_bz2
Definition: femparameters.h:68
FEMParameters::extrarefinements
unsigned int extrarefinements
Definition: femparameters.h:51
FEMParameters::initial_extra_quadrature
unsigned int initial_extra_quadrature
Definition: femparameters.h:140
FEMParameters::mesh_redistribute_func
std::string mesh_redistribute_func
Definition: femparameters.h:52
FEMParameters::output_xdr
bool output_xdr
Definition: femparameters.h:68
FEMParameters::print_element_jacobians
bool print_element_jacobians
Definition: femparameters.h:118
FEMParameters::coarsegridy
unsigned int coarsegridy
Definition: femparameters.h:50
FEMParameters::absolute_residual_tolerance
libMesh::Real absolute_residual_tolerance
Definition: femparameters.h:132
FEMParameters::output_nemesis
bool output_nemesis
Definition: femparameters.h:68
FEMParameters::relative_residual_tolerance
libMesh::Real relative_residual_tolerance
Definition: femparameters.h:132
FEMParameters::write_exodus_error
bool write_exodus_error
Definition: femparameters.h:68
FEMParameters::relative_step_tolerance
libMesh::Real relative_step_tolerance
Definition: femparameters.h:132
FEMParameters::minimum_linear_tolerance
double minimum_linear_tolerance
Definition: femparameters.h:134
FEMParameters::other_interior_functions
std::map< int, std::map< libMesh::subdomain_id_type, libMesh::FunctionBase< libMesh::Number > * > > other_interior_functions
Definition: femparameters.h:97
FEMParameters::print_element_solutions
bool print_element_solutions
Definition: femparameters.h:118
FEMParameters::domain_edge_height
libMesh::Real domain_edge_height
Definition: femparameters.h:49
FEMParameters::write_interval
unsigned int write_interval
Definition: femparameters.h:67
libMesh::PeriodicBoundary
The definition of a periodic boundary.
Definition: periodic_boundary.h:44
FEMParameters::other_boundary_functions
std::map< int, std::map< libMesh::boundary_id_type, libMesh::FunctionBase< libMesh::Number > * > > other_boundary_functions
Definition: femparameters.h:100
FEMParameters::timesolver_theta
libMesh::Real timesolver_theta
Definition: femparameters.h:38
FEMParameters::global_tolerance
libMesh::Real global_tolerance
Definition: femparameters.h:60
FEMParameters::elementtype
std::string elementtype
Definition: femparameters.h:46
FEMParameters::min_step_length
libMesh::Real min_step_length
Definition: femparameters.h:130
FEMParameters::read
void read(GetPot &input, const std::vector< std::string > *other_variable_names=nullptr)
Definition: femparameters.C:154
FEMParameters::extra_quadrature_order
int extra_quadrature_order
Definition: femparameters.h:110
FEMParameters::patch_reuse
bool patch_reuse
Definition: femparameters.h:146
FEMParameters::coarsegridz
unsigned int coarsegridz
Definition: femparameters.h:50
FEMParameters::~FEMParameters
~FEMParameters()
Definition: femparameters.C:101
FEMParameters::deltat
libMesh::Real deltat
Definition: femparameters.h:38
std
Definition: float128_shims.h:27
libMesh::L2
Definition: enum_norm_type.h:36
FEMParameters::initial_conditions
std::map< libMesh::subdomain_id_type, libMesh::FunctionBase< libMesh::Number > * > initial_conditions
Definition: femparameters.h:88
FEMParameters::numerical_jacobian_h
libMesh::Real numerical_jacobian_h
Definition: femparameters.h:116
FEMParameters::initial_sobolev_order
unsigned int initial_sobolev_order
Definition: femparameters.h:139
FEMParameters::domaintype
std::string domaintype
Definition: femparameters.h:46
FEMParameters::analytic_jacobians
bool analytic_jacobians
Definition: femparameters.h:114
FEMParameters::coarsegridx
unsigned int coarsegridx
Definition: femparameters.h:50
FEMParameters::verify_analytic_jacobians
libMesh::Real verify_analytic_jacobians
Definition: femparameters.h:115
FEMParameters::print_jacobian_norms
bool print_jacobian_norms
Definition: femparameters.h:118
FEMParameters::system_types
std::vector< std::string > system_types
Definition: femparameters.h:77
FEMParameters::FEMParameters
FEMParameters(const libMesh::Parallel::Communicator &comm_in)
FEMParameters::domain_edge_width
libMesh::Real domain_edge_width
Definition: femparameters.h:49
FEMParameters::write_tecplot_error
bool write_tecplot_error
Definition: femparameters.h:68
libMesh::err
OStreamProxy err
libMesh::TestClass
Definition: id_types.h:33
libMesh::H1
Definition: enum_norm_type.h:37
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
FEMParameters::neumann_condition_variables
std::map< libMesh::boundary_id_type, std::vector< unsigned int > > neumann_condition_variables
Definition: femparameters.h:93
libMesh::ParallelObject
An object whose state is distributed along a set of processors.
Definition: parallel_object.h:55
FEMParameters::deltat_reductions
unsigned int deltat_reductions
Definition: femparameters.h:36
FEMParameters::nelem_target
unsigned int nelem_target
Definition: femparameters.h:59
FEMParameters::linear_tolerance_multiplier
double linear_tolerance_multiplier
Definition: femparameters.h:134
FEMParameters::dirichlet_condition_variables
std::map< libMesh::boundary_id_type, std::vector< unsigned int > > dirichlet_condition_variables
Definition: femparameters.h:93
FEMParameters::output_tecplot
bool output_tecplot
Definition: femparameters.h:68
libMesh::out
OStreamProxy out
FEMParameters::output_gz
bool output_gz
Definition: femparameters.h:68
FEMParameters::timesolver_upper_tolerance
libMesh::Real timesolver_upper_tolerance
Definition: femparameters.h:38
FEMParameters::print_element_residuals
bool print_element_residuals
Definition: femparameters.h:118
FEMParameters::time_solver_quiet
bool time_solver_quiet
Definition: femparameters.h:128
libMesh::Quality::name
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
FEMParameters::initial_timestep
unsigned int initial_timestep
Definition: femparameters.h:34