The source file assembly.h with comments:

        #ifndef __assembly_h__
        #define __assembly_h__
        
        #if defined(LIBMESH_HAVE_SLEPC) && defined(LIBMESH_HAVE_GLPK)
        
        #include "libmesh/sparse_matrix.h"
        #include "libmesh/numeric_vector.h"
        #include "libmesh/dense_matrix.h"
        #include "libmesh/dense_vector.h"
        #include "libmesh/fe.h"
        #include "libmesh/fe_interface.h"
        #include "libmesh/fe_base.h"
        #include "libmesh/elem_assembly.h"
        #include "libmesh/quadrature_gauss.h"
        
rbOOmit includes
        #include "libmesh/rb_theta.h"
        #include "libmesh/rb_assembly_expansion.h"
        
Bring in bits from the libMesh namespace. Just the bits we're using, since this is a header.
        using libMesh::ElemAssembly;
        using libMesh::FEInterface;
        using libMesh::FEMContext;
        using libMesh::Number;
        using libMesh::Point;
        using libMesh::RBTheta;
        using libMesh::Real;
        using libMesh::RealGradient;
        
Functors for the parameter-dependent part of the affine decomposition of the PDE The RHS and outputs just require a constant value of 1, so use a default RBTheta object there
        struct ThetaA0 : RBTheta { virtual Number evaluate(const RBParameters& mu) { return mu.get_value("mu_0"); } };
        struct ThetaA1 : RBTheta { virtual Number evaluate(const RBParameters& mu) { return mu.get_value("mu_1"); } };
        struct ThetaA2 : RBTheta { virtual Number evaluate(const RBParameters& mu) { return mu.get_value("mu_2"); } };
        
        struct B : ElemAssembly
        {
Assemble the H1 inner product. This will be used as the inner product for this problem.
          virtual void interior_assembly(FEMContext &c)
          {
            const unsigned int u_var = 0;
        
            const std::vector<Real> &JxW =
              c.element_fe_var[u_var]->get_JxW();
        
            const std::vector<std::vector<Real> >& phi =
              c.element_fe_var[u_var]->get_phi();
        
The velocity shape function gradients at interior quadrature points.
            const std::vector<std::vector<RealGradient> >& dphi =
              c.element_fe_var[u_var]->get_dphi();
        
The number of local degrees of freedom in each variable
            const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
        
Now we will build the affine operator
            unsigned int n_qpoints = (c.get_element_qrule())->n_points();
        
            for (unsigned int qp=0; qp != n_qpoints; qp++)
              for (unsigned int i=0; i != n_u_dofs; i++)
                for (unsigned int j=0; j != n_u_dofs; j++)
                  c.get_elem_jacobian()(i,j) += JxW[qp] * (phi[j][qp]*phi[i][qp] + dphi[j][qp]*dphi[i][qp]);
          }
        };
        
        
        struct A0 : ElemAssembly
        {
Assemble the Laplacian operator
          virtual void interior_assembly(FEMContext &c)
          {
            const unsigned int u_var = 0;
        
            const std::vector<Real> &JxW =
              c.element_fe_var[u_var]->get_JxW();
        
The velocity shape function gradients at interior quadrature points.
            const std::vector<std::vector<RealGradient> >& dphi =
              c.element_fe_var[u_var]->get_dphi();
        
The number of local degrees of freedom in each variable
            const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
        
Now we will build the affine operator
            unsigned int n_qpoints = (c.get_element_qrule())->n_points();
        
            Real min_x=0.,max_x=0.5;
        
            Point centroid = c.elem->centroid();
            if( (min_x <= centroid(0)) && (centroid(0) < max_x) )
            for (unsigned int qp=0; qp != n_qpoints; qp++)
              for (unsigned int i=0; i != n_u_dofs; i++)
                for (unsigned int j=0; j != n_u_dofs; j++)
                  c.get_elem_jacobian()(i,j) += JxW[qp] * dphi[j][qp]*dphi[i][qp];
          }
        };
        
        struct A1 : ElemAssembly
        {
Assemble the Laplacian operator
          virtual void interior_assembly(FEMContext &c)
          {
            const unsigned int u_var = 0;
        
            const std::vector<Real> &JxW =
              c.element_fe_var[u_var]->get_JxW();
        
The velocity shape function gradients at interior quadrature points.
            const std::vector<std::vector<RealGradient> >& dphi =
              c.element_fe_var[u_var]->get_dphi();
        
The number of local degrees of freedom in each variable
            const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
        
Now we will build the affine operator
            unsigned int n_qpoints = (c.get_element_qrule())->n_points();
        
            Real min_x=0.5, max_x=1.;
        
            Point centroid = c.elem->centroid();
            if( (min_x <= centroid(0)) && (centroid(0) <= max_x) )
            for (unsigned int qp=0; qp != n_qpoints; qp++)
              for (unsigned int i=0; i != n_u_dofs; i++)
                for (unsigned int j=0; j != n_u_dofs; j++)
                  c.get_elem_jacobian()(i,j) += JxW[qp] * dphi[j][qp]*dphi[i][qp];
          }
        };
        
        struct A2 : ElemAssembly
        {
Convection in the x-direction
          virtual void interior_assembly(FEMContext &c)
          {
            const unsigned int u_var = 0;
        
            const std::vector<Real> &JxW =
              c.element_fe_var[u_var]->get_JxW();
        
            const std::vector<std::vector<Real> >& phi =
              c.element_fe_var[u_var]->get_phi();
        
            const std::vector<std::vector<RealGradient> >& dphi =
              c.element_fe_var[u_var]->get_dphi();
        
The number of local degrees of freedom in each variable
            const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
        
Now we will build the affine operator
            unsigned int n_qpoints = (c.get_element_qrule())->n_points();
        
            for (unsigned int qp=0; qp != n_qpoints; qp++)
              for (unsigned int i=0; i != n_u_dofs; i++)
                for (unsigned int j=0; j != n_u_dofs; j++)
                  c.get_elem_jacobian()(i,j) -= JxW[qp] *dphi[i][qp](0)*phi[j][qp];
          }
        };
        
        
        struct F0 : ElemAssembly
        {
Source term, 1 throughout the domain
          virtual void interior_assembly(FEMContext &c)
          {
            const unsigned int u_var = 0;
        
            const std::vector<Real> &JxW =
              c.element_fe_var[u_var]->get_JxW();
        
            const std::vector<std::vector<Real> >& phi =
              c.element_fe_var[u_var]->get_phi();
        
The number of local degrees of freedom in each variable
            const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
        
Now we will build the affine operator
            unsigned int n_qpoints = (c.get_element_qrule())->n_points();
        
            for (unsigned int qp=0; qp != n_qpoints; qp++)
              for (unsigned int i=0; i != n_u_dofs; i++)
                c.get_elem_residual()(i) += JxW[qp] * ( 1.*phi[i][qp] );
          }
        };
        
        struct OutputAssembly : ElemAssembly
        {
          OutputAssembly(Real min_x_in, Real max_x_in,
                         Real min_y_in, Real max_y_in)
                        :
                        min_x(min_x_in),
                        max_x(max_x_in),
                        min_y(min_y_in),
                        max_y(max_y_in)
          {}
        
Output: Average value over the region [min_x,max_x]x[min_y,max_y]
          virtual void interior_assembly(FEMContext &c)
          {
            const unsigned int u_var = 0;
        
            const std::vector<Real> &JxW =
              c.element_fe_var[u_var]->get_JxW();
        
            const std::vector<std::vector<Real> >& phi =
              c.element_fe_var[u_var]->get_phi();
        
The number of local degrees of freedom in each variable
            const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
        
Now we will build the affine operator
            unsigned int n_qpoints = (c.get_element_qrule())->n_points();
            
            Real output_area = (max_x-min_x) * (max_y-min_y);
        
            Point centroid = c.elem->centroid();
            if( (min_x <= centroid(0)) && (centroid(0) <= max_x) &&
                (min_y <= centroid(1)) && (centroid(1) <= max_y) )
              for (unsigned int qp=0; qp != n_qpoints; qp++)
                for (unsigned int i=0; i != n_u_dofs; i++)
                  c.get_elem_residual()(i) += JxW[qp] * ( 1.*phi[i][qp] ) / output_area;
          }
          
Member variables that define the output region in 2D
          Real min_x, max_x, min_y, max_y;
        };
        
Define an RBThetaExpansion class for this PDE
        struct Ex02RBThetaExpansion : RBThetaExpansion
        {
        
          /**
           * Constructor.
           */
          Ex02RBThetaExpansion()
          {
set up the RBThetaExpansion object
            attach_A_theta(&theta_a_0);   // Attach the lhs theta
            attach_A_theta(&theta_a_1);
            attach_A_theta(&theta_a_2);
        
            attach_F_theta(&rb_theta);    // Attach the rhs theta
        
            attach_output_theta(&rb_theta); // Attach output 0 theta
            attach_output_theta(&rb_theta); // Attach output 1 theta
            attach_output_theta(&rb_theta); // Attach output 2 theta
            attach_output_theta(&rb_theta); // Attach output 3 theta
          }
        
The RBTheta member variables
          ThetaA0 theta_a_0;
          ThetaA1 theta_a_1;
          ThetaA2 theta_a_2;
          RBTheta rb_theta; // Default RBTheta object, just returns 1.
        };
        
Define an RBAssemblyExpansion class for this PDE
        struct Ex02RBAssemblyExpansion : RBAssemblyExpansion
        {
        
          /**
           * Constructor.
           */
          Ex02RBAssemblyExpansion()
            :
            L0(0.72,0.88,0.72,0.88), // We make sure these output regions conform to the mesh
            L1(0.12,0.28,0.72,0.88),
            L2(0.12,0.28,0.12,0.28),
            L3(0.72,0.88,0.12,0.28)
          {
And set up the RBAssemblyExpansion object
            attach_A_assembly(&A0_assembly); // Attach the lhs assembly
            attach_A_assembly(&A1_assembly);
            attach_A_assembly(&A2_assembly);
            
            attach_F_assembly(&F0_assembly); // Attach the rhs assembly
            
            attach_output_assembly(&L0);       // Attach output 0 assembly
            attach_output_assembly(&L1);       // Attach output 1 assembly
            attach_output_assembly(&L2);       // Attach output 2 assembly
            attach_output_assembly(&L3);       // Attach output 3 assembly
          }
        
The ElemAssembly objects
          B B_assembly;
          A0 A0_assembly;
          A1 A1_assembly;
          A2 A2_assembly;
          F0 F0_assembly;
          OutputAssembly L0;
          OutputAssembly L1;
          OutputAssembly L2;
          OutputAssembly L3;
        };
        
        #endif // LIBMESH_HAVE_SLEPC && LIBMESH_HAVE_GLPK
        
        #endif
        



The source file rb_classes.h with comments:



rbOOmit is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.

You should have received a copy of the GNU Lesser General Public License along with this library; if not, write to the Free Software Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA

        #ifndef __rb_classes_h__
        #define __rb_classes_h__
        
        #if defined(LIBMESH_HAVE_SLEPC) && defined(LIBMESH_HAVE_GLPK)
        
        #include "libmesh/rb_construction.h"
        #include "libmesh/rb_scm_construction.h"
        #include "libmesh/fe_base.h"
        
Bring in bits from the libMesh namespace. Just the bits we're using, since this is a header.
        using libMesh::EquationSystems;
        using libMesh::FEMContext;
        using libMesh::RBConstruction;
        using libMesh::RBSCMConstruction;
        using libMesh::RBEvaluation;
        using libMesh::RBSCMEvaluation;
        using libMesh::RBParameters;
        using libMesh::RBThetaExpansion;
        using libMesh::RBAssemblyExpansion;
        using libMesh::AutoPtr;
        using libMesh::DirichletBoundary;
        using libMesh::Real;
        
local include
        #include "assembly.h"
        
A simple subclass of RBEvaluation. We also store the theta expansion object for the affine expansion of the PDE as a member variable.
        class SimpleRBEvaluation : public RBEvaluation
        {
        public:
        
          /**
           * Constructor. Just set the theta expansion.
           */
          SimpleRBEvaluation()
          {
            set_rb_theta_expansion(ex02_rb_theta_expansion);
          }
        
          /**
           * We override get_stability_lower_bound so that it calls rb_scm_eval to return
           * a parameter-dependent lower bound for the coercivity constant.
           */
          virtual Real get_stability_lower_bound()
          { 
            rb_scm_eval->set_parameters( get_parameters() );
            return rb_scm_eval->get_SCM_LB() ;
          }
        
          /**
           * Pointer to the SCM object that will provide our coercivity constant lower bound.
           */
          RBSCMEvaluation* rb_scm_eval;
          
          /**
           * The object that stores the "theta" expansion of the parameter dependent PDE,
           * i.e. the set of parameter-dependent functions in the affine expansion of the PDE.
           */
          Ex02RBThetaExpansion ex02_rb_theta_expansion;
        
        };
        
A simple subclass of RBConstruction, which initializes libMesh-related data such as the number of variables and their finite element type. We also store the objects that define the affine expansion of the PDE as member variables.
        class SimpleRBConstruction : public RBConstruction
        {
        public:
        
          SimpleRBConstruction (EquationSystems& es,
                                const std::string& name,
                                const unsigned int number)
          : Parent(es, name, number)
          {}
        
          /**
           * Destructor.
           */
          virtual ~SimpleRBConstruction () {}
        
          /**
           * The type of system.
           */
          typedef SimpleRBConstruction sys_type;
        
          /**
           * The type of the parent.
           */
          typedef RBConstruction Parent;
        
          /**
           * Initialize data structures.
           */
          virtual void init_data()
          {
            u_var = this->add_variable ("u", FIRST);
        
Generate a DirichletBoundary object
            dirichlet_bc = build_zero_dirichlet_boundary_object();
            
Set the Dirichet boundary IDs and the Dirichlet boundary variable numbers
            dirichlet_bc->b.insert(3);
            dirichlet_bc->variables.push_back(u_var);
            
Attach dirichlet_bc (must do this _before_ Parent::init_data)
            get_dof_map().add_dirichlet_boundary(*dirichlet_bc);
        
            Parent::init_data();
            
Set the rb_assembly_expansion for this Construction object.
            set_rb_assembly_expansion(ex02_rb_assembly_expansion);
        
We need to define an inner product matrix for this problem
            set_inner_product_assembly(ex02_rb_assembly_expansion.B_assembly);
          }
        
          /**
           * Pre-request all relevant element data. (This is not essential, but it
           * allows libMesh to cache data and hence run faster.)
           */
          virtual void init_context(FEMContext &c)
          {
For efficiency, we should prerequest all the data we will need to build the linear system before doing an element loop.
            c.element_fe_var[u_var]->get_JxW();
            c.element_fe_var[u_var]->get_phi();
            c.element_fe_var[u_var]->get_dphi();
          }
        
          /**
           * Variable number for u.
           */
          unsigned int u_var;
          
          /**
           * The object that stores the "assembly" expansion of the parameter dependent PDE,
           * i.e. the objects that define how to assemble the set of parameter-independent
           * operators in the affine expansion of the PDE.
           */
          Ex02RBAssemblyExpansion ex02_rb_assembly_expansion;
        
          /**
           * The object that defines which degrees of freedom are on a Dirichlet boundary.
           */
          AutoPtr<DirichletBoundary> dirichlet_bc;
        
        };
        
        #endif // LIBMESH_HAVE_SLEPC && LIBMESH_HAVE_GLPK
        
        #endif



The source file reduced_basis_ex2.C with comments:



        /* rbOOmit is distributed in the hope that it will be useful, */
        /* but WITHOUT ANY WARRANTY; without even the implied warranty of */
        /* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU */
        /* Lesser General Public License for more details. */
          
        /* You should have received a copy of the GNU Lesser General Public */
        /* License along with this library; if not, write to the Free Software */
        /* Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA */
        

Reduced Basis Example 2 - Successive Constraint Method



In this example we extend reduced_basis_ex1 to solve a steady convection-diffusion problem on the unit square via the Reduced Basis Method. In this case, we modify the PDE so that it no longer has a parameter-independent coercivity constant. Therefore, in order to obtain an error bound, we need to employ the Successive Constraint Method (SCM) implemented in RBSCMConstruction/RBSCMEvaluation to obtain a parameter-dependent lower bound for the coercivity constant.

The PDE being solved is div(k*grad(u)) + Beta*grad(u) = f k is the diffusion coefficient : - constant in the domain 0<=x<0.5 , its value is given by the first parameter mu[0] - constant in the domain 0.5<=x<=1 , its value is given by the second parameter mu[1] Beta is the convection velocity : - constant in the whole domain - equal to zero in the y-direction - its value in the x-direction is given by the third (and last) parameter mu[2] Boundary conditions : - dyu=0 on top and bottom - u=0 on the left side - dxu + Beta*u = 0 on the right side

C++ include files that we need
        #include <iostream>
        #include <algorithm>
        #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
        #include <cmath>
        #include <set>
        
Basic include file needed for the mesh functionality.
        #include "libmesh/libmesh.h"
        #include "libmesh/mesh.h"
        #include "libmesh/mesh_generation.h"
        #include "libmesh/exodusII_io.h"
        #include "libmesh/equation_systems.h"
        #include "libmesh/dof_map.h"
        #include "libmesh/getpot.h"
        #include "libmesh/elem.h"
        
local includes
        #include "rb_classes.h"
        #include "assembly.h"
        
Bring in everything from the libMesh namespace
        using namespace libMesh;
        
The main program.
        int main (int argc, char** argv)
        {
Initialize libMesh.
          LibMeshInit init (argc, argv);
        
This example requires SLEPc and GLPK
        #if !defined(LIBMESH_HAVE_SLEPC) || !defined(LIBMESH_HAVE_GLPK)
          libmesh_example_assert(false, "--enable-slepc --enable-glpk");
        #else
        
        #if !defined(LIBMESH_HAVE_XDR)
We need XDR support to write out reduced bases
          libmesh_example_assert(false, "--enable-xdr");
        #elif defined(LIBMESH_DEFAULT_SINGLE_PRECISION)
XDR binary support requires double precision
          libmesh_example_assert(false, "--disable-singleprecision");
        #endif
FIXME: This example currently segfaults with Trilinos?
          libmesh_example_assert(libMesh::default_solver_package() == PETSC_SOLVERS, "--enable-petsc");
        
Skip this 2D example if libMesh was compiled as 1D-only.
          libmesh_example_assert(2 <= LIBMESH_DIM, "2D support");
          
Parse the input file (reduced_basis_ex2.in) using GetPot
          std::string parameters_filename = "reduced_basis_ex2.in";
          GetPot infile(parameters_filename);
        
          unsigned int n_elem = infile("n_elem", 1);       // Determines the number of elements in the "truth" mesh
          const unsigned int dim = 2;                      // The number of spatial dimensions
        
          bool store_basis_functions = infile("store_basis_functions", true); // Do we write the RB basis functions to disk?
        
Read the "online_mode" flag from the command line
          GetPot command_line (argc, argv);
          int online_mode = 0;
          if ( command_line.search(1, "-online_mode") )
            online_mode = command_line.next(online_mode);
        
Build a mesh.
          Mesh mesh (dim);
          MeshTools::Generation::build_square (mesh,
                                               n_elem, n_elem,
                                               0., 1.,
                                               0., 1.,
                                               QUAD4);
        
Create an equation systems object.
          EquationSystems equation_systems (mesh);
        
We override RBConstruction with SimpleRBConstruction in order to specialize a few functions for this particular problem.
          SimpleRBConstruction & rb_con =
            equation_systems.add_system<SimpleRBConstruction> ("RBConvectionDiffusion");
        
Initialize the SCM Construction object
          RBSCMConstruction & rb_scm_con =
            equation_systems.add_system<RBSCMConstruction> ("RBSCMConvectionDiffusion");
          rb_scm_con.set_RB_system_name("RBConvectionDiffusion");
          rb_scm_con.add_variable("p", FIRST);
        
Initialize the data structures for the equation system.
          equation_systems.init ();
        
Print out some information about the "truth" discretization
          equation_systems.print_info();
          mesh.print_info();
        
Set parameters for the eigenvalue problems that will be solved by rb_scm_con
          equation_systems.parameters.set<unsigned int>("eigenpairs")    = 1;
          equation_systems.parameters.set<unsigned int>("basis vectors") = 3;
          equation_systems.parameters.set<unsigned int>
            ("linear solver maximum iterations") = 1000;
        
Build a new RBEvaluation object which will be used to perform Reduced Basis calculations. This is required in both the "Offline" and "Online" stages.
          SimpleRBEvaluation rb_eval;
        
We need to give the RBConstruction object a pointer to our RBEvaluation object
          rb_con.set_rb_evaluation(rb_eval);
        
We also need a SCM evaluation object to perform SCM calculations
          RBSCMEvaluation rb_scm_eval;
          rb_scm_eval.set_rb_theta_expansion( rb_eval.get_rb_theta_expansion() );
        
Tell rb_eval about rb_scm_eval
          rb_eval.rb_scm_eval = &rb_scm_eval;
          
Finally, need to give rb_scm_con and rb_eval a pointer to the SCM evaluation object, rb_scm_eval
          rb_scm_con.set_rb_scm_evaluation(rb_scm_eval);
          
          if(!online_mode) // Perform the Offline stage of the RB method
          {
Read in the data that defines this problem from the specified text file
            rb_con.process_parameters_file(parameters_filename);
            rb_scm_con.process_parameters_file(parameters_filename);
        
Print out info that describes the current setup of rb_con
            rb_con.print_info();
            rb_scm_con.print_info();
        
Prepare rb_con for the Construction stage of the RB method. This sets up the necessary data structures and performs initial assembly of the "truth" affine expansion of the PDE.
            rb_con.initialize_rb_construction();
            
Perform the SCM Greedy algorithm to derive the data required for rb_scm_eval to provide a coercivity lower bound.
            rb_scm_con.perform_SCM_greedy();
        
Compute the reduced basis space by computing "snapshots", i.e. "truth" solves, at well-chosen parameter values and employing these snapshots as basis functions.
            rb_con.train_reduced_basis();
            
Write out the data that will subsequently be required for the Evaluation stage
            rb_con.get_rb_evaluation().write_offline_data_to_files("rb_data");
            rb_scm_con.get_rb_scm_evaluation().write_offline_data_to_files("scm_data");
            
If requested, write out the RB basis functions for visualization purposes
            if(store_basis_functions)
            {
Write out the basis functions
              rb_con.get_rb_evaluation().write_out_basis_functions(rb_con,"rb_data");
            }
          }
          else // Perform the Online stage of the RB method
          {
        
Read in the reduced basis data
            rb_eval.read_offline_data_from_files("rb_data");
            rb_scm_eval.read_offline_data_from_files("scm_data");
        
Read in online_N and initialize online parameters
            unsigned int online_N = infile("online_N",1);
            Real online_mu_0 = infile("online_mu_0", 0.);
            Real online_mu_1 = infile("online_mu_1", 0.);
            Real online_mu_2 = infile("online_mu_2", 0.);
            RBParameters online_mu;
            online_mu.set_value("mu_0", online_mu_0);
            online_mu.set_value("mu_1", online_mu_1);
            online_mu.set_value("mu_2", online_mu_2);
            rb_eval.set_parameters(online_mu);
            rb_eval.print_parameters();
         
Now do the Online solve using the precomputed reduced basis
            rb_eval.rb_solve(online_N);
        
Print out outputs as well as the corresponding output error bounds.
            std::cout << "output 1, value = " << rb_eval.RB_outputs[0]
                      << ", bound = " << rb_eval.RB_output_error_bounds[0]
                      << std::endl;
            std::cout << "output 2, value = " << rb_eval.RB_outputs[1]
                      << ", bound = " << rb_eval.RB_output_error_bounds[1]
                      << std::endl;
            std::cout << "output 3, value = " << rb_eval.RB_outputs[2]
                      << ", bound = " << rb_eval.RB_output_error_bounds[2]
                      << std::endl;
            std::cout << "output 4, value = " << rb_eval.RB_outputs[3]
                      << ", bound = " << rb_eval.RB_output_error_bounds[3]
                      << std::endl << std::endl;
        
            if(store_basis_functions)
            {
Read in the basis functions
              rb_eval.read_in_basis_functions(rb_con,"rb_data");
              
Plot the solution
              rb_con.load_rb_solution();
        #ifdef LIBMESH_HAVE_EXODUS_API
              ExodusII_IO(mesh).write_equation_systems ("RB_sol.e",equation_systems);
        #endif
              
Plot the first basis function that was generated from the train_reduced_basis call in the Offline stage
              rb_con.load_basis_function(0);
        #ifdef LIBMESH_HAVE_EXODUS_API
              ExodusII_IO(mesh).write_equation_systems ("bf0.e",equation_systems);
        #endif
            }
          }
        
          return 0;
        
        #endif // LIBMESH_HAVE_SLEPC && LIBMESH_HAVE_GLPK
        }
        



