LCOV - code coverage report
Current view: top level - include/systems - equation_systems.h (source / functions) Hit Total Coverage
Test: libMesh/libmesh: #4229 (6a9aeb) with base 727f46 Lines: 48 55 87.3 %
Date: 2025-08-19 19:27:09 Functions: 19 45 42.2 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : // The libMesh Finite Element Library.
       2             : // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
       3             : 
       4             : // This library is free software; you can redistribute it and/or
       5             : // modify it under the terms of the GNU Lesser General Public
       6             : // License as published by the Free Software Foundation; either
       7             : // version 2.1 of the License, or (at your option) any later version.
       8             : 
       9             : // This library is distributed in the hope that it will be useful,
      10             : // but WITHOUT ANY WARRANTY; without even the implied warranty of
      11             : // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
      12             : // Lesser General Public License for more details.
      13             : 
      14             : // You should have received a copy of the GNU Lesser General Public
      15             : // License along with this library; if not, write to the Free Software
      16             : // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
      17             : 
      18             : 
      19             : 
      20             : #ifndef LIBMESH_EQUATION_SYSTEMS_H
      21             : #define LIBMESH_EQUATION_SYSTEMS_H
      22             : 
      23             : // Local Includes
      24             : #include "libmesh/libmesh_common.h"
      25             : #include "libmesh/parameters.h"
      26             : #include "libmesh/system.h"
      27             : #include "libmesh/parallel_object.h"
      28             : 
      29             : // HP aCC needs these for some reason
      30             : #ifdef __HP_aCC
      31             : # include "libmesh/frequency_system.h"
      32             : # include "libmesh/transient_system.h"
      33             : # include "libmesh/newmark_system.h"
      34             : # include "libmesh/steady_system.h"
      35             : #endif
      36             : 
      37             : // C++ includes
      38             : #include <cstddef>
      39             : #include <map>
      40             : #include <set>
      41             : #include <string>
      42             : #include <string_view>
      43             : #include <vector>
      44             : #include <memory>
      45             : 
      46             : namespace libMesh
      47             : {
      48             : 
      49             : // Forward Declarations
      50             : class Elem;
      51             : class MeshBase;
      52             : enum XdrMODE : int;
      53             : 
      54             : /**
      55             :  * This is the \p EquationSystems class.  It is in charge
      56             :  * of handling all the various equation systems defined
      57             :  * for a \p MeshBase.  It may have multiple systems, which may
      58             :  * be active or inactive, so that at different solution
      59             :  * stages only a sub-set may be solved for.  Also, through
      60             :  * the templated access, different types of systems
      61             :  * may be handled.  Also other features, like flags,
      62             :  * parameters, I/O etc are provided.
      63             :  *
      64             :  * \author Benjamin S. Kirk
      65             :  * \date 2002-2007
      66             :  * \brief Manages multiples systems of equations.
      67             :  */
      68       24256 : class EquationSystems : public ReferenceCountedObject<EquationSystems>,
      69             :                         public ParallelObject
      70             : 
      71             : {
      72             : public:
      73             : 
      74             :   /**
      75             :    * Define enumeration to set properties in EquationSystems::read()
      76             :    */
      77             :   enum ReadFlags { READ_HEADER           = 1,
      78             :                    READ_DATA             = 2,
      79             :                    READ_ADDITIONAL_DATA  = 4,
      80             :                    READ_LEGACY_FORMAT    = 8,
      81             :                    TRY_READ_IFEMS        = 16,
      82             :                    READ_BASIC_ONLY       = 32 };
      83             : 
      84             :   /**
      85             :    * Define enumeration to set properties in EquationSystems::write()
      86             :    */
      87             :   enum WriteFlags { WRITE_DATA             = 1,
      88             :                     WRITE_ADDITIONAL_DATA  = 2,
      89             :                     WRITE_PARALLEL_FILES   = 4,
      90             :                     WRITE_SERIAL_FILES     = 8 };
      91             : 
      92             :   /**
      93             :    * Constructor.
      94             :    */
      95             :   EquationSystems (MeshBase & mesh);
      96             : 
      97             :   /**
      98             :    * Destructor.  Should be virtual, since the user may want to derive
      99             :    * subclasses of EquationSystems.
     100             :    */
     101             :   virtual ~EquationSystems ();
     102             : 
     103             :   /**
     104             :    * Restores the data structure to a pristine state.
     105             :    */
     106             :   virtual void clear ();
     107             : 
     108             :   /**
     109             :    * Initialize all the systems
     110             :    */
     111             :   virtual void init ();
     112             : 
     113             :   /**
     114             :    * Handle any mesh changes and reinitialize all the systems on the
     115             :    * updated mesh
     116             :    */
     117             :   virtual void reinit ();
     118             : 
     119             :   /**
     120             :    * Handle the association of a completely *new* mesh with the
     121             :    * EquationSystem and all the Systems assigned to it.
     122             :    */
     123             :   virtual void reinit_mesh ();
     124             : 
     125             :   /**
     126             :    * Enable or disable default ghosting functors on the Mesh and on
     127             :    * all Systems.  Standard ghosting is enabled by default.  If
     128             :    * disabled, default ghosting will also be disabled on any later
     129             :    * added systems.
     130             :    *
     131             :    * Unless other equivalent ghosting functors have been added,
     132             :    * removing the default coupling functor is only safe for explicit
     133             :    * solves, and removing the default algebraic ghosting functor is
     134             :    * only safe for codes where no evaluations on neighbor cells (e.g.
     135             :    * no jump error estimators) are done.
     136             :    */
     137             :   virtual void enable_default_ghosting (bool enable);
     138             : 
     139             :   /**
     140             :    * Updates local values for all the systems
     141             :    */
     142             :   void update ();
     143             : 
     144             :   /**
     145             :    * \returns The number of equation systems.
     146             :    */
     147             :   unsigned int n_systems() const;
     148             : 
     149             :   /**
     150             :    * \returns \p true if the system named \p name exists within
     151             :    * this EquationSystems object.
     152             :    */
     153             :   bool has_system (std::string_view name) const;
     154             : 
     155             :   /**
     156             :    * \returns A constant reference to the system named \p name.
     157             :    * The template argument defines the return type.  For example,
     158             :    * const SteadySystem & sys = eq.get_system<SteadySystem> ("sys");
     159             :    * is an example of how the method might be used
     160             :    */
     161             :   template <typename T_sys>
     162             :   const T_sys & get_system (std::string_view name) const;
     163             : 
     164             :   /**
     165             :    * \returns A writable reference to the system named \p name.
     166             :    * The template argument defines the return type.  For example,
     167             :    * const SteadySystem & sys = eq.get_system<SteadySystem> ("sys");
     168             :    * is an example of how the method might be used
     169             :    */
     170             :   template <typename T_sys>
     171             :   T_sys & get_system (std::string_view name);
     172             : 
     173             :   /**
     174             :    * \returns A constant reference to system number \p num.
     175             :    * The template argument defines the return type.  For example,
     176             :    * const SteadySystem & sys = eq.get_system<SteadySystem> (0);
     177             :    * is an example of how the method might be used
     178             :    */
     179             :   template <typename T_sys>
     180             :   const T_sys & get_system (const unsigned int num) const;
     181             : 
     182             :   /**
     183             :    * \returns A writable reference to the system number \p num.
     184             :    * The template argument defines the return type.  For example,
     185             :    * const SteadySystem & sys = eq.get_system<SteadySystem> (0);
     186             :    * is an example of how the method might be used
     187             :    */
     188             :   template <typename T_sys>
     189             :   T_sys & get_system (const unsigned int num);
     190             : 
     191             :   /**
     192             :    * \returns A constant reference to the system named \p name.
     193             :    */
     194             :   const System & get_system (std::string_view name) const;
     195             : 
     196             :   /**
     197             :    * \returns A writable reference to the system named \p name.
     198             :    */
     199             :   System & get_system (std::string_view name);
     200             : 
     201             :   /**
     202             :    * \returns A constant reference to system number \p num.
     203             :    */
     204             :   const System & get_system (const unsigned int num) const;
     205             : 
     206             :   /**
     207             :    * \returns A writable reference to the system number \p num.
     208             :    */
     209             :   System & get_system (const unsigned int num);
     210             : 
     211             :   /**
     212             :    * Add the system of type \p system_type named \p name to the
     213             :    * systems array.
     214             :    */
     215             :   virtual System & add_system (std::string_view system_type,
     216             :                                std::string_view name);
     217             : 
     218             :   /**
     219             :    * Add the system named \p name to the systems array.
     220             :    */
     221             :   template <typename T_sys>
     222             :   T_sys & add_system (std::string_view name);
     223             : 
     224             :   /**
     225             :    * \returns The total number of variables in all
     226             :    * systems.
     227             :    */
     228             :   unsigned int n_vars () const;
     229             : 
     230             :   /**
     231             :    * \returns The total number of degrees of freedom
     232             :    * in all systems.
     233             :    */
     234             :   std::size_t n_dofs () const;
     235             : 
     236             :   /**
     237             :    * \returns The number of active degrees of freedom
     238             :    * for the EquationSystems object.
     239             :    */
     240             :   std::size_t n_active_dofs() const;
     241             : 
     242             :   /**
     243             :    * Call \p solve on all the individual equation systems.
     244             :    *
     245             :    * By default this function solves each equation system once,
     246             :    * in the order they were added.  For more sophisticated decoupled
     247             :    * problems the user may with to override this behavior in a derived
     248             :    * class.
     249             :    */
     250             :   virtual void solve ();
     251             : 
     252             :   /**
     253             :    * Call \p adjoint_solve on all the individual equation systems.
     254             :    *
     255             :    * By default this function solves each system's adjoint once,
     256             :    * in the reverse order from that in which they were added.  For
     257             :    * more sophisticated decoupled problems the user may with to
     258             :    * override this behavior in a derived class.
     259             :    */
     260             :   virtual void adjoint_solve (const QoISet & qoi_indices = QoISet());
     261             : 
     262             :   /**
     263             :    * Call \p sensitivity_solve on all the individual equation systems.
     264             :    *
     265             :    * By default this function solves each sensitivity system once,
     266             :    * in the order in which in which they were added.  For
     267             :    * more sophisticated decoupled problems the user may with to
     268             :    * override this behavior in a derived class.
     269             :    */
     270             :   virtual void sensitivity_solve (const ParameterVector & parameters);
     271             : 
     272             :   /**
     273             :    * Fill the input vector \p var_names with the names
     274             :    * of the variables for each system. If \p type is passed,
     275             :    * only variables of the specified type will be populated.
     276             :    * If systems_names!=nullptr, only include names from the
     277             :    * specified systems.
     278             :    */
     279             :   void build_variable_names (std::vector<std::string> & var_names,
     280             :                              const FEType * type=nullptr,
     281             :                              const std::set<std::string> * system_names=nullptr) const;
     282             : 
     283             :   /**
     284             :    * Fill the input vector \p soln with the solution values for the
     285             :    * system named \p name.
     286             :    *
     287             :    * \note The input vector \p soln will only be assembled on
     288             :    * processor 0, so this method is only applicable to outputting plot
     289             :    * files from processor 0.
     290             :    */
     291             :   void build_solution_vector (std::vector<Number> & soln,
     292             :                               std::string_view system_name,
     293             :                               std::string_view variable_name = "all_vars") const;
     294             : 
     295             :   /**
     296             :    * Fill the input vector \p soln with solution values.  The
     297             :    * entries will be in variable-major format (corresponding to
     298             :    * the names from \p build_variable_names()).
     299             :    *
     300             :    * If systems_names!=nullptr, only include data from the
     301             :    * specified systems.
     302             :    *
     303             :    * If \p add_sides is true, append data for plotting on "side
     304             :    * elements" too.
     305             :    */
     306             :   void build_solution_vector (std::vector<Number> & soln,
     307             :                               const std::set<std::string> * system_names=nullptr,
     308             :                               bool add_sides=false) const;
     309             : 
     310             :   /**
     311             :    * A version of build_solution_vector which is appropriate for
     312             :    * "parallel" output formats like Nemesis.
     313             :    *
     314             :    * \returns A std::unique_ptr to a node-major NumericVector of total
     315             :    * length n_nodes*n_vars that various I/O classes can then use to
     316             :    * get the local values they need to write on each processor.
     317             :    */
     318             :   std::unique_ptr<NumericVector<Number>>
     319             :   build_parallel_solution_vector(const std::set<std::string> * system_names=nullptr,
     320             :                                  bool add_sides=false) const;
     321             : 
     322             :   /**
     323             :    * Retrieve \p vars_active_subdomains, which indicates the active
     324             :    * subdomains for each variable in \p names.
     325             :    */
     326             :   void get_vars_active_subdomains(const std::vector<std::string> & names,
     327             :                                   std::vector<std::set<subdomain_id_type>> & vars_active_subdomains) const;
     328             : 
     329             : #ifdef LIBMESH_ENABLE_DEPRECATED
     330             :   /**
     331             :    * Retrieve the solution data for CONSTANT MONOMIALs and/or components of
     332             :    * CONSTANT MONOMIAL_VECs. If 'names' is populated, only the variables
     333             :    * corresponding to those names will be retrieved. This can be used to
     334             :    * filter which variables are retrieved.
     335             :    *
     336             :    * \deprecated Call the more appropriately-named build_elemental_solution_vector()
     337             :    * instead.
     338             :    */
     339             :   void get_solution (std::vector<Number> & soln,
     340             :                      std::vector<std::string> & names) const;
     341             : #endif // LIBMESH_ENABLE_DEPRECATED
     342             : 
     343             :   /**
     344             :    * Retrieve the solution data for CONSTANT MONOMIALs and/or components of
     345             :    * CONSTANT MONOMIAL_VECs. If 'names' is populated, only the variables
     346             :    * corresponding to those names will be retrieved. This can be used to
     347             :    * filter which variables are retrieved.
     348             :    *
     349             :    * This is the more appropriately-named replacement for the get_solution()
     350             :    * function defined above.
     351             :    */
     352             :   void build_elemental_solution_vector (std::vector<Number> & soln,
     353             :                                         std::vector<std::string> & names) const;
     354             : 
     355             :   /**
     356             :    * Finds system and variable numbers for any variables of 'type' or of
     357             :    * 'types' corresponding to the entries in the input 'names' vector.
     358             :    * If 'names' is empty, this returns all variables of the type.
     359             :    * The names of vector variables are decomposed into individual ones
     360             :    * suffixed with their cartesian component, but there will still be
     361             :    * a single pair of system numbers for such vector variables. Thus,
     362             :    * the size of the 'names' vector modified by this function may not
     363             :    * be equal to that of the returned vector of pairs. Nevertheless, both
     364             :    * should be sorted in accordance with ExodusII format, and so the
     365             :    * developer just needs to know to separate dof_indices when
     366             :    * accessing the system solution for vector variables.
     367             :    *
     368             :    * This function is designed to work for either a single type or a
     369             :    * vector of types, but not both. This is because it can't simply be
     370             :    * called a second time with another type as it filters (deletes) the
     371             :    * names of those on the first call that used a different type. Thus,
     372             :    * the 'types' argument is for the case where variables of multiple
     373             :    * types are allowed to pass through.
     374             :    *
     375             :    * TODO: find a more generic way to handle this whole procedure.
     376             :    */
     377             :   std::vector<std::pair<unsigned int, unsigned int>>
     378             :   find_variable_numbers (std::vector<std::string> & names,
     379             :                          const FEType * type=nullptr,
     380             :                          const std::vector<FEType> * types=nullptr) const;
     381             : 
     382             :   /**
     383             :    * Builds a parallel vector of CONSTANT MONOMIAL and/or components of
     384             :    * CONSTANT MONOMIAL_VEC solution values corresponding to the entries
     385             :    * in the input 'names' vector. This vector is approximately uniformly
     386             :    * distributed across all of the available processors.
     387             :    *
     388             :    * The related function build_elemental_solution_vector() is
     389             :    * implemented by calling this function and then calling
     390             :    * localize_to_one() on the resulting vector.
     391             :    *
     392             :    * Returns a nullptr if no CONSTANT, MONOMIAL/MONOMIAL_VEC variables
     393             :    * exist in the 'names' vector (if it is empty, then it will return all
     394             :    * variables in the system of this type if any) or a std::unique_ptr to
     395             :    * a var-major numeric vector of total length n_elem * n_vars, where
     396             :    * n_vars includes all components of vectors, ordered according to:
     397             :    * [u0, u1, ... uN, v0, v1, ... vN, w0, w1, ... wN] for constant monomial
     398             :    * variables (u, v, w) on a mesh with N elements.
     399             :    */
     400             :   std::unique_ptr<NumericVector<Number>>
     401             :   build_parallel_elemental_solution_vector (std::vector<std::string> & names) const;
     402             : 
     403             :   /**
     404             :    * Fill the input vector \p soln with solution values.  The
     405             :    * entries will be in variable-major format (corresponding to
     406             :    * the names from \p build_variable_names()).
     407             : 
     408             :    * If systems_names!=nullptr, only include data from the
     409             :    * specified systems.
     410             : 
     411             :    * If vertices_only == true, then for higher-order elements only
     412             :    * the solution at the vertices is computed. This can be useful
     413             :    * when plotting discontinuous solutions on higher-order elements
     414             :    * and only a lower-order representation is required.
     415             :    *
     416             :    * If add_sides == true, append data for plotting on "side
     417             :    * elements" too.
     418             :    */
     419             :   void build_discontinuous_solution_vector
     420             :   (std::vector<Number> & soln,
     421             :    const std::set<std::string> * system_names = nullptr,
     422             :    const std::vector<std::string> * var_names = nullptr,
     423             :    bool vertices_only = false,
     424             :    bool add_sides = false) const;
     425             : 
     426             :   /*
     427             :    * Returns true iff the given side of the given element is *never*
     428             :    * added to output from that element, because it is considered to be
     429             :    * redundant with respect to the same data added from the
     430             :    * neighboring element sharing that side.  This helper function is
     431             :    * used with the add_sides option when building solution vectors
     432             :    * here and when outputting solution vectors in MeshOutput
     433             :    * (currently just Exodus) I/O.
     434             :    */
     435             :   static bool redundant_added_side(const Elem & elem, unsigned int side);
     436             : 
     437             : 
     438             :   /**
     439             :    * Read & initialize the systems from disk using the XDR data format.
     440             :    * This format allows for machine-independent binary output.
     441             :    *
     442             :    * Set which sections of the file to read by bitwise OR'ing the
     443             :    * EquationSystems::ReadFlags enumeration together. For example, to
     444             :    * read all sections of the file, set read_flags to:
     445             :    * (READ_HEADER | READ_DATA | READ_ADDITIONAL_DATA)
     446             :    *
     447             :    * \note The equation system can be defined without initializing
     448             :    * the data vectors to any solution values.  This can be done
     449             :    * by omitting READ_DATA in the read_flags parameter.
     450             :    *
     451             :    * If XdrMODE is omitted, it will be inferred as READ for filenames
     452             :    * containing .xda or as DECODE for filenames containing .xdr
     453             :    *
     454             :    * \param name Name of the file to be read.
     455             :    * \param read_flags Single flag created by bitwise-OR'ing several flags together.
     456             :    * \param mode Controls whether reading is done in binary or ascii mode.
     457             :    * \param partition_agnostic If true then the mesh and degrees of freedom
     458             :    * will be temporarily renumbered in a partition agnostic way so that
     459             :    * files written using "n" mpi processes can be re-read on "m" mpi
     460             :    * processes.  This renumbering is not compatible with meshes
     461             :    * that have two nodes in exactly the same position!
     462             :    */
     463             :   template <typename InValType = Number>
     464             :   void read (std::string_view name,
     465             :              const XdrMODE,
     466             :              const unsigned int read_flags=(READ_HEADER | READ_DATA),
     467             :              bool partition_agnostic = true);
     468             : 
     469             :   template <typename InValType = Number>
     470             :   void read (std::string_view name,
     471             :              const unsigned int read_flags=(READ_HEADER | READ_DATA),
     472             :              bool partition_agnostic = true);
     473             : 
     474             :   template <typename InValType = Number>
     475             :   void read (Xdr & io,
     476             :              std::function<std::unique_ptr<Xdr>()> & local_io_functor,
     477             :              const unsigned int read_flags=(READ_HEADER | READ_DATA),
     478             :              bool partition_agnostic = true);
     479             : 
     480             :   /**
     481             :    * Write the systems to disk using the XDR data format.
     482             :    * This format allows for machine-independent binary output.
     483             :    *
     484             :    * Set the writing properties using the EquationSystems::WriteFlags
     485             :    * enumeration. Set which sections to write out by bitwise OR'ing
     486             :    * the enumeration values. Write everything by setting write_flags to:
     487             :    * (WRITE_DATA | WRITE_ADDITIONAL_DATA)
     488             :    *
     489             :    * \note The solution data can be omitted by calling
     490             :    * this routine with WRITE_DATA omitted in the write_flags argument.
     491             :    *
     492             :    * If XdrMODE is omitted, it will be inferred as WRITE for filenames
     493             :    * containing .xda or as ENCODE for filenames containing .xdr
     494             :    *
     495             :    * \param name Name of the file to be read.
     496             :    * \param write_flags Single flag created by bitwise-OR'ing several flags together.
     497             :    * \param mode Controls whether reading is done in binary or ascii mode.
     498             :    * \param partition_agnostic If true then the mesh and degrees of freedom
     499             :    * will be temporarily renumbered in a partition agnostic way so that
     500             :    * files written using "n" mpi processes can be re-read on "m" mpi
     501             :    * processes.  This renumbering is not compatible with meshes
     502             :    * that have two nodes in exactly the same position!
     503             :    */
     504             :   void write (std::string_view name,
     505             :               const XdrMODE,
     506             :               const unsigned int write_flags=(WRITE_DATA),
     507             :               bool partition_agnostic = true) const;
     508             : 
     509             :   void write (std::string_view name,
     510             :               const unsigned int write_flags=(WRITE_DATA),
     511             :               bool partition_agnostic = true) const;
     512             : 
     513             :   void write (std::ostream name,
     514             :               const unsigned int write_flags=(WRITE_DATA),
     515             :               bool partition_agnostic = true) const;
     516             : 
     517             :   void write (Xdr & io,
     518             :               const unsigned int write_flags=(WRITE_DATA),
     519             :               bool partition_agnostic = true,
     520             :               Xdr * const local_io = nullptr) const;
     521             : 
     522             :   /**
     523             :    * \returns \p true when this equation system contains
     524             :    * identical data, up to the given threshold.  Delegates
     525             :    * most of the comparisons to perform to the responsible
     526             :    * systems
     527             :    */
     528             :   virtual bool compare (const EquationSystems & other_es,
     529             :                         const Real threshold,
     530             :                         const bool verbose) const;
     531             : 
     532             :   /**
     533             :    * \returns A string containing information about the
     534             :    * systems, flags, and parameters.
     535             :    */
     536             :   virtual std::string get_info() const;
     537             : 
     538             :   /**
     539             :    * Prints information about the equation systems, by default to
     540             :    * libMesh::out.
     541             :    */
     542             :   void print_info (std::ostream & os=libMesh::out) const;
     543             : 
     544             :   /**
     545             :    * Same as above, but allows you to also use stream syntax.
     546             :    */
     547             :   friend std::ostream & operator << (std::ostream & os,
     548             :                                      const EquationSystems & es);
     549             : 
     550             :   /**
     551             :    * \returns A constant reference to the mesh
     552             :    */
     553             :   const MeshBase & get_mesh() const;
     554             : 
     555             :   /**
     556             :    * \returns A reference to the mesh
     557             :    */
     558             :   MeshBase & get_mesh();
     559             : 
     560             :   /**
     561             :    * Serializes a distributed mesh and its associated
     562             :    * degree of freedom numbering for all systems
     563             :    **/
     564             :   void allgather ();
     565             : 
     566             :   /**
     567             :    * Calls to reinit() will also do two-step coarsen-then-refine
     568             :    **/
     569             :   void enable_refine_in_reinit() { this->_refine_in_reinit = true; }
     570             : 
     571             :   /**
     572             :    * Calls to reinit() will not try to coarsen or refine the mesh
     573             :    **/
     574             :   void disable_refine_in_reinit() { this->_refine_in_reinit = false; }
     575             : 
     576             :   /**
     577             :    * \returns Whether or not calls to reinit() will try to coarsen/refine the mesh
     578             :    **/
     579             :   bool refine_in_reinit_flag() { return this->_refine_in_reinit; }
     580             : 
     581             :   /**
     582             :    * Handle any mesh changes and project any solutions onto the
     583             :    * updated mesh.
     584             :    *
     585             :    * \returns Whether or not the mesh may have changed.
     586             :    */
     587             :   bool reinit_solutions ();
     588             : 
     589             :   /**
     590             :    * Reinitialize all systems on the current mesh.
     591             :    */
     592             :   virtual void reinit_systems ();
     593             : 
     594             :   /**
     595             :    * Data structure holding arbitrary parameters.
     596             :    */
     597             :   Parameters parameters;
     598             : 
     599             : 
     600             : protected:
     601             : 
     602             : 
     603             :   /**
     604             :    * The mesh data structure
     605             :    */
     606             :   MeshBase & _mesh;
     607             : 
     608             :   /**
     609             :    * Data structure holding the systems.
     610             :    */
     611             :   std::map<std::string, std::unique_ptr<System>, std::less<>> _systems;
     612             : 
     613             :   /**
     614             :    * Flag for whether to call coarsen/refine in reinit().
     615             :    * Default value: true
     616             :    */
     617             :   bool _refine_in_reinit;
     618             : 
     619             :   /**
     620             :    * Flag for whether to enable default ghosting on newly added Systems.
     621             :    * Default value: true
     622             :    */
     623             :   bool _enable_default_ghosting;
     624             : 
     625             : private:
     626             :   /**
     627             :    * This function is used in the implementation of add_system,
     628             :    * it loops over the nodes and elements of the Mesh, adding the
     629             :    * system to each one.  The main reason to separate this part
     630             :    * is to avoid coupling this header file to mesh.h, and elem.h.
     631             :    */
     632             :   void _add_system_to_nodes_and_elems();
     633             : 
     634             :   /**
     635             :    * This just calls DofMap::remove_default_ghosting() but using a
     636             :    * shim lets us forward-declare DofMap.
     637             :    */
     638             :   void _remove_default_ghosting(unsigned int sys_num);
     639             : };
     640             : 
     641             : 
     642             : 
     643             : // ------------------------------------------------------------
     644             : // EquationSystems inline methods
     645             : inline
     646        6292 : const MeshBase & EquationSystems::get_mesh () const
     647             : {
     648      122072 :   return _mesh;
     649             : }
     650             : 
     651             : 
     652             : 
     653             : inline
     654         266 : MeshBase & EquationSystems::get_mesh ()
     655             : {
     656      270504 :   return _mesh;
     657             : }
     658             : 
     659             : 
     660             : inline
     661           0 : unsigned int EquationSystems::n_systems () const
     662             : {
     663           0 :   return cast_int<unsigned int>(_systems.size());
     664             : }
     665             : 
     666             : 
     667             : 
     668             : 
     669             : template <typename T_sys>
     670             : inline
     671       13131 : T_sys & EquationSystems::add_system (std::string_view name)
     672             : {
     673       13131 :   if (!_systems.count(name))
     674             :     {
     675       12701 :       const unsigned int sys_num = this->n_systems();
     676             : 
     677        1290 :       auto result = _systems.emplace
     678       24542 :         (name, std::make_unique<T_sys>(*this, std::string(name),
     679             :                                        sys_num));
     680             : 
     681       13131 :       if (!_enable_default_ghosting)
     682           0 :         this->_remove_default_ghosting(sys_num);
     683             : 
     684             :       // Tell all the \p DofObject entities to add a system.
     685       13131 :       this->_add_system_to_nodes_and_elems();
     686             : 
     687             :       // Return reference to newly added item
     688         430 :       auto it = result.first;
     689         430 :       auto & sys_ptr = it->second;
     690         430 :       return cast_ref<T_sys &>(*sys_ptr);
     691             :     }
     692             :   else
     693             :     {
     694             :       // We now allow redundant add_system calls, to make it
     695             :       // easier to load data from files for user-derived system
     696             :       // subclasses
     697           0 :       return this->get_system<T_sys>(name);
     698             :     }
     699             : }
     700             : 
     701             : 
     702             : 
     703             : inline
     704         960 : bool EquationSystems::has_system (std::string_view name) const
     705             : {
     706         960 :   if (_systems.find(name) == _systems.end())
     707           0 :     return false;
     708         960 :   return true;
     709             : }
     710             : 
     711             : 
     712             : 
     713             : 
     714             : template <typename T_sys>
     715             : inline
     716       25148 : const T_sys & EquationSystems::get_system (const unsigned int num) const
     717             : {
     718         710 :   libmesh_assert_less (num, this->n_systems());
     719             : 
     720       26564 :   for (auto & pr : _systems)
     721             :     {
     722         750 :       const auto & sys_ptr = pr.second;
     723       26564 :       if (sys_ptr->number() == num)
     724       25148 :         return cast_ref<const T_sys &>(*sys_ptr);
     725             :     }
     726             :   // Error if we made it here
     727           0 :   libmesh_error_msg("ERROR: no system number " << num << " found!");
     728             : }
     729             : 
     730             : 
     731             : 
     732             : 
     733             : template <typename T_sys>
     734             : inline
     735      348500 : T_sys & EquationSystems::get_system (const unsigned int num)
     736             : {
     737       16952 :   libmesh_assert_less (num, this->n_systems());
     738             : 
     739      355584 :   for (auto & pr : _systems)
     740             :     {
     741       17188 :       auto & sys_ptr = pr.second;
     742      355584 :       if (sys_ptr->number() == num)
     743      348500 :         return cast_ref<T_sys &>(*sys_ptr);
     744             :     }
     745             : 
     746             :   // Error if we made it here
     747           0 :   libmesh_error_msg("ERROR: no system number " << num << " found!");
     748             : }
     749             : 
     750             : 
     751             : 
     752             : 
     753             : 
     754             : 
     755             : template <typename T_sys>
     756             : inline
     757       88826 : const T_sys & EquationSystems::get_system (std::string_view name) const
     758             : {
     759        2690 :   auto pos = _systems.find(name);
     760             : 
     761             :   // Check for errors
     762       88826 :   libmesh_error_msg_if(pos == _systems.end(), "ERROR: no system named \"" << name << "\" found!");
     763             : 
     764             :   // Attempt dynamic cast
     765        2690 :   const auto & sys_ptr = pos->second;
     766       91516 :   return cast_ref<const T_sys &>(*sys_ptr);
     767             : }
     768             : 
     769             : 
     770             : 
     771             : 
     772             : 
     773             : 
     774             : template <typename T_sys>
     775             : inline
     776       10791 : T_sys & EquationSystems::get_system (std::string_view name)
     777             : {
     778         310 :   auto pos = _systems.find(name);
     779             : 
     780             :   // Check for errors
     781       10791 :   libmesh_error_msg_if(pos == _systems.end(), "ERROR: no system named " << name << " found!");
     782             : 
     783             :   // Attempt dynamic cast
     784         310 :   auto & sys_ptr = pos->second;
     785       11101 :   return cast_ref<T_sys &>(*sys_ptr);
     786             : }
     787             : 
     788             : 
     789             : 
     790             : 
     791             : 
     792             : 
     793             : 
     794             : inline
     795        1922 : const System & EquationSystems::get_system (std::string_view name) const
     796             : {
     797       61577 :   return this->get_system<System>(name);
     798             : }
     799             : 
     800             : 
     801             : 
     802             : inline
     803         310 : System & EquationSystems::get_system (std::string_view name)
     804             : {
     805       12224 :   return this->get_system<System>(name);
     806             : }
     807             : 
     808             : 
     809             : 
     810             : inline
     811         710 : const System & EquationSystems::get_system (const unsigned int num) const
     812             : {
     813       24872 :   return this->get_system<System>(num);
     814             : }
     815             : 
     816             : 
     817             : 
     818             : inline
     819       16952 : System & EquationSystems::get_system (const unsigned int num)
     820             : {
     821      348500 :   return this->get_system<System>(num);
     822             : }
     823             : 
     824             : 
     825             : } // namespace libMesh
     826             : 
     827             : 
     828             : #endif // LIBMESH_EQUATION_SYSTEMS_H

Generated by: LCOV version 1.14