Example 4 - Solving a 1D, 2D or 3D Poisson Problem in Parallel



This is the fourth example program. It builds on the third example program by showing how to formulate the code in a dimension-independent way. Very minor changes to the example will allow the problem to be solved in one, two or three dimensions.

This example will also introduce the PerfLog class as a way to monitor your code's performance. We will use it to instrument the matrix assembly code and look for bottlenecks where we should focus optimization efforts.

This example also shows how to extend example 3 to run in parallel. Notice how litte has changed! The significant differences are marked with "PARALLEL CHANGE".



C++ include files that we need
        #include <iostream>
        #include <algorithm>
        #include <math.h>
        
Basic include file needed for the mesh functionality.
        #include "libmesh.h"
        #include "mesh.h"
        #include "mesh_generation.h"
        #include "gmv_io.h"
        #include "gnuplot_io.h"
        #include "linear_implicit_system.h"
        #include "equation_systems.h"
        
Define the Finite Element object.
        #include "fe.h"
        
Define Gauss quadrature rules.
        #include "quadrature_gauss.h"
        
Define the DofMap, which handles degree of freedom indexing.
        #include "dof_map.h"
        
Define useful datatypes for finite element matrix and vector components.
        #include "sparse_matrix.h"
        #include "numeric_vector.h"
        #include "dense_matrix.h"
        #include "dense_vector.h"
        
Define the PerfLog, a performance logging utility. It is useful for timing events in a code and giving you an idea where bottlenecks lie.
        #include "perf_log.h"
        
The definition of a geometric element
        #include "elem.h"
        
        #include "string_to_enum.h"
        #include "getpot.h"
        
        
        
Function prototype. This is the function that will assemble the linear system for our Poisson problem. Note that the function will take the \p EquationSystems object and the name of the system we are assembling as input. From the \p EquationSystems object we have acess to the \p Mesh and other objects we might need.
        void assemble_poisson(EquationSystems& es,
                              const std::string& system_name);
        
Exact solution function prototype.
        Real exact_solution (const Real x,
        		     const Real y = 0.,
        		     const Real z = 0.);
        
