Example 10 - Solving a Transient System with Adaptive Mesh Refinement



This example shows how a simple, linear transient system can be solved in parallel. The system is simple scalar convection-diffusion with a specified external velocity. The initial condition is given, and the solution is advanced in time with a standard Crank-Nicolson time-stepping strategy.

Also, we use this example to demonstrate writing out and reading in of solutions. We do 25 time steps, then save the solution and do another 25 time steps starting from the saved solution.

C++ include files that we need
        #include <iostream>
        #include <algorithm>
        #include <cmath>
        
Basic include file needed for the mesh functionality.
        #include "libmesh.h"
        #include "mesh.h"
        #include "mesh_refinement.h"
        #include "gmv_io.h"
        #include "equation_systems.h"
        #include "fe.h"
        #include "quadrature_gauss.h"
        #include "dof_map.h"
        #include "sparse_matrix.h"
        #include "numeric_vector.h"
        #include "dense_matrix.h"
        #include "dense_vector.h"
        
        #include "getpot.h"
        
Some (older) compilers do not offer full stream functionality, \p OStringStream works around this. Check example 9 for details.
        #include "o_string_stream.h"
        
This example will solve a linear transient system, so we need to include the \p TransientLinearImplicitSystem definition.
        #include "transient_system.h"
        #include "linear_implicit_system.h"
        #include "vector_value.h"
        
To refine the mesh we need an \p ErrorEstimator object to figure out which elements to refine.
        #include "error_vector.h"
        #include "kelly_error_estimator.h"
        
The definition of a geometric element
        #include "elem.h"
        
Bring in everything from the libMesh namespace
        using namespace libMesh;
        
Function prototype. This function will assemble the system matrix and right-hand-side at each time step. Note that since the system is linear we technically do not need to assmeble the matrix at each time step, but we will anyway. In subsequent examples we will employ adaptive mesh refinement, and with a changing mesh it will be necessary to rebuild the system matrix.
        void assemble_cd (EquationSystems& es,
                          const std::string& system_name);
        
Function prototype. This function will initialize the system. Initialization functions are optional for systems. They allow you to specify the initial values of the solution. If an initialization function is not provided then the default (0) solution is provided.
        void init_cd (EquationSystems& es,
                      const std::string& system_name);
        
Exact solution function prototype. This gives the exact solution as a function of space and time. In this case the initial condition will be taken as the exact solution at time 0, as will the Dirichlet boundary conditions at time t.
        Real exact_solution (const Real x,
                             const Real y,
                             const Real t);
        
        Number exact_value (const Point& p,
                            const Parameters& parameters,
                            const std::string&,
                            const std::string&)
        {
          return exact_solution(p(0), p(1), parameters.get<Real> ("time"));
        }
        
        
        
Begin the main program. Note that the first statement in the program throws an error if you are in complex number mode, since this example is only intended to work with real numbers.
        int main (int argc, char** argv)
        {
Initialize libMesh.
          LibMeshInit init (argc, argv);
        
        #ifndef LIBMESH_ENABLE_AMR
          libmesh_example_assert(false, "--enable-amr");
        #else
        
Our Trilinos interface does not yet support adaptive transient problems
          libmesh_example_assert(libMesh::default_solver_package() != TRILINOS_SOLVERS, "--enable-petsc");
        
Brief message to the user regarding the program name and command line arguments.

Use commandline parameter to specify if we are to read in an initial solution or generate it ourself
          std::cout << "Usage:\n"
            <<"\t " << argv[0] << " -init_timestep 0\n"
            << "OR\n"
            <<"\t " << argv[0] << " -read_solution -init_timestep 26\n"
            << std::endl;
        
          std::cout << "Running: " << argv[0];
        
          for (int i=1; i<argc; i++)
            std::cout << " " << argv[i];
        
          std::cout << std::endl << std::endl;
        
Create a GetPot object to parse the command line
          GetPot command_line (argc, argv);
        
        
This boolean value is obtained from the command line, it is true if the flag "-read_solution" is present, false otherwise. It indicates whether we are going to read in the mesh and solution files "saved_mesh.xda" and "saved_solution.xda" or whether we are going to start from scratch by just reading "mesh.xda"
          const bool read_solution   = command_line.search("-read_solution");
        
This value is also obtained from the commandline and it specifies the initial value for the t_step looping variable. We must distinguish between the two cases here, whether we read in the solution or we started from scratch, so that we do not overwrite the gmv output files.
          unsigned int init_timestep = 0;
          
Search the command line for the "init_timestep" flag and if it is present, set init_timestep accordingly.
          if(command_line.search("-init_timestep"))
            init_timestep = command_line.next(0);
          else
            {
              if (libMesh::processor_id() == 0)
                std::cerr << "ERROR: Initial timestep not specified\n" << std::endl;
        
This handy function will print the file name, line number, and then abort. Currrently the library does not use C++ exception handling.
              libmesh_error();
            }
        
This value is also obtained from the command line, and specifies the number of time steps to take.
          unsigned int n_timesteps = 0;
        
Again do a search on the command line for the argument
          if(command_line.search("-n_timesteps"))
            n_timesteps = command_line.next(0);
          else
            {
              std::cout << "ERROR: Number of timesteps not specified\n" << std::endl;
              libmesh_error();
            }
        
        
Skip this 2D example if libMesh was compiled as 1D-only.
          libmesh_example_assert(2 <= LIBMESH_DIM, "2D support");
        
Create a new mesh.
          Mesh mesh;
        
Create an equation systems object.
          EquationSystems equation_systems (mesh);
          MeshRefinement mesh_refinement (mesh);
        
First we process the case where we do not read in the solution
          if(!read_solution)
            {
Read the mesh from file.
              mesh.read ("mesh.xda");
        
Again do a search on the command line for an argument
              unsigned int n_refinements = 5;
              if(command_line.search("-n_refinements"))
                n_refinements = command_line.next(0);
        
Uniformly refine the mesh 5 times
              if(!read_solution)
                mesh_refinement.uniformly_refine (n_refinements);
        
Print information about the mesh to the screen.
              mesh.print_info();
              
              
Declare the system and its variables. Begin by creating a transient system named "Convection-Diffusion".
              TransientLinearImplicitSystem & system = 
                equation_systems.add_system<TransientLinearImplicitSystem> 
                ("Convection-Diffusion");
        
Adds the variable "u" to "Convection-Diffusion". "u" will be approximated using first-order approximation.
              system.add_variable ("u", FIRST);
        
Give the system a pointer to the matrix assembly and initialization functions.
              system.attach_assemble_function (assemble_cd);
              system.attach_init_function (init_cd);
        
Initialize the data structures for the equation system.
              equation_systems.init ();
            }
Otherwise we read in the solution and mesh
          else 
            {
Read in the mesh stored in "saved_mesh.xda"
              mesh.read("saved_mesh.xda");
        
Print information about the mesh to the screen.
              mesh.print_info();
        
Read in the solution stored in "saved_solution.xda"
              equation_systems.read("saved_solution.xda", libMeshEnums::READ);
        
Get a reference to the system so that we can call update() on it
              TransientLinearImplicitSystem & system = 
                equation_systems.get_system<TransientLinearImplicitSystem> 
                ("Convection-Diffusion");
        
We need to call update to put system in a consistent state with the solution that was read in
              system.update();
        
Attach the same matrix assembly function as above. Note, we do not have to attach an init() function since we are initializing the system by reading in "saved_solution.xda"
              system.attach_assemble_function (assemble_cd);
        
Print out the H1 norm of the saved solution, for verification purposes:
              Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
        
              std::cout << "Initial H1 norm = " << H1norm << std::endl << std::endl;
            }
        
Prints information about the system to the screen.
          equation_systems.print_info();
            
          equation_systems.parameters.set<unsigned int>
            ("linear solver maximum iterations") = 250;
          equation_systems.parameters.set<Real>
            ("linear solver tolerance") = TOLERANCE;
            
          if(!read_solution)
Write out the initial condition
            GMVIO(mesh).write_equation_systems ("out.gmv.000",
                                                equation_systems);
          else
Write out the solution that was read in
            GMVIO(mesh).write_equation_systems ("solution_read_in.gmv",
                                                equation_systems);
        
The Convection-Diffusion system requires that we specify the flow velocity. We will specify it as a RealVectorValue data type and then use the Parameters object to pass it to the assemble function.
          equation_systems.parameters.set<RealVectorValue>("velocity") = 
            RealVectorValue (0.8, 0.8);
        
The Convection-Diffusion system also requires a specified diffusivity. We use an isotropic (hence Real) value.
          equation_systems.parameters.set<Real>("diffusivity") = 0.01;
            
Solve the system "Convection-Diffusion". This will be done by looping over the specified time interval and calling the \p solve() member at each time step. This will assemble the system and call the linear solver.
          const Real dt = 0.025;
          Real time     = init_timestep*dt;
          
We do 25 timesteps both before and after writing out the intermediate solution
          for(unsigned int t_step=init_timestep; 
                           t_step<(init_timestep+n_timesteps); 
                           t_step++)
            {
Increment the time counter, set the time and the time step size as parameters in the EquationSystem.
              time += dt;
        
              equation_systems.parameters.set<Real> ("time") = time;
              equation_systems.parameters.set<Real> ("dt")   = dt;
        
A pretty update message
              std::cout << " Solving time step ";
              
As already seen in example 9, use a work-around for missing stream functionality (of older compilers). Add a set of scope braces to enforce data locality.
              {
                OStringStream out;
        
                OSSInt(out,2,t_step);
                out << ", time=";
                OSSRealzeroleft(out,6,3,time);
                out <<  "...";
                std::cout << out.str() << std::endl;
              }
              
At this point we need to update the old solution vector. The old solution vector will be the current solution vector from the previous time step. We will do this by extracting the system from the \p EquationSystems object and using vector assignment. Since only \p TransientLinearImplicitSystems (and systems derived from them) contain old solutions we need to specify the system type when we ask for it.
              TransientLinearImplicitSystem &  system =
                equation_systems.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
        
              *system.old_local_solution = *system.current_local_solution;
              
The number of refinement steps per time step.
              const unsigned int max_r_steps = 2;
              
A refinement loop.
              for (unsigned int r_step=0; r_step<max_r_steps; r_step++)
                {
Assemble & solve the linear system
                  system.solve();
        
Print out the H1 norm, for verification purposes:
                  Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
        
                  std::cout << "H1 norm = " << H1norm << std::endl;
                  
Possibly refine the mesh
                  if (r_step+1 != max_r_steps)
                    {
                      std::cout << "  Refining the mesh..." << std::endl;
        
The \p ErrorVector is a particular \p StatisticsVector for computing error information on a finite element mesh.
                      ErrorVector error;
        
The \p ErrorEstimator class interrogates a finite element solution and assigns to each element a positive error value. This value is used for deciding which elements to refine and which to coarsen. ErrorEstimator* error_estimator = new KellyErrorEstimator;
                      KellyErrorEstimator error_estimator;
                      
Compute the error for each active element using the provided \p flux_jump indicator. Note in general you will need to provide an error estimator specifically designed for your application.
                      error_estimator.estimate_error (system,
                                                      error);
                      
This takes the error in \p error and decides which elements will be coarsened or refined. Any element within 20% of the maximum error on any element will be refined, and any element within 7% of the minimum error on any element might be coarsened. Note that the elements flagged for refinement will be refined, but those flagged for coarsening _might_ be coarsened.
                      mesh_refinement.refine_fraction() = 0.80;
                      mesh_refinement.coarsen_fraction() = 0.07;
                      mesh_refinement.max_h_level() = 5;
                      mesh_refinement.flag_elements_by_error_fraction (error);
                      
This call actually refines and coarsens the flagged elements.
                      mesh_refinement.refine_and_coarsen_elements();
                      
This call reinitializes the \p EquationSystems object for the newly refined mesh. One of the steps in the reinitialization is projecting the \p solution, \p old_solution, etc... vectors from the old mesh to the current one.
                      equation_systems.reinit ();
                    }            
                }
                
Again do a search on the command line for an argument
              unsigned int output_freq = 10;
              if(command_line.search("-output_freq"))
                output_freq = command_line.next(0);
        
Output every 10 timesteps to file.
              if ( (t_step+1)%output_freq == 0)
                {
                  OStringStream file_name;
        
                  file_name << "out.gmv.";
                  OSSRealzeroright(file_name,3,0,t_step+1);
        
                  GMVIO(mesh).write_equation_systems (file_name.str(),
                                                      equation_systems);
                }
            }
        
          if(!read_solution)
            {
Print out the H1 norm of the saved solution, for verification purposes:
              TransientLinearImplicitSystem& system =
        	equation_systems.get_system<TransientLinearImplicitSystem>
                  ("Convection-Diffusion");
              Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
        
              std::cout << "Final H1 norm = " << H1norm << std::endl << std::endl;
        
              mesh.write("saved_mesh.xda");
              equation_systems.write("saved_solution.xda", libMeshEnums::WRITE);
              GMVIO(mesh).write_equation_systems ("saved_solution.gmv",
                                                  equation_systems);
            }
        #endif // #ifndef LIBMESH_ENABLE_AMR
          
          return 0;
        }
        
Here we define the initialization routine for the Convection-Diffusion system. This routine is responsible for applying the initial conditions to the system.
        void init_cd (EquationSystems& es,
                      const std::string& system_name)
        {
It is a good idea to make sure we are initializing the proper system.
          libmesh_assert (system_name == "Convection-Diffusion");
        
Get a reference to the Convection-Diffusion system object.
          TransientLinearImplicitSystem & system =
            es.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
        
Project initial conditions at time 0
          es.parameters.set<Real> ("time") = 0;
          
          system.project_solution(exact_value, NULL, es.parameters);
        }
        
        
        
This function defines the assembly routine which will be called at each time step. It is responsible for computing the proper matrix entries for the element stiffness matrices and right-hand sides.
        void assemble_cd (EquationSystems& es,
                          const std::string& system_name)
        {
        #ifdef LIBMESH_ENABLE_AMR
It is a good idea to make sure we are assembling the proper system.
          libmesh_assert (system_name == "Convection-Diffusion");
          
Get a constant reference to the mesh object.
          const MeshBase& mesh = es.get_mesh();
          
The dimension that we are running
          const unsigned int dim = mesh.mesh_dimension();
          
Get a reference to the Convection-Diffusion system object.
          TransientLinearImplicitSystem & system =
            es.get_system<TransientLinearImplicitSystem> ("Convection-Diffusion");
          
Get the Finite Element type for the first (and only) variable in the system.
          FEType fe_type = system.variable_type(0);
          
Build a Finite Element object of the specified type. Since the \p FEBase::build() member dynamically creates memory we will store the object as an \p AutoPtr. This can be thought of as a pointer that will clean up after itself.
          AutoPtr<FEBase> fe      (FEBase::build(dim, fe_type));
          AutoPtr<FEBase> fe_face (FEBase::build(dim, fe_type));
          
A Gauss quadrature rule for numerical integration. Let the \p FEType object decide what order rule is appropriate.
          QGauss qrule (dim,   fe_type.default_quadrature_order());
          QGauss qface (dim-1, fe_type.default_quadrature_order());
        
Tell the finite element object to use our quadrature rule.
          fe->attach_quadrature_rule      (&qrule);
          fe_face->attach_quadrature_rule (&qface);
        
Here we define some references to cell-specific data that will be used to assemble the linear system. We will start with the element Jacobian * quadrature weight at each integration point.
          const std::vector<Real>& JxW      = fe->get_JxW();
          const std::vector<Real>& JxW_face = fe_face->get_JxW();
          
The element shape functions evaluated at the quadrature points.
          const std::vector<std::vector<Real> >& phi = fe->get_phi();
          const std::vector<std::vector<Real> >& psi = fe_face->get_phi();
        
The element shape function gradients evaluated at the quadrature points.
          const std::vector<std::vector<RealGradient> >& dphi = fe->get_dphi();
        
The XY locations of the quadrature points used for face integration
          const std::vector<Point>& qface_points = fe_face->get_xyz();
            
A reference to the \p DofMap object for this system. The \p DofMap object handles the index translation from node and element numbers to degree of freedom numbers. We will talk more about the \p DofMap in future examples.
          const DofMap& dof_map = system.get_dof_map();
          
Define data structures to contain the element matrix and right-hand-side vector contribution. Following basic finite element terminology we will denote these "Ke" and "Fe".
          DenseMatrix<Number> Ke;
          DenseVector<Number> Fe;
          
This vector will hold the degree of freedom indices for the element. These define where in the global system the element degrees of freedom get mapped.
          std::vector<unsigned int> dof_indices;
        
Here we extract the velocity & parameters that we put in the EquationSystems object.
          const RealVectorValue velocity =
            es.parameters.get<RealVectorValue> ("velocity");
        
          const Real diffusivity =
            es.parameters.get<Real> ("diffusivity");
        
          const Real dt = es.parameters.get<Real>   ("dt");
          const Real time = es.parameters.get<Real> ("time");
        
Now we will loop over all the elements in the mesh that live on the local processor. We will compute the element matrix and right-hand-side contribution. Since the mesh will be refined we want to only consider the ACTIVE elements, hence we use a variant of the \p active_elem_iterator.
          MeshBase::const_element_iterator       el     = mesh.active_local_elements_begin();
          const MeshBase::const_element_iterator end_el = mesh.active_local_elements_end(); 
          
          for ( ; el != end_el; ++el)
            {    
Store a pointer to the element we are currently working on. This allows for nicer syntax later.
              const Elem* elem = *el;
              
Get the degree of freedom indices for the current element. These define where in the global matrix and right-hand-side this element will contribute to.
              dof_map.dof_indices (elem, dof_indices);
        
Compute the element-specific data for the current element. This involves computing the location of the quadrature points (q_point) and the shape functions (phi, dphi) for the current element.
              fe->reinit (elem);
              
Zero the element matrix and right-hand side before summing them. We use the resize member here because the number of degrees of freedom might have changed from the last element. Note that this will be the case if the element type is different (i.e. the last element was a triangle, now we are on a quadrilateral).
              Ke.resize (dof_indices.size(),
                         dof_indices.size());
        
              Fe.resize (dof_indices.size());
              
Now we will build the element matrix and right-hand-side. Constructing the RHS requires the solution and its gradient from the previous timestep. This myst be calculated at each quadrature point by summing the solution degree-of-freedom values by the appropriate weight functions.
              for (unsigned int qp=0; qp<qrule.n_points(); qp++)
                {
Values to hold the old solution & its gradient.
                  Number   u_old = 0.;
                  Gradient grad_u_old;
                  
Compute the old solution & its gradient.
                  for (unsigned int l=0; l<phi.size(); l++)
                    {
                      u_old      += phi[l][qp]*system.old_solution  (dof_indices[l]);
                      
This will work, grad_u_old += dphi[l][qp]*system.old_solution (dof_indices[l]); but we can do it without creating a temporary like this:
                      grad_u_old.add_scaled (dphi[l][qp],system.old_solution (dof_indices[l]));
                    }
                  
Now compute the element matrix and RHS contributions.
                  for (unsigned int i=0; i<phi.size(); i++)
                    {
The RHS contribution
                      Fe(i) += JxW[qp]*(
Mass matrix term
                                        u_old*phi[i][qp] + 
                                        -.5*dt*(
Convection term (grad_u_old may be complex, so the order here is important!)
                                                (grad_u_old*velocity)*phi[i][qp] +
                                                
Diffusion term
                                                diffusivity*(grad_u_old*dphi[i][qp]))     
                                        );
                      
                      for (unsigned int j=0; j<phi.size(); j++)
                        {
The matrix contribution
                          Ke(i,j) += JxW[qp]*(
Mass-matrix
                                              phi[i][qp]*phi[j][qp] + 
                                              .5*dt*(
Convection term
                                                     (velocity*dphi[j][qp])*phi[i][qp] +
Diffusion term
                                                     diffusivity*(dphi[i][qp]*dphi[j][qp]))      
                                              );
                        }
                    } 
                } 
        
At this point the interior element integration has been completed. However, we have not yet addressed boundary conditions. For this example we will only consider simple Dirichlet boundary conditions imposed via the penalty method.

The following loops over the sides of the element. If the element has no neighbor on a side then that side MUST live on a boundary of the domain.
              {
The penalty value.
                const Real penalty = 1.e10;
        
The following loops over the sides of the element. If the element has no neighbor on a side then that side MUST live on a boundary of the domain.
                for (unsigned int s=0; s<elem->n_sides(); s++)
                  if (elem->neighbor(s) == NULL)
                    {
                      fe_face->reinit(elem,s);
                      
                      for (unsigned int qp=0; qp<qface.n_points(); qp++)
                        {
                          const Number value = exact_solution (qface_points[qp](0),
                                                               qface_points[qp](1),
                                                               time);
                                                               
RHS contribution
                          for (unsigned int i=0; i<psi.size(); i++)
                            Fe(i) += penalty*JxW_face[qp]*value*psi[i][qp];
        
Matrix contribution
                          for (unsigned int i=0; i<psi.size(); i++)
                            for (unsigned int j=0; j<psi.size(); j++)
                              Ke(i,j) += penalty*JxW_face[qp]*psi[i][qp]*psi[j][qp];
                        }
                    } 
              } 
        
              
We have now built the element matrix and RHS vector in terms of the element degrees of freedom. However, it is possible that some of the element DOFs are constrained to enforce solution continuity, i.e. they are not really "free". We need to constrain those DOFs in terms of non-constrained DOFs to ensure a continuous solution. The \p DofMap::constrain_element_matrix_and_vector() method does just that.
              dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
              
The element matrix and right-hand-side are now built for this element. Add them to the global matrix and right-hand-side vector. The \p SparseMatrix::add_matrix() and \p NumericVector::add_vector() members do this for us.
              system.matrix->add_matrix (Ke, dof_indices);
              system.rhs->add_vector    (Fe, dof_indices);
              
            }
Finished computing the sytem matrix and right-hand side.
        #endif // #ifdef LIBMESH_ENABLE_AMR
        }



The program without comments:

 
   
  #include <iostream>
  #include <algorithm>
  #include <cmath>
  
  #include "libmesh.h"
  #include "mesh.h"
  #include "mesh_refinement.h"
  #include "gmv_io.h"
  #include "equation_systems.h"
  #include "fe.h"
  #include "quadrature_gauss.h"
  #include "dof_map.h"
  #include "sparse_matrix.h"
  #include "numeric_vector.h"
  #include "dense_matrix.h"
  #include "dense_vector.h"
  
  #include "getpot.h"
  
  #include "o_string_stream.h"
  
  #include "transient_system.h"
  #include "linear_implicit_system.h"
  #include "vector_value.h"
  
  #include "error_vector.h"
  #include "kelly_error_estimator.h"
  
  #include "elem.h"
  
  using namespace libMesh;
  
  void assemble_cd (EquationSystems& es,
                    const std::string& system_name);
  
  void init_cd (EquationSystems& es,
                const std::string& system_name);
  
  Real exact_solution (const Real x,
                       const Real y,
                       const Real t);
  
  Number exact_value (const Point& p,
                      const Parameters& parameters,
                      const std::string&,
                      const std::string&)
  {
    return exact_solution(p(0), p(1), parameters.get<Real> ("time"));
  }
  
  
  
  int main (int argc, char** argv)
  {
    LibMeshInit init (argc, argv);
  
  #ifndef LIBMESH_ENABLE_AMR
    libmesh_example_assert(false, "--enable-amr");
  #else
  
    libmesh_example_assert(libMesh::default_solver_package() != TRILINOS_SOLVERS, "--enable-petsc");
  
  
    std::cout << "Usage:\n"
      <<"\t " << argv[0] << " -init_timestep 0\n"
      << "OR\n"
      <<"\t " << argv[0] << " -read_solution -init_timestep 26\n"
      << std::endl;
  
    std::cout << "Running: " << argv[0];
  
    for (int i=1; i<argc; i++)
      std::cout << " " << argv[i];
  
    std::cout << std::endl << std::endl;
  
    GetPot command_line (argc, argv);
  
  
    const bool read_solution   = command_line.search("-read_solution");
  
    unsigned int init_timestep = 0;
    
    if(command_line.search("-init_timestep"))
      init_timestep = command_line.next(0);
    else
      {
        if (libMesh::processor_id() == 0)
          std::cerr << "ERROR: Initial timestep not specified\n" << std::endl;
  
        libmesh_error();
      }
  
    unsigned int n_timesteps = 0;
  
    if(command_line.search("-n_timesteps"))
      n_timesteps = command_line.next(0);
    else
      {
        std::cout << "ERROR: Number of timesteps not specified\n" << std::endl;
        libmesh_error();
      }
  
  
    libmesh_example_assert(2 <= LIBMESH_DIM, "2D support");
  
    Mesh mesh;
  
    EquationSystems equation_systems (mesh);
    MeshRefinement mesh_refinement (mesh);
  
    if(!read_solution)
      {
        mesh.read ("mesh.xda");
  
        unsigned int n_refinements = 5;
        if(command_line.search("-n_refinements"))
          n_refinements = command_line.next(0);
  
        if(!read_solution)
          mesh_refinement.uniformly_refine (n_refinements);
  
        mesh.print_info();
        
        
        TransientLinearImplicitSystem & system = 
          equation_systems.add_system<TransientLinearImplicitSystem> 
          ("Convection-Diffusion");
  
        system.add_variable ("u", FIRST);
  
        system.attach_assemble_function (assemble_cd);
        system.attach_init_function (init_cd);
  
        equation_systems.init ();
      }
    else 
      {
        mesh.read("saved_mesh.xda");
  
        mesh.print_info();
  
        equation_systems.read("saved_solution.xda", libMeshEnums::READ);
  
        TransientLinearImplicitSystem & system = 
          equation_systems.get_system<TransientLinearImplicitSystem> 
          ("Convection-Diffusion");
  
        system.update();
  
        system.attach_assemble_function (assemble_cd);
  
        Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
  
        std::cout << "Initial H1 norm = " << H1norm << std::endl << std::endl;
      }
  
    equation_systems.print_info();
      
    equation_systems.parameters.set<unsigned int>
      ("linear solver maximum iterations") = 250;
    equation_systems.parameters.set<Real>
      ("linear solver tolerance") = TOLERANCE;
      
    if(!read_solution)
      GMVIO(mesh).write_equation_systems ("out.gmv.000",
                                          equation_systems);
    else
      GMVIO(mesh).write_equation_systems ("solution_read_in.gmv",
                                          equation_systems);
  
    equation_systems.parameters.set<RealVectorValue>("velocity") = 
      RealVectorValue (0.8, 0.8);
  
    equation_systems.parameters.set<Real>("diffusivity") = 0.01;
      
    const Real dt = 0.025;
    Real time     = init_timestep*dt;
    
    for(unsigned int t_step=init_timestep; 
                     t_step<(init_timestep+n_timesteps); 
                     t_step++)
      {
        time += dt;
  
        equation_systems.parameters.set<Real> ("time") = time;
        equation_systems.parameters.set<Real> ("dt")   = dt;
  
        std::cout << " Solving time step ";
        
        {
          OStringStream out;
  
          OSSInt(out,2,t_step);
          out << ", time=";
          OSSRealzeroleft(out,6,3,time);
          out <<  "...";
          std::cout << out.str() << std::endl;
        }
        
        TransientLinearImplicitSystem &  system =
          equation_systems.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
  
        *system.old_local_solution = *system.current_local_solution;
        
        const unsigned int max_r_steps = 2;
        
        for (unsigned int r_step=0; r_step<max_r_steps; r_step++)
          {
            system.solve();
  
            Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
  
            std::cout << "H1 norm = " << H1norm << std::endl;
            
            if (r_step+1 != max_r_steps)
              {
                std::cout << "  Refining the mesh..." << std::endl;
  
                ErrorVector error;
  
                KellyErrorEstimator error_estimator;
                
                error_estimator.estimate_error (system,
                                                error);
                
                mesh_refinement.refine_fraction() = 0.80;
                mesh_refinement.coarsen_fraction() = 0.07;
                mesh_refinement.max_h_level() = 5;
                mesh_refinement.flag_elements_by_error_fraction (error);
                
                mesh_refinement.refine_and_coarsen_elements();
                
                equation_systems.reinit ();
              }            
          }
          
        unsigned int output_freq = 10;
        if(command_line.search("-output_freq"))
          output_freq = command_line.next(0);
  
        if ( (t_step+1)%output_freq == 0)
          {
            OStringStream file_name;
  
            file_name << "out.gmv.";
            OSSRealzeroright(file_name,3,0,t_step+1);
  
            GMVIO(mesh).write_equation_systems (file_name.str(),
                                                equation_systems);
          }
      }
  
    if(!read_solution)
      {
        TransientLinearImplicitSystem& system =
  	equation_systems.get_system<TransientLinearImplicitSystem>
            ("Convection-Diffusion");
        Real H1norm = system.calculate_norm(*system.solution, SystemNorm(H1));
  
        std::cout << "Final H1 norm = " << H1norm << std::endl << std::endl;
  
        mesh.write("saved_mesh.xda");
        equation_systems.write("saved_solution.xda", libMeshEnums::WRITE);
        GMVIO(mesh).write_equation_systems ("saved_solution.gmv",
                                            equation_systems);
      }
  #endif // #ifndef LIBMESH_ENABLE_AMR
    
    return 0;
  }
  
  void init_cd (EquationSystems& es,
                const std::string& system_name)
  {
    libmesh_assert (system_name == "Convection-Diffusion");
  
    TransientLinearImplicitSystem & system =
      es.get_system<TransientLinearImplicitSystem>("Convection-Diffusion");
  
    es.parameters.set<Real> ("time") = 0;
    
    system.project_solution(exact_value, NULL, es.parameters);
  }
  
  
  
  void assemble_cd (EquationSystems& es,
                    const std::string& system_name)
  {
  #ifdef LIBMESH_ENABLE_AMR
    libmesh_assert (system_name == "Convection-Diffusion");
    
    const MeshBase& mesh = es.get_mesh();
    
    const unsigned int dim = mesh.mesh_dimension();
    
    TransientLinearImplicitSystem & system =
      es.get_system<TransientLinearImplicitSystem> ("Convection-Diffusion");
    
    FEType fe_type = system.variable_type(0);
    
    AutoPtr<FEBase> fe      (FEBase::build(dim, fe_type));
    AutoPtr<FEBase> fe_face (FEBase::build(dim, fe_type));
    
    QGauss qrule (dim,   fe_type.default_quadrature_order());
    QGauss qface (dim-1, fe_type.default_quadrature_order());
  
    fe->attach_quadrature_rule      (&qrule);
    fe_face->attach_quadrature_rule (&qface);
  
    const std::vector<Real>& JxW      = fe->get_JxW();
    const std::vector<Real>& JxW_face = fe_face->get_JxW();
    
    const std::vector<std::vector<Real> >& phi = fe->get_phi();
    const std::vector<std::vector<Real> >& psi = fe_face->get_phi();
  
    const std::vector<std::vector<RealGradient> >& dphi = fe->get_dphi();
  
    const std::vector<Point>& qface_points = fe_face->get_xyz();
      
    const DofMap& dof_map = system.get_dof_map();
    
    DenseMatrix<Number> Ke;
    DenseVector<Number> Fe;
    
    std::vector<unsigned int> dof_indices;
  
    const RealVectorValue velocity =
      es.parameters.get<RealVectorValue> ("velocity");
  
    const Real diffusivity =
      es.parameters.get<Real> ("diffusivity");
  
    const Real dt = es.parameters.get<Real>   ("dt");
    const Real time = es.parameters.get<Real> ("time");
  
    MeshBase::const_element_iterator       el     = mesh.active_local_elements_begin();
    const MeshBase::const_element_iterator end_el = mesh.active_local_elements_end(); 
    
    for ( ; el != end_el; ++el)
      {    
        const Elem* elem = *el;
        
        dof_map.dof_indices (elem, dof_indices);
  
        fe->reinit (elem);
        
        Ke.resize (dof_indices.size(),
                   dof_indices.size());
  
        Fe.resize (dof_indices.size());
        
        for (unsigned int qp=0; qp<qrule.n_points(); qp++)
          {
            Number   u_old = 0.;
            Gradient grad_u_old;
            
            for (unsigned int l=0; l<phi.size(); l++)
              {
                u_old      += phi[l][qp]*system.old_solution  (dof_indices[l]);
                
                grad_u_old.add_scaled (dphi[l][qp],system.old_solution (dof_indices[l]));
              }
            
            for (unsigned int i=0; i<phi.size(); i++)
              {
                Fe(i) += JxW[qp]*(
                                  u_old*phi[i][qp] + 
                                  -.5*dt*(
                                          (grad_u_old*velocity)*phi[i][qp] +
                                          
                                          diffusivity*(grad_u_old*dphi[i][qp]))     
                                  );
                
                for (unsigned int j=0; j<phi.size(); j++)
                  {
                    Ke(i,j) += JxW[qp]*(
                                        phi[i][qp]*phi[j][qp] + 
                                        .5*dt*(
                                               (velocity*dphi[j][qp])*phi[i][qp] +
                                               diffusivity*(dphi[i][qp]*dphi[j][qp]))      
                                        );
                  }
              } 
          } 
  
        {
          const Real penalty = 1.e10;
  
          for (unsigned int s=0; s<elem->n_sides(); s++)
            if (elem->neighbor(s) == NULL)
              {
                fe_face->reinit(elem,s);
                
                for (unsigned int qp=0; qp<qface.n_points(); qp++)
                  {
                    const Number value = exact_solution (qface_points[qp](0),
                                                         qface_points[qp](1),
                                                         time);
                                                         
                    for (unsigned int i=0; i<psi.size(); i++)
                      Fe(i) += penalty*JxW_face[qp]*value*psi[i][qp];
  
                    for (unsigned int i=0; i<psi.size(); i++)
                      for (unsigned int j=0; j<psi.size(); j++)
                        Ke(i,j) += penalty*JxW_face[qp]*psi[i][qp]*psi[j][qp];
                  }
              } 
        } 
  
        
        dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
        
        system.matrix->add_matrix (Ke, dof_indices);
        system.rhs->add_vector    (Fe, dof_indices);
        
      }
  #endif // #ifdef LIBMESH_ENABLE_AMR
  }



