-
-
Notifications
You must be signed in to change notification settings - Fork 13
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Unify parameter interface for solid solver #51
Merged
Merged
Changes from 7 commits
Commits
Show all changes
11 commits
Select commit
Hold shift + click to select a range
bf653bd
Move and unify parameter header file into include path
davidscn c3fe006
Move and unify parameter sources into include path
davidscn 4057bbe
Configure CMake for new source file and adjust includes
davidscn 8f8369c
Simplify parameter handling using add_parameter
davidscn bb7ed54
Remove specific parameter files and add unified one
davidscn d691fa2
Add parameter option for model selection
davidscn 8997e9f
Add an entry for the parameter file location
davidscn 5476c49
Fix Hook Hooke typo
davidscn 81363a6
Apply suggestions from code review
davidscn 005535a
Fix 'Force' data check for non-linear model
davidscn 54893cc
Enfore compatibility to tutorial parameter file
davidscn File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,198 @@ | ||
#include <adapter/parameters.h> | ||
|
||
namespace Parameters | ||
{ | ||
void | ||
Time::add_output_parameters(ParameterHandler &prm) | ||
{ | ||
prm.enter_subsection("Time"); | ||
{ | ||
prm.add_parameter("End time", end_time, "End time", Patterns::Double()); | ||
|
||
prm.add_parameter("Time step size", | ||
delta_t, | ||
"Time step size", | ||
Patterns::Double()); | ||
|
||
prm.add_parameter("Output interval", | ||
output_interval, | ||
"Write results every x timesteps", | ||
Patterns::Integer(0)); | ||
} | ||
prm.leave_subsection(); | ||
} | ||
|
||
|
||
void | ||
System::add_output_parameters(ParameterHandler &prm) | ||
{ | ||
prm.enter_subsection("System properties"); | ||
{ | ||
prm.add_parameter("Shear modulus", | ||
mu, | ||
"Shear modulus", | ||
Patterns::Double()); | ||
|
||
prm.add_parameter("Poisson's ratio", | ||
nu, | ||
"Poisson's ratio", | ||
Patterns::Double(-1.0, 0.5)); | ||
|
||
prm.add_parameter("rho", rho, "density", Patterns::Double(0.0)); | ||
|
||
prm.add_parameter("body forces", | ||
body_force, | ||
"body forces x,y,z", | ||
Patterns::List(Patterns::Double())); | ||
} | ||
prm.leave_subsection(); | ||
|
||
lambda = 2 * mu * nu / (1 - 2 * nu); | ||
} | ||
|
||
|
||
void | ||
Solver::add_output_parameters(ParameterHandler &prm) | ||
{ | ||
prm.enter_subsection("Solver"); | ||
{ | ||
prm.add_parameter("Model", | ||
model, | ||
"Structural model to be used: linear or neo-Hook", | ||
Patterns::Selection("linear|neo-Hook")); | ||
|
||
prm.add_parameter("Solver type", | ||
type_lin, | ||
"Linear solver: CG or Direct", | ||
Patterns::Selection("CG|Direct")); | ||
|
||
prm.add_parameter("Residual", | ||
tol_lin, | ||
"Linear solver residual (scaled by residual norm)", | ||
Patterns::Double(0.0)); | ||
|
||
prm.add_parameter( | ||
"Max iteration multiplier", | ||
max_iterations_lin, | ||
"Linear solver iterations (multiples of the system matrix size)", | ||
Patterns::Double(0.0)); | ||
|
||
prm.add_parameter("Max iterations Newton-Raphson", | ||
max_iterations_NR, | ||
"Number of Newton-Raphson iterations allowed", | ||
Patterns::Integer(0)); | ||
|
||
prm.add_parameter("Tolerance force", | ||
tol_f, | ||
"Force residual tolerance", | ||
Patterns::Double(0.0)); | ||
|
||
prm.add_parameter("Tolerance displacement", | ||
tol_u, | ||
"Displacement error tolerance", | ||
Patterns::Double(0.0)); | ||
} | ||
prm.leave_subsection(); | ||
} | ||
|
||
|
||
void | ||
Discretization::add_output_parameters(ParameterHandler &prm) | ||
{ | ||
prm.enter_subsection("Discretization"); | ||
{ | ||
prm.add_parameter("Polynomial degree", | ||
poly_degree, | ||
"Polynomial degree of the FE system", | ||
Patterns::Integer(0)); | ||
|
||
prm.add_parameter("theta", | ||
theta, | ||
"Time integration scheme", | ||
Patterns::Double(0, 1)); | ||
|
||
prm.add_parameter("beta", beta, "Newmark beta", Patterns::Double(0, 0.5)); | ||
|
||
prm.add_parameter("gamma", | ||
gamma, | ||
"Newmark gamma", | ||
Patterns::Double(0, 1)); | ||
} | ||
prm.leave_subsection(); | ||
} | ||
|
||
|
||
void | ||
PreciceAdapterConfiguration::add_output_parameters(ParameterHandler &prm) | ||
{ | ||
prm.enter_subsection("precice configuration"); | ||
{ | ||
prm.add_parameter("Scenario", | ||
scenario, | ||
"Cases: FSI3 or PF for perpendicular flap", | ||
Patterns::Selection("FSI3|PF")); | ||
|
||
prm.add_parameter("precice config-file", | ||
config_file, | ||
"Name of the precice configuration file", | ||
Patterns::Anything()); | ||
|
||
prm.add_parameter( | ||
"Participant name", | ||
participant_name, | ||
"Name of the participant in the precice-config.xml file", | ||
Patterns::Anything()); | ||
|
||
prm.add_parameter( | ||
"Mesh name", | ||
mesh_name, | ||
"Name of the coupling mesh in the precice-config.xml file", | ||
Patterns::Anything()); | ||
|
||
prm.add_parameter("Read data name", | ||
read_data_name, | ||
"Name of the read data in the precice-config.xml file", | ||
Patterns::Anything()); | ||
|
||
prm.add_parameter("Write data name", | ||
write_data_name, | ||
"Name of the write data in the precice-config.xml file", | ||
Patterns::Anything()); | ||
|
||
prm.add_parameter("Flap location", | ||
flap_location, | ||
"PF x-location", | ||
Patterns::Double(-3, 3)); | ||
} | ||
prm.leave_subsection(); | ||
|
||
// Look at the specific type of read data | ||
if ((read_data_name.find("Stress") == 0)) | ||
data_consistent = true; | ||
else if ((read_data_name.find("Force") == 0)) | ||
data_consistent = false; | ||
else | ||
AssertThrow( | ||
false, | ||
ExcMessage( | ||
"Unknown read data type. Please use 'Force' or 'Stress' in the read data naming.")); | ||
} | ||
|
||
|
||
AllParameters::AllParameters(const std::string &input_file) | ||
{ | ||
ParameterHandler prm; | ||
|
||
Solver::add_output_parameters(prm); | ||
Discretization::add_output_parameters(prm); | ||
System::add_output_parameters(prm); | ||
Time::add_output_parameters(prm); | ||
PreciceAdapterConfiguration::add_output_parameters(prm); | ||
|
||
prm.parse_input(input_file); | ||
|
||
// Optional, if we want to print all parameters in the beginning of the | ||
// simulation | ||
// prm.print_parameters(std::cout,ParameterHandler::Text); | ||
} | ||
} // namespace Parameters |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,112 @@ | ||
#ifndef PARAMETERS_H | ||
#define PARAMETERS_H | ||
|
||
#include <deal.II/base/parameter_handler.h> | ||
|
||
/** | ||
* This class declares all parameters, which can be specified in the parameter | ||
* file. The subsection abut preCICE configurations is directly interlinked | ||
* to the Adapter class. | ||
*/ | ||
namespace Parameters | ||
{ | ||
using namespace dealii; | ||
/** | ||
* @brief Time: Specifies simulation time properties | ||
*/ | ||
struct Time | ||
{ | ||
double end_time = 1; | ||
double delta_t = 0.1; | ||
int output_interval = 1; | ||
|
||
void | ||
add_output_parameters(ParameterHandler &prm); | ||
}; | ||
|
||
/** | ||
* @brief The System struct keeps material properties and body force contributions | ||
*/ | ||
struct System | ||
{ | ||
double nu = 0.3; | ||
double mu = 1538462; | ||
double lambda = 2307692; | ||
double rho = 1000; | ||
Tensor<1, 3, double> body_force; | ||
|
||
void | ||
add_output_parameters(ParameterHandler &prm); | ||
}; | ||
|
||
|
||
/** | ||
* @brief LinearSolver: Specifies linear solver properties | ||
*/ | ||
struct Solver | ||
{ | ||
std::string model = "linear"; | ||
std::string type_lin = "Direct"; | ||
double tol_lin = 1e-6; | ||
double max_iterations_lin = 1; | ||
unsigned int max_iterations_NR = 10; | ||
double tol_f = 1e-9; | ||
double tol_u = 1e-6; | ||
|
||
void | ||
add_output_parameters(ParameterHandler &prm); | ||
}; | ||
|
||
|
||
|
||
/** | ||
* @brief Discretization: Specifies parameters for time integration by a | ||
* theta scheme and polynomial degree of the FE system | ||
*/ | ||
struct Discretization | ||
{ | ||
unsigned int poly_degree = 3; | ||
// For the linear elastic model (theta-scheme) | ||
double theta = 0.5; | ||
// For the nonlinear elastic model (Newmark) | ||
double beta = 0.25; | ||
double gamma = 0.5; | ||
|
||
void | ||
add_output_parameters(ParameterHandler &prm); | ||
}; | ||
|
||
|
||
/** | ||
* @brief PreciceAdapterConfiguration: Specifies preCICE related information. | ||
* A lot of these information need to be consistent with the | ||
* precice-config.xml file. | ||
*/ | ||
struct PreciceAdapterConfiguration | ||
{ | ||
std::string scenario = "FSI3"; | ||
std::string config_file = "precice-config.xml"; | ||
std::string participant_name = "dealiisolver"; | ||
std::string mesh_name = "dealii-mesh"; | ||
std::string read_data_name = "Stress"; | ||
std::string write_data_name = "Displacement"; | ||
double flap_location = 0.0; | ||
bool data_consistent = true; | ||
|
||
void | ||
add_output_parameters(ParameterHandler &prm); | ||
}; | ||
|
||
|
||
struct AllParameters : public Solver, | ||
public Discretization, | ||
public System, | ||
public Time, | ||
public PreciceAdapterConfiguration | ||
|
||
{ | ||
AllParameters(const std::string &input_file); | ||
}; | ||
} // namespace Parameters | ||
|
||
#endif // PARAMETERS_H |
Oops, something went wrong.
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Done. I added now an option here, where the user can select between the linear model and the neo-Hookean model. I hope this is intuitive enough, otherwise let me know.