The source file assembly.h without comments:

 
  #ifndef __assembly_h__
  #define __assembly_h__
  
  #if defined(LIBMESH_HAVE_SLEPC) && defined(LIBMESH_HAVE_GLPK)
  
  #include "libmesh/sparse_matrix.h"
  #include "libmesh/numeric_vector.h"
  #include "libmesh/dense_matrix.h"
  #include "libmesh/dense_vector.h"
  #include "libmesh/fe.h"
  #include "libmesh/fe_interface.h"
  #include "libmesh/fe_base.h"
  #include "libmesh/elem_assembly.h"
  #include "libmesh/quadrature_gauss.h"
  
  #include "libmesh/rb_theta.h"
  #include "libmesh/rb_assembly_expansion.h"
  
  using libMesh::ElemAssembly;
  using libMesh::FEInterface;
  using libMesh::FEMContext;
  using libMesh::Number;
  using libMesh::Point;
  using libMesh::RBTheta;
  using libMesh::Real;
  using libMesh::RealGradient;
  
  struct ThetaA0 : RBTheta { virtual Number evaluate(const RBParameters& mu) { return mu.get_value("mu_0"); } };
  struct ThetaA1 : RBTheta { virtual Number evaluate(const RBParameters& mu) { return mu.get_value("mu_1"); } };
  struct ThetaA2 : RBTheta { virtual Number evaluate(const RBParameters& mu) { return mu.get_value("mu_2"); } };
  
  struct B : ElemAssembly
  {
    virtual void interior_assembly(FEMContext &c)
    {
      const unsigned int u_var = 0;
  
      const std::vector<Real> &JxW =
        c.element_fe_var[u_var]->get_JxW();
  
      const std::vector<std::vector<Real> >& phi =
        c.element_fe_var[u_var]->get_phi();
  
      const std::vector<std::vector<RealGradient> >& dphi =
        c.element_fe_var[u_var]->get_dphi();
  
      const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
  
      unsigned int n_qpoints = (c.get_element_qrule())->n_points();
  
      for (unsigned int qp=0; qp != n_qpoints; qp++)
        for (unsigned int i=0; i != n_u_dofs; i++)
          for (unsigned int j=0; j != n_u_dofs; j++)
            c.get_elem_jacobian()(i,j) += JxW[qp] * (phi[j][qp]*phi[i][qp] + dphi[j][qp]*dphi[i][qp]);
    }
  };
  
  
  struct A0 : ElemAssembly
  {
    virtual void interior_assembly(FEMContext &c)
    {
      const unsigned int u_var = 0;
  
      const std::vector<Real> &JxW =
        c.element_fe_var[u_var]->get_JxW();
  
      const std::vector<std::vector<RealGradient> >& dphi =
        c.element_fe_var[u_var]->get_dphi();
  
      const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
  
      unsigned int n_qpoints = (c.get_element_qrule())->n_points();
  
      Real min_x=0.,max_x=0.5;
  
      Point centroid = c.elem->centroid();
      if( (min_x <= centroid(0)) && (centroid(0) < max_x) )
      for (unsigned int qp=0; qp != n_qpoints; qp++)
        for (unsigned int i=0; i != n_u_dofs; i++)
          for (unsigned int j=0; j != n_u_dofs; j++)
            c.get_elem_jacobian()(i,j) += JxW[qp] * dphi[j][qp]*dphi[i][qp];
    }
  };
  
  struct A1 : ElemAssembly
  {
    virtual void interior_assembly(FEMContext &c)
    {
      const unsigned int u_var = 0;
  
      const std::vector<Real> &JxW =
        c.element_fe_var[u_var]->get_JxW();
  
      const std::vector<std::vector<RealGradient> >& dphi =
        c.element_fe_var[u_var]->get_dphi();
  
      const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
  
      unsigned int n_qpoints = (c.get_element_qrule())->n_points();
  
      Real min_x=0.5, max_x=1.;
  
      Point centroid = c.elem->centroid();
      if( (min_x <= centroid(0)) && (centroid(0) <= max_x) )
      for (unsigned int qp=0; qp != n_qpoints; qp++)
        for (unsigned int i=0; i != n_u_dofs; i++)
          for (unsigned int j=0; j != n_u_dofs; j++)
            c.get_elem_jacobian()(i,j) += JxW[qp] * dphi[j][qp]*dphi[i][qp];
    }
  };
  
  struct A2 : ElemAssembly
  {
    virtual void interior_assembly(FEMContext &c)
    {
      const unsigned int u_var = 0;
  
      const std::vector<Real> &JxW =
        c.element_fe_var[u_var]->get_JxW();
  
      const std::vector<std::vector<Real> >& phi =
        c.element_fe_var[u_var]->get_phi();
  
      const std::vector<std::vector<RealGradient> >& dphi =
        c.element_fe_var[u_var]->get_dphi();
  
      const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
  
      unsigned int n_qpoints = (c.get_element_qrule())->n_points();
  
      for (unsigned int qp=0; qp != n_qpoints; qp++)
        for (unsigned int i=0; i != n_u_dofs; i++)
          for (unsigned int j=0; j != n_u_dofs; j++)
            c.get_elem_jacobian()(i,j) -= JxW[qp] *dphi[i][qp](0)*phi[j][qp];
    }
  };
  
  
  struct F0 : ElemAssembly
  {
    virtual void interior_assembly(FEMContext &c)
    {
      const unsigned int u_var = 0;
  
      const std::vector<Real> &JxW =
        c.element_fe_var[u_var]->get_JxW();
  
      const std::vector<std::vector<Real> >& phi =
        c.element_fe_var[u_var]->get_phi();
  
      const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
  
      unsigned int n_qpoints = (c.get_element_qrule())->n_points();
  
      for (unsigned int qp=0; qp != n_qpoints; qp++)
        for (unsigned int i=0; i != n_u_dofs; i++)
          c.get_elem_residual()(i) += JxW[qp] * ( 1.*phi[i][qp] );
    }
  };
  
  struct OutputAssembly : ElemAssembly
  {
    OutputAssembly(Real min_x_in, Real max_x_in,
                   Real min_y_in, Real max_y_in)
                  :
                  min_x(min_x_in),
                  max_x(max_x_in),
                  min_y(min_y_in),
                  max_y(max_y_in)
    {}
  
    virtual void interior_assembly(FEMContext &c)
    {
      const unsigned int u_var = 0;
  
      const std::vector<Real> &JxW =
        c.element_fe_var[u_var]->get_JxW();
  
      const std::vector<std::vector<Real> >& phi =
        c.element_fe_var[u_var]->get_phi();
  
      const unsigned int n_u_dofs = c.dof_indices_var[u_var].size();
  
      unsigned int n_qpoints = (c.get_element_qrule())->n_points();
      
      Real output_area = (max_x-min_x) * (max_y-min_y);
  
      Point centroid = c.elem->centroid();
      if( (min_x <= centroid(0)) && (centroid(0) <= max_x) &&
          (min_y <= centroid(1)) && (centroid(1) <= max_y) )
        for (unsigned int qp=0; qp != n_qpoints; qp++)
          for (unsigned int i=0; i != n_u_dofs; i++)
            c.get_elem_residual()(i) += JxW[qp] * ( 1.*phi[i][qp] ) / output_area;
    }
    
    Real min_x, max_x, min_y, max_y;
  };
  
  struct Ex02RBThetaExpansion : RBThetaExpansion
  {
  
    /**
     * Constructor.
     */
    Ex02RBThetaExpansion()
    {
      attach_A_theta(&theta_a_0);   // Attach the lhs theta
      attach_A_theta(&theta_a_1);
      attach_A_theta(&theta_a_2);
  
      attach_F_theta(&rb_theta);    // Attach the rhs theta
  
      attach_output_theta(&rb_theta); // Attach output 0 theta
      attach_output_theta(&rb_theta); // Attach output 1 theta
      attach_output_theta(&rb_theta); // Attach output 2 theta
      attach_output_theta(&rb_theta); // Attach output 3 theta
    }
  
    ThetaA0 theta_a_0;
    ThetaA1 theta_a_1;
    ThetaA2 theta_a_2;
    RBTheta rb_theta; // Default RBTheta object, just returns 1.
  };
  
  struct Ex02RBAssemblyExpansion : RBAssemblyExpansion
  {
  
    /**
     * Constructor.
     */
    Ex02RBAssemblyExpansion()
      :
      L0(0.72,0.88,0.72,0.88), // We make sure these output regions conform to the mesh
      L1(0.12,0.28,0.72,0.88),
      L2(0.12,0.28,0.12,0.28),
      L3(0.72,0.88,0.12,0.28)
    {
      attach_A_assembly(&A0_assembly); // Attach the lhs assembly
      attach_A_assembly(&A1_assembly);
      attach_A_assembly(&A2_assembly);
      
      attach_F_assembly(&F0_assembly); // Attach the rhs assembly
      
      attach_output_assembly(&L0);       // Attach output 0 assembly
      attach_output_assembly(&L1);       // Attach output 1 assembly
      attach_output_assembly(&L2);       // Attach output 2 assembly
      attach_output_assembly(&L3);       // Attach output 3 assembly
    }
  
    B B_assembly;
    A0 A0_assembly;
    A1 A1_assembly;
    A2 A2_assembly;
    F0 F0_assembly;
    OutputAssembly L0;
    OutputAssembly L1;
    OutputAssembly L2;
    OutputAssembly L3;
  };
  
  #endif // LIBMESH_HAVE_SLEPC && LIBMESH_HAVE_GLPK
  
  #endif
  



