This page was generated from jupyter-files/unit-1.1-poisson/poisson.ipynb.

1.1 First NGSolve example

Let us solve the Poisson problem.

\[\text{find: } u \in H_{0,D}^1 \quad \int_\Omega \nabla u \nabla v = \int_\Omega f v, \quad \text{ for all } v \in H_{0,D}^1.\]

Quick steps to solution:

1. Import NGSolve and Netgen Python modules:

In [1]:
import netgen.gui
%gui tk
from ngsolve import *
from netgen.geom2d import unit_square

2. Generate an unstructured mesh

In [2]:
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
mesh.nv, mesh.ne   # number of vertices & elements
Out[2]:
(39, 56)
  • Here we prescribed a maximal mesh-size of 0.2 using the maxh flag.
  • The mesh can be viewed by switching to the Mesh tab in the Netgen GUI.

3. Declare a finite element space:

In [3]:
fes = H1(mesh, order=2, dirichlet="bottom|right")
fes.ndof  # number of unknowns in this space
Out[3]:
133

Python’s help system displays further documentation.

In [4]:
help(fes)
Help on H1 in module ngsolve.comp object:

class H1(FESpace)
 |   Keyword arguments can be:
 |  order: int = 1
 |    order of finite element space
 |  complex: bool = False
 |  dirichlet: regexpr
 |    Regular expression string defining the dirichlet boundary.
 |    More than one boundary can be combined by the | operator,
 |    i.e.: dirichlet = 'top|right'
 |  definedon: Region or regexpr
 |    FESpace is only defined on specific Region, created with mesh.Materials('regexpr')
 |    or mesh.Boundaries('regexpr'). If given a regexpr, the region is assumed to be
 |    mesh.Materials('regexpr').
 |  dim: int = 1
 |    Create multi dimensional FESpace (i.e. [H1]^3)
 |  dgjumps: bool = False
 |    Enable discontinuous space for DG methods, this flag is needed for DG methods,
 |    since the dofs have a different coupling then and this changes the sparsity
 |    pattern of matrices.
 |
 |  Method resolution order:
 |      H1
 |      FESpace
 |      pybind11_builtins.pybind11_object
 |      builtins.object
 |
 |  Methods defined here:
 |
 |  __getstate__(...)
 |      __getstate__(self: ngsolve.comp.FESpace) -> tuple
 |
 |  __init__(...)
 |      __init__(self: ngsolve.comp.H1, mesh: ngsolve.comp.Mesh, **kwargs) -> None
 |
 |  __setstate__(...)
 |      __setstate__(self: ngsolve.comp.H1, arg0: tuple) -> None
 |
 |  ----------------------------------------------------------------------
 |  Data descriptors defined here:
 |
 |  __dict__
 |
 |  ----------------------------------------------------------------------
 |  Methods inherited from FESpace:
 |
 |  ApplyM(...)
 |      ApplyM(self: ngsolve.comp.FESpace, vec: ngsolve.la.BaseVector, rho: ngsolve.fem.CoefficientFunction=None) -> None
 |
 |      Apply mass-matrix. Available only for L2-like spaces
 |
 |  CouplingType(...)
 |      CouplingType(self: ngsolve.comp.FESpace, dofnr: int) -> ngsolve.comp.COUPLING_TYPE
 |
 |      get coupling type of a degree of freedom
 |
 |  Elements(...)
 |      Elements(self: ngsolve.comp.FESpace, VOL_or_BND: ngsolve.comp.VorB=VorB.VOL) -> ngsolve.comp.FESpaceElementRange
 |
 |  FinalizeUpdate(...)
 |      FinalizeUpdate(self: ngsolve.comp.FESpace) -> None
 |
 |      finalize update
 |
 |  FreeDofs(...)
 |      FreeDofs(self: ngsolve.comp.FESpace, coupling: bool=False) -> ngsolve.ngstd.BitArray
 |
 |      Return BitArray of free (non-Dirichlet) dofs
 |      coupling=False ... all free dofs including local dofs
 |      coupling=True .... only element-boundary free dofs
 |
 |  GetDofNrs(...)
 |      GetDofNrs(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. GetDofNrs(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.ElementId) -> tuple
 |
 |      2. GetDofNrs(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.NodeId) -> tuple
 |
 |  GetDofs(...)
 |      GetDofs(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.Region) -> ngsolve.ngstd.BitArray
 |
 |  GetFE(...)
 |      GetFE(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.ElementId) -> object
 |
 |  ParallelDofs(...)
 |      ParallelDofs(self: ngsolve.comp.FESpace) -> ngsolve.la.ParallelDofs
 |
 |      Return dof-identification for MPI-distributed meshes
 |
 |  Range(...)
 |      Range(self: ngsolve.comp.FESpace, component: int) -> slice
 |
 |      Return interval of dofs of a component of a product space
 |
 |  SetCouplingType(...)
 |      SetCouplingType(self: ngsolve.comp.FESpace, dofnr: int, coupling_type: ngsolve.comp.COUPLING_TYPE) -> None
 |
 |      set coupling type of a degree of freedom
 |
 |  SetDefinedOn(...)
 |      SetDefinedOn(self: ngsolve.comp.FESpace, Region: ngsolve.comp.Region) -> None
 |
 |  SetOrder(...)
 |      SetOrder(*args, **kwargs)
 |      Overloaded function.
 |
 |      1. SetOrder(self: ngsolve.comp.FESpace, element_type: ngsolve.fem.ET, order: int) -> None
 |
 |      2. SetOrder(self: ngsolve.comp.FESpace, nodeid: ngsolve.comp.NodeId, order: int) -> None
 |
 |  SolveM(...)
 |      SolveM(self: ngsolve.comp.FESpace, vec: ngsolve.la.BaseVector, rho: ngsolve.fem.CoefficientFunction=None) -> None
 |
 |      Solve with the mass-matrix. Available only for L2-like spaces
 |
 |  TestFunction(...)
 |      TestFunction(self: ngsolve.comp.FESpace) -> object
 |
 |      Return a proxy to be used as a testfunction for Symbolic Integrators
 |
 |  TnT(...)
 |      TnT(self: ngsolve.comp.FESpace) -> Tuple[object, object]
 |
 |      Return a tuple of trial and testfunction
 |
 |  TrialFunction(...)
 |      TrialFunction(self: ngsolve.comp.FESpace) -> object
 |
 |      Return a proxy to be used as a trialfunction in Symbolic Integrators
 |
 |  Update(...)
 |      Update(self: ngsolve.comp.FESpace) -> None
 |
 |      update space after mesh-refinement
 |
 |  UpdateDofTables(...)
 |      UpdateDofTables(self: ngsolve.comp.FESpace) -> None
 |
 |      update dof-tables after changing polynomial order distribution
 |
 |  __eq__(...)
 |      __eq__(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.FESpace) -> bool
 |
 |  __flags_doc__(...) from builtins.PyCapsule
 |      __flags_doc__() -> dict
 |
 |  __special_treated_flags__(...) from builtins.PyCapsule
 |      __special_treated_flags__() -> dict
 |
 |  __str__(...)
 |      __str__(self: ngsolve.comp.FESpace) -> str
 |
 |  __timing__(...)
 |      __timing__(self: ngsolve.comp.FESpace) -> object
 |
 |  ----------------------------------------------------------------------
 |  Data descriptors inherited from FESpace:
 |
 |  components
 |      Return a list of the components of a product space
 |
 |  globalorder
 |      query global order of space
 |
 |  is_complex
 |
 |  lospace
 |
 |  mesh
 |
 |  ndof
 |      number of degrees of freedom
 |
 |  ndofglobal
 |      global number of dofs on MPI-distributed mesh
 |
 |  type
 |      type of finite element space
 |
 |  ----------------------------------------------------------------------
 |  Methods inherited from pybind11_builtins.pybind11_object:
 |
 |  __new__(*args, **kwargs) from pybind11_builtins.pybind11_type
 |      Create and return a new object.  See help(type) for accurate signature.

4. Declare test function, trial function, and grid function

  • Test and trial function are symbolic objects - called ProxyFunctions - that help you construct bilinear forms (and have no space to hold solutions).
  • GridFunctions, on the other hand, represent functions in the finite element space and contains memory to hold coefficient vectors.
In [5]:
u = fes.TrialFunction()  # symbolic object
v = fes.TestFunction()   # symbolic object
gfu = GridFunction(fes)  # solution

A shorter command for obtaining at once both trial and test variables is now available:

In [6]:
u, v = fes.TnT()

5. Define and assemble linear and bilinear forms:

In [7]:
a = BilinearForm(fes, symmetric=True)
a += SymbolicBFI(grad(u)*grad(v))
a.Assemble()

f = LinearForm(fes)
f += SymbolicLFI(x*v)
f.Assemble()

You can examine the linear system in more detail:

In [8]:
print(f.vec)
 0.000333333
 0.0090505
 0.00633333
 0.000770891
 0.00419003
 0.00873851
 0.0110113
 0.0118831
 0.0154585
 0.0172657
 0.0153513
 0.0180248
 0.0203006
 0.00930857
 0.00722848
 0.00380471
 0.000687974
 0.000703175
 0.00144642
 0.00174641
 0.0130156
 0.0196218
 0.0235625
 0.0209991
 0.032979
 0.0275741
 0.0284488
 0.0170589
 0.0176338
 0.0124153
 0.0051304
 0.00284308
 0.00630915
 0.0220354
 0.0157804
 0.0210096
 0.0100578
 0.0181327
 0.0217548
 -6.66667e-05
 -3.33333e-05
 -0.000536217
 -0.000585765
 -0.00111116
 -0.0008
 -0.000766667
 -8.25863e-05
 -2.35728e-05
 -0.000125108
 -0.000290565
 -0.000209771
 -0.000451004
 -0.00040564
 -0.000661732
 -0.000800436
 -0.000504522
 -0.000864384
 -0.000967242
 -0.000947681
 -0.000970631
 -0.000782608
 -0.00113705
 -0.00131313
 -0.000661412
 -0.00146165
 -0.00135486
 -0.000642205
 -0.00125381
 -0.0012366
 -0.00157121
 -0.00143302
 -0.000438352
 -0.00133842
 -0.000991598
 -0.000377208
 -0.000778967
 -0.000731847
 -0.000237359
 -0.000666891
 -0.000536209
 -0.00039775
 -0.000251178
 -1.37451e-05
 -9.64747e-05
 -7.25997e-05
 -2.40468e-05
 -7.40193e-05
 -9.91414e-05
 -4.37045e-05
 -0.000208445
 -0.000157729
 -0.000237114
 -0.000737581
 -0.000370323
 -0.000591075
 -0.00090079
 -0.000832185
 -0.000759186
 -0.000954243
 -0.00110341
 -0.00101272
 -0.00107714
 -0.00121445
 -0.00106391
 -0.00100948
 -0.00115025
 -0.000971716
 -0.000950052
 -0.000983723
 -0.000968857
 -0.000711981
 -0.000799923
 -0.000661911
 -0.000816822
 -0.00081628
 -0.000348442
 -0.000484742
 -0.000677524
 -0.00016251
 -0.000294359
 -0.000169603
 -0.000232291
 -0.000459975
 -0.000328541
 -0.000828049
 -0.000903769
 -0.000868363
 -0.000543998
 -0.000763656
 -0.000844269
 -0.000848631
 -0.000633923
 -0.0008243


In [9]:
print(a.mat)
Row 0:   0: 1
Row 1:   1: 0.829127
Row 2:   2: 1
Row 3:   3: 0.880753
Row 4:   0: -0.5   4: 1.90755
Row 5:   4: -0.559177   5: 1.91948
Row 6:   5: -0.37742   6: 1.7799
Row 7:   1: -0.222678   6: -0.298266   7: 1.83885
Row 8:   1: -0.223343   8: 1.90243
Row 9:   8: -0.375377   9: 1.75223
Row 10:   9: -0.283322   10: 1.78161
Row 11:   2: -0.5   10: -0.268388   11: 1.99954
Row 12:   2: -0.5   11: -0.098129   12: 1.84322
Row 13:   12: -0.211616   13: 1.75354
Row 14:   13: -0.319115   14: 1.77273
Row 15:   3: -0.359363   14: -0.356559   15: 1.85161
Row 16:   3: -0.354152   16: 1.80192
Row 17:   16: -0.146234   17: 1.85831
Row 18:   17: -0.280268   18: 1.8809
Row 19:   0: -0.5   4: -0.159147   18: -0.497592   19: 1.89108
Row 20:   4: -0.689223   5: -0.134026   18: -0.213492   19: -0.734339   20: 3.5107
Row 21:   5: -0.848855   6: -0.338846   20: -0.607596   21: 3.57058
Row 22:   6: -0.76537   7: -0.325948   21: -0.590849   22: 3.53156
Row 23:   1: -0.383106   7: -0.991955   8: -1.08169   22: -0.711434   23: 3.78819
Row 24:   8: -0.222018   9: -0.654805   22: -0.588928   23: -0.62   24: 3.44463
Row 25:   9: -0.438728   10: -0.852694   24: -0.516043   25: 3.57605
Row 26:   10: -0.377206   11: -1.13302   12: -0.22963   25: -0.541968   26: 3.73524
Row 27:   12: -0.803846   13: -0.666276   26: -0.773272   27: 3.68753
Row 28:   13: -0.556538   14: -0.353124   27: -0.781766   28: 3.60218
Row 29:   14: -0.74393   15: -0.202714   28: -0.943134   29: 3.72717
Row 30:   3: -0.167237   15: -0.932971   16: -0.826096   29: -0.840918   30: 3.7686
Row 31:   16: -0.475439   17: -1.04297   30: -0.848015   31: 3.7525
Row 32:   17: -0.388835   18: -0.889547   20: -0.459182   31: -0.532046   32: 3.57532
Row 33:   21: -0.763143   22: -0.54903   24: -0.350551   33: 3.66745
Row 34:   20: -0.672842   21: -0.421292   32: -0.676988   33: -0.704385   34: 3.51434
Row 35:   25: -0.347729   26: -0.680143   27: -0.662371   28: -0.447405   35: 3.54016
Row 36:   29: -0.830595   30: -0.153364   31: -0.854028   32: -0.628726   34: -0.40922   36: 3.64742
Row 37:   24: -0.492286   25: -0.878889   33: -1.07273   35: -0.76226   37: 3.75396
Row 38:   28: -0.520216   29: -0.165884   33: -0.227614   34: -0.629613   35: -0.640254   36: -0.771484   37: -0.547794   38: 3.50286
Row 39:   0: -0.0833333   4: 0   19: 0.0833333   39: 0.0416667
Row 40:   0: -0.0833333   4: 0.0833333   19: 0   39: 1.73472e-18   40: 0.0416667
Row 41:   1: -0.0318638   7: -0.0835338   23: 0.115398   41: 0.0381276
Row 42:   1: -0.0319873   8: -0.0831334   23: 0.115121   42: 0.0380861
Row 43:   1: -0.0743368   7: 0.120647   8: 0.120357   23: -0.166667   41: -0.0208834   42: -0.0207833   43: 0.0762138
Row 44:   2: -0.0833333   11: 0   12: 0.0833333   44: 0.0416667
Row 45:   2: -0.0833333   11: 0.0833333   12: 0   44: 0   45: 0.0416667
Row 46:   3: -0.0151987   15: -0.0803559   30: 0.0955546   46: 0.0388621
Row 47:   3: -0.0126741   16: -0.0864211   30: 0.0990952   47: 0.0395302
Row 48:   3: -0.118919   15: 0.14025   16: 0.145447   30: -0.166777   46: -0.020089   47: -0.0216053   48: 0.0783923
Row 49:   4: -0.0583984   5: -0.00990762   20: 0.068306   49: 0.0403755
Row 50:   0: 0.166667   4: -0.139806   19: -0.148957   20: 0.122096   39: -0.0208333   40: -0.0208333   50: 0.0788218
Row 51:   4: -0.119721   5: 0.103104   19: 0.0921482   20: -0.0755312   49: -0.0024769   50: -0.0164059   51: 0.0775306
Row 52:   5: -0.0550429   6: -0.0295225   21: 0.0845653   52: 0.0368672
Row 53:   4: 0.151595   5: -0.179629   20: -0.117774   21: 0.145809   49: -0.0145996   51: -0.023299   53: 0.0799353
Row 54:   5: -0.0753334   6: 0.0924258   20: 0.0718059   21: -0.0888984   52: -0.00738062   53: -0.014844   54: 0.0764269
Row 55:   6: -0.0582582   7: -0.0374957   22: 0.0957539   55: 0.0363662
Row 56:   5: 0.117946   6: -0.132207   21: -0.107784   22: 0.122044   52: -0.0137607   54: -0.0157258   56: 0.0741162
Row 57:   6: -0.0766627   7: 0.0872066   21: 0.0796928   22: -0.0902367   55: -0.00937392   56: -0.0131853   57: 0.0736152
Row 58:   6: 0.107969   7: -0.131503   22: -0.114716   23: 0.13825   55: -0.0145646   57: -0.0124277   58: 0.075136
Row 59:   1: 0.0689767   7: -0.053942   22: 0.0732872   23: -0.0883219   41: -0.00796594   43: -0.00927824   58: -0.0141145   59: 0.0768974
Row 60:   8: -0.0282254   9: -0.0570403   24: 0.0852656   60: 0.0369571
Row 61:   1: 0.0692111   8: -0.0460015   23: -0.0894962   24: 0.0662866   42: -0.00799682   43: -0.00930597   61: 0.0789449
Row 62:   8: -0.159711   9: 0.119603   23: 0.154658   24: -0.114549   60: -0.0142601   61: -0.0143772   62: 0.0778159
Row 63:   9: -0.0300697   10: -0.071478   25: 0.101548   63: 0.037192
Row 64:   8: 0.0907882   9: -0.105614   24: -0.0776415   25: 0.0924677   60: -0.00705634   62: -0.0156407   64: 0.0730975
Row 65:   9: -0.0993143   10: 0.118698   24: 0.10151   25: -0.120894   63: -0.0178695   64: -0.012354   65: 0.0733324
Row 66:   10: -0.0304682   11: -0.0742233   26: 0.104691   66: 0.0373557
Row 67:   9: 0.0772901   10: -0.0796198   25: -0.0752556   26: 0.0775854   63: -0.00751743   65: -0.0118051   67: 0.0742478
Row 68:   10: -0.115369   11: 0.118955   25: 0.115823   26: -0.119409   66: -0.0185558   67: -0.0112965   68: 0.0744115
Row 69:   2: 0.166667   11: -0.197946   12: -0.122045   26: 0.153325   44: -0.0208333   45: -0.0208333   69: 0.0840865
Row 70:   10: 0.0751996   11: -0.0610862   12: 0.0550663   26: -0.0691797   66: -0.00761705   68: -0.0111829   69: -0.00967788   70: 0.0797756
Row 71:   12: -0.0509965   13: -0.0596508   27: 0.110647   71: 0.0364792
Row 72:   11: 0.130968   12: -0.0993327   26: -0.199191   27: 0.167556   69: -0.0286533   70: -0.00408871   72: 0.084199
Row 73:   12: -0.0348296   13: 0.0949202   26: 0.0841386   27: -0.144229   71: -0.0149127   72: -0.0211446   73: 0.0782583
Row 74:   13: -0.0597076   14: -0.0333842   28: 0.0930917   74: 0.0365694
Row 75:   12: 0.0862659   13: -0.0683181   27: -0.113119   28: 0.0951715   71: -0.0127491   73: -0.00881734   75: 0.0731209
Row 76:   13: -0.104581   14: 0.08657   27: 0.113518   28: -0.095507   74: -0.00834605   75: -0.0155307   76: 0.0732111
Row 77:   14: -0.0676137   15: -0.0230352   29: 0.0906489   77: 0.0375189
Row 78:   13: 0.112893   14: -0.10956   28: -0.127013   29: 0.12368   74: -0.0149269   76: -0.0132964   78: 0.0738569
Row 79:   14: -0.0848962   15: 0.0824617   28: 0.0927755   29: -0.090341   77: -0.00575881   78: -0.0168264   79: 0.0748064
Row 80:   14: 0.12704   15: -0.134566   29: -0.139062   30: 0.146588   77: -0.0169034   79: -0.0148566   80: 0.0768533
Row 81:   3: 0.0750926   15: -0.0706443   29: 0.0821986   30: -0.086647   46: -0.00379968   48: -0.0149735   80: -0.0178621   81: 0.0781966
Row 82:   16: -0.0441816   17: -0.0855914   31: 0.129773   82: 0.0385364
Row 83:   3: 0.0716995   16: -0.0940836   30: -0.0723049   31: 0.0946889   47: -0.00316853   48: -0.0147563   83: 0.0760178
Row 84:   16: -0.0756339   17: 0.109964   30: 0.110892   31: -0.145222   82: -0.0213979   83: -0.0149077   84: 0.075024
Row 85:   17: -0.042337   18: -0.0557767   32: 0.0981137   85: 0.0362063
Row 86:   16: 0.068554   17: -0.0468412   31: -0.0890017   32: 0.0672889   82: -0.0110454   84: -0.00609309   86: 0.0774179
Row 87:   17: -0.134949   18: 0.102488   31: 0.133057   32: -0.100597   85: -0.0139442   86: -0.011205   87: 0.0750878
Row 88:   18: -0.016011   19: -0.0567662   20: 0.0727772   88: 0.0389273
Row 89:   18: -0.175413   19: 0.139698   20: -0.102588   32: 0.138303   88: -0.0141915   89: 0.0783959
Row 90:   17: 0.0890484   18: -0.0662824   20: 0.0653934   32: -0.0881593   85: -0.0105842   87: -0.0116778   89: -0.0114556   90: 0.0756749
Row 91:   4: 0.0829967   18: 0.0989431   19: -0.109457   20: -0.0724832   50: -0.014118   51: -0.00663113   88: -0.00400276   89: -0.020733   91: 0.0760824
Row 92:   5: 0.098863   20: -0.0714262   21: -0.130771   34: 0.103334   53: -0.0216082   54: -0.0031075   92: 0.0758658
Row 93:   18: 0.112052   20: -0.0727152   32: -0.155837   34: 0.1165   89: -0.0231203   90: -0.00489276   93: 0.0762705
Row 94:   20: -0.072598   21: 0.086228   32: 0.0940636   34: -0.107694   92: -0.0110845   93: -0.0158389   94: 0.073108
Row 95:   6: 0.0962552   21: -0.0854376   22: -0.110271   33: 0.0994538   56: -0.0173258   57: -0.00673795   95: 0.0735459
Row 96:   21: -0.0716112   22: 0.0867017   33: -0.095593   34: 0.0805024   95: -0.010242   96: 0.0735987
Row 97:   20: 0.100886   21: -0.110595   33: 0.12333   34: -0.113621   92: -0.014749   94: -0.0104725   96: -0.0136563   97: 0.0736078
Row 98:   7: 0.0986211   22: -0.0547959   23: -0.127616   24: 0.0837913   58: -0.020448   59: -0.00420727   98: 0.0752462
Row 99:   22: -0.112651   23: 0.107939   24: -0.0810715   33: 0.0857841   98: -0.0114561   99: 0.0729694
Row 100:   21: 0.10422   22: -0.105922   24: 0.0954349   33: -0.093733   95: -0.0146215   96: -0.0114335   99: -0.00881178   100: 0.0727899
Row 101:   8: 0.105926   22: 0.100081   23: -0.159263   24: -0.0467445   61: -0.00219441   62: -0.0242872   98: -0.00949172   99: -0.0155286   101: 0.0773352
Row 102:   9: 0.0951455   24: -0.0879729   25: -0.120761   37: 0.113588   64: -0.0107629   65: -0.0130235   102: 0.0736852
Row 103:   22: 0.110725   24: -0.106357   33: -0.135247   37: 0.130879   99: -0.0126343   100: -0.015047   103: 0.0750071
Row 104:   24: -0.059769   25: 0.1143   33: 0.107888   37: -0.162419   102: -0.0194273   103: -0.0211775   104: 0.076059
Row 105:   10: 0.103037   25: -0.106824   26: -0.097702   35: 0.101489   67: -0.00809987   68: -0.0176594   105: 0.0737134
Row 106:   25: -0.113914   26: 0.110445   35: -0.125467   37: 0.128937   105: -0.0163256   106: 0.074334
Row 107:   24: 0.07247   25: -0.0583599   35: 0.0819333   37: -0.0960434   102: -0.00896975   104: -0.00914776   106: -0.0150411   107: 0.0752212
Row 108:   12: 0.0825381   26: -0.0476148   27: -0.13512   35: 0.100197   72: -0.0207445   73: 0.000109951   108: 0.0779035
Row 109:   25: 0.0813281   26: -0.0894424   27: 0.0964427   35: -0.0883284   105: -0.0090465   106: -0.0112855   108: -0.0130356   109: 0.072782
Row 110:   13: 0.084444   27: -0.109648   28: -0.056567   35: 0.0817709   75: -0.00826218   76: -0.0128488   110: 0.0741273
Row 111:   26: 0.0923549   27: -0.112472   28: 0.0916897   35: -0.0715729   108: -0.0120137   109: -0.0110751   110: -0.00587956   111: 0.0736099
Row 112:   14: 0.0818443   28: -0.0773065   29: -0.0724993   38: 0.0679614   78: -0.0140936   79: -0.00636744   112: 0.0767487
Row 113:   27: 0.126424   28: -0.103037   35: -0.118365   38: 0.0949784   110: -0.0145632   111: -0.0170429   113: 0.0739925
Row 114:   28: -0.140932   29: 0.106008   35: 0.111162   38: -0.0762371   112: -0.00403117   113: -0.0150281   114: 0.0759681
Row 115:   15: 0.0858897   29: -0.0709417   30: -0.0969322   36: 0.0819842   80: -0.0187848   81: -0.00268759   115: 0.0770067
Row 116:   29: -0.0802273   30: 0.0904976   36: -0.0891129   38: 0.0788427   115: -0.00544822   116: 0.0769432
Row 117:   28: 0.14172   29: -0.168124   36: 0.145561   38: -0.119157   112: -0.0129592   114: -0.0224708   116: -0.01683   117: 0.0787321
Row 118:   16: 0.0863197   30: -0.0550293   31: -0.112704   36: 0.0814134   83: -0.00876454   84: -0.0128154   118: 0.0772672
Row 119:   29: 0.128896   30: -0.15041   31: 0.159351   36: -0.137837   115: -0.0150478   116: -0.0171762   118: -0.0194114   119: 0.0784519
Row 120:   17: 0.110706   31: -0.15293   32: -0.0603088   36: 0.102532   86: -0.0056172   87: -0.0220593   120: 0.0754782
Row 121:   30: 0.085473   31: -0.125559   32: 0.0816942   36: -0.0416078   118: -0.000941958   119: -0.0204263   120: -0.00945999   121: 0.0773763
Row 122:   20: 0.083852   32: -0.0976556   34: -0.084342   36: 0.0981455   93: -0.013286   94: -0.00767698   122: 0.0737072
Row 123:   31: 0.108547   32: -0.09333   34: 0.0806736   36: -0.0958903   120: -0.0161731   121: -0.0109636   122: -0.00779947   123: 0.073502
Row 124:   21: 0.094582   33: -0.0965577   34: -0.0832128   38: 0.0851884   96: -0.00646935   97: -0.0171762   124: 0.074292
Row 125:   24: 0.0693467   33: -0.0332607   37: -0.103733   38: 0.0676465   103: -0.0115422   104: -0.0057945   125: 0.0789452
Row 126:   33: -0.15685   34: 0.120108   37: 0.151642   38: -0.114899   124: -0.0143338   125: -0.014391   126: 0.0774213
Row 127:   32: 0.116423   34: -0.0970758   36: -0.128208   38: 0.108861   122: -0.0167369   123: -0.0123689   127: 0.0733718
Row 128:   33: 0.0906254   34: -0.0997778   36: 0.098266   38: -0.0891136   124: -0.00696326   126: -0.0156931   127: -0.0153151   128: 0.0734567
Row 129:   25: 0.0905411   35: -0.0683654   37: -0.102507   38: 0.0803318   106: -0.0171931   107: -0.00544221   129: 0.074479
Row 130:   28: 0.0859151   35: -0.117928   37: 0.100614   38: -0.0686011   113: -0.00871649   114: -0.0127623   129: -0.0084338   130: 0.0733096
Row 131:   29: 0.0897638   34: 0.0846056   36: -0.115247   38: -0.0591227   116: -0.00288065   117: -0.0195603   127: -0.0119   128: -0.00925136   131: 0.0757375
Row 132:   33: 0.104161   35: 0.113475   37: -0.160957   38: -0.0566793   125: -0.00252067   126: -0.0235195   129: -0.0116491   130: -0.0167197   132: 0.0772338

6. Solve the system:

In [10]:
gfu.vec.data = \
    a.mat.Inverse(freedofs=fes.FreeDofs()) * f.vec
Draw(gfu)

The Dirichlet boundary condition constrains some degrees of freedom. The argument fes.FreeDofs() indicates that only the remaining “free” degrees of freedom should participate in the linear solve.

You can examine the coefficient vector of solution if needed:

In [11]:
print(gfu.vec)
       0
       0
       0
 0.0923019
       0
       0
       0
       0
       0
       0
       0
       0
 0.0578968
 0.086337
 0.0954017
 0.0944824
 0.0888267
 0.078059
 0.0595638
 0.0331316
 0.0428238
 0.0391783
 0.0329405
 0.0190228
 0.0396205
 0.0437688
 0.0470417
 0.0751353
 0.0891838
 0.0930866
 0.0914133
 0.0837669
 0.0706131
 0.0597334
 0.0673558
 0.0730976
 0.0844105
 0.0647765
 0.0814921
       0
 -0.00563511
       0
       0
 0.0213178
       0
 -0.0350904
 0.00259849
 -0.00728873
 -0.00280466
       0
 -0.00967023
 -0.0123906
       0
 -0.0240987
 -0.0135412
       0
 -0.0171313
 -0.0200732
 -0.00676928
 -0.0212666
       0
 -0.0254979
 -0.00995785
       0
 -0.0390824
 -0.013713
       0
 -0.0269058
 -0.0176991
 -0.0333808
 -0.0238911
 -0.0243047
 -0.00442431
 -0.0102753
 -0.0146405
 -0.006598
 -0.00679179
 -0.00570818
 -0.00931926
 -0.00609179
 -0.00771503
 -0.00532019
 -0.00755367
 0.00262541
 -0.00112314
 -0.00818072
 0.00105745
 0.000896631
 -0.00820643
 -0.00077952
 0.00155068
 0.00436689
 -0.00490164
 -0.01261
 -0.00903107
 -0.00751952
 -0.0136094
 -0.013109
 -0.00783553
 -0.0252117
 -0.00819149
 -0.00330874
 -0.00919931
 -0.0212406
 -0.00570119
 -0.00323099
 -0.0135482
 -0.0185787
 -0.0142031
 -0.0209225
 -0.013664
 -0.00405249
 -0.00922214
 -0.0103465
 -0.00890607
 -0.00277651
 -0.00675828
 -0.0105432
 -0.00452866
 -0.00654881
 -0.00618555
 -0.000547291
 -0.00262347
 -0.00559157
 -0.00883386
 -0.0116287
 -0.010581
 -0.00921447
 -0.0126657
 -0.00446274
 -0.0171414
 -0.00785433
 -0.0118158


Ways to interact with NGSolve

  • A jupyter notebook (like this one) gives you one way to interact with NGSolve. When you have a complex sequence of tasks to perform, the notebook may not be adequate.

  • You can write an entire python script in a text editor and call python on the command line. (A script of the above is provided in poisson.py.)

    python3 poisson.py
    
  • If you want the Netgen GUI, then use netgen on the command line:

    netgen poisson.py
    

    You can then ask for a python shell from the GUI’s menu options (Solve -> Python shell).