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 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282
|
// MFEM Example 13 - Parallel Version
//
// Compile with: make ex13p
//
// Sample runs: mpirun -np 4 ex13p -m ../data/star.mesh
// mpirun -np 4 ex13p -m ../data/square-disc.mesh -o 2 -n 4
// mpirun -np 4 ex13p -m ../data/beam-tet.mesh
// mpirun -np 4 ex13p -m ../data/beam-tet.mesh -nc -o 2 -rs 1
// mpirun -np 4 ex13p -m ../data/beam-hex.mesh
// mpirun -np 4 ex13p -m ../data/escher.mesh
// mpirun -np 4 ex13p -m ../data/fichera.mesh
// mpirun -np 4 ex13p -m ../data/fichera-q2.vtk
// mpirun -np 4 ex13p -m ../data/fichera-q3.mesh
// mpirun -np 4 ex13p -m ../data/square-disc-nurbs.mesh
// mpirun -np 4 ex13p -m ../data/beam-hex-nurbs.mesh
// mpirun -np 4 ex13p -m ../data/amr-quad.mesh -o 2
// mpirun -np 4 ex13p -m ../data/amr-hex.mesh
// mpirun -np 4 ex13p -m ../data/mobius-strip.mesh -n 8 -o 2
// mpirun -np 4 ex13p -m ../data/klein-bottle.mesh -n 10 -o 2
//
// Description: This example code solves the Maxwell (electromagnetic)
// eigenvalue problem curl curl E = lambda E with homogeneous
// Dirichlet boundary conditions E x n = 0.
//
// We compute a number of the lowest nonzero eigenmodes by
// discretizing the curl curl operator using a Nedelec FE space of
// the specified order in 2D or 3D.
//
// The example highlights the use of the AME subspace eigenvalue
// solver from HYPRE, which uses LOBPCG and AMS internally.
// Reusing a single GLVis visualization window for multiple
// eigenfunctions is also illustrated.
//
// We recommend viewing examples 3 and 11 before viewing this
// example.
#include "mfem.hpp"
#include <fstream>
#include <iostream>
using namespace std;
using namespace mfem;
int main(int argc, char *argv[])
{
// 1. Initialize MPI and HYPRE.
Mpi::Init(argc, argv);
int num_procs = Mpi::WorldSize();
int myid = Mpi::WorldRank();
Hypre::Init();
// 2. Parse command-line options.
const char *mesh_file = "../data/beam-tet.mesh";
int ser_ref_levels = 2;
int par_ref_levels = 1;
int order = 1;
int nev = 5;
bool nc = 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(&ser_ref_levels, "-rs", "--refine-serial",
"Number of times to refine the mesh uniformly in serial.");
args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
"Number of times to refine the mesh uniformly in parallel.");
args.AddOption(&order, "-o", "--order",
"Finite element order (polynomial degree) or -1 for"
" isoparametric space.");
args.AddOption(&nev, "-n", "--num-eigs",
"Number of desired eigenmodes.");
args.AddOption(&nc, "-nc", "--non-conforming", "-c",
"--conforming",
"Mark the mesh as nonconforming before partitioning.");
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())
{
if (myid == 0)
{
args.PrintUsage(cout);
}
return 1;
}
if (myid == 0)
{
args.PrintOptions(cout);
}
// 3. 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);
if (myid == 0) { device.Print(); }
// 4. Read the (serial) mesh from the given mesh file on all processors. We
// can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
// and volume meshes with the same code.
Mesh *mesh = new Mesh(mesh_file, 1, 1);
int dim = mesh->Dimension();
if (nc)
{
mesh->EnsureNCMesh(true);
}
// 5. Refine the serial mesh on all processors to increase the resolution. In
// this example we do 'ref_levels' of uniform refinement (2 by default, or
// specified on the command line with -rs).
for (int lev = 0; lev < ser_ref_levels; lev++)
{
mesh->UniformRefinement();
}
// 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
// this mesh further in parallel to increase the resolution (1 time by
// default, or specified on the command line with -rp). Once the parallel
// mesh is defined, the serial mesh can be deleted.
ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
delete mesh;
for (int lev = 0; lev < par_ref_levels; lev++)
{
pmesh->UniformRefinement();
}
// 7. Define a parallel finite element space on the parallel mesh. Here we
// use the Nedelec finite elements of the specified order.
FiniteElementCollection *fec = new ND_FECollection(order, dim);
ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
HYPRE_BigInt size = fespace->GlobalTrueVSize();
if (myid == 0)
{
cout << "Number of unknowns: " << size << endl;
}
// 8. Set up the parallel bilinear forms a(.,.) and m(.,.) on the finite
// element space. The first corresponds to the curl curl, while the second
// is a simple mass matrix needed on the right hand side of the
// generalized eigenvalue problem below. The boundary conditions are
// implemented by marking all the boundary attributes from the mesh as
// essential. The corresponding degrees of freedom are eliminated with
// special values on the diagonal to shift the Dirichlet eigenvalues out
// of the computational range. After serial and parallel assembly we
// extract the corresponding parallel matrices A and M.
ConstantCoefficient one(1.0);
Array<int> ess_bdr;
if (pmesh->bdr_attributes.Size())
{
ess_bdr.SetSize(pmesh->bdr_attributes.Max());
ess_bdr = 1;
}
ParBilinearForm *a = new ParBilinearForm(fespace);
a->AddDomainIntegrator(new CurlCurlIntegrator(one));
if (pmesh->bdr_attributes.Size() == 0)
{
// Add a mass term if the mesh has no boundary, e.g. periodic mesh or
// closed surface.
a->AddDomainIntegrator(new VectorFEMassIntegrator(one));
}
a->Assemble();
a->EliminateEssentialBCDiag(ess_bdr, 1.0);
a->Finalize();
ParBilinearForm *m = new ParBilinearForm(fespace);
m->AddDomainIntegrator(new VectorFEMassIntegrator(one));
m->Assemble();
// shift the eigenvalue corresponding to eliminated dofs to a large value
m->EliminateEssentialBCDiag(ess_bdr, numeric_limits<real_t>::min());
m->Finalize();
HypreParMatrix *A = a->ParallelAssemble();
HypreParMatrix *M = m->ParallelAssemble();
delete a;
delete m;
// 9. Define and configure the AME eigensolver and the AMS preconditioner for
// A to be used within the solver. Set the matrices which define the
// generalized eigenproblem A x = lambda M x.
HypreAMS *ams = new HypreAMS(*A,fespace);
ams->SetPrintLevel(0);
ams->SetSingularProblem();
HypreAME *ame = new HypreAME(MPI_COMM_WORLD);
ame->SetNumModes(nev);
ame->SetPreconditioner(*ams);
ame->SetMaxIter(100);
ame->SetTol(1e-8);
ame->SetPrintLevel(1);
ame->SetMassMatrix(*M);
ame->SetOperator(*A);
// 10. Compute the eigenmodes and extract the array of eigenvalues. Define a
// parallel grid function to represent each of the eigenmodes returned by
// the solver.
Array<real_t> eigenvalues;
ame->Solve();
ame->GetEigenvalues(eigenvalues);
ParGridFunction x(fespace);
// 11. Save the refined mesh and the modes in parallel. This output can be
// viewed later using GLVis: "glvis -np <np> -m mesh -g mode".
{
ostringstream mesh_name, mode_name;
mesh_name << "mesh." << setfill('0') << setw(6) << myid;
ofstream mesh_ofs(mesh_name.str().c_str());
mesh_ofs.precision(8);
pmesh->Print(mesh_ofs);
for (int i=0; i<nev; i++)
{
// convert eigenvector from HypreParVector to ParGridFunction
x = ame->GetEigenvector(i);
mode_name << "mode_" << setfill('0') << setw(2) << i << "."
<< setfill('0') << setw(6) << myid;
ofstream mode_ofs(mode_name.str().c_str());
mode_ofs.precision(8);
x.Save(mode_ofs);
mode_name.str("");
}
}
// 12. Send the solution by socket to a GLVis server.
if (visualization)
{
char vishost[] = "localhost";
int visport = 19916;
socketstream mode_sock(vishost, visport);
mode_sock.precision(8);
for (int i=0; i<nev; i++)
{
if ( myid == 0 )
{
cout << "Eigenmode " << i+1 << '/' << nev
<< ", Lambda = " << eigenvalues[i] << endl;
}
// convert eigenvector from HypreParVector to ParGridFunction
x = ame->GetEigenvector(i);
mode_sock << "parallel " << num_procs << " " << myid << "\n"
<< "solution\n" << *pmesh << x << flush
<< "window_title 'Eigenmode " << i+1 << '/' << nev
<< ", Lambda = " << eigenvalues[i] << "'" << endl;
char c;
if (myid == 0)
{
cout << "press (q)uit or (c)ontinue --> " << flush;
cin >> c;
}
MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
if (c != 'c')
{
break;
}
}
mode_sock.close();
}
// 13. Free the used memory.
delete ame;
delete ams;
delete M;
delete A;
delete fespace;
delete fec;
delete pmesh;
return 0;
}
|