The source file rb_classes.h without comments:

 
    
    
  
  #ifndef __rb_classes_h__
  #define __rb_classes_h__
  
  #if defined(LIBMESH_HAVE_SLEPC) && defined(LIBMESH_HAVE_GLPK)
  
  #include "libmesh/rb_construction.h"
  #include "libmesh/rb_scm_construction.h"
  #include "libmesh/fe_base.h"
  
  using libMesh::EquationSystems;
  using libMesh::FEMContext;
  using libMesh::RBConstruction;
  using libMesh::RBSCMConstruction;
  using libMesh::RBEvaluation;
  using libMesh::RBSCMEvaluation;
  using libMesh::RBParameters;
  using libMesh::RBThetaExpansion;
  using libMesh::RBAssemblyExpansion;
  using libMesh::AutoPtr;
  using libMesh::DirichletBoundary;
  using libMesh::Real;
  
  #include "assembly.h"
  
  class SimpleRBEvaluation : public RBEvaluation
  {
  public:
  
    /**
     * Constructor. Just set the theta expansion.
     */
    SimpleRBEvaluation()
    {
      set_rb_theta_expansion(ex02_rb_theta_expansion);
    }
  
    /**
     * We override get_stability_lower_bound so that it calls rb_scm_eval to return
     * a parameter-dependent lower bound for the coercivity constant.
     */
    virtual Real get_stability_lower_bound()
    { 
      rb_scm_eval->set_parameters( get_parameters() );
      return rb_scm_eval->get_SCM_LB() ;
    }
  
    /**
     * Pointer to the SCM object that will provide our coercivity constant lower bound.
     */
    RBSCMEvaluation* rb_scm_eval;
    
    /**
     * The object that stores the "theta" expansion of the parameter dependent PDE,
     * i.e. the set of parameter-dependent functions in the affine expansion of the PDE.
     */
    Ex02RBThetaExpansion ex02_rb_theta_expansion;
  
  };
  
  class SimpleRBConstruction : public RBConstruction
  {
  public:
  
    SimpleRBConstruction (EquationSystems& es,
                          const std::string& name,
                          const unsigned int number)
    : Parent(es, name, number)
    {}
  
    /**
     * Destructor.
     */
    virtual ~SimpleRBConstruction () {}
  
    /**
     * The type of system.
     */
    typedef SimpleRBConstruction sys_type;
  
    /**
     * The type of the parent.
     */
    typedef RBConstruction Parent;
  
    /**
     * Initialize data structures.
     */
    virtual void init_data()
    {
      u_var = this->add_variable ("u", FIRST);
  
      dirichlet_bc = build_zero_dirichlet_boundary_object();
      
      dirichlet_bc->b.insert(3);
      dirichlet_bc->variables.push_back(u_var);
      
      get_dof_map().add_dirichlet_boundary(*dirichlet_bc);
  
      Parent::init_data();
      
      set_rb_assembly_expansion(ex02_rb_assembly_expansion);
  
      set_inner_product_assembly(ex02_rb_assembly_expansion.B_assembly);
    }
  
    /**
     * Pre-request all relevant element data. (This is not essential, but it
     * allows libMesh to cache data and hence run faster.)
     */
    virtual void init_context(FEMContext &c)
    {
      c.element_fe_var[u_var]->get_JxW();
      c.element_fe_var[u_var]->get_phi();
      c.element_fe_var[u_var]->get_dphi();
    }
  
    /**
     * Variable number for u.
     */
    unsigned int u_var;
    
    /**
     * The object that stores the "assembly" expansion of the parameter dependent PDE,
     * i.e. the objects that define how to assemble the set of parameter-independent
     * operators in the affine expansion of the PDE.
     */
    Ex02RBAssemblyExpansion ex02_rb_assembly_expansion;
  
    /**
     * The object that defines which degrees of freedom are on a Dirichlet boundary.
     */
    AutoPtr<DirichletBoundary> dirichlet_bc;
  
  };
  
  #endif // LIBMESH_HAVE_SLEPC && LIBMESH_HAVE_GLPK
  
  #endif



The source file reduced_basis_ex2.C without comments:

 
    
  /* rbOOmit is distributed in the hope that it will be useful, */
  /* but WITHOUT ANY WARRANTY; without even the implied warranty of */
  /* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU */
  /* Lesser General Public License for more details. */
    
  /* You should have received a copy of the GNU Lesser General Public */
  /* License along with this library; if not, write to the Free Software */
  /* Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA */
  
  
  
  
  #include <iostream>
  #include <algorithm>
  #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
  #include <cmath>
  #include <set>
  
  #include "libmesh/libmesh.h"
  #include "libmesh/mesh.h"
  #include "libmesh/mesh_generation.h"
  #include "libmesh/exodusII_io.h"
  #include "libmesh/equation_systems.h"
  #include "libmesh/dof_map.h"
  #include "libmesh/getpot.h"
  #include "libmesh/elem.h"
  
  #include "rb_classes.h"
  #include "assembly.h"
  
  using namespace libMesh;
  
  int main (int argc, char** argv)
  {
    LibMeshInit init (argc, argv);
  
  #if !defined(LIBMESH_HAVE_SLEPC) || !defined(LIBMESH_HAVE_GLPK)
    libmesh_example_assert(false, "--enable-slepc --enable-glpk");
  #else
  
  #if !defined(LIBMESH_HAVE_XDR)
    libmesh_example_assert(false, "--enable-xdr");
  #elif defined(LIBMESH_DEFAULT_SINGLE_PRECISION)
    libmesh_example_assert(false, "--disable-singleprecision");
  #endif
    libmesh_example_assert(libMesh::default_solver_package() == PETSC_SOLVERS, "--enable-petsc");
  
    libmesh_example_assert(2 <= LIBMESH_DIM, "2D support");
    
    std::string parameters_filename = "reduced_basis_ex2.in";
    GetPot infile(parameters_filename);
  
    unsigned int n_elem = infile("n_elem", 1);       // Determines the number of elements in the "truth" mesh
    const unsigned int dim = 2;                      // The number of spatial dimensions
  
    bool store_basis_functions = infile("store_basis_functions", true); // Do we write the RB basis functions to disk?
  
    GetPot command_line (argc, argv);
    int online_mode = 0;
    if ( command_line.search(1, "-online_mode") )
      online_mode = command_line.next(online_mode);
  
    Mesh mesh (dim);
    MeshTools::Generation::build_square (mesh,
                                         n_elem, n_elem,
                                         0., 1.,
                                         0., 1.,
                                         QUAD4);
  
    EquationSystems equation_systems (mesh);
  
    SimpleRBConstruction & rb_con =
      equation_systems.add_system<SimpleRBConstruction> ("RBConvectionDiffusion");
  
    RBSCMConstruction & rb_scm_con =
      equation_systems.add_system<RBSCMConstruction> ("RBSCMConvectionDiffusion");
    rb_scm_con.set_RB_system_name("RBConvectionDiffusion");
    rb_scm_con.add_variable("p", FIRST);
  
    equation_systems.init ();
  
    equation_systems.print_info();
    mesh.print_info();
  
    equation_systems.parameters.set<unsigned int>("eigenpairs")    = 1;
    equation_systems.parameters.set<unsigned int>("basis vectors") = 3;
    equation_systems.parameters.set<unsigned int>
      ("linear solver maximum iterations") = 1000;
  
    SimpleRBEvaluation rb_eval;
  
    rb_con.set_rb_evaluation(rb_eval);
  
    RBSCMEvaluation rb_scm_eval;
    rb_scm_eval.set_rb_theta_expansion( rb_eval.get_rb_theta_expansion() );
  
    rb_eval.rb_scm_eval = &rb_scm_eval;
    
    rb_scm_con.set_rb_scm_evaluation(rb_scm_eval);
    
    if(!online_mode) // Perform the Offline stage of the RB method
    {
      rb_con.process_parameters_file(parameters_filename);
      rb_scm_con.process_parameters_file(parameters_filename);
  
      rb_con.print_info();
      rb_scm_con.print_info();
  
      rb_con.initialize_rb_construction();
      
      rb_scm_con.perform_SCM_greedy();
  
      rb_con.train_reduced_basis();
      
      rb_con.get_rb_evaluation().write_offline_data_to_files("rb_data");
      rb_scm_con.get_rb_scm_evaluation().write_offline_data_to_files("scm_data");
      
      if(store_basis_functions)
      {
        rb_con.get_rb_evaluation().write_out_basis_functions(rb_con,"rb_data");
      }
    }
    else // Perform the Online stage of the RB method
    {
  
      rb_eval.read_offline_data_from_files("rb_data");
      rb_scm_eval.read_offline_data_from_files("scm_data");
  
      unsigned int online_N = infile("online_N",1);
      Real online_mu_0 = infile("online_mu_0", 0.);
      Real online_mu_1 = infile("online_mu_1", 0.);
      Real online_mu_2 = infile("online_mu_2", 0.);
      RBParameters online_mu;
      online_mu.set_value("mu_0", online_mu_0);
      online_mu.set_value("mu_1", online_mu_1);
      online_mu.set_value("mu_2", online_mu_2);
      rb_eval.set_parameters(online_mu);
      rb_eval.print_parameters();
   
      rb_eval.rb_solve(online_N);
  
      std::cout << "output 1, value = " << rb_eval.RB_outputs[0]
                << ", bound = " << rb_eval.RB_output_error_bounds[0]
                << std::endl;
      std::cout << "output 2, value = " << rb_eval.RB_outputs[1]
                << ", bound = " << rb_eval.RB_output_error_bounds[1]
                << std::endl;
      std::cout << "output 3, value = " << rb_eval.RB_outputs[2]
                << ", bound = " << rb_eval.RB_output_error_bounds[2]
                << std::endl;
      std::cout << "output 4, value = " << rb_eval.RB_outputs[3]
                << ", bound = " << rb_eval.RB_output_error_bounds[3]
                << std::endl << std::endl;
  
      if(store_basis_functions)
      {
        rb_eval.read_in_basis_functions(rb_con,"rb_data");
        
        rb_con.load_rb_solution();
  #ifdef LIBMESH_HAVE_EXODUS_API
        ExodusII_IO(mesh).write_equation_systems ("RB_sol.e",equation_systems);
  #endif
        
        rb_con.load_basis_function(0);
  #ifdef LIBMESH_HAVE_EXODUS_API
        ExodusII_IO(mesh).write_equation_systems ("bf0.e",equation_systems);
  #endif
      }
    }
  
    return 0;
  
  #endif // LIBMESH_HAVE_SLEPC && LIBMESH_HAVE_GLPK
  }
  



The console output of the program:

***************************************************************
* Running Example reduced_basis_ex2:
*  mpirun -np 2 example-devel -online_mode 0 -eps_type lapack -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
 
*** Warning, This code is untested, experimental, or likely to see future API changes: src/reduced_basis/rb_parametrized.C, line 41, compiled Feb  5 2013 at 13:19:15 ***
 EquationSystems
  n_systems()=2
   System #0, "RBConvectionDiffusion"
    Type "RBConstruction"
    Variables="u" 
    Finite Element Types="LAGRANGE", "JACOBI_20_00" 
    Infinite Element Mapping="CARTESIAN" 
    Approximation Orders="FIRST", "THIRD" 
    n_dofs()=676
    n_local_dofs()=354
    n_constrained_dofs()=26
    n_local_constrained_dofs()=26
    n_vectors()=1
    n_matrices()=1
    DofMap Sparsity
      Average  On-Processor Bandwidth <= 8.44083
      Average Off-Processor Bandwidth <= 0.242604
      Maximum  On-Processor Bandwidth <= 11
      Maximum Off-Processor Bandwidth <= 4
    DofMap Constraints
      Number of DoF Constraints = 26
      Average DoF Constraint Length= 0
      Number of Node Constraints = 0
   System #1, "RBSCMConvectionDiffusion"
    Type "Eigen"
    Variables="p" 
    Finite Element Types="LAGRANGE", "JACOBI_20_00" 
    Infinite Element Mapping="CARTESIAN" 
    Approximation Orders="FIRST", "THIRD" 
    n_dofs()=676
    n_local_dofs()=354
    n_constrained_dofs()=0
    n_local_constrained_dofs()=0
    n_vectors()=0
    n_matrices()=2
    DofMap Sparsity
      Average  On-Processor Bandwidth <= 8.44083
      Average Off-Processor Bandwidth <= 0.242604
      Maximum  On-Processor Bandwidth <= 11
      Maximum Off-Processor Bandwidth <= 4
    DofMap Constraints
      Number of DoF Constraints = 0
      Number of Node Constraints = 0

 Mesh Information:
  mesh_dimension()=2
  spatial_dimension()=3
  n_nodes()=676
    n_local_nodes()=354
  n_elem()=625
    n_local_elem()=313
    n_active_elem()=625
  n_subdomains()=1
  n_partitions()=2
  n_processors()=2
  n_threads()=1
  processor_id()=0

Initializing training parameters with random training set...
Parameter mu_0: log scaling = 1
Parameter mu_1: log scaling = 1
Parameter mu_2: log scaling = 1

Initializing training parameters with random training set...
Parameter mu_0: log scaling = 1
Parameter mu_1: log scaling = 1
Parameter mu_2: log scaling = 1


RBConstruction parameters:
system name: RBConvectionDiffusion
constrained_problem: 0
Nmax: 10
Basis training error tolerance: 1e-05
Aq operators attached: 3
Fq functions attached: 1
n_outputs: 4
output 0, Q_l = 1
output 1, Q_l = 1
output 2, Q_l = 1
output 3, Q_l = 1
Number of parameters: 3
Parameter mu_0: Min = 0.1, Max = 1, value = 0.2
Parameter mu_1: Min = 0.1, Max = 1, value = 0.7
Parameter mu_2: Min = 0.01, Max = 0.1, value = 0.1
n_training_samples: 100
single-matrix mode? 0
reuse preconditioner? 1
use a relative error bound in greedy? 0
write out data during basis training? 0
quiet mode? 1


RBSCMConstruction parameters:
system name: RBSCMConvectionDiffusion
SCM Greedy tolerance: 0.1
A_q operators attached: 3
Number of parameters: 3
Parameter mu_0: Min = 0.1, Max = 1, value = 0.2
Parameter mu_1: Min = 0.1, Max = 1, value = 0.7
Parameter mu_2: Min = 0.01, Max = 0.1, value = 0.1
n_training_samples: 100


B_min(0) = -1.56354e-15
B_max(0) = 0.999932

B_min(1) = -9.4918e-16
B_max(1) = 0.999933

B_min(2) = -0.380786
B_max(2) = 4.96715e-17

SCM: Added mu = (0.69213,0.335734,0.0303732)

Stability constant for C_J(0) = 0.28917

SCM iteration 0, max_SCM_error = 0.963081

SCM: Added mu = (0.143443,0.852727,0.0965499)

-----------------------------------


Stability constant for C_J(1) = 0.0891787

SCM iteration 1, max_SCM_error = 0.611379

SCM: Added mu = (0.126756,0.166705,0.088367)

-----------------------------------


Stability constant for C_J(2) = 0.0686695

SCM iteration 2, max_SCM_error = 0.309361

SCM: Added mu = (0.81592,0.776663,0.0140495)

-----------------------------------


Stability constant for C_J(3) = 0.569835

SCM iteration 3, max_SCM_error = 0.145772

SCM: Added mu = (0.69586,0.132675,0.0642938)

-----------------------------------


Stability constant for C_J(4) = 0.104137

SCM iteration 4, max_SCM_error = 0.135439

SCM: Added mu = (0.109383,0.235457,0.0903173)

-----------------------------------


Stability constant for C_J(5) = 0.0611163

SCM iteration 5, max_SCM_error = 0.12401

SCM: Added mu = (0.628695,0.460497,0.0832009)

-----------------------------------


Stability constant for C_J(6) = 0.36118

SCM iteration 6, max_SCM_error = 0.066033

SCM tolerance of 0.1 reached.

Compute output dual inner products
output_dual_innerprods[0][0] = 0.839698
output_dual_innerprods[1][0] = 0.318298
output_dual_innerprods[2][0] = 0.318298
output_dual_innerprods[3][0] = 0.839698

---- Performing Greedy basis enrichment ----

---- Basis dimension: 0 ----
Performing RB solves on training set
Maximum (absolute) error bound is 7.98785

Performing truth solve at parameter:
mu_0: 1.093829e-01
mu_1: 2.354570e-01
mu_2: 9.031732e-02

Enriching the RB space
Updating RB matrices
Updating RB residual terms

---- Basis dimension: 1 ----
Performing RB solves on training set
Maximum (absolute) error bound is 6.12771

Performing truth solve at parameter:
mu_0: 1.434429e-01
mu_1: 8.527266e-01
mu_2: 9.654986e-02

Enriching the RB space
Updating RB matrices
Updating RB residual terms

---- Basis dimension: 2 ----
Performing RB solves on training set
Maximum (absolute) error bound is 2.21489

Performing truth solve at parameter:
mu_0: 5.864481e-01
mu_1: 1.047184e-01
mu_2: 2.703668e-02

Enriching the RB space
Updating RB matrices
Updating RB residual terms

---- Basis dimension: 3 ----
Performing RB solves on training set
Maximum (absolute) error bound is 0.689732

Performing truth solve at parameter:
mu_0: 1.038247e-01
mu_1: 3.400998e-01
mu_2: 1.268152e-02

Enriching the RB space
Updating RB matrices
Updating RB residual terms

---- Basis dimension: 4 ----
Performing RB solves on training set
Maximum (absolute) error bound is 0.0271113

Performing truth solve at parameter:
mu_0: 3.867963e-01
mu_1: 1.086049e-01
mu_2: 5.303670e-02

Enriching the RB space
Updating RB matrices
Updating RB residual terms

---- Basis dimension: 5 ----
Performing RB solves on training set
Maximum (absolute) error bound is 0.015854

Performing truth solve at parameter:
mu_0: 1.267559e-01
mu_1: 1.667054e-01
mu_2: 8.836701e-02

Enriching the RB space
Updating RB matrices
Updating RB residual terms

---- Basis dimension: 6 ----
Performing RB solves on training set
Maximum (absolute) error bound is 0.000187333

Performing truth solve at parameter:
mu_0: 1.174236e-01
mu_1: 4.833019e-01
mu_2: 4.534858e-02

Enriching the RB space
Updating RB matrices
Updating RB residual terms

---- Basis dimension: 7 ----
Performing RB solves on training set
Maximum (absolute) error bound is 0.000170426

Performing truth solve at parameter:
mu_0: 6.367141e-01
mu_1: 1.410916e-01
mu_2: 1.142793e-02

Enriching the RB space
Updating RB matrices
Updating RB residual terms

---- Basis dimension: 8 ----
Performing RB solves on training set
Maximum (absolute) error bound is 8.34799e-07

Specified error tolerance reached.
************************************************************************************************************************
***             WIDEN YOUR WINDOW TO 120 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document            ***
************************************************************************************************************************

---------------------------------------------- PETSc Performance Summary: ----------------------------------------------

/workspace/libmesh/examples/reduced_basis/reduced_basis_ex2/.libs/lt-example-devel on a intel-12. named hbar.ices.utexas.edu with 2 processors, by benkirk Tue Feb  5 13:46:46 2013
Using Petsc Release Version 3.3.0, Patch 2, Fri Jul 13 15:42:00 CDT 2012 

                         Max       Max/Min        Avg      Total 
Time (sec):           4.743e+00      1.00000   4.743e+00
Objects:              9.549e+03      1.00000   9.549e+03
Flops:                5.530e+07      1.11724   5.240e+07  1.048e+08
Flops/sec:            1.166e+07      1.11724   1.105e+07  2.210e+07
MPI Messages:         2.104e+03      1.00000   2.104e+03  4.209e+03
MPI Message Lengths:  5.985e+05      1.00000   2.844e+02  1.197e+06
MPI Reductions:       2.446e+04      1.00000

Flop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)
                            e.g., VecAXPY() for real vectors of length N --> 2N flops
                            and VecAXPY() for complex vectors of length N --> 8N flops

Summary of Stages:   ----- Time ------  ----- Flops -----  --- Messages ---  -- Message Lengths --  -- Reductions --
                        Avg     %Total     Avg     %Total   counts   %Total     Avg         %Total   counts   %Total 
 0:      Main Stage: 4.7426e+00 100.0%  1.0480e+08 100.0%  4.209e+03 100.0%  2.844e+02      100.0%  2.445e+04 100.0% 

------------------------------------------------------------------------------------------------------------------------
See the 'Profiling' chapter of the users' manual for details on interpreting output.
Phase summary info:
   Count: number of times phase was executed
   Time and Flops: Max - maximum over all processors
                   Ratio - ratio of maximum to minimum over all processors
   Mess: number of messages sent
   Avg. len: average message length
   Reduct: number of global reductions
   Global: entire computation
   Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().
      %T - percent time in this phase         %f - percent flops in this phase
      %M - percent messages in this phase     %L - percent message lengths in this phase
      %R - percent reductions in this phase
   Total Mflop/s: 10e-6 * (sum of flops over all processors)/(max time over all processors)
------------------------------------------------------------------------------------------------------------------------
Event                Count      Time (sec)     Flops                             --- Global ---  --- Stage ---   Total
                   Max Ratio  Max     Ratio   Max  Ratio  Mess   Avg len Reduct  %T %f %M %L %R  %T %f %M %L %R Mflop/s
------------------------------------------------------------------------------------------------------------------------

--- Event Stage 0: Main Stage

STSetUp               13 1.0 2.1243e-03 2.7 0.00e+00 0.0 0.0e+00 0.0e+00 2.6e+01  0  0  0  0  0   0  0  0  0  0     0
EPSSetUp              13 1.0 1.0178e-01 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 1.7e+04  2  0  0  0 71   2  0  0  0 71     0
EPSSolve              13 1.0 1.5630e+00 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 33  0  0  0  0  33  0  0  0  0     0
DSSolve               13 1.0 1.5382e+00 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 32  0  0  0  0  32  0  0  0  0     0
DSVectors             13 1.0 5.3697e-03 1.2 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
DSOther               13 1.0 1.1403e-02 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecDot               749 1.0 2.3134e-03 1.4 5.30e+05 1.1 0.0e+00 0.0e+00 7.5e+02  0  1  0  0  3   0  1  0  0  3   437
VecMDot             1113 1.0 6.8488e-03 1.3 1.21e+07 1.1 0.0e+00 0.0e+00 1.1e+03  0 22  0  0  5   0 22  0  0  5  3375
VecNorm             1192 1.0 2.4030e-03 1.1 8.44e+05 1.1 0.0e+00 0.0e+00 1.2e+03  0  2  0  0  5   0  2  0  0  5   671
VecScale            1195 1.0 4.6182e-04 1.0 4.20e+05 1.1 0.0e+00 0.0e+00 0.0e+00  0  1  0  0  0   0  1  0  0  0  1738
VecCopy              178 1.0 7.4148e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecSet              2360 1.0 7.2479e-04 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecAXPY              125 1.0 9.2030e-05 1.3 8.85e+04 1.1 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0  1836
VecMAXPY            1155 1.0 3.5090e-03 1.1 1.29e+07 1.1 0.0e+00 0.0e+00 0.0e+00  0 24  0  0  0   0 24  0  0  0  7026
VecAssemblyBegin     244 1.0 1.9443e-03 1.1 0.00e+00 0.0 1.0e+01 3.7e+02 7.3e+02  0  0  0  0  3   0  0  0  0  3     0
VecAssemblyEnd       244 1.0 1.1063e-04 1.2 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecScatterBegin     1890 1.0 2.0957e-03 1.1 0.00e+00 0.0 3.8e+03 2.6e+02 0.0e+00  0  0 90 82  0   0  0 90 82  0     0
VecScatterEnd       1890 1.0 2.3348e-03 1.5 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecNormalize        1155 1.0 3.0959e-03 1.1 1.23e+06 1.1 0.0e+00 0.0e+00 1.2e+03  0  2  0  0  5   0  2  0  0  5   757
MatMult             1155 1.0 7.5083e-03 1.0 6.60e+06 1.1 2.3e+03 2.6e+02 0.0e+00  0 12 55 49  0   0 12 55 49  0  1673
MatMultAdd           685 1.0 4.0443e-03 1.1 4.15e+06 1.1 1.4e+03 2.6e+02 0.0e+00  0  8 33 29  0   0  8 33 29  0  1957
MatSolve            1192 1.0 7.7748e-03 1.1 1.65e+07 1.2 0.0e+00 0.0e+00 0.0e+00  0 29  0  0  0   0 29  0  0  0  3951
MatLUFactorNum        18 1.0 1.8251e-03 1.1 1.18e+06 1.1 0.0e+00 0.0e+00 0.0e+00  0  2  0  0  0   0  2  0  0  0  1229
MatILUFactorSym       18 1.0 3.8581e-03 1.2 0.00e+00 0.0 0.0e+00 0.0e+00 5.4e+01  0  0  0  0  0   0  0  0  0  0     0
MatConvert            26 1.0 1.8203e-02 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 5.2e+01  0  0  0  0  0   0  0  0  0  0     0
MatAssemblyBegin    1044 1.0 3.1349e-02 1.1 0.00e+00 0.0 1.5e+02 1.3e+03 1.9e+03  1  0  3 16  8   1  0  3 16  8     0
MatAssemblyEnd      1044 1.0 1.0674e-02 1.1 0.00e+00 0.0 2.6e+02 6.6e+01 5.4e+02  0  0  6  1  2   0  0  6  1  2     0
MatGetRow          41680 1.1 3.0260e-03 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
MatGetRowIJ           18 1.0 8.1062e-06 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
MatGetSubMatrice      26 1.0 4.1554e-03 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 2.1e+02  0  0  0  0  1   0  0  0  0  1     0
MatGetOrdering        18 1.0 3.2020e-04 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 3.6e+01  0  0  0  0  0   0  0  0  0  0     0
MatZeroEntries        68 1.0 2.0003e-04 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
MatAXPY               35 1.0 1.5570e-02 1.0 0.00e+00 0.0 1.4e+02 6.6e+01 5.6e+02  0  0  3  1  2   0  0  3  1  2     0
KSPGMRESOrthog      1113 1.0 1.0626e-02 1.2 2.42e+07 1.1 0.0e+00 0.0e+00 1.1e+03  0 44  0  0  5   0 44  0  0  5  4354
KSPSetUp              68 1.0 1.9646e-04 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
KSPSolve              37 1.0 4.3604e-02 1.0 5.06e+07 1.1 2.3e+03 2.6e+02 2.4e+03  1 91 55 49 10   1 91 55 49 10  2197
PCSetUp               49 1.0 7.9837e-03 1.1 1.18e+06 1.1 0.0e+00 0.0e+00 9.2e+01  0  2  0  0  0   0  2  0  0  0   281
PCSetUpOnBlocks       37 1.0 7.3440e-03 1.1 1.18e+06 1.1 0.0e+00 0.0e+00 9.0e+01  0  2  0  0  0   0  2  0  0  0   305
PCApply             1192 1.0 1.4040e-02 1.1 1.65e+07 1.2 0.0e+00 0.0e+00 0.0e+00  0 29  0  0  0   0 29  0  0  0  2188
------------------------------------------------------------------------------------------------------------------------

Memory usage is given in bytes:

Object Type          Creations   Destructions     Memory  Descendants' Mem.
Reports information only for process 0.

--- Event Stage 0: Main Stage

           Container    13             13         7124     0
  Spectral Transform     1              1          736     0
 Eigenproblem Solver     1              1       305340     0
       Inner product     1              1          624     0
       Direct solver     1              1    264249620     0
              Vector  8802           8802     13717920     0
      Vector Scatter    70             70        72520     0
           Index Set   324            324       274272     0
   IS L to G Mapping     2              2         1128     0
              Matrix   326            326     96402744     0
         PetscRandom     1              1          608     0
       Krylov Solver     3              3        20432     0
      Preconditioner     3              3         2640     0
              Viewer     1              0            0     0
========================================================================================================================
Average time to get PetscTime(): 0
Average time for MPI_Barrier(): 8.10623e-07
Average time for zero size MPI_Send(): 1.5974e-05
#PETSc Option Table entries:
-eps_type lapack
-ksp_right_pc
-log_summary
-online_mode 0
-pc_type bjacobi
-sub_pc_factor_levels 4
-sub_pc_factor_zeropivot 0
-sub_pc_type ilu
#End of PETSc Option Table entries
Compiled without FORTRAN kernels
Compiled with full precision matrices (default)
sizeof(short) 2 sizeof(int) 4 sizeof(long) 8 sizeof(void*) 8 sizeof(PetscScalar) 8 sizeof(PetscInt) 4
Configure run at: Thu Nov  8 11:21:02 2012
Configure options: --with-debugging=false --COPTFLAGS=-O3 --CXXOPTFLAGS=-O3 --FOPTFLAGS=-O3 --with-clanguage=C++ --with-shared-libraries=1 --with-mpi-dir=/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1 --with-mumps=true --download-mumps=1 --with-metis=true --download-metis=1 --with-parmetis=true --download-parmetis=1 --with-superlu=true --download-superlu=1 --with-superludir=true --download-superlu_dist=1 --with-blacs=true --download-blacs=1 --with-scalapack=true --download-scalapack=1 --with-hypre=true --download-hypre=1 --with-blas-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_intel_lp64.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_sequential.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_core.so]" --with-lapack-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_lapack95_lp64.a]"
-----------------------------------------
Libraries compiled on Thu Nov  8 11:21:02 2012 on daedalus.ices.utexas.edu 
Machine characteristics: Linux-2.6.32-279.1.1.el6.x86_64-x86_64-with-redhat-6.3-Carbon
Using PETSc directory: /opt/apps/ossw/libraries/petsc/petsc-3.3-p2
Using PETSc arch: intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt
-----------------------------------------

Using C compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx  -wd1572 -O3   -fPIC   ${COPTFLAGS} ${CFLAGS}
Using Fortran compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90  -fPIC -O3   ${FOPTFLAGS} ${FFLAGS} 
-----------------------------------------

Using include paths: -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/include
-----------------------------------------

Using C linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx
Using Fortran linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90
Using libraries: -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lpetsc -lX11 -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lcmumps -ldmumps -lsmumps -lzmumps -lmumps_common -lpord -lHYPRE -lpthread -lsuperlu_dist_3.0 -lparmetis -lmetis -lscalapack -lblacs -lsuperlu_4.3 -Wl,-rpath,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -L/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -lmkl_lapack95_lp64 -lmkl_intel_lp64 -lmkl_sequential -lmkl_core -Wl,-rpath,/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -L/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -Wl,-rpath,/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -L/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -Wl,-rpath,/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -L/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -lmpichf90 -lifport -lifcore -lm -lm -lmpichcxx -ldl -lmpich -lopa -lmpl -lrt -lpthread -limf -lsvml -lipgo -ldecimal -lcilkrts -lstdc++ -lgcc_s -lirc -lirc_s -ldl 
-----------------------------------------


 ----------------------------------------------------------------------------------------------------------------------
| Processor id:   0                                                                                                    |
| Num Processors: 2                                                                                                    |
| Time:           Tue Feb  5 13:46:46 2013                                                                             |
| OS:             Linux                                                                                                |
| HostName:       hbar.ices.utexas.edu                                                                                 |
| OS Release:     2.6.32-279.1.1.el6.x86_64                                                                            |
| OS Version:     #1 SMP Tue Jul 10 11:24:23 CDT 2012                                                                  |
| Machine:        x86_64                                                                                               |
| Username:       benkirk                                                                                              |
| Configuration:  ./configure  '--enable-everything'                                                                   |
|  '--prefix=/workspace/libmesh/install'                                                                               |
|  'CXX=mpicxx'                                                                                                        |
|  'CC=mpicc'                                                                                                          |
|  'F77=mpif77'                                                                                                        |
|  'FC=mpif90'                                                                                                         |
|  'PETSC_DIR=/opt/apps/ossw/libraries/petsc/petsc-3.3-p2'                                                             |
|  'PETSC_ARCH=intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt'                                                |
|  'SLEPC_DIR=/opt/apps/ossw/libraries/slepc/slepc-3.3-p2-petsc-3.3-p2-cxx-opt'                                        |
|  'TRILINOS_DIR=/opt/apps/ossw/libraries/trilinos/trilinos-10.12.2/sl6/intel-12.1/mpich2-1.4.1p1/mkl-intel-10.3.12.361'|
|  'VTK_DIR=/opt/apps/ossw/libraries/vtk/vtk-5.10.0/sl6/intel-12.1'                                                    |
 ----------------------------------------------------------------------------------------------------------------------
 --------------------------------------------------------------------------------------------------------------------
| libMesh Performance: Alive time=4.75649, Active time=4.69616                                                       |
 --------------------------------------------------------------------------------------------------------------------
| Event                                  nCalls    Total Time  Avg Time    Total Time  Avg Time    % of Active Time  |
|                                                  w/o Sub     w/o Sub     With Sub    With Sub    w/o S    With S   |
|--------------------------------------------------------------------------------------------------------------------|
|                                                                                                                    |
|                                                                                                                    |
| CondensedEigenSystem                                                                                               |
|   get_eigenpair()                      13        0.0020      0.000156    0.0247      0.001901    0.04     0.53     |
|   solve()                              13        0.0136      0.001048    1.6729      0.128686    0.29     35.62    |
|                                                                                                                    |
| DofMap                                                                                                             |
|   add_neighbors_to_send_list()         2         0.0144      0.007188    0.0177      0.008855    0.31     0.38     |
|   build_constraint_matrix()            16902     0.0630      0.000004    0.0630      0.000004    1.34     1.34     |
|   build_sparsity()                     2         0.0092      0.004588    0.0274      0.013708    0.20     0.58     |
|   cnstrn_elem_mat_vec()                16902     0.0411      0.000002    0.0411      0.000002    0.88     0.88     |
|   create_dof_constraints()             2         0.0065      0.003253    0.0240      0.011996    0.14     0.51     |
|   distribute_dofs()                    2         0.0135      0.006746    0.0369      0.018429    0.29     0.78     |
|   dof_indices()                        35179     0.9677      0.000028    0.9677      0.000028    20.61    20.61    |
|   prepare_send_list()                  2         0.0001      0.000029    0.0001      0.000029    0.00     0.00     |
|   reinit()                             2         0.0221      0.011039    0.0221      0.011039    0.47     0.47     |
|                                                                                                                    |
| FE                                                                                                                 |
|   compute_shape_functions()            39096     0.5215      0.000013    0.5215      0.000013    11.10    11.10    |
|   init_shape_functions()               5400      0.0730      0.000014    0.0730      0.000014    1.55     1.55     |
|   inverse_map()                        10584     0.0531      0.000005    0.0531      0.000005    1.13     1.13     |
|                                                                                                                    |
| FEMap                                                                                                              |
|   compute_affine_map()                 39096     0.2107      0.000005    0.2107      0.000005    4.49     4.49     |
|   compute_face_map()                   5292      0.0539      0.000010    0.1085      0.000021    1.15     2.31     |
|   init_face_shape_functions()          108       0.0010      0.000009    0.0010      0.000009    0.02     0.02     |
|   init_reference_to_physical_map()     5400      0.0433      0.000008    0.0433      0.000008    0.92     0.92     |
|                                                                                                                    |
| Mesh                                                                                                               |
|   find_neighbors()                     1         0.0072      0.007192    0.0077      0.007680    0.15     0.16     |
|   renumber_nodes_and_elem()            2         0.0003      0.000160    0.0003      0.000160    0.01     0.01     |
|                                                                                                                    |
| MeshCommunication                                                                                                  |
|   assign_global_indices()              1         0.0162      0.016173    0.0166      0.016637    0.34     0.35     |
|   compute_hilbert_indices()            2         0.0057      0.002859    0.0057      0.002859    0.12     0.12     |
|   find_global_indices()                2         0.0021      0.001055    0.0091      0.004566    0.04     0.19     |
|   parallel_sort()                      2         0.0009      0.000442    0.0010      0.000500    0.02     0.02     |
|                                                                                                                    |
| MeshTools::Generation                                                                                              |
|   build_cube()                         1         0.0015      0.001527    0.0015      0.001527    0.03     0.03     |
|                                                                                                                    |
| MetisPartitioner                                                                                                   |
|   partition()                          1         0.0065      0.006492    0.0107      0.010708    0.14     0.23     |
|                                                                                                                    |
| Parallel                                                                                                           |
|   allgather()                          17        0.0012      0.000068    0.0012      0.000070    0.02     0.03     |
|   barrier()                            1         0.0000      0.000014    0.0000      0.000014    0.00     0.00     |
|   broadcast()                          34        0.0002      0.000005    0.0002      0.000005    0.00     0.00     |
|   max(bool)                            7         0.0000      0.000002    0.0000      0.000002    0.00     0.00     |
|   max(scalar)                          167       0.0005      0.000003    0.0005      0.000003    0.01     0.01     |
|   max(vector)                          37        0.0002      0.000007    0.0005      0.000014    0.01     0.01     |
|   maxloc(scalar)                       16        0.0012      0.000077    0.0012      0.000077    0.03     0.03     |
|   min(bool)                            183       0.0004      0.000002    0.0004      0.000002    0.01     0.01     |
|   min(scalar)                          149       0.0022      0.000014    0.0022      0.000014    0.05     0.05     |
|   min(vector)                          37        0.0003      0.000008    0.0006      0.000016    0.01     0.01     |
|   probe()                              38        0.0003      0.000007    0.0003      0.000007    0.01     0.01     |
|   receive()                            30        0.0002      0.000007    0.0005      0.000015    0.00     0.01     |
|   send()                               26        0.0001      0.000003    0.0001      0.000003    0.00     0.00     |
|   send_receive()                       30        0.0003      0.000010    0.0008      0.000027    0.01     0.02     |
|   sum()                                56        0.0229      0.000410    0.0230      0.000411    0.49     0.49     |
|                                                                                                                    |
| Parallel::Request                                                                                                  |
|   wait()                               26        0.0001      0.000002    0.0001      0.000002    0.00     0.00     |
|                                                                                                                    |
| Partitioner                                                                                                        |
|   set_node_processor_ids()             1         0.0008      0.000788    0.0009      0.000865    0.02     0.02     |
|   set_parent_processor_ids()           1         0.0005      0.000542    0.0005      0.000542    0.01     0.01     |
|                                                                                                                    |
| PetscLinearSolver                                                                                                  |
|   solve()                              37        0.0457      0.001236    0.0457      0.001236    0.97     0.97     |
|                                                                                                                    |
| RBConstruction                                                                                                     |
|   add_scaled_Aq()                      45        0.0032      0.000071    2.2012      0.048915    0.07     46.87    |
|   add_scaled_matrix_and_vector()       54        0.6126      0.011345    2.6285      0.048676    13.04    55.97    |
|   clear()                              1         0.0005      0.000500    0.0005      0.000500    0.01     0.01     |
|   compute_Fq_representor_innerprods()  1         0.0006      0.000554    0.0020      0.001964    0.01     0.04     |
|   compute_max_error_bound()            9         0.0074      0.000818    0.0943      0.010475    0.16     2.01     |
|   compute_output_dual_innerprods()     1         0.0012      0.001242    0.0067      0.006653    0.03     0.14     |
|   enrich_RB_space()                    8         0.0018      0.000228    0.0018      0.000228    0.04     0.04     |
|   train_reduced_basis()                1         0.0009      0.000871    0.1764      0.176362    0.02     3.76     |
|   truth_assembly()                     8         0.0116      0.001448    0.0116      0.001448    0.25     0.25     |
|   truth_solve()                        8         0.0006      0.000070    0.0236      0.002950    0.01     0.50     |
|   update_RB_system_matrices()          8         0.0051      0.000632    0.0051      0.000632    0.11     0.11     |
|   update_residual_terms()              8         0.0143      0.001792    0.0418      0.005223    0.31     0.89     |
|                                                                                                                    |
| RBEvaluation                                                                                                       |
|   clear()                              1         0.0000      0.000050    0.0000      0.000050    0.00     0.00     |
|   compute_residual_dual_norm()         450       0.0436      0.000097    0.0436      0.000097    0.93     0.93     |
|   rb_solve()                           450       0.0149      0.000033    0.0863      0.000192    0.32     1.84     |
|   resize_data_structures()             1         0.0003      0.000321    0.0003      0.000321    0.01     0.01     |
|   write_offline_data_to_files()        1         0.0008      0.000807    0.0008      0.000807    0.02     0.02     |
|   write_out_basis_functions()          1         0.0001      0.000109    0.0513      0.051253    0.00     1.09     |
|   write_out_vectors()                  1         0.0343      0.034307    0.0511      0.051142    0.73     1.09     |
|                                                                                                                    |
| RBSCMConstruction                                                                                                  |
|   add_scaled_symm_Aq()                 45        0.0002      0.000005    2.2014      0.048920    0.00     46.88    |
|   compute_SCM_bounding_box()           1         0.0004      0.000387    0.8470      0.847050    0.01     18.04    |
|   compute_SCM_bounds_on_training_set() 7         0.0054      0.000767    0.0248      0.003549    0.11     0.53     |
|   enrich_C_J()                         7         0.0007      0.000093    0.0008      0.000113    0.01     0.02     |
|   evaluate_stability_constant()        7         0.0019      0.000275    3.0543      0.436328    0.04     65.04    |
|   perform_SCM_greedy()                 1         0.0012      0.001193    3.9282      3.928174    0.03     83.65    |
|                                                                                                                    |
| RBSCMEvaluation                                                                                                    |
|   get_SCM_LB()                         800       0.0433      0.000054    0.0433      0.000054    0.92     0.92     |
|   get_SCM_UB()                         350       0.0029      0.000008    0.0029      0.000008    0.06     0.06     |
|   write_offline_data_to_files()        1         0.0002      0.000153    0.0002      0.000153    0.00     0.00     |
|                                                                                                                    |
| SlepcEigenSolver                                                                                                   |
|   solve_generalized()                  13        1.6593      0.127638    1.6593      0.127638    35.33    35.33    |
 --------------------------------------------------------------------------------------------------------------------
| Totals:                                177191    4.6962                                          100.00            |
 --------------------------------------------------------------------------------------------------------------------

 
***************************************************************
* Done Running Example reduced_basis_ex2:
*  mpirun -np 2 example-devel -online_mode 0 -eps_type lapack -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
***************************************************************
* Running Example reduced_basis_ex2:
*  mpirun -np 2 example-devel -online_mode 1 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
 
*** Warning, This code is untested, experimental, or likely to see future API changes: src/reduced_basis/rb_parametrized.C, line 41, compiled Feb  5 2013 at 13:19:15 ***
 EquationSystems
  n_systems()=2
   System #0, "RBConvectionDiffusion"
    Type "RBConstruction"
    Variables="u" 
    Finite Element Types="LAGRANGE", "JACOBI_20_00" 
    Infinite Element Mapping="CARTESIAN" 
    Approximation Orders="FIRST", "THIRD" 
    n_dofs()=676
    n_local_dofs()=354
    n_constrained_dofs()=26
    n_local_constrained_dofs()=26
    n_vectors()=1
    n_matrices()=1
    DofMap Sparsity
      Average  On-Processor Bandwidth <= 8.44083
      Average Off-Processor Bandwidth <= 0.242604
      Maximum  On-Processor Bandwidth <= 11
      Maximum Off-Processor Bandwidth <= 4
    DofMap Constraints
      Number of DoF Constraints = 26
      Average DoF Constraint Length= 0
      Number of Node Constraints = 0
   System #1, "RBSCMConvectionDiffusion"
    Type "Eigen"
    Variables="p" 
    Finite Element Types="LAGRANGE", "JACOBI_20_00" 
    Infinite Element Mapping="CARTESIAN" 
    Approximation Orders="FIRST", "THIRD" 
    n_dofs()=676
    n_local_dofs()=354
    n_constrained_dofs()=0
    n_local_constrained_dofs()=0
    n_vectors()=0
    n_matrices()=2
    DofMap Sparsity
      Average  On-Processor Bandwidth <= 8.44083
      Average Off-Processor Bandwidth <= 0.242604
      Maximum  On-Processor Bandwidth <= 11
      Maximum Off-Processor Bandwidth <= 4
    DofMap Constraints
      Number of DoF Constraints = 0
      Number of Node Constraints = 0

 Mesh Information:
  mesh_dimension()=2
  spatial_dimension()=3
  n_nodes()=676
    n_local_nodes()=354
  n_elem()=625
    n_local_elem()=313
    n_active_elem()=625
  n_subdomains()=1
  n_partitions()=2
  n_processors()=2
  n_threads()=1
  processor_id()=0

mu_0: 2.000000e-01
mu_1: 7.000000e-01
mu_2: 1.000000e-01

output 1, value = 2.35239, bound = 7.99528e-05
output 2, value = 0.944936, bound = 4.92254e-05
output 3, value = 0.944936, bound = 4.92254e-05
output 4, value = 2.35239, bound = 7.99528e-05

************************************************************************************************************************
***             WIDEN YOUR WINDOW TO 120 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document            ***
************************************************************************************************************************

---------------------------------------------- PETSc Performance Summary: ----------------------------------------------

/workspace/libmesh/examples/reduced_basis/reduced_basis_ex2/.libs/lt-example-devel on a intel-12. named hbar.ices.utexas.edu with 2 processors, by benkirk Tue Feb  5 13:46:47 2013
Using Petsc Release Version 3.3.0, Patch 2, Fri Jul 13 15:42:00 CDT 2012 

                         Max       Max/Min        Avg      Total 
Time (sec):           2.355e-01      1.00000   2.355e-01
Objects:              3.500e+01      1.00000   3.500e+01
Flops:                4.248e+03      1.09938   4.056e+03  8.112e+03
Flops/sec:            1.804e+04      1.09938   1.722e+04  3.445e+04
MPI Messages:         6.000e+00      1.00000   6.000e+00  1.200e+01
MPI Message Lengths:  1.160e+03      1.00000   1.933e+02  2.320e+03
MPI Reductions:       8.100e+01      1.00000

Flop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)
                            e.g., VecAXPY() for real vectors of length N --> 2N flops
                            and VecAXPY() for complex vectors of length N --> 8N flops

Summary of Stages:   ----- Time ------  ----- Flops -----  --- Messages ---  -- Message Lengths --  -- Reductions --
                        Avg     %Total     Avg     %Total   counts   %Total     Avg         %Total   counts   %Total 
 0:      Main Stage: 2.3547e-01 100.0%  8.1120e+03 100.0%  1.200e+01 100.0%  1.933e+02      100.0%  8.000e+01  98.8% 

------------------------------------------------------------------------------------------------------------------------
See the 'Profiling' chapter of the users' manual for details on interpreting output.
Phase summary info:
   Count: number of times phase was executed
   Time and Flops: Max - maximum over all processors
                   Ratio - ratio of maximum to minimum over all processors
   Mess: number of messages sent
   Avg. len: average message length
   Reduct: number of global reductions
   Global: entire computation
   Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().
      %T - percent time in this phase         %f - percent flops in this phase
      %M - percent messages in this phase     %L - percent message lengths in this phase
      %R - percent reductions in this phase
   Total Mflop/s: 10e-6 * (sum of flops over all processors)/(max time over all processors)
------------------------------------------------------------------------------------------------------------------------
Event                Count      Time (sec)     Flops                             --- Global ---  --- Stage ---   Total
                   Max Ratio  Max     Ratio   Max  Ratio  Mess   Avg len Reduct  %T %f %M %L %R  %T %f %M %L %R Mflop/s
------------------------------------------------------------------------------------------------------------------------

--- Event Stage 0: Main Stage

VecCopy                3 1.0 3.8147e-06 2.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecSet                16 1.0 1.0252e-05 1.3 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecAXPY                6 1.0 1.5259e-05 1.3 4.25e+03 1.1 0.0e+00 0.0e+00 0.0e+00  0100  0  0  0   0100  0  0  0   532
VecAssemblyBegin      11 1.0 3.6809e-0340.6 0.00e+00 0.0 0.0e+00 0.0e+00 3.3e+01  1  0  0  0 41   1  0  0  0 41     0
VecAssemblyEnd        11 1.0 1.7881e-05 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
VecScatterBegin        2 1.0 3.2902e-05 1.0 0.00e+00 0.0 4.0e+00 3.8e+02 0.0e+00  0  0 33 66  0   0  0 33 66  0     0
VecScatterEnd          2 1.0 7.1526e-06 1.2 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
MatZeroEntries         6 1.0 3.1948e-05 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0     0
------------------------------------------------------------------------------------------------------------------------

Memory usage is given in bytes:

Object Type          Creations   Destructions     Memory  Descendants' Mem.
Reports information only for process 0.

--- Event Stage 0: Main Stage

              Vector    17             17        68296     0
      Vector Scatter     2              2         2072     0
           Index Set     4              4         3216     0
   IS L to G Mapping     2              2         1128     0
              Matrix     9              9       160656     0
              Viewer     1              0            0     0
========================================================================================================================
Average time to get PetscTime(): 9.53674e-08
Average time for MPI_Barrier(): 2.00272e-06
Average time for zero size MPI_Send(): 1.34706e-05
#PETSc Option Table entries:
-ksp_right_pc
-log_summary
-online_mode 1
-pc_type bjacobi
-sub_pc_factor_levels 4
-sub_pc_factor_zeropivot 0
-sub_pc_type ilu
#End of PETSc Option Table entries
Compiled without FORTRAN kernels
Compiled with full precision matrices (default)
sizeof(short) 2 sizeof(int) 4 sizeof(long) 8 sizeof(void*) 8 sizeof(PetscScalar) 8 sizeof(PetscInt) 4
Configure run at: Thu Nov  8 11:21:02 2012
Configure options: --with-debugging=false --COPTFLAGS=-O3 --CXXOPTFLAGS=-O3 --FOPTFLAGS=-O3 --with-clanguage=C++ --with-shared-libraries=1 --with-mpi-dir=/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1 --with-mumps=true --download-mumps=1 --with-metis=true --download-metis=1 --with-parmetis=true --download-parmetis=1 --with-superlu=true --download-superlu=1 --with-superludir=true --download-superlu_dist=1 --with-blacs=true --download-blacs=1 --with-scalapack=true --download-scalapack=1 --with-hypre=true --download-hypre=1 --with-blas-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_intel_lp64.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_sequential.so,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_core.so]" --with-lapack-lib="[/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64/libmkl_lapack95_lp64.a]"
-----------------------------------------
Libraries compiled on Thu Nov  8 11:21:02 2012 on daedalus.ices.utexas.edu 
Machine characteristics: Linux-2.6.32-279.1.1.el6.x86_64-x86_64-with-redhat-6.3-Carbon
Using PETSc directory: /opt/apps/ossw/libraries/petsc/petsc-3.3-p2
Using PETSc arch: intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt
-----------------------------------------

Using C compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx  -wd1572 -O3   -fPIC   ${COPTFLAGS} ${CFLAGS}
Using Fortran compiler: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90  -fPIC -O3   ${FOPTFLAGS} ${FFLAGS} 
-----------------------------------------

Using include paths: -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/include -I/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/include -I/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/include
-----------------------------------------

Using C linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpicxx
Using Fortran linker: /opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/bin/mpif90
Using libraries: -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lpetsc -lX11 -Wl,-rpath,/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -L/opt/apps/ossw/libraries/petsc/petsc-3.3-p2/intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt/lib -lcmumps -ldmumps -lsmumps -lzmumps -lmumps_common -lpord -lHYPRE -lpthread -lsuperlu_dist_3.0 -lparmetis -lmetis -lscalapack -lblacs -lsuperlu_4.3 -Wl,-rpath,/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -L/opt/apps/sysnet/intel/12.1/mkl/10.3.12.361/lib/intel64 -lmkl_lapack95_lp64 -lmkl_intel_lp64 -lmkl_sequential -lmkl_core -Wl,-rpath,/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -L/opt/apps/ossw/libraries/mpich2/mpich2-1.4.1p1/sl6/intel-12.1/lib -Wl,-rpath,/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -L/opt/apps/sysnet/intel/12.1/composer_xe_2011_sp1.7.256/compiler/lib/intel64 -Wl,-rpath,/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -L/usr/lib/gcc/x86_64-redhat-linux/4.4.6 -lmpichf90 -lifport -lifcore -lm -lm -lmpichcxx -ldl -lmpich -lopa -lmpl -lrt -lpthread -limf -lsvml -lipgo -ldecimal -lcilkrts -lstdc++ -lgcc_s -lirc -lirc_s -ldl 
-----------------------------------------


 ----------------------------------------------------------------------------------------------------------------------
| Processor id:   0                                                                                                    |
| Num Processors: 2                                                                                                    |
| Time:           Tue Feb  5 13:46:47 2013                                                                             |
| OS:             Linux                                                                                                |
| HostName:       hbar.ices.utexas.edu                                                                                 |
| OS Release:     2.6.32-279.1.1.el6.x86_64                                                                            |
| OS Version:     #1 SMP Tue Jul 10 11:24:23 CDT 2012                                                                  |
| Machine:        x86_64                                                                                               |
| Username:       benkirk                                                                                              |
| Configuration:  ./configure  '--enable-everything'                                                                   |
|  '--prefix=/workspace/libmesh/install'                                                                               |
|  'CXX=mpicxx'                                                                                                        |
|  'CC=mpicc'                                                                                                          |
|  'F77=mpif77'                                                                                                        |
|  'FC=mpif90'                                                                                                         |
|  'PETSC_DIR=/opt/apps/ossw/libraries/petsc/petsc-3.3-p2'                                                             |
|  'PETSC_ARCH=intel-12.1-mkl-intel-10.3.12.361-mpich2-1.4.1p1-cxx-opt'                                                |
|  'SLEPC_DIR=/opt/apps/ossw/libraries/slepc/slepc-3.3-p2-petsc-3.3-p2-cxx-opt'                                        |
|  'TRILINOS_DIR=/opt/apps/ossw/libraries/trilinos/trilinos-10.12.2/sl6/intel-12.1/mpich2-1.4.1p1/mkl-intel-10.3.12.361'|
|  'VTK_DIR=/opt/apps/ossw/libraries/vtk/vtk-5.10.0/sl6/intel-12.1'                                                    |
 ----------------------------------------------------------------------------------------------------------------------
 --------------------------------------------------------------------------------------------------------------
| libMesh Performance: Alive time=0.244959, Active time=0.221636                                               |
 --------------------------------------------------------------------------------------------------------------
| Event                            nCalls    Total Time  Avg Time    Total Time  Avg Time    % of Active Time  |
|                                            w/o Sub     w/o Sub     With Sub    With Sub    w/o S    With S   |
|--------------------------------------------------------------------------------------------------------------|
|                                                                                                              |
|                                                                                                              |
| DofMap                                                                                                       |
|   add_neighbors_to_send_list()   2         0.0144      0.007192    0.0177      0.008861    6.49     8.00     |
|   build_sparsity()               2         0.0092      0.004618    0.0270      0.013520    4.17     12.20    |
|   create_dof_constraints()       2         0.0062      0.003111    0.0231      0.011531    2.81     10.41    |
|   distribute_dofs()              2         0.0135      0.006758    0.0361      0.018061    6.10     16.30    |
|   dof_indices()                  2627      0.0704      0.000027    0.0704      0.000027    31.77    31.77    |
|   prepare_send_list()            2         0.0001      0.000030    0.0001      0.000030    0.03     0.03     |
|   reinit()                       2         0.0216      0.010791    0.0216      0.010791    9.74     9.74     |
|                                                                                                              |
| EquationSystems                                                                                              |
|   build_solution_vector()        2         0.0044      0.002221    0.0385      0.019248    2.00     17.37    |
|                                                                                                              |
| ExodusII_IO                                                                                                  |
|   write_nodal_data()             2         0.0067      0.003327    0.0067      0.003327    3.00     3.00     |
|                                                                                                              |
| Mesh                                                                                                         |
|   find_neighbors()               1         0.0075      0.007506    0.0076      0.007558    3.39     3.41     |
|   renumber_nodes_and_elem()      2         0.0003      0.000152    0.0003      0.000152    0.14     0.14     |
|                                                                                                              |
| MeshCommunication                                                                                            |
|   assign_global_indices()        1         0.0169      0.016885    0.0174      0.017399    7.62     7.85     |
|   compute_hilbert_indices()      2         0.0058      0.002890    0.0058      0.002890    2.61     2.61     |
|   find_global_indices()          2         0.0022      0.001077    0.0092      0.004616    0.97     4.17     |
|   parallel_sort()                2         0.0009      0.000444    0.0010      0.000487    0.40     0.44     |
|                                                                                                              |
| MeshOutput                                                                                                   |
|   write_equation_systems()       2         0.0001      0.000054    0.0454      0.022677    0.05     20.46    |
|                                                                                                              |
| MeshTools::Generation                                                                                        |
|   build_cube()                   1         0.0016      0.001598    0.0016      0.001598    0.72     0.72     |
|                                                                                                              |
| MetisPartitioner                                                                                             |
|   partition()                    1         0.0064      0.006417    0.0107      0.010711    2.90     4.83     |
|                                                                                                              |
| Parallel                                                                                                     |
|   allgather()                    17        0.0009      0.000055    0.0010      0.000056    0.42     0.43     |
|   barrier()                      1         0.0000      0.000016    0.0000      0.000016    0.01     0.01     |
|   broadcast()                    13        0.0002      0.000015    0.0002      0.000012    0.09     0.07     |
|   max(bool)                      3         0.0000      0.000002    0.0000      0.000002    0.00     0.00     |
|   max(scalar)                    186       0.0005      0.000002    0.0005      0.000002    0.21     0.21     |
|   max(vector)                    41        0.0003      0.000007    0.0006      0.000014    0.12     0.25     |
|   min(bool)                      213       0.0005      0.000002    0.0005      0.000002    0.21     0.21     |
|   min(scalar)                    175       0.0013      0.000008    0.0013      0.000008    0.61     0.61     |
|   min(vector)                    41        0.0003      0.000008    0.0006      0.000016    0.14     0.29     |
|   probe()                        30        0.0002      0.000008    0.0002      0.000008    0.11     0.11     |
|   receive()                      28        0.0002      0.000008    0.0005      0.000016    0.10     0.21     |
|   send()                         28        0.0001      0.000004    0.0001      0.000004    0.05     0.05     |
|   send_receive()                 30        0.0003      0.000010    0.0008      0.000028    0.14     0.38     |
|   sum()                          44        0.0003      0.000007    0.0005      0.000011    0.13     0.22     |
|                                                                                                              |
| Parallel::Request                                                                                            |
|   wait()                         28        0.0001      0.000002    0.0001      0.000002    0.03     0.03     |
|                                                                                                              |
| Partitioner                                                                                                  |
|   set_node_processor_ids()       1         0.0008      0.000801    0.0009      0.000875    0.36     0.39     |
|   set_parent_processor_ids()     1         0.0005      0.000498    0.0005      0.000498    0.22     0.22     |
|                                                                                                              |
| RBConstruction                                                                                               |
|   clear()                        1         0.0002      0.000235    0.0002      0.000235    0.11     0.11     |
|   load_basis_function()          1         0.0000      0.000044    0.0000      0.000044    0.02     0.02     |
|   load_rb_solution()             1         0.0001      0.000139    0.0001      0.000139    0.06     0.06     |
|                                                                                                              |
| RBEvaluation                                                                                                 |
|   clear()                        1         0.0001      0.000057    0.0001      0.000057    0.03     0.03     |
|   compute_residual_dual_norm()   1         0.0001      0.000146    0.0001      0.000146    0.07     0.07     |
|   rb_solve()                     1         0.0064      0.006448    0.0069      0.006859    2.91     3.09     |
|   read_in_basis_functions()      1         0.0000      0.000038    0.0363      0.036254    0.02     16.36    |
|   read_in_vectors()              1         0.0184      0.018400    0.0362      0.036216    8.30     16.34    |
|   read_offline_data_from_files() 1         0.0007      0.000661    0.0010      0.001001    0.30     0.45     |
|   resize_data_structures()       1         0.0003      0.000340    0.0003      0.000340    0.15     0.15     |
|                                                                                                              |
| RBSCMEvaluation                                                                                              |
|   get_SCM_LB()                   1         0.0003      0.000263    0.0003      0.000263    0.12     0.12     |
|   read_offline_data_from_files() 1         0.0002      0.000174    0.0002      0.000174    0.08     0.08     |
 --------------------------------------------------------------------------------------------------------------
| Totals:                          3549      0.2216                                          100.00            |
 --------------------------------------------------------------------------------------------------------------

 
***************************************************************
* Done Running Example reduced_basis_ex2:
*  mpirun -np 2 example-devel -online_mode 1 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************

Site Created By: libMesh Developers
Last modified: February 05 2013 19:46:47 UTC

Hosted By:
SourceForge.net Logo