Permalink
Fetching contributors…
Cannot retrieve contributors at this time
266 lines (235 sloc) 9.21 KB
// MFEM Example 4
//
// Compile with: make ex4
//
// Sample runs: ex4 -m ../data/square-disc.mesh
// ex4 -m ../data/star.mesh
// ex4 -m ../data/beam-tet.mesh
// ex4 -m ../data/beam-hex.mesh
// ex4 -m ../data/escher.mesh
// ex4 -m ../data/fichera.mesh -o 2 -hb
// ex4 -m ../data/fichera-q2.vtk
// ex4 -m ../data/fichera-q3.mesh -o 2 -sc
// ex4 -m ../data/square-disc-nurbs.mesh
// ex4 -m ../data/beam-hex-nurbs.mesh
// ex4 -m ../data/periodic-square.mesh -no-bc
// ex4 -m ../data/periodic-cube.mesh -no-bc
// ex4 -m ../data/amr-quad.mesh
// ex4 -m ../data/amr-hex.mesh
// ex4 -m ../data/amr-hex.mesh -o 2 -hb
// ex4 -m ../data/fichera-amr.mesh -o 2 -sc
// ex4 -m ../data/star-surf.mesh -o 1
//
// Description: This example code solves a simple 2D/3D H(div) diffusion
// problem corresponding to the second order definite equation
// -grad(alpha div F) + beta F = f with boundary condition F dot n
// = <given normal field>. Here, we use a given exact solution F
// and compute the corresponding r.h.s. f. We discretize with
// Raviart-Thomas finite elements.
//
// The example demonstrates the use of H(div) finite element
// spaces with the grad-div and H(div) vector finite element mass
// bilinear form, as well as the computation of discretization
// error when the exact solution is known. Bilinear form
// hybridization and static condensation are also illustrated.
//
// We recommend viewing examples 1-3 before viewing this example.
#include "mfem.hpp"
#include <fstream>
#include <iostream>
using namespace std;
using namespace mfem;
// Exact solution, F, and r.h.s., f. See below for implementation.
void F_exact(const Vector &, Vector &);
void f_exact(const Vector &, Vector &);
double freq = 1.0, kappa;
int main(int argc, char *argv[])
{
// 1. Parse command-line options.
const char *mesh_file = "../data/star.mesh";
int order = 1;
bool set_bc = true;
bool static_cond = false;
bool hybridization = false;
bool visualization = 1;
OptionsParser args(argc, argv);
args.AddOption(&mesh_file, "-m", "--mesh",
"Mesh file to use.");
args.AddOption(&order, "-o", "--order",
"Finite element order (polynomial degree).");
args.AddOption(&set_bc, "-bc", "--impose-bc", "-no-bc", "--dont-impose-bc",
"Impose or not essential boundary conditions.");
args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
" solution.");
args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
"--no-static-condensation", "Enable static condensation.");
args.AddOption(&hybridization, "-hb", "--hybridization", "-no-hb",
"--no-hybridization", "Enable hybridization.");
args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
"--no-visualization",
"Enable or disable GLVis visualization.");
args.Parse();
if (!args.Good())
{
args.PrintUsage(cout);
return 1;
}
args.PrintOptions(cout);
kappa = freq * M_PI;
// 2. Read the mesh from the given mesh file. We can handle triangular,
// quadrilateral, tetrahedral, hexahedral, surface and volume, as well as
// periodic meshes with the same code.
Mesh *mesh = new Mesh(mesh_file, 1, 1);
int dim = mesh->Dimension();
int sdim = mesh->SpaceDimension();
// 3. Refine the mesh to increase the resolution. In this example we do
// 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
// largest number that gives a final mesh with no more than 25,000
// elements.
{
int ref_levels =
(int)floor(log(25000./mesh->GetNE())/log(2.)/dim);
for (int l = 0; l < ref_levels; l++)
{
mesh->UniformRefinement();
}
}
// 4. Define a finite element space on the mesh. Here we use the
// Raviart-Thomas finite elements of the specified order.
FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
cout << "Number of finite element unknowns: "
<< fespace->GetTrueVSize() << endl;
// 5. Determine the list of true (i.e. conforming) essential boundary dofs.
// In this example, the boundary conditions are defined by marking all
// the boundary attributes from the mesh as essential (Dirichlet) and
// converting them to a list of true dofs.
Array<int> ess_tdof_list;
if (mesh->bdr_attributes.Size())
{
Array<int> ess_bdr(mesh->bdr_attributes.Max());
ess_bdr = set_bc ? 1 : 0;
fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
}
// 6. Set up the linear form b(.) which corresponds to the right-hand side
// of the FEM linear system, which in this case is (f,phi_i) where f is
// given by the function f_exact and phi_i are the basis functions in the
// finite element fespace.
VectorFunctionCoefficient f(sdim, f_exact);
LinearForm *b = new LinearForm(fespace);
b->AddDomainIntegrator(new VectorFEDomainLFIntegrator(f));
b->Assemble();
// 7. Define the solution vector x as a finite element grid function
// corresponding to fespace. Initialize x by projecting the exact
// solution. Note that only values from the boundary faces will be used
// when eliminating the non-homogeneous boundary condition to modify the
// r.h.s. vector b.
GridFunction x(fespace);
VectorFunctionCoefficient F(sdim, F_exact);
x.ProjectCoefficient(F);
// 8. Set up the bilinear form corresponding to the H(div) diffusion operator
// grad alpha div + beta I, by adding the div-div and the mass domain
// integrators.
Coefficient *alpha = new ConstantCoefficient(1.0);
Coefficient *beta = new ConstantCoefficient(1.0);
BilinearForm *a = new BilinearForm(fespace);
a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
a->AddDomainIntegrator(new VectorFEMassIntegrator(*beta));
// 9. Assemble the bilinear form and the corresponding linear system,
// applying any necessary transformations such as: eliminating boundary
// conditions, applying conforming constraints for non-conforming AMR,
// static condensation, hybridization, etc.
FiniteElementCollection *hfec = NULL;
FiniteElementSpace *hfes = NULL;
if (static_cond)
{
a->EnableStaticCondensation();
}
else if (hybridization)
{
hfec = new DG_Interface_FECollection(order-1, dim);
hfes = new FiniteElementSpace(mesh, hfec);
a->EnableHybridization(hfes, new NormalTraceJumpIntegrator(),
ess_tdof_list);
}
a->Assemble();
SparseMatrix A;
Vector B, X;
a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
cout << "Size of linear system: " << A.Height() << endl;
#ifndef MFEM_USE_SUITESPARSE
// 10. Define a simple symmetric Gauss-Seidel preconditioner and use it to
// solve the system A X = B with PCG.
GSSmoother M(A);
PCG(A, M, B, X, 1, 10000, 1e-20, 0.0);
#else
// 10. If compiled with SuiteSparse support, use UMFPACK to solve the system.
UMFPackSolver umf_solver;
umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
umf_solver.SetOperator(A);
umf_solver.Mult(B, X);
#endif
// 11. Recover the solution as a finite element grid function.
a->RecoverFEMSolution(X, *b, x);
// 12. Compute and print the L^2 norm of the error.
cout << "\n|| F_h - F ||_{L^2} = " << x.ComputeL2Error(F) << '\n' << endl;
// 13. Save the refined mesh and the solution. This output can be viewed
// later using GLVis: "glvis -m refined.mesh -g sol.gf".
{
ofstream mesh_ofs("refined.mesh");
mesh_ofs.precision(8);
mesh->Print(mesh_ofs);
ofstream sol_ofs("sol.gf");
sol_ofs.precision(8);
x.Save(sol_ofs);
}
// 14. Send the solution by socket to a GLVis server.
if (visualization)
{
char vishost[] = "localhost";
int visport = 19916;
socketstream sol_sock(vishost, visport);
sol_sock.precision(8);
sol_sock << "solution\n" << *mesh << x << flush;
}
// 15. Free the used memory.
delete hfes;
delete hfec;
delete a;
delete alpha;
delete beta;
delete b;
delete fespace;
delete fec;
delete mesh;
return 0;
}
// The exact solution (for non-surface meshes)
void F_exact(const Vector &p, Vector &F)
{
int dim = p.Size();
double x = p(0);
double y = p(1);
// double z = (dim == 3) ? p(2) : 0.0;
F(0) = cos(kappa*x)*sin(kappa*y);
F(1) = cos(kappa*y)*sin(kappa*x);
if (dim == 3)
{
F(2) = 0.0;
}
}
// The right hand side
void f_exact(const Vector &p, Vector &f)
{
int dim = p.Size();
double x = p(0);
double y = p(1);
// double z = (dim == 3) ? p(2) : 0.0;
double temp = 1 + 2*kappa*kappa;
f(0) = temp*cos(kappa*x)*sin(kappa*y);
f(1) = temp*cos(kappa*y)*sin(kappa*x);
if (dim == 3)
{
f(2) = 0;
}
}