Skip to content
Snippets Groups Projects
Commit c3252338 authored by François Henrotte's avatar François Henrotte
Browse files

reworked Gelerkin {} description and moved tuto1 into a new subdirectory

parent fd69ab25
No related branches found
No related tags found
No related merge requests found
Pipeline #
/* -------------------------------------------------------------------
File "mStrip.geo"
File "microstrip.geo"
This file is the geometrical description used by GMSH to produce
the file "mStrip.msh".
the file "microstrip.msh".
------------------------------------------------------------------- */
/* Definition of some parameters for geometrical dimensions, i.e.
......@@ -43,11 +43,12 @@ Line(11) = {5,7};
Line Loop(12) = {8,-2,-1,-7,9}; Plane Surface(13) = {12};
Line Loop(14) = {10,-11,8,3,4,5}; Plane Surface(15) = {14};
/* Definition of Physical entities (surfaces, lines). The Physical
entities tell GMSH the elements and their associated region numbers
to save in the file 'mStrip.msh'. For example, the Region
111 is made of elements of surface 13, while the Region 121 is
made of elements of lines 9, 10 and 11 */
/* Definition of Physical entities (surfaces, lines).
The definition of Physical entities (Surfaces and Lines)
tells GMSH the elements and associated region numbers
that have to be saved in the mesh file 'microstrip.msh'.
For example, Region 111 is made of the triangle elements of the geometric surface 13,
whereas Region 121 is made of line elements of the geometric lines 9, 10 and 11. */
Physical Surface (101) = {15} ; /* Air */
Physical Surface (111) = {13} ; /* Diel1 */
......
File moved
......@@ -33,7 +33,8 @@ Group {
/* We now define abstract regions to be used below
in the definition of the scalar electric potential formulation:
Vol_Dielectric_Ele : dielectric volume regions where "Div ( epsr[] Grad v)" is solved
Vol_Dielectric_Ele : dielectric volume regions where
"Div ( epsr[] Grad v)" is solved
Sur_Dir_Ele : Dirichlet boundary condition (v imposed)
Sur_Neu_Ele : Neumann bondary condition ( epsr[] n.Grad v = 0 )
......@@ -55,7 +56,8 @@ Function {
}
Constraint {
/* As for material laws, the Dirichlet boundary condition is defined piecewise.
/* As for material laws, the Dirichlet boundary condition
is defined piecewise.
The constraint "Dirichlet_Ele" is invoked in the FunctionSpace below */
{ Name Dirichlet_Ele; Type Assign;
......@@ -72,14 +74,16 @@ Group{
Contrary to Vol_xxx and Sur_xxx regions,
which contain only volume or surface regions, resp.,
domains of definition Dom_xxx regions may contain both volume and surface regions.
domains of definition Dom_xxx regions may contain
both volume and surface regions.
Hence the use of the prefixes Vol_, Sur_ and Dom_ to avoid confusions.*/
Dom_Hgrad_v_Ele = Region[ {Vol_Dielectric_Ele, Sur_Dir_Ele, Sur_Neu_Ele} ];
}
FunctionSpace {
/* The function space in which we shall pick the electric scalar potential "v" solution
/* The function space in which we shall pick
the electric scalar potential "v" solution
is definied by
- a domain of definition ("Dom_Hgrad_v_Ele")
- a type ("Form0" means scalar field)
......@@ -88,12 +92,15 @@ FunctionSpace {
The FE representation of the unknown field "v" reads
v = Sum_k sn_k vn_k
v = Sum_k vn_k sn_k
where the "vn_k" are the nodal values and "sn_k" the basis functions.
Not all nodal values are unknowns of the FE problem, due to the "Constraint"
which assigns particular values to the nodes of the region "Sur_Dir_Ele".
GetDP deals with that automatically on basis of the definition of the FunctionSpace.
where the "vn_k" are the nodal values (connectors)
and "sn_k" the basis functions.
Not all connectors are unknowns of the FE problem,
due to the "Constraint", which assigns particular values
to the nodes of the region "Sur_Dir_Ele".
GetDP deals with that automatically
on basis of the definition of the FunctionSpace.
*/
{ Name Hgrad_v_Ele; Type Form0;
......@@ -130,37 +137,73 @@ Formulation {
/* The syntax of the Formulation{} section is a hard nut to crack.
So let's deal with it carefully.
The weak form of the electrostatic problem is particularly simple in this model,
as it has only one term:
The weak form of the electrostatic problem is particularly simple
in this model, as it has only one term:
(epsr[] Grad v, Grad vp)_Vol_Dielectric_Ele = 0 for all vp in S0(Vol_Dielectric_Ele).
(epsr[] Grad v, Grad vp)_Vol_Dielectric_Ele = 0
for all vp in S0(Vol_Dielectric_Ele).
The corresponding Euler-Lagrange are:
The corresponding Euler-Lagrange equations are:
* Div ( epsr[] Grad v) = 0 on Vol_Dielectric_Ele
* epsr[] n.Grad v = 0 on Sur_Neu_Ele
The Galerkin{} statement is a symbolic representation of this weak formulation term.
The first argument represents the density to be integrated over finite elements,
the result of which being the bilinear form, i.e., the element-based matrix
that will be assembled in the linear system of FE equations.
What is a bit confusing is that the two comma-separated arguments of the bracket []
are not exactly interpreted the same way.
The symbol "d" represents the exterior derivative and is a synonym of "Grad".
The expression "d v" stands then for the gradient of the v field,
The Galerkin{} statement is a symbolic representation
of this weak formulation term.
It has got 4 semicolon separated arguments:
* the density [,] to be integrated,
* the domain of integration,
* the jacobian of the transformation reference element -> real element,
* the integration method to be used.
The symbol "d" represents the exterior derivative,
and it is a synonym of "Grad" when applied to a scalar function.
The expression "d v" stands thus for the gradient of the v field,
i.e., for the electric field up to a minus sign.
The Dof{} operator, finally, returns the array of basis functions of its argument.
As the Galerkin method uses as trial (test) functions the basis functions
of the unknown field, the density of the bilinear form should be
What is a bit confusing is that the two comma-separated terms
of the bracket [,], in the first argument,
are not interpreted the same way. Let us unravel this in detail.
As the Galerkin method uses as trial (test) functions
the basis functions "sn_k" of the unknown field "v",
the density should be something like this
[ epsr[] * {d v} , basis_functions_of {d v} ].
Since the second argument is devoted to the trial functions,
the operator "basis_functions_of" would always be there.
It can therefore be made implicit, and,
according to the GetDP syntax, it is omitted.
So, one writes simply "{d v}".
The first argument, on the other hand, can contain
a much wider variety of expressions than the second one.
In this problem, it contains
epsr[] * {d v} = Sum_k vn_k epsr[]*{d sn_k}
which is the eletric displacement vector, up to a sign.
Here, we have two valid syntaxes, with very different meanings.
[ epsr[] * Dof{d v} , Dof{d v} ]
The first argument can be expressed in terms of
the FE expansion of "v" at the present system solution.
This is indicated by invoking the Dof{} operator:
to highlight the symmetry.
But the second argument of the density is always reserved for the trial functions
and the Dof{} operator should therefore always be there.
It can thus be made implicit and, according to the GetDP syntax, it must be omitted.
[ epsr[] * Dof{d v} , {d v} ].
Another option, which would not work here,
is to evaluate the first argument
with the last available already computed solution.
For this the Dof{} operator is omitted:
[ epsr[] * {d v} , {d v} ].
Both choices are commonly used in different contexts,
and we shall come back on this often in subsequent tutorials.
To express the stiffness matrix of the electrostatic problem at hand,
we have to take the first option.
Hence the final expression of the density below.
*/
{ Name Electrostatics_v; Type FemEquation;
Quantity {
......@@ -193,9 +236,9 @@ eps0 = 8.854187818e-12; // permittivity of empty space
which itself refers to the FunctionSpace,
a number of quantities that can be evaluated at the postprocessing stage.
The three quantities defined here are :
- the scalar vector potential
- the electric fied
- the electric displacement
- the scalar vector potential,
- the electric field,
- the electric displacement.
The second part consists in defining groups of post-processing operations,
which can be invoked separately.
The first group is invoked by default when Gmsh is run interactively.
......@@ -204,7 +247,8 @@ eps0 = 8.854187818e-12; // permittivity of empty space
- the region on which the evaluation is done,
- the name of the output file.
The generated post-processing files are automatically displayed by Gmsh
if the "Merge result automatically" option is enabled (which is the default). */
if the "Merge result automatically" option is enabled
(which is the default). */
PostProcessing {
{ Name EleSta_v; NameOfFormulation Electrostatics_v;
......@@ -237,11 +281,11 @@ h = 5.e-3; // vertical position of the cut
PostOperation {
{ Name Map; NameOfPostProcessing EleSta_v;
Operation {
Print [ v, OnElementsOf Dom_Hgrad_v_Ele, File "mStrip_v.pos" ];
Echo[ StrCat["l=PostProcessing.NbViews-1;",
"View[l].IntervalsType = 1;",
"View[l].NbIso = 40;"],
File "tmp.geo", LastTimeStepOnly] ;
Print [ v, OnElementsOf Dom_Hgrad_v_Ele, File "mStrip_v.pos" ];
Print [ e, OnLine {{e,h,0}{14.e-3,h,0}}{60}, File "Cut_e.pos" ];
}
}
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment