Let us solve the Poisson problem.
$$\text{find: } u \in H_{0,D}^1 \quad \int_\Omega \nabla u \nabla v = \int_\Omega f v, \quad \text{ for all } v \in H_{0,D}^1.$$import netgen.gui
%gui tk
from ngsolve import *
from netgen.geom2d import unit_square
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
mesh.nv, mesh.ne # number of vertices & elements
(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.
fes = H1(mesh, order=2, dirichlet="bottom|right")
fes.ndof # number of unknowns in this space
125
Python's help system displays further documentation.
help(fes)
Help on FESpace in module ngsolve.comp object: class FESpace(pybind11_builtins.pybind11_object_56) | Finite Element Space | | Provides the functionality for finite element calculations. Use | the finite element space generator functions to construct the space, | this can sometimes provide additional functionality (e.g HCurl). | When createing a FESpace with a generator function the set parameters | are passed to the FESpace constructor with and the type parameter is | set. If the space has additional functionality it is added to the space. | | Available generator functions | | H1 | HCurl | HDiv | L2 | FacetFESpace | HDivDiv | | 2 __init__ overloads: | 1) To create a registered FESpace | 2) To create a compound FESpace from multiple created FESpaces | | 1) | | Parameters | | type : string | Type of the finite element space. This parameter is automatically | set if the space is constructed with a generator function. | | mesh : ngsolve.Mesh | Mesh on which the finite element space is defined on. | | flags : dict | Provide additional flags for the finite element space, possible options | are: | dgjumps : bool | Enable DG functionality | print : bool | Write additional debug information to testout file. This | file must be set by ngsolve.SetTestoutFile. | | order : int | Order of the finite element space | | is_complex : bool | Set to true if you want to specify complex (bi-)linearforms on the | FESpace. | | dirichlet : regexpr | Regular expression string defining the dirichlet boundary. | More than one boundary can be combined by the | operator, | Example: dirichlet = "dirichlet1|dirichlet2" | | definedon : list of bits | Define FESpace only on the given domain numbers. Must be list of | 0s and 1s, 0 for not defined, 1 for defined. | | dim : int | Create multi dimensional FESpace (i.e. [H1]^3) | | 2) | | Parameters: | | spaces : list of ngsolve.FESpace | List of the spaces for the compound finite element space | | flags : dict | Additional flags for the compound FESpace | | Method resolution order: | FESpace | pybind11_builtins.pybind11_object_56 | builtins.object | | Methods defined here: | | CouplingType(...) | CouplingType(self: ngsolve.comp.FESpace, dofnr: int) -> ngsolve.comp.COUPLING_TYPE | | get coupling type of a degree of freedom | | Elements(...) | Elements(*args, **kwargs) | Overloaded function. | | 1. Elements(self: ngsolve.comp.FESpace, VOL_or_BND: ngsolve.comp.VorB=VorB.VOL, heapsize: int=10000) -> ngsolve.comp.FESpaceElementRange | | 2. Elements(self: ngsolve.comp.FESpace, VOL_or_BND: ngsolve.comp.VorB=VorB.VOL, heap: ngsolve.ngstd.LocalHeap) -> ngsolve.comp.FESpaceElementRange | | FreeDofs(...) | FreeDofs(self: ngsolve.comp.FESpace, coupling: bool=False) -> ngsolve.ngstd.BitArray | | GetDofNrs(...) | GetDofNrs(*args, **kwargs) | Overloaded function. | | 1. GetDofNrs(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.ElementId) -> tuple | | 2. GetDofNrs(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.NodeId) -> tuple | | GetFE(...) | GetFE(*args, **kwargs) | Overloaded function. | | 1. GetFE(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.ElementId) -> object | | 2. GetFE(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.ElementId, arg1: ngsolve.ngstd.LocalHeap) -> ngsolve.fem.FiniteElement | | Range(...) | Range(self: ngsolve.comp.FESpace, arg0: int) -> slice | | SetCouplingType(...) | SetCouplingType(self: ngsolve.comp.FESpace, dofnr: int, coupling_type: ngsolve.comp.COUPLING_TYPE) -> None | | set coupling type of a degree of freedom | | SetOrder(...) | SetOrder(*args, **kwargs) | Overloaded function. | | 1. SetOrder(self: ngsolve.comp.FESpace, element_type: ngsolve.fem.ET, order: object=<ngsolve.ngstd.DummyArgument>, order_left: object=<ngsolve.ngstd.DummyArgument>, order_right: object=<ngsolve.ngstd.DummyArgument>) -> None | | 2. SetOrder(self: ngsolve.comp.FESpace, nodeid: ngsolve.comp.NodeId, order: int) -> None | | SolveM(...) | SolveM(self: ngsolve.comp.FESpace, rho: ngsolve.fem.CoefficientFunction, vec: ngsolve.la.BaseVector, heapsize: int=1000000) -> None | | TestFunction(...) | TestFunction(self: ngsolve.comp.FESpace) -> object | | Gives a proxy to be used as a testfunction for Symbolic Integrators | | TrialFunction(...) | TrialFunction(self: ngsolve.comp.FESpace) -> object | | Gives a proxy to be used as a trialfunction in Symbolic Integrators | | Update(...) | Update(self: ngsolve.comp.FESpace, heapsize: int=1000000) -> None | | update space after mesh-refinement | | __eq__(...) | __eq__(self: ngsolve.comp.FESpace, arg0: ngsolve.comp.FESpace) -> bool | | __init__(self, /, *args, **kwargs) | Initialize self. See help(type(self)) for accurate signature. | | __new__ = patched_new(class_t, *args, **kwargs) | | __ngsid__(...) | __ngsid__(self: ngsolve.comp.FESpace) -> int | | __reduce__(...) | __reduce__(self: object) -> tuple | | __setstate__(...) | __setstate__(self: object, arg0: tuple) -> None | | __str__(...) | __str__(self: ngsolve.comp.FESpace) -> str | | __timing__(...) | __timing__(self: ngsolve.comp.FESpace) -> object | | ---------------------------------------------------------------------- | Data descriptors defined here: | | __dict__ | | components | list of gridfunctions for compound gridfunction | | globalorder | query global order of space | | mesh | | ndof | number of degrees of freedom | | ndofglobal | global number of dofs on MPI-distributed mesh | | order | proxy to set order for individual nodes | | type | type of finite element space
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.
u = fes.TrialFunction() # symbolic object
v = fes.TestFunction() # symbolic object
gfu = GridFunction(fes) # solution
a = BilinearForm(fes, symmetric=True)
a += SymbolicBFI(grad(u)*grad(v))
a.Assemble()
f = LinearForm(fes)
f += SymbolicLFI(x*v)
f.Assemble()
You can examine the linear system in more detail:
print(f.vec)
0.000333333 0.0103012 0.00633333 0.000333333 0.00387671 0.006844 0.00946082 0.0131276 0.00932238 0.0231198 0.0160023 0.024313 0.0168221 0.0111949 0.00812513 0.00410081 0.00276968 0.000493909 0.00137444 0.00297784 0.0117697 0.0200805 0.0283719 0.0282328 0.0279711 0.0134645 0.00752298 0.0155354 0.0257406 0.0228613 0.0214477 0.00993661 0.0137531 0.0156134 0.0216362 0.0237827 0.0210528 -6.66667e-05 -3.33333e-05 -0.000687449 -0.000587035 -0.00126541 -0.0008 -0.000766667 -6.66667e-05 -3.33333e-05 -0.000228362 -0.000220805 -0.00042906 -0.00033644 -0.000514805 -0.000625037 -0.000462964 -0.000842943 -0.000719098 -0.000956807 -0.00116943 -0.000587035 -0.00113898 -0.000776995 -0.00136243 -0.00146798 -0.00149124 -0.000668081 -0.00130141 -0.00120446 -0.00164251 -0.00152774 -0.00132933 -0.00049277 -0.00133031 -0.000399239 -0.00100556 -0.000902792 -0.000250598 -0.00060796 -0.000769314 -0.000227209 -0.000456957 -2.46955e-05 -0.000350467 -0.0001952 -2.46955e-05 -9.87818e-05 -3.89193e-05 -0.000155955 -0.000192762 -0.000355936 -0.000244359 -0.000594642 -0.000462583 -0.000560307 -0.000951433 -0.000705779 -0.000849162 -0.00071472 -0.00122166 -0.00121405 -0.00104716 -0.0010076 -0.00102318 -0.00114876 -0.000998108 -0.000993348 -0.000471213 -0.000771239 -0.000693185 -0.000388595 -0.000563781 -0.000620059 -0.00061514 -0.000777933 -0.000893067 -0.000986579 -0.000911771 -0.00100144 -0.000578666 -0.000611576 -0.000615209 -0.000710176 -0.000769118 -0.000760015 -0.00096259 -0.000889339 -0.000981825
print(a.mat)
Row 0: 0: 1 Row 1: 1: 0.870166 Row 2: 2: 1 Row 3: 3: 1 Row 4: 0: -0.5 4: 1.89594 Row 5: 4: -0.328601 5: 1.75373 Row 6: 5: -0.252529 6: 1.74397 Row 7: 1: -0.355392 6: -0.246158 7: 1.75039 Row 8: 1: -0.315451 8: 2.11144 Row 9: 8: -0.417618 9: 1.74594 Row 10: 9: -0.35985 10: 1.7563 Row 11: 2: -0.5 10: -0.267891 11: 1.8234 Row 12: 2: -0.5 11: -0.157475 12: 1.9414 Row 13: 12: -0.279511 13: 1.73971 Row 14: 13: -0.347754 14: 1.75594 Row 15: 3: -0.5 14: -0.389763 15: 1.87732 Row 16: 3: -0.5 15: -0.203686 16: 1.80816 Row 17: 16: -0.436369 17: 2.02265 Row 18: 17: -0.424487 18: 1.88532 Row 19: 0: -0.5 4: -0.19922 18: -0.505055 19: 1.8246 Row 20: 4: -0.868115 5: -0.442585 19: -0.376723 20: 3.59693 Row 21: 6: -0.713275 7: -0.495693 21: 3.52299 Row 22: 1: -0.199323 7: -0.653151 8: -1.37837 9: -0.102645 21: -0.73776 22: 3.94083 Row 23: 9: -0.484277 10: -0.45587 23: 3.57987 Row 24: 11: -0.1789 12: -1.00442 13: -0.54338 24: 3.70855 Row 25: 14: -0.437784 15: -0.783873 16: -0.40993 25: 3.54731 Row 26: 16: -0.258171 17: -1.1618 18: -0.149049 25: -0.770254 26: 3.78149 Row 27: 5: -0.730015 6: -0.532007 20: -0.57912 21: -0.447455 27: 3.51184 Row 28: 9: -0.381555 21: -0.393059 22: -0.869584 23: -0.884276 28: 3.62025 Row 29: 10: -0.672687 11: -0.719133 23: -0.86098 24: -0.905134 29: 3.7132 Row 30: 13: -0.569066 14: -0.580637 24: -0.482194 25: -0.65918 30: 3.49528 Row 31: 18: -0.806729 19: -0.243604 20: -0.822289 26: -0.783458 31: 3.67336 Row 32: 21: -0.735752 27: -0.733473 28: -0.670829 32: 3.63948 Row 33: 20: -0.508096 27: -0.489774 31: -0.764287 32: -0.695649 33: 3.53766 Row 34: 25: -0.486293 26: -0.658754 30: -0.483077 31: -0.252992 33: -0.601314 34: 3.41839 Row 35: 23: -0.419347 24: -0.594525 29: -0.555266 30: -0.721121 34: -0.401719 35: 3.54681 Row 36: 23: -0.475121 28: -0.420952 32: -0.80378 33: -0.478544 34: -0.534244 35: -0.854834 36: 3.56748 Row 37: 0: -0.0833333 4: 0 19: 0.0833333 37: 0.0416667 Row 38: 0: -0.0833333 4: 0.0833333 19: -3.00685e-18 37: 1.73472e-18 38: 0.0416667 Row 39: 1: -0.0247066 7: -0.065298 22: 0.0900047 39: 0.0373091 Row 40: 1: -0.00851394 8: -0.10635 22: 0.114864 40: 0.0418598 Row 41: 1: -0.111807 7: 0.12453 8: 0.158925 22: -0.171648 39: -0.0163245 40: -0.0265875 41: 0.0791689 Row 42: 2: -0.0833333 11: 2.65991e-18 12: 0.0833333 42: 0.0416667 Row 43: 2: -0.0833333 11: 0.0833333 12: 0 42: 0 43: 0.0416667 Row 44: 3: -0.0833333 15: 2.65991e-18 16: 0.0833333 44: 0.0416667 Row 45: 3: -0.0833333 15: 0.0833333 16: 0 44: 0 45: 0.0416667 Row 46: 4: -0.0517721 5: -0.0385685 20: 0.0903406 46: 0.0362769 Row 47: 0: 0.166667 4: -0.176247 19: -0.113935 20: 0.123515 37: -0.0208333 38: -0.0208333 47: 0.0808464 Row 48: 4: -0.0879702 5: 0.0933354 19: 0.063805 20: -0.0691702 46: -0.00964212 47: -0.00765043 48: 0.0754566 Row 49: 5: -0.0563741 6: -0.0464315 27: 0.102806 49: 0.0362235 Row 50: 4: 0.106539 5: -0.120062 20: -0.0980085 27: 0.111532 46: -0.012943 48: -0.0136917 50: 0.0729587 Row 51: 5: -0.0772839 6: 0.0885198 20: 0.0814321 27: -0.092668 49: -0.0116079 50: -0.0115591 51: 0.0729053 Row 52: 6: -0.0608303 7: -0.0436771 21: 0.104507 52: 0.0363834 Row 53: 6: -0.0832626 7: 0.0847034 21: -0.088476 27: 0.0870353 52: -0.0109193 53: 0.0726545 Row 54: 5: 0.0984623 6: -0.100137 21: 0.102848 27: -0.101173 49: -0.0140935 51: -0.0105221 53: -0.0111997 54: 0.0724945 Row 55: 6: 0.101857 7: -0.0845868 21: -0.124446 22: 0.107177 52: -0.0152076 53: -0.0102566 55: 0.0729122 Row 56: 1: 0.0839386 7: -0.0981704 21: 0.102555 22: -0.0883228 39: -0.00617666 41: -0.014808 55: -0.015904 56: 0.0738379 Row 57: 8: -0.123378 9: 0.00851394 22: 0.114864 57: 0.0461167 Row 58: 1: 0.0610891 8: -0.122178 9: 0.0610891 22: 3.46945e-18 40: -0.00212848 41: -0.0131438 57: 0.00212848 58: 0.0879765 Row 59: 9: -0.0356793 10: -0.0502285 23: 0.0859078 59: 0.0364707 Row 60: 8: 0.192981 9: -0.107927 22: -0.216622 28: 0.131568 57: -0.0308445 58: -0.0174008 60: 0.0854142 Row 61: 9: -0.0852437 10: 0.110204 23: -0.132822 28: 0.107862 59: -0.0125571 61: 0.0746947 Row 62: 9: -0.0706549 22: 0.118866 23: 0.127627 28: -0.175838 60: -0.0233112 61: -0.0206484 62: 0.0775214 Row 63: 10: -0.0584739 11: -0.0420246 29: 0.100499 63: 0.0362867 Row 64: 9: 0.0956543 10: -0.113616 23: -0.105753 29: 0.123715 59: -0.00891982 61: -0.0149938 64: 0.0738369 Row 65: 10: -0.0703982 11: 0.086673 23: 0.0958239 29: -0.112099 63: -0.0105061 64: -0.0175186 65: 0.0736529 Row 66: 2: 0.166667 11: -0.115053 12: -0.188774 24: 0.137161 42: -0.0208333 43: -0.0208333 66: 0.0825183 Row 67: 11: -0.104077 12: 0.131687 24: -0.198853 29: 0.171242 66: -0.0263603 67: 0.0831865 Row 68: 10: 0.103122 11: -0.0427455 24: 0.0915086 29: -0.151886 63: -0.0146185 65: -0.0111621 67: -0.0233529 68: 0.0786216 Row 69: 12: -0.0619621 13: -0.037384 24: 0.0993461 69: 0.0364828 Row 70: 11: 0.0579654 12: -0.0728309 13: 0.0839692 24: -0.0691036 66: -0.00792991 67: -0.00656144 69: -0.009346 70: 0.0773344 Row 71: 13: -0.0456943 14: -0.0414464 30: 0.0871407 71: 0.0362749 Row 72: 12: 0.108547 13: -0.0957353 24: -0.104283 30: 0.091471 69: -0.0154905 70: -0.0116463 72: 0.0726454 Row 73: 13: -0.111138 14: 0.0994054 24: 0.0955002 30: -0.0837673 71: -0.0103616 72: -0.0105802 73: 0.0724375 Row 74: 14: -0.0328517 15: -0.0491798 25: 0.0820314 74: 0.036748 Row 75: 14: -0.120287 15: 0.11414 25: -0.0986898 30: 0.104836 74: -0.0122949 75: 0.0729852 Row 76: 13: 0.103653 14: -0.0980712 25: 0.0896223 30: -0.0952042 71: -0.0114236 73: -0.0144897 75: -0.0123775 76: 0.0725121 Row 77: 3: 0.166667 15: -0.164799 16: -0.119541 25: 0.117674 44: -0.0208333 45: -0.0208333 77: 0.079572 Row 78: 14: 0.0978122 15: -0.0989082 16: 0.0701555 25: -0.0690595 74: -0.00821292 75: -0.0162401 77: -0.00905195 78: 0.0746533 Row 79: 16: 0.000990177 17: -0.0978067 26: 0.0968166 79: 0.0423862 Row 80: 15: 0.115413 16: -0.0779664 25: -0.154113 26: 0.116666 77: -0.0203664 78: -0.00848692 80: 0.0751003 Row 81: 16: -0.104842 17: 0.170535 25: 0.104761 26: -0.170454 79: -0.0244517 80: -0.0181619 81: 0.0795812 Row 82: 17: -0.0958264 18: -0.000990177 26: 0.0968166 82: 0.0418911 Row 83: 16: 0.071738 17: -0.143476 18: 0.071738 26: -4.98733e-18 79: 0.000247544 81: -0.018182 82: -0.000247544 83: 0.0842773 Row 84: 18: -0.0707206 19: -0.00640091 31: 0.0771215 84: 0.0403243 Row 85: 17: 0.166574 18: -0.134482 26: -0.157758 31: 0.125666 82: -0.0239566 83: -0.017687 85: 0.0792703 Row 86: 18: -0.108027 19: 0.0905768 26: 0.0857827 31: -0.0683322 84: -0.00160023 85: -0.0154828 86: 0.0777036 Row 87: 4: 0.126117 19: -0.067403 20: -0.180941 31: 0.122227 47: -0.0232284 48: -0.00830083 87: 0.0777828 Row 88: 18: 0.154896 19: -0.116361 20: 0.120213 31: -0.158748 84: -0.0176801 86: -0.021044 87: -0.0220068 88: 0.0789275 Row 89: 5: 0.100491 20: -0.0760107 27: -0.118997 33: 0.0945165 50: -0.0163238 51: -0.00879891 89: 0.0728818 Row 90: 19: 0.0663852 20: -0.076053 31: -0.0858102 33: 0.095478 87: -0.00854993 88: -0.00804637 90: 0.0747279 Row 91: 20: -0.0993044 27: 0.103985 31: 0.100631 33: -0.105312 89: -0.0134253 90: -0.0129026 91: 0.0723247 Row 92: 7: 0.082499 21: -0.0738587 22: -0.0952466 28: 0.0866064 55: -0.0108901 56: -0.00973461 92: 0.0730163 Row 93: 6: 0.100285 21: -0.117804 27: -0.0999266 32: 0.117445 53: -0.0105591 54: -0.0145122 93: 0.0730766 Row 94: 21: -0.122214 22: 0.11103 28: -0.105413 32: 0.116597 92: -0.0129215 94: 0.0732841 Row 95: 21: -0.0603663 27: 0.0874671 28: 0.0843162 32: -0.111417 93: -0.0144226 94: -0.0134317 95: 0.0736022 Row 96: 9: 0.0639451 21: 0.0942641 22: -0.0849651 28: -0.0732441 60: -0.00958094 62: -0.00640534 92: -0.00873008 94: -0.0148359 96: 0.075785 Row 97: 9: 0.0703022 23: -0.0823257 28: -0.0696312 36: 0.0816547 61: -0.00631716 62: -0.0112584 97: 0.0748341 Row 98: 10: 0.0793903 23: -0.0549424 29: -0.100427 35: 0.0759796 64: -0.0134101 65: -0.00643743 98: 0.0747166 Row 99: 23: -0.115317 29: 0.120209 35: -0.110223 36: 0.105331 98: -0.0116967 99: 0.0738579 Row 100: 23: -0.105484 28: 0.109148 35: 0.104135 36: -0.107799 97: -0.0110906 99: -0.0158591 100: 0.0731175 Row 101: 11: 0.075928 24: -0.0399176 29: -0.123588 35: 0.087578 67: -0.0194577 68: 0.000475734 101: 0.0785904 Row 102: 13: 0.102329 24: -0.110446 30: -0.0991518 35: 0.107269 72: -0.0122875 73: -0.0132948 102: 0.0724909 Row 103: 24: -0.0954889 29: 0.103202 30: 0.0880464 35: -0.0957591 101: -0.0114394 102: -0.0125004 103: 0.0725839 Row 104: 16: 0.0761325 25: -0.0704368 26: -0.095148 34: 0.0894523 80: -0.0110047 81: -0.00802846 104: 0.0734901 Row 105: 14: 0.0954388 25: -0.0828381 30: -0.0976803 34: 0.0850796 75: -0.0138316 76: -0.0100281 105: 0.0725954 Row 106: 25: -0.116082 26: 0.106858 30: 0.102707 34: -0.0934831 104: -0.0127823 105: -0.0105885 106: 0.0726534 Row 107: 18: 0.0875857 26: -0.0825145 31: -0.0866513 34: 0.0815801 85: -0.0159336 86: -0.00596284 107: 0.0749355 Row 108: 25: 0.0940512 26: -0.124373 31: 0.091562 34: -0.06124 104: -0.00958075 106: -0.0139321 107: -0.00572925 108: 0.0738514 Row 109: 21: 0.089532 27: -0.0577045 32: -0.11535 33: 0.0835224 93: -0.0149388 95: -0.00744421 109: 0.073825 Row 110: 20: 0.0910986 27: -0.114839 32: 0.12015 33: -0.0964098 89: -0.0102038 91: -0.0125709 109: -0.0138987 110: 0.0732194 Row 111: 21: 0.0934596 28: -0.0563856 32: -0.127804 36: 0.0907298 94: -0.0157175 95: -0.00764736 111: 0.0739985 Row 112: 23: 0.102078 28: -0.122864 32: 0.123012 36: -0.102226 97: -0.00932305 100: -0.0161965 111: -0.0162334 112: 0.0738121 Row 113: 23: 0.102615 24: 0.0992647 29: -0.130867 35: -0.0710132 98: -0.00729817 99: -0.0183556 101: -0.0104551 103: -0.014361 113: 0.0736061 Row 114: 25: 0.103079 30: -0.130539 34: -0.0821025 35: 0.109562 105: -0.0106815 106: -0.0150883 114: 0.0732885 Row 115: 24: 0.0953117 30: -0.0762038 34: 0.0775358 35: -0.0966436 102: -0.0143167 103: -0.00951119 114: -0.00984418 115: 0.0732586 Row 116: 20: 0.0928885 31: -0.0682693 33: -0.101603 34: 0.076984 90: -0.0109669 91: -0.0122552 116: 0.0743134 Row 117: 26: 0.127308 31: -0.144416 33: 0.133506 34: -0.116399 107: -0.0146658 108: -0.0171612 116: -0.0144339 117: 0.0757449 Row 118: 27: 0.092483 32: -0.133585 33: -0.0494054 36: 0.0905075 109: -0.00698191 110: -0.0161388 118: 0.074733 Row 119: 28: 0.0838742 32: -0.118425 33: 0.0818246 36: -0.0472739 111: -0.00644903 112: -0.0145195 118: -0.00536944 119: 0.0749155 Row 120: 31: 0.095019 33: -0.13405 34: -0.0635953 36: 0.102626 116: -0.00481209 117: -0.0189427 120: 0.0744661 Row 121: 32: 0.129377 33: -0.10283 34: 0.0868304 36: -0.113377 118: -0.0172574 119: -0.0150867 120: -0.0110867 121: 0.0739911 Row 122: 30: 0.108344 34: -0.0828528 35: -0.149221 36: 0.12373 114: -0.0175463 115: -0.00953977 122: 0.0747568 Row 123: 33: 0.100763 34: -0.0700598 35: 0.106613 36: -0.137316 120: -0.0145699 121: -0.0106209 122: -0.019759 123: 0.074104 Row 124: 23: 0.0825931 34: 0.0722701 35: -0.0682748 36: -0.0865885 99: -0.0104737 100: -0.0101746 122: -0.0111734 123: -0.0068941 124: 0.0743339
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:
print(gfu.vec)
0 0 0 0.0924091 0 0 0 0 0 0 0 0 0.0579249 0.0863395 0.095412 0.0945061 0.0887126 0.0780983 0.0595882 0.0330527 0.0347201 0.0306088 0.0226997 0.0485255 0.0702498 0.0920839 0.0799272 0.0334133 0.0451654 0.0475731 0.089964 0.0608074 0.0521306 0.0604622 0.0816294 0.0755969 0.067892 0 -0.00580903 0 0 0.0247054 0 -0.0349026 0.00339392 -0.00739309 0 -0.0100181 -0.00505004 0 -0.0147034 -0.00938597 0 -0.017922 -0.0121489 -0.00529069 -0.0292704 0 -0.0183813 0 -0.0461399 -0.0215322 -0.05256 0 -0.0378931 -0.0206654 -0.0330051 -0.0606934 -0.0276814 -0.024247 -0.00671107 -0.0145564 -0.0107763 -0.00927653 -0.00544311 -0.00971665 -0.00812781 -0.00479834 -0.00711323 -0.00781592 0.00251074 -0.00483204 -0.00805012 0.00347202 -0.0083864 -0.0049253 0.00455957 0.00401028 -0.00241395 -0.00346299 -0.00961795 -0.00949918 -0.0143154 -0.00652048 -0.0197189 -0.00467488 -0.000876779 -0.00958371 -0.00431468 -0.0129541 -0.0262452 -0.00799141 -0.0244375 -0.00989479 -0.00826822 -0.0118755 -0.00901224 -0.00914765 -0.0047194 -0.0107081 -0.0105955 -0.0151661 -0.00862588 -0.0250001 -0.0135096 -0.00898515 -0.00363723 -0.0139331 -0.00600571 -0.00733149 -0.0111797 -0.0153338 -0.0227022 -0.011174 -0.00783562
You can write an entire python script 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
).