From: bangerth Date: Mon, 19 May 2008 16:23:05 +0000 (+0000) Subject: Now happy with the parameters section. X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=e52c7fa87409d94ddb862fd22a83cdbe7341457e;p=dealii-svn.git Now happy with the parameters section. git-svn-id: https://svn.dealii.org/trunk@16122 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/examples/step-33/step-33.cc b/deal.II/examples/step-33/step-33.cc index 9ad0f071f5..7aaeb15511 100644 --- a/deal.II/examples/step-33/step-33.cc +++ b/deal.II/examples/step-33/step-33.cc @@ -1020,49 +1020,51 @@ namespace Parameters // being specified for each component and // each boundary indicator separately. // - // The data structure used to store the - // boundary indicators is a bit + // The data structure used to store + // the boundary indicators is a bit // complicated. It is an array of - // max_n_boundaries elements - // indicating the range of boundary - // indicators that will be accepted. For - // each entry in this array, we store a - // pair of data: first, an array of size - // n_components that for each - // component of the solution vector - // indicates whether it is an inflow, - // outflow, or other kind of boundary, and - // second a FunctionParser object that - // describes all components of the solution - // vector for this boundary id at once. + // max_n_boundaries + // elements indicating the range of + // boundary indicators that will be + // accepted. For each entry in this + // array, we store a pair of data + // in the + // BoundaryCondition + // structure: first, an array of + // size n_components + // that for each component of the + // solution vector indicates + // whether it is an inflow, + // outflow, or other kind of + // boundary, and second a + // FunctionParser object that + // describes all components of the + // solution vector for this + // boundary id at once. // - // The data structure is made a bit more - // inconvenient by the fact that C++ has no - // way to initialize arrays in - // constructors. Now, we need to tell the - // function parser object at construction - // time how many vector components it is to - // describe. Since this can't be done in - // the constructor of this class, what we - // do is not to have the second part of the - // array elements be a FunctionParser - // object, but a pointer to one, and - // initialize the pointer in the - // parser_parameters() - // function below. In order to avoid - // writing a destructor for this class that - // later releases this memory again, we use - // the boost::shared_ptr class - // instead that will make sure that memory - // is released whenever the object pointed - // to is not used anywhere any more. + // The + // BoundaryCondition + // structure requires a constructor + // since we need to tell the + // function parser object at + // construction time how many + // vector components it is to + // describe. This initialization + // can therefore not wait till we + // actually set the formulas the + // FunctionParser object represents + // later in + // AllParameters::parse_parameters() // - // For the same reason of having to tell - // Function objects their vector size at - // construction time, we have to have a - // constructor of this class that at least - // initializes the other FunctionParser - // object, i.e. the one describing initial + // For the same reason of having to + // tell Function objects their + // vector size at construction + // time, we have to have a + // constructor of the + // AllParameters class + // that at least initializes the + // other FunctionParser object, + // i.e. the one describing initial // conditions. template struct AllParameters : public Solver, @@ -1080,6 +1082,15 @@ namespace Parameters pressure_boundary }; + struct BoundaryConditions + { + BoundaryKind kind[EulerEquations::n_components]; + FunctionParser values; + + BoundaryConditions (); + }; + + AllParameters (); double diffusion_power; @@ -1091,10 +1102,7 @@ namespace Parameters std::string mesh_filename; FunctionParser initial_conditions; - - std::pair::n_components], - boost::shared_ptr > > - boundary_conditions[max_n_boundaries]; + BoundaryConditions boundary_conditions[max_n_boundaries]; static void declare_parameters (ParameterHandler &prm); void parse_parameters (ParameterHandler &prm); @@ -1102,6 +1110,13 @@ namespace Parameters + template + AllParameters::BoundaryConditions::BoundaryConditions () + : + values (EulerEquations::n_components) + {} + + template AllParameters::AllParameters () : @@ -1144,7 +1159,7 @@ namespace Parameters { prm.declare_entry("no penetration", "false", Patterns::Bool(), - "Whether the names boundary is allows gas to " + "whether the named boundary allows gas to " "penetrate or is a rigid wall"); for (unsigned int di=0; di::n_components; ++di) @@ -1210,54 +1225,39 @@ namespace Parameters prm.enter_subsection("boundary_" + Utilities::int_to_string(boundary_id)); { - boundary_conditions[boundary_id].second - = boost::shared_ptr > ( - new FunctionParser(EulerEquations::n_components)); - std::vector expressions(EulerEquations::n_components, "0.0"); - const bool nopen = prm.get_bool("no penetration"); + const bool no_penetration = prm.get_bool("no penetration"); - // Determine how each component is - // handled. for (unsigned int di=0; di::n_components; ++di) { const std::string boundary_type = prm.get("w_" + Utilities::int_to_string(di)); - const std::string var_value - = prm.get("w_" + Utilities::int_to_string(di) + - " value"); - if (di < dim && nopen) - boundary_conditions[boundary_id].first[di] = no_penetration_boundary; + if ((di < dim) && (no_penetration == true)) + boundary_conditions[boundary_id].kind[di] = no_penetration_boundary; else if (boundary_type == "inflow") - { - boundary_conditions[boundary_id].first[di] = inflow_boundary; - expressions[di] = var_value; - } + boundary_conditions[boundary_id].kind[di] = inflow_boundary; else if (boundary_type == "pressure") - { - boundary_conditions[boundary_id].first[di] = pressure_boundary; - expressions[di] = var_value; - } + boundary_conditions[boundary_id].kind[di] = pressure_boundary; else if (boundary_type == "outflow") - boundary_conditions[boundary_id].first[di] = outflow_boundary; + boundary_conditions[boundary_id].kind[di] = outflow_boundary; else AssertThrow (false, ExcNotImplemented()); - } + expressions[di] = prm.get("w_" + Utilities::int_to_string(di) + + " value"); + } - // Add the boundary condition to the - // law. - boundary_conditions[boundary_id].second->initialize (FunctionParser::default_variable_names(), - expressions, - std::map()); + boundary_conditions[boundary_id].values + .initialize (FunctionParser::default_variable_names(), + expressions, + std::map()); } prm.leave_subsection(); } - // Initial conditions. prm.enter_subsection("initial condition"); { std::vector expressions (EulerEquations::n_components, @@ -1275,9 +1275,7 @@ namespace Parameters Parameters::Refinement::parse_parameters (prm); Parameters::Flux::parse_parameters (prm); Parameters::Output::parse_parameters (prm); - } - - + } } @@ -1683,7 +1681,7 @@ void ConsLaw::assemble_face_term( // not prescribed, the values evaluate to // zero and are ignored, below. std::vector > bvals(n_q_points, Vector(EulerEquations::n_components)); - parameters.boundary_conditions[boundary_id].second->vector_value_list(fe_v.get_quadrature_points(), bvals); + parameters.boundary_conditions[boundary_id].values.vector_value_list(fe_v.get_quadrature_points(), bvals); // We loop the quadrature points, and we treat each // component individualy. @@ -1691,9 +1689,9 @@ void ConsLaw::assemble_face_term( for (unsigned int di = 0; di < EulerEquations::n_components; di++) { // An inflow/dirichlet type of boundary condition - if (parameters.boundary_conditions[boundary_id].first[di] == Parameters::AllParameters::inflow_boundary) { + if (parameters.boundary_conditions[boundary_id].kind[di] == Parameters::AllParameters::inflow_boundary) { Wminus[q][di] = bvals[q](di); - } else if (parameters.boundary_conditions[boundary_id].first[di] == Parameters::AllParameters::pressure_boundary) { + } else if (parameters.boundary_conditions[boundary_id].kind[di] == Parameters::AllParameters::pressure_boundary) { // A prescribed pressure boundary // condition. This boundary // condition is complicated by the @@ -1711,11 +1709,11 @@ void ConsLaw::assemble_face_term( Sacado::Fad::DFad rho_vel_sqr = 0; Sacado::Fad::DFad dens; - dens = parameters.boundary_conditions[boundary_id].first[EulerEquations::density_component] == Parameters::AllParameters::inflow_boundary ? bvals[q](EulerEquations::density_component) : + dens = parameters.boundary_conditions[boundary_id].kind[EulerEquations::density_component] == Parameters::AllParameters::inflow_boundary ? bvals[q](EulerEquations::density_component) : Wplus[q][EulerEquations::density_component]; for (unsigned int d=0; d < dim; d++) { - if (parameters.boundary_conditions[boundary_id].first[d] == Parameters::AllParameters::inflow_boundary) + if (parameters.boundary_conditions[boundary_id].kind[d] == Parameters::AllParameters::inflow_boundary) rho_vel_sqr += bvals[q](d)*bvals[q](d); else rho_vel_sqr += Wplus[q][d]*Wplus[q][d]; @@ -1726,7 +1724,7 @@ void ConsLaw::assemble_face_term( Wminus[q][di] = bvals[q](di)/(EulerEquations::gas_gamma-1.0) + 0.5*rho_vel_sqr; - } else if (parameters.boundary_conditions[boundary_id].first[di] == Parameters::AllParameters::outflow_boundary) { + } else if (parameters.boundary_conditions[boundary_id].kind[di] == Parameters::AllParameters::outflow_boundary) { // A free/outflow boundary, very simple. Wminus[q][di] = Wplus[q][di];