{ "cells": [ { "cell_type": "markdown", "id": "owned-extraction", "metadata": {}, "source": [ "## Load the MFEM library\n", "\n", "Any non-default libraries must be loaded before you can `#include` files that use them. For more info see the [xeus-cling help](https://xeus-cling.readthedocs.io/en/latest/build_options.html)." ] }, { "cell_type": "code", "execution_count": null, "id": "waiting-portrait", "metadata": {}, "outputs": [], "source": [ "#pragma cling load(\"mfem\")" ] }, { "cell_type": "markdown", "id": "foreign-recycling", "metadata": {}, "source": [ "## MFEM Example 1" ] }, { "cell_type": "markdown", "id": "public-white", "metadata": {}, "source": [ "This is the simplest MFEM example and a good starting point for new users. The example demonstrates the use of MFEM to define and solve an $H^1$ finite element discretization of the Laplace problem\n", "\n", "$$\n", "-\\Delta u = 1\n", "$$\n", "\n", "with homogeneous Dirichlet boundary conditions $u=0$.\n", "\n", "The example illustrates the use of the basic MFEM classes for defining the mesh, finite element space, as well as linear and bilinear forms corresponding to the left-hand side and right-hand side of the discrete linear system.\n", "\n", "Compare with MFEM's [ex1.cpp](https://github.com/mfem/mfem/blob/master/examples/ex1.cpp) and PyMFEM's [ex1.py](https://github.com/mfem/PyMFEM/blob/master/examples/ex1.py)." ] }, { "cell_type": "code", "execution_count": null, "id": "protective-darkness", "metadata": {}, "outputs": [], "source": [ "#include \n", "#include \n", "#include \n", "\n", "#include \n", "#include " ] }, { "cell_type": "code", "execution_count": null, "id": "falling-monkey", "metadata": {}, "outputs": [], "source": [ "using namespace std;\n", "using namespace mfem;\n", "\n", "Mesh mesh = Mesh::MakeCartesian2D(5, 5, Element::TRIANGLE);\n", "mesh.UniformRefinement();\n", "\n", "H1_FECollection fec(2, mesh.Dimension());\n", "\n", "FiniteElementSpace fespace(&mesh, &fec);\n", "cout << \"Number of finite element unknowns: \" << fespace.GetTrueVSize() << endl;\n", "\n", "Array ess_tdof_list;\n", "if (mesh.bdr_attributes.Size())\n", "{\n", " Array ess_bdr(mesh.bdr_attributes.Max());\n", " ess_bdr = 1;\n", " fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);\n", "}\n", "\n", "LinearForm b(&fespace);\n", "ConstantCoefficient one(1.0);\n", "b.AddDomainIntegrator(new DomainLFIntegrator(one));\n", "b.Assemble();\n", "\n", "GridFunction x(&fespace);\n", "x = 0.0;\n", "\n", "BilinearForm a(&fespace);\n", "a.AddDomainIntegrator(new DiffusionIntegrator(one));\n", "a.Assemble();\n", "\n", "OperatorPtr A;\n", "Vector B, X;\n", "a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);\n", "\n", "cout << \"Size of linear system: \" << A->Height() << endl;\n", "\n", "GSSmoother M((SparseMatrix&)(*A));\n", "PCG(*A, M, B, X, 1, 200, 1e-12, 0.0);\n", "a.RecoverFEMSolution(X, b, x);" ] }, { "cell_type": "markdown", "id": "hawaiian-republican", "metadata": {}, "source": [ "## GLVis Visualization\n", "\n", "For now we save the computational mesh and finite element solution in a string and pass that to the glvis widget, see https://github.com/glvis/xeus-glvis for the widget backend and https://github.com/GLVis/pyglvis/tree/master/js for the widget frontend." ] }, { "cell_type": "code", "execution_count": null, "id": "ordinary-equation", "metadata": {}, "outputs": [], "source": [ "std::stringstream ss;\n", "ss << \"solution\\n\" << mesh << x << flush;\n", "\n", "auto glv = glvis::glvis();\n", "glv.plot(ss.str() + \"keys Rjml\"); // the `+ \"keys ....\"' is optional\n", "glv" ] } ], "metadata": { "kernelspec": { "display_name": "C++14", "language": "C++14", "name": "xcpp14" }, "language_info": { "codemirror_mode": "text/x-c++src", "file_extension": ".cpp", "mimetype": "text/x-c++src", "name": "c++", "version": "14" } }, "nbformat": 4, "nbformat_minor": 5 }