36 #include "libmesh/composite_function.h" 
   37 #include "libmesh/getpot.h" 
   38 #include "libmesh/parameter_multiaccessor.h" 
   44                                           const std::string& name,
 
   45                                           const unsigned int number )
 
   46     : FEMSystem(es, name, number),
 
   47       _use_numerical_jacobians_only(false)
 
   61     this->verify_analytic_jacobians = input(
"linear-nonlinear-solver/verify_analytic_jacobians", 0.0 );
 
   62     this->print_solution_norms = input(
"screen-options/print_solution_norms", 
false );
 
   63     this->print_solutions = input(
"screen-options/print_solutions", 
false );
 
   64     this->print_residual_norms = input(
"screen-options/print_residual_norms", 
false );
 
   68     this->print_residuals = input(
"screen-options/print_residual", 
false );
 
   69     if (this->print_residuals)
 
   72     this->print_residuals = input(
"screen-options/print_residuals", this->print_residuals );
 
   73     this->print_jacobian_norms = input(
"screen-options/print_jacobian_norms", 
false );
 
   74     this->print_jacobians = input(
"screen-options/print_jacobians", 
false );
 
   75     this->print_element_solutions = input(
"screen-options/print_element_solutions", 
false );
 
   76     this->print_element_residuals = input(
"screen-options/print_element_residuals", 
false );
 
   77     this->print_element_jacobians = input(
"screen-options/print_element_jacobians", 
false );
 
   81     numerical_jacobian_h =
 
   82       input(
"linear-nonlinear-solver/numerical_jacobian_h",
 
   83             numerical_jacobian_h);
 
   85     const unsigned int n_numerical_jacobian_h_values =
 
   86       input.vector_variable_size
 
   87         (
"linear-nonlinear-solver/numerical_jacobian_h_values");
 
   89     if (n_numerical_jacobian_h_values !=
 
   90         input.vector_variable_size
 
   91           (
"linear-nonlinear-solver/numerical_jacobian_h_variables"))
 
   93         std::cerr << 
"Error: found " << n_numerical_jacobian_h_values
 
   94                   << 
" numerical_jacobian_h_values" << std::endl;
 
   96                   << input.vector_variable_size
 
   97                        (
"linear-nonlinear-solver/numerical_jacobian_h_variables")
 
   98                 << 
" numerical_jacobian_h_variables" << std::endl;
 
  104     for (
unsigned int i=0; i != n_numerical_jacobian_h_values; ++i)
 
  107           input(
"linear-nonlinear-solver/numerical_jacobian_h_variables",
 
  110           input(
"linear-nonlinear-solver/numerical_jacobian_h_values",
 
  111                 libMesh::Real(0), i);
 
  120     use_fixed_solution = 
true;
 
  131         (physics_iter->second)->init_variables( 
this );
 
  141         unsigned int var_num =
 
  143         this->set_numerical_jacobian_h_for_var
 
  152         (physics_iter->second)->set_time_evolving_vars( 
this );
 
  158       (
_physics_list.begin()->second)->set_is_steady((this->time_solver)->is_steady());
 
  162     libMesh::FEMSystem::init_data();
 
  172         (physics_iter->second)->init_ics( 
this, ic_function );
 
  175     if (ic_function.n_subfunctions())
 
  177         this->project_solution(&ic_function);
 
  185         (physics_iter->second)->auxiliary_init( *
this );
 
  195     libMesh::UniquePtr<libMesh::DiffContext> ap(context);
 
  197     libMesh::DifferentiablePhysics* phys = libMesh::FEMSystem::get_physics();
 
  199     libmesh_assert(phys);
 
  202     context->set_mesh_system(phys->get_mesh_system());
 
  203     context->set_mesh_x_var(phys->get_mesh_x_var());
 
  204     context->set_mesh_y_var(phys->get_mesh_y_var());
 
  205     context->set_mesh_z_var(phys->get_mesh_z_var());
 
  207     ap->set_deltat_pointer( &deltat );
 
  210     ap->is_adjoint() = this->get_time_solver().is_adjoint();
 
  229       ( 
const std::string & param_name,
 
  234          physics_iter != _physics_list.end();
 
  237         (physics_iter->second)->register_parameter( param_name,
 
  246     AssemblyContext& c = libMesh::libmesh_cast_ref<AssemblyContext&>(context);
 
  261                                               libMesh::DiffContext& context,
 
  263                                               CacheFuncType cachefunc)
 
  265     AssemblyContext& c = libMesh::libmesh_cast_ref<AssemblyContext&>(context);
 
  267     bool compute_jacobian = 
true;
 
  278         ((*(physics_iter->second)).*cachefunc)( c, cache );
 
  288             if( (physics_iter->second)->enabled_on_elem( &c.get_elem() ) )
 
  290                 ((*(physics_iter->second)).*resfunc)( compute_jacobian, c, cache );
 
  295             ((*(physics_iter->second)).*resfunc)( compute_jacobian, c, cache );
 
  301     return compute_jacobian;
 
  305                                                     libMesh::DiffContext& context )
 
  315                                                  libMesh::DiffContext& context )
 
  320     jacobian_computed = jacobian_computed &&
 
  327     return jacobian_computed;
 
  331                                                      libMesh::DiffContext& context )
 
  341                                                libMesh::DiffContext& context )
 
  351                                             libMesh::DiffContext& context )
 
  361                                                 libMesh::DiffContext& context )
 
  371                                              libMesh::DiffContext& context )
 
  381                                           libMesh::DiffContext& context )
 
  391                                                    libMesh::DiffContext& context )
 
  404         std::cerr << 
