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]:
(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:¶
[3]:
fes = H1(mesh, order=2, dirichlet="bottom|right")
fes.ndof # number of unknowns in this space
[3]:
133
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'
| dirichlet_bbnd: regexpr
| Regular expression string defining the dirichlet bboundary,
| i.e. points in 2D and edges in 3D.
| More than one boundary can be combined by the | operator,
| i.e.: dirichlet_bbnd = '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
|
| loembedding
|
| 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.00888574
0.00633333
0.00074628
0.00399953
0.00769352
0.0104008
0.0115294
0.0150217
0.0167909
0.015258
0.0182768
0.0211091
0.0106045
0.00711336
0.00320035
0.000735854
0.000876249
0.000529403
0.00282218
0.0124495
0.0182096
0.0224435
0.0204234
0.0317338
0.0270543
0.0291603
0.0187131
0.0217393
0.00854715
0.00505989
0.0038006
0.00883394
0.0205505
0.0163797
0.0213793
0.0150785
0.0170279
0.0191555
-6.66667e-05
-3.33333e-05
-0.000526369
-0.000575033
-0.00109143
-0.0008
-0.000766667
-8.08556e-05
-2.24397e-05
-0.000120589
-0.000249038
-0.000218716
-0.000440093
-0.000380259
-0.000572665
-0.000711098
-0.000487162
-0.00080946
-0.000916707
-0.000914728
-0.000944161
-0.000759221
-0.00110619
-0.00127157
-0.000646306
-0.0014197
-0.00131979
-0.000648788
-0.00123697
-0.0012392
-0.0015957
-0.00146355
-0.000469258
-0.0014082
-0.00106589
-0.000437144
-0.000880349
-0.000856527
-0.000190932
-0.000699086
-0.000461922
-0.000315546
-0.000220954
-1.64401e-05
-9.82846e-05
-8.35917e-05
-2.64702e-05
-9.73025e-05
-0.000122662
-2.64702e-05
-0.000105881
-0.000353264
-0.00021487
-0.00066398
-0.000468351
-0.00059592
-0.000843865
-0.000796408
-0.000722285
-0.000922806
-0.00105484
-0.000955488
-0.00104202
-0.00117171
-0.00100951
-0.000963748
-0.00114996
-0.000965717
-0.00092025
-0.00104152
-0.00098591
-0.000825307
-0.000865806
-0.000574259
-0.000880372
-0.000710639
-0.000885225
-0.000308496
-0.000477574
-0.000191583
-0.00032424
-0.000250767
-0.000328592
-0.000559968
-0.00047995
-0.000773918
-0.00082872
-0.000774775
-0.000675938
-0.000761632
-0.000801479
-0.000841885
-0.000775647
-0.000743459
[9]:
print(a.mat)
Row 0: 0: 1
Row 1: 1: 0.828441
Row 2: 2: 1
Row 3: 3: 0.870927
Row 4: 0: -0.5 4: 1.88347
Row 5: 4: -0.397227 5: 1.77564
Row 6: 5: -0.326749 6: 1.7517
Row 7: 1: -0.208968 6: -0.273913 7: 1.82559
Row 8: 1: -0.209694 8: 1.8809
Row 9: 8: -0.350939 9: 1.74429
Row 10: 9: -0.262838 10: 1.77715
Row 11: 2: -0.5 10: -0.267845 11: 1.9829
Row 12: 2: -0.5 11: -0.115438 12: 1.81266
Row 13: 12: -0.253469 13: 1.79007
Row 14: 13: -0.404897 14: 1.83152
Row 15: 3: -0.344377 14: -0.218833 15: 1.78023
Row 16: 3: -0.336878 16: 1.81347
Row 17: 16: -0.20561 17: 1.94288
Row 18: 17: -0.368349 18: 2.04376
Row 19: 0: -0.5 4: -0.177276 18: -0.553239 19: 1.87497
Row 20: 4: -0.808965 5: -0.345183 19: -0.545937 20: 3.54246
Row 21: 5: -0.706484 6: -0.432506 20: -0.647227 21: 3.50609
Row 22: 6: -0.718531 7: -0.374026 21: -0.571505 22: 3.51428
Row 23: 1: -0.40978 7: -0.968685 8: -1.05698 22: -0.708085 23: 3.76596
Row 24: 8: -0.263291 9: -0.64648 22: -0.536887 23: -0.62243 24: 3.433
Row 25: 9: -0.484029 10: -0.84644 24: -0.459506 25: 3.57375
Row 26: 10: -0.400025 11: -1.09961 12: -0.239874 25: -0.517151 26: 3.71385
Row 27: 12: -0.703881 13: -0.775929 26: -0.788273 27: 3.68527
Row 28: 13: -0.355771 14: -0.309609 27: -0.733256 28: 3.48151
Row 29: 14: -0.898185 15: -0.432178 28: -0.790531 29: 3.74275
Row 30: 3: -0.189673 15: -0.784846 16: -0.908731 29: -0.98463 30: 3.81689
Row 31: 16: -0.36225 17: -1.16439 30: -0.851669 31: 3.85451
Row 32: 17: -0.204537 18: -1.12218 19: -0.098519 20: -0.584699 31: -0.517495 32: 3.69368
Row 33: 21: -0.714154 22: -0.60525 24: -0.354764 33: 3.57925
Row 34: 20: -0.610448 21: -0.43421 32: -0.800343 33: -0.624553 34: 3.5401
Row 35: 25: -0.379385 26: -0.66891 27: -0.683932 28: -0.449095 35: 3.53739
Row 36: 28: -0.434461 29: -0.637229 30: -0.0973365 31: -0.958711 32: -0.365907 34: -0.341534 36: 3.57061
Row 37: 24: -0.549644 25: -0.887236 33: -0.928035 35: -0.674877 37: 3.68681
Row 38: 28: -0.408789 33: -0.352491 34: -0.729009 35: -0.68119 36: -0.735427 37: -0.64702 38: 3.55392
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.0340869 7: -0.0835417 23: 0.117629 41: 0.0381141
Row 42: 1: -0.0342097 8: -0.0831255 23: 0.117335 42: 0.038071
Row 43: 1: -0.069777 7: 0.11837 8: 0.118074 23: -0.166667 41: -0.0208854 42: -0.0207814 43: 0.0761852
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.0173948 15: -0.0795022 30: 0.096897 46: 0.0385733
Row 47: 3: -0.0142173 16: -0.0873491 30: 0.101566 47: 0.0394282
Row 48: 3: -0.113542 15: 0.136898 16: 0.143495 30: -0.166851 46: -0.0198756 47: -0.0218373 48: 0.0780015
Row 49: 4: -0.0535736 5: -0.028366 20: 0.0819397 49: 0.037036
Row 50: 0: 0.166667 4: -0.164587 19: -0.124342 20: 0.122262 39: -0.0208333 40: -0.0208333 50: 0.0796187
Row 51: 4: -0.0957505 5: 0.0945705 19: 0.0705543 20: -0.0693744 49: -0.00709151 50: -0.0102521 51: 0.0749881
Row 52: 5: -0.0526905 6: -0.0380314 21: 0.0907219 52: 0.036295
Row 53: 4: 0.119778 5: -0.131261 20: -0.10714 21: 0.118623 49: -0.0133934 51: -0.0165511 53: 0.073983
Row 54: 5: -0.0836226 6: 0.0924896 20: 0.082731 21: -0.091598 52: -0.00950784 53: -0.0133916 54: 0.073242
Row 55: 6: -0.0573097 7: -0.0420363 22: 0.099346 55: 0.0362495
Row 56: 5: 0.107149 6: -0.116904 21: -0.102619 22: 0.112374 52: -0.0131726 54: -0.0136146 56: 0.0729017
Row 57: 6: -0.0797051 7: 0.0876885 21: 0.0839811 22: -0.0919644 55: -0.0105091 56: -0.012482 57: 0.0728562
Row 58: 6: 0.102962 7: -0.123558 22: -0.111917 23: 0.132513 55: -0.0143274 57: -0.011413 58: 0.0744532
Row 59: 1: 0.0689149 7: -0.0551294 22: 0.0749089 23: -0.0886944 41: -0.00852172 43: -0.008707 58: -0.0136519 59: 0.0763178
Row 60: 8: -0.0322109 9: -0.0557926 24: 0.0880035 60: 0.0366233
Row 61: 1: 0.0691587 8: -0.0466199 23: -0.0901615 24: 0.0676228 42: -0.00855243 43: -0.00873724 61: 0.078236
Row 62: 8: -0.151527 9: 0.114282 23: 0.148989 24: -0.111744 60: -0.0139482 61: -0.013988 62: 0.0767883
Row 63: 9: -0.0342301 10: -0.0697744 25: 0.104004 63: 0.0369527
Row 64: 8: 0.0907008 9: -0.104931 24: -0.0782661 25: 0.0924966 60: -0.00805273 62: -0.0146225 64: 0.072736
Row 65: 9: -0.0957605 10: 0.113581 24: 0.0980093 25: -0.11583 63: -0.0174436 64: -0.0115138 65: 0.0730654
Row 66: 10: -0.0327463 11: -0.0708467 26: 0.103593 66: 0.0370585
Row 67: 9: 0.0780365 10: -0.0777309 25: -0.0772401 26: 0.0769345 63: -0.00855753 65: -0.0109516 67: 0.0740111
Row 68: 10: -0.11594 11: 0.115488 25: 0.114309 26: -0.113857 66: -0.0177117 67: -0.0107525 68: 0.0741168
Row 69: 2: 0.166667 11: -0.195755 12: -0.11965 26: 0.148738 44: -0.0208333 45: -0.0208333 69: 0.0836612
Row 70: 10: 0.0773871 11: -0.0638806 12: 0.0555561 26: -0.0690626 66: -0.00818658 68: -0.0111602 69: -0.00907908 70: 0.079053
Row 71: 12: -0.0432149 13: -0.0598976 27: 0.103112 71: 0.0363393
Row 72: 11: 0.131662 12: -0.0933383 26: -0.198236 27: 0.159913 69: -0.0281055 70: -0.00480994 72: 0.0828885
Row 73: 12: -0.0459076 13: 0.102142 26: 0.0894771 27: -0.145712 71: -0.0149744 72: -0.0214536 73: 0.0772333
Row 74: 13: -0.0406718 14: -0.0388314 28: 0.0795032 74: 0.0367465
Row 75: 12: 0.0854598 13: -0.0608683 27: -0.107402 28: 0.0828109 71: -0.0108037 73: -0.0105612 75: 0.074398
Row 76: 13: -0.136907 14: 0.106314 27: 0.133611 28: -0.103019 74: -0.00970785 75: -0.0160469 76: 0.0748052
Row 77: 14: -0.0683139 15: -0.0424951 29: 0.110809 77: 0.0368203
Row 78: 13: 0.108155 14: -0.148866 28: -0.10339 29: 0.144102 74: -0.010168 76: -0.0168707 78: 0.0759645
Row 79: 14: -0.0492423 15: 0.0789672 28: 0.0754885 29: -0.105213 77: -0.0106238 78: -0.0156796 79: 0.0760383
Row 80: 14: 0.104786 15: -0.0877775 29: -0.135473 30: 0.118465 77: -0.0170785 79: -0.00911803 80: 0.0738201
Row 81: 3: 0.074791 15: -0.0869308 29: 0.0966939 30: -0.0845541 46: -0.00434871 48: -0.014349 80: -0.0167898 81: 0.0755731
Row 82: 16: -0.033994 17: -0.0846663 31: 0.11866 82: 0.0382322
Row 83: 3: 0.0703636 16: -0.0825273 30: -0.0722726 31: 0.0844363 47: -0.00355433 48: -0.0140366 83: 0.0765638
Row 84: 16: -0.0983745 17: 0.118935 30: 0.122161 31: -0.142722 82: -0.0211666 83: -0.0145138 84: 0.0753678
Row 85: 17: -0.0154075 18: -0.0781072 32: 0.0935147 85: 0.0387265
Row 86: 16: 0.0682624 17: -0.0529503 31: -0.0722565 32: 0.0569444 82: -0.0084985 84: -0.00856709 86: 0.0798179
Row 87: 17: -0.17079 18: 0.139499 31: 0.147661 32: -0.11637 85: -0.0195268 86: -0.00956562 87: 0.0803122
Row 88: 18: -0.108922 19: 0.0154075 32: 0.0935147 88: 0.0464303
Row 89: 17: 0.076799 18: -0.153598 19: 0.076799 32: 0 85: -0.00385189 87: -0.0153479 88: 0.00385189 89: 0.0851569
Row 90: 4: 0.1108 19: -0.0613734 20: -0.146695 32: 0.0972688 50: -0.0203135 51: -0.0073865 90: 0.0747646
Row 91: 18: 0.201129 19: -0.142188 20: 0.115423 32: -0.174364 88: -0.0272306 89: -0.0230516 90: -0.0163604 91: 0.0832428
Row 92: 5: 0.0942211 20: -0.077413 21: -0.107223 34: 0.0904152 53: -0.0162642 54: -0.00729109 92: 0.0731269
Row 93: 19: 0.0818086 20: -0.103474 32: -0.0930223 34: 0.114688 90: -0.00795685 91: -0.0124953 93: 0.0734865
Row 94: 20: -0.0863131 21: 0.0964712 32: 0.0932033 34: -0.103361 92: -0.0105417 93: -0.0152987 94: 0.0728539
Row 95: 6: 0.0964985 21: -0.0880292 22: -0.107744 33: 0.0992748 56: -0.0156114 57: -0.00851324 95: 0.072756
Row 96: 21: -0.0755245 22: 0.0906212 33: -0.0975796 34: 0.0824828 95: -0.0113246 96: 0.0730324
Row 97: 20: 0.102553 21: -0.119354 33: 0.11733 34: -0.10053 92: -0.0120622 94: -0.0135762 96: -0.0130703 97: 0.073063
Row 98: 7: 0.0982072 22: -0.0555055 23: -0.125692 24: 0.0829906 58: -0.0194765 59: -0.00507535 98: 0.074803
Row 99: 22: -0.118983 23: 0.111193 24: -0.0835425 33: 0.0913324 98: -0.0119466 99: 0.0730017
Row 100: 21: 0.0992989 22: -0.0995997 24: 0.0900331 33: -0.0897323 95: -0.0134941 96: -0.0113307 99: -0.00893901 100: 0.0725517
Row 101: 8: 0.104708 22: 0.0986109 23: -0.156444 24: -0.0468751 61: -0.00291774 62: -0.0232593 98: -0.00880103 99: -0.0158517 101: 0.0767643
Row 102: 9: 0.0983955 24: -0.0942456 25: -0.124075 37: 0.119925 64: -0.0116104 65: -0.0129885 102: 0.0737262
Row 103: 22: 0.109853 24: -0.103593 33: -0.135259 37: 0.128998 99: -0.0138941 100: -0.0135693 103: 0.0744947
Row 104: 24: -0.0539003 25: 0.108163 33: 0.103054 37: -0.157316 102: -0.0194084 103: -0.0199206 104: 0.0757059
Row 105: 10: 0.105223 25: -0.111385 26: -0.0965358 35: 0.102697 67: -0.00848113 68: -0.0178247 105: 0.073528
Row 106: 25: -0.113421 26: 0.105793 35: -0.119567 37: 0.127194 105: -0.0156528 106: 0.0740546
Row 107: 24: 0.0728206 25: -0.0536744 35: 0.0801006 37: -0.0992469 102: -0.0105729 104: -0.00763227 106: -0.0142388 107: 0.0751985
Row 108: 12: 0.0777613 26: -0.0525365 27: -0.124052 35: 0.0988275 72: -0.0185246 73: -0.000915684 108: 0.0769919
Row 109: 25: 0.0832674 26: -0.0887462 27: 0.0955182 35: -0.0900393 105: -0.0100214 106: -0.0107954 108: -0.0124884 109: 0.0725676
Row 110: 13: 0.0880473 27: -0.133459 28: -0.0450785 35: 0.09049 75: -0.00465584 76: -0.017356 110: 0.0751867
Row 111: 26: 0.0944382 27: -0.103586 28: 0.0844769 35: -0.0753288 108: -0.0122184 109: -0.0113911 110: -0.00661378 111: 0.0732259
Row 112: 14: 0.0941536 28: -0.0400872 29: -0.133883 36: 0.079817 78: -0.0203459 79: -0.00319254 112: 0.0764314
Row 113: 27: 0.122057 28: -0.0863546 35: -0.136673 38: 0.100971 110: -0.0160087 111: -0.0145055 113: 0.0744692
Row 114: 28: -0.108835 29: 0.121537 36: -0.113163 38: 0.100462 112: -0.013125 114: 0.0736021
Row 115: 28: -0.0934871 35: 0.121032 36: 0.105756 38: -0.133301 113: -0.0181596 114: -0.0151658 115: 0.07373
Row 116: 15: 0.0808401 29: -0.0832396 30: -0.0628419 36: 0.0652415 80: -0.0128264 81: -0.00738366 116: 0.0775466
Row 117: 28: 0.0963537 29: -0.165982 30: 0.108482 36: -0.0388535 112: -0.00682926 114: -0.0172592 116: -0.00288412 117: 0.0777602
Row 118: 16: 0.0904872 30: -0.0687924 31: -0.100344 36: 0.078649 83: -0.00659526 84: -0.0160265 118: 0.0777703
Row 119: 29: 0.150651 30: -0.180835 31: 0.157852 36: -0.127668 116: -0.0134262 117: -0.0242364 118: -0.0184907 119: 0.0811814
Row 120: 17: 0.12808 31: -0.195221 32: -0.0398024 36: 0.106943 86: -0.00467048 87: -0.0273496 120: 0.0803181
Row 121: 30: 0.0885758 31: -0.131876 32: 0.0691071 36: -0.0258068 118: -0.00117156 119: -0.0209724 120: -0.00528012 121: 0.0793671
Row 122: 20: 0.0855011 32: -0.0718724 34: -0.0897809 36: 0.0761522 93: -0.0133732 94: -0.00800211 122: 0.0737609
Row 123: 31: 0.133809 32: -0.120182 34: 0.108484 36: -0.122111 120: -0.0214556 121: -0.0119967 122: -0.00907206 123: 0.0758193
Row 124: 21: 0.0952513 33: -0.105243 34: -0.083047 38: 0.0930384 96: -0.00755045 97: -0.0162624 124: 0.0730955
Row 125: 24: 0.072687 33: -0.0419266 37: -0.108677 38: 0.0779166 103: -0.0123289 104: -0.00584281 125: 0.076319
Row 126: 33: -0.126801 34: 0.104656 37: 0.134351 38: -0.112206 124: -0.0132113 125: -0.0148403 126: 0.0744391
Row 127: 32: 0.11206 34: -0.140852 36: -0.101772 38: 0.130564 122: -0.00996599 123: -0.0180489 127: 0.0748865
Row 128: 33: 0.0920044 34: -0.0724453 36: 0.0825421 38: -0.102101 124: -0.0100483 126: -0.0129528 127: -0.015477 128: 0.074012
Row 129: 25: 0.0933844 35: -0.0640386 37: -0.118715 38: 0.0893687 106: -0.0175598 107: -0.00578633 129: 0.0738082
Row 130: 28: 0.0767268 35: -0.103918 37: 0.104 38: -0.0768082 113: -0.00708319 115: -0.0120985 129: -0.0121189 130: 0.0735645
Row 131: 28: 0.0848918 34: 0.0892903 36: -0.0657272 38: -0.108455 114: -0.00994967 115: -0.0112733 127: -0.0171641 128: -0.00515852 131: 0.0741883
Row 132: 33: 0.0935455 35: 0.0964175 37: -0.130514 38: -0.0594486 125: -0.00463885 126: -0.0187475 129: -0.0102233 130: -0.0138811 132: 0.0744499
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.0923044
0
0
0
0
0
0
0
0
0.0578972
0.0863393
0.0954128
0.0944891
0.0888235
0.0780439
0.059628
0.0330628
0.0374521
0.0370446
0.0317995
0.0185123
0.0383632
0.0427969
0.047351
0.0760518
0.0895602
0.0939297
0.0914143
0.0833171
0.0654761
0.0574569
0.0647326
0.0723792
0.0852186
0.0626708
0.0776805
0
-0.00576774
0
0
0.0208108
0
-0.0350847
0.00263379
-0.00729372
-0.0027615
0
-0.00993549
-0.00711876
0
-0.0174982
-0.0120617
0
-0.0150936
-0.0195578
-0.00638541
-0.0213554
0
-0.0254715
-0.00955261
0
-0.0381111
-0.0132127
0
-0.0262138
-0.0182002
-0.0333693
-0.0245117
-0.0242838
-0.00451911
-0.0119297
-0.014536
-0.00612786
-0.0105674
-0.00558557
-0.00891829
-0.00493114
-0.0051195
-0.00519767
-0.007602
0.00247896
-0.00133557
-0.00800354
0.00157024
-0.00177994
-0.00829596
0.00374357
0.00436653
-0.00759115
-0.00373981
-0.0120787
-0.0104187
-0.00708368
-0.0136752
-0.0120743
-0.00778966
-0.0256197
-0.0072408
-0.00331602
-0.00900117
-0.0206949
-0.0050238
-0.00300073
-0.0140407
-0.0179806
-0.0153171
-0.020859
-0.0164647
-0.00404961
-0.00810738
-0.0130797
-0.0135373
-0.00892006
-0.00206078
-0.00747051
-0.00532125
-0.00675915
-0.00966483
-0.00121964
-0.00276467
-0.0111685
-0.00868831
-0.0118313
-0.00814322
-0.0124053
-0.0120342
-0.00417412
-0.0154546
-0.00866456
-0.010265
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
).