Skip to content
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

Add a parameter object on the system, use it in derived classes #4007

Open
wants to merge 5 commits into
base: devel
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
30 changes: 15 additions & 15 deletions examples/miscellaneous/miscellaneous_ex7/biharmonic_jr.C
Original file line number Diff line number Diff line change
Expand Up @@ -99,18 +99,18 @@ void Biharmonic::JR::initialize()
if (_biharmonic._verbose)
libMesh::out << ">>> Initializing Biharmonic::JR\n";

Parameters parameters;
parameters.set<Point>("center") = _biharmonic._initialCenter;
parameters.set<Real>("width") = _biharmonic._initialWidth;
Parameters params;
params.set<Point>("center") = _biharmonic._initialCenter;
params.set<Real>("width") = _biharmonic._initialWidth;

if (_biharmonic._initialState == Biharmonic::BALL)
project_solution(Biharmonic::JR::InitialDensityBall, Biharmonic::JR::InitialGradientZero, parameters);
project_solution(Biharmonic::JR::InitialDensityBall, Biharmonic::JR::InitialGradientZero, params);

if (_biharmonic._initialState == Biharmonic::ROD)
project_solution(Biharmonic::JR::InitialDensityRod, Biharmonic::JR::InitialGradientZero, parameters);
project_solution(Biharmonic::JR::InitialDensityRod, Biharmonic::JR::InitialGradientZero, params);

if (_biharmonic._initialState == Biharmonic::STRIP)
project_solution(Biharmonic::JR::InitialDensityStrip, Biharmonic::JR::InitialGradientZero, parameters);
project_solution(Biharmonic::JR::InitialDensityStrip, Biharmonic::JR::InitialGradientZero, params);

// both states are equal
*(old_local_solution) = *(current_local_solution);
Expand All @@ -125,13 +125,13 @@ void Biharmonic::JR::initialize()


Number Biharmonic::JR::InitialDensityBall(const Point & p,
const Parameters & parameters,
const Parameters & params,
const std::string &,
const std::string &)
{
// Initialize with a ball in the middle, which is a segment in 1D, a disk in 2D and a ball in 3D.
Point center = parameters.get<Point>("center");
Real width = parameters.get<Real>("width");
Point center = params.get<Point>("center");
Real width = params.get<Real>("width");
Point pc = p-center;
Real r = pc.norm();
return (r < width) ? 1.0 : -0.5;
Expand All @@ -141,13 +141,13 @@ Number Biharmonic::JR::InitialDensityBall(const Point & p,


Number Biharmonic::JR::InitialDensityRod(const Point & p,
const Parameters & parameters,
const Parameters & params,
const std::string &,
const std::string &)
{
// Initialize with a rod in the middle so that we have a z-homogeneous system to model the 2D disk.
Point center = parameters.get<Point>("center");
Real width = parameters.get<Real>("width");
Point center = params.get<Point>("center");
Real width = params.get<Real>("width");
Point pc = p-center;
#if LIBMESH_DIM > 2
pc(2) = 0;
Expand All @@ -161,13 +161,13 @@ Number Biharmonic::JR::InitialDensityRod(const Point & p,


Number Biharmonic::JR::InitialDensityStrip(const Point & p,
const Parameters & parameters,
const Parameters & params,
const std::string &,
const std::string &)
{
// Initialize with a wide strip in the middle so that we have a yz-homogeneous system to model the 1D.
Point center = parameters.get<Point>("center");
Real width = parameters.get<Real>("width");
Point center = params.get<Point>("center");
Real width = params.get<Real>("width");
Real r = sqrt((p(0)-center(0))*(p(0)-center(0)));
return (r < width) ? 1.0 : -0.5;
}
Expand Down
6 changes: 3 additions & 3 deletions examples/miscellaneous/miscellaneous_ex7/biharmonic_jr.h
Original file line number Diff line number Diff line change
Expand Up @@ -46,17 +46,17 @@ class Biharmonic::JR : public TransientNonlinearImplicitSystem,
* Static functions to be used for initialization
*/
static Number InitialDensityBall(const Point & p,
const Parameters & parameters,
const Parameters & params,
const std::string &,
const std::string &);

static Number InitialDensityRod(const Point & p,
const Parameters & parameters,
const Parameters & params,
const std::string &,
const std::string &);

static Number InitialDensityStrip(const Point & p,
const Parameters & parameters,
const Parameters & params,
const std::string &,
const std::string &);

Expand Down
17 changes: 4 additions & 13 deletions include/systems/system.h
Original file line number Diff line number Diff line change
Expand Up @@ -27,6 +27,7 @@
#include "libmesh/fem_function_base.h"
#include "libmesh/libmesh_common.h"
#include "libmesh/parallel_object.h"
#include "libmesh/parameters.h"
#include "libmesh/qoi_set.h"
#include "libmesh/reference_counted_object.h"
#include "libmesh/tensor_value.h" // For point_hessian
Expand Down Expand Up @@ -60,7 +61,6 @@ class MeshBase;
class Xdr;
class DofMap;
template <typename Output> class FunctionBase;
class Parameters;
class ParameterVector;
class Point;
class SensitivityData;
Expand Down Expand Up @@ -510,8 +510,6 @@ class System : public ReferenceCountedObject<System>,
* user-provided cloneable functors.
* A gradient \p g is only required/used for projecting onto finite
* element spaces with continuous derivatives.
* If non-default \p Parameters are to be used, they can be provided
* in the \p parameters argument.
*/
void project_solution (FunctionBase<Number> * f,
FunctionBase<Gradient> * g = nullptr) const;
Expand All @@ -522,8 +520,6 @@ class System : public ReferenceCountedObject<System>,
* user-provided cloneable functors.
* A gradient \p g is only required/used for projecting onto finite
* element spaces with continuous derivatives.
* If non-default \p Parameters are to be used, they can be provided
* in the \p parameters argument.
*/
void project_solution (FEMFunctionBase<Number> * f,
FEMFunctionBase<Gradient> * g = nullptr) const;
Expand Down Expand Up @@ -554,8 +550,6 @@ class System : public ReferenceCountedObject<System>,
* user-provided cloneable functors.
* A gradient \p g is only required/used for projecting onto finite
* element spaces with continuous derivatives.
* If non-default \p Parameters are to be used, they can be provided
* in the \p parameters argument.
*
* Constrain the new vector using the requested adjoint rather than
* primal constraints if is_adjoint is non-negative.
Expand All @@ -572,8 +566,6 @@ class System : public ReferenceCountedObject<System>,
* user-provided cloneable functors.
* A gradient \p g is only required/used for projecting onto finite
* element spaces with continuous derivatives.
* If non-default \p Parameters are to be used, they can be provided
* in the \p parameters argument.
*
* Constrain the new vector using the requested adjoint rather than
* primal constraints if is_adjoint is non-negative.
Expand Down Expand Up @@ -611,8 +603,6 @@ class System : public ReferenceCountedObject<System>,
* user-provided cloneable functors.
* A gradient \p g is only required/used for projecting onto finite
* element spaces with continuous derivatives.
* If non-default \p Parameters are to be used, they can be provided
* in the \p parameters argument.
*/
void boundary_project_solution (const std::set<boundary_id_type> & b,
const std::vector<unsigned int> & variables,
Expand Down Expand Up @@ -648,8 +638,6 @@ class System : public ReferenceCountedObject<System>,
* user-provided cloneable functors.
* A gradient \p g is only required/used for projecting onto finite
* element spaces with continuous derivatives.
* If non-default \p Parameters are to be used, they can be provided
* in the \p parameters argument.
*
* Constrain the new vector using the requested adjoint rather than
* primal constraints if is_adjoint is non-negative.
Expand Down Expand Up @@ -1531,6 +1519,9 @@ class System : public ReferenceCountedObject<System>,
*/
virtual void prolong_vectors ();

/// Parameters for the system. If a parameter is not provided, it should be retrieved from the EquationSystems
Parameters parameters;

/**
* Flag which tells the system to whether or not to
* call the user assembly function during each call to solve().
Expand Down
12 changes: 8 additions & 4 deletions src/systems/eigen_system.C
Original file line number Diff line number Diff line change
Expand Up @@ -222,16 +222,20 @@ EigenSystem::solve_helper(SparseMatrix<Number> * const A,
// Get the tolerance for the solver and the maximum
// number of iterations. Here, we simply adopt the linear solver
// specific parameters.
const double tol =
const double tol = parameters.have_parameter<Real>("linear solver tolerance") ?
double(parameters.get<Real>("linear solver tolerance")) :
double(es.parameters.get<Real>("linear solver tolerance"));

const unsigned int maxits =
const unsigned int maxits = parameters.have_parameter<unsigned int>("linear solver maximum iterations") ?
parameters.get<unsigned int>("linear solver maximum iterations") :
es.parameters.get<unsigned int>("linear solver maximum iterations");

const unsigned int nev =
const unsigned int nev = parameters.have_parameter<unsigned int>("eigenpairs") ?
parameters.get<unsigned int>("eigenpairs") :
es.parameters.get<unsigned int>("eigenpairs");

const unsigned int ncv =
const unsigned int ncv = parameters.have_parameter<unsigned int>("basis vectors") ?
parameters.get<unsigned int>("basis vectors") :
es.parameters.get<unsigned int>("basis vectors");

std::pair<unsigned int, unsigned int> solve_data;
Expand Down
9 changes: 5 additions & 4 deletions src/systems/frequency_system.C
Original file line number Diff line number Diff line change
Expand Up @@ -112,6 +112,10 @@ void FrequencySystem::init_data ()
// make sure we have frequencies to solve for
if (!_finished_set_frequencies)
{
// not supported for now
if (parameters.have_parameter<unsigned int> ("n_frequencies"))
libmesh_not_implemented_msg("ERROR: Setting the n_frequencies parameter on the system is not supported");

// when this system was read from file, check
// if this has a "n_frequencies" parameter,
// and initialize us with these.
Expand Down Expand Up @@ -339,10 +343,7 @@ void FrequencySystem::solve (const unsigned int n_start,
// Get the user-specified linear solver tolerance,
// the user-specified maximum # of linear solver iterations,
// the user-specified wave speed
const Real tol =
es.parameters.get<Real>("linear solver tolerance");
const unsigned int maxits =
es.parameters.get<unsigned int>("linear solver maximum iterations");
const auto [maxits, tol] = this->get_linear_solve_parameters();

// start solver loop
for (unsigned int n=n_start; n<= n_stop; n++)
Expand Down
Loading