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 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310
|
// MFEM Example 26 - Parallel Version
//
// Compile with: make ex26p
//
// Sample runs: mpirun -np 4 ex26p -m ../data/star.mesh
// mpirun -np 4 ex26p -m ../data/fichera.mesh
// mpirun -np 4 ex26p -m ../data/beam-hex.mesh
//
// Device sample runs:
// mpirun -np 4 ex26p -d cuda
// mpirun -np 4 ex26p -d occa-cuda
// mpirun -np 4 ex26p -d raja-omp
// mpirun -np 4 ex26p -d ceed-cpu
// mpirun -np 4 ex26p -d ceed-cuda
//
// Description: This example code demonstrates the use of MFEM to define a
// simple finite element discretization of the Laplace problem
// -Delta u = 1 with homogeneous Dirichlet boundary conditions
// as in Example 1.
//
// It highlights on the creation of a hierarchy of discretization
// spaces with partial assembly and the construction of an
// efficient multigrid preconditioner for the iterative solver.
//
// We recommend viewing Example 1 before viewing this example.
#include "mfem.hpp"
#include <fstream>
#include <iostream>
using namespace std;
using namespace mfem;
// Class for constructing a multigrid preconditioner for the diffusion operator.
// This example multigrid preconditioner class demonstrates the creation of the
// parallel diffusion bilinear forms and operators using partial assembly for
// all spaces except the coarsest one in the ParFiniteElementSpaceHierarchy.
// The multigrid uses a PCG solver preconditioned with AMG on the coarsest level
// and second order Chebyshev accelerated smoothers on the other levels.
class DiffusionMultigrid : public GeometricMultigrid
{
private:
ConstantCoefficient coeff;
HypreBoomerAMG* amg;
public:
// Constructs a diffusion multigrid for the ParFiniteElementSpaceHierarchy
// and the array of essential boundaries
DiffusionMultigrid(ParFiniteElementSpaceHierarchy& fespaces,
Array<int>& ess_bdr)
: GeometricMultigrid(fespaces, ess_bdr), coeff(1.0)
{
ConstructCoarseOperatorAndSolver(fespaces.GetFESpaceAtLevel(0));
for (int level = 1; level < fespaces.GetNumLevels(); ++level)
{
ConstructOperatorAndSmoother(fespaces.GetFESpaceAtLevel(level), level);
}
}
virtual ~DiffusionMultigrid()
{
delete amg;
}
private:
void ConstructBilinearForm(ParFiniteElementSpace& fespace,
bool partial_assembly)
{
ParBilinearForm* form = new ParBilinearForm(&fespace);
if (partial_assembly)
{
form->SetAssemblyLevel(AssemblyLevel::PARTIAL);
}
form->AddDomainIntegrator(new DiffusionIntegrator(coeff));
form->Assemble();
bfs.Append(form);
}
void ConstructCoarseOperatorAndSolver(ParFiniteElementSpace& coarse_fespace)
{
ConstructBilinearForm(coarse_fespace, false);
HypreParMatrix* hypreCoarseMat = new HypreParMatrix();
bfs[0]->FormSystemMatrix(*essentialTrueDofs[0], *hypreCoarseMat);
amg = new HypreBoomerAMG(*hypreCoarseMat);
amg->SetPrintLevel(-1);
CGSolver* pcg = new CGSolver(MPI_COMM_WORLD);
pcg->SetPrintLevel(-1);
pcg->SetMaxIter(10);
pcg->SetRelTol(sqrt(1e-4));
pcg->SetAbsTol(0.0);
pcg->SetOperator(*hypreCoarseMat);
pcg->SetPreconditioner(*amg);
AddLevel(hypreCoarseMat, pcg, true, true);
}
void ConstructOperatorAndSmoother(ParFiniteElementSpace& fespace, int level)
{
const Array<int> &ess_tdof_list = *essentialTrueDofs[level];
ConstructBilinearForm(fespace, true);
OperatorPtr opr;
opr.SetType(Operator::ANY_TYPE);
bfs.Last()->FormSystemMatrix(ess_tdof_list, opr);
opr.SetOperatorOwner(false);
Vector diag(fespace.GetTrueVSize());
bfs.Last()->AssembleDiagonal(diag);
Solver* smoother = new OperatorChebyshevSmoother(
*opr, diag, ess_tdof_list, 2, fespace.GetParMesh()->GetComm());
AddLevel(opr.Ptr(), smoother, true, true);
}
};
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/star.mesh";
int geometric_refinements = 0;
int order_refinements = 2;
const char *device_config = "cpu";
bool visualization = true;
OptionsParser args(argc, argv);
args.AddOption(&mesh_file, "-m", "--mesh",
"Mesh file to use.");
args.AddOption(&geometric_refinements, "-gr", "--geometric-refinements",
"Number of geometric refinements done prior to order refinements.");
args.AddOption(&order_refinements, "-or", "--order-refinements",
"Number of order refinements. Finest level in the hierarchy has order 2^{or}.");
args.AddOption(&device_config, "-d", "--device",
"Device configuration string, see Device::Configure().");
args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
"--no-visualization",
"Enable or disable GLVis visualization.");
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();
// 5. Refine the serial mesh on all processors 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 1,000 elements.
{
int ref_levels =
(int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
for (int l = 0; l < ref_levels; l++)
{
mesh->UniformRefinement();
}
}
// 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
// this mesh further in parallel to increase the resolution. Once the
// parallel mesh is defined, the serial mesh can be deleted.
ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
delete mesh;
{
int par_ref_levels = 2;
for (int l = 0; l < par_ref_levels; l++)
{
pmesh->UniformRefinement();
}
}
// 7. Define a parallel finite element space hierarchy on the parallel mesh.
// Here we use continuous Lagrange finite elements. We start with order 1
// on the coarse level and geometrically refine the spaces by the specified
// amount. Afterwards, we increase the order of the finite elements by a
// factor of 2 for each additional level.
FiniteElementCollection *fec = new H1_FECollection(1, dim);
ParFiniteElementSpace *coarse_fespace = new ParFiniteElementSpace(pmesh, fec);
Array<FiniteElementCollection*> collections;
collections.Append(fec);
ParFiniteElementSpaceHierarchy* fespaces = new ParFiniteElementSpaceHierarchy(
pmesh, coarse_fespace, true, true);
for (int level = 0; level < geometric_refinements; ++level)
{
fespaces->AddUniformlyRefinedLevel();
}
for (int level = 0; level < order_refinements; ++level)
{
collections.Append(new H1_FECollection((int)std::pow(2, level+1), dim));
fespaces->AddOrderRefinedLevel(collections.Last());
}
HYPRE_BigInt size = fespaces->GetFinestFESpace().GlobalTrueVSize();
if (myid == 0)
{
cout << "Number of finite element unknowns: " << size << endl;
}
// 8. Set up the parallel linear form b(.) which corresponds to the
// right-hand side of the FEM linear system, which in this case is
// (1,phi_i) where phi_i are the basis functions in fespace.
ParLinearForm *b = new ParLinearForm(&fespaces->GetFinestFESpace());
ConstantCoefficient one(1.0);
b->AddDomainIntegrator(new DomainLFIntegrator(one));
b->Assemble();
// 9. Define the solution vector x as a parallel finite element grid function
// corresponding to fespace. Initialize x with initial guess of zero,
// which satisfies the boundary conditions.
ParGridFunction x(&fespaces->GetFinestFESpace());
x = 0.0;
// 10. Create the multigrid operator using the previously created parallel
// FiniteElementSpaceHierarchy and additional boundary information. This
// operator is then used to create the MultigridSolver as preconditioner
// in the iterative solver.
Array<int> ess_bdr(pmesh->bdr_attributes.Max());
if (pmesh->bdr_attributes.Size())
{
ess_bdr = 1;
}
DiffusionMultigrid* M = new DiffusionMultigrid(*fespaces, ess_bdr);
M->SetCycleType(Multigrid::CycleType::VCYCLE, 1, 1);
OperatorPtr A;
Vector X, B;
M->FormFineLinearSystem(x, *b, A, X, B);
// 11. Solve the linear system A X = B.
CGSolver cg(MPI_COMM_WORLD);
cg.SetRelTol(1e-12);
cg.SetMaxIter(2000);
cg.SetPrintLevel(1);
cg.SetOperator(*A);
cg.SetPreconditioner(*M);
cg.Mult(B, X);
// 12. Recover the parallel grid function corresponding to X. This is the
// local finite element solution on each processor.
M->RecoverFineFEMSolution(X, *b, x);
// 13. Save the refined mesh and the solution in parallel. This output can be
// viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
{
ostringstream mesh_name, sol_name;
mesh_name << "mesh." << setfill('0') << setw(6) << myid;
sol_name << "sol." << setfill('0') << setw(6) << myid;
ofstream mesh_ofs(mesh_name.str().c_str());
mesh_ofs.precision(8);
fespaces->GetFinestFESpace().GetParMesh()->Print(mesh_ofs);
ofstream sol_ofs(sol_name.str().c_str());
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 << "parallel " << num_procs << " " << myid << "\n";
sol_sock.precision(8);
sol_sock << "solution\n" << *fespaces->GetFinestFESpace().GetParMesh()
<< x << flush;
}
// 15. Free the used memory.
delete M;
delete b;
delete fespaces;
for (int level = 0; level < collections.Size(); ++level)
{
delete collections[level];
}
return 0;
}
|