Begin the main program.
        int main (int argc, char** argv)
        {
Initialize libMesh and any dependent libaries, like in example 2.
          libMesh::init (argc, argv);
        
Declare a performance log for the main program PerfLog perf_main("Main Program");

Braces are used to force object scope, like in example 2
          {
Create a GetPot object to parse the command line
            GetPot command_line (argc, argv);
            
Check for proper calling arguments.
            if (argc < 3)
              {
        	std::cerr << "Usage:\n"
        		  <<"\t " << argv[0] << " -d 2(3)" << " -n 15"
        		  << std::endl;
        
This handy function will print the file name, line number, and then abort. Currrently the library does not use C++ exception handling.
                error();
              }
            
Brief message to the user regarding the program name and command line arguments.
            else 
              {
        	std::cout << "Running " << argv[0];
        	
        	for (int i=1; i<argc; i++)
        	  std::cout << " " << argv[i];
        	
        	std::cout << std::endl << std::endl;
              }
            
        
Read problem dimension from command line. Use int instead of unsigned since the GetPot overload is ambiguous otherwise.
            int dim = 2;
            if ( command_line.search(1, "-d") )
              dim = command_line.next(dim);
            
Read number of elements from command line
            int ps = 15;
            if ( command_line.search(1, "-n") )
              ps = command_line.next(ps);
            
Read FE order from command line
            std::string order = "SECOND"; 
            if ( command_line.search(2, "-Order", "-o") )
              order = command_line.next(order);
        
Read FE Family from command line
            std::string family = "LAGRANGE"; 
            if ( command_line.search(2, "-FEFamily", "-f") )
              family = command_line.next(family);
            
Cannot use dicontinuous basis.
            if ((family == "MONOMIAL") || (family == "XYZ"))
              {
        	std::cerr << "ex4 currently requires a C^0 (or higher) FE basis." << std::endl;
        	error();
              }
              
Create a mesh with user-defined dimension.
            Mesh mesh (dim);
            
        
Use the MeshTools::Generation mesh generator to create a uniform grid on the square [-1,1]^D. We instruct the mesh generator to build a mesh of 8x8 \p Quad9 elements in 2D, or \p Hex27 elements in 3D. Building these higher-order elements allows us to use higher-order approximation, as in example 3.
            if ((family == "LAGRANGE") && (order == "FIRST"))
              {
No reason to use high-order geometric elements if we are solving with low-order finite elements.
                MeshTools::Generation::build_cube (mesh,
        					   ps, ps, ps,
        					   -1., 1.,
        					   -1., 1.,
        					   -1., 1.,
        					   (dim==1)    ? EDGE2 : 
        					   ((dim == 2) ? QUAD4 : HEX8));
              }
            
            else
              {
        	MeshTools::Generation::build_cube (mesh,
        					   ps, ps, ps,
        					   -1., 1.,
        					   -1., 1.,
        					   -1., 1.,
        					   (dim==1)    ? EDGE3 : 
        					   ((dim == 2) ? QUAD9 : HEX27));
              }
        
            
Print information about the mesh to the screen.
            mesh.print_info();
            
            
Create an equation systems object.
            EquationSystems equation_systems (mesh);
            
Declare the system and its variables.
            {
Creates a system named "Poisson"
              LinearImplicitSystem& system =
        	equation_systems.add_system<LinearImplicitSystem> ("Poisson");
        
              
Adds the variable "u" to "Poisson". "u" will be approximated using second-order approximation.
              system.add_variable("u",
        			  Utility::string_to_enum<Order>   (order),
        			  Utility::string_to_enum<FEFamily>(family));
        
Give the system a pointer to the matrix assembly function.
              system.attach_assemble_function (assemble_poisson);
              
Initialize the data structures for the equation system.
              equation_systems.init();
        
Prints information about the system to the screen.
              equation_systems.print_info();
            }
        
Solve the system "Poisson", just like example 2.
            equation_systems.get_system("Poisson").solve();
        
After solving the system write the solution to a GMV-formatted plot file.
            if(dim == 1)
            {        
              GnuPlotIO plot(mesh,"Example 4, 1D",GnuPlotIO::GRID_ON);
              plot.write_equation_systems("out_1",equation_systems);
            }
            else
            {
              GMVIO (mesh).write_equation_systems ((dim == 3) ? 
                "out_3.gmv" : "out_2.gmv",equation_systems);
            }
          }
          
All done.
          return libMesh::close ();
        }
        
        
        
        






We now define the matrix assembly function for the Poisson system. We need to first compute element matrices and right-hand sides, and then take into account the boundary conditions, which will be handled via a penalty method.
        void assemble_poisson(EquationSystems& es,
                              const std::string& system_name)
        {
It is a good idea to make sure we are assembling the proper system.
          assert (system_name == "Poisson");
        
        
Declare a performance log. Give it a descriptive string to identify what part of the code we are logging, since there may be many PerfLogs in an application.
          PerfLog perf_log ("Matrix Assembly");
          
Get a constant reference to the mesh object.
          const Mesh& mesh = es.get_mesh();
        
The dimension that we are running
          const unsigned int dim = mesh.mesh_dimension();
        
Get a reference to the LinearImplicitSystem we are solving
          LinearImplicitSystem& system = es.get_system<LinearImplicitSystem>("Poisson");
          
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();
        
Get a constant reference to the Finite Element type for the first (and only) variable in the system.
          FEType fe_type = dof_map.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));
          
A 5th order Gauss quadrature rule for numerical integration.
          QGauss qrule (dim, FIFTH);
        
Tell the finite element object to use our quadrature rule.
          fe->attach_quadrature_rule (&qrule);
        
Declare a special finite element object for boundary integration.
          AutoPtr<FEBase> fe_face (FEBase::build(dim, fe_type));
        	      
Boundary integration requires one quadraure rule, with dimensionality one less than the dimensionality of the element.
          QGauss qface(dim-1, FIFTH);
          
Tell the finte element object to use our quadrature rule.
          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 begin with the element Jacobian * quadrature weight at each integration point.
          const std::vector<Real>& JxW = fe->get_JxW();
        
The physical XY locations of the quadrature points on the element. These might be useful for evaluating spatially varying material properties at the quadrature points.
          const std::vector<Point>& q_point = fe->get_xyz();
        
The element shape functions evaluated at the quadrature points.
          const std::vector<std::vector<Real> >& phi = fe->get_phi();
        
The element shape function gradients evaluated at the quadrature points.
          const std::vector<std::vector<RealGradient> >& dphi = fe->get_dphi();
        
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". More detail is in example 3.
          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;
        
Now we will loop over all the elements in the mesh. We will compute the element matrix and right-hand-side contribution. See example 3 for a discussion of the element iterators. Here we use the \p const_local_elem_iterator to indicate we only want to loop over elements that are assigned to the local processor. This allows each processor to compute its components of the global matrix.

"PARALLEL CHANGE"
          MeshBase::const_element_iterator       el     = mesh.local_elements_begin();
          const MeshBase::const_element_iterator end_el = mesh.local_elements_end();
        
          for ( ; el != end_el; ++el)
            {
Start logging the shape function initialization. This is done through a simple function call with the name of the event to log.
              perf_log.start_event("elem init");      
        
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());
        
Stop logging the shape function initialization. If you forget to stop logging an event the PerfLog object will probably catch the error and abort.
              perf_log.stop_event("elem init");      
        
Now we will build the element matrix. This involves a double loop to integrate the test funcions (i) against the trial functions (j).

We have split the numeric integration into two loops so that we can log the matrix and right-hand-side computation seperately.

Now start logging the element matrix computation
              perf_log.start_event ("Ke");
        
              for (unsigned int qp=0; qp<qrule.n_points(); qp++)
        	for (unsigned int i=0; i<phi.size(); i++)
        	  for (unsigned int j=0; j<phi.size(); j++)
        	    Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
        	    
        
Stop logging the matrix computation
              perf_log.stop_event ("Ke");
        
Now we build the element right-hand-side contribution. This involves a single loop in which we integrate the "forcing function" in the PDE against the test functions.

Start logging the right-hand-side computation
              perf_log.start_event ("Fe");
              
              for (unsigned int qp=0; qp<qrule.n_points(); qp++)
        	{
fxy is the forcing function for the Poisson equation. In this case we set fxy to be a finite difference Laplacian approximation to the (known) exact solution.

We will use the second-order accurate FD Laplacian approximation, which in 2D on a structured grid is

u_xx + u_yy = (u(i-1,j) + u(i+1,j) + u(i,j-1) + u(i,j+1) + -4*u(i,j))/h^2

Since the value of the forcing function depends only on the location of the quadrature point (q_point[qp]) we will compute it here, outside of the i-loop
                  const Real x = q_point[qp](0);
        	  const Real y = q_point[qp](1);
        	  const Real z = q_point[qp](2);
        	  const Real eps = 1.e-3;
        
        	  const Real uxx = (exact_solution(x-eps,y,z) +
        			    exact_solution(x+eps,y,z) +
        			    -2.*exact_solution(x,y,z))/eps/eps;
        	      
        	  const Real uyy = (exact_solution(x,y-eps,z) +
        			    exact_solution(x,y+eps,z) +
        			    -2.*exact_solution(x,y,z))/eps/eps;
        	  
        	  const Real uzz = (exact_solution(x,y,z-eps) +
        			    exact_solution(x,y,z+eps) +
        			    -2.*exact_solution(x,y,z))/eps/eps;
        
                  Real fxy;
                  if(dim==1)
                  {
In 1D, compute the rhs by differentiating the exact solution twice.
                    const Real pi = libMesh::pi;
                    fxy = (0.25*pi*pi)*sin(.5*pi*x);
                  }
                  else
                  {
        	    fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
                  } 
        
Add the RHS contribution
                  for (unsigned int i=0; i<phi.size(); i++)
        	    Fe(i) += JxW[qp]*fxy*phi[i][qp];	  
        	}
              
Stop logging the right-hand-side computation
              perf_log.stop_event ("Fe");
        
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. This is discussed at length in example 3.
              {
        	
Start logging the boundary condition computation
                perf_log.start_event ("BCs");
        
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 side=0; side<elem->n_sides(); side++)
        	  if (elem->neighbor(side) == NULL)
        	    {
                    
The penalty value. \frac{1}{\epsilon} in the discussion above.
                      const Real penalty = 1.e10;
        
The value of the shape functions at the quadrature points.
                      const std::vector<std::vector<Real> >&  phi_face = fe_face->get_phi();
        
The Jacobian * Quadrature Weight at the quadrature points on the face.
                      const std::vector<Real>& JxW_face = fe_face->get_JxW();
        
The XYZ locations (in physical space) of the quadrature points on the face. This is where we will interpolate the boundary value function.
                      const std::vector<Point >& qface_point = fe_face->get_xyz();
        
Compute the shape function values on the element face.
                      fe_face->reinit(elem, side);
        
Loop over the face quadrature points for integration.
                      for (unsigned int qp=0; qp<qface.n_points(); qp++)
                      {
The location on the boundary of the current face quadrature point.
                        const Real xf = qface_point[qp](0);
                        const Real yf = qface_point[qp](1);
                        const Real zf = qface_point[qp](2);
        
        
The boundary value.
                        const Real value = exact_solution(xf, yf, zf);
        
Matrix contribution of the L2 projection.
                        for (unsigned int i=0; i<phi_face.size(); i++)
                          for (unsigned int j=0; j<phi_face.size(); j++)
                            Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
        
Right-hand-side contribution of the L2 projection.
                        for (unsigned int i=0; i<phi_face.size(); i++)
                          Fe(i) += JxW_face[qp]*penalty*value*phi_face[i][qp];
                      } 
                    }
                    
        	
Stop logging the boundary condition computation
                perf_log.stop_event ("BCs");
              } 
              
        
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 PetscMatrix::add_matrix() and \p PetscVector::add_vector() members do this for us. Start logging the insertion of the local (element) matrix and vector into the global matrix and vector
              perf_log.start_event ("matrix insertion");
              
              system.matrix->add_matrix (Ke, dof_indices);
              system.rhs->add_vector    (Fe, dof_indices);
        
Start logging the insertion of the local (element) matrix and vector into the global matrix and vector
              perf_log.stop_event ("matrix insertion");
            }
        
That's it. We don't need to do anything else to the PerfLog. When it goes out of scope (at this function return) it will print its log to the screen. Pretty easy, huh?
        }



The program without comments:

 
  
  
  #include <iostream>
  #include <algorithm>
  #include <math.h>
  
  #include "libmesh.h"
  #include "mesh.h"
  #include "mesh_generation.h"
  #include "gmv_io.h"
  #include "gnuplot_io.h"
  #include "linear_implicit_system.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 "perf_log.h"
  
  #include "elem.h"
  
  #include "string_to_enum.h"
  #include "getpot.h"
  
  
  
  void assemble_poisson(EquationSystems& es,
                        const std::string& system_name);
  
  Real exact_solution (const Real x,
  		     const Real y = 0.,
  		     const Real z = 0.);
  
  int main (int argc, char** argv)
  {
    libMesh::init (argc, argv);
  
    
    {
      GetPot command_line (argc, argv);
      
      if (argc < 3)
        {
  	std::cerr << "Usage:\n"
  		  <<"\t " << argv[0] << " -d 2(3)" << " -n 15"
  		  << std::endl;
  
  	error();
        }
      
      else 
        {
  	std::cout << "Running " << argv[0];
  	
  	for (int i=1; i<argc; i++)
  	  std::cout << " " << argv[i];
  	
  	std::cout << std::endl << std::endl;
        }
      
  
      int dim = 2;
      if ( command_line.search(1, "-d") )
        dim = command_line.next(dim);
      
      int ps = 15;
      if ( command_line.search(1, "-n") )
        ps = command_line.next(ps);
      
      std::string order = "SECOND"; 
      if ( command_line.search(2, "-Order", "-o") )
        order = command_line.next(order);
  
      std::string family = "LAGRANGE"; 
      if ( command_line.search(2, "-FEFamily", "-f") )
        family = command_line.next(family);
      
      if ((family == "MONOMIAL") || (family == "XYZ"))
        {
  	std::cerr << "ex4 currently requires a C^0 (or higher) FE basis." << std::endl;
  	error();
        }
        
      Mesh mesh (dim);
      
  
      if ((family == "LAGRANGE") && (order == "FIRST"))
        {
  	MeshTools::Generation::build_cube (mesh,
  					   ps, ps, ps,
  					   -1., 1.,
  					   -1., 1.,
  					   -1., 1.,
  					   (dim==1)    ? EDGE2 : 
  					   ((dim == 2) ? QUAD4 : HEX8));
        }
      
      else
        {
  	MeshTools::Generation::build_cube (mesh,
  					   ps, ps, ps,
  					   -1., 1.,
  					   -1., 1.,
  					   -1., 1.,
  					   (dim==1)    ? EDGE3 : 
  					   ((dim == 2) ? QUAD9 : HEX27));
        }
  
      
      mesh.print_info();
      
      
      EquationSystems equation_systems (mesh);
      
      {
        LinearImplicitSystem& system =
  	equation_systems.add_system<LinearImplicitSystem> ("Poisson");
  
        
        system.add_variable("u",
  			  Utility::string_to_enum<Order>   (order),
  			  Utility::string_to_enum<FEFamily>(family));
  
        system.attach_assemble_function (assemble_poisson);
        
        equation_systems.init();
  
        equation_systems.print_info();
      }
  
      equation_systems.get_system("Poisson").solve();
  
      if(dim == 1)
      {        
        GnuPlotIO plot(mesh,"Example 4, 1D",GnuPlotIO::GRID_ON);
        plot.write_equation_systems("out_1",equation_systems);
      }
      else
      {
        GMVIO (mesh).write_equation_systems ((dim == 3) ? 
          "out_3.gmv" : "out_2.gmv",equation_systems);
      }
    }
    
    return libMesh::close ();
  }
  
  
  
  
  void assemble_poisson(EquationSystems& es,
                        const std::string& system_name)
  {
    assert (system_name == "Poisson");
  
  
    PerfLog perf_log ("Matrix Assembly");
    
    const Mesh& mesh = es.get_mesh();
  
    const unsigned int dim = mesh.mesh_dimension();
  
    LinearImplicitSystem& system = es.get_system<LinearImplicitSystem>("Poisson");
    
    const DofMap& dof_map = system.get_dof_map();
  
    FEType fe_type = dof_map.variable_type(0);
  
    AutoPtr<FEBase> fe (FEBase::build(dim, fe_type));
    
    QGauss qrule (dim, FIFTH);
  
    fe->attach_quadrature_rule (&qrule);
  
    AutoPtr<FEBase> fe_face (FEBase::build(dim, fe_type));
  	      
    QGauss qface(dim-1, FIFTH);
    
    fe_face->attach_quadrature_rule (&qface);
  
    const std::vector<Real>& JxW = fe->get_JxW();
  
    const std::vector<Point>& q_point = fe->get_xyz();
  
    const std::vector<std::vector<Real> >& phi = fe->get_phi();
  
    const std::vector<std::vector<RealGradient> >& dphi = fe->get_dphi();
  
    DenseMatrix<Number> Ke;
    DenseVector<Number> Fe;
  
    std::vector<unsigned int> dof_indices;
  
    MeshBase::const_element_iterator       el     = mesh.local_elements_begin();
    const MeshBase::const_element_iterator end_el = mesh.local_elements_end();
  
    for ( ; el != end_el; ++el)
      {
        perf_log.start_event("elem init");      
  
        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());
  
        perf_log.stop_event("elem init");      
  
        perf_log.start_event ("Ke");
  
        for (unsigned int qp=0; qp<qrule.n_points(); qp++)
  	for (unsigned int i=0; i<phi.size(); i++)
  	  for (unsigned int j=0; j<phi.size(); j++)
  	    Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
  	    
  
        perf_log.stop_event ("Ke");
  
        perf_log.start_event ("Fe");
        
        for (unsigned int qp=0; qp<qrule.n_points(); qp++)
  	{
  	  const Real x = q_point[qp](0);
  	  const Real y = q_point[qp](1);
  	  const Real z = q_point[qp](2);
  	  const Real eps = 1.e-3;
  
  	  const Real uxx = (exact_solution(x-eps,y,z) +
  			    exact_solution(x+eps,y,z) +
  			    -2.*exact_solution(x,y,z))/eps/eps;
  	      
  	  const Real uyy = (exact_solution(x,y-eps,z) +
  			    exact_solution(x,y+eps,z) +
  			    -2.*exact_solution(x,y,z))/eps/eps;
  	  
  	  const Real uzz = (exact_solution(x,y,z-eps) +
  			    exact_solution(x,y,z+eps) +
  			    -2.*exact_solution(x,y,z))/eps/eps;
  
            Real fxy;
            if(dim==1)
            {
              const Real pi = libMesh::pi;
              fxy = (0.25*pi*pi)*sin(.5*pi*x);
            }
            else
            {
  	    fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
            } 
  
  	  for (unsigned int i=0; i<phi.size(); i++)
  	    Fe(i) += JxW[qp]*fxy*phi[i][qp];	  
  	}
        
        perf_log.stop_event ("Fe");
  
        {
  	
  	perf_log.start_event ("BCs");
  
  	for (unsigned int side=0; side<elem->n_sides(); side++)
  	  if (elem->neighbor(side) == NULL)
  	    {
              
  	      const Real penalty = 1.e10;
  
                const std::vector<std::vector<Real> >&  phi_face = fe_face->get_phi();
  
                const std::vector<Real>& JxW_face = fe_face->get_JxW();
  
                const std::vector<Point >& qface_point = fe_face->get_xyz();
  
                fe_face->reinit(elem, side);
  
                for (unsigned int qp=0; qp<qface.n_points(); qp++)
                {
                  const Real xf = qface_point[qp](0);
                  const Real yf = qface_point[qp](1);
                  const Real zf = qface_point[qp](2);
  
  
                  const Real value = exact_solution(xf, yf, zf);
  
                  for (unsigned int i=0; i<phi_face.size(); i++)
                    for (unsigned int j=0; j<phi_face.size(); j++)
                      Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
  
                  for (unsigned int i=0; i<phi_face.size(); i++)
                    Fe(i) += JxW_face[qp]*penalty*value*phi_face[i][qp];
                } 
              }
              
  	
  	perf_log.stop_event ("BCs");
        } 
        
  
        perf_log.start_event ("matrix insertion");
        
        system.matrix->add_matrix (Ke, dof_indices);
        system.rhs->add_vector    (Fe, dof_indices);
  
        perf_log.stop_event ("matrix insertion");
      }
  
  }



