4 #include <deal.II/base/utilities.h>     6 #include "parameters_grid_refinement.h"    10 namespace Parameters {
    16         prm.declare_entry(
"refinement_steps", 
"3",
    17                           dealii::Patterns::Integer(),
    18                           "Number of iterations to be performed");
    20         prm.declare_entry(
"refinement_method", 
"uniform",
    21                            dealii::Patterns::Selection(
    25                           "Enum of refinement methods."    31         prm.declare_entry(
"refinement_type", 
"h",
    32                           dealii::Patterns::Selection(
    36                           "Enum of refinement types."    42         prm.declare_entry(
"anisotropic", 
"false",
    43                           dealii::Patterns::Bool(),
    44                           "Inidcates whether the refinement should be done anisotropically.");
    46         prm.declare_entry(
"anisotropic_ratio_max", 
"1.0e+12",
    47                           dealii::Patterns::Double(1.0, dealii::Patterns::Double::max_double_value),
    48                           "maximum anisotropic ratio for continuous size field targets.");
    50         prm.declare_entry(
"anisotropic_ratio_min", 
"1.0e-12",
    51                           dealii::Patterns::Double(dealii::Patterns::Double::min_double_value, 1.0),
    52                           "miniumum anistropic ratio for continuous size field targets.");
    54         prm.declare_entry(
"anisotropic_threshold_ratio", 
"3.0",
    55                           dealii::Patterns::Double(1.0, dealii::Patterns::Double::max_double_value),
    56                           "Threshold for flagging cells with anisotropic refinement.");
    58         prm.declare_entry(
"anisotropic_indicator", 
"jump_based",
    59                           dealii::Patterns::Selection(
    61                           " reconstruction_based"),
    62                           "Enum of anisotropic indicators (unused for isotropic refinement)."    65                           "  reconstruction_based>.");
    67         prm.declare_entry(
"error_indicator", 
"error_based",
    68                           dealii::Patterns::Selection(
    73                           "Enum of error indicators (unused for uniform refinement)."    80         prm.declare_entry(
"output_type", 
"msh_out",
    81                           dealii::Patterns::Selection(
    84                           "Enum of output data types (for interface with mesh generators)."    89         prm.declare_entry(
"output_data_type", 
"size_field",
    90                           dealii::Patterns::Selection(
    94                           "Enum of output data types of refinement indicator (used for msh_out only currently)."   100         prm.declare_entry(
"norm_Lq", 
"2.0",
   101                           dealii::Patterns::Double(1.0, dealii::Patterns::Double::max_double_value),
   102                           "Degree of q for use in the Lq norm of some indicators.");
   104         prm.declare_entry(
"refinement_fraction", 
"0.3",
   105                           dealii::Patterns::Double(0.0, 1.0),
   106                           "Fraction of elements to undergo refinement for fixed_fraction method.");
   108         prm.declare_entry(
"coarsening_fraction", 
"0.03",
   109                           dealii::Patterns::Double(0.0, 1.0),
   110                           "Fraction of elements to undergo coarsening for fixed_fraction method.");
   112         prm.declare_entry(
"r_max", 
"20",
   113                           dealii::Patterns::Double(1.0, dealii::Patterns::Double::max_double_value),
   114                           "Maximum refinement factor for adjoint-based size-field (from log DWR).");
   116         prm.declare_entry(
"c_max", 
"4",
   117                           dealii::Patterns::Double(1.0, dealii::Patterns::Double::max_double_value),
   118                           "Maximum coarsening factor for adjoint-based size-field (from log DWR).");
   120         prm.declare_entry(
"complexity_scale", 
"2.0",
   121                           dealii::Patterns::Double(),
   122                           "Scaling factor multiplying previous complexity.");
   124         prm.declare_entry(
"complexity_add", 
"0.0",
   125                           dealii::Patterns::Double(),
   126                           "Constant added to the complexity at each step.");
   128         prm.declare_entry(
"complexity_vector", 
"[]",
   129                           dealii::Patterns::Anything(),
   130                           "Stores an initial vector of values for complexity targets. "   131                           "Will iterate over and then switch to scaling and adding. "   132                           "Formatted in square brackets and seperated by commas, eg. \"[1000,2000]\"");
   134         prm.declare_entry(
"exit_after_refine", 
"false",
   135                           dealii::Patterns::Bool(),
   136                           "Option to exit after call to the grid refinement (for debugging mesh write).");
   148         const std::string refinement_method_string = prm.get(
"refinement_method");
   149         if(refinement_method_string == 
"uniform")            {
refinement_method = RefinementMethod::uniform;}
   150         else if(refinement_method_string == 
"fixed_fraction"){
refinement_method = RefinementMethod::fixed_fraction;}
   151         else if(refinement_method_string == 
"continuous")    {
refinement_method = RefinementMethod::continuous;}
   153         const std::string refinement_type_string = prm.get(
"refinement_type");
   154         if(refinement_type_string == 
"h")      {
refinement_type = RefinementType::h;}
   155         else if(refinement_type_string == 
"p") {
refinement_type = RefinementType::p;}
   156         else if(refinement_type_string == 
"hp"){
refinement_type = RefinementType::hp;}
   165         const std::string anisotropic_indicator_string = prm.get(
"anisotropic_indicator");
   166         if(anisotropic_indicator_string == 
"jump_based")               {
anisotropic_indicator = AnisoIndicator::jump_based;}
   167         else if(anisotropic_indicator_string == 
"reconstruction_based"){
anisotropic_indicator = AnisoIndicator::reconstruction_based;}
   169         const std::string error_indicator_string = prm.get(
"error_indicator");
   170         if(error_indicator_string == 
"error_based")        {
error_indicator = ErrorIndicator::error_based;}
   171         else if(error_indicator_string == 
"hessian_based") {
error_indicator = ErrorIndicator::hessian_based;}
   172         else if(error_indicator_string == 
"residual_based"){
error_indicator = ErrorIndicator::residual_based;}
   173         else if(error_indicator_string == 
"adjoint_based") {
error_indicator = ErrorIndicator::adjoint_based;}
   175         const std::string output_type_string = prm.get(
"output_type");
   176         if(output_type_string == 
"gmsh_out")     {
output_type = OutputType::gmsh_out;}
   177         else if(output_type_string == 
"msh_out") {
output_type = OutputType::msh_out;}
   179         const std::string output_data_type_string = prm.get(
"output_data_type");
   180         if(output_data_type_string == 
"size_field")        {
output_data_type = OutputDataType::size_field;}
   181         else if(output_data_type_string == 
"frame_field")  {
output_data_type = OutputDataType::frame_field;}
   182         else if(output_data_type_string == 
"metric_field") {
output_data_type = OutputDataType::metric_field;}
   184         norm_Lq             = prm.get_double(
"norm_Lq");
   188         r_max = prm.get_double(
"r_max");
   189         c_max = prm.get_double(
"c_max");
   194         std::string complexity_string = prm.get(
"complexity_vector");
   197         std::string removeChars = 
"[]";
   198         for(
unsigned int i = 0; i < removeChars.length(); ++i)
   199             complexity_string.erase(std::remove(complexity_string.begin(), complexity_string.end(), removeChars.at(i)), complexity_string.end());
   200         std::vector<std::string> complexity_string_vector = dealii::Utilities::split_string_list(complexity_string);
   203         for(
auto entry : complexity_string_vector)
 RefinementType refinement_type
Selected type of refinement to be performed. 
double r_max
refinement factor for log DWR size field 
Files for the baseline physics. 
AnisoIndicator anisotropic_indicator
Selected anisotropic splitting indicator. 
double coarsening_fraction
coarsening fraction for fixed-fraction methods 
ErrorIndicator error_indicator
Selected error indicator type. 
bool anisotropic
Flag for performing anisotropic refinement. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
double anisotropic_threshold_ratio
threshold value in anisotropic indicator to enable anisotropic splitting 
RefinementMethod refinement_method
Selected method of refinement. 
double norm_Lq
Lq norm exponent selection. 
unsigned int refinement_steps
Number of refinement steps to be performed. 
double complexity_add
additive constant to complexity between grid refinement iterations 
bool exit_after_refine
Flag to exit after call to refinement. 
OutputDataType output_data_type
Selected data storage type. 
double anisotropic_ratio_max
Maximum anisotropic ratio for continuous size field targets. 
double refinement_fraction
refinement fraction for fixed-fraction methods 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
double anisotropic_ratio_min
Minimum anisotropic ratio for continuous zie field targets. 
std::vector< double > complexity_vector
Vector of complexities to be used for initial continuous grid refinement iterations. 
double complexity_scale
multiplier to complexity between grid refinement iterations 
OutputType output_type
Selected file output type. 
double c_max
coarsening factor for log DWR size field