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]:
from ngsolve import *
from ngsolve.webgui import Draw
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
| element-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'
| dirichlet_bbbnd: regexpr
| Regular expression string defining the dirichlet bbboundary,
| i.e. points in 3D.
| More than one boundary can be combined by the | operator,
| i.e.: dirichlet_bbbnd = '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,
| VARIABLE ... 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) -> 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
|
| CreateDirectSolverCluster(...)
| CreateDirectSolverCluster(self: ngsolve.comp.FESpace, **kwargs) -> list
|
| Elements(...)
| Elements(self: ngsolve.comp.FESpace, VOL_or_BND: ngsolve.comp.VorB = <VorB.VOL: 0>) -> 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) -> BaseMatrix
|
| Mass(...)
| Mass(self: ngsolve.comp.FESpace, rho: ngsolve.fem.CoefficientFunction = None, definedon: Optional[ngsolve.comp.Region] = None) -> 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, arg0: int) -> ngsolve.la.DofRange
|
| deprecated, will be only available for ProductSpace
|
| 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) -> 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
|
| __mul__(...)
| __mul__(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.FESpace) -> ngcomp::CompoundFESpace
|
| __pow__(...)
| __pow__(self: ngsolve.comp.FESpace, arg0: int) -> ngcomp::CompoundFESpaceAllSame
|
| __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
|
| ----------------------------------------------------------------------
| Readonly properties inherited from FESpace:
|
| components
| deprecated, will be only available for ProductSpace
|
| 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 and other attributes inherited from FESpace:
|
| __hash__ = None
|
| ----------------------------------------------------------------------
| Readonly properties inherited from NGS_Object:
|
| __memory__
|
| ----------------------------------------------------------------------
| Data descriptors inherited from NGS_Object:
|
| name
|
| ----------------------------------------------------------------------
| Static 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()
[7]:
<ngsolve.comp.LinearForm at 0x7f79dcaa2c70>
You can examine the linear system in more detail:
[8]:
print(f.vec)
0.000333333
0.00888667
0.00633333
0.000772552
0.00388886
0.00745593
0.0102799
0.0115007
0.0150668
0.0169176
0.0154171
0.0183849
0.0212923
0.0109064
0.00750505
0.00352697
0.000376105
0.00262832
0.00051516
0.00266077
0.0117876
0.0178168
0.0223235
0.0204184
0.0318376
0.0272957
0.029483
0.0190557
0.0227156
0.00905798
0.00529671
0.00779569
0.020463
0.0157845
0.0216524
0.0157801
0.0172011
0.019586
-6.66667e-05
-3.33333e-05
-0.000525777
-0.000575773
-0.00109151
-0.0008
-0.000766667
-9.18288e-05
-1.88053e-05
-0.000121132
-0.000238755
-0.000215067
-0.000426062
-0.000372909
-0.000551384
-0.000691744
-0.000485002
-0.000796725
-0.000909071
-0.000911696
-0.000942778
-0.000763411
-0.00110769
-0.0012761
-0.000653931
-0.00142763
-0.00133086
-0.000654795
-0.00125029
-0.00125106
-0.0016032
-0.00147642
-0.000475024
-0.00142506
-0.0010807
-0.00045191
-0.000903122
-0.000886228
-0.000204734
-0.00073304
-0.000494315
-0.00034316
-0.000246364
-1.88053e-05
-7.52211e-05
-2.5758e-05
-0.000173161
-0.000213128
-0.000357643
-2.5758e-05
-0.000103032
-0.000329228
-0.000194844
-0.00064017
-0.000422578
-0.000560041
-0.000834752
-0.000785222
-0.00070094
-0.000921418
-0.00105231
-0.000949261
-0.00104196
-0.00117992
-0.00101029
-0.000970793
-0.00116359
-0.000974281
-0.000925954
-0.00105527
-0.00099741
-0.000847599
-0.00087727
-0.000615988
-0.000891989
-0.000762229
-0.000931528
-0.000320916
-0.000490049
-0.0003756
-0.000511811
-0.000491896
-0.000762552
-0.000835836
-0.000773654
-0.000642566
-0.00076271
-0.00081234
-0.000855256
-0.000808305
-0.000755427
[9]:
print(a.mat)
Row 0: 0: 1 4: -0.5 19: -0.5 38: -0.0833333 39: -0.0833333 49: 0.166667
Row 1: 1: 0.828444 7: -0.208707 8: -0.210132 23: -0.409604 40: -0.034013 41: -0.0342543 42: -0.0698066 58: 0.0687976 60: 0.0692764
Row 2: 2: 1 11: -0.5 12: -0.5 43: -0.0833333 44: -0.0833333 68: 0.166667
Row 3: 3: 0.907169 15: -0.407338 16: -0.376463 30: -0.123368 45: -0.020713 46: 0.000151649 47: -0.130634 80: 0.0886026 82: 0.0625923
Row 4: 0: -0.5 4: 1.8918 5: -0.369703 19: -0.170695 20: -0.851406 38: 0 39: 0.0833333 48: -0.0553087 49: -0.169926 50: -0.0900664 52: 0.116926 89: 0.115041
Row 5: 4: -0.369703 5: 1.76479 6: -0.315355 20: -0.373272 21: -0.70646 48: -0.0302454 50: 0.0918626 51: -0.0537908 52: -0.12557 53: -0.0845258 55: 0.10635 91: 0.0959192
Row 6: 5: -0.315355 6: 1.74934 7: -0.271622 21: -0.444807 22: -0.71756 51: -0.0387141 53: 0.0912733 54: -0.0577312 55: -0.114421 56: -0.0806907 57: 0.103002 94: 0.0972825
Row 7: 1: -0.208707 6: -0.271622 7: 1.82575 22: -0.375475 23: -0.969947 40: -0.0837431 42: 0.118528 54: -0.0420473 56: 0.0873176 57: -0.123185 58: -0.0553164 97: 0.0984466
Row 8: 1: -0.210132 8: 1.88178 9: -0.354044 23: -1.05651 24: -0.261086 41: -0.0829256 42: 0.117948 59: -0.0322254 60: -0.0463111 61: -0.152168 63: 0.0912327 100: 0.104449
Row 9: 8: -0.354044 9: 1.7442 10: -0.268795 24: -0.640588 25: -0.480773 59: -0.0552752 61: 0.114283 62: -0.0341102 63: -0.105026 64: -0.0962885 66: 0.0789093 101: 0.0975081
Row 10: 9: -0.268795 10: 1.77359 11: -0.272339 25: -0.832083 26: -0.400373 62: -0.0686401 64: 0.113439 65: -0.0327251 66: -0.0788029 67: -0.11543 69: 0.078115 104: 0.104044
Row 11: 2: -0.5 10: -0.272339 11: 1.97625 12: -0.118504 26: -1.0854 43: 0 44: 0.0833333 65: -0.0698849 67: 0.115275 68: -0.194349 69: -0.0651405 71: 0.130766
Row 12: 2: -0.5 11: -0.118504 12: 1.80573 13: -0.26228 26: -0.248259 27: -0.676685 43: 0.0833333 44: 0 68: -0.119671 69: 0.0560888 70: -0.0411842 71: -0.0913473 72: -0.0487517 74: 0.0848976 107: 0.0766349
Row 13: 12: -0.26228 13: 1.81175 14: -0.428891 27: -0.813095 28: -0.307486 70: -0.0605924 72: 0.104306 73: -0.0361696 74: -0.0587915 75: -0.146405 77: 0.107651 109: 0.0900017
Row 14: 13: -0.428891 14: 1.79855 15: -0.254586 28: -0.343721 29: -0.771353 73: -0.0404916 75: 0.111973 76: -0.0578574 77: -0.142183 78: -0.0592263 79: 0.100288 111: 0.0874967
Row 15: 3: -0.407338 14: -0.254586 15: 1.75967 29: -0.52416 30: -0.573581 45: -0.0625066 47: 0.130396 76: -0.0447659 78: 0.0871969 79: -0.0755213 80: -0.110484 115: 0.0756844
Row 16: 3: -0.376463 16: 2.08248 17: -0.374644 30: -1.33137 46: -0.111099 47: 0.173843 81: -0.110796 82: -0.125185 84: 0.173237
Row 17: 16: -0.374644 17: 1.78398 18: -0.43722 30: -0.336404 31: -0.30274 35: -0.332973 81: -0.000151649 82: 0.0625923 83: -0.000386213 84: -0.0860908 85: -0.104715 86: -0.105986 88: 0.0732563 117: 0.0795659 119: 0.0819158
Row 18: 17: -0.43722 18: 2.01666 19: -0.441855 31: -1.13758 83: -0.0944124 85: 0.167282 87: -0.0951848 88: -0.146513 90: 0.168827
Row 19: 0: -0.5 4: -0.170695 18: -0.441855 19: 1.80929 20: -0.457008 31: -0.239727 38: 0.0833333 39: 0 49: -0.122284 50: 0.0674001 87: 0.000386213 88: 0.0732563 89: -0.0687899 90: -0.11086 92: 0.0775578
Row 20: 4: -0.851406 5: -0.373272 19: -0.457008 20: 3.55115 21: -0.621608 31: -0.660832 33: -0.587027 48: 0.0855542 49: 0.125543 50: -0.0691963 52: -0.106394 53: 0.0830523 89: -0.156773 90: 0.107398 91: -0.0811964 92: -0.0858253 93: -0.0924734 96: 0.101746 118: 0.0885658
Row 21: 5: -0.70646 6: -0.444807 20: -0.621608 21: 3.49846 22: -0.562802 32: -0.707916 33: -0.454863 51: 0.0925049 52: 0.115038 53: -0.0897997 55: -0.102651 56: 0.0842808 91: -0.106796 93: 0.095359 94: -0.0895218 95: -0.0779072 96: -0.1164 99: 0.0990414 120: 0.0968518
Row 22: 6: -0.71756 7: -0.375475 21: -0.562802 22: 3.51385 23: -0.705204 24: -0.539682 32: -0.613124 54: 0.0997784 55: 0.110723 56: -0.0909077 57: -0.112022 58: 0.0748223 94: -0.10743 95: 0.090508 97: -0.055892 98: -0.119863 99: -0.0995268 100: 0.0986036 102: 0.111206
Row 23: 1: -0.409604 7: -0.969947 8: -1.05651 22: -0.705204 23: 3.76626 24: -0.624994 40: 0.117756 41: 0.11718 42: -0.166669 57: 0.132205 58: -0.0883035 60: -0.0906718 61: 0.149578 97: -0.125663 98: 0.110992 100: -0.156404
Row 24: 8: -0.261086 9: -0.640588 22: -0.539682 23: -0.624994 24: 3.4335 25: -0.471228 32: -0.347567 36: -0.548354 59: 0.0875006 60: 0.0677065 61: -0.111693 63: -0.0791442 64: 0.0984083 97: 0.0831084 98: -0.0824026 99: 0.0892412 100: -0.0466491 101: -0.0934947 102: -0.103974 103: -0.0548927 106: 0.0736245 121: 0.0726606
Row 25: 9: -0.480773 10: -0.832083 24: -0.471228 25: 3.56549 26: -0.53029 34: -0.379179 36: -0.871933 62: 0.10275 63: 0.0929375 64: -0.115559 66: -0.0779648 67: 0.113895 101: -0.122301 103: 0.107902 104: -0.110436 105: -0.113567 106: -0.0544203 108: 0.0849225 125: 0.0918406
Row 26: 10: -0.400373 11: -1.0854 12: -0.248259 25: -0.53029 26: 3.70321 27: -0.7952 34: -0.643685 65: 0.10261 66: 0.0778583 67: -0.113739 68: 0.147354 69: -0.0690633 71: -0.196926 72: 0.0909485 104: -0.0945973 105: 0.105121 107: -0.0517256 108: -0.0911502 110: 0.0933104
Row 27: 12: -0.676685 13: -0.813095 26: -0.7952 27: 3.69033 28: -0.707702 34: -0.697653 70: 0.101777 71: 0.157507 72: -0.146503 74: -0.105791 75: 0.139531 107: -0.122692 108: 0.0977186 109: -0.140104 110: -0.0999664 112: 0.118523
Row 28: 13: -0.307486 14: -0.343721 27: -0.707702 28: 3.52468 29: -0.895438 34: -0.473338 35: -0.287796 37: -0.509201 73: 0.0766612 74: 0.0796851 75: -0.105099 77: -0.101966 78: 0.082592 109: -0.0443342 110: 0.0825994 111: -0.0325616 112: -0.0842096 113: -0.137443 114: -0.0818333 116: 0.0992092 126: 0.0804999 127: 0.0862002
Row 29: 14: -0.771353 15: -0.52416 28: -0.895438 29: 3.72354 30: -0.968571 35: -0.564014 76: 0.102623 77: 0.136498 78: -0.110563 79: -0.13099 80: 0.115727 111: -0.127439 113: 0.14018 115: -0.0798592 116: -0.171739 117: 0.125561
Row 30: 3: -0.123368 15: -0.573581 16: -1.33137 17: -0.336404 29: -0.968571 30: 3.93186 35: -0.598564 45: 0.0832196 46: 0.110948 47: -0.173606 79: 0.106223 80: -0.0938455 81: 0.110948 82: 0 84: -0.181455 86: 0.126574 115: -0.0621923 116: 0.117398 117: -0.144212
Row 31: 17: -0.30274 18: -1.13758 19: -0.239727 20: -0.660832 31: 3.75805 33: -0.767936 35: -0.649229 83: 0.0947986 85: -0.159723 86: 0.115381 87: 0.0947986 88: 0 89: 0.110522 90: -0.165366 92: -0.0968204 93: 0.0964375 118: -0.0828523 119: -0.12158 123: 0.114404
Row 32: 21: -0.707916 22: -0.613124 24: -0.347567 32: 3.57504 33: -0.609107 36: -0.927505 37: -0.369818 94: 0.0996696 95: -0.0955243 96: 0.113841 98: 0.0912736 99: -0.0887558 102: -0.136374 103: 0.103028 120: -0.104609 121: -0.0441856 122: -0.126391 124: 0.0922859 128: 0.0957414
Row 33: 20: -0.587027 21: -0.454863 31: -0.767936 32: -0.609107 33: 3.5243 35: -0.372794 37: -0.732577 91: 0.092073 92: 0.105088 93: -0.099323 95: 0.0829234 96: -0.0991859 118: -0.0824973 119: 0.105399 120: -0.0854633 122: 0.104058 123: -0.14111 124: -0.0798048 127: 0.0978433
Row 34: 25: -0.379179 26: -0.643685 27: -0.697653 28: -0.473338 34: 3.53245 36: -0.710183 37: -0.628416 104: 0.100989 105: -0.119067 106: 0.0812748 107: 0.0977827 108: -0.0914909 109: 0.0944361 110: -0.0759434 112: -0.132415 114: 0.116869 125: -0.0603016 126: -0.109524 128: 0.0973906
Row 35: 17: -0.332973 28: -0.287796 29: -0.564014 30: -0.598564 31: -0.649229 33: -0.372794 35: 3.47302 37: -0.667649 84: 0.0943089 85: 0.0971555 86: -0.135969 111: 0.0725036 113: -0.117127 114: 0.0925899 115: 0.0663671 116: -0.0448683 117: -0.0609153 118: 0.0767838 119: -0.0657345 123: -0.0937792 124: 0.0791278 127: -0.0604428
Row 36: 24: -0.548354 25: -0.871933 32: -0.927505 34: -0.710183 36: 3.68457 37: -0.626599 101: 0.118288 102: 0.129142 103: -0.156037 105: 0.127513 106: -0.100479 121: -0.105576 122: 0.131018 125: -0.117284 126: 0.108134 128: -0.134719
Row 37: 28: -0.509201 32: -0.369818 33: -0.732577 34: -0.628416 35: -0.667649 36: -0.626599 37: 3.53426 112: 0.0981016 113: 0.114391 114: -0.127626 120: 0.0932201 121: 0.0771011 122: -0.108685 123: 0.120485 124: -0.0916088 125: 0.0857449 126: -0.0791105 127: -0.123601 128: -0.0584127
Row 38: 0: -0.0833333 4: 0 19: 0.0833333 38: 0.0416667 39: 0 49: -0.0208333
Row 39: 0: -0.0833333 4: 0.0833333 19: 0 38: 0 39: 0.0416667 49: -0.0208333
Row 40: 1: -0.034013 7: -0.0837431 23: 0.117756 40: 0.0381352 42: -0.0209358 58: -0.00850325
Row 41: 1: -0.0342543 8: -0.0829256 23: 0.11718 41: 0.0380505 42: -0.0207314 60: -0.00856358
Row 42: 1: -0.0698066 7: 0.118528 8: 0.117948 23: -0.166669 40: -0.0209358 41: -0.0207314 42: 0.0761857 58: -0.00869614 60: -0.00875552
Row 43: 2: -0.0833333 11: 0 12: 0.0833333 43: 0.0416667 44: 1.73472e-18 68: -0.0208333
Row 44: 2: -0.0833333 11: 0.0833333 12: 0 43: 1.73472e-18 44: 0.0416667 68: -0.0208333
Row 45: 3: -0.020713 15: -0.0625066 30: 0.0832196 45: 0.0377773 47: -0.0156267 80: -0.00517824
Row 46: 3: 0.000151649 16: -0.111099 30: 0.110948 46: 0.0434229 47: -0.0277748 82: 3.79121e-05
Row 47: 3: -0.130634 15: 0.130396 16: 0.173843 30: -0.173606 45: -0.0156267 46: -0.0277748 47: 0.0812002 80: -0.0169724 82: -0.015686
Row 48: 4: -0.0553087 5: -0.0302454 20: 0.0855542 48: 0.0367928 50: -0.00756135 52: -0.0138272
Row 49: 0: 0.166667 4: -0.169926 19: -0.122284 20: 0.125543 38: -0.0208333 39: -0.0208333 49: 0.0801647 50: -0.00973772 89: -0.021648
Row 50: 4: -0.0900664 5: 0.0918626 19: 0.0674001 20: -0.0691963 48: -0.00756135 49: -0.00973772 50: 0.0752909 52: -0.0154043 89: -0.0071123
Row 51: 5: -0.0537908 6: -0.0387141 21: 0.0925049 51: 0.036266 53: -0.00967852 55: -0.0134477
Row 52: 4: 0.116926 5: -0.12557 20: -0.106394 21: 0.115038 48: -0.0138272 50: -0.0154043 52: 0.0735441 53: -0.0127714 91: -0.0159881
Row 53: 5: -0.0845258 6: 0.0912733 20: 0.0830523 21: -0.0897997 51: -0.00967852 52: -0.0127714 53: 0.0730172 55: -0.0131398 91: -0.00799166
Row 54: 6: -0.0577312 7: -0.0420473 22: 0.0997784 54: 0.0362622 56: -0.0105118 57: -0.0144328
Row 55: 5: 0.10635 6: -0.114421 21: -0.102651 22: 0.110723 51: -0.0134477 53: -0.0131398 55: 0.0728018 56: -0.0122151 94: -0.0154655
Row 56: 6: -0.0806907 7: 0.0873176 21: 0.0842808 22: -0.0909077 54: -0.0105118 55: -0.0122151 56: 0.0727979 57: -0.0113176 94: -0.00885509
Row 57: 6: 0.103002 7: -0.123185 22: -0.112022 23: 0.132205 54: -0.0144328 56: -0.0113176 57: 0.0744465 58: -0.0135726 97: -0.0194787
Row 58: 1: 0.0687976 7: -0.0553164 22: 0.0748223 23: -0.0883035 40: -0.00850325 42: -0.00869614 57: -0.0135726 58: 0.0763194 97: -0.00513297
Row 59: 8: -0.0322254 9: -0.0552752 24: 0.0875006 59: 0.036627 61: -0.0138188 63: -0.00805635
Row 60: 1: 0.0692764 8: -0.0463111 23: -0.0906718 24: 0.0677065 41: -0.00856358 42: -0.00875552 60: 0.0782672 61: -0.0141044 100: -0.00282225
Row 61: 8: -0.152168 9: 0.114283 23: 0.149578 24: -0.111693 59: -0.0138188 60: -0.0141044 61: 0.0768437 63: -0.0147518 100: -0.0232901
Row 62: 9: -0.0341102 10: -0.0686401 25: 0.10275 62: 0.0368873 64: -0.01716 66: -0.00852755
Row 63: 8: 0.0912327 9: -0.105026 24: -0.0791442 25: 0.0929375 59: -0.00805635 61: -0.0147518 63: 0.0727337 64: -0.0117297 101: -0.0115047
Row 64: 9: -0.0962885 10: 0.113439 24: 0.0984083 25: -0.115559 62: -0.01716 63: -0.0117297 64: 0.0729941 66: -0.0111998 101: -0.0128724
Row 65: 10: -0.0327251 11: -0.0698849 26: 0.10261 65: 0.037 67: -0.0174712 69: -0.00818128
Row 66: 9: 0.0789093 10: -0.0788029 25: -0.0779648 26: 0.0778583 62: -0.00852755 64: -0.0111998 66: 0.073862 67: -0.0109636 104: -0.00850094
Row 67: 10: -0.11543 11: 0.115275 25: 0.113895 26: -0.113739 65: -0.0174712 66: -0.0109636 67: 0.0739746 69: -0.0113475 104: -0.0175101
Row 68: 2: 0.166667 11: -0.194349 12: -0.119671 26: 0.147354 43: -0.0208333 44: -0.0208333 68: 0.0834428 69: -0.00908454 71: -0.0277539
Row 69: 10: 0.078115 11: -0.0651405 12: 0.0560888 26: -0.0690633 65: -0.00818128 67: -0.0113475 68: -0.00908454 69: 0.0787761 71: -0.00493766
Row 70: 12: -0.0411842 13: -0.0605924 27: 0.101777 70: 0.0363725 72: -0.0151481 74: -0.010296
Row 71: 11: 0.130766 12: -0.0913473 26: -0.196926 27: 0.157507 68: -0.0277539 69: -0.00493766 71: 0.0824124 72: -0.0214775 107: -0.0178992
Row 72: 12: -0.0487517 13: 0.104306 26: 0.0909485 27: -0.146503 70: -0.0151481 71: -0.0214775 72: 0.0770088 74: -0.0109283 107: -0.00125957
Row 73: 13: -0.0361696 14: -0.0404916 28: 0.0766612 73: 0.0370358 75: -0.0101229 77: -0.00904239
Row 74: 12: 0.0848976 13: -0.0587915 27: -0.105791 28: 0.0796851 70: -0.010296 72: -0.0109283 74: 0.0750247 75: -0.0161518 109: -0.00376953
Row 75: 13: -0.146405 14: 0.111973 27: 0.139531 28: -0.105099 73: -0.0101229 74: -0.0161518 75: 0.0756879 77: -0.0178705 109: -0.0187309
Row 76: 14: -0.0578574 15: -0.0447659 29: 0.102623 76: 0.0362636 78: -0.0111915 79: -0.0144643
Row 77: 13: 0.107651 14: -0.142183 28: -0.101966 29: 0.136498 73: -0.00904239 75: -0.0178705 77: 0.0753591 78: -0.0164492 111: -0.0176754
Row 78: 14: -0.0592263 15: 0.0871969 28: 0.082592 29: -0.110563 76: -0.0111915 77: -0.0164492 78: 0.0745869 79: -0.0106078 111: -0.00419882
Row 79: 14: 0.100288 15: -0.0755213 29: -0.13099 30: 0.106223 76: -0.0144643 78: -0.0106078 79: 0.0734678 80: -0.0182831 115: -0.00827257
Row 80: 3: 0.0886026 15: -0.110484 29: 0.115727 30: -0.0938455 45: -0.00517824 47: -0.0169724 79: -0.0182831 80: 0.0749815 115: -0.0106485
Row 81: 16: -0.110796 17: -0.000151649 30: 0.110948 81: 0.0433471 82: -3.79121e-05 84: -0.027699
Row 82: 3: 0.0625923 16: -0.125185 17: 0.0625923 30: 0 46: 3.79121e-05 47: -0.015686 81: -3.79121e-05 82: 0.08677 84: -0.0156102
Row 83: 17: -0.000386213 18: -0.0944124 31: 0.0947986 83: 0.0419172 85: -0.0236031 88: -9.65532e-05
Row 84: 16: 0.173237 17: -0.0860908 30: -0.181455 35: 0.0943089 81: -0.027699 82: -0.0156102 84: 0.0809032 86: -0.0176647 117: -0.00591256
Row 85: 17: -0.104715 18: 0.167282 31: -0.159723 35: 0.0971555 83: -0.0236031 85: 0.0787237 86: -0.0163275 88: -0.0182175 119: -0.00796132
Row 86: 17: -0.105986 30: 0.126574 31: 0.115381 35: -0.135969 84: -0.0176647 85: -0.0163275 86: 0.0743626 117: -0.0139789 119: -0.0125176
Row 87: 18: -0.0951848 19: 0.000386213 31: 0.0947986 87: 0.0421103 88: 9.65532e-05 90: -0.0237962
Row 88: 17: 0.0732563 18: -0.146513 19: 0.0732563 31: 0 83: -9.65532e-05 85: -0.0182175 87: 9.65532e-05 88: 0.0840274 90: -0.0184106
Row 89: 4: 0.115041 19: -0.0687899 20: -0.156773 31: 0.110522 49: -0.021648 50: -0.0071123 89: 0.0754328 90: -0.0175452 92: -0.0100852
Row 90: 18: 0.168827 19: -0.11086 20: 0.107398 31: -0.165366 87: -0.0237962 88: -0.0184106 89: -0.0175452 90: 0.079045 92: -0.00930429
Row 91: 5: 0.0959192 20: -0.0811964 21: -0.106796 33: 0.092073 52: -0.0159881 53: -0.00799166 91: 0.0728984 93: -0.0107108 96: -0.0123074
Row 92: 19: 0.0775578 20: -0.0858253 31: -0.0968204 33: 0.105088 89: -0.0100852 90: -0.00930429 92: 0.0731961 93: -0.0141199 118: -0.012152
Row 93: 20: -0.0924734 21: 0.095359 31: 0.0964375 33: -0.099323 91: -0.0107108 92: -0.0141199 93: 0.0724086 96: -0.0131289 118: -0.00998942
Row 94: 6: 0.0972825 21: -0.0895218 22: -0.10743 32: 0.0996696 55: -0.0154655 56: -0.00885509 94: 0.0726881 95: -0.011392 99: -0.0135254
Row 95: 21: -0.0779072 22: 0.090508 32: -0.0955243 33: 0.0829234 94: -0.011392 95: 0.0728544 96: -0.012489 99: -0.011235 120: -0.00824182
Row 96: 20: 0.101746 21: -0.1164 32: 0.113841 33: -0.0991859 91: -0.0123074 93: -0.0131289 95: -0.012489 96: 0.0728492 120: -0.0159711
Row 97: 7: 0.0984466 22: -0.055892 23: -0.125663 24: 0.0831084 57: -0.0194787 58: -0.00513297 97: 0.0747722 98: -0.0119371 100: -0.00884002
Row 98: 22: -0.119863 23: 0.110992 24: -0.0824026 32: 0.0912736 97: -0.0119371 98: 0.0730531 99: -0.00866359 100: -0.0158109 102: -0.0141548
Row 99: 21: 0.0990414 22: -0.0995268 24: 0.0892412 32: -0.0887558 94: -0.0135254 95: -0.011235 98: -0.00866359 99: 0.0726175 102: -0.0136467
Row 100: 8: 0.104449 22: 0.0986036 23: -0.156404 24: -0.0466491 60: -0.00282225 61: -0.0232901 97: -0.00884002 98: -0.0158109 100: 0.0768046
Row 101: 9: 0.0975081 24: -0.0934947 25: -0.122301 36: 0.118288 63: -0.0115047 64: -0.0128724 101: 0.0735835 103: -0.0190706 106: -0.0105013
Row 102: 22: 0.111206 24: -0.103974 32: -0.136374 36: 0.129142 98: -0.0141548 99: -0.0136467 102: 0.074569 103: -0.0199387 121: -0.0123468
Row 103: 24: -0.0548927 25: 0.107902 32: 0.103028 36: -0.156037 101: -0.0190706 102: -0.0199387 103: 0.0755806 106: -0.00790481 121: -0.00581837
Row 104: 10: 0.104044 25: -0.110436 26: -0.0945973 34: 0.100989 66: -0.00850094 67: -0.0175101 104: 0.0733537 105: -0.0151484 108: -0.0100989
Row 105: 25: -0.113567 26: 0.105121 34: -0.119067 36: 0.127513 104: -0.0151484 105: 0.0739576 106: -0.0146185 108: -0.0111318 125: -0.0172599
Row 106: 24: 0.0736245 25: -0.0544203 34: 0.0812748 36: -0.100479 101: -0.0105013 103: -0.00790481 105: -0.0146185 106: 0.0750554 125: -0.00570025
Row 107: 12: 0.0766349 26: -0.0517256 27: -0.122692 34: 0.0977827 71: -0.0178992 72: -0.00125957 107: 0.0767377 108: -0.0127739 110: -0.0116718
Row 108: 25: 0.0849225 26: -0.0911502 27: 0.0977186 34: -0.0914909 104: -0.0100989 105: -0.0111318 107: -0.0127739 108: 0.0724805 110: -0.0116558
Row 109: 13: 0.0900017 27: -0.140104 28: -0.0443342 34: 0.0944361 74: -0.00376953 75: -0.0187309 109: 0.075597 110: -0.00731402 112: -0.016295
Row 110: 26: 0.0933104 27: -0.0999664 28: 0.0825994 34: -0.0759434 107: -0.0116718 108: -0.0116558 109: -0.00731402 110: 0.0730463 112: -0.0133358
Row 111: 14: 0.0874967 28: -0.0325616 29: -0.127439 35: 0.0725036 77: -0.0176754 78: -0.00419882 111: 0.07731 113: -0.0141843 116: -0.00394157
Row 112: 27: 0.118523 28: -0.0842096 34: -0.132415 37: 0.0981016 109: -0.016295 110: -0.0133358 112: 0.0738787 114: -0.0168088 126: -0.00771657
Row 113: 28: -0.137443 29: 0.14018 35: -0.117127 37: 0.114391 111: -0.0141843 113: 0.0756342 114: -0.0150976 116: -0.0208607 127: -0.0135001
Row 114: 28: -0.0818333 34: 0.116869 35: 0.0925899 37: -0.127626 112: -0.0168088 113: -0.0150976 114: 0.0735814 126: -0.0124084 127: -0.00804992
Row 115: 15: 0.0756844 29: -0.0798592 30: -0.0621923 35: 0.0663671 79: -0.00827257 80: -0.0106485 115: 0.07587 116: -0.00727552 117: -0.00931627
Row 116: 28: 0.0992092 29: -0.171739 30: 0.117398 35: -0.0448683 111: -0.00394157 113: -0.0208607 115: -0.00727552 116: 0.0776524 117: -0.022074
Row 117: 17: 0.0795659 29: 0.125561 30: -0.144212 35: -0.0609153 84: -0.00591256 86: -0.0139789 115: -0.00931627 116: -0.022074 117: 0.0762219
Row 118: 20: 0.0885658 31: -0.0828523 33: -0.0824973 35: 0.0767838 92: -0.012152 93: -0.00998942 118: 0.0733347 119: -0.0084723 123: -0.0107237
Row 119: 17: 0.0819158 31: -0.12158 33: 0.105399 35: -0.0657345 85: -0.00796132 86: -0.0125176 118: -0.0084723 119: 0.0738798 123: -0.0178774
Row 120: 21: 0.0968518 32: -0.104609 33: -0.0854633 37: 0.0932201 95: -0.00824182 96: -0.0159711 120: 0.0728975 122: -0.013124 124: -0.010181
Row 121: 24: 0.0726606 32: -0.0441856 36: -0.105576 37: 0.0771011 102: -0.0123468 103: -0.00581837 121: 0.0760864 122: -0.0140472 128: -0.00522804
Row 122: 32: -0.126391 33: 0.104058 36: 0.131018 37: -0.108685 120: -0.013124 121: -0.0140472 122: 0.074178 124: -0.0128904 128: -0.0187073
Row 123: 31: 0.114404 33: -0.14111 35: -0.0937792 37: 0.120485 118: -0.0107237 119: -0.0178774 123: 0.0742553 124: -0.0127212 127: -0.0174
Row 124: 32: 0.0922859 33: -0.0798048 35: 0.0791278 37: -0.0916088 120: -0.010181 122: -0.0128904 123: -0.0127212 124: 0.0733774 127: -0.00706078
Row 125: 25: 0.0918406 34: -0.0603016 36: -0.117284 37: 0.0857449 105: -0.0172599 106: -0.00570025 125: 0.0739873 126: -0.0120611 128: -0.00937515
Row 126: 28: 0.0804999 34: -0.109524 36: 0.108134 37: -0.0791105 112: -0.00771657 114: -0.0124084 125: -0.0120611 126: 0.0733425 128: -0.0149725
Row 127: 28: 0.0862002 33: 0.0978433 35: -0.0604428 37: -0.123601 113: -0.0135001 114: -0.00804992 123: -0.0174 124: -0.00706078 127: 0.0738296
Row 128: 32: 0.0957414 34: 0.0973906 36: -0.134719 37: -0.0584127 121: -0.00522804 122: -0.0187073 125: -0.00937515 126: -0.0149725 128: 0.0743913
6. Solve the system:¶
[10]:
gfu.vec.data = \
a.mat.Inverse(freedofs=fes.FreeDofs()) * f.vec
Draw(gfu)
[10]:
BaseWebGuiScene
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.0923051
0
0
0
0
0
0
0
0
0.0578976
0.0863397
0.095413
0.0944941
0.088859
0.0779511
0.059639
0.0330669
0.0361131
0.0364559
0.0317036
0.0185143
0.0384776
0.0431703
0.0476901
0.0763835
0.0896404
0.0935756
0.0902834
0.0619423
0.05721
0.0634411
0.072785
0.0832449
0.0629249
0.0775277
0
-0.00576212
0
0
0.0208108
0
-0.0350834
0.00265432
-0.00721308
-0.00488766
0
-0.00992425
-0.00630809
0
-0.0163826
-0.0115972
0
-0.0146117
-0.019397
-0.00629918
-0.021289
0
-0.0255422
-0.00978494
0
-0.0385158
-0.0136376
0
-0.026807
-0.018579
-0.0333669
-0.024961
-0.0242793
-0.00449407
-0.0124615
-0.0145086
-0.00594865
-0.0117808
-0.00555001
-0.00890522
-0.00588952
-0.00539555
-0.00706629
-0.00787172
0.0025744
-0.00832323
-0.00587822
-0.00475479
0.00342728
-0.00825162
0.00390131
0.00421041
-0.00475943
-0.00353907
-0.00993744
-0.009953
-0.00708995
-0.0135307
-0.0114551
-0.0078331
-0.0254921
-0.0071086
-0.00325064
-0.0089176
-0.0210183
-0.00522333
-0.0030572
-0.0144307
-0.0182813
-0.0155115
-0.0215164
-0.0172692
-0.00414894
-0.00783019
-0.0132532
-0.016378
-0.00903592
-0.00194529
-0.00888682
-0.00574139
-0.00227548
-0.0105216
-0.00871044
-0.0118078
-0.00808258
-0.0113092
-0.0118318
-0.00416845
-0.0163461
-0.00951852
-0.0109196
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
).
[ ]: