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
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.
[3]:
Draw(mesh);
3. Declare a finite element space:¶
[4]:
fes = H1(mesh, order=2, dirichlet="bottom|right")
fes.ndof # number of unknowns in this space
[4]:
129
Python’s help system displays further documentation.
[5]:
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.
| hoprolongation: bool = False
| Create high order prolongation operators,
| only available for H1 and L2 on simplicial meshes
| 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
| print: bool = False
| (historic) print some output into file set by 'SetTestoutFile'
| 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
| hoprolongation: bool = false
| (experimental, only trigs) creates high order prolongation,
| and switches off low-order space
|
| 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__(...) method of builtins.PyCapsule instance
| __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
|
| CreateSmoothingBlocks(...)
| CreateSmoothingBlocks(self: ngsolve.comp.FESpace, **kwargs) -> pyngcore.pyngcore.Table_I
|
|
| Create table of smoothing blocks for block-Jacobi/block-Gauss-Seidel preconditioners.
|
| Every table entry describes the set of dofs belonging a Jacobi/Gauss-Seidel block.
|
| Paramters:
|
| blocktype: string | [ string ] | int
| describes blocktype.
| string form ["vertex", "edge", "face", "facet", "vertexedge", ....]
| or list of strings for combining multiple blocktypes
| int is for backward compatibility with old style blocktypes
|
| condense: bool = False
| exclude dofs eliminated by static condensation
|
| Elements(...)
| Elements(*args, **kwargs)
| Overloaded function.
|
| 1. 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,...
|
|
|
| 2. Elements(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.Region) -> Iterator[ngsolve.comp.FESpaceElement]
|
| 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__(...) method of builtins.PyCapsule instance
| __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) class method of pybind11_builtins.pybind11_object
| 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.
[6]:
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:
[7]:
u, v = fes.TnT()
5. Define and assemble linear and bilinear forms:¶
[8]:
a = BilinearForm(fes)
a += grad(u)*grad(v)*dx
a.Assemble()
f = LinearForm(fes)
f += x*v*dx
f.Assemble();
Alternately, we can do one-liners:
[9]:
a = BilinearForm(grad(u)*grad(v)*dx).Assemble()
f = LinearForm(x*v*dx).Assemble()
You can examine the linear system in more detail:
[10]:
print(f.vec)
0.000333333
0.00633333
0.00633333
0.000529288
0.00423206
0.00879483
0.0115687
0.0176999
0.0273014
0.0166436
0.0151442
0.0182045
0.0210073
0.0107111
0.00781528
0.00306235
0.00067471
0.00103167
0.00166604
0.00181674
0.0137015
0.0205046
0.0304556
0.0291057
0.0260001
0.0287722
0.0186625
0.0219344
0.0100352
0.00291583
0.00518005
0.00777967
0.0205471
0.0158929
0.0207949
0.0129617
0.0158859
0.0179665
-6.66667e-05
-3.33333e-05
-0.000766667
-0.0008
-0.0008
-0.000766667
-6.06392e-05
-1.49949e-05
-8.31522e-05
-0.000289108
-0.000214565
-0.000459295
-0.000398603
-0.000676002
-0.000812886
-0.000522023
-0.000916349
-0.00105268
-0.00172646
-0.00143407
-0.000777223
-0.00178112
-0.00160576
-0.000646778
-0.00139999
-0.00128445
-0.000648431
-0.00122331
-0.00122479
-0.00158722
-0.00145468
-0.000467351
-0.00139751
-0.00106135
-0.000445345
-0.000885342
-0.000870842
-0.000215685
-0.00076362
-0.000534472
-0.000320341
-0.000177293
-2.19785e-05
-7.42829e-05
-9.11566e-05
-2.9417e-05
-0.000119731
-0.000138374
-4.68667e-05
-0.000234449
-0.00018908
-0.000250255
-0.000763332
-0.000417469
-0.00062207
-0.000994616
-0.000858919
-0.000770058
-0.00135916
-0.0010122
-0.00108461
-0.00096346
-0.000872767
-0.0011199
-0.000934733
-0.000855051
-0.00103163
-0.000963124
-0.000829911
-0.00085758
-0.00063321
-0.000865444
-0.000687957
-0.000829446
-0.00023004
-0.000301583
-0.000494185
-0.0001632
-0.000257852
-0.000363996
-0.000519229
-0.000426959
-0.000777698
-0.000773376
-0.000750942
-0.00057903
-0.000701523
-0.000764691
-0.000809458
-0.000691204
-0.000705732
[11]:
print(a.mat)
Row 0: 0: 1 4: -0.5 19: -0.5 38: -0.0833333 39: -0.0833333 48: 0.166667
Row 1: 1: 1 7: -0.5 8: -0.5 40: -0.0833333 41: -0.0833333 56: 0.166667
Row 2: 2: 1 11: -0.5 12: -0.5 42: -0.0833333 43: -0.0833333 67: 0.166667
Row 3: 3: 0.831314 15: -0.17846 16: -0.173448 29: -0.479406 44: -0.040286 45: -0.039615 46: -0.0586512 79: 0.0700292 81: 0.068523
Row 4: 0: -0.5 4: 1.89829 5: -0.540848 19: -0.168769 20: -0.688673 38: 5.83717e-17 39: 0.0833333 47: -0.0562171 48: -0.141895 49: -0.118269 51: 0.146359 89: 0.0866899
Row 5: 4: -0.540848 5: 1.88645 6: -0.354514 20: -0.188208 21: -0.802875 47: -0.0128244 49: 0.102966 50: -0.0499561 51: -0.173998 52: -0.0776293 54: 0.109042 90: 0.1024
Row 6: 5: -0.354514 6: 1.74447 7: -0.307659 21: -0.489892 22: -0.592403 50: -0.0366167 52: 0.0957024 53: -0.0414522 54: -0.116367 55: -0.0963085 57: 0.0927287 93: 0.102314
Row 7: 1: -0.5 6: -0.307659 7: 1.90404 8: -0.20374 22: -0.892641 40: 0 41: 0.0833333 53: -0.051968 55: 0.103244 56: -0.180139 57: -0.0852332 59: 0.130762
Row 8: 1: -0.5 7: -0.20374 8: 1.79311 9: -0.389613 22: -0.247724 23: -0.452028 40: 0.0833333 41: 0 56: -0.111302 57: 0.0619254 58: -0.0211183 59: -0.0881764 60: -0.078254 62: 0.0860538 96: 0.0675384
Row 9: 8: -0.389613 9: 1.76257 10: -0.271361 23: -0.692568 24: -0.409029 58: -0.0647632 60: 0.129699 61: -0.03061 62: -0.102497 63: -0.0958916 65: 0.0758369 98: 0.0882263
Row 10: 9: -0.271361 10: 1.78763 11: -0.270137 24: -0.87395 25: -0.372178 61: -0.0733159 63: 0.118543 64: -0.031622 65: -0.0756344 66: -0.117365 68: 0.0766449 101: 0.10275
Row 11: 2: -0.5 10: -0.270137 11: 1.98752 12: -0.108449 25: -1.10893 42: 0 43: 0.0833333 64: -0.07203 66: 0.117053 67: -0.196125 68: -0.0630977 70: 0.130867
Row 12: 2: -0.5 11: -0.108449 12: 1.81026 13: -0.251506 25: -0.253079 26: -0.697226 42: 0.0833333 43: 0 67: -0.12082 68: 0.0555613 69: -0.0431575 70: -0.0911215 71: -0.046611 73: 0.0850753 104: 0.07774
Row 13: 12: -0.251506 13: 1.80179 14: -0.417484 26: -0.798894 27: -0.333901 69: -0.0603628 71: 0.10228 72: -0.0389519 73: -0.058616 74: -0.142367 76: 0.108533 106: 0.0894845
Row 14: 13: -0.417484 14: 1.76811 15: -0.286119 27: -0.397757 28: -0.666748 72: -0.0390126 74: 0.108593 75: -0.0488037 76: -0.131902 77: -0.0749666 78: 0.0964902 108: 0.0896011
Row 15: 3: -0.17846 14: -0.286119 15: 1.82432 28: -0.407682 29: -0.952064 44: -0.0820544 46: 0.111798 75: -0.0478511 77: 0.0955376 78: -0.124309 79: -0.0498391 112: 0.0967187
Row 16: 3: -0.173448 16: 1.8265 17: -0.252776 29: -0.969793 30: -0.430482 45: -0.0846322 46: 0.11354 80: -0.0478231 81: -0.0528319 82: -0.119129 84: 0.0899525 115: 0.100924
Row 17: 16: -0.252776 17: 1.75225 18: -0.338143 30: -0.71401 31: -0.447317 80: -0.0548032 82: 0.0969326 83: -0.039654 84: -0.0770282 85: -0.120556 88: 0.0960112 116: 0.0990973
Row 18: 17: -0.338143 18: 1.87019 19: -0.523655 20: -0.213918 31: -0.794476 83: -0.0490532 85: 0.10541 86: -0.014548 87: -0.170635 88: -0.0774621 89: 0.101824 91: 0.104464
Row 19: 0: -0.5 4: -0.168769 18: -0.523655 19: 1.89344 20: -0.701018 38: 0.0833333 39: 6.245e-17 48: -0.144439 49: 0.0892334 86: -0.055731 87: 0.143007 89: -0.115404
Row 20: 4: -0.688673 5: -0.188208 18: -0.213918 19: -0.701018 20: 3.49193 21: -0.526678 31: -0.482287 33: -0.691148 47: 0.0690415 48: 0.119667 49: -0.0739297 51: -0.108848 52: 0.0711749 86: 0.0702791 87: -0.105367 88: 0.0707404 89: -0.0731098 90: -0.0775275 91: -0.0773124 92: -0.0658942 95: 0.0941323 118: 0.0869532
Row 21: 5: -0.802875 6: -0.489892 20: -0.526678 21: 3.54814 22: -0.455722 32: -0.781753 33: -0.491222 50: 0.0865728 51: 0.136488 52: -0.0892479 54: -0.0926185 55: 0.0876944 90: -0.135605 92: 0.0868975 93: -0.109228 94: -0.0634125 95: -0.101245 97: 0.0974871 120: 0.0962176
Row 22: 6: -0.592403 7: -0.892641 8: -0.247724 21: -0.455722 22: 3.65035 23: -0.897438 32: -0.56442 53: 0.0934201 54: 0.099944 55: -0.0946304 56: 0.124774 57: -0.0694209 59: -0.188936 60: 0.105449 93: -0.106594 94: 0.0826032 96: -0.0580765 97: -0.0907341 99: 0.102201
Row 23: 8: -0.452028 9: -0.692568 22: -0.897438 23: 3.58584 24: -0.546045 32: -0.391235 36: -0.606525 58: 0.0858815 59: 0.14635 60: -0.156893 62: -0.07826 63: 0.107807 96: -0.0970122 97: 0.100235 98: -0.0961596 99: -0.113036 100: -0.0562792 103: 0.0793605 121: 0.0780061
Row 24: 9: -0.409029 10: -0.87395 23: -0.546045 24: 3.5706 25: -0.544913 34: -0.374111 36: -0.822549 61: 0.103926 62: 0.0947032 63: -0.130458 65: -0.0767869 66: 0.118519 98: -0.105652 100: 0.101957 101: -0.112246 102: -0.113643 103: -0.056314 105: 0.0845456 125: 0.091449
Row 25: 10: -0.372178 11: -1.10893 12: -0.253079 24: -0.544913 25: 3.71425 26: -0.776042 34: -0.659104 64: 0.103652 65: 0.0765845 66: -0.118207 67: 0.150279 68: -0.0691085 70: -0.197711 71: 0.0896124 101: -0.0914761 102: 0.10571 104: -0.0534755 105: -0.0890633 107: 0.0932034
Row 26: 12: -0.697226 13: -0.798894 25: -0.776042 26: 3.68738 27: -0.711151 34: -0.704067 69: 0.10352 70: 0.157966 71: -0.145282 73: -0.10718 74: 0.136809 104: -0.124305 105: 0.0956799 106: -0.138872 107: -0.0989237 109: 0.120588
Row 27: 13: -0.333901 14: -0.397757 26: -0.711151 27: 3.43912 28: -0.688614 34: -0.446754 35: -0.481784 37: -0.379158 72: 0.0779646 73: 0.080721 74: -0.103035 76: -0.0974815 77: 0.0858098 106: -0.0444174 107: 0.0822215 108: -0.0642451 109: -0.0838915 110: -0.0900435 111: -0.0900722 114: 0.0932045 126: 0.0761289 127: 0.0771363
Row 28: 14: -0.666748 15: -0.407682 27: -0.688614 28: 3.51458 29: -0.681592 30: -0.520895 35: -0.549052 75: 0.0966549 76: 0.120851 77: -0.106381 78: -0.104684 79: 0.0759759 108: -0.114524 110: 0.108442 112: -0.053174 113: -0.0970242 114: -0.109977 115: 0.0907968 117: 0.0930433
Row 29: 3: -0.479406 15: -0.952064 16: -0.969793 28: -0.681592 29: 3.71893 30: -0.636072 44: 0.12234 45: 0.124247 46: -0.166687 78: 0.132503 79: -0.0961661 81: -0.0901711 82: 0.127556 112: -0.132079 113: 0.113174 115: -0.134718
Row 30: 16: -0.430482 17: -0.71401 28: -0.520895 29: -0.636072 30: 3.50313 31: -0.557898 35: -0.643775 80: 0.102626 81: 0.0744799 82: -0.105359 84: -0.0952916 85: 0.111667 112: 0.0885341 113: -0.107391 114: 0.105673 115: -0.0570021 116: -0.119559 117: -0.0992521 119: 0.100875
Row 31: 17: -0.447317 18: -0.794476 20: -0.482287 30: -0.557898 31: 3.52545 33: -0.639524 35: -0.603946 83: 0.0887072 84: 0.0823673 85: -0.0965217 87: 0.132995 88: -0.0892895 91: -0.143349 92: 0.0907355 116: -0.0787626 117: 0.0893782 118: -0.0875396 119: -0.0921121 123: 0.103391
Row 32: 21: -0.781753 22: -0.56442 23: -0.391235 32: 3.56434 33: -0.5832 36: -0.866831 37: -0.376902 93: 0.113508 94: -0.100795 95: 0.117579 96: 0.0875503 97: -0.106988 99: -0.117809 100: 0.0954644 120: -0.109837 121: -0.0414894 122: -0.117138 124: 0.0894585 128: 0.0904968
Row 33: 20: -0.691148 21: -0.491222 31: -0.639524 32: -0.5832 33: 3.4985 35: -0.434135 37: -0.659269 90: 0.110733 91: 0.116197 92: -0.111739 94: 0.081604 95: -0.110466 118: -0.0977775 119: 0.0881677 120: -0.0853952 122: 0.100991 123: -0.101202 124: -0.0765031 127: 0.0853902
Row 34: 24: -0.374111 25: -0.659104 26: -0.704067 27: -0.446754 34: 3.54115 36: -0.686395 37: -0.670721 101: 0.100972 102: -0.120069 103: 0.0814484 104: 0.100041 105: -0.0911622 106: 0.093805 107: -0.0765012 109: -0.139262 111: 0.119916 125: -0.0610592 126: -0.102139 128: 0.0940099
Row 35: 27: -0.481784 28: -0.549052 30: -0.643775 31: -0.603946 33: -0.434135 35: 3.50876 37: -0.796067 108: 0.089168 110: -0.123241 111: 0.114371 113: 0.0912412 114: -0.0889005 116: 0.0992241 117: -0.0831694 118: 0.0983639 119: -0.0969304 123: -0.118433 124: 0.0924253 127: -0.0741181
Row 36: 23: -0.606525 24: -0.822549 32: -0.866831 34: -0.686395 36: 3.65909 37: -0.67679 98: 0.113586 99: 0.128644 100: -0.141142 102: 0.128001 103: -0.104495 121: -0.117275 122: 0.133103 125: -0.120117 126: 0.106515 128: -0.12682
Row 37: 27: -0.379158 32: -0.376902 33: -0.659269 34: -0.670721 35: -0.796067 36: -0.67679 37: 3.55891 109: 0.102565 110: 0.104842 111: -0.144214 120: 0.099015 121: 0.080758 122: -0.116956 123: 0.116244 124: -0.105381 125: 0.0897272 126: -0.0805053 127: -0.0884084 128: -0.0576867
Row 38: 0: -0.0833333 4: 5.83717e-17 19: 0.0833333 38: 0.0416667 39: 1.73472e-17 48: -0.0208333
Row 39: 0: -0.0833333 4: 0.0833333 19: 6.245e-17 38: 1.73472e-17 39: 0.0416667 48: -0.0208333
Row 40: 1: -0.0833333 7: 0 8: 0.0833333 40: 0.0416667 41: 0 56: -0.0208333
Row 41: 1: -0.0833333 7: 0.0833333 8: 0 40: 0 41: 0.0416667 56: -0.0208333
Row 42: 2: -0.0833333 11: 0 12: 0.0833333 42: 0.0416667 43: 0 67: -0.0208333
Row 43: 2: -0.0833333 11: 0.0833333 12: 0 42: 0 43: 0.0416667 67: -0.0208333
Row 44: 3: -0.040286 15: -0.0820544 29: 0.12234 44: 0.0380209 46: -0.0205136 79: -0.0100715
Row 45: 3: -0.039615 16: -0.0846322 29: 0.124247 45: 0.0382888 46: -0.0211581 81: -0.00990376
Row 46: 3: -0.0586512 15: 0.111798 16: 0.11354 29: -0.166687 44: -0.0205136 45: -0.0211581 46: 0.0763097 79: -0.00743582 81: -0.007227
Row 47: 4: -0.0562171 5: -0.0128244 20: 0.0690415 47: 0.0397957 49: -0.00320609 51: -0.0140543
Row 48: 0: 0.166667 4: -0.141895 19: -0.144439 20: 0.119667 38: -0.0208333 39: -0.0208333 48: 0.0786155 49: -0.0152763 89: -0.0146404
Row 49: 4: -0.118269 5: 0.102966 19: 0.0892334 20: -0.0739297 47: -0.00320609 48: -0.0152763 49: 0.0767445 51: -0.0225353 89: -0.00703203
Row 50: 5: -0.0499561 6: -0.0366167 21: 0.0865728 50: 0.0364146 52: -0.00915418 54: -0.012489
Row 51: 4: 0.146359 5: -0.173998 20: -0.108848 21: 0.136488 47: -0.0140543 49: -0.0225353 51: 0.0785535 52: -0.0131578 90: -0.0209641
Row 52: 5: -0.0776293 6: 0.0957024 20: 0.0711749 21: -0.0892479 50: -0.00915418 51: -0.0131578 52: 0.0751725 54: -0.0147714 90: -0.00463591
Row 53: 6: -0.0414522 7: -0.051968 22: 0.0934201 53: 0.0361742 55: -0.012992 57: -0.010363
Row 54: 5: 0.109042 6: -0.116367 21: -0.0926185 22: 0.099944 50: -0.012489 52: -0.0147714 54: 0.0726586 55: -0.0106656 93: -0.0143204
Row 55: 6: -0.0963085 7: 0.103244 21: 0.0876944 22: -0.0946304 53: -0.012992 54: -0.0106656 55: 0.0724182 57: -0.0128191 93: -0.011258
Row 56: 1: 0.166667 7: -0.180139 8: -0.111302 22: 0.124774 40: -0.0208333 41: -0.0208333 56: 0.0813494 57: -0.00699219 59: -0.0242014
Row 57: 6: 0.0927287 7: -0.0852332 8: 0.0619254 22: -0.0694209 53: -0.010363 55: -0.0128191 56: -0.00699219 57: 0.0758569 59: -0.00848916
Row 58: 8: -0.0211183 9: -0.0647632 23: 0.0858815 58: 0.0377042 60: -0.0161908 62: -0.00527958
Row 59: 7: 0.130762 8: -0.0881764 22: -0.188936 23: 0.14635 56: -0.0242014 57: -0.00848916 59: 0.0795998 60: -0.0230325 96: -0.0135549
Row 60: 8: -0.078254 9: 0.129699 22: 0.105449 23: -0.156893 58: -0.0161908 59: -0.0230325 60: 0.0776213 62: -0.0162339 96: -0.00332965
Row 61: 9: -0.03061 10: -0.0733159 24: 0.103926 61: 0.0372882 63: -0.018329 65: -0.00765251
Row 62: 8: 0.0860538 9: -0.102497 23: -0.07826 24: 0.0947032 58: -0.00527958 60: -0.0162339 62: 0.0740462 63: -0.0142854 98: -0.00939036
Row 63: 9: -0.0958916 10: 0.118543 23: 0.107807 24: -0.130458 61: -0.018329 62: -0.0142854 63: 0.0736302 65: -0.0113067 98: -0.0126662
Row 64: 10: -0.031622 11: -0.07203 25: 0.103652 64: 0.0371687 66: -0.0180075 68: -0.00790549
Row 65: 9: 0.0758369 10: -0.0756344 24: -0.0767869 25: 0.0765845 61: -0.00765251 63: -0.0113067 65: 0.0745199 66: -0.0115442 101: -0.00760191
Row 66: 10: -0.117365 11: 0.117053 24: 0.118519 25: -0.118207 64: -0.0180075 65: -0.0115442 66: 0.0744005 68: -0.0112557 101: -0.0180856
Row 67: 2: 0.166667 11: -0.196125 12: -0.12082 25: 0.150279 42: -0.0208333 43: -0.0208333 67: 0.083755 68: -0.00937164 70: -0.028198
Row 68: 10: 0.0766449 11: -0.0630977 12: 0.0555613 25: -0.0691085 64: -0.00790549 66: -0.0112557 67: -0.00937164 68: 0.0792571 70: -0.00451869
Row 69: 12: -0.0431575 13: -0.0603628 26: 0.10352 69: 0.0363595 71: -0.0150907 73: -0.0107894
Row 70: 11: 0.130867 12: -0.0911215 25: -0.197711 26: 0.157966 67: -0.028198 68: -0.00451869 70: 0.0827532 71: -0.0212298 104: -0.0182617
Row 71: 12: -0.046611 13: 0.10228 25: 0.0896124 26: -0.145282 69: -0.0150907 70: -0.0212298 71: 0.0770243 73: -0.0104794 104: -0.00117332
Row 72: 13: -0.0389519 14: -0.0390126 27: 0.0779646 72: 0.0368863 74: -0.00975316 76: -0.00973798
Row 73: 12: 0.0850753 13: -0.058616 26: -0.10718 27: 0.080721 69: -0.0107894 71: -0.0104794 73: 0.0747363 74: -0.0160057 106: -0.00417457
Row 74: 13: -0.142367 14: 0.108593 26: 0.136809 27: -0.103035 72: -0.00975316 73: -0.0160057 74: 0.0752631 76: -0.0173952 106: -0.0181966
Row 75: 14: -0.0488037 15: -0.0478511 28: 0.0966549 75: 0.0360853 77: -0.0119628 78: -0.0122009
Row 76: 13: 0.108533 14: -0.131902 27: -0.0974815 28: 0.120851 72: -0.00973798 74: -0.0173952 76: 0.073919 77: -0.0146324 108: -0.0155802
Row 77: 14: -0.0749666 15: 0.0955376 27: 0.0858098 28: -0.106381 75: -0.0119628 76: -0.0146324 77: 0.073118 78: -0.0119216 108: -0.00682004
Row 78: 14: 0.0964902 15: -0.124309 28: -0.104684 29: 0.132503 75: -0.0122009 77: -0.0119216 78: 0.074235 79: -0.01397 112: -0.0191557
Row 79: 3: 0.0700292 15: -0.0498391 28: 0.0759759 29: -0.0961661 44: -0.0100715 46: -0.00743582 78: -0.01397 79: 0.0761706 112: -0.00502395
Row 80: 16: -0.0478231 17: -0.0548032 30: 0.102626 80: 0.0361889 82: -0.0137008 84: -0.0119558
Row 81: 3: 0.068523 16: -0.0528319 29: -0.0901711 30: 0.0744799 45: -0.00990376 46: -0.007227 81: 0.0761588 82: -0.012639 115: -0.00598097
Row 82: 16: -0.119129 17: 0.0969326 29: 0.127556 30: -0.105359 80: -0.0137008 81: -0.012639 82: 0.0740589 84: -0.0105324 115: -0.01925
Row 83: 17: -0.039654 18: -0.0490532 31: 0.0887072 83: 0.0362661 85: -0.0122633 88: -0.0099135
Row 84: 16: 0.0899525 17: -0.0770282 30: -0.0952916 31: 0.0823673 80: -0.0119558 82: -0.0105324 84: 0.0728304 85: -0.0118671 116: -0.0087247
Row 85: 17: -0.120556 18: 0.10541 30: 0.111667 31: -0.0965217 83: -0.0122633 84: -0.0118671 85: 0.0729075 88: -0.0140893 116: -0.0160496
Row 86: 18: -0.014548 19: -0.055731 20: 0.0702791 86: 0.0393887 87: -0.0139328 89: -0.00363701
Row 87: 18: -0.170635 19: 0.143007 20: -0.105367 31: 0.132995 86: -0.0139328 87: 0.0779137 88: -0.0124089 89: -0.021819 91: -0.0208399
Row 88: 17: 0.0960112 18: -0.0774621 20: 0.0707404 31: -0.0892895 83: -0.0099135 85: -0.0140893 87: -0.0124089 88: 0.0747911 91: -0.00527623
Row 89: 4: 0.0866899 18: 0.101824 19: -0.115404 20: -0.0731098 48: -0.0146404 49: -0.00703203 86: -0.00363701 87: -0.021819 89: 0.0763375
Row 90: 5: 0.1024 20: -0.0775275 21: -0.135605 33: 0.110733 51: -0.0209641 52: -0.00463591 90: 0.0752281 92: -0.0129372 95: -0.014746
Row 91: 18: 0.104464 20: -0.0773124 31: -0.143349 33: 0.116197 87: -0.0208399 88: -0.00527623 91: 0.0752607 92: -0.0149975 118: -0.0140519
Row 92: 20: -0.0658942 21: 0.0868975 31: 0.0907355 33: -0.111739 90: -0.0129372 91: -0.0149975 92: 0.0732061 95: -0.00878712 118: -0.00768643
Row 93: 6: 0.102314 21: -0.109228 22: -0.106594 32: 0.113508 54: -0.0143204 55: -0.011258 93: 0.0729438 94: -0.012328 97: -0.016049
Row 94: 21: -0.0634125 22: 0.0826032 32: -0.100795 33: 0.081604 93: -0.012328 94: 0.0736248 95: -0.0128706 97: -0.0083228 120: -0.00753034
Row 95: 20: 0.0941323 21: -0.101245 32: 0.117579 33: -0.110466 90: -0.014746 92: -0.00878712 94: -0.0128706 95: 0.0733954 120: -0.0165241
Row 96: 8: 0.0675384 22: -0.0580765 23: -0.0970122 32: 0.0875503 59: -0.0135549 60: -0.00332965 96: 0.0761654 97: -0.0106981 99: -0.0111895
Row 97: 21: 0.0974871 22: -0.0907341 23: 0.100235 32: -0.106988 93: -0.016049 94: -0.0083228 96: -0.0106981 97: 0.0729481 99: -0.0143607
Row 98: 9: 0.0882263 23: -0.0961596 24: -0.105652 36: 0.113586 62: -0.00939036 63: -0.0126662 98: 0.0732049 100: -0.0170227 103: -0.0113737
Row 99: 22: 0.102201 23: -0.113036 32: -0.117809 36: 0.128644 96: -0.0111895 97: -0.0143607 99: 0.0740126 100: -0.0182627 121: -0.0138982
Row 100: 23: -0.0562792 24: 0.101957 32: 0.0954644 36: -0.141142 98: -0.0170227 99: -0.0182627 100: 0.0746271 103: -0.00846643 121: -0.00560337
Row 101: 10: 0.10275 24: -0.112246 25: -0.0914761 34: 0.100972 65: -0.00760191 66: -0.0180856 101: 0.0736353 102: -0.0152671 105: -0.0099759
Row 102: 24: -0.113643 25: 0.10571 34: -0.120069 36: 0.128001 101: -0.0152671 102: 0.0740158 103: -0.01475 105: -0.0111605 125: -0.0172502
Row 103: 23: 0.0793605 24: -0.056314 34: 0.0814484 36: -0.104495 98: -0.0113737 100: -0.00846643 102: -0.01475 103: 0.0744751 125: -0.00561208
Row 104: 12: 0.07774 25: -0.0534755 26: -0.124305 34: 0.100041 70: -0.0182617 71: -0.00117332 104: 0.0767803 105: -0.0128147 107: -0.0121955
Row 105: 24: 0.0845456 25: -0.0890633 26: 0.0956799 34: -0.0911622 101: -0.0099759 102: -0.0111605 104: -0.0128147 105: 0.072519 107: -0.0111053
Row 106: 13: 0.0894845 26: -0.138872 27: -0.0444174 34: 0.093805 73: -0.00417457 74: -0.0181966 106: 0.0754537 107: -0.00692977 109: -0.0165215
Row 107: 25: 0.0932034 26: -0.0989237 27: 0.0822215 34: -0.0765012 104: -0.0121955 105: -0.0111053 106: -0.00692977 107: 0.0731924 109: -0.0136256
Row 108: 14: 0.0896011 27: -0.0642451 28: -0.114524 35: 0.089168 76: -0.0155802 77: -0.00682004 108: 0.0733846 110: -0.0130507 114: -0.00924125
Row 109: 26: 0.120588 27: -0.0838915 34: -0.139262 37: 0.102565 106: -0.0165215 107: -0.0136256 109: 0.0744031 111: -0.018294 126: -0.00734725
Row 110: 27: -0.0900435 28: 0.108442 35: -0.123241 37: 0.104842 108: -0.0130507 110: 0.0733955 111: -0.0177596 114: -0.0140599 127: -0.008451
Row 111: 27: -0.0900722 34: 0.119916 35: 0.114371 37: -0.144214 109: -0.018294 110: -0.0177596 111: 0.0743699 126: -0.011685 127: -0.0108331
Row 112: 15: 0.0967187 28: -0.053174 29: -0.132079 30: 0.0885341 78: -0.0191557 79: -0.00502395 112: 0.0747129 113: -0.013864 115: -0.00826956
Row 113: 28: -0.0970242 29: 0.113174 30: -0.107391 35: 0.0912412 112: -0.013864 113: 0.0728079 114: -0.0129839 115: -0.0144296 117: -0.00982642
Row 114: 27: 0.0932045 28: -0.109977 30: 0.105673 35: -0.0889005 108: -0.00924125 110: -0.0140599 113: -0.0129839 114: 0.0725966 117: -0.0134344
Row 115: 16: 0.100924 28: 0.0907968 29: -0.134718 30: -0.0570021 81: -0.00598097 82: -0.01925 112: -0.00826956 113: -0.0144296 115: 0.0744331
Row 116: 17: 0.0990973 30: -0.119559 31: -0.0787626 35: 0.0992241 84: -0.0087247 85: -0.0160496 116: 0.0728261 117: -0.0109659 119: -0.0138401
Row 117: 28: 0.0930433 30: -0.0992521 31: 0.0893782 35: -0.0831694 113: -0.00982642 114: -0.0134344 116: -0.0109659 117: 0.0724293 119: -0.0113786
Row 118: 20: 0.0869532 31: -0.0875396 33: -0.0977775 35: 0.0983639 91: -0.0140519 92: -0.00768643 118: 0.0729761 119: -0.0103925 123: -0.0141985
Row 119: 30: 0.100875 31: -0.0921121 33: 0.0881677 35: -0.0969304 116: -0.0138401 117: -0.0113786 118: -0.0103925 119: 0.072425 123: -0.0116494
Row 120: 21: 0.0962176 32: -0.109837 33: -0.0853952 37: 0.099015 94: -0.00753034 95: -0.0165241 120: 0.0731081 122: -0.0138185 124: -0.0109353
Row 121: 23: 0.0780061 32: -0.0414894 36: -0.117275 37: 0.080758 99: -0.0138982 100: -0.00560337 121: 0.075809 122: -0.0154205 128: -0.00476897
Row 122: 32: -0.117138 33: 0.100991 36: 0.133103 37: -0.116956 120: -0.0138185 121: -0.0154205 122: 0.0742278 124: -0.0114293 128: -0.0178552
Row 123: 31: 0.103391 33: -0.101202 35: -0.118433 37: 0.116244 118: -0.0141985 119: -0.0116494 123: 0.0729978 124: -0.0154099 127: -0.0136511
Row 124: 32: 0.0894585 33: -0.0765031 35: 0.0924253 37: -0.105381 120: -0.0109353 122: -0.0114293 123: -0.0154099 124: 0.0729405 127: -0.00769644
Row 125: 24: 0.091449 34: -0.0610592 36: -0.120117 37: 0.0897272 102: -0.0172502 103: -0.00561208 125: 0.0738938 126: -0.0127791 128: -0.00965271
Row 126: 27: 0.0761289 34: -0.102139 36: 0.106515 37: -0.0805053 109: -0.00734725 111: -0.011685 125: -0.0127791 126: 0.0736078 128: -0.0138498
Row 127: 27: 0.0771363 33: 0.0853902 35: -0.0741181 37: -0.0884084 110: -0.008451 111: -0.0108331 123: -0.0136511 124: -0.00769644 127: 0.0738011
Row 128: 32: 0.0904968 34: 0.0940099 36: -0.12682 37: -0.0576867 121: -0.00476897 122: -0.0178552 125: -0.00965271 126: -0.0138498 128: 0.0743263
6. Solve the system:¶
[12]:
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:
[13]:
print(gfu.vec)
0
0
0
0.0923181
0
0
0
0
0
0
0
0
0.0578964
0.0863391
0.0954146
0.0945008
0.0888215
0.0780489
0.0595585
0.0331336
0.0426142
0.0384925
0.0322658
0.0404939
0.0430216
0.0473617
0.076118
0.0895333
0.0932199
0.0919131
0.0857457
0.0712628
0.0581085
0.0661322
0.0723769
0.0842824
0.0630681
0.0775393
0
-0.00562995
0
0
0
-0.035088
0.002793
-0.00729426
-0.00144059
0
-0.00965991
-0.0118565
0
-0.023484
-0.0130763
0
-0.016122
-0.0242459
-0.0872423
-0.00749163
0
-0.053399
-0.0070418
0
-0.0361924
-0.0128913
0
-0.0258153
-0.0181385
-0.033376
-0.0244207
-0.0242838
-0.00446541
-0.0119724
-0.01453
-0.00607083
-0.0110517
-0.00550777
-0.0090215
-0.00690385
-0.00564518
-0.00320331
-0.00759005
0.00119038
0.000777561
-0.0082326
0.000951081
0.00144642
-0.00822658
-0.00086348
0.0018314
0.00459756
-0.00544951
-0.0133503
-0.0087701
-0.010031
-0.013985
-0.0122591
-0.0222892
-0.0067813
-0.00811706
-0.0208728
-0.00480275
-0.00284346
-0.0139733
-0.0176085
-0.015435
-0.0205775
-0.0168732
-0.00387324
-0.00933138
-0.0134445
-0.0124145
-0.00837241
-0.00130944
-0.00480625
-0.00707163
-0.00509836
-0.00790602
-0.00259686
-0.00365637
-0.00637113
-0.0091458
-0.0114363
-0.00754309
-0.00951891
-0.011657
-0.00386109
-0.0148107
-0.00755826
-0.00980606
You can see the zeros coming from the zero boundary conditions.
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
).