The console output of the program:

***************************************************************
* Running Example  ./ex4-devel
***************************************************************
 
Running ./ex4-devel -d 1 -n 20

 Mesh Information:
  mesh_dimension()=1
  spatial_dimension()=3
  n_nodes()=41
  n_elem()=20
   n_local_elem()=20
   n_active_elem()=20
  n_subdomains()=1
  n_processors()=1
  processor_id()=0

 EquationSystems
  n_systems()=1
   System "Poisson"
    Type "LinearImplicit"
    Variables="u" 
    Finite Element Types="LAGRANGE" 
    Approximation Orders="SECOND" 
    n_dofs()=41
    n_local_dofs()=41
    n_constrained_dofs()=0
    n_vectors()=1


-------------------------------------------------------
| Time:           Wed Jun  6 12:18:07 2007             |
| OS:             Linux                                |
| HostName:       orville                              |
| OS Release:     2.6.21-1.3194.fc7PAE                 |
| OS Version:     #1 SMP Wed May 23 22:27:31 EDT 2007  |
| Machine:        i686                                 |
| Username:       peterson                             |
-------------------------------------------------------
 ------------------------------------------------------------------------------
| Matrix Assembly Performance: Alive time=0.001533, Active time=0.000777       |
 ------------------------------------------------------------------------------
| Event                         nCalls    Total       Avg         Percent of   |
|                                         Time        Time        Active Time  |
|------------------------------------------------------------------------------|
|                                                                              |
| BCs                           20        0.0001      0.000005    13.26        |
| Fe                            20        0.0002      0.000010    26.25        |
| Ke                            20        0.0000      0.000002    5.53         |
| elem init                     20        0.0002      0.000011    28.44        |
| matrix insertion              20        0.0002      0.000010    26.51        |
 ------------------------------------------------------------------------------
| Totals:                       100       0.0008                  100.00       |
 ------------------------------------------------------------------------------

Running ./ex4-devel -d 2 -n 15

 Mesh Information:
  mesh_dimension()=2
  spatial_dimension()=3
  n_nodes()=961
  n_elem()=225
   n_local_elem()=225
   n_active_elem()=225
  n_subdomains()=1
  n_processors()=1
  processor_id()=0

 EquationSystems
  n_systems()=1
   System "Poisson"
    Type "LinearImplicit"
    Variables="u" 
    Finite Element Types="LAGRANGE" 
    Approximation Orders="SECOND" 
    n_dofs()=961
    n_local_dofs()=961
    n_constrained_dofs()=0
    n_vectors()=1


-------------------------------------------------------
| Time:           Wed Jun  6 12:18:07 2007             |
| OS:             Linux                                |
| HostName:       orville                              |
| OS Release:     2.6.21-1.3194.fc7PAE                 |
| OS Version:     #1 SMP Wed May 23 22:27:31 EDT 2007  |
| Machine:        i686                                 |
| Username:       peterson                             |
-------------------------------------------------------
 ------------------------------------------------------------------------------
| Matrix Assembly Performance: Alive time=0.02087, Active time=0.018112        |
 ------------------------------------------------------------------------------
| Event                         nCalls    Total       Avg         Percent of   |
|                                         Time        Time        Active Time  |
|------------------------------------------------------------------------------|
|                                                                              |
| BCs                           225       0.0042      0.000019    23.18        |
| Fe                            225       0.0046      0.000021    25.52        |
| Ke                            225       0.0038      0.000017    21.18        |
| elem init                     225       0.0037      0.000016    20.30        |
| matrix insertion              225       0.0018      0.000008    9.82         |
 ------------------------------------------------------------------------------
| Totals:                       1125      0.0181                  100.00       |
 ------------------------------------------------------------------------------

Running ./ex4-devel -d 3 -n 6

 Mesh Information:
  mesh_dimension()=3
  spatial_dimension()=3
  n_nodes()=2197
  n_elem()=216
   n_local_elem()=216
   n_active_elem()=216
  n_subdomains()=1
  n_processors()=1
  processor_id()=0

 EquationSystems
  n_systems()=1
   System "Poisson"
    Type "LinearImplicit"
    Variables="u" 
    Finite Element Types="LAGRANGE" 
    Approximation Orders="SECOND" 
    n_dofs()=2197
    n_local_dofs()=2197
    n_constrained_dofs()=0
    n_vectors()=1


-------------------------------------------------------
| Time:           Wed Jun  6 12:18:07 2007             |
| OS:             Linux                                |
| HostName:       orville                              |
| OS Release:     2.6.21-1.3194.fc7PAE                 |
| OS Version:     #1 SMP Wed May 23 22:27:31 EDT 2007  |
| Machine:        i686                                 |
| Username:       peterson                             |
-------------------------------------------------------
 ------------------------------------------------------------------------------
| Matrix Assembly Performance: Alive time=0.240115, Active time=0.237185       |
 ------------------------------------------------------------------------------
| Event                         nCalls    Total       Avg         Percent of   |
|                                         Time        Time        Active Time  |
|------------------------------------------------------------------------------|
|                                                                              |
| BCs                           216       0.1129      0.000523    47.59        |
| Fe                            216       0.0128      0.000059    5.41         |
| Ke                            216       0.0812      0.000376    34.22        |
| elem init                     216       0.0169      0.000078    7.12         |
| matrix insertion              216       0.0134      0.000062    5.66         |
 ------------------------------------------------------------------------------
| Totals:                       1080      0.2372                  100.00       |
 ------------------------------------------------------------------------------

 
***************************************************************
* Done Running Example  ./ex4-devel
***************************************************************

Site Created By: libMesh Developers
Last modified: October 22 2008 00:23:47.

Hosted By:
SourceForge.net Logo