"Error: Could not find physics " << physics_name << std::endl;
 
  423                                                            const libMesh::Point& point,
 
  424                                                            libMesh::Real& value )
 
  431         if( (physics_iter->second)->enabled_on_elem( &context.get_elem() ) )
 
  440                                                    const std::vector<SharedPtr<NeumannBCContainer> >& neumann_bcs,
 
  441                                                    std::vector<SharedPtr<NeumannBCContainer> >& active_neumann_bcs )
 
  444     for( std::vector<SharedPtr<NeumannBCContainer> >::const_iterator it = neumann_bcs.begin();
 
  445          it < neumann_bcs.end(); ++it )
 
  446       if( (*it)->has_bc_id( bc_id ) )
 
  447           active_neumann_bcs.push_back( *it );
 
  451                                               libMesh::DiffContext& context )
 
  454       libMesh::libmesh_cast_ref<AssemblyContext&>( context );
 
  456     std::vector<BoundaryID> ids = assembly_context.side_boundary_ids();
 
  458     bool compute_jacobian = request_jacobian;
 
  461     for( std::vector<BoundaryID>::const_iterator it = ids.begin();
 
  462          it != ids.end(); it++ )
 
  466         libmesh_assert_not_equal_to(bc_id, libMesh::BoundaryInfo::invalid_id);
 
  469         std::vector<SharedPtr<NeumannBCContainer> > active_neumann_bcs;
 
  472         if( !active_neumann_bcs.empty() )
 
  474             typedef std::vector<SharedPtr<NeumannBCContainer> >::iterator BCIt;
 
  476             for( BCIt container = active_neumann_bcs.begin(); container < active_neumann_bcs.end(); ++container )
 
  478                 SharedPtr<NeumannBCAbstract>& func = (*container)->get_func();
 
  482                 func->eval_flux( compute_jacobian, assembly_context,
 
  488     return compute_jacobian;
 
  491 #ifdef GRINS_USE_GRVY_TIMERS 
  492   void MultiphysicsSystem::attach_grvy_timer( GRVY::GRVY_Timer_Class* grvy_timer )
 
  501         (physics_iter->second)->attach_grvy_timer( grvy_timer );
 
static bool is_axisymmetric()
 
bool _use_numerical_jacobians_only
 
bool has_physics(const std::string physics_name) const 
Query to check if a particular physics has been enabled. 
 
static void build_boundary_conditions(const GetPot &input, MultiphysicsSystem &system, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
 
void attach_physics_list(PhysicsList physics_list)
PhysicsList gets built by GRINS::PhysicsFactory and attached here. 
 
std::vector< std::string > _numerical_jacobian_h_variables
 
std::vector< SharedPtr< NeumannBCContainer > > _neumann_bcs
Neumann boundary conditions. 
 
virtual void init_data()
System initialization. Calls each physics implementation of init_variables() 
 
void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Each Physics will register their postprocessed quantities with this call. 
 
virtual bool side_constraint(bool request_jacobian, libMesh::DiffContext &context)
Boundary contributions to  which do not have time varying components. 
 
libMesh::boundary_id_type BoundaryID
More descriptive name of the type used for boundary ids. 
 
virtual void damping_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
 
virtual bool damping_residual(bool request_jacobian, libMesh::DiffContext &context)
Contributions to . 
 
virtual void nonlocal_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for scalar variables. 
 
SharedPtr< GRINS::Physics > get_physics(const std::string physics_name)
 
PhysicsList _physics_list
Container of pointers to GRINS::Physics classes requested at runtime. 
 
virtual void compute_mass_residual_cache(const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_element_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
 
virtual bool nonlocal_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Contributions to  on SCALAR variables which have time varying components. 
 
virtual void element_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for element interiors. 
 
virtual void compute_side_constraint_cache(const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
libMesh::Real neumann_bc_sign() const 
 
void register_parameter(const std::string ¶m_name, libMesh::ParameterMultiAccessor< libMesh::Number > ¶m_pointer)
Each Physics will register its copy(s) of an independent variable. 
 
virtual bool nonlocal_mass_residual(bool request_jacobian, libMesh::DiffContext &context)
Contributions to  on SCALAR variables. 
 
virtual bool element_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Element interior contributions to  which have time varying components. 
 
virtual void compute_element_constraint_cache(const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
 
virtual libMesh::UniquePtr< libMesh::DiffContext > build_context()
Override FEMSystem::build_context in order to use our own AssemblyContext. 
 
virtual bool nonlocal_constraint(bool request_jacobian, libMesh::DiffContext &context)
Contributions to  on SCALAR variables which do not have time varying components. 
 
virtual void side_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for boundaries of elements on the domain boundary. 
 
virtual void nonlocal_mass_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Mass matrix part(s) for scalar variables. 
 
void get_active_neumann_bcs(BoundaryID bc_id, const std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs, std::vector< SharedPtr< NeumannBCContainer > > &active_neumann_bcs)
Extract the bcs from neumann_bcs that are active on bc_id and return them in active_neumann_bcs. 
 
virtual void compute_damping_residual_cache(const AssemblyContext &context, CachedValues &cache)
 
std::map< std::string, SharedPtr< GRINS::Physics > > PhysicsList
Container for GRINS::Physics object pointers. 
 
std::vector< libMesh::Real > _numerical_jacobian_h_values
 
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors. 
 
bool _general_residual(bool request_jacobian, libMesh::DiffContext &context, ResFuncType resfunc, CacheFuncType cachefunc)
 
virtual void side_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for boundaries of elements on the domain boundary. 
 
virtual bool element_constraint(bool request_jacobian, libMesh::DiffContext &context)
 
virtual bool mass_residual(bool request_jacobian, libMesh::DiffContext &context)
Contributions to . 
 
virtual void compute_nonlocal_mass_residual_cache(const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_time_derivative_cache(const AssemblyContext &context, CachedValues &cache)
 
bool apply_neumann_bcs(bool request_jacobian, libMesh::DiffContext &context)
Applies the subset of _neumann_bcs that are active on the current element side. 
 
virtual void compute_nonlocal_constraint_cache(const AssemblyContext &context, CachedValues &cache)
 
MultiphysicsSystem(libMesh::EquationSystems &es, const std::string &name, const unsigned int number)
Constructor. Will be called by libMesh only. 
 
virtual bool side_time_derivative(bool request_jacobian, libMesh::DiffContext &context)
Boundary contributions to  which have time varying components. 
 
std::map< std::string, SharedPtr< GRINS::Physics > >::const_iterator PhysicsListIter
Iterator for PhysicsList. 
 
virtual void init_context(libMesh::DiffContext &context)
Context initialization. Calls each physics implementation of init_context() 
 
virtual void nonlocal_constraint(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Constraint part(s) of physics for scalar variables. 
 
virtual void read_input_options(const GetPot &input)
Reads input options for this class and all physics that are enabled. 
 
const GetPot * _input
Cached for helping build boundary conditions. 
 
virtual void mass_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...