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]:
(37, 52)
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]:
125
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.
| autoupdate: bool = False
| Automatically update on a change to the mesh.
| 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, **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.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.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:
|
| autoupdate
|
| 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__
|
| flags
|
| ----------------------------------------------------------------------
| 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 0x7f31b46890b0>
You can examine the linear system in more detail:
[8]:
print(f.vec)
0.000333333
0.00633333
0.00633333
0.000800951
0.00388508
0.0076445
0.0113284
0.0178159
0.0278275
0.0173112
0.0154685
0.0181138
0.0204443
0.00955428
0.00775661
0.0042821
0.000386462
0.00215946
0.000477049
0.00259004
0.0118211
0.0192953
0.0307721
0.0302149
0.0269091
0.0285724
0.0173381
0.0183376
0.013492
0.00550585
0.00686885
0.0217628
0.0157415
0.0211763
0.0113116
0.0174036
0.0226306
-6.66667e-05
-3.33333e-05
-0.000766667
-0.0008
-0.0008
-0.000766667
-9.52746e-05
-1.93231e-05
-0.000125688
-0.000238492
-0.000214885
-0.000425538
-0.000384163
-0.000563154
-0.000717295
-0.000526134
-0.000884261
-0.00103647
-0.0017373
-0.0014478
-0.000807656
-0.00181551
-0.00165701
-0.000673671
-0.00145372
-0.00133575
-0.000648933
-0.00126099
-0.00123885
-0.00157569
-0.0014435
-0.000443027
-0.00135145
-0.00100443
-0.000389019
-0.00079759
-0.00075629
-0.000261064
-0.000708227
-0.000587395
-0.000442862
-0.000285407
-1.93231e-05
-7.72924e-05
-2.38524e-05
-0.000154323
-0.000172093
-0.000278246
-2.38524e-05
-9.54097e-05
-0.00032181
-0.000183132
-0.000663312
-0.00040581
-0.000559774
-0.000983527
-0.000855246
-0.000718479
-0.00139785
-0.00103282
-0.00113578
-0.00100047
-0.000920078
-0.00114204
-0.000959829
-0.000896019
-0.000996606
-0.000969011
-0.00072904
-0.000808183
-0.000701041
-0.000840338
-0.000845828
-0.000375439
-0.000528432
-0.000739237
-0.000345317
-0.000459317
-0.000387811
-0.000819442
-0.000864813
-0.000865653
-0.000578816
-0.000797631
-0.000841733
-0.000870701
-0.000703596
-0.000824523
[9]:
print(a.mat)
Row 0: 0: 1 4: -0.5 19: -0.5 37: -0.0833333 38: -0.0833333 47: 0.166667
Row 1: 1: 1 7: -0.5 8: -0.5 39: -0.0833333 40: -0.0833333 55: 0.166667
Row 2: 2: 1 11: -0.5 12: -0.5 41: -0.0833333 42: -0.0833333 66: 0.166667
Row 3: 3: 0.92654 15: -0.429776 16: -0.404531 29: -0.0922334 43: -0.0192115 44: 0.00383926 45: -0.139051 78: 0.0908408 80: 0.0635826
Row 4: 0: -0.5 4: 1.89214 5: -0.369141 19: -0.170235 20: -0.852766 37: 0 38: 0.0833333 46: -0.0554044 47: -0.170057 48: -0.089896 50: 0.116928 87: 0.115096
Row 5: 4: -0.369141 5: 1.75902 6: -0.330249 20: -0.393103 21: -0.666525 46: -0.0302389 48: 0.0917624 49: -0.0503815 50: -0.12223 51: -0.0903197 53: 0.105423 89: 0.0959843
Row 6: 5: -0.330249 6: 1.7383 7: -0.312052 21: -0.525655 22: -0.570346 49: -0.039568 51: 0.0946095 52: -0.0418649 53: -0.108234 54: -0.10005 56: 0.0938736 92: 0.101234
Row 7: 1: -0.5 6: -0.312052 7: 1.90043 8: -0.209289 22: -0.879085 39: 0 40: 0.0833333 52: -0.0507822 54: 0.102791 55: -0.179065 56: -0.0868901 58: 0.130613
Row 8: 1: -0.5 7: -0.209289 8: 1.79706 9: -0.419471 22: -0.261905 23: -0.406398 39: 0.0833333 40: 0 55: -0.110935 56: 0.0624833 57: -0.0183232 58: -0.0842912 59: -0.0859608 61: 0.088235 95: 0.0654588
Row 9: 8: -0.419471 9: 1.77344 10: -0.299906 23: -0.681965 24: -0.372103 57: -0.0641858 59: 0.134098 60: -0.0270683 61: -0.104861 62: -0.0994593 64: 0.0770526 97: 0.0844239
Row 10: 9: -0.299906 10: 1.78385 11: -0.274651 24: -0.849309 25: -0.359984 60: -0.0725667 62: 0.122551 63: -0.0299194 64: -0.0800624 65: -0.11476 67: 0.0756945 100: 0.0990629
Row 11: 2: -0.5 10: -0.274651 11: 1.99359 12: -0.0983666 25: -1.12058 41: 0 42: 0.0833333 63: -0.0736497 65: 0.119425 66: -0.196446 67: -0.0621696 69: 0.129507
Row 12: 2: -0.5 11: -0.0983666 12: 1.83136 13: -0.218895 25: -0.244791 26: -0.769308 41: 0.0833333 42: 0 66: -0.122636 67: 0.0556974 68: -0.0486254 69: -0.0959869 70: -0.037978 72: 0.085108 103: 0.081088
Row 13: 12: -0.218895 13: 1.75747 14: -0.333769 26: -0.704669 27: -0.500134 68: -0.0607518 70: 0.0972344 71: -0.0542939 72: -0.0655443 73: -0.112321 75: 0.109922 105: 0.0857548
Row 14: 13: -0.333769 14: 1.75676 15: -0.403014 27: -0.429757 28: -0.590218 71: -0.0356996 73: 0.0913278 74: -0.0553925 75: -0.0986054 76: -0.103096 77: 0.122561 107: 0.0789038
Row 15: 3: -0.429776 14: -0.403014 15: 1.80299 28: -0.281759 29: -0.688445 43: -0.0612978 45: 0.132927 74: -0.0263034 76: 0.0934724 77: -0.120612 78: -0.0922857 110: 0.0740994
Row 16: 3: -0.404531 16: 2.0764 17: -0.35846 29: -1.31341 44: -0.11329 45: 0.180712 79: -0.105612 80: -0.127165 82: 0.165355
Row 17: 16: -0.35846 17: 1.76572 18: -0.441688 29: -0.2915 30: -0.172197 34: -0.501877 79: -0.00383926 80: 0.0635826 81: 0.0030286 82: -0.0992307 83: -0.117773 84: -0.0764721 86: 0.0705861 113: 0.0842313 115: 0.0758868
Row 18: 17: -0.441688 18: 2.02918 19: -0.405345 30: -1.18215 81: -0.101541 83: 0.175156 85: -0.095484 86: -0.141172 88: 0.163041
Row 19: 0: -0.5 4: -0.170235 18: -0.405345 19: 1.79579 20: -0.434821 30: -0.285392 37: 0.0833333 38: 0 47: -0.122291 48: 0.0673304 85: -0.0030286 86: 0.0705861 87: -0.0729092 88: -0.10107 90: 0.078049
Row 20: 4: -0.852766 5: -0.393103 19: -0.434821 20: 3.55048 21: -0.611678 30: -0.667364 32: -0.590751 46: 0.0856433 47: 0.125681 48: -0.0691968 50: -0.105442 51: 0.085316 87: -0.156576 88: 0.103365 89: -0.0870282 90: -0.0802209 91: -0.0932833 94: 0.103658 114: 0.0880834
Row 21: 5: -0.666525 6: -0.525655 20: -0.611678 21: 3.49993 22: -0.49992 31: -0.763067 32: -0.433089 49: 0.0899494 50: 0.110744 51: -0.0896057 53: -0.0937366 54: 0.0913963 89: -0.101785 91: 0.0929875 92: -0.111805 93: -0.0710673 94: -0.115323 96: 0.103729 116: 0.0945164
Row 22: 6: -0.570346 7: -0.879085 8: -0.261905 21: -0.49992 22: 3.64558 23: -0.933015 31: -0.501309 52: 0.0926471 53: 0.0965478 54: -0.0941373 55: 0.123334 56: -0.0694667 58: -0.189707 59: 0.110024 92: -0.0955738 93: 0.0823459 95: -0.0572194 96: -0.101493 98: 0.102698
Row 23: 8: -0.406398 9: -0.681965 22: -0.933015 23: 3.59109 24: -0.614411 31: -0.407342 35: -0.547958 57: 0.082509 58: 0.143384 59: -0.158161 61: -0.080099 62: 0.111251 95: -0.092364 96: 0.104482 97: -0.0898208 98: -0.112508 99: -0.0655623 102: 0.0809719 117: 0.0759167
Row 24: 9: -0.372103 10: -0.849309 23: -0.614411 24: 3.56414 25: -0.584965 33: -0.356721 35: -0.786629 60: 0.099635 61: 0.0967247 62: -0.134342 64: -0.0762241 65: 0.118141 97: -0.10047 99: 0.106147 100: -0.106083 101: -0.113922 102: -0.0629812 104: 0.0854369 121: 0.0879389
Row 25: 10: -0.359984 11: -1.12058 12: -0.244791 24: -0.584965 25: 3.72246 26: -0.780094 33: -0.632047 63: 0.103569 64: 0.0792338 65: -0.122805 66: 0.152416 67: -0.0692223 69: -0.197286 70: 0.0856684 100: -0.09037 101: 0.10863 103: -0.0465446 104: -0.0941812 106: 0.0908919
Row 26: 12: -0.769308 13: -0.704669 25: -0.780094 26: 3.67759 27: -0.749926 33: -0.673588 68: 0.109377 69: 0.163765 70: -0.144925 72: -0.110393 73: 0.11846 103: -0.133275 104: 0.0995248 105: -0.115276 106: -0.109063 108: 0.121803
Row 27: 13: -0.500134 14: -0.429757 26: -0.749926 27: 3.57447 28: -0.916705 33: -0.508712 36: -0.469238 71: 0.0899935 72: 0.0908291 73: -0.097467 75: -0.122737 76: 0.10437 105: -0.055669 106: 0.0898276 107: -0.0812277 108: -0.0961255 109: -0.142519 112: 0.129642 122: 0.0910832
Row 28: 14: -0.590218 15: -0.281759 27: -0.916705 28: 3.63466 29: -0.913179 34: -0.598184 36: -0.334611 74: 0.0816959 75: 0.11142 76: -0.0947465 77: -0.134212 78: 0.0994762 107: -0.067072 109: 0.108436 110: -0.0602816 111: -0.10505 112: -0.144413 113: 0.113002 123: 0.0917459
Row 29: 3: -0.0922334 15: -0.688445 16: -1.31341 17: -0.2915 28: -0.913179 29: 3.88193 34: -0.583158 43: 0.0805093 44: 0.109451 45: -0.174588 77: 0.132263 78: -0.0980313 79: 0.109451 80: 0 82: -0.167081 84: 0.106213 110: -0.089167 111: 0.109101 113: -0.118121
Row 30: 17: -0.172197 18: -1.18215 19: -0.285392 20: -0.667364 30: 3.78248 32: -0.710649 34: -0.764723 81: 0.0985126 83: -0.174589 84: 0.104776 85: 0.0985126 86: 1.38778e-17 87: 0.114389 88: -0.165337 90: -0.101436 91: 0.0982741 114: -0.0957807 115: -0.0932702 119: 0.115948
Row 31: 21: -0.763067 22: -0.501309 23: -0.407342 31: 3.63771 32: -0.715588 35: -1.00537 36: -0.245036 92: 0.106145 93: -0.0949869 94: 0.11602 95: 0.0841246 96: -0.106718 98: -0.115899 99: 0.0996652 116: -0.0929254 117: -0.0362641 118: -0.159491 120: 0.0961702 124: 0.10416
Row 32: 20: -0.590751 21: -0.433089 30: -0.710649 31: -0.715588 32: 3.50593 34: -0.429218 36: -0.62664 89: 0.0928289 90: 0.103608 91: -0.0979783 93: 0.0837082 94: -0.104356 114: -0.0777241 115: 0.0925576 116: -0.082857 118: 0.118413 119: -0.114227 120: -0.10718 123: 0.0932062
Row 33: 24: -0.356721 25: -0.632047 26: -0.673588 27: -0.508712 33: 3.54391 35: -0.812646 36: -0.5602 100: 0.0973905 101: -0.122589 102: 0.0846519 103: 0.0987313 104: -0.0907806 105: 0.0851898 106: -0.0716564 108: -0.113808 109: 0.113403 121: -0.0604965 122: -0.131322 124: 0.111286
Row 34: 17: -0.501877 28: -0.598184 29: -0.583158 30: -0.764723 32: -0.429218 34: 3.51649 36: -0.639329 82: 0.100956 83: 0.117207 84: -0.134517 110: 0.0753492 111: -0.0906772 112: 0.115025 113: -0.0791125 114: 0.0854215 115: -0.0751743 119: -0.106008 120: 0.0921225 123: -0.100593
Row 35: 23: -0.547958 24: -0.786629 31: -1.00537 33: -0.812646 35: 3.71605 36: -0.563445 97: 0.105867 98: 0.12571 99: -0.14025 101: 0.127881 102: -0.102643 117: -0.107555 118: 0.149407 121: -0.102917 122: 0.110477 124: -0.165977
Row 36: 27: -0.469238 28: -0.334611 31: -0.245036 32: -0.62664 33: -0.5602 34: -0.639329 35: -0.563445 36: 3.4385 107: 0.0693959 108: 0.0881303 109: -0.0793199 111: 0.0866269 112: -0.100254 116: 0.081266 117: 0.0679026 118: -0.108329 119: 0.104287 120: -0.081113 121: 0.0754744 122: -0.0702381 123: -0.084359 124: -0.0494695
Row 37: 0: -0.0833333 4: 0 19: 0.0833333 37: 0.0416667 38: 0 47: -0.0208333
Row 38: 0: -0.0833333 4: 0.0833333 19: 0 37: 0 38: 0.0416667 47: -0.0208333
Row 39: 1: -0.0833333 7: 0 8: 0.0833333 39: 0.0416667 40: 0 55: -0.0208333
Row 40: 1: -0.0833333 7: 0.0833333 8: 0 39: 0 40: 0.0416667 55: -0.0208333
Row 41: 2: -0.0833333 11: 0 12: 0.0833333 41: 0.0416667 42: 1.73472e-18 66: -0.0208333
Row 42: 2: -0.0833333 11: 0.0833333 12: 0 41: 1.73472e-18 42: 0.0416667 66: -0.0208333
Row 43: 3: -0.0192115 15: -0.0612978 29: 0.0805093 43: 0.0380346 45: -0.0153244 78: -0.00480287
Row 44: 3: 0.00383926 16: -0.11329 29: 0.109451 44: 0.0442182 45: -0.0283226 80: 0.000959815
Row 45: 3: -0.139051 15: 0.132927 16: 0.180712 29: -0.174588 43: -0.0153244 44: -0.0283226 45: 0.0822529 78: -0.0179073 80: -0.0168555
Row 46: 4: -0.0554044 5: -0.0302389 20: 0.0856433 46: 0.0367917 48: -0.00755972 50: -0.0138511
Row 47: 0: 0.166667 4: -0.170057 19: -0.122291 20: 0.125681 37: -0.0208333 38: -0.0208333 47: 0.0801801 48: -0.00973948 87: -0.0216808
Row 48: 4: -0.089896 5: 0.0917624 19: 0.0673304 20: -0.0691968 46: -0.00755972 47: -0.00973948 48: 0.0753051 50: -0.0153809 87: -0.00709313
Row 49: 5: -0.0503815 6: -0.039568 21: 0.0899494 49: 0.0362477 51: -0.00989199 53: -0.0125954
Row 50: 4: 0.116928 5: -0.12223 20: -0.105442 21: 0.110744 46: -0.0138511 48: -0.0153809 50: 0.0732972 51: -0.0125094 89: -0.0151765
Row 51: 5: -0.0903197 6: 0.0946095 20: 0.085316 21: -0.0896057 49: -0.00989199 50: -0.0125094 51: 0.0727533 53: -0.0137604 89: -0.00881956
Row 52: 6: -0.0418649 7: -0.0507822 22: 0.0926471 52: 0.0361639 54: -0.0126956 56: -0.0104662
Row 53: 5: 0.105423 6: -0.108234 21: -0.0937366 22: 0.0965478 49: -0.0125954 51: -0.0137604 53: 0.072395 54: -0.0108388 92: -0.0132982
Row 54: 6: -0.10005 7: 0.102791 21: 0.0913963 22: -0.0941373 52: -0.0126956 53: -0.0108388 54: 0.0723112 56: -0.0130022 92: -0.0120103
Row 55: 1: 0.166667 7: -0.179065 8: -0.110935 22: 0.123334 39: -0.0208333 40: -0.0208333 55: 0.0812205 56: -0.00690047 58: -0.023933
Row 56: 6: 0.0938736 7: -0.0868901 8: 0.0624833 22: -0.0694667 52: -0.0104662 54: -0.0130022 55: -0.00690047 56: 0.0757178 58: -0.00872036
Row 57: 8: -0.0183232 9: -0.0641858 23: 0.082509 57: 0.0381052 59: -0.0160465 61: -0.0045808
Row 58: 7: 0.130613 8: -0.0842912 22: -0.189707 23: 0.143384 55: -0.023933 56: -0.00872036 58: 0.0794122 59: -0.0234937 95: -0.0123524
Row 59: 8: -0.0859608 9: 0.134098 22: 0.110024 23: -0.158161 57: -0.0160465 58: -0.0234937 59: 0.0779636 61: -0.0174779 95: -0.00401225
Row 60: 9: -0.0270683 10: -0.0725667 24: 0.099635 60: 0.0374048 62: -0.0181417 64: -0.00676708
Row 61: 8: 0.088235 9: -0.104861 23: -0.080099 24: 0.0967247 57: -0.0045808 59: -0.0174779 61: 0.0746551 62: -0.0154439 97: -0.00873723
Row 62: 9: -0.0994593 10: 0.122551 23: 0.111251 24: -0.134342 60: -0.0181417 61: -0.0154439 62: 0.0739547 64: -0.0124961 97: -0.0123687
Row 63: 10: -0.0299194 11: -0.0736497 25: 0.103569 63: 0.0373361 65: -0.0184124 67: -0.00747984
Row 64: 9: 0.0770526 10: -0.0800624 24: -0.0762241 25: 0.0792338 60: -0.00676708 62: -0.0124961 64: 0.0744595 65: -0.0122889 100: -0.00751951
Row 65: 10: -0.11476 11: 0.119425 24: 0.118141 25: -0.122805 63: -0.0184124 64: -0.0122889 65: 0.0743907 67: -0.0114438 100: -0.0172462
Row 66: 2: 0.166667 11: -0.196446 12: -0.122636 25: 0.152416 41: -0.0208333 42: -0.0208333 66: 0.0838693 67: -0.00982573 69: -0.0282782
Row 67: 10: 0.0756945 11: -0.0621696 12: 0.0556974 25: -0.0692223 63: -0.00747984 65: -0.0114438 66: -0.00982573 67: 0.0795386 69: -0.00409861
Row 68: 12: -0.0486254 13: -0.0607518 26: 0.109377 68: 0.036465 70: -0.015188 72: -0.0121564
Row 69: 11: 0.129507 12: -0.0959869 25: -0.197286 26: 0.163765 66: -0.0282782 67: -0.00409861 69: 0.0835178 70: -0.0210432 103: -0.0198981
Row 70: 12: -0.037978 13: 0.0972344 25: 0.0856684 26: -0.144925 68: -0.015188 69: -0.0210432 70: 0.0777802 72: -0.00912064 103: -0.000373871
Row 71: 13: -0.0542939 14: -0.0356996 27: 0.0899935 71: 0.0364054 73: -0.0089249 75: -0.0135735
Row 72: 12: 0.085108 13: -0.0655443 26: -0.110393 27: 0.0908291 68: -0.0121564 70: -0.00912064 72: 0.0733455 73: -0.0154418 105: -0.00726543
Row 73: 13: -0.112321 14: 0.0913278 26: 0.11846 27: -0.097467 71: -0.0089249 72: -0.0154418 73: 0.073286 75: -0.013907 105: -0.0141733
Row 74: 14: -0.0553925 15: -0.0263034 28: 0.0816959 74: 0.0372162 76: -0.00657586 77: -0.0138481
Row 75: 13: 0.109922 14: -0.0986054 27: -0.122737 28: 0.11142 71: -0.0135735 73: -0.013907 75: 0.0732421 76: -0.0171108 107: -0.0107443
Row 76: 14: -0.103096 15: 0.0934724 27: 0.10437 28: -0.0947465 74: -0.00657586 75: -0.0171108 76: 0.0740529 77: -0.0167923 107: -0.00898164
Row 77: 14: 0.122561 15: -0.120612 28: -0.134212 29: 0.132263 74: -0.0138481 76: -0.0167923 77: 0.075446 78: -0.019705 110: -0.0133608
Row 78: 3: 0.0908408 15: -0.0922857 28: 0.0994762 29: -0.0980313 43: -0.00480287 45: -0.0179073 77: -0.019705 78: 0.0762645 110: -0.00516409
Row 79: 16: -0.105612 17: -0.00383926 29: 0.109451 79: 0.0422986 80: -0.000959815 82: -0.0264029
Row 80: 3: 0.0635826 16: -0.127165 17: 0.0635826 29: 0 44: 0.000959815 45: -0.0168555 79: -0.000959815 80: 0.0865168 82: -0.0149358
Row 81: 17: 0.0030286 18: -0.101541 30: 0.0985126 81: 0.0430318 83: -0.0253853 86: 0.00075715
Row 82: 16: 0.165355 17: -0.0992307 29: -0.167081 34: 0.100956 79: -0.0264029 80: -0.0149358 82: 0.0787237 84: -0.0153673 113: -0.00987183
Row 83: 17: -0.117773 18: 0.175156 30: -0.174589 34: 0.117207 81: -0.0253853 83: 0.0802655 84: -0.018262 86: -0.0184037 115: -0.0110397
Row 84: 17: -0.0764721 29: 0.106213 30: 0.104776 34: -0.134517 82: -0.0153673 83: -0.018262 84: 0.0736588 113: -0.011186 115: -0.00793202
Row 85: 18: -0.095484 19: -0.0030286 30: 0.0985126 85: 0.0415175 86: -0.00075715 88: -0.023871
Row 86: 17: 0.0705861 18: -0.141172 19: 0.0705861 30: 1.38778e-17 81: 0.00075715 83: -0.0184037 85: -0.00075715 86: 0.0845493 88: -0.0168894
Row 87: 4: 0.115096 19: -0.0729092 20: -0.156576 30: 0.114389 47: -0.0216808 48: -0.00709313 87: 0.0754888 88: -0.0174631 90: -0.0111342
Row 88: 18: 0.163041 19: -0.10107 20: 0.103365 30: -0.165337 85: -0.023871 86: -0.0168894 87: -0.0174631 88: 0.0784929 90: -0.00837806
Row 89: 5: 0.0959843 20: -0.0870282 21: -0.101785 32: 0.0928289 50: -0.0151765 51: -0.00881956 89: 0.0726899 91: -0.0102697 94: -0.0129375
Row 90: 19: 0.078049 20: -0.0802209 30: -0.101436 32: 0.103608 87: -0.0111342 88: -0.00837806 90: 0.0732211 91: -0.0142248 114: -0.0116772
Row 91: 20: -0.0932833 21: 0.0929875 30: 0.0982741 32: -0.0979783 89: -0.0102697 90: -0.0142248 91: 0.07243 94: -0.0129771 114: -0.0103437
Row 92: 6: 0.101234 21: -0.111805 22: -0.0955738 31: 0.106145 53: -0.0132982 54: -0.0120103 92: 0.0726747 93: -0.0105953 96: -0.015941
Row 93: 21: -0.0710673 22: 0.0823459 31: -0.0949869 32: 0.0837082 92: -0.0105953 93: 0.073308 94: -0.0131515 96: -0.00999121 116: -0.00777561
Row 94: 20: 0.103658 21: -0.115323 31: 0.11602 32: -0.104356 89: -0.0129375 91: -0.0129771 93: -0.0131515 94: 0.0729649 116: -0.0158535
Row 95: 8: 0.0654588 22: -0.0572194 23: -0.092364 31: 0.0841246 58: -0.0123524 59: -0.00401225 95: 0.0762715 96: -0.0107386 98: -0.0102926
Row 96: 21: 0.103729 22: -0.101493 23: 0.104482 31: -0.106718 92: -0.015941 93: -0.00999121 95: -0.0107386 96: 0.0729405 98: -0.015382
Row 97: 9: 0.0844239 23: -0.0898208 24: -0.10047 35: 0.105867 61: -0.00873723 62: -0.0123687 97: 0.0731731 99: -0.0163802 102: -0.0100864
Row 98: 22: 0.102698 23: -0.112508 31: -0.115899 35: 0.12571 95: -0.0102926 96: -0.015382 98: 0.0740746 99: -0.0186823 117: -0.0127451
Row 99: 23: -0.0655623 24: 0.106147 31: 0.0996652 35: -0.14025 97: -0.0163802 98: -0.0186823 99: 0.0742846 102: -0.0101565 117: -0.00623403
Row 100: 10: 0.0990629 24: -0.106083 25: -0.09037 33: 0.0973905 64: -0.00751951 65: -0.0172462 100: 0.0734869 101: -0.015073 104: -0.00927461
Row 101: 24: -0.113922 25: 0.10863 33: -0.122589 35: 0.127881 100: -0.015073 101: 0.0739912 102: -0.0155742 104: -0.0120846 121: -0.016396
Row 102: 23: 0.0809719 24: -0.0629812 33: 0.0846519 35: -0.102643 97: -0.0100864 99: -0.0101565 101: -0.0155742 102: 0.0741821 121: -0.00558875
Row 103: 12: 0.081088 25: -0.0465446 26: -0.133275 33: 0.0987313 69: -0.0198981 70: -0.000373871 103: 0.0774587 104: -0.0134205 106: -0.0112623
Row 104: 24: 0.0854369 25: -0.0941812 26: 0.0995248 33: -0.0907806 100: -0.00927461 101: -0.0120846 103: -0.0134205 104: 0.0725757 106: -0.0114607
Row 105: 13: 0.0857548 26: -0.115276 27: -0.055669 33: 0.0851898 72: -0.00726543 73: -0.0141733 105: 0.0739831 106: -0.00665181 108: -0.0146456
Row 106: 25: 0.0908919 26: -0.109063 27: 0.0898276 33: -0.0716564 103: -0.0112623 104: -0.0114607 105: -0.00665181 106: 0.073246 108: -0.0158051
Row 107: 14: 0.0789038 27: -0.0812277 28: -0.067072 36: 0.0693959 75: -0.0107443 76: -0.00898164 107: 0.075271 109: -0.00602369 112: -0.0113253
Row 108: 26: 0.121803 27: -0.0961255 33: -0.113808 36: 0.0881303 105: -0.0146456 106: -0.0158051 108: 0.0736796 109: -0.0138063 122: -0.00822629
Row 109: 27: -0.142519 28: 0.108436 33: 0.113403 36: -0.0793199 107: -0.00602369 108: -0.0138063 109: 0.0750113 112: -0.0210853 122: -0.0145445
Row 110: 15: 0.0740994 28: -0.0602816 29: -0.089167 34: 0.0753492 77: -0.0133608 78: -0.00516409 110: 0.0754113 111: -0.00893099 113: -0.0099063
Row 111: 28: -0.10505 29: 0.109101 34: -0.0906772 36: 0.0866269 110: -0.00893099 111: 0.0738563 112: -0.0137383 113: -0.0183442 123: -0.00791843
Row 112: 27: 0.129642 28: -0.144413 34: 0.115025 36: -0.100254 107: -0.0113253 109: -0.0210853 111: -0.0137383 112: 0.075109 123: -0.0150181
Row 113: 17: 0.0842313 28: 0.113002 29: -0.118121 34: -0.0791125 82: -0.00987183 84: -0.011186 110: -0.0099063 111: -0.0183442 113: 0.0736066
Row 114: 20: 0.0880834 30: -0.0957807 32: -0.0777241 34: 0.0854215 90: -0.0116772 91: -0.0103437 114: 0.0729866 115: -0.00775388 119: -0.0136015
Row 115: 17: 0.0758868 30: -0.0932702 32: 0.0925576 34: -0.0751743 83: -0.0110397 84: -0.00793202 114: -0.00775388 115: 0.0739746 119: -0.0153855
Row 116: 21: 0.0945164 31: -0.0929254 32: -0.082857 36: 0.081266 93: -0.00777561 94: -0.0158535 116: 0.0737618 118: -0.0129387 120: -0.00737784
Row 117: 23: 0.0759167 31: -0.0362641 35: -0.107555 36: 0.0679026 98: -0.0127451 99: -0.00623403 117: 0.0778452 118: -0.0141437 124: -0.002832
Row 118: 31: -0.159491 32: 0.118413 35: 0.149407 36: -0.108329 116: -0.0129387 117: -0.0141437 118: 0.077165 120: -0.0166647 124: -0.0232081
Row 119: 30: 0.115948 32: -0.114227 34: -0.106008 36: 0.104287 114: -0.0136015 115: -0.0153855 119: 0.0729428 120: -0.0129004 123: -0.0131713
Row 120: 31: 0.0961702 32: -0.10718 34: 0.0921225 36: -0.081113 116: -0.00737784 118: -0.0166647 119: -0.0129004 120: 0.0731832 123: -0.0101302
Row 121: 24: 0.0879389 33: -0.0604965 35: -0.102917 36: 0.0754744 101: -0.016396 102: -0.00558875 121: 0.0747136 122: -0.00933324 124: -0.00953537
Row 122: 27: 0.0910832 33: -0.131322 35: 0.110477 36: -0.0702381 108: -0.00822629 109: -0.0145445 121: -0.00933324 122: 0.0737317 124: -0.018286
Row 123: 28: 0.0917459 32: 0.0932062 34: -0.100593 36: -0.084359 111: -0.00791843 112: -0.0150181 119: -0.0131713 120: -0.0101302 123: 0.0728767
Row 124: 31: 0.10416 33: 0.111286 35: -0.165977 36: -0.0494695 117: -0.002832 118: -0.0232081 121: -0.00953537 122: -0.018286 124: 0.0773384
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.0923003
0
0
0
0
0
0
0
0
0.0578963
0.0863378
0.0954043
0.0944885
0.0888514
0.0779852
0.0596341
0.0330681
0.0360787
0.037316
0.0325369
0.0418652
0.0445357
0.0474396
0.0755392
0.0894102
0.0924353
0.0900805
0.0610257
0.0590117
0.0639013
0.0736543
0.0809406
0.0652939
0.0809321
0
-0.00575858
0
0
0
-0.0350907
0.00263511
-0.00722838
-0.00529619
0
-0.00991716
-0.00628959
0
-0.0163492
-0.0123264
0
-0.0151224
-0.0244265
-0.0872104
-0.00757595
0
-0.0535158
-0.00798349
0
-0.0380047
-0.0140984
0
-0.0275011
-0.0180975
-0.0333814
-0.0243562
-0.0242978
-0.00435856
-0.0108133
-0.014613
-0.00632623
-0.00765274
-0.00565419
-0.00897655
-0.00777747
-0.00842581
-0.00753412
-0.00774811
0.00263938
-0.00823644
-0.00510169
-0.00532211
0.00395315
-0.00825199
0.00370076
0.00421692
-0.00435681
-0.00368699
-0.0092949
-0.0102197
-0.00997628
-0.0140086
-0.0115245
-0.0213015
-0.00731494
-0.00823074
-0.0213212
-0.00565518
-0.00319583
-0.0139782
-0.0188197
-0.0144226
-0.0217368
-0.0142698
-0.00416346
-0.00959546
-0.0103297
-0.00989169
-0.00229709
-0.0091171
-0.0100331
-0.00695374
-0.00194445
-0.00844938
-0.00903942
-0.0111037
-0.0102999
-0.0090134
-0.0131175
-0.00427374
-0.0185946
-0.00885348
-0.0126483
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
).
[ ]: