This page was generated from unit-1.1-poisson/poisson.ipynb.
1.1 First NGSolve example¶
Let us solve the Poisson problem of finding \(u\) satisfying
Quick steps to solution:¶
1. Import NGSolve and Netgen Python modules:¶
[1]:
import netgen.gui
from ngsolve import *
from netgen.geom2d import unit_square
2. Generate an unstructured mesh¶
[2]:
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
mesh.nv, mesh.ne # number of vertices & elements
[2]:
(38, 54)
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:¶
[3]:
fes = H1(mesh, order=2, dirichlet="bottom|right")
fes.ndof # number of unknowns in this space
[3]:
129
Python’s help system displays further documentation.
[4]:
help(fes)
Help on H1 in module ngsolve.comp object:
class H1(FESpace)
| An H1-conforming finite element space.
|
| The H1 finite element space consists of continuous and
| elemenet-wise polynomial functions. It uses a hierarchical (=modal)
| basis built from integrated Legendre polynomials on tensor-product elements,
| and Jaboci polynomials on simplicial elements.
|
| Boundary values are well defined. The function can be used directly on the
| boundary, using the trace operator is optional.
|
| The H1 space supports variable order, which can be set individually for edges,
| faces and cells.
|
| Internal degrees of freedom are declared as local dofs and are eliminated
| if static condensation is on.
|
| The wirebasket consists of all vertex dofs. Optionally, one can include the
| first (the quadratic bubble) edge basis function, or all edge basis functions
| into the wirebasket.
|
| Keyword arguments can be:
| order: int = 1
| order of finite element space
| complex: bool = False
| Set if FESpace should be complex
| 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.
| low_order_space: bool = True
| Generate a lowest order space together with the high-order space,
| needed for some preconditioners.
| order_policy: ORDER_POLICY = ORDER_POLICY.OLDSTYLE
| CONSTANT .. use the same fixed order for all elements,
| NODAL ..... use the same order for nodes of same shape,
| VARIBLE ... use an individual order for each edge, face and cell,
| OLDSTYLE .. as it used to be for the last decade
| wb_withedges: bool = true(3D) / false(2D)
| use lowest-order edge dofs for BDDC wirebasket
| wb_fulledges: bool = false
| use all edge dofs for BDDC wirebasket
|
| Method resolution order:
| H1
| FESpace
| NGS_Object
| 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, autoupdate: bool = False, **kwargs) -> None
|
| __setstate__(...)
| __setstate__(self: ngsolve.comp.H1, arg0: tuple) -> None
|
| ----------------------------------------------------------------------
| Static methods defined here:
|
| __flags_doc__(...) from builtins.PyCapsule
| __flags_doc__() -> dict
|
| ----------------------------------------------------------------------
| Data descriptors defined here:
|
| __dict__
|
| ----------------------------------------------------------------------
| Methods inherited from FESpace:
|
| ApplyM(...)
| ApplyM(self: ngsolve.comp.FESpace, vec: ngsolve.la.BaseVector, rho: ngsolve.fem.CoefficientFunction = None, definedon: ngsolve.comp.Region = None) -> None
|
| Apply mass-matrix. Available only for L2-like spaces
|
| ConvertL2Operator(...)
| ConvertL2Operator(self: ngsolve.comp.FESpace, l2space: ngsolve.comp.FESpace) -> ngsolve.la.BaseMatrix
|
| CouplingType(...)
| CouplingType(self: ngsolve.comp.FESpace, dofnr: int) -> ngsolve.comp.COUPLING_TYPE
|
|
| Get coupling type of a degree of freedom.
|
| Parameters:
|
| dofnr : int
| input dof number
|
| Elements(...)
| Elements(self: ngsolve.comp.FESpace, VOL_or_BND: ngsolve.comp.VorB = VorB.VOL) -> ngsolve.comp.FESpaceElementRange
|
|
| Returns an iterable range of elements.
|
| Parameters:
|
| VOL_or_BND : ngsolve.comp.VorB
| input VOL, BND, BBND,...
|
| FinalizeUpdate(...)
| FinalizeUpdate(self: ngsolve.comp.FESpace) -> None
|
| finalize update
|
| FreeDofs(...)
| FreeDofs(self: ngsolve.comp.FESpace, coupling: bool = False) -> pyngcore.BitArray
|
|
|
| Return BitArray of free (non-Dirichlet) dofs\n
| coupling=False ... all free dofs including local dofs\n
| coupling=True .... only element-boundary free dofs
|
| Parameters:
|
| coupling : bool
| input coupling
|
| GetDofNrs(...)
| GetDofNrs(*args, **kwargs)
| Overloaded function.
|
| 1. GetDofNrs(self: ngsolve.comp.FESpace, ei: ngsolve.comp.ElementId) -> tuple
|
|
|
| Parameters:
|
| ei : ngsolve.comp.ElementId
| input element id
|
|
|
| 2. GetDofNrs(self: ngsolve.comp.FESpace, ni: ngsolve.comp.NodeId) -> tuple
|
|
|
| Parameters:
|
| ni : ngsolve.comp.NodeId
| input node id
|
| GetDofs(...)
| GetDofs(self: ngsolve.comp.FESpace, region: ngsolve.comp.Region) -> pyngcore.BitArray
|
|
| Returns all degrees of freedom in given region.
|
| Parameters:
|
| region : ngsolve.comp.Region
| input region
|
| GetFE(...)
| GetFE(self: ngsolve.comp.FESpace, ei: ngsolve.comp.ElementId) -> object
|
|
| Get the finite element to corresponding element id.
|
| Parameters:
|
| ei : ngsolve.comp.ElementId
| input element id
|
| GetOrder(...)
| GetOrder(self: ngsolve.comp.FESpace, nodeid: ngsolve.comp.NodeId) -> int
|
| return order of node.
| by now, only isotropic order is supported here
|
| GetTrace(...)
| GetTrace(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.FESpace, arg1: ngsolve.la.BaseVector, arg2: ngsolve.la.BaseVector, arg3: bool) -> None
|
| GetTraceTrans(...)
| GetTraceTrans(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.FESpace, arg1: ngsolve.la.BaseVector, arg2: ngsolve.la.BaseVector, arg3: bool) -> None
|
| HideAllDofs(...)
| HideAllDofs(self: ngsolve.comp.FESpace, component: object = <ngsolve.ngstd.DummyArgument>) -> None
|
| set all visible coupling types to HIDDEN_DOFs (will be overwritten by any Update())
|
| InvM(...)
| InvM(self: ngsolve.comp.FESpace, rho: ngsolve.fem.CoefficientFunction = None) -> ngsolve.la.BaseMatrix
|
| Mass(...)
| Mass(self: ngsolve.comp.FESpace, rho: ngsolve.fem.CoefficientFunction = None, definedon: Optional[ngsolve.comp.Region] = None) -> ngsolve.la.BaseMatrix
|
| ParallelDofs(...)
| ParallelDofs(self: ngsolve.comp.FESpace) -> ngsolve.la.ParallelDofs
|
| Return dof-identification for MPI-distributed meshes
|
| Prolongation(...)
| Prolongation(self: ngsolve.comp.FESpace) -> ngmg::Prolongation
|
| Return prolongation operator for use in multi-grid
|
| Range(...)
| Range(self: ngsolve.comp.FESpace, component: int) -> ngsolve.ngstd.IntRange
|
|
| Return interval of dofs of a component of a product space.
|
| Parameters:
|
| component : int
| input component
|
| SetCouplingType(...)
| SetCouplingType(*args, **kwargs)
| Overloaded function.
|
| 1. SetCouplingType(self: ngsolve.comp.FESpace, dofnr: int, coupling_type: ngsolve.comp.COUPLING_TYPE) -> None
|
|
| Set coupling type of a degree of freedom.
|
| Parameters:
|
| dofnr : int
| input dof number
|
| coupling_type : ngsolve.comp.COUPLING_TYPE
| input coupling type
|
|
|
| 2. SetCouplingType(self: ngsolve.comp.FESpace, dofnrs: ngsolve.ngstd.IntRange, coupling_type: ngsolve.comp.COUPLING_TYPE) -> None
|
|
| Set coupling type for interval of dofs.
|
| Parameters:
|
| dofnrs : Range
| range of dofs
|
| coupling_type : ngsolve.comp.COUPLING_TYPE
| input coupling type
|
| SetDefinedOn(...)
| SetDefinedOn(self: ngsolve.comp.FESpace, region: ngsolve.comp.Region) -> None
|
|
| Set the regions on which the FESpace is defined.
|
| Parameters:
|
| region : ngsolve.comp.Region
| input region
|
| SetOrder(...)
| SetOrder(*args, **kwargs)
| Overloaded function.
|
| 1. SetOrder(self: ngsolve.comp.FESpace, element_type: ngsolve.fem.ET, order: int) -> None
|
|
|
| Parameters:
|
| element_type : ngsolve.fem.ET
| input element type
|
| order : object
| input polynomial order
|
|
| 2. SetOrder(self: ngsolve.comp.FESpace, nodeid: ngsolve.comp.NodeId, order: int) -> None
|
|
|
| Parameters:
|
| nodeid : ngsolve.comp.NodeId
| input node id
|
| order : int
| input polynomial order
|
| SolveM(...)
| SolveM(self: ngsolve.comp.FESpace, vec: ngsolve.la.BaseVector, rho: ngsolve.fem.CoefficientFunction = None, definedon: ngsolve.comp.Region = None) -> None
|
|
| Solve with the mass-matrix. Available only for L2-like spaces.
|
| Parameters:
|
| vec : ngsolve.la.BaseVector
| input right hand side vector
|
| rho : ngsolve.fem.CoefficientFunction
| input CF
|
| TestFunction(...)
| TestFunction(self: ngsolve.comp.FESpace) -> object
|
| Return a proxy to be used as a testfunction for :any:`Symbolic Integrators<symbolic-integrators>`
|
| TnT(...)
| TnT(self: ngsolve.comp.FESpace) -> Tuple[object, object]
|
| Return a tuple of trial and testfunction
|
| TraceOperator(...)
| TraceOperator(self: ngsolve.comp.FESpace, tracespace: ngsolve.comp.FESpace, average: bool) -> ngsolve.la.BaseMatrix
|
| TrialFunction(...)
| TrialFunction(self: ngsolve.comp.FESpace) -> object
|
| Return a proxy to be used as a trialfunction in :any:`Symbolic Integrators<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, space: ngsolve.comp.FESpace) -> bool
|
| __str__(...)
| __str__(self: ngsolve.comp.FESpace) -> str
|
| __timing__(...)
| __timing__(self: ngsolve.comp.FESpace) -> object
|
| ----------------------------------------------------------------------
| Static methods inherited from FESpace:
|
| __special_treated_flags__(...) from builtins.PyCapsule
| __special_treated_flags__() -> dict
|
| ----------------------------------------------------------------------
| Data descriptors inherited from FESpace:
|
| components
| Return a list of the components of a product space
|
| couplingtype
|
| dim
| multi-dim of FESpace
|
| globalorder
| query global order of space
|
| is_complex
|
| lospace
|
| mesh
| mesh on which the FESpace is created
|
| ndof
| number of degrees of freedom
|
| ndofglobal
| global number of dofs on MPI-distributed mesh
|
| type
| type of finite element space
|
| ----------------------------------------------------------------------
| Data descriptors inherited from NGS_Object:
|
| __memory__
|
| name
|
| ----------------------------------------------------------------------
| 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.
[5]:
u = fes.TrialFunction() # symbolic object
v = fes.TestFunction() # symbolic object
gfu = GridFunction(fes) # solution
Alternately, you can get both the trial and test variables at once:
[6]:
u, v = fes.TnT()
5. Define and assemble linear and bilinear forms:¶
[7]:
a = BilinearForm(fes, symmetric=True)
a += grad(u)*grad(v)*dx
a.Assemble()
f = LinearForm(fes)
f += x*v*dx
f.Assemble()
You can examine the linear system in more detail:
[8]:
print(f.vec)
0.000333333
0.00896813
0.00633333
0.000801205
0.00384743
0.00745274
0.010402
0.0116458
0.0153532
0.0173012
0.0155111
0.0182062
0.0206364
0.00967366
0.00782633
0.00429241
0.000385723
0.00214091
0.000471409
0.00255107
0.0114933
0.0179248
0.0228051
0.0207208
0.0328187
0.0278923
0.0289917
0.0175791
0.0185541
0.0135822
0.0054952
0.00674838
0.0216669
0.0154523
0.0216449
0.0112586
0.0184036
0.0228345
-6.66667e-05
-3.33333e-05
-0.000530412
-0.000581324
-0.00110125
-0.0008
-0.000766667
-9.53866e-05
-1.92861e-05
-0.000125689
-0.000236921
-0.000212375
-0.000419804
-0.000376692
-0.000547145
-0.00069263
-0.000491432
-0.000804894
-0.000920196
-0.000924592
-0.000954328
-0.000780634
-0.0011262
-0.00130465
-0.000667812
-0.00146005
-0.00136067
-0.000648932
-0.00126658
-0.00124969
-0.00158695
-0.00145486
-0.000449393
-0.0013689
-0.00101734
-0.000393375
-0.000808082
-0.000765005
-0.000262255
-0.000715978
-0.00059172
-0.000443981
-0.000285588
-1.92861e-05
-7.71446e-05
-2.35704e-05
-0.000153523
-0.000170177
-0.000275715
-2.35704e-05
-9.42817e-05
-0.000315732
-0.000180311
-0.000634821
-0.000397348
-0.00054224
-0.000851315
-0.000796378
-0.000696661
-0.000936411
-0.00108382
-0.000981779
-0.00106254
-0.0012218
-0.00105316
-0.0010168
-0.0011697
-0.00099108
-0.000964265
-0.00101041
-0.000992376
-0.000739054
-0.000821161
-0.000709278
-0.000850196
-0.000854996
-0.000375107
-0.00053044
-0.000747141
-0.000343662
-0.000448901
-0.000381786
-0.000808027
-0.000909199
-0.000867983
-0.000572495
-0.000794308
-0.00086756
-0.000885566
-0.000705928
-0.000841344
[9]:
print(a.mat)
Row 0: 0: 1
Row 1: 1: 0.828659
Row 2: 2: 1
Row 3: 3: 0.927121
Row 4: 0: -0.5 4: 1.89564
Row 5: 4: -0.368222 5: 1.76675
Row 6: 5: -0.323989 6: 1.75371
Row 7: 1: -0.215377 6: -0.28107 7: 1.83079
Row 8: 1: -0.216995 8: 1.89633
Row 9: 8: -0.370903 9: 1.74923
Row 10: 9: -0.28503 10: 1.77374
Row 11: 2: -0.5 10: -0.2708 11: 1.98742
Row 12: 2: -0.5 11: -0.107953 12: 1.8314
Row 13: 12: -0.226615 13: 1.75649
Row 14: 13: -0.340475 14: 1.75841
Row 15: 3: -0.430464 14: -0.405649 15: 1.80397
Row 16: 3: -0.40524 16: 2.07708
Row 17: 16: -0.357165 17: 1.76574
Row 18: 17: -0.4411 18: 2.03144
Row 19: 0: -0.5 4: -0.162559 18: -0.401142 19: 1.79475
Row 20: 4: -0.864858 5: -0.363931 19: -0.445914 20: 3.55639
Row 21: 5: -0.710604 6: -0.421297 20: -0.63885 21: 3.50636
Row 22: 6: -0.727351 7: -0.354504 21: -0.590012 22: 3.51941
Row 23: 1: -0.396287 7: -0.979835 8: -1.07358 22: -0.708758 23: 3.77912
Row 24: 8: -0.234858 9: -0.642955 22: -0.575333 23: -0.620659 24: 3.43908
Row 25: 9: -0.450342 10: -0.827115 24: -0.519524 25: 3.55949
Row 26: 10: -0.390798 11: -1.10866 12: -0.229033 25: -0.551848 26: 3.72126
Row 27: 12: -0.767794 13: -0.694967 26: -0.795614 27: 3.68173
Row 28: 13: -0.494437 14: -0.423252 27: -0.768236 28: 3.58147
Row 29: 14: -0.58903 15: -0.280159 28: -0.925237 29: 3.64107
Row 30: 3: -0.0914169 15: -0.687702 16: -1.31467 17: -0.291354 29: -0.918013 30: 3.88347
Row 31: 17: -0.167863 18: -1.1892 19: -0.285134 20: -0.656359 31: 3.78758
Row 32: 21: -0.721448 22: -0.56345 24: -0.348528 32: 3.65403
Row 33: 20: -0.586481 21: -0.424154 31: -0.727336 32: -0.714696 33: 3.50962
Row 34: 25: -0.359375 26: -0.645303 27: -0.655116 28: -0.505535 34: 3.53735
Row 35: 17: -0.508263 29: -0.60229 30: -0.580311 31: -0.761684 33: -0.427885 35: 3.51679
Row 36: 24: -0.497222 25: -0.851288 32: -1.05767 34: -0.787214 36: 3.74238
Row 37: 28: -0.464771 29: -0.326346 32: -0.248246 33: -0.629064 34: -0.584808 35: -0.636353 36: -0.54899 37: 3.43858
Row 38: 0: -0.0833333 4: 0 19: 0.0833333 38: 0.0416667
Row 39: 0: -0.0833333 4: 0.0833333 19: 0 38: 1.73472e-18 39: 0.0416667
Row 40: 1: -0.0328806 7: -0.0838095 23: 0.11669 40: 0.0381466
Row 41: 1: -0.0331673 8: -0.0828599 23: 0.116027 41: 0.0380482
Row 42: 1: -0.072062 7: 0.119706 8: 0.119026 23: -0.166669 40: -0.0209524 41: -0.020715 42: 0.0761948
Row 43: 2: -0.0833333 11: 0 12: 0.0833333 43: 0.0416667
Row 44: 2: -0.0833333 11: 0.0833333 12: 0 43: 0 44: 0.0416667
Row 45: 3: -0.0192424 15: -0.0611511 30: 0.0803935 45: 0.0380344
Row 46: 3: 0.00400625 16: -0.113562 30: 0.109556 46: 0.044274
Row 47: 3: -0.139284 15: 0.132895 16: 0.181102 30: -0.174713 45: -0.0152878 46: -0.0283905 47: 0.0823083
Row 48: 4: -0.0567722 5: -0.0292894 20: 0.0860616 48: 0.036858
Row 49: 0: 0.166667 4: -0.170704 19: -0.123322 20: 0.12736 38: -0.0208333 39: -0.0208333 49: 0.0802799
Row 50: 4: -0.0884634 5: 0.0906596 19: 0.0670821 20: -0.0692783 48: -0.00732235 49: -0.00999723 50: 0.0754712
Row 51: 5: -0.0548556 6: -0.0365843 21: 0.0914399 51: 0.0363595
Row 52: 4: 0.118142 5: -0.124949 20: -0.108911 21: 0.115717 48: -0.014193 50: -0.0153426 52: 0.0736287
Row 53: 5: -0.0853641 6: 0.0905825 20: 0.0835044 21: -0.0887229 51: -0.00914608 52: -0.0130346 53: 0.0731302
Row 54: 6: -0.0583949 7: -0.0399938 22: 0.0983887 54: 0.0363084
Row 55: 5: 0.108854 6: -0.116829 21: -0.104941 22: 0.112916 51: -0.0137139 53: -0.0134996 55: 0.0729964
Row 56: 6: -0.0804768 7: 0.0868387 21: 0.0837172 22: -0.0900791 54: -0.00999844 55: -0.0125213 56: 0.0729453
Row 57: 6: 0.10524 7: -0.126341 22: -0.113441 23: 0.134543 54: -0.0145987 56: -0.0117112 57: 0.0747167
Row 58: 1: 0.0687769 7: -0.0549866 22: 0.0741367 23: -0.087927 40: -0.00822015 42: -0.00897406 57: -0.0137616 58: 0.0765548
Row 59: 8: -0.0298036 9: -0.0556868 24: 0.0854904 59: 0.0368269
Row 60: 1: 0.069333 8: -0.0455052 23: -0.0905362 24: 0.0667084 41: -0.00829182 42: -0.00904144 60: 0.0787427
Row 61: 8: -0.157887 9: 0.117504 23: 0.153439 24: -0.113056 59: -0.0139217 60: -0.0143422 61: 0.0775214
Row 62: 9: -0.0313276 10: -0.0692128 25: 0.10054 62: 0.0370114
Row 63: 8: 0.0916208 9: -0.105547 24: -0.079105 25: 0.0930309 59: -0.0074509 61: -0.0154543 63: 0.0729527
Row 64: 9: -0.0989773 10: 0.116718 24: 0.100774 25: -0.118514 62: -0.0173032 63: -0.0123254 64: 0.0731372
Row 65: 10: -0.0314963 11: -0.0720729 26: 0.103569 65: 0.0371756
Row 66: 9: 0.0788326 10: -0.0811417 25: -0.0766521 26: 0.0789612 62: -0.00783189 64: -0.0118763 66: 0.0739116
Row 67: 10: -0.113773 11: 0.117206 25: 0.113964 26: -0.117397 65: -0.0180182 66: -0.0113311 67: 0.0740759
Row 68: 2: 0.166667 11: -0.196038 12: -0.120952 26: 0.150323 43: -0.0208333 44: -0.0208333 68: 0.0837455
Row 69: 10: 0.0766295 11: -0.0631255 12: 0.0556111 26: -0.0691151 65: -0.00787407 67: -0.0112833 68: -0.0094047 69: 0.0792544
Row 70: 12: -0.04799 13: -0.0598409 27: 0.107831 70: 0.0364
Row 71: 11: 0.130696 12: -0.0979679 26: -0.19839 27: 0.165661 68: -0.028176 69: -0.00449806 71: 0.0836325
Row 72: 12: -0.0383225 13: 0.09761 26: 0.0862391 27: -0.145527 70: -0.0149602 71: -0.0214214 72: 0.0779537
Row 73: 13: -0.0538257 14: -0.0351814 28: 0.0890071 73: 0.0364382
Row 74: 12: 0.0857591 13: -0.0663496 27: -0.111186 28: 0.0917764 70: -0.0119975 72: -0.00944228 74: 0.0733408
Row 75: 13: -0.112733 14: 0.0919272 27: 0.119183 28: -0.0983773 73: -0.00879534 74: -0.015799 75: 0.0733791
Row 76: 14: -0.0549708 15: -0.0263338 29: 0.0813045 76: 0.0372282
Row 77: 13: 0.110572 14: -0.0999467 28: -0.122776 29: 0.112151 73: -0.0134564 75: -0.0141865 77: 0.0733162
Row 78: 14: -0.102969 15: 0.093942 28: 0.104311 29: -0.095284 76: -0.00658344 77: -0.0172375 78: 0.0741061
Row 79: 14: 0.122579 15: -0.121074 29: -0.134292 30: 0.132787 76: -0.0137427 78: -0.0169021 79: 0.0755148
Row 80: 3: 0.0909864 15: -0.0921033 29: 0.0996805 30: -0.0985636 45: -0.0048106 47: -0.017936 79: -0.0198303 80: 0.076321
Row 81: 16: -0.10555 17: -0.00400625 30: 0.109556 81: 0.0422708
Row 82: 3: 0.0635338 16: -0.127068 17: 0.0635338 30: 6.93889e-18 46: 0.00100156 47: -0.016885 81: -0.00100156 82: 0.0865448
Row 83: 17: 0.00332986 18: -0.10243 31: 0.0991002 83: 0.0431542
Row 84: 16: 0.165077 17: -0.0995766 30: -0.166543 35: 0.101042 81: -0.0263874 82: -0.0148819 84: 0.0786696
Row 85: 17: -0.118178 18: 0.175947 31: -0.175437 35: 0.117668 83: -0.0256075 85: 0.0803981
Row 86: 17: -0.0758597 30: 0.105546 31: 0.104314 35: -0.134 84: -0.0152483 85: -0.0182518 86: 0.0736427
Row 87: 18: -0.0957703 19: -0.00332986 31: 0.0991002 87: 0.0414893
Row 88: 17: 0.0701869 18: -0.140374 19: 0.0701869 31: 0 83: 0.000832464 85: -0.0183792 87: -0.000832464 88: 0.0846435
Row 89: 4: 0.114464 19: -0.0712856 20: -0.156489 31: 0.11331 49: -0.0218427 50: -0.0067733 89: 0.0755233
Row 90: 18: 0.162627 19: -0.101187 20: 0.103448 31: -0.164888 87: -0.0239426 88: -0.0167143 89: -0.0172795 90: 0.0783994
Row 91: 5: 0.0949442 20: -0.0823254 21: -0.103233 33: 0.0906143 52: -0.0158946 53: -0.00784146 91: 0.0730084
Row 92: 19: 0.0785225 20: -0.0811174 31: -0.102313 33: 0.104907 89: -0.0110481 90: -0.00858251 92: 0.0732058
Row 93: 20: -0.0946116 21: 0.0939912 31: 0.0983953 33: -0.0977749 91: -0.00991368 92: -0.01453 93: 0.0725333
Row 94: 6: 0.0964622 21: -0.0896614 22: -0.1035 32: 0.0966993 55: -0.0157076 56: -0.00840796 94: 0.0728742
Row 95: 21: -0.0792876 22: 0.0889197 32: -0.0926539 33: 0.0830219 94: -0.0101674 95: 0.0730457
Row 96: 20: 0.105296 21: -0.118548 32: 0.116196 33: -0.102944 91: -0.0127399 93: -0.0135841 95: -0.0129961 96: 0.0730461
Row 97: 7: 0.0985866 22: -0.0560847 23: -0.125571 24: 0.0830686 57: -0.0198741 58: -0.00477258 97: 0.0749454
Row 98: 22: -0.116319 23: 0.109154 24: -0.0800426 32: 0.0872071 97: -0.0115186 98: 0.0730625
Row 99: 21: 0.10428 22: -0.107144 24: 0.0928628 32: -0.0899979 94: -0.0140074 95: -0.0120625 98: -0.00849209 99: 0.0727627
Row 100: 8: 0.105409 22: 0.100074 23: -0.15915 24: -0.0463337 60: -0.00233485 61: -0.0240174 97: -0.00924858 98: -0.01577 100: 0.0772316
Row 101: 9: 0.0952018 24: -0.0883787 25: -0.118783 36: 0.111959 63: -0.0109324 64: -0.0128681 101: 0.0734371
Row 102: 22: 0.112133 24: -0.104858 32: -0.136508 36: 0.129233 98: -0.0133097 99: -0.0147236 102: 0.0748636
Row 103: 24: -0.0614056 25: 0.112339 32: 0.107389 36: -0.158322 101: -0.0187633 102: -0.0208173 103: 0.0756495
Row 104: 10: 0.102276 25: -0.105873 26: -0.0957173 34: 0.0993136 66: -0.00840918 67: -0.0171599 104: 0.0733912
Row 105: 25: -0.113478 26: 0.108731 34: -0.122756 36: 0.127504 104: -0.0155202 105: 0.0740326
Row 106: 24: 0.0741923 25: -0.0599487 34: 0.0833384 36: -0.0975819 101: -0.00922659 103: -0.00932148 105: -0.0151689 106: 0.074853
Row 107: 12: 0.080529 26: -0.0460231 27: -0.132052 34: 0.0975463 71: -0.0199939 72: -0.000138337 107: 0.0776694
Row 108: 25: 0.0838832 26: -0.0935667 27: 0.0989931 34: -0.0893095 104: -0.00930825 105: -0.0116625 107: -0.0130191 108: 0.0726067
Row 109: 13: 0.0845674 27: -0.113096 28: -0.0551585 34: 0.0836875 74: -0.00714512 75: -0.0139967 109: 0.0740736
Row 110: 26: 0.0923864 27: -0.11176 28: 0.0914214 34: -0.0720478 107: -0.0113674 108: -0.0117291 109: -0.0066445 110: 0.0732485
Row 111: 14: 0.0785615 28: -0.0806626 29: -0.0668086 37: 0.0689097 77: -0.0108002 78: -0.00884017 111: 0.0754194
Row 112: 27: 0.121953 28: -0.0970033 34: -0.113762 37: 0.0888123 109: -0.0142774 110: -0.0162109 112: 0.0737553
Row 113: 28: -0.142934 29: 0.108864 34: 0.11433 37: -0.0802602 111: -0.00590194 112: -0.0141631 113: 0.075164
Row 114: 15: 0.0738253 29: -0.0597738 30: -0.0891911 35: 0.0751396 79: -0.0133665 80: -0.00508983 114: 0.0754918
Row 115: 29: -0.104464 30: 0.109406 35: -0.0909585 37: 0.0860166 114: -0.00893127 115: 0.0739511
Row 116: 28: 0.130558 29: -0.146223 35: 0.116201 37: -0.100535 111: -0.0113255 113: -0.021314 115: -0.0138084 116: 0.0752874
Row 117: 17: 0.0846018 29: 0.113095 30: -0.118234 35: -0.0794636 84: -0.0100123 86: -0.0111382 114: -0.00985363 115: -0.0184202 117: 0.0736039
Row 118: 20: 0.0870624 31: -0.0942154 33: -0.0770384 35: 0.0841913 92: -0.0116968 93: -0.0100688 118: 0.0731191
Row 119: 17: 0.0759684 31: -0.0944096 33: 0.0933537 35: -0.0749125 85: -0.0111653 86: -0.00782675 118: -0.00756278 119: 0.0740674
Row 120: 21: 0.0952493 32: -0.0932157 33: -0.0830199 37: 0.0809863 95: -0.00775942 96: -0.0160529 120: 0.0738379
Row 121: 24: 0.0700837 32: -0.03649 36: -0.100946 37: 0.0673519 102: -0.011491 103: -0.00602993 121: 0.0784283
Row 122: 32: -0.16014 33: 0.119114 36: 0.14799 37: -0.106964 120: -0.0129955 121: -0.0137454 122: 0.0771196
Row 123: 31: 0.117043 33: -0.115964 35: -0.104766 37: 0.103687 118: -0.0134851 119: -0.0157757 123: 0.0730111
Row 124: 32: 0.0961358 33: -0.108195 35: 0.0918885 37: -0.0798295 120: -0.00725102 122: -0.0167829 123: -0.0127064 124: 0.0732171
Row 125: 25: 0.0894911 34: -0.0634784 36: -0.103345 37: 0.0773323 105: -0.0167071 106: -0.00566571 125: 0.0745064
Row 126: 28: 0.0898378 34: -0.128205 36: 0.107043 37: -0.0686766 112: -0.00803997 113: -0.0144195 125: -0.00912918 126: 0.0735873
Row 127: 29: 0.0917503 33: 0.0939249 35: -0.10203 37: -0.083645 115: -0.0076958 116: -0.0152418 123: -0.0132155 124: -0.0102658 127: 0.0729335
Row 128: 32: 0.105379 34: 0.111342 36: -0.163535 37: -0.0531859 121: -0.00309257 122: -0.0232522 125: -0.0102039 126: -0.0176317 128: 0.0770549
6. Solve the system:¶
[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:
[11]:
print(gfu.vec)
0
0
0
0.0923002
0
0
0
0
0
0
0
0
0.0578974
0.0863383
0.0954046
0.0944886
0.088851
0.0779859
0.0596329
0.0330695
0.0358778
0.0367923
0.0321534
0.0187666
0.039323
0.0440006
0.0473955
0.0754942
0.089345
0.0924018
0.0900675
0.0609401
0.0584553
0.0636311
0.0735012
0.0808489
0.0646725
0.0806824
0
-0.005753
0
0
0.0210629
0
-0.035086
0.00263586
-0.00722906
-0.0053199
0
-0.00990601
-0.00627918
0
-0.016323
-0.0118361
0
-0.0149549
-0.0195579
-0.00642189
-0.0212263
0
-0.025571
-0.0101959
0
-0.0393962
-0.0141925
0
-0.0275764
-0.0181556
-0.033372
-0.0244469
-0.0242966
-0.00443999
-0.0108775
-0.0146141
-0.0063948
-0.00774334
-0.00565315
-0.00909429
-0.00785633
-0.00847728
-0.00755551
-0.00774445
0.00263616
-0.00823406
-0.00507479
-0.00539107
0.00395556
-0.00824646
0.00366889
0.00416289
-0.00432812
-0.00340091
-0.00931158
-0.0100027
-0.00716045
-0.0135339
-0.0116171
-0.00788447
-0.0251341
-0.00768871
-0.00317827
-0.00903235
-0.0216162
-0.00597072
-0.0032929
-0.0141649
-0.0192129
-0.014382
-0.0219207
-0.0142459
-0.0043319
-0.00964124
-0.0103001
-0.0100999
-0.00227504
-0.00914321
-0.0102146
-0.00700589
-0.00179493
-0.00840332
-0.00865067
-0.0115504
-0.0105298
-0.00906982
-0.0131517
-0.00454235
-0.01852
-0.00882419
-0.0126031
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 module 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
).