1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217
|
// MFEM Example 14
//
// Compile with: make ex14
//
// Sample runs: ex14 -m ../data/inline-quad.mesh -o 0
// ex14 -m ../data/star.mesh -r 4 -o 2
// ex14 -m ../data/star-mixed.mesh -r 4 -o 2
// ex14 -m ../data/star-mixed.mesh -r 2 -o 2 -k 0 -e 1
// ex14 -m ../data/escher.mesh -s 1
// ex14 -m ../data/fichera.mesh -s 1 -k 1
// ex14 -m ../data/fichera-mixed.mesh -s 1 -k 1
// ex14 -m ../data/square-disc-p2.vtk -r 3 -o 2
// ex14 -m ../data/square-disc-p3.mesh -r 2 -o 3
// ex14 -m ../data/square-disc-nurbs.mesh -o 1
// ex14 -m ../data/disc-nurbs.mesh -r 3 -o 2 -s 1 -k 0
// ex14 -m ../data/pipe-nurbs.mesh -o 1
// ex14 -m ../data/inline-segment.mesh -r 5
// ex14 -m ../data/amr-quad.mesh -r 3
// ex14 -m ../data/amr-hex.mesh
// ex14 -m ../data/fichera-amr.mesh
// ex14 -pa -r 1 -o 3
// ex14 -pa -r 1 -o 3 -m ../data/fichera.mesh
//
// Device sample runs:
// ex14 -pa -r 2 -d cuda -o 3
// ex14 -pa -r 2 -d cuda -o 3 -m ../data/fichera.mesh
//
// Description: This example code demonstrates the use of MFEM to define a
// discontinuous Galerkin (DG) finite element discretization of
// the Laplace problem -Delta u = 1 with homogeneous Dirichlet
// boundary conditions. Finite element spaces of any order,
// including zero on regular grids, are supported. The example
// highlights the use of discontinuous spaces and DG-specific face
// integrators.
//
// We recommend viewing examples 1 and 9 before viewing this
// example.
#include "mfem.hpp"
#include <fstream>
#include <iostream>
using namespace std;
using namespace mfem;
int main(int argc, char *argv[])
{
// 1. Parse command-line options.
const char *mesh_file = "../data/star.mesh";
int ref_levels = -1;
int order = 1;
real_t sigma = -1.0;
real_t kappa = -1.0;
real_t eta = 0.0;
bool pa = false;
bool visualization = 1;
const char *device_config = "cpu";
OptionsParser args(argc, argv);
args.AddOption(&mesh_file, "-m", "--mesh",
"Mesh file to use.");
args.AddOption(&ref_levels, "-r", "--refine",
"Number of times to refine the mesh uniformly, -1 for auto.");
args.AddOption(&order, "-o", "--order",
"Finite element order (polynomial degree) >= 0.");
args.AddOption(&sigma, "-s", "--sigma",
"One of the three DG penalty parameters, typically +1/-1."
" See the documentation of class DGDiffusionIntegrator.");
args.AddOption(&kappa, "-k", "--kappa",
"One of the three DG penalty parameters, should be positive."
" Negative values are replaced with (order+1)^2.");
args.AddOption(&eta, "-e", "--eta", "BR2 penalty parameter.");
args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
"--no-partial-assembly", "Enable Partial Assembly.");
args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
"--no-visualization",
"Enable or disable GLVis visualization.");
args.AddOption(&device_config, "-d", "--device",
"Device configuration string, see Device::Configure().");
args.Parse();
if (!args.Good())
{
args.PrintUsage(cout);
return 1;
}
if (kappa < 0)
{
kappa = (order+1)*(order+1);
}
args.PrintOptions(cout);
// 2. Enable hardware devices such as GPUs, and programming models such as
// CUDA, OCCA, RAJA and OpenMP based on command line options.
Device device(device_config);
device.Print();
// 3. Read the mesh from the given mesh file. We can handle triangular,
// quadrilateral, tetrahedral and hexahedral meshes with the same code.
// NURBS meshes are projected to second order meshes.
Mesh mesh(mesh_file);
const int dim = mesh.Dimension();
// 4. Refine the mesh to increase the resolution. In this example we do
// 'ref_levels' of uniform refinement. By default, or if ref_levels < 0,
// we choose it to be the largest number that gives a final mesh with no
// more than 50,000 elements.
{
if (ref_levels < 0)
{
ref_levels = (int)floor(log(50000./mesh.GetNE())/log(2.)/dim);
}
for (int l = 0; l < ref_levels; l++)
{
mesh.UniformRefinement();
}
}
if (mesh.NURBSext)
{
mesh.SetCurvature(max(order, 1));
}
// 5. Define a finite element space on the mesh. Here we use discontinuous
// finite elements of the specified order >= 0.
const auto bt = pa ? BasisType::GaussLobatto : BasisType::GaussLegendre;
DG_FECollection fec(order, dim, bt);
FiniteElementSpace fespace(&mesh, &fec);
cout << "Number of unknowns: " << fespace.GetVSize() << endl;
// 6. Set up the linear form b(.) which corresponds to the right-hand side of
// the FEM linear system.
LinearForm b(&fespace);
ConstantCoefficient one(1.0);
ConstantCoefficient zero(0.0);
b.AddDomainIntegrator(new DomainLFIntegrator(one));
b.AddBdrFaceIntegrator(
new DGDirichletLFIntegrator(zero, one, sigma, kappa));
b.Assemble();
// 7. Define the solution vector x as a finite element grid function
// corresponding to fespace. Initialize x with initial guess of zero.
GridFunction x(&fespace);
x = 0.0;
// 8. Set up the bilinear form a(.,.) on the finite element space
// corresponding to the Laplacian operator -Delta, by adding the Diffusion
// domain integrator and the interior and boundary DG face integrators.
// Note that boundary conditions are imposed weakly in the form, so there
// is no need for dof elimination. After assembly and finalizing we
// extract the corresponding sparse matrix A.
BilinearForm a(&fespace);
a.AddDomainIntegrator(new DiffusionIntegrator(one));
a.AddInteriorFaceIntegrator(new DGDiffusionIntegrator(one, sigma, kappa));
a.AddBdrFaceIntegrator(new DGDiffusionIntegrator(one, sigma, kappa));
if (eta > 0)
{
MFEM_VERIFY(!pa, "BR2 not yet compatible with partial assembly.");
a.AddInteriorFaceIntegrator(new DGDiffusionBR2Integrator(fespace, eta));
a.AddBdrFaceIntegrator(new DGDiffusionBR2Integrator(fespace, eta));
}
if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
a.Assemble();
a.Finalize();
// 9. Define a simple symmetric Gauss-Seidel preconditioner and use it to
// solve the system Ax=b with PCG in the symmetric case, and GMRES in the
// non-symmetric one. (Note that tolerances are squared: 1e-12 corresponds
// to a relative tolerance of 1e-6).
//
// If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
if (pa)
{
MFEM_VERIFY(sigma == -1.0,
"The case of PA with sigma != -1 is not yet supported.");
CG(a, b, x, 1, 500, 1e-12, 0.0);
}
else
{
const SparseMatrix &A = a.SpMat();
#ifndef MFEM_USE_SUITESPARSE
GSSmoother M(A);
if (sigma == -1.0)
{
PCG(A, M, b, x, 1, 500, 1e-12, 0.0);
}
else
{
GMRES(A, M, b, x, 1, 500, 10, 1e-12, 0.0);
}
#else
UMFPackSolver umf_solver;
umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
umf_solver.SetOperator(A);
umf_solver.Mult(b, x);
#endif
}
// 10. 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);
// 11. 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;
}
return 0;
}
|