The console output of the program:

Compiling C++ (in optimized mode) ex10.C...
Compiling C++ (in optimized mode) ../ex9/exact_solution.C...
Linking ex10-opt...
***************************************************************
* Running Example  mpirun -np 2 ./ex10-opt [-read_solution] -n_timesteps 25 -n_refinements 5 -init_timestep [0|25] -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary
***************************************************************
 
Usage:
	 ./ex10-opt -init_timestep 0
OR
	 ./ex10-opt -read_solution -init_timestep 26

Running: ./ex10-opt -n_timesteps 25 -n_refinements 5 -output_freq 10 -init_timestep 0 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary

 Mesh Information:
  mesh_dimension()=2
  spatial_dimension()=3
  n_nodes()=6273
    n_local_nodes()=3171
  n_elem()=13650
    n_local_elem()=6839
    n_active_elem()=10240
  n_subdomains()=1
  n_processors()=2
  processor_id()=0

 EquationSystems
  n_systems()=1
   System "Convection-Diffusion"
    Type "TransientLinearImplicit"
    Variables="u" 
    Finite Element Types="LAGRANGE" 
    Approximation Orders="FIRST" 
    n_dofs()=6273
    n_local_dofs()=3171
    n_constrained_dofs()=0
    n_vectors()=3

 Solving time step  0, time=0.0250...
H1 norm = 1.58843
  Refining the mesh...
H1 norm = 1.58839
 Solving time step  1, time=0.0500...
H1 norm = 1.46061
  Refining the mesh...
H1 norm = 1.45992
 Solving time step  2, time=0.0750...
H1 norm = 1.35107
  Refining the mesh...
H1 norm = 1.35069
 Solving time step  3, time=0.1000...
H1 norm = 1.25698
  Refining the mesh...
H1 norm = 1.25636
 Solving time step  4, time=0.1250...
H1 norm = 1.17465
  Refining the mesh...
H1 norm = 1.1744
 Solving time step  5, time=0.1500...
H1 norm = 1.10273
  Refining the mesh...
H1 norm = 1.10224
 Solving time step  6, time=0.1750...
H1 norm = 1.03876
  Refining the mesh...
H1 norm = 1.03853
 Solving time step  7, time=0.2000...
H1 norm = 0.981437
  Refining the mesh...
H1 norm = 0.981583
 Solving time step  8, time=0.2250...
H1 norm = 0.930545
  Refining the mesh...
H1 norm = 0.930546
 Solving time step  9, time=0.2500...
H1 norm = 0.884244
  Refining the mesh...
H1 norm = 0.884219
 Solving time step 10, time=0.2750...
H1 norm = 0.842798
  Refining the mesh...
H1 norm = 0.842573
 Solving time step 11, time=0.3000...
H1 norm = 0.8047
  Refining the mesh...
H1 norm = 0.804681
 Solving time step 12, time=0.3250...
H1 norm = 0.770235
  Refining the mesh...
H1 norm = 0.770234
 Solving time step 13, time=0.3500...
H1 norm = 0.73834
  Refining the mesh...
H1 norm = 0.738326
 Solving time step 14, time=0.3750...
H1 norm = 0.709153
  Refining the mesh...
H1 norm = 0.709047
 Solving time step 15, time=0.4000...
H1 norm = 0.682057
  Refining the mesh...
H1 norm = 0.682047
 Solving time step 16, time=0.4250...
H1 norm = 0.657029
  Refining the mesh...
H1 norm = 0.657062
 Solving time step 17, time=0.4500...
H1 norm = 0.633825
  Refining the mesh...
H1 norm = 0.633816
 Solving time step 18, time=0.4750...
H1 norm = 0.612257
  Refining the mesh...
H1 norm = 0.61217
 Solving time step 19, time=0.5000...
H1 norm = 0.59211
  Refining the mesh...
H1 norm = 0.59199
 Solving time step 20, time=0.5250...
H1 norm = 0.573032
  Refining the mesh...
H1 norm = 0.573031
 Solving time step 21, time=0.5500...
H1 norm = 0.555293
  Refining the mesh...
H1 norm = 0.555287
 Solving time step 22, time=0.5750...
H1 norm = 0.538735
  Refining the mesh...
H1 norm = 0.53873
 Solving time step 23, time=0.6000...
H1 norm = 0.523123
  Refining the mesh...
H1 norm = 0.523115
 Solving time step 24, time=0.6250...
H1 norm = 0.508358
  Refining the mesh...
H1 norm = 0.508359
Final H1 norm = 0.508359

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

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

./ex10-opt on a gcc-4.5-l named daedalus with 2 processors, by roystgnr Tue Feb 22 12:19:31 2011
Using Petsc Release Version 3.1.0, Patch 5, Mon Sep 27 11:51:54 CDT 2010

                         Max       Max/Min        Avg      Total 
Time (sec):           9.241e-01      1.00223   9.231e-01
Objects:              2.766e+03      1.00000   2.766e+03
Flops:                1.872e+07      1.20341   1.714e+07  3.428e+07
Flops/sec:            2.026e+07      1.20073   1.857e+07  3.713e+07
MPI Messages:         2.297e+03      1.00000   2.297e+03  4.594e+03
MPI Message Lengths:  1.901e+06      1.03837   8.124e+02  3.732e+06
MPI Reductions:       5.780e+03      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: 9.2303e-01 100.0%  3.4281e+07 100.0%  4.594e+03 100.0%  8.124e+02      100.0%  4.859e+03  84.1% 

------------------------------------------------------------------------------------------------------------------------
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

VecMDot              121 1.0 8.7690e-04 1.5 1.60e+05 1.1 0.0e+00 0.0e+00 1.2e+02  0  1  0  0  2   0  1  0  0  2   344
VecNorm              221 1.0 1.0788e-03 1.1 1.47e+05 1.1 0.0e+00 0.0e+00 2.2e+02  0  1  0  0  4   0  1  0  0  5   259
VecScale             171 1.0 9.2030e-05 1.1 5.54e+04 1.1 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0  1141
VecCopy              238 1.0 1.2040e-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
VecSet               479 1.0 1.8167e-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               88 1.0 3.7370e-03 1.0 6.51e+04 1.1 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0    33
VecMAXPY             159 1.0 1.0157e-04 1.3 2.35e+05 1.1 0.0e+00 0.0e+00 0.0e+00  0  1  0  0  0   0  1  0  0  0  4368
VecAssemblyBegin     702 1.0 4.7653e-03 1.1 0.00e+00 0.0 2.0e+02 2.9e+02 1.9e+03  0  0  4  2 33   0  0  4  2 39     0
VecAssemblyEnd       702 1.0 2.8825e-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
VecScatterBegin      924 1.0 1.2133e-03 1.1 0.00e+00 0.0 1.3e+03 8.3e+02 0.0e+00  0  0 29 30  0   0  0 29 30  0     0
VecScatterEnd        924 1.0 1.5054e-02 3.1 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  1  0  0  0  0   1  0  0  0  0     0
VecNormalize         171 1.0 1.1020e-03 1.1 1.66e+05 1.1 0.0e+00 0.0e+00 1.7e+02  0  1  0  0  3   0  1  0  0  4   286
MatMult              171 1.0 1.5715e-02 2.8 9.36e+05 1.1 3.4e+02 3.1e+02 0.0e+00  1  5  7  3  0   1  5  7  3  0   112
MatSolve             159 1.0 3.5260e-03 1.3 3.94e+06 1.2 0.0e+00 0.0e+00 0.0e+00  0 21  0  0  0   0 21  0  0  0  2053
MatLUFactorNum        50 1.0 1.6851e-02 1.2 1.32e+07 1.2 0.0e+00 0.0e+00 0.0e+00  2 70  0  0  0   2 70  0  0  0  1426
MatILUFactorSym       50 1.0 4.4662e-02 1.2 0.00e+00 0.0 0.0e+00 0.0e+00 5.0e+01  4  0  0  0  1   4  0  0  0  1     0
MatAssemblyBegin     100 1.0 2.7528e-02 5.9 0.00e+00 0.0 2.2e+02 1.2e+03 2.0e+02  2  0  5  7  3   2  0  5  7  4     0
MatAssemblyEnd       100 1.0 3.1877e-03 1.0 0.00e+00 0.0 1.0e+02 8.3e+01 2.6e+02  0  0  2  0  4   0  0  2  0  5     0
MatGetRowIJ           50 1.0 9.0599e-06 4.8 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
MatGetOrdering        50 1.0 3.5000e-04 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 1.0e+02  0  0  0  0  2   0  0  0  0  2     0
MatZeroEntries       102 1.0 1.7715e-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
KSPGMRESOrthog       121 1.0 1.0164e-03 1.4 3.21e+05 1.1 0.0e+00 0.0e+00 1.2e+02  0  2  0  0  2   0  2  0  0  2   595
KSPSetup             100 1.0 3.3569e-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              50 1.0 8.0929e-02 1.0 1.87e+07 1.2 3.4e+02 3.1e+02 4.9e+02  9100  7  3  9   9100  7  3 10   424
PCSetUp              100 1.0 6.3950e-02 1.2 1.32e+07 1.2 0.0e+00 0.0e+00 1.5e+02  6 70  0  0  3   6 70  0  0  3   376
PCSetUpOnBlocks       50 1.0 6.2954e-02 1.2 1.32e+07 1.2 0.0e+00 0.0e+00 1.5e+02  6 70  0  0  3   6 70  0  0  3   382
PCApply              159 1.0 4.4537e-03 1.3 3.94e+06 1.2 0.0e+00 0.0e+00 0.0e+00  0 21  0  0  0   0 21  0  0  0  1626
------------------------------------------------------------------------------------------------------------------------

Memory usage is given in bytes:

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

--- Event Stage 0: Main Stage

                 Vec  1090           1090      4322464     0
         Vec Scatter   506            506       439208     0
           Index Set   810            810       710416     0
   IS L to G Mapping   128            128       333688     0
              Matrix   128            128      9863268     0
       Krylov Solver    52             52       490880     0
      Preconditioner    52             52        36608     0
========================================================================================================================
Average time to get PetscTime(): 0
Average time for MPI_Barrier(): 1.19209e-06
Average time for zero size MPI_Send(): 5.48363e-06
#PETSc Option Table entries:
-init_timestep 0
-ksp_right_pc
-log_summary
-n_refinements 5
-n_timesteps 25
-output_freq 10
-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
Configure run at: Fri Oct 15 13:01:23 2010
Configure options: --with-debugging=false --COPTFLAGS=-O3 --CXXOPTFLAGS=-O3 --FOPTFLAGS=-O3 --with-clanguage=C++ --with-shared=1 --with-mpi-dir=/org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid --with-mumps=true --download-mumps=ifneeded --with-parmetis=true --download-parmetis=ifneeded --with-superlu=true --download-superlu=ifneeded --with-superludir=true --download-superlu_dist=ifneeded --with-blacs=true --download-blacs=ifneeded --with-scalapack=true --download-scalapack=ifneeded --with-hypre=true --download-hypre=ifneeded --with-blas-lib="[/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t/libmkl_intel_lp64.so,/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t/libmkl_sequential.so,/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t/libmkl_core.so]" --with-lapack-lib=/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t/libmkl_solver_lp64_sequential.a
-----------------------------------------
Libraries compiled on Fri Oct 15 13:01:23 CDT 2010 on atreides 
Machine characteristics: Linux atreides 2.6.32-25-generic #44-Ubuntu SMP Fri Sep 17 20:05:27 UTC 2010 x86_64 GNU/Linux 
Using PETSc directory: /org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5
Using PETSc arch: gcc-4.5-lucid-mpich2-1.2.1-cxx-opt
-----------------------------------------
Using C compiler: /org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/bin/mpicxx -Wall -Wwrite-strings -Wno-strict-aliasing -O3   -fPIC   
Using Fortran compiler: /org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/bin/mpif90 -fPIC -Wall -Wno-unused-variable -O3    
-----------------------------------------
Using include paths: -I/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/include -I/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/include -I/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/include -I/org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/include  
------------------------------------------
Using C linker: /org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/bin/mpicxx -Wall -Wwrite-strings -Wno-strict-aliasing -O3 
Using Fortran linker: /org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/bin/mpif90 -fPIC -Wall -Wno-unused-variable -O3  
Using libraries: -Wl,-rpath,/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/lib -L/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/lib -lpetsc       -lX11 -Wl,-rpath,/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/lib -L/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/lib -lHYPRE -lsuperlu_dist_2.4 -lcmumps -ldmumps -lsmumps -lzmumps -lmumps_common -lpord -lparmetis -lmetis -lscalapack -lblacs -lsuperlu_4.0 -Wl,-rpath,/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t -L/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t -lmkl_solver_lp64_sequential -lmkl_intel_lp64 -lmkl_sequential -lmkl_core -lm -Wl,-rpath,/org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/lib -L/org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/lib -Wl,-rpath,/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib/gcc/x86_64-unknown-linux-gnu/4.5.1 -L/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib/gcc/x86_64-unknown-linux-gnu/4.5.1 -Wl,-rpath,/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib64 -L/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib64 -Wl,-rpath,/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib -L/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib -ldl -lmpich -lopa -lpthread -lrt -lgcc_s -lmpichf90 -lgfortran -lm -lm -lmpichcxx -lstdc++ -ldl -lmpich -lopa -lpthread -lrt -lgcc_s -ldl  
------------------------------------------
 
***** Finished first 25 steps, now read in saved solution and continue *****
 
Usage:
	 ./ex10-opt -init_timestep 0
OR
	 ./ex10-opt -read_solution -init_timestep 26

Running: ./ex10-opt -read_solution -n_timesteps 25 -output_freq 10 -init_timestep 25 -pc_type bjacobi -sub_pc_type ilu -sub_pc_factor_levels 4 -sub_pc_factor_zeropivot 0 -ksp_right_pc -log_summary

 Mesh Information:
  mesh_dimension()=2
  spatial_dimension()=3
  n_nodes()=721
    n_local_nodes()=396
  n_elem()=1030
    n_local_elem()=500
    n_active_elem()=775
  n_subdomains()=1
  n_processors()=2
  processor_id()=0

Initial H1 norm = 0.508359

 EquationSystems
  n_systems()=1
   System "Convection-Diffusion"
    Type "TransientLinearImplicit"
    Variables="u" 
    Finite Element Types="LAGRANGE" 
    Approximation Orders="FIRST" 
    n_dofs()=721
    n_local_dofs()=396
    n_constrained_dofs()=123
    n_vectors()=3

 Solving time step 25, time=0.6500...
H1 norm = 0.494397
  Refining the mesh...
H1 norm = 0.494393
 Solving time step 26, time=0.6750...
H1 norm = 0.481245
  Refining the mesh...
H1 norm = 0.481242
 Solving time step 27, time=0.7000...
H1 norm = 0.468821
  Refining the mesh...
H1 norm = 0.468719
 Solving time step 28, time=0.7250...
H1 norm = 0.456941
  Refining the mesh...
H1 norm = 0.456854
 Solving time step 29, time=0.7500...
H1 norm = 0.445682
  Refining the mesh...
H1 norm = 0.445682
 Solving time step 30, time=0.7750...
H1 norm = 0.435013
  Refining the mesh...
H1 norm = 0.435006
 Solving time step 31, time=0.8000...
H1 norm = 0.424893
  Refining the mesh...
H1 norm = 0.424892
 Solving time step 32, time=0.8250...
H1 norm = 0.415206
  Refining the mesh...
H1 norm = 0.415202
 Solving time step 33, time=0.8500...
H1 norm = 0.405957
  Refining the mesh...
H1 norm = 0.405957
 Solving time step 34, time=0.8750...
H1 norm = 0.397088
  Refining the mesh...
H1 norm = 0.397044
 Solving time step 35, time=0.9000...
H1 norm = 0.388558
  Refining the mesh...
H1 norm = 0.388557
 Solving time step 36, time=0.9250...
H1 norm = 0.380394
  Refining the mesh...
H1 norm = 0.380392
 Solving time step 37, time=0.9500...
H1 norm = 0.372562
  Refining the mesh...
H1 norm = 0.37254
 Solving time step 38, time=0.9750...
H1 norm = 0.365056
  Refining the mesh...
H1 norm = 0.365054
 Solving time step 39, time=1.0000...
H1 norm = 0.357873
  Refining the mesh...
H1 norm = 0.357836
 Solving time step 40, time=1.0300...
H1 norm = 0.350969
  Refining the mesh...
H1 norm = 0.350967
 Solving time step 41, time=1.0500...
H1 norm = 0.344366
  Refining the mesh...
H1 norm = 0.344366
 Solving time step 42, time=1.0700...
H1 norm = 0.338035
  Refining the mesh...
H1 norm = 0.338033
 Solving time step 43, time=1.1000...
H1 norm = 0.33193
  Refining the mesh...
H1 norm = 0.331907
 Solving time step 44, time=1.1200...
H1 norm = 0.326024
  Refining the mesh...
H1 norm = 0.326002
 Solving time step 45, time=1.1500...
H1 norm = 0.32032
  Refining the mesh...
H1 norm = 0.320298
 Solving time step 46, time=1.1700...
H1 norm = 0.314798
  Refining the mesh...
H1 norm = 0.31478
 Solving time step 47, time=1.2000...
H1 norm = 0.309464
  Refining the mesh...
H1 norm = 0.309462
 Solving time step 48, time=1.2200...
H1 norm = 0.30432
  Refining the mesh...
H1 norm = 0.30432
 Solving time step 49, time=1.2500...
H1 norm = 0.299319
  Refining the mesh...
H1 norm = 0.299323
************************************************************************************************************************
***             WIDEN YOUR WINDOW TO 120 CHARACTERS.  Use 'enscript -r -fCourier9' to print this document            ***
************************************************************************************************************************

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

./ex10-opt on a gcc-4.5-l named daedalus with 2 processors, by roystgnr Tue Feb 22 12:19:33 2011
Using Petsc Release Version 3.1.0, Patch 5, Mon Sep 27 11:51:54 CDT 2010

                         Max       Max/Min        Avg      Total 
Time (sec):           1.384e+00      1.00302   1.382e+00
Objects:              2.808e+03      1.00000   2.808e+03
Flops:                4.433e+07      1.06635   4.295e+07  8.591e+07
Flops/sec:            3.203e+07      1.06315   3.108e+07  6.215e+07
MPI Messages:         2.360e+03      1.00000   2.360e+03  4.719e+03
MPI Message Lengths:  2.598e+06      1.02201   1.089e+03  5.140e+06
MPI Reductions:       5.866e+03      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: 1.3821e+00 100.0%  8.5905e+07 100.0%  4.719e+03 100.0%  1.089e+03      100.0%  4.945e+03  84.3% 

------------------------------------------------------------------------------------------------------------------------
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

VecMDot              161 1.0 1.3199e-03 1.6 5.83e+05 1.1 0.0e+00 0.0e+00 1.6e+02  0  1  0  0  3   0  1  0  0  3   861
VecNorm              261 1.0 1.7569e-03 1.1 3.20e+05 1.1 0.0e+00 0.0e+00 2.6e+02  0  1  0  0  4   0  1  0  0  5   354
VecScale             211 1.0 1.2040e-04 1.0 1.29e+05 1.1 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0  2095
VecCopy              235 1.0 1.2803e-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
VecSet               513 1.0 1.9908e-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               85 1.0 9.5367e-05 1.0 1.04e+05 1.1 0.0e+00 0.0e+00 0.0e+00  0  0  0  0  0   0  0  0  0  0  2118
VecMAXPY             196 1.0 2.4939e-04 1.1 7.82e+05 1.1 0.0e+00 0.0e+00 0.0e+00  0  2  0  0  0   0  2  0  0  0  6108
VecAssemblyBegin     703 1.0 5.2691e-03 1.1 0.00e+00 0.0 2.1e+02 4.4e+02 1.9e+03  0  0  4  2 32   0  0  4  2 38     0
VecAssemblyEnd       703 1.0 3.2187e-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
VecScatterBegin      966 1.0 1.3852e-03 1.1 0.00e+00 0.0 1.4e+03 1.0e+03 0.0e+00  0  0 30 28  0   0  0 30 28  0     0
VecScatterEnd        966 1.0 1.5233e-02 2.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00  1  0  0  0  0   1  0  0  0  0     0
VecNormalize         211 1.0 1.8098e-03 1.1 3.88e+05 1.1 0.0e+00 0.0e+00 2.1e+02  0  1  0  0  4   0  1  0  0  4   418
MatMult              211 1.0 1.6989e-02 1.8 2.03e+06 1.0 4.2e+02 4.3e+02 0.0e+00  1  5  9  4  0   1  5  9  4  0   234
MatSolve             196 1.0 7.1387e-03 1.1 1.02e+07 1.1 0.0e+00 0.0e+00 0.0e+00  0 23  0  0  0   0 23  0  0  0  2778
MatLUFactorNum        50 1.0 3.5021e-02 1.1 3.02e+07 1.1 0.0e+00 0.0e+00 0.0e+00  2 68  0  0  0   2 68  0  0  0  1666
MatILUFactorSym       50 1.0 9.1278e-02 1.1 0.00e+00 0.0 0.0e+00 0.0e+00 5.0e+01  6  0  0  0  1   6  0  0  0  1     0
MatAssemblyBegin     100 1.0 1.1788e-02 2.2 0.00e+00 0.0 2.4e+02 1.6e+03 2.0e+02  1  0  5  8  3   1  0  5  8  4     0
MatAssemblyEnd       100 1.0 4.2770e-03 1.1 0.00e+00 0.0 1.0e+02 1.1e+02 2.6e+02  0  0  2  0  4   0  0  2  0  5     0
MatGetRowIJ           50 1.0 1.1444e-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
MatGetOrdering        50 1.0 3.8862e-04 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 1.0e+02  0  0  0  0  2   0  0  0  0  2     0
MatZeroEntries       102 1.0 2.0480e-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
KSPGMRESOrthog       161 1.0 1.6124e-03 1.4 1.17e+06 1.1 0.0e+00 0.0e+00 1.6e+02  0  3  0  0  3   0  3  0  0  3  1410
KSPSetup             100 1.0 3.2258e-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              50 1.0 1.5128e-01 1.0 4.43e+07 1.1 4.2e+02 4.3e+02 5.7e+02 11100  9  4 10  11100  9  4 12   568
PCSetUp              100 1.0 1.2879e-01 1.1 3.02e+07 1.1 0.0e+00 0.0e+00 1.5e+02  9 68  0  0  3   9 68  0  0  3   453
PCSetUpOnBlocks       50 1.0 1.2778e-01 1.1 3.02e+07 1.1 0.0e+00 0.0e+00 1.5e+02  9 68  0  0  3   9 68  0  0  3   457
PCApply              196 1.0 8.3489e-03 1.1 1.02e+07 1.1 0.0e+00 0.0e+00 0.0e+00  1 23  0  0  0   1 23  0  0  0  2375
------------------------------------------------------------------------------------------------------------------------

Memory usage is given in bytes:

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

--- Event Stage 0: Main Stage

                 Vec  1130           1130      6354688     0
         Vec Scatter   507            507       440076     0
           Index Set   811            811       864976     0
   IS L to G Mapping   128            128       434808     0
              Matrix   128            128     17887820     0
       Krylov Solver    52             52       490880     0
      Preconditioner    52             52        36608     0
========================================================================================================================
Average time to get PetscTime(): 0
Average time for MPI_Barrier(): 5.72205e-07
Average time for zero size MPI_Send(): 6.07967e-06
#PETSc Option Table entries:
-init_timestep 25
-ksp_right_pc
-log_summary
-n_timesteps 25
-output_freq 10
-pc_type bjacobi
-read_solution
-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
Configure run at: Fri Oct 15 13:01:23 2010
Configure options: --with-debugging=false --COPTFLAGS=-O3 --CXXOPTFLAGS=-O3 --FOPTFLAGS=-O3 --with-clanguage=C++ --with-shared=1 --with-mpi-dir=/org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid --with-mumps=true --download-mumps=ifneeded --with-parmetis=true --download-parmetis=ifneeded --with-superlu=true --download-superlu=ifneeded --with-superludir=true --download-superlu_dist=ifneeded --with-blacs=true --download-blacs=ifneeded --with-scalapack=true --download-scalapack=ifneeded --with-hypre=true --download-hypre=ifneeded --with-blas-lib="[/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t/libmkl_intel_lp64.so,/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t/libmkl_sequential.so,/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t/libmkl_core.so]" --with-lapack-lib=/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t/libmkl_solver_lp64_sequential.a
-----------------------------------------
Libraries compiled on Fri Oct 15 13:01:23 CDT 2010 on atreides 
Machine characteristics: Linux atreides 2.6.32-25-generic #44-Ubuntu SMP Fri Sep 17 20:05:27 UTC 2010 x86_64 GNU/Linux 
Using PETSc directory: /org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5
Using PETSc arch: gcc-4.5-lucid-mpich2-1.2.1-cxx-opt
-----------------------------------------
Using C compiler: /org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/bin/mpicxx -Wall -Wwrite-strings -Wno-strict-aliasing -O3   -fPIC   
Using Fortran compiler: /org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/bin/mpif90 -fPIC -Wall -Wno-unused-variable -O3    
-----------------------------------------
Using include paths: -I/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/include -I/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/include -I/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/include -I/org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/include  
------------------------------------------
Using C linker: /org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/bin/mpicxx -Wall -Wwrite-strings -Wno-strict-aliasing -O3 
Using Fortran linker: /org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/bin/mpif90 -fPIC -Wall -Wno-unused-variable -O3  
Using libraries: -Wl,-rpath,/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/lib -L/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/lib -lpetsc       -lX11 -Wl,-rpath,/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/lib -L/org/centers/pecos/LIBRARIES/PETSC3/petsc-3.1-p5/gcc-4.5-lucid-mpich2-1.2.1-cxx-opt/lib -lHYPRE -lsuperlu_dist_2.4 -lcmumps -ldmumps -lsmumps -lzmumps -lmumps_common -lpord -lparmetis -lmetis -lscalapack -lblacs -lsuperlu_4.0 -Wl,-rpath,/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t -L/org/centers/pecos/LIBRARIES/MKL/mkl-10.0.3.020-gcc-4.5-lucid/lib/em64t -lmkl_solver_lp64_sequential -lmkl_intel_lp64 -lmkl_sequential -lmkl_core -lm -Wl,-rpath,/org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/lib -L/org/centers/pecos/LIBRARIES/MPICH2/mpich2-1.2.1-gcc-4.5-lucid/lib -Wl,-rpath,/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib/gcc/x86_64-unknown-linux-gnu/4.5.1 -L/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib/gcc/x86_64-unknown-linux-gnu/4.5.1 -Wl,-rpath,/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib64 -L/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib64 -Wl,-rpath,/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib -L/org/centers/pecos/LIBRARIES/GCC/gcc-4.5.1-lucid/lib -ldl -lmpich -lopa -lpthread -lrt -lgcc_s -lmpichf90 -lgfortran -lm -lm -lmpichcxx -lstdc++ -ldl -lmpich -lopa -lpthread -lrt -lgcc_s -ldl  
------------------------------------------
 
***************************************************************
* Done Running Example  mpirun -np 2 ./ex10-opt [-read_solution] -n_timesteps 25 -init_timestep [0|25] -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: December 02 2011 21:19:02 UTC

Hosted By:
SourceForge.net Logo