MFEM v2.0
ex1.cpp
Go to the documentation of this file.
00001 //                                MFEM Example 1
00002 //
00003 // Compile with: make ex1
00004 //
00005 // Sample runs:  ex1 ../data/square-disc.mesh
00006 //               ex1 ../data/star.mesh
00007 //               ex1 ../data/escher.mesh
00008 //               ex1 ../data/fichera.mesh
00009 //               ex1 ../data/square-disc-p2.vtk
00010 //               ex1 ../data/square-disc-p3.mesh
00011 //               ex1 ../data/square-disc-nurbs.mesh
00012 //               ex1 ../data/disc-nurbs.mesh
00013 //               ex1 ../data/pipe-nurbs.mesh
00014 //
00015 // Description:  This example code demonstrates the use of MFEM to define a
00016 //               simple isoparametric finite element discretization of the
00017 //               Laplace problem -Delta u = 1 with homogeneous Dirichlet
00018 //               boundary conditions. Specifically, we discretize with the
00019 //               FE space coming from the mesh (linear by default, quadratic
00020 //               for quadratic curvilinear mesh, NURBS for NURBS mesh, etc.)
00021 //
00022 //               The example highlights the use of mesh refinement, finite
00023 //               element grid functions, as well as linear and bilinear forms
00024 //               corresponding to the left-hand side and right-hand side of the
00025 //               discrete linear system. We also cover the explicit elimination
00026 //               of boundary conditions on all boundary edges, and the optional
00027 //               connection to the GLVis tool for visualization.
00028 
00029 #include <fstream>
00030 #include "mfem.hpp"
00031 
00032 int main (int argc, char *argv[])
00033 {
00034    Mesh *mesh;
00035 
00036    if (argc == 1)
00037    {
00038       cout << "\nUsage: ex1 <mesh_file>\n" << endl;
00039       return 1;
00040    }
00041 
00042    // 1. Read the mesh from the given mesh file. We can handle triangular,
00043    //    quadrilateral, tetrahedral or hexahedral elements with the same code.
00044    ifstream imesh(argv[1]);
00045    if (!imesh)
00046    {
00047       cerr << "\nCan not open mesh file: " << argv[1] << '\n' << endl;
00048       return 2;
00049    }
00050    mesh = new Mesh(imesh, 1, 1);
00051    imesh.close();
00052 
00053    // 2. Refine the mesh to increase the resolution. In this example we do
00054    //    'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
00055    //    largest number that gives a final mesh with no more than 50,000
00056    //    elements.
00057    {
00058       int ref_levels =
00059          (int)floor(log(50000./mesh->GetNE())/log(2.)/mesh->Dimension());
00060       for (int l = 0; l < ref_levels; l++)
00061          mesh->UniformRefinement();
00062    }
00063 
00064    // 3. Define a finite element space on the mesh. Here we use isoparametric
00065    //    finite elements coming from the mesh nodes (linear by default).
00066    FiniteElementCollection *fec;
00067    if (mesh->GetNodes())
00068       fec = mesh->GetNodes()->OwnFEC();
00069    else
00070       fec = new LinearFECollection;
00071    FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
00072    cout << "Number of unknowns: " << fespace->GetVSize() << endl;
00073 
00074    // 4. Set up the linear form b(.) which corresponds to the right-hand side of
00075    //    the FEM linear system, which in this case is (1,phi_i) where phi_i are
00076    //    the basis functions in the finite element fespace.
00077    LinearForm *b = new LinearForm(fespace);
00078    ConstantCoefficient one(1.0);
00079    b->AddDomainIntegrator(new DomainLFIntegrator(one));
00080    b->Assemble();
00081 
00082    // 5. Define the solution vector x as a finite element grid function
00083    //    corresponding to fespace. Initialize x with initial guess of zero,
00084    //    which satisfies the boundary conditions.
00085    GridFunction x(fespace);
00086    x = 0.0;
00087 
00088    // 6. Set up the bilinear form a(.,.) on the finite element space
00089    //    corresponding to the Laplacian operator -Delta, by adding the Diffusion
00090    //    domain integrator and imposing homogeneous Dirichlet boundary
00091    //    conditions. The boundary conditions are implemented by marking all the
00092    //    boundary attributes from the mesh as essential (Dirichlet). After
00093    //    assembly and finalizing we extract the corresponding sparse matrix A.
00094    BilinearForm *a = new BilinearForm(fespace);
00095    a->AddDomainIntegrator(new DiffusionIntegrator(one));
00096    a->Assemble();
00097    Array<int> ess_bdr(mesh->bdr_attributes.Max());
00098    ess_bdr = 1;
00099    a->EliminateEssentialBC(ess_bdr, x, *b);
00100    a->Finalize();
00101    const SparseMatrix &A = a->SpMat();
00102 
00103    // 7. Define a simple symmetric Gauss-Seidel preconditioner and use it to
00104    //    solve the system Ax=b with PCG.
00105    GSSmoother M(A);
00106    PCG(A, M, *b, x, 1, 200, 1e-12, 0.0);
00107 
00108    // 8. Save the refined mesh and the solution. This output can be viewed later
00109    //    using GLVis: "glvis -m refined.mesh -g sol.gf".
00110    {
00111       ofstream mesh_ofs("refined.mesh");
00112       mesh_ofs.precision(8);
00113       mesh->Print(mesh_ofs);
00114       ofstream sol_ofs("sol.gf");
00115       sol_ofs.precision(8);
00116       x.Save(sol_ofs);
00117    }
00118 
00119    // 9. (Optional) Send the solution by socket to a GLVis server.
00120    char vishost[] = "localhost";
00121    int  visport   = 19916;
00122    osockstream sol_sock(visport, vishost);
00123    sol_sock << "solution\n";
00124    sol_sock.precision(8);
00125    mesh->Print(sol_sock);
00126    x.Save(sol_sock);
00127    sol_sock.send();
00128 
00129    // 10. Free the used memory.
00130    delete a;
00131    delete b;
00132    delete fespace;
00133    if (!mesh->GetNodes())
00134       delete fec;
00135    delete mesh;
00136 
00137    return 0;
00138 }
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines