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

minor improvements

parent 6b2ff6d4
No related branches found
No related tags found
No related merge requests found
/* -------------------------------------------------------------------
Tutorial 1 : electrostatic field of a microstrip
Features:
- Physical and abstract regions
- Scalar FunctionSpace with Dirichlet constraint
- Galerkin term for stiffness matrix
To compute the solution in a terminal:
getdp microstrip -solve EleSta_v
getdp microstrip -pos Map
......@@ -11,12 +16,9 @@
Run (button at the bottom of the left panel)
The microstrip.pro file describes the finite element (FE) model.
------------------------------------------------------------------- */
Group { // CG utiliser le mot-cle Regions {} ?
Group {
/* One starts by giving explicit meaningful names to
the Physical regions defined in the "microstrip.msh" mesh file.
There are 2 volume regions and 3 surface regions in this model. */
......@@ -65,30 +67,34 @@ Constraint {
}
Group{
/* The domain of definition lists all regions on which the unknown field "v" is defined.
/* The domain of definition of a FunctionSpace lists all regions
on which a field is defined.
Contrary Vol_xxx and Sur_xxx regions, which contain only volume or surface regions, resp.,
Dom_xxx regions may contain both volume and surface regions.
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.
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)
- a set of scalar basis functions ("BF_Node" means nodal basis functions)
- a constraint (here the Dirichlet boundary conditions)
The FE representation of the unknown fields now reads
The FE representation of the unknown field "v" reads
v = Sum_k sn_k vn_k
where the "vn_k" are the nodal values and "sn_k" the basis functions.
Not all nodal values are unknowns of the problem, due to the "Constraint"
which assigns particular values to the nodes of the region "Sur_Dir_Ele". */
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.
*/
{ Name Hgrad_v_Ele; Type Form0;
BasisFunction {
......@@ -111,7 +117,7 @@ Jacobian {
}
Integration {
{ Name GradGrad ;
{ Name Int ;
Case { {Type Gauss ;
Case { { GeoElement Triangle ; NumberOfPoints 4 ; }
{ GeoElement Quadrangle ; NumberOfPoints 4 ; } }
......@@ -135,7 +141,7 @@ Formulation {
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 matrix
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 []
......@@ -152,7 +158,7 @@ Formulation {
to highlight the symmetry.
But the second argument of the density is always reserved for the trial functions
and the the Dof operator should therefore be there systematically.
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.
Hence the final expression of the density below.
*/
......@@ -163,7 +169,7 @@ Formulation {
Equation {
Galerkin { [ epsr[] * Dof{d v} , {d v} ];
In Vol_Dielectric_Ele;
Jacobian Vol; Integration GradGrad; }
Jacobian Vol; Integration Int; }
}
}
}
......@@ -182,6 +188,23 @@ Resolution {
eps0 = 8.854187818e-12; // permittivity of empty space
/* Post-processing is done in two parts.
The first part defines, in terms of the Formulation,
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 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.
Each Operation specifies
- a quantity to be eveluated,
- 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). */
PostProcessing {
{ Name EleSta_v; NameOfFormulation Electrostatics_v;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment