[Top] [Contents] [Index] [ ? ]

GetDP 1.2

Patrick Dular and Christophe Geuzaine

GetDP is a scientific computation software for the numerical solution of integro-differential equations, using finite element and integral type methods. This is the GetDP Reference Manual for GetDP 1.2 (November, 29 2010).

Copying conditions  Terms and conditions of use.
Introduction  What is GetDP?
1. Overview  Quick overview of the general philosophy of GetDP.
2. Expressions  Definition of basic expressions in GetDP.
3. Objects  Definition of the 10 GetDP objects.
4. Types for objects  Definition of all available types for the 10 objects.
5. Short examples  Simple object examples.
6. Complete examples  Some simple complete examples.
7. Running GetDP  How to run GetDP on your operating system.
8. File formats  Input and output file formats.
9. Bugs, versions and credits  Contact information and ChangeLog
A. Tips and tricks  Some tips to make your life easier with GetDP.
B. Frequently asked questions  The GetDP FAQ
C. Gmsh examples  Sample Gmsh input files.
D. License  Complete copy of the license.
Concept index  Index of concepts.
Metasyntactic variable index  Index of metasyntactic variables used in this manual.
Syntax index  Index of reserved keywords in the GetDP language.

 -- The Detailed Node Listing ---

Introduction

Research and collaboration activities  
Education and training activities  
Industrial studies  
How to read this manual  

Overview

1.1 Numerical tools as objects  
1.2 Syntactic rules used in this document  
1.3 Comments  
1.4 Includes  
1.5 Which problems can GetDP actually solve?  

Expressions 

2.1 Definition  
2.2 Constants  
2.3 Operators  
2.4 Functions  
2.5 Current values  
2.6 Arguments  
2.7 Registers  
2.8 Fields  
2.9 Loops and conditionals  

Definition

2.3 Operators  
2.2 Constants  
2.4 Functions  
2.5 Current values  
2.8 Fields  

Operators

2.3.1 Operator types  
2.3.2 Evaluation order  

Objects

3.1 Group: defining topological entities  
3.2 Function: defining global and piecewise expressions  
3.3 Constraint: specifying constraints on function spaces and formulations  
3.4 FunctionSpace: building function spaces  
3.5 Jacobian: defining jacobian methods  
3.6 Integration: defining integration methods  
3.7 Formulation: building equations  
3.8 Resolution: solving systems of equations  
3.9 PostProcessing: exploiting computational results  
3.10 PostOperation: exporting results  

Types for objects

4.1 Types for Group  
4.2 Types for Function  
4.3 Types for Constraint  
4.4 Types for FunctionSpace  
4.5 Types for Jacobian  
4.6 Types for Integration  
4.7 Types for Formulation  
4.8 Types for Resolution  
4.9 Types for PostProcessing  
4.10 Types for PostOperation  

Types for Function

4.2.1 Math functions  
4.2.2 Extended math functions  
4.2.3 Green functions  
4.2.4 Type manipulation functions  
4.2.5 Coordinate functions  
4.2.6 Miscellaneous functions  

Short examples

5.1 Constant expression examples  
5.2 Group examples  
5.3 Function examples  
5.4 Constraint examples  
5.5 FunctionSpace examples  
5.6 Jacobian examples  
5.7 Integration examples  
5.8 Formulation examples  
5.9 Resolution examples  
5.10 PostProcessing examples  
5.11 PostOperation examples  

FunctionSpace examples

5.5.1 Nodal finite element spaces  
5.5.2 High order nodal finite element space  
5.5.3 Nodal finite element space with floating potentials  
5.5.4 Edge finite element space  
5.5.5 Edge finite element space with gauge condition  
5.5.6 Coupled edge and nodal finite element spaces  
5.5.7 Coupled edge and nodal finite element spaces for multiply connected domains  

Formulation examples

5.8.1 Electrostatic scalar potential formulation  
5.8.2 Electrostatic scalar potential formulation with floating potentials and electric charges  
5.8.3 Magnetostatic 3D vector potential formulation  
5.8.4 Magnetodynamic 3D or 2D magnetic field and magnetic scalar potential formulation  
5.8.5 Nonlinearities, Mixed formulations, ...  

Resolution examples

5.9.1 Static resolution (electrostatic problem)  
5.9.2 Frequency domain resolution (magnetodynamic problem)  
5.9.3 Time domain resolution (magnetodynamic problem)  
5.9.4 Nonlinear time domain resolution (magnetodynamic problem)  
5.9.5 Coupled formulations  

Complete examples

6.1 Electrostatic problem  
6.2 Magnetostatic problem  
6.3 Magnetodynamic problem  

File formats

8.1 Input file format  
8.2 Output file format  

Output file format

8.2.1 File `.pre'  
8.2.2 File `.res'  

Bugs, versions and credits

9.1 Bugs  
9.2 Versions  
9.3 Credits  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

Copying conditions

GetDP is "free software"; this means that everyone is free to use it and to redistribute it on a free basis. GetDP is not in the public domain; it is copyrighted and there are restrictions on its distribution, but these restrictions are designed to permit everything that a good cooperating citizen would want to do. What is not allowed is to try to prevent others from further sharing any version of GetDP that they might get from you.

Specifically, we want to make sure that you have the right to give away copies of GetDP, that you receive source code or else can get it if you want it, that you can change GetDP or use pieces of GetDP in new free programs, and that you know you can do these things.

To make sure that everyone has such rights, we have to forbid you to deprive anyone else of these rights. For example, if you distribute copies of GetDP, you must give the recipients all the rights that you have. You must make sure that they, too, receive or can get the source code. And you must tell them their rights.

Also, for our own protection, we must make certain that everyone finds out that there is no warranty for GetDP. If GetDP is modified by someone else and passed on, we want their recipients to know that what they have is not what we distributed, so that any problems introduced by others will not reflect on our reputation.

The precise conditions of the license for GetDP are found in the General Public License that accompanies the source code (see section D. License). Further information about this license is available from the GNU Project webpage http://www.gnu.org/copyleft/gpl-faq.html. Detailed copyright information can be found in 9.3 Credits.

The source code and various pre-compiled versions of GetDP (for Unix, Windows and Mac OS) can be downloaded from the web site http://www.geuz.org/getdp/.

If you use GetDP, we would appreciate that you mention it in your work. References and the latest news about GetDP are always available on http://www.geuz.org/getdp/. Please send all GetDP-related questions to the public GetDP mailing list at getdp@geuz.org.

If you want to integrate GetDP into a closed-source software, or want to sell a modified closed-source version of GetDP, please contact one of the authors. You can purchase a version of GetDP under a different license, with "no strings attached" (for example allowing you to take parts of GetDP and integrate them into your own proprietary code).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

Introduction

GetDP (a "General environment for the treatment of Discrete Problems") is a scientific software environment for the numerical solution of integro-differential equations, open to the coupling of physical problems (electromagnetic, thermal, etc.) as well as of numerical methods (finite element method, integral methods, etc.). It can deal with such problems of various dimensions (1D, 2D or 3D) and time states (static, transient or harmonic).

The main feature of GetDP is the closeness between its internal structure (written in C), the organization of data defining discrete problems (written by the user in ASCII data files) and the symbolic mathematical expressions of these problems. Its aim is to be welcoming and of easy use for both development and application levels: it consists of a working environment in which the definition of any problem makes use of a limited number of objects, which makes the environment structured and concise. It therefore gives researchers advanced developing tools and a large freedom in adding new functionalities.

The modeling tools provided by GetDP can be tackled at various levels of complexity: this opens the software to a wide range of activities, such as research, collaboration, education, training and industrial studies.

Research and collaboration activities  
Education and training activities  
Industrial studies  
How to read this manual  

Research and collaboration activities

The internal structure of the software is very close to the structure used to define discrete problems in the input data files. As a result, a unicity and conciseness of both development and application levels is obtained without any interface between them, which facilitates work of any kind--particularly validations. GetDP permits to rapidly develop tools for the comparison of methods and the exchange of solutions between research teams, which is one of the main aims of collaboration. Moreover, after a short training, GetDP could be used by teams as a basis for their own developments, while allowing a large freedom in the latter.

Education and training activities

The software environment consists of modeling tools applicable to various physical problems. Education and training can therefore be offered in various domains and at different levels. The closeness between the definition of discrete problems and their symbolic mathematical expressions entails that the theory and mathematical bases of numerical methods, which are essential to anyone who wants to tackle the solving of discrete problems, can be directly followed by their practical applications. These applications can be evolutionary, in the sense that the offered tools are of various levels of complexity. Everyone can tackle, step by step, according to his apprenticeship level, the tools adapted to the solving of more and more complex problems, as well as various methods for the solving of the same problem.

Industrial studies

The treatment of industrial problems can also be facilitated because GetDP is adapted to the study of a wide range of physical problems using various numerical methods. In particular, adapted made to measure and ready-to-use software could be rapidly developed for given problems.

How to read this manual

After reading 1. Overview, which depicts the general philosophy of GetDP, you might want to skip 2. Expressions, 3. Objects and 4. Types for objects and directly run the demo files bundled in the distribution on your computer (see section 7. Running GetDP). You should then open these examples with a text editor and compare their structure with the examples given in 5. Short examples and 6. Complete examples. For each new syntax element that you fall onto, you can then go back to 2. Expressions, 3. Objects, and 4. Types for objects, and find in these chapters the detailed description of the syntactic rules as well as all the available options.

Indexes for many concepts (see section Concept index) and for all the syntax elements (see section Syntax index) are available at the end of this manual.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

1. Overview

1.1 Numerical tools as objects  
1.2 Syntactic rules used in this document  
1.3 Comments  
1.4 Includes  
1.5 Which problems can GetDP actually solve?  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

1.1 Numerical tools as objects

An assembly of computational tools (or objects) in GetDP leads to a problem definition structure, which is a transcription of the mathematical expression of the problem, and forms a text data file: the equations describing a phenomenon, written in a mathematical form adapted to a chosen numerical method, directly constitute data for GetDP.

The resolution of a discrete problem with GetDP requires the definition, in a text data file, of the GetDP objects listed (together with their dependencies) in the following figure and table.

objects-wrap

 
Group           ---
Function        Group
Constraint      Group, Function, (Resolution)
FunctionSpace   Group, Constraint, (Formulation), (Resolution)
Jacobian        Group
Integration     ---
Formulation     Group, Function, (Constraint), FunctionSpace,
                Jacobian, Integration
Resolution      Function, Formulation
PostProcessing  Group, Function, Jacobian, Integration, 
                Formulation, Resolution
PostOperation   Group, PostProcessing

The gathering of all these objects constitutes the problem definition structure, which is a copy of the formal mathematical formulation of the problem. Reading the first column of the table from top to bottom pictures the working philosophy and the linking of operations peculiar to GetDP, from group definition to results visualization. The decomposition highlighted in the figure points out the separation between the objects defining the method of resolution, which may be isolated in a "black box" (bottom) and those defining the data peculiar to a given problem (top).

The computational tools which are in the center of a problem definition structure are formulations (Formulation) and function spaces (FunctionSpace). Formulations define systems of equations that have to be built and solved, while function spaces contain all the quantities, i.e., functions, fields of vectors or covectors, known or not, involved in formulations.

Each object of a problem definition structure must be defined before being referred to by others. A linking which always respects this property is the following: it first contains the objects defining particular data of a problem, such as geometry, physical characteristics and boundary conditions (i.e., Group, Function and Constraint) followed by those defining a resolution method, such as unknowns, equations and related objects (i.e., Jacobian, Integration, FunctionSpace, Formulation, Resolution and PostProcessing). The processing cycle ends with the presentation of the results (i.e., lists of numbers in various formats), defined in PostOperation fields. This decomposition points out the possibility of building black boxes, containing objects of the second group, adapted to treatment of general classes of problems that share the same resolution methods.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

1.2 Syntactic rules used in this document

Here are the rules we tried to follow when writing this user's guide. Note that metasyntactic variable definitions stay valid throughout all the manual (and not only in the sections where the definitions appear). See Metasyntactic variable index, for an index of all metasyntactic variables.

  1. Keywords and literal symbols are printed like this.
  2. Metasyntactic variables (i.e., text bits that are not part of the syntax, but stand for other text bits) are printed like this.
  3. A colon (:) after a metasyntactic variable separates the variable from its definition.
  4. Optional rules are enclosed in < > pairs.
  5. Multiple choices are separated by |.
  6. Three dots (...) indicate a possible repetition of the preceding rule.
  7. For conciseness, the notation rule <, rule > ... is replaced by rule <,...>.
  8. The etc symbol replaces nonlisted rules.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

1.3 Comments

Both C and C++ style comments are supported and can be used in the input data files to comment selected text regions:

  1. the text region comprised between /* and */ pairs is ignored;
  2. the rest of a line after a double slash // is ignored.

Comments cannot be used inside double quotes or inside GetDP keywords.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

1.4 Includes

An input data file can be included in another input data file by placing one of the following commands (expression-char represents a file name) on a separate line, outside the GetDP objects. Any text placed after an include command on the same line is ignored.

 
Include expression-char 
#include expression-char 

See 2.2 Constants, for the definition of the character expression expression-char.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

1.5 Which problems can GetDP actually solve?

The preceding explanations may seem very (too) general. Which are the problems that GetDP can actually solve? To answer this question, here is a list of methods that we have considered and coupled until now:

Numerical methods
finite element method
boundary element method (experimental, undocumented)
volume integral methods (experimental, undocumented)
Geometrical models
one-dimensional models (1D)
two-dimensional models (2D), plane and axisymmetric
three-dimensional models (3D)
Time states
static states
sinusoidal and harmonic states
transient states
eigenvalue problems

These methods have been successfully applied to build coupled physical models involving electromagnetic phenomena (magnetostatics, magnetodynamics, electrostatics, electrokinetics, electrodynamics, wave propagation, lumped electric circuits), acoustic phenomena, thermal phenomena and mechanical phenomena (elasticity, rigid body movement).

As can be guessed from the preceding list, GetDP has been initially developed in the field of computational electromagnetics, which fully uses all the offered coupling features. We believe that this does not interfere with the expected generality of the software because a particular modeling forms a problem definition structure which is totally external to the software: GetDP offers computational tools; the user freely applies them to define and solve his problem.

Nevertheless, specific numerical tools will always need to be implemented to solve specific problems in areas other than those mentionned above. If you think the general phisosophy of GetDP is right for you and your problem, but you discover that GetDP lacks the tools necessary to handle it, let us know: we would love to discuss it with you. For example, at the time of this writing, many areas of GetDP would need to be improved to make GetDP as useful for computational mechanics or computational fluid dynamics as it is for computational electromagnetics... So if you have the skills and some free time, feel free to join the project: we gladly accept all code contributions!


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2. Expressions

This chapter and the next two describe in a rather formal way all the commands that can be used in the ASCII text input files. If you are just beginning to use GetDP, or just want to see what GetDP is all about, you should skip this chapter and the next two for now, have a quick look at 7. Running GetDP, and run the demo problems bundled in the distribution on your computer. You should then open the `.pro' files in a text editor and compare their structure with the examples given in 5. Short examples and 6. Complete examples. Once you have a general idea of how the files are organized, you might want to come back here to learn more about the specific syntax of all the objects, and all the available options.

2.1 Definition  
2.2 Constants  
2.3 Operators  
2.4 Functions  
2.5 Current values  
2.6 Arguments  
2.7 Registers  
2.8 Fields  
2.9 Loops and conditionals  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.1 Definition

Expressions are the basic tool of GetDP. They cover a wide range of functional expressions, from constants to formal expressions containing functions (built-in or user-defined, depending on space and time, etc.), arguments, discrete quantities and their associated differential operators, etc. Note that `white space' (spaces, tabs, new line characters) is ignored inside expressions (as well as inside all GetDP objects).

Expressions are denoted by the metasyntactic variable expression (remember the definition of the syntactic rules in 1.2 Syntactic rules used in this document):

 
expression:
  ( expression ) |
  integer |
  real |
  constant-id |
  quantity |
  argument |
  current-value |
  register-value-set |
  register-value-get |
  operator-unary expression |
  expression operator-binary expression |
  expression operator-ternary-left expression operator-ternary-right expression |
  built-in-function-id [ < expression-list > ] < { expression-cst-list } > |
  function-id [ < expression-list > ] |
  < Real | Complex > [ expression ] |
  Dt [ expression ] |
  AtAnteriorTimeStep [ expression, integer ] |
  Order [ quantity ] |
  Trace [ expression, group-id ] |
  expression ##integer

The following sections introduce the quantities that can appear in expressions, i.e., constant terminals (integer, real) and constant expression identifiers (constant-id, expression-cst-list), discretized fields (quantity), arguments (argument), current values (current-value), register values (register-value-set, register-value-get), operators (operator-unary, operator-binary, operator-ternary-left, operator-ternary-right) and built-in or user-defined functions (built-in-function-id, function-id). The last seven cases in this definition permit to cast an expression as real or complex, get the time derivative or evaluate an expression at an anterior time step, retrieve the interpolation order of a discretized quantity, evaluate the trace of an expression, and print the value of an expression for debugging purposes.

List of expressions are defined as:

 
expression-list: 
  expression <,...>

2.3 Operators  
2.2 Constants  
2.4 Functions  
2.5 Current values  
2.8 Fields  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.2 Constants

The three constant types used in GetDP are integer, real and string. These types have the same meaning and syntax as in the C or C++ programming languages. Besides general expressions (expression), purely constant expressions, denoted by the metasyntactic variable expression-cst, are also used:

 
expression-cst:
  ( expression-cst ) |
  integer |
  real |
  constant-id |
  operator-unary expression-cst |
  expression-cst operator-binary expression-cst |
  expression-cst operator-ternary-left expression-cst operator-ternary-right 
      expression-cst |
  math-function-id [ < expression-cst-list > ]

List of constant expressions are defined as:

 
expression-cst-list:
  expression-cst-list-item <,...>

with

 
expression-cst-list-item:
  expression-cst |
  expression-cst : expression-cst |
  expression-cst : expression-cst : expression-cst |
  constant-id {} |
  constant-id { expression-cst-list } |
  List[ constant-id ] |
  ListAlt[ constant-id, constant-id ] |
  LinSpace[ expression-cst, expression-cst, expression-cst ] |
  LogSpace[ expression-cst, expression-cst, expression-cst ]

The second case in this last definition permits to create a list containing the range of numbers comprised between the two expression-cst, with a unit incrementation step. The third case also permits to create a list containing the range of numbers comprised between the two expression-cst, but with a positive or negative incrementation step equal to the third expression-cst. The fourth and fifth cases permit to reference constant identifiers (constant-ids) of lists of constants and constant identifiers of sublists of constants (see below for the definition of constant identifiers) . The sixth case is a synonym for the fourth. The seventh case permits to create alternate lists: the arguments of ListAlt must be constant-ids of lists of constants of the same dimension. The result is an alternate list of these constants: first constant of argument 1, first constant of argument 2, second constant of argument 1, etc. These kinds of lists of constants are for example often used for function parameters (see section 2.4 Functions). The last two cases permit to create linear and logarithmic lists of numbers, respectively.

Contrary to a general expression which is evaluated at runtime (thanks to an internal stack mechanism), an expression-cst is completely evaluated during the syntactic analysis of the problem (when GetDP reads the `.pro' file). The definition of such constants or lists of constants with identifiers can be made outside or inside any GetDP object. The syntax for the definition of constants is:

 
affectation:
  DefineConstant [ constant-id < = expression-cst > 
                   string-id < = "string" > <,...> ]; |
  constant-id = constant-def; |
  string-id = string-def; |
  Printf("string"); |
  Printf("string", expression-cst-list); |
  Read(constant-id); |
  Read(constant-id)[expression-cst];

with

 
constant-id:
  string |
  string ~ { expression-cst }

constant-def:
  expression-cst-list-item |
  { expression-cst-list } |
  ListFromFile [ expression-char ]

string-id:
  string |
  string ~ { expression-cst }

string-def:
  "string" |
  Str[ expression-char ] |
  StrCat[ expression-char, expression-char ]

Notes:

  1. Five constants are predefined in GetDP: Pi (3.1415926535897932), 0D (0), 1D (1), 2D (2) and 3D (3).
  2. When ~{expression-cst} is appended to a string string, the result is a new string formed by the concatenation of string, _ (an underscore) and the value of the expression-cst. This is most useful in loops (see section 2.9 Loops and conditionals), where it permits to define unique strings automatically. For example,
     
    For i In {1:3}
      x~{i} = i;
    EndFor
    
    is the same as
     
    x_1 = 1;
    x_2 = 2;
    x_3 = 3;
    
  3. The assignment in DefineConstant (zero if no expression-cst is given) is performed only if constant-id has not yet been defined. This kind of explicit default definition mechanism is most useful in general problem definition structures making use of a large number of generic constants, functions or groups. When exploiting only a part of a complex problem definition structure, the default definition mechanism allows to define the quantities of interest only, the others being assigned a default value (that will not be used during the processing but that avoids the error messages produced when references to undefined quantities are made).

See 5.1 Constant expression examples, as well as 5.3 Function examples, for some examples.

Character expressions are defined as follows:

 
expression-char:
  "string" |
  string-id |
  StrCat[ expression-char , expression-char ] |
  Sprintf( expression-char ) |
  Sprintf( expression-char, expression-cst-list ) |
  Date

The third case in this definition permits to concatenate two character expressions; the next two cases permit to print the value of variables using standard C formatting; the last case permits to access the current date.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.3 Operators

2.3.1 Operator types  
2.3.2 Evaluation order  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.3.1 Operator types

The operators in GetDP are similar to the corresponding operators in the C or C++ programming languages.

operator-unary:

-
Unary minus.
!
Logical not.

operator-binary:

^
Exponentiation. The evaluation of the both arguments must result in a scalar value.
*
Multiplication or scalar product, depending on the type of the arguments.
/\
Cross product. The evaluation of both arguments must result in vectors.
/
Division.
%
Modulo. The evaluation of the second argument must result in a scalar value.
+
Addition.
-
Subtraction.
==
Equality.
!=
Inequality.
>
Greater. The evaluation of both arguments must result in scalar values.
>=
Greater or equality. The evaluation of both arguments must result in scalar values.
<
Less. The evaluation of both arguments must result in scalar values.
<=
Less or equality. The evaluation of both arguments must result in scalar values.
&&
Logical `and'. The evaluation of both arguments must result in scalar values.
||
Logical `or'. The evaluation of both arguments must result in floating point values. Warning: the logical `or' always (unlike in C or C++) implies the evaluation of both arguments. That is, the second operand of || is evaluated even if the first one is true.

operator-ternary-left:

?
operator-ternary-right:
:
The only ternary operator, formed by operator-ternary-left and operator-ternary-right is defined as in the C or C++ programming languages. The ternary operator first evaluates its first argument (the expression-cst located before the ?), which must result in a scalar value. If it is true (non-zero) the second argument (located between ? and :) is evaluated and returned; otherwise the third argument (located after :) is evaluated and returned.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.3.2 Evaluation order

The evaluation priorities are summarized below (from stronger to weaker, i.e., ^ has the highest evaluation priority). Parentheses () may be used anywhere to change the order of evaluation.

^
- (unary), !
/\
*, /, %
+, -
<, >, <=, >=
!=, ==
&&, ||
?:


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.4 Functions

Two types of functions coexist in GetDP: user-defined functions (function-id, see 3.2 Function: defining global and piecewise expressions) and built-in functions (built-in-function-id, defined in this section).

Both types of functions are always followed by a pair of brackets [] that can possibly contain arguments (see section 2.6 Arguments). This makes it simple to distinguish a function-id or a built-in-function-id from a constant-id. As shown below, built-in functions might also have parameters, given between braces {}, and which are completely evaluated during the analysis of the syntax (since they are of expression-cst-list type):

 
built-in-function-id [ < expression-list > ] < { expression-cst-list } >

with

 
built-in-function-id:
  math-function-id |
  extended-math-function-id |
  green-function-id |
  type-function-id |
  coord-function-id |
  misc-function-id

Notes:

  1. All possible values for built-in-function-id are listed in 4.2 Types for Function.
  2. Classical mathematical functions (see section 4.2.1 Math functions) are the only functions allowed in a constant definition (see the definition of expression-cst in 2.2 Constants).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.5 Current values

Current values are a special kind of arguments (see section 2.6 Arguments) which return the current integer or floating point value of an internal GetDP variable:

$Time
Value of the current time. This value is set to zero for non time dependent analyses.
$DTime
Value of the current time increment used in a time stepping algorithm.
$Theta
Current theta value in a theta time stepping algorithm.
$TimeStep
Number of the current time step in a time stepping algorithm.
$Iteration
Number of the current iteration in a nonlinear loop.
$EigenvalueReal
Real part of the current eigenvalue.
$EigenvalueImag
Imaginary part of the current eigenvalue.
$X, $XS
Value of the current (destination or source) X-coordinate.
$Y, $YS
Value of the current (destination or source) Y-coordinate.
$Z, $ZS
Value of the current (destination or source) Z-coordinate.
$A, $B, $C
Value of the current parametric coordinates used in the parametric OnGrid PostOperation (see section 4.10 Types for PostOperation).

Note:

  1. The current X, Y and Z coordinates refer to the `physical world' coordinates, i.e., coordinates in which the mesh is expressed.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.6 Arguments

Function arguments can be used in expressions and have the following syntax (integer indicates the position of the argument in the expression-list of the function, starting from 1):

 
argument:
  $integer

See 3.2 Function: defining global and piecewise expressions, and 5.3 Function examples, for more details.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.7 Registers

In many situations, identical parts of expressions are used more than once. If this is not a problem with constant expressions (since expression-csts are evaluated only once during the analysis of the problem definition structure, cf. 2.2 Constants), it may introduce some important overhead while evaluating complex expressions (which are evaluated at runtime, thanks to an internal stack mechanism). In order to circumvent this problem, the evaluation result of any part of an expression can be saved in a register: a memory location where this partial result will be accessible without any costly reevaluation of the partial expression.

Registers have the following syntax:
 
register-value-set:
  expression#integer

register-value-get:
  #integer

Thus, to store any part of an expression in the register 5, one should add #5 directly after the expression. To reuse the value stored in this register, one simply uses #5 instead of the expression it should replace.

See 5.3 Function examples, for an example.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.8 Fields

A discretized quantity (defined in a function space, cf. 3.4 FunctionSpace: building function spaces) is represented between braces {}, and can only appear in well-defined expressions in Formulation (see section 3.7 Formulation: building equations) and PostProcessing (see section 3.9 PostProcessing: exploiting computational results) objects:

 
quantity:
  < quantity-dof > { < quantity-operator > quantity-id } |
  { < quantity-operator > quantity-id } [ expression-cst-list ]

with

 
quantity-id:
  string |
  string ~ { expression-cst }

and

quantity-dof:

Dof
Defines a vector of discrete quantities (vector of Degrees of freedom), to be used only in Equation terms of formulations to define (elementary) matrices. Roughly said, the Dof symbol in front of a discrete quantity indicates that this quantity is an unknown quantity, and should therefore not be considered as already computed.

BF
Indicates that only a basis function will be used (only valid with basis functions associated with regions).

quantity-operator:

d
Exterior derivative (d): applied to a p-form, gives a (p+1)-form.

Grad
Gradient: applied to a scalar field, gives a vector.

Curl
Rot
Curl: applied to a vector field, gives a vector.

Div
Divergence (div): applied to a vector field, gives a scalar.

dInv
d^(-1): applied to a p-form, gives a (p-1)-form.

GradInv
Inverse grad: applied to a gradient field, gives a scalar.

CurlInv
RotInv
Inverse curl: applied to a curl field, gives a vector.

DivInv
Inverse div: applied to a divergence field.

Notes:

  1. While the operators Grad, Curl and Div can be applied to 0, 1 and 2-forms respectively, the exterior derivative operator d is usually preferred with such fields.
  2. The second case permits to evaluate a discretized quantity at a certain position X, Y, Z (when expression-cst-list contains three items) or at a specific time, N time steps ago (when expression-cst-list contains a single item).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

2.9 Loops and conditionals

Loops and conditionals are defined as follows, and can be imbricated:

loop:

For ( expression-cst : expression-cst )
Iterates from the value of the first expression-cst to the value of the second expression-cst, with a unit incrementation step. At each iteration, the commands comprised between `For ( expression-cst : expression-cst )' and the matching EndFor are executed.

For ( expression-cst : expression-cst : expression-cst )
Iterates from the value of the first expression-cst to the value of the second expression-cst, with a positive or negative incrementation step equal to the third expression-cst. At each iteration, the commands comprised between `For ( expression-cst : expression-cst : expression-cst )' and the matching EndFor are executed.

For string In { expression-cst : expression-cst }
Iterates from the value of the first expression-cst to the value of the second expression-cst, with a unit incrementation step. At each iteration, the value of the iterate is affected to an expression named string, and the commands comprised between `For string In { expression-cst : expression-cst }' and the matching EndFor are executed.

For string In { expression-cst : expression-cst : expression-cst }
Iterates from the value of the first expression-cst to the value of the second expression-cst, with a positive or negative incrementation step equal to the third expression-cst. At each iteration, the value of the iterate is affected to an expression named string, and the commands comprised between `For string In { expression-cst : expression-cst : expression-cst }' and the matching EndFor are executed.

EndFor
Ends a matching For command.

If ( expression-cst )
The body enclosed between `If ( expression-cst )' and the matching Endif is evaluated if expression-cst is non-zero.

EndIf
Ends a matching If command.

Loops and conditionals can be used in any of the following objects: Group, Function, Constraint (as well as in a contraint-case), FunctionSpace, Formulation (as well as in the quantity and equation defintions), Resolution (as well as resolution-term, system defintion and operations), PostProcessing (in the definition of the PostQuantities) and PostOperation (as well as in the operation list).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3. Objects

This chapter presents the formal definition of the ten GetDP objects mentioned in 1. Overview. To be concise, all the possible parameters for these objects are not given here (cf. the etc syntactic rule defined in 1.2 Syntactic rules used in this document). Please refer to 4. Types for objects, for the list of all available options.

3.1 Group: defining topological entities  
3.2 Function: defining global and piecewise expressions  
3.3 Constraint: specifying constraints on function spaces and formulations  
3.4 FunctionSpace: building function spaces  
3.5 Jacobian: defining jacobian methods  
3.6 Integration: defining integration methods  
3.7 Formulation: building equations  
3.8 Resolution: solving systems of equations  
3.9 PostProcessing: exploiting computational results  
3.10 PostOperation: exporting results  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3.1 Group: defining topological entities

Meshes (grids) constitute the input data of GetDP. All that is needed by GetDP as a mesh is a file containing a list of nodes (with their coordinates) and a list of geometrical elements with, for each one, a number characterizing its geometrical type (i.e., line, triangle, quadrangle, tetrahedron, hexahedron, prism, etc.), a number characterizing the physical region to which it belongs and the list of its nodes. This minimal input set should be easy to extract from most of the classical mesh file formats (see section 8.1 Input file format, for a complete description of the mesh file format read by GetDP).

Groups of geometrical entities of various types can be considered and are used in many objects. There are region groups, of which the entities are regions, and function groups, with nodes, edges, facets, volumes, groups of nodes, edges of tree, facets of tree, ... of regions.

Amongst region groups, elementary and global groups can be distinguished: elementary groups are relative to single regions (e.g., physical regions in which piecewise defined functions or constraints can be defined) while global groups are relative to sets of regions for which given treatments have to be performed (e.g., domain of integration, support of a function space, etc.).

Groups of function type contain lists of entities built on some region groups (e.g., nodes for nodal elements, edges for edge elements, edges of tree for gauge conditions, groups of nodes for floating potentials, elements on one side of a surface for cuts, etc.).

A definition of initially empty groups can be obtained thanks to a DefineGroup command, so that their identifiers exist and can be referred to in other objects, even if these groups are not explicitly defined. This procedure is similar to the DefineConstant procedure introduced for constants in 2.2 Constants.

The syntax for the definition of groups is:

 
Group {
  < DefineGroup [ group-id <{integer}> <,...> ]; > ...
  < group-id = group-def; > ...
  < group-id += group-def; > ...
  < affectation > ...
  < loop > ...
}

with

 
group-id:
  string |
  string ~ { expression-cst }

group-def:
  group-type [ group-list <, group-sub-type group-list > ] |
  group-id <{<integer>}> |
  #group-list

group-type: 
  Region | Global | NodesOf | EdgesOf | etc

group-list:
  All | group-list-item | { group-list-item <,...> }

group-list-item:
  integer | 
  integer : integer | 
  integer : integer : integer |
  group-id <{<integer>}>

group-sub-type: 
  Not | StartingOn | OnOneSideOf | etc

Notes:

  1. integer as a group-list-item is the only interface with the mesh; with each element is associated a region number, being this integer, and a geometrical type (see section 8.1 Input file format). Ranges of integers can be specified in the same way as ranges of constant expressions in an expression-cst-list-item (see section 2.2 Constants). For example, i:j replaces the list of consecutive integers i, i+1, ..., j-1, j.
  2. Array of groups: DefineGroup[group-id{n}] defines the empty groups group-id{i}, i=1, ..., n. Such a definition is optional, i.e., each group-id{i} can be separately defined, in any order.
  3. #group-list is an abbreviation of Region[group-list].

See 4.1 Types for Group, for the complete list of options and 5.2 Group examples, for some examples.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3.2 Function: defining global and piecewise expressions

A user-defined function can be global in space or piecewise defined in region groups. A physical characteristic is an example of a piecewise defined function (e.g., magnetic permeability, electric conductivity, etc.) and can be simply a constant, for linear materials, or a function of one or several arguments for nonlinear materials. Such functions can of course depend on space coordinates or time, which can be needed to express complex constraints.

A definition of initially empty functions can be made thanks to the DefineFunction command so that their identifiers exist and can be referred to (but cannot be used) in other objects. The syntax for the definition of functions is:

 
Function {
  < DefineFunction [ function-id <,...> ]; > ...
  < function-id [ < group-def > ] = expression; > ...
  < affectation > ...
  < loop > ...
}

with

 
function-id:
  string

Note:

  1. The optional group-def in brackets must be of Region type, and indicates on which region the (piecewise) function is defined. Warning: it is incorrect to write f[reg1]=1; g[reg2]=f[]+1; since the domains of definition of f[] and g[] don't match.

See 4.2 Types for Function, for the complete list of built-in functions and 5.3 Function examples, for some examples.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3.3 Constraint: specifying constraints on function spaces and formulations

Constraints can be referred to in FunctionSpace objects to be used for boundary conditions, to impose global quantities or to initialize quantities. These constraints can be expressed with functions or be imposed by the pre-resolution of another discrete problem. Other constraints can also be defined, e.g., constraints of network type for the definition of circuit connections, to be used in Formulation objects.

The syntax for the definition of constraints is:

 
Constraint {
  { Name constraint-id; Type constraint-type;
    Case {
      { Region group-def; < Type constraint-type; >
        < SubRegion group-def; > < TimeFunction expression; > 
        < RegionRef group-def; > < SubRegionRef group-def; > 
        < Coefficient expression; > < Function expression; >
        < Filter expression; > 
        constraint-val; } ...
      < loop > ...
    } 
  | Case constraint-case-id { 
      { Region group-def; < Type constraint-type; >
        constraint-case-val; } ...
      < loop > ...
    } ...
  } ...
  < affectation > ...
  < loop > ...
}

with

 
constraint-id:
constraint-case-id:
  string |
  string ~ { expression-cst }

constraint-type: 
  Assign | Init | Network | Link | etc

constraint-val:
  Value expression | NameOfResolution resolution-id | etc

constraint-case-val:
  Branch { integer, integer } | etc

Notes:

  1. The constraint type constraint-type defined outside the Case fields is applied to all the cases of the constraint, unless other types are explicitly given in these cases. The default type is Assign.
  2. The region type Region group-def will be the main group-list argument of the group-def to be built for the constraints of FunctionSpaces. The optional region type SubRegion group-def will be the argument of the associated group-sub-type.
  3. expression in Value of constraint-val cannot be time dependent ($Time) because it is evaluated only once during the pre-processing (for efficiency reasons). Time dependences must be defined in TimeFunction expression.

See 4.3 Types for Constraint, for the complete list of options and 5.4 Constraint examples, for some examples.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3.4 FunctionSpace: building function spaces

A FunctionSpace is characterized by the type of its interpolated fields, one or several basis functions and optional constraints (in space and time). Subspaces of a function space can be defined (e.g., for the use with hierarchical elements), as well as direct associations of global quantities (e.g., floating potential, electric charge, current, voltage, magnetomotive force, etc.).

A key point is that basis functions are defined by any number of subsets of functions, being added. Each subset is characterized by associated built-in functions for evaluation, a support of definition and a set of associated supporting geometrical entities (e.g., nodes, edges, facets, volumes, groups of nodes, edges incident to a node, etc.). The freedom in defining various kinds of basis functions associated with different geometrical entities to interpolate a field permits to build made-to-measure function spaces adapted to a wide variety of field approximations (see section 5.5 FunctionSpace examples).

The syntax for the definition of function spaces is:

 
FunctionSpace {
  { Name function-space-id;
    Type function-space-type;
    BasisFunction { 
     { Name basis-function-id; NameOfCoef coef-id; 
       Function basis-function-type
         < { Quantity quantity-id;
             Formulation formulation-id {#integer}; 
             Group group-def; Resolution resolution-id {} } >;
       Support group-def; Entity group-def; } ...
    }
  < SubSpace { 
     { Name sub-space-id; 
       NameOfBasisFunction basis-function-list; } ...
    } >
  < GlobalQuantity { 
     { Name global-quantity-id; Type global-quantity-type; 
       NameOfCoef coef-id; } ...
    } >
  < Constraint { 
     { NameOfCoef coef-id;
       EntityType group-type; < EntitySubType group-sub-type; >
       NameOfConstraint constraint-id <{}>; } ...
    } >
  } ...
  < affectation > ...
  < loop > ...
}

with

 
function-space-id: 
formulation-id:
resolution-id:
  string |
  string ~ { expression-cst }

basis-function-id:
coef-id:
sub-space-id:
global-quantity-id: 
  string

function-space-type:   
  Scalar | Vector | Form0 | Form1 | etc 

basis-function-type:
  BF_Node | BF_Edge | etc 

basis-function-list:
  basis-function-id | { basis-function-id <,...> } 

global-quantity-type:
  AliasOf | AssociatedWith

Notes:

  1. When the definition region of a function type group used as an Entity of a BasisFunction is the same as that of the associated Support, it is replaced by All for more efficient treatments during the computation process (this prevents the construction and the analysis of a list of geometrical entities).
  2. Piecewise defined basis functions: the same Name for several BasisFunction fields permits to define piecewise basis functions; separate NameOfCoefs must be defined for those fields.
  3. Constraint: a constraint is associated with geometrical entities defined by an automatically created Group of type group-type, using the Region defined in a Constraint object as its main argument, and the optional SubRegion in the same object as a group-sub-type argument.
  4. Function: a global basis function (BF_Global or BF_dGlobal) needs parameters, i.e., it is given by the quantity (quantity-id) pre-computed from multiresolutions performed on multiformulations.

See 4.4 Types for FunctionSpace, for the complete list of options and 5.5 FunctionSpace examples, for some examples.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3.5 Jacobian: defining jacobian methods

Jacobian methods can be referred to in Formulation and PostProcessing objects to be used in the computation of integral terms and for changes of coordinates. They are based on Group objects and define the geometrical transformations applied to the reference elements (i.e., lines, triangles, quadrangles, tetrahedra, prisms, hexahedra, etc.). Besides the classical lineic, surfacic and volume Jacobians, the Jacobian object allows the construction of various transformation methods (e.g., infinite transformations for unbounded domains) thanks to dedicated jacobian methods.

The syntax for the definition of Jacobian methods is:

 
Jacobian {
  { Name jacobian-id;
    Case { 
      { Region group-def | All; 
        Jacobian jacobian-type < { expression-cst-list } >; } ...
    } 
  } ...
}

with

 
jacobian-id:
  string

jacobian-type:
  Vol | Sur | VolAxi | etc

Note:

  1. The default case of a Jacobian object is defined by Region All and must follow all the other cases.

See 4.5 Types for Jacobian, for the complete list of options and 5.6 Jacobian examples, for some examples.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3.6 Integration: defining integration methods

Various numerical or analytical integration methods can be referred to in Formulation and PostProcessing objects to be used in the computation of integral terms, each with a set of particular options (number of integration points for quadrature methods--which can be linked to an error criterion for adaptative methods, definition of transformations for singular integrations, etc.). Moreover, a choice can be made between several integration methods according to a criterion (e.g., on the proximity between the source and computation points in integral formulations).

The syntax for the definition of integration methods is:

 
Integration {
  { Name integration-id; < Criterion expression; >
    Case { 
    < { Type integration-type; 
        Case { 
          { GeoElement element-type; NumberOfPoints expression-cst } ... 
        } 
      } ... >
    < { Type Analytic; } ... >
    } 
  } ... 
}

with

 
integration-id:
  string

integration-type:
  Gauss | etc

element-type:
  Line | Triangle | Tetrahedron etc 

See 4.6 Types for Integration, for the complete list of options and 5.7 Integration examples, for some examples.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3.7 Formulation: building equations

The Formulation tool permits to deal with volume, surface and line integrals with many kinds of densities to integrate, written in a form that is similar to their symbolic expressions (it uses the same expression syntax as elsewhere in GetDP), which therefore permits to directly take into account various kinds of elementary matrices (e.g., with scalar or cross products, anisotropies, nonlinearities, time derivatives, various test functions, etc.). In case nonlinear physical characteristics are considered, arguments are used for associated functions. In that way, many formulations can be directly written in the data file, as they are written symbolically. Fields involved in each formulation are declared as belonging to beforehand defined function spaces. The uncoupling between formulations and function spaces allows to maintain a generality in both their definitions.

A Formulation is characterized by its type, the involved quantities (of local, global or integral type) and a list of equation terms. Global equations can also be considered, e.g., for the coupling with network relations.

The syntax for the definition of formulations is:

 
Formulation {
  { Name formulation-id; Type formulation-type; 
    Quantity { 
      { Name quantity-id; Type quantity-type; 
        NameOfSpace function-space-id <{}>
                  < [ sub-space-id | global-quantity-id ] >;
        < Symmetry expression-cst; >
        < [ expression ]; In group-def;
          Jacobian jacobian-id; Integration integration-id; >
        < IndexOfSystem integer; >  } ...
    }
    Equation { 
     < local-term-type 
         { < term-op-type > [ expression, expression ]; 
           In group-def; Jacobian jacobian-id;
           Integration integration-id; } > ...
     < GlobalTerm 
         { < term-op-type > [ expression, expression ]; 
           In group-def; } > ...
     < GlobalEquation 
         { Type Network; NameOfConstraint constraint-id;
           { Node expression; Loop expression; Equation expression;
             In group-def; } ...
         } > ...
     < affectation > ...
     < loop > ...
    }
  } ...
  < affectation > ...
  < loop > ...
}

with

 
formulation-id:
  string |
  string ~ { expression-cst }

formulation-type:
  FemEquation | etc

local-term-type:
  Galerkin | deRham

quantity-type:
  Local | Global | Integral

term-op-type:
  Dt | DtDt | JacNL | etc

Note:

  1. IndexOfSystem permits to resolve ambiguous cases when several quantities belong to the same function space, but to different systems of equations. The integer parameter then specifies the index in the list of an OriginSystem command (see section 3.8 Resolution: solving systems of equations).
  2. A GlobalTerm defines a term to be assembled in an equation associated with a global quantity. This equation is a finite element equation if that global quantity is linked with local quantities.
  3. A GlobalEquation defines a global equation to be assembled in the matrix of the system.

See 4.7 Types for Formulation, for the complete list of options and 5.8 Formulation examples, for some examples.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3.8 Resolution: solving systems of equations

The operations available in a Resolution include: the generation of a linear system, its solving with various kinds of linear solvers, the saving of the solution or its transfer to another system, the definition of various time stepping methods, the construction of iterative loops for nonlinear problems (Newton-Raphson and fixed point methods), etc. Multi-harmonic resolutions, coupled problems (e.g., magneto-thermal) or linked problems (e.g., pre-computations of source fields) are thus easily defined in GetDP.

The Resolution object is characterized by a list of systems to build and their associated formulations, using time or frequency domain, and a list of elementary operations:

 
Resolution {
  { Name resolution-id; 
    System { 
      { Name system-id; NameOfFormulation formulation-list; 
        < Type system-type; >
        < Frequency expression-cst-list-item | 
          Frequency { expression-cst-list }; >
        < DestinationSystem system-id; > 
        < OriginSystem system-id; | OriginSystem { system-id <,...> }; >
        < NameOfMesh expression-char > < Solver expression-char > 
        < loop > } ...
      < loop > ...
    }
    Operation { 
      < resolution-op; > ... 
      < loop > ...
    }
  } ...
  < affectation > ...
  < loop > ...
}

with

 
resolution-id:
system-id:
  string |
  string ~ { expression-cst }

formulation-list:
  formulation-id <{}> | { formulation-id <{}> <,...> }

system-type:
  Real | Complex

resolution-op:
  Generate[system-id] | Solve[system-id] | etc

Notes:

  1. The default type for a system of equations is Real. A frequency domain analysis is defined through the definition of one or several frequencies (Frequency expression-cst-list-item | Frequency { expression-cst-list }). Complex systems of equations with no predefined list of frequencies (e.g., in modal analyses) can be explicitely defined with Type Complex.
  2. NameOfMesh permits to explicitely specify the mesh to be used for the construction of the system of equations.
  3. Solver permits to explicitely specify the name of the solver parameter file to use for the solving of the system of equations. This is ony valid if GetDP was compiled against the default solver library (it is the case if you downloaded a pre-compiled copy of GetDP from the internet).
  4. DestinationSystem permits to specify the destination system of a TransferSolution operation (see section 4.8 Types for Resolution).
  5. OriginSystem permits to specify the systems from which ambiguous quantity definitions can be solved (see section 3.7 Formulation: building equations).

See 4.8 Types for Resolution, for the complete list of options and 5.9 Resolution examples, for some examples.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3.9 PostProcessing: exploiting computational results

The PostProcessing object is based on the quantities defined in a Formulation and permits the construction (thanks to the expression syntax) of any useful piecewise defined quantity of interest:

 
PostProcessing {
  { Name post-processing-id; 
    NameOfFormulation formulation-id <{}>; < NameOfSystem system-id; >
    Quantity {
      { Name post-quantity-id; Value { post-value ... } } ...
      < loop > ...
    }
  } ...
  < affectation > ...
  < loop > ...
}

with

 
post-processing-id:
post-quantity-id:
  string |
  string ~ { expression-cst }

post-value:
  Local { local-value } | Integral { integral-value }

local-value:
  [ expression ]; In group-def; Jacobian jacobian-id;

integral-value:
  [ expression ]; In group-def; 
  Integration integration-id; Jacobian jacobian-id;

Notes:

  1. The quantity defined with integral-value is piecewise defined over the elements of the mesh of group-def, and takes, in each element, the value of the integration of expression over this element. The global integral of expression over a whole region (being either group-def or a subset of group-def) has to be defined in the PostOperation with the post-quantity-id[group-def] command (see section 3.10 PostOperation: exporting results).
  2. If NameOfSystem system-id is not given, the system is automatically selected as the one to which the first quantity listed in the Quantity field of formulation-id is associated.

See 4.9 Types for PostProcessing, for the complete list of options and 5.10 PostProcessing examples, for some examples.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

3.10 PostOperation: exporting results

The PostOperation is the bridge between results obtained with GetDP and the external world. It defines several elementary operations on PostProcessing quantities (e.g., plot on a region, section on a user-defined plane, etc.), and outputs the results in several file formats.

 
PostOperation {
  { Name post-operation-id; NameOfPostProcessing post-processing-id; 
    < Format post-operation-fmt; > < Append expression-char; >
    Operation { 
      < post-operation-op; > ... 
    }
  } ...
  < affectation > ...
  < loop > ...
} |
PostOperation post-operation-id UsingPost post-processing-id {
  < post-operation-op; > ...
} ...

with

 
post-operation-id:
  string |
  string ~ { expression-cst }

post-operation-op:
  Print[ post-quantity-term, print-support <,print-option> ... ] |
  Print[ "string", expression <,print-option> ... ] |
  Print[ "string", Str[ expression-char ] <,print-option> ... ] |
  Echo[ "string" <,print-option> ... ] |
  PrintGroup[ group-id, print-support <,print-option> ... ] |
  < loop > ...
  etc

post-quantity-term:
  post-quantity-id <[group-def]> |
  post-quantity-id post-quantity-op post-quantity-id[group-def] |
  post-quantity-id[group-def] post-quantity-op post-quantity-id

post-quantity-op:
  + | - | * | /

print-support:
  OnElementsOf group-def | OnRegion group-def | OnGlobal | etc

print-option:
  File expression-char | Format post-operation-fmt | etc

post-operation-fmt:
  Table | TimeTable | etc

Notes:

  1. Both PostOperation syntaxes are equivalent. The first one conforms to the overall interface, but the second one is more concise.
  2. The format post-operation-fmt defined outside the Operation field is applied to all the post-processing operations, unless other formats are explicitly given in these operations with the Format option (see section 4.10 Types for PostOperation). The default format is Gmsh.
  3. The optional argument [group-def] of the post-quantity-id can only be used when this quantity has been defined as an integral-value (see section 3.9 PostProcessing: exploiting computational results). In this case, the sum of all elementary integrals is performed over the region group-def.
  4. The post-quantity-op allows the simple combination of space-dependent quantities (post-quantity-id) with global integral quantities (post-quantity-id[group-def]).

See 4.10 Types for PostOperation, for the complete list of options and 5.11 PostOperation examples, for some examples.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4. Types for objects

This chapter presents the complete list of choices associated with metasyntactic variables introduced for the ten GetDP objects.

4.1 Types for Group  
4.2 Types for Function  
4.3 Types for Constraint  
4.4 Types for FunctionSpace  
4.5 Types for Jacobian  
4.6 Types for Integration  
4.7 Types for Formulation  
4.8 Types for Resolution  
4.9 Types for PostProcessing  
4.10 Types for PostOperation  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.1 Types for Group

Types in

 
group-type [ R1 <, group-sub-type R2 > ]

group-type < group-sub-type >:

Region
Regions in R1.

Global
Regions in R1 (variant of Region used with global BasisFunctions BF_Global and BF_dGlobal).

NodesOf
Nodes of elements of R1

< Not: but not those of R2 >.

EdgesOf
Edges of elements of R1

< Not: but not those of R2 >.

FacetsOf
Facets of elements of R1

< Not: but not those of R2 >.

VolumesOf
Volumes of elements of R1

< Not: but not those of R2 >.

ElementsOf
Elements of regions in R1

< OnOneSideOf: only elements on one side of R2) >.

GroupsOfNodesOf
Groups of nodes of elements of R1 (a group is associated with each region).

GroupsOfEdgesOf
Groups of edges of elements of R1 (a group is associated with each region).

< InSupport: in a support R2 being a group of type ElementOf, i.e., containing elements >.

GroupsOfEdgesOnNodesOf
Groups of edges incident to nodes of elements of R1 (a group is associated with each node).

< Not: but not those of R2) >.

EdgesOfTreeIn
Edges of a tree of edges of R1

< StartingOn: a complete tree is first built on R2 >.

FacetsOfTreeIn
Facets of a tree of facets of R1

< StartingOn: a complete tree is first built on R2 >.

DualNodesOf
Dual nodes of elements of R1.

DualEdgesOf
Dual edges of elements of R1.

DualFacetsOf
Dual facets of elements of R1.

DualVolumesOf
Dual volumes of elements of R1.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.2 Types for Function

4.2.1 Math functions  
4.2.2 Extended math functions  
4.2.3 Green functions  
4.2.4 Type manipulation functions  
4.2.5 Coordinate functions  
4.2.6 Miscellaneous functions  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.2.1 Math functions

The following functions are the equivalent of the functions of the C math library, and always return real-valued expressions. These are the only functions allowed in constant expressions (expression-cst, see 2.2 Constants).

math-function-id:

Exp
[expression]

Exponential function: e^expression.

Log
[expression]

Natural logarithm: ln(expression), expression>0.

Log10
[expression]

Base 10 logarithm: log10(expression), expression>0.

Sqrt
[expression]

Square root, expression>=0.

Sin
[expression]

Sine of expression.

Asin
[expression]

Arc sine (inverse sine) of expression in [-Pi/2,Pi/2], expression in [-1,1].

Cos
[expression]

Cosine of expression.

Acos
[expression]

Arc cosine (inverse cosine) of expression in [0,Pi], expression in [-1,1].

Tan
[expression]

Tangent of expression.

Atan
[expression]

Arc tangent (inverse tangent) of expression in [-Pi/2,Pi/2].

Atan2
[expression,expression]

Arc tangent (inverse tangent) of the first expression divided by the second, in [-Pi,Pi].

Sinh
[expression]

Hyperbolic sine of expression.

Cosh
[expression]

Hyperbolic cosine of expression.

Tanh
[expression]

Hyperbolic tangent of expression.

Fabs
[expression]

Absolute value of expression.

Fmod
[expression,expression]

Remainder of the division of the first expression by the second, with the sign of the first.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.2.2 Extended math functions

extended-math-function-id:

Cross
[expression,expression]

Cross product of the two arguments; expression must be a vector.

Hypot
[expression,expression]

Square root of the sum of the squares of its arguments.

Norm
[expression]

Absolute value if expression is a scalar; euclidian norm if expression is a vector.

SquNorm
[expression]

Square norm: Norm[expression]^2.

Unit
[expression]

Normalization: expression/Norm[expression]. Returns 0 if the norm is smaller than 1.e-30.

Transpose
[expression]

Transposition; expression must be a tensor.

TTrace
[expression]

Trace; expression must be a tensor.

F_Cos_wt_p
[]{expression-cst,expression-cst}

The first parameter represents the angular frequency and the second represents the phase. If the type of the current system is Real, F_Cos_wt_p[]{w,p} is identical to Cos[w*$Time+p]. If the type of the current system is Complex, it is identical to Complex[Cos[w],Sin[w]].

F_Sin_wt_p
[]{expression-cst,expression-cst}

The first parameter represents the angular frequency and the second represents the phase. If the type of the current system is Real, F_Sin_wt_p[]{w,p} is identical to Sin[w*$Time+p]. If the type of the current system is Complex, it is identical to Complex[Sin[w],-Cos[w]].

F_Period
[expression]{expression-cst}

Fmod[expression,expression-cst] + (expression<0 ? expression-cst : 0); the result is always in [0,expression-cst[.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.2.3 Green functions

The Green functions are only used in integral quantities (see section 3.7 Formulation: building equations). The first parameter represents the dimension of the problem:

The triplets of values given in the definitions below correspond to the 1D, 2D and 3D cases.

green-function-id:

Laplace
[]{expression-cst}

r/2, 1/(2*Pi)*ln(1/r), 1/(4*Pi*r).

GradLaplace
[]{expression-cst}

Gradient of Laplace relative to the destination point ($X, $Y, $Z).

Helmholtz
[]{expression-cst, expression-cst}

exp(j*k0*r)/(4*Pi*r), where k0 is given by the second parameter.

GradHelmholtz
[]{expression-cst, expression-cst}

Gradient of Helmholtz relative to the destination point ($X, $Y, $Z).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.2.4 Type manipulation functions

type-function-id:

Complex
[expression-list]

Creates a (multi-harmonic) complex expression from an number of real-valued expressions. The number of expressions in expression-list must be even.

Re
[expression]

Takes the real part of a complex-valued expression.

Im
[expression]

Takes the imaginary part of a complex-valued expression.

Vector
[expression,expression,expression]

Creates a vector from 3 scalars.

Tensor
[expression,expression,expression,expression,expression,expression,

expression,expression,expression]

Creates a second-rank tensor of order 3 from 9 scalars.

TensorV
[expression,expression,expression]

Creates a second-rank tensor of order 3 from 3 vectors.

TensorSym
[expression,expression,expression,expression,expression,expression]

Creates a symmetrical second-rank tensor of order 3 from 6 scalars.

TensorDiag
[expression,expression,expression]

Creates a diagonal second-rank tensor of order 3 from 3 scalars.

CompX
[expression]

Gets the X component of a vector.

CompY
[expression]

Gets the Y component of a vector.

CompZ
[expression]

Gets the Z component of a vector.

CompXX
[expression]

Gets the XX component of a tensor.

CompXY
[expression]

Gets the XY component of a tensor.

CompXZ
[expression]

Gets the XZ component of a tensor.

CompYX
[expression]

Gets the YX component of a tensor.

CompYY
[expression]

Gets the YY component of a tensor.

CompYZ
[expression]

Gets the YZ component of a tensor.

CompZX
[expression]

Gets the ZX component of a tensor.

CompZY
[expression]

Gets the ZY component of a tensor.

CompZZ
[expression]

Gets the ZZ component of a tensor.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.2.5 Coordinate functions

coord-function-id:

X
[]

Gets the X coordinate.

Y
[]

Gets the Y coordinate.

Z
[]

Gets the Z coordinate.

XYZ
[]

Gets X, Y and Z in a vector.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.2.6 Miscellaneous functions

misc-function-id:

Printf
[expression]

Prints the value of expression when evaluated.

Normal
[]

Computes the normal to the element.

NormalSource
[]

Computes the normal to the source element (only valid in a quantity of Integral type).

F_CompElementNum
[]

Returns 0 if the current element and the current source element are identical.

InterpolationLinear
[]{expression-cst-list}

Linear interpolation of points. The number of constant expressions in expression-cst-list must be even.

dInterpolationLinear
[]{expression-cst-list}

Derivative of linear interpolation of points. The number of constant expressions in expression-cst-list must be even.

InterpolationAkima
[]{expression-cst-list}

Akima interpolation of points. The number of constant expressions in expression-cst-list must be even.

dInterpolationAkima
[]{expression-cst-list}

Derivative of Akima interpolation of points. The number of constant expressions in expression-cst-list must be even.

Order
[quantity]

Returns the interpolation order of the quantity.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.3 Types for Constraint

constraint-type:

Assign
To assign a value (e.g., for boundary condition).

Init
To give an initial value (e.g., initial value in a time analysis).

AssignFromResolution
To assign a value to be computed by a pre-resolution.

InitFromResolution
To give an initial value to be computed by a pre-resolution.

Network
To describe the node connections of branches in a network.

Link
To define links between degrees of freedom in the constrained region with degrees of freedom in a "reference" region, with some coefficient. For example, to link the degrees of freedom in the contrained region Left with the degrees of freedom in the reference region Right, located 1 unit to the right of the region Left along the X-axis, with the coeficient -1, one could write:

 
{ Name periodic;
  Case {
    { Region Left; Type Link ; RegionRef Right;
      Coefficient -1; Function Vector[$X+Pi,$Y,$Z] ;
    }
  }
}

In this example, Function defines the mapping that translates the geometrical elements in the region Left by 1 along the X-axis, so that they correspond with the elements in the region Right. For this mapping to work, the meshes of Left and Right must be identical.

LinkCplx
To define complex-valued links between degrees of freedom. The syntax is the same as for constraints of type Link, but Coeficient can be complex.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.4 Types for FunctionSpace

function-space-type:

Form0
0-form, i.e., scalar field of potential type.

Form1
1-form, i.e., curl-conform field (associated with a curl).

Form2
2-form, i.e., div-conform field (associated with a divergence).

Form3
3-form, i.e., scalar field of density type.

Form1P
1-form perpendicular to the z=0 plane, i.e., perpendicular curl-conform field (associated with a curl).

Form2P
2-form in the z=0 plane, i.e., parallel div-conform field (associated with a divergence).

Scalar
Scalar field.

Vector
Vector field.

basis-function-type:

BF_Node
Nodal function (on NodesOf, value Form0).

BF_Edge
Edge function (on EdgesOf, value Form1).

BF_Facet
Facet function (on FacetsOf, value Form2).

BF_Volume
Volume function (on VolumesOf, value Form3).

BF_GradNode
Gradient of nodal function (on NodesOf, value Form1).

BF_CurlEdge
Curl of edge function (on EdgesOf, value Form2).

BF_DivFacet
Divergence of facet function (on FacetsOf, value Form3).

BF_GroupOfNodes
Sum of nodal functions (on GroupsOfNodesOf, value Form0).

BF_GradGroupOfNodes
Gradient of sum of nodal functions (on GroupsOfNodesOf, value Form1).

BF_GroupOfEdges
Sum of edge functions (on GroupsOfEdgesOf, value Form1).

BF_CurlGroupOfEdges
Curl of sum of edge functions (on GroupsOfEdgesOf, value Form2).

BF_PerpendicularEdge
1-form (0, 0, BF_Node) (on NodesOf, value Form1P).

BF_CurlPerpendicularEdge
Curl of 1-form (0, 0, BF_Node) (on NodesOf, value Form2P).

BF_GroupOfPerpendicularEdge
Sum of 1-forms (0, 0, BF_Node) (on NodesOf, value Form1P).

BF_CurlGroupOfPerpendicularEdge
Curl of sum of 1-forms (0, 0, BF_Node) (on NodesOf, value Form2P).

BF_PerpendicularFacet
2-form (90 degree rotation of BF_Edge) (on EdgesOf, value Form2P).

BF_DivPerpendicularFacet
Div of 2-form (90 degree rotation of BF_Edge) (on EdgesOf, value Form3).

BF_Region
Unit value 1 (on Region, value Scalar).

BF_RegionX
Unit vector (1, 0, 0) (on Region, value Vector).

BF_RegionY
Unit vector (0, 1, 0) (on Region, value Vector).

BF_RegionZ
Unit vector (0, 0, 1) (on Region, value Vector).

BF_Global
Global pre-computed quantity (on Global, value depends on parameters).

BF_dGlobal
Exterior derivative of global pre-computed quantity (on Global, value depends on parameters).

BF_NodeX
Vector (BF_Node, 0, 0) (on NodesOf, value Vector).

BF_NodeY
Vector (0, BF_Node, 0) (on NodesOf, value Vector).

BF_NodeZ
Vector (0, 0, BF_Node) (on NodesOf, value Vector).

BF_Zero
Zero value 0 (on all regions, value Scalar).

BF_One
Unit value 1 (on all regions, value Scalar).

global-quantity-type:

AliasOf
Another name for a name of coefficient of basis function.

AssociatedWith
A global quantity associated with a name of coefficient of basis function, and therefore with this basis function.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.5 Types for Jacobian

jacobian-type:

Vol
Volume Jacobian, for n-D regions in n-D geometries, n = 1, 2 or 3.

Sur
Surface Jacobian, for (n-1)-D regions in n-D geometries, n = 1, 2 or 3.

Lin
Line Jacobian, for (n-2)-D regions in n-D geometries, n = 2 or 3.

VolAxi
Axisymmetrical volume Jacobian (1st type: r), for 2-D regions in axisymmetrical geometries.

SurAxi
Axisymmetrical surface Jacobian (1st type: r), for 1-D regions in axisymmetrical geometries.

VolAxiSqu
Axisymmetrical volume Jacobian (2nd type: r^2), for 2-D regions in axisymmetrical geometries.

VolSphShell
Volume Jacobian with spherical shell transformation, for n-D regions in n-D geometries, n = 2 or 3.

Parameters: radius-internal, radius-external <, center-X, center-Y, center-Z, power, 1/infinity >.

VolAxiSphShell
Same as VolAxi, but with spherical shell transformation.

Parameters: radius-internal, radius-external <, center-X, center-Y, center-Z, power, 1/infinity >.

VolAxiSquSphShell
Same as VolAxiSqu, but with spherical shell transformation.

Parameters: radius-internal, radius-external <, center-X, center-Y, center-Z, power, 1/infinity >.

VolRectShell
Volume Jacobian with rectangular shell transformation, for n-D regions in n-D geometries, n = 2 or 3.

Parameters: radius-internal, radius-external <, direction, center-X, center-Y, center-Z, power, 1/infinity >.

VolAxiRectShell
Same as VolAxi, but with rectangular shell transformation.

Parameters: radius-internal, radius-external <, direction, center-X, center-Y, center-Z, power, 1/infinity >.

VolAxiSquRectShell
Same as VolAxiSqu, but with rectangular shell transformation.

Parameters: radius-internal, radius-external <, direction, center-X, center-Y, center-Z, power, 1/infinity >.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.6 Types for Integration

integration-type:

Gauss
Numerical Gauss integration.

GaussLegendre
Numerical Gauss integration obtained by application of a multiplicative rule on the one-dimensional Gauss integration.

element-type:

Line
Line (2 nodes, 1 edge, 1 volume) (#1).

Triangle
Triangle (3 nodes, 3 edges, 1 facet, 1 volume) (#2).

Quadrangle
Quadrangle (4 nodes, 4 edges, 1 facet, 1 volume) (#3).

Tetrahedron
Tetrahedron (4 nodes, 6 edges, 4 facets, 1 volume) (#4).

Hexahedron
Hexahedron (8 nodes, 12 edges, 6 facets, 1 volume) (#5).

Prism
Prism (6 nodes, 9 edges, 5 facets, 1 volume) (#6).

Pyramid
Pyramid (5 nodes, 8 edges, 5 facets, 1 volume) (#7).

Point
Point (1 node) (#15).

Note:

  1. n in (#n) is the type number of the element (see section 8.1 Input file format).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.7 Types for Formulation

formulation-type:

FemEquation
Finite element method formulation (all methods of moments, integral methods).

local-term-type:

Galerkin
Integral of Galerkin type.

deRham
deRham projection (collocation).

quantity-type:

Local
Local quantity defining a field in a function space. In case a subspace is considered, its identifier has to be given between the brackets following the NameOfSpace function-space-id.

Global
Global quantity defining a global quantity from a function space. The identifier of this quantity has to be given between the brackets following the NameOfSpace function-space-id.

Integral
Integral quantity obtained by the integration of a LocalQuantity before its use in an Equation term.

term-op-type:

Dt
Time derivative applied to the whole term of the equation.

DtDof
Time derivative applied only to the Dof{} term of the equation.

DtDt
Time derivative of 2nd order applied to the whole term of the equation.

DtDtDof
Time derivative of 2nd order applied only to the Dof{} term of the equation.

JacNL
Jacobian term to be assembled in the Jacobian matrix for nonlinear analysis.

NeverDt
No time scheme applied to the term (e.g., Theta is always 1 even if a theta scheme is applied).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.8 Types for Resolution

resolution-op:

Generate
[system-id]

Generate the system of equations system-id.

Solve
[system-id]

Solve the system of equations system-id.

GenerateJac
[system-id]

Generate the system of equations system-id using a jacobian matrix (of which the unknowns are corrections dx of the current solution x).

SolveJac
[system-id]

Solve the system of equations system-id using a jacobian matrix (of which the unknowns are corrections dx of the current solution x). Then, Increment the solution (x=x+dx) and compute the relative error dx/x.

GenerateSeparate
[system-id]

Generate iteration matrices separately for system system-id. It is destined to be used with Update in order to create more efficiently the actual system to solve (this is only useful in linear transient problems with one single excitation) or with EigenSolve in order to generate the matrices of a (generalized) eigenvalue problem.

GenerateOnly
[system-id, expression-cst-list]

Not documented yet.

GenerateOnlyJac
[system-id, expression-cst-list]

Not documented yet.

Update
[system-id, expression]

Update the system of equations system-id (built from iteration matrices generated separately with GenerateSeparate) with expression

UpdateConstraint
[system-id, group-id, constraint-type]

Not documented yet.

InitSolution
[system-id]

Initialize the solution of system-id to zero (default) or to the values given in a Constraint of Init type.

SaveSolution
[system-id]

Save the solution of the system of equations system-id.

SaveSolutions
[system-id]

Save all the solutions available for the system of equations system-id. This should be used with algorithms that generate more than one solution at once, e.g., EigenSolve or FourierTransform.

TransferSolution
[system-id]

Transfer the solution of system system-id, as an Assign constraint, to the system of equations defined with a DestinationSystem command. This is used with the AssignFromResolution constraint type (see section 4.3 Types for Constraint).

TransferInitSolution
[system-id]

Transfer the solution of system system-id, as an Init constraint, to the system of equations defined with a DestinationSystem command. This is used with the InitFromResolution constraint type (see section 4.3 Types for Constraint).

Evaluate
[expression]

Evaluate expression.

SetTime
[expression]

Change the current time.

SetFrequency
[system-id, expression]

Change the frequency of system system-id.

SystemCommand
[expression-char]

Execute the system command given by expression-char.

If
[expression] { resolution-op }

If expression is true (nonzero), perform the operations in resolution-op.

If
[expression] { resolution-op }
Else
{ resolution-op }

If expression is true (nonzero), perform the operations in the first resolution-op, else perform the operations in the second resolution-op.

Break

Not implemented yet.

Print
[ { expression-list }, < File expression-char > ]

Print the expressions listed in expression-list.

Print
[ system-id <, File expression-char > <, { expression-cst-list } >
<, TimeStep { expression-cst-list } >]

Print the system system-id. If the expression-cst-list is given, print only the values of the degrees of freedom given in that list. If the TimeStep option is present, limit the printing to the selected time steps.

EigenSolve
[system-id, expression-cst, expression-cst, expression-cst]

Eigenvalue/eigenvector computation using Arpack. The parameters are: the system (which has to be generated with GenerateSeparate[]), the number of eigenvalues/eigenvectors to compute and the real and imaginary spectral shift (around which to look for eigenvalues).

Lanczos
[system-id, expression-cst, { expression-cst-list } , expression-cst]

Eigenvalue/eigenvector computation using the Lanczos algorithm. The parameters are: the system (which has to be generated with GenerateSeparate[]), the size of the Lanczos space, the indices of the eigenvalues/eigenvectors to store, the spectral shift. This routine is deprecated: use EigenSolve instead.

FourierTransform
[system-id, system-id, { expression-cst-list }]

On-the-fly computation of a discrete Fourier transform. The parameters are: the (time domain) system, the destination system in which the result of the Fourier tranform is to be saved (it should be declared with Type Complex), the list of frequencies to consider in the discrete Fourier transform.

TimeLoopTheta
[expression-cst,expression-cst,expression,expression-cst] { resolution-op }

Time loop of a theta scheme. The parameters are: the initial time, the end time, the time step and the theta parameter (e.g., 1 for implicit Euler, 0.5 for Crank-Nicholson).

TimeLoopNewmark
[expression-cst,expression-cst,expression,expression-cst,expression-cst]
{ resolution-op }

Time loop of a Newmark scheme. The parameters are: the initial time, the end time, the time step, the beta and the gamma parameter.

IterativeLoop
[expression-cst,expression,expression-cst<,expression-cst>] { resolution-op }

Iterative loop for nonlinear analysis. The parameters are: the maximum number of iterations (if no convergence), the relaxation factor (multiplies the iterative correction dx) and the relative error to achieve. The optional parameter is a flag for testing purposes.

PostOperation
[post-operation-id]

Perform the specified PostOperation.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.9 Types for PostProcessing

post-value:

Local
{ local-value }

To compute a local quantity.

Integral
{ integral-value }

To integrate the expression over each element.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

4.10 Types for PostOperation

print-support:

OnElementsOf
group-def

To compute a quantity on the elements belonging to the region group-def, where the solution was computed during the processing stage.

OnRegion
group-def

To compute a global quantity associated with the region group-def.

OnGlobal
To compute a global integral quantity, with no associated region.

OnSection
{ { expression-cst-list } { expression-cst-list } { expression-cst-list } }

To compute a quantity on a section of the mesh defined by three points (i.e., on the intersection of the mesh with a cutting a plane, specified by three points). Each expression-cst-list must contain exactly three elements (the coordinates of the points).

OnGrid
group-def

To compute a quantity in elements of a mesh which differs from the real support of the solution. OnGrid group-def differs from OnElementsOf group-def by the reinterpolation that must be performed.

OnGrid
{ expression, expression, expression }
{ expression-cst-list-item | { expression-cst-list } ,
expression-cst-list-item | { expression-cst-list } ,
expression-cst-list-item | { expression-cst-list } }

To compute a quantity on a parametric grid. The three expressions represent the three cartesian coordinates x, y and z, and can be functions of the current values $A, $B and $C. The values for $A, $B and $C are specified by each expression-cst-list-item or expression-cst-list. For example, OnGrid {Cos[$A], Sin[$A], 0} { 0:2*Pi:Pi/180, 0, 0 } will compute the quantity on 360 points equally distributed on a circle in the z=0 plane, and centered on the origin.

OnPoint
{ expression-cst-list }

To compute a quantity at a point. The expression-cst-list must contain exactly three elements (the coordinates of the point).

OnLine
{ { expression-cst-list } { expression-cst-list } } { expression-cst }

To compute a quantity along a line (given by its two end points), with an associated number of divisions equal to expression-cst. The interpolation points on the line are equidistant. Each expression-cst-list must contain exactly three elements (the coordinates of the points).

OnPlane
{ { expression-cst-list } { expression-cst-list } { expression-cst-list } }
{ expression-cst, expression-cst }

To compute a quantity on a plane (specified by three points), with an associated number of divisions equal to each expression-cst along both generating directions. Each expression-cst-list must contain exactly three elements (the coordinates of the points).

OnBox
{ { expression-cst-list } { expression-cst-list } { expression-cst-list }
{ expression-cst-list } } { expression-cst, expression-cst, expression-cst }

To compute a quantity in a box (specified by four points), with an associated number of divisions equal to each expression-cst along the three generating directions. Each expression-cst-list must contain exactly three elements (the coordinates of the points).

print-option:

File
expression-char

Outputs the result in a file named expression-char.

File
> expression-char

Same as File expression-char, except that, if several File > expression-char options appear in the same PostOperation, the results are concatenated in the file expression-char.

File
>> expression-char

Appends the result to a file named expression-char.

Depth
expression-cst

Recursive division of the elements if expression-cst is greater than zero, derefinement if expression-cst is smaller than zero. If expression-cst is equal to zero, evaluation at the barycenter of the elements.

Skin
Computes the result on the boundary of the region.

Smoothing
Smoothes the solution at the nodes.

HarmonicToTime
expression-cst

Converts a harmonic solution into a time-dependent one (with expression-cst steps).

Dimension
expression-cst

Forces the dimension of the elements to consider in an element search. Specifies the problem dimension during an adaptation (h- or p-refinement).

TimeStep
expression-cst-list-item | { expression-cst-list }

Outputs results for the specified time steps only.

LastTimeStepOnly

Outputs results for the last time step only (useful when calling a PostOperation directly in a Resolution, for example).

Frequency
expression-cst-list-item | { expression-cst-list }

Outputs results for the specified frequencies only.

Format
post-operation-fmt

Outputs results in the specified format.

Adapt
P1 | H1 | H2

Performs p- or h-refinement on the post-processing result, considered as an error map.

Target
expression-cst

Specifies the target for the optimizer during adaptation (error for P1|H1, number of elements for H2).

Value
expression-cst-list-item | { expression-cst-list }

Specifies acceptable output values for discrete optimization (e.g. the available interpolation orders with Adapt P1).

Sort
Position | Connection

Sorts the output by position (x, y, z) or by connection (for LINE elements only).

Iso
expression-cst

Outputs directly contour prints (with expression-cst values) instead of elementary values.

Iso
{ expression-cst-list }

Outputs directly contour prints for the values specified in the expression-cst-list instead of elementary values.

NoNewLine
Suppresses the new lines in the output when printing global quantities (i.e., with Print OnRegion or Print OnGlobal).

ChangeOfCoordinates
{ expression, expression, expression }

Changes the coordinates of the results according to the three expressions given in argument. The three expressions represent the three new cartesian coordinates x, y and z, and can be functions of the current values of the cartesian coordinates $X, $Y and $Z.

ChangeOfValues
{ expression-list }

Changes the values of the results according to the expressions given in argument. The expressions represent the new values (x-compoment, y-component, etc.), and can be functions of the current values of the solution ($Val0, $Val1, etc.).

DecomposeInSimplex

Decomposes all output elements in simplices (points, lines, triangles or tetrahedra).

Store
expression-cst

Stores the result of an OnRegion post-processing operation in the register expression-cst.

TimeLegend
< { expression, expression, expression } >

Includes a time legend in Gmsh plots. If the three optional expressions giving the position of the legend are not specified, the legend is centered on top of the plot.

FrequencyLegend
< { expression, expression, expression } >

Includes a frequency legend in Gmsh plots. If the three optional expressions giving the position of the legend are not specified, the legend is centered on top of the plot.

EigenvalueLegend
< { expression, expression, expression } >

Includes an eigenvalue legend in Gmsh plots. If the three optional expressions giving the position of the legend are not specified, the legend is centered on top of the plot.

post-operation-fmt:

Gmsh
GmshParsed
Gmsh output. See the documentation of Gmsh (http://www.geuz.org/gmsh/) for a description of the file formats.

Table
Space oriented column output, e.g., suitable for Gnuplot, Excel, Caleida Graph, etc. The columns are: element-type element-index x-coord y-coord z-coord <x-coord y-coord z-coord> ... real real real values. The three real numbers preceding the values contain context-dependent information, depending on the type of plot: curvilinear abscissa for OnLine plots, normal to the plane for OnPlane plots, parametric coordinates for parametric OnGrid plots, etc.

TimeTable
Time oriented column output, e.g., suitable for Gnuplot, Excel, Caleida Graph, etc. The columns are: time-step time x-coord y-coord z-coord <x-coord y-coord z-coord> ... value.

Gnuplot
Space oriented column output similar to the Table format, except that a new line is created for each node of each element, with a repetition of the first node if the number of nodes in the element is greater than 2. This permits to draw unstructured meshes and nice three-dimensional elevation plots in Gnuplot. The columns are: element-type element-index x-coord y-coord z-coord real real real values. The three real numbers preceding the values contain context-dependent information, depending on the type of plot: curvilinear abscissa for OnLine plots, normal to the plane for OnPlane plots, parametric coordinates for parametric OnGrid plots, etc.

Adaptation
Adaptation map, suitable for the GetDP -adapt command line option.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5. Short examples

5.1 Constant expression examples  
5.2 Group examples  
5.3 Function examples  
5.4 Constraint examples  
5.5 FunctionSpace examples  
5.6 Jacobian examples  
5.7 Integration examples  
5.8 Formulation examples  
5.9 Resolution examples  
5.10 PostProcessing examples  
5.11 PostOperation examples  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.1 Constant expression examples

The simplest constant expression consists of an integer or a real number as in

 
21
-3

or

 
-3.1415
27e3
-290.53e-12

Using operators and the classic math functions, constant-ids can be defined:

 
c1 = Sin[2/3*3.1415] * 5000^2;
c2 = -1/c1;


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.2 Group examples

Let us assume that some elements in the input mesh have the region numbers 1000, 2000 and 3000. In the definitions

 
Group {
  Air = Region[1000]; Core = Region[2000]; Inductor = Region[3000];
  NonConductingDomain = Region[{Air, Core}];
  ConductingDomain    = Region[{Inductor}];
}

Air, Core, Inductor are identifiers of elementary region groups while NonConductingDomain and ConductingDomain are global region groups.

Groups of function type contain lists of entities built on the region groups appearing in their arguments. For example,

 
NodesOf[NonConductingDomain]

represents the group of nodes of geometrical elements belonging to the regions in NonConductingDomain and

 
EdgesOf[DomainC, Not SkinDomainC]

represents the group of edges of geometrical elements belonging to the regions in DomainC but not to those of SkinDomainC.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.3 Function examples

A physical characteristic is a piecewise defined function. The magnetic permeability mu[] can for example be defined in the considered regions by

 
Function {
  mu[Air] = 4.e-7*Pi; 
  mu[Core] = 1000.*4.e-7*Pi;
}

A nonlinear characteristic can be defined through an expression with arguments, e.g.,

 
Function {
  mu0 = 4.e-7*Pi;
  a1 = 1000.; b1 = 100.; // Constants
  mu[NonlinearCore] = mu0 + 1./(a1+b1*Norm[$1]^6);
}

where function mu[] in region NonLinearCore has one argument $1 which has to be the magnetic flux density. This function is actually called when writing the equations of a formulation, which permits to directly extend it to a nonlinear form by adding only the necessary arguments. For example, in a magnetic vector potential formulation, one may write mu[{Curl a}] instead of mu[] in Equation terms (see section 5.8 Formulation examples). Multiple arguments can be specified in a similar way: writing mu[{Curl a},{T}] in an Equation term will provide the function mu[] with two usable arguments, $1 (the magnetic flux density) and $2 (the temperature).

It is also possible to directly interpolate one-dimensional functions from tabulated data. In the following example, the function f(x) as well as its derivative f'(x) are interpolated from the (x,f(x)) couples (0,0.65), (1,0.72), (2,0.98) and (3,1.12):

 
Function {
  couples = {0, 0.65 , 1, 0.72 , 2, 0.98 , 3, 1.12};
  f[] = InterpolationLinear[$1]{List[couples]};
  dfdx[] = dInterpolationLinear[$1]{List[couples]};
}

The function f[] may then be called in an Equation term of a Formulation with one argument, x. Notice how the list of constants List[couples] is supplied as a list of parameters to the built-in function InterpolationLinear (see section 2.2 Constants, as well as 2.4 Functions). In order to facilitate the construction of such interpolations, the couples can also be specified in two separate lists, merged with the alternate list ListAlt command (see section 2.2 Constants):

 
Function {
  data_x = {0, 1, 2, 3};
  data_f = {0.65, 0.72, 0.98, 1.12};
  f[] = InterpolationLinear[$1]{ListAlt[data_x, data_f]};
  dfdx[] = dInterpolationLinear[$1]{ListAlt[data_x, data_f]};
}

In order to optimize the evaluation time of complex expressions, registers may be used (see section 2.7 Registers). For example, the evaluation of g[] = f[$1]*Sin[f[$1]^2] would require two (costly) linear interpolations. But the result of the evaluation of f[] may be stored in a register (for example the register 0) with

 
g[] = f[$1]#0 * Sin[#0^2];

thus reducing the number of evaluations of f[] (and of the argument $1) to one.

A function can also be time dependent, e.g.,

 
Function {
  Freq = 50.; Phase = 30./180.*Pi; // Constants
  TimeFct_Sin[] = Sin [ 2.*Pi*Freq * $Time + Phase ];
  TimeFct_Exp[] = Exp [ - $Time / 0.0119 ];
  TimeFct_ExtSin[] = F_Sin_wt_p [] {2.*Pi*Freq, Phase};
}

Note that TimeFct_ExtSin[] is identical to TimeFct_Sin[] in a time domain analysis, but also permits to define phasors implicitely in the case of harmonic analyses.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.4 Constraint examples

Constraints are referred to in FunctionSpaces and are usually used for boundary conditions (Assign type). For example, essential conditions on two surface regions, Surf0 and Surf1, will be first defined by

 
Constraint {
  { Name DirichletBoundaryCondition1; Type Assign;
    Case {
      { Region Surf0; Value 0.; }
      { Region Surf1; Value 1.; }
    }
  }
}

The way the Values are associated with Regions (with their nodes, their edges, their global regions, ...) is defined in the FunctionSpaces which use the Constraint. In other words, a Constraint defines data but does not define the method to process them. A time dependent essential boundary condition on Surf1 would be introduced as (cf. 5.3 Function examples for the definition of TimeFct_Exp[]):

 
      { Region Surf1; Value 1.; TimeFunction 3*TimeFct_Exp[] }

It is important to notice that the time dependence cannot be introduced in the Value field, since the Value is only evaluated once during the pre-processing.

Other constraints can be referred to in Formulations. It is the case of those defining electrical circuit connections (Network type), e.g.,

 
Constraint {
  { Name ElectricalCircuit; Type Network; 
    Case Circuit1 {
      { Region VoltageSource; Branch {1,2}; }
      { Region PrimaryCoil; Branch {1,2}; }
    }
    Case Circuit2 {
      { Region SecondaryCoil; Branch {1,2}; }
      { Region Charge; Branch {1,2}; }
    }
  }
}

which defines two non-connected circuits (Circuit1 and Circuit2), with an independent numbering of nodes: region VoltageSource is connected in parallel with region PrimaryCoil, and region SecondaryCoil is connected in parallel with region Charge.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.5 FunctionSpace examples

Various discrete function spaces can be defined in the frame of the finite element method.

5.5.1 Nodal finite element spaces  
5.5.2 High order nodal finite element space  
5.5.3 Nodal finite element space with floating potentials  
5.5.4 Edge finite element space  
5.5.5 Edge finite element space with gauge condition  
5.5.6 Coupled edge and nodal finite element spaces  
5.5.7 Coupled edge and nodal finite element spaces for multiply connected domains  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.5.1 Nodal finite element spaces

The most elementary function space is the nodal finite element space, defined on a mesh of a domain W and denoted S0(W) (associated finite elements can be of various geometries), and associated with essential boundary conditions (Dirichlet conditions). It contains 0-forms, i.e., scalar fields of potential type:

v = Sum [ vn * sn, for all n in N ], v in S0(W)

where N is the set of nodes of W, sn is the nodal basis function associated with node n and vn is the value of v at node n. It is defined by

 
FunctionSpace {
  { Name Hgrad_v; Type Form0;
    BasisFunction {
      { Name sn; NameOfCoef vn; Function BF_Node;
        Support Domain; Entity NodesOf[All]; }
    }
    Constraint {
      { NameOfCoef vn; EntityType NodesOf;
        NameOfConstraint DirichletBoundaryCondition1; }
    }
  }
}

Function sn is the built-in basis function BF_Node associated with all nodes (NodesOf) in the mesh of W (Domain). Previously defined Constraint DirichletBoundaryCondition1 (see section 5.4 Constraint examples) is used as boundary condition.

In the example above, Entity NodesOf[All] is preferred to Entity NodesOf[Domain]. In this way, the list of all the nodes of Domain will not have to be generated. All the nodes of each geometrical element in Support Domain will be directly taken into account.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.5.2 High order nodal finite element space

Higher order finite elements can be directly taken into account by BF_Node. Hierarchical finite elements for 0-forms can be used by simply adding other basis functions (associated with other geometrical entities, e.g., edges and facets) to BasisFunction, e.g.,

 
    ...
    BasisFunction {
      { Name sn; NameOfCoef vn; Function BF_Node;
        Support Domain; Entity NodesOf[All]; }
      { Name s2; NameOfCoef v2; Function BF_Node_2E;
        Support Domain; Entity EdgesOf[All]; }
    }
    ...


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.5.3 Nodal finite element space with floating potentials

A scalar potential with floating values vf on certain boundaries Gf, f in Cf, e.g., for electrostatic problems, can be expressed as

v = Sum [ vn * sn, for all n in Nv ] + Sum [ vf * sf, for all f in Cf ], v in S0(W)

where Nv is the set of inner nodes of W and each function sf is associated with the group of nodes of boundary Gf, f in Cf (SkinDomainC); sf is the sum of the nodal basis functions of all the nodes of Cf. Its function space is defined by

 
FunctionSpace {
  { Name Hgrad_v_floating; Type Form0;
    BasisFunction {
      { Name sn; NameOfCoef vn; Function BF_Node;
        Support Domain; Entity NodesOf[All, Not SkinDomainC]; }
      { Name sf; NameOfCoef vf; Function BF_GroupOfNodes; 
        Support Domain; Entity GroupsOfNodesOf[SkinDomainC]; }
    }
    GlobalQuantity {
      { Name GlobalElectricPotential; Type AliasOf; NameOfCoef vf; }
      { Name GlobalElectricCharge; Type AssociatedWith; 
        NameOfCoef vf; }
    }
    Constraint { ... }
  }
}

Two global quantities have been associated with this space: the electric potential GlobalElectricPotential, being an alias of coefficient vf, and the electric charge GlobalElectricCharge, being associated with coefficient vf.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.5.4 Edge finite element space

Another space is the edge finite element space, denoted S1(W), containing 1-forms, i.e., curl-conform fields:

h = Sum [ he * se, for all e in E ], h in S1(W)

where E is the set of edges of W, se is the edge basis function for edge e and he is the circulation of h along edge e. It is defined by

 
FunctionSpace {
  { Name Hcurl_h; Type Form1;
    BasisFunction {
      { Name se; NameOfCoef he; Function BF_Edge;
        Support Domain; Entity EdgesOf[All]; }
    }
    Constraint { ... }
  }
}


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.5.5 Edge finite element space with gauge condition

A 1-form function space containing vector potentials can be associated with a gauge condition, which can be defined as a constraint, e.g., a zero value is fixed for all circulations along edges of a tree (EdgesOfTreeIn) built in the mesh (Domain), having to be complete on certain boundaries (StartingOn Surf):

 
Constraint {
  { Name GaugeCondition_a_Mag_3D; Type Assign;
    Case {
      { Region Domain; SubRegion Surf; Value 0.; }
    }
  }
}

FunctionSpace {
  { Name Hcurl_a_Gauge; Type Form1;
    BasisFunction {
      { Name se; NameOfCoef ae; Function BF_Edge;
        Support Domain; Entity EdgesOf[All]; }
    }
    Constraint {
      { NameOfCoef ae;
        EntityType EdgesOfTreeIn; EntitySubType StartingOn;
        NameOfConstraint GaugeCondition_a_Mag_3D; }
      ...
    }
  }
}


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.5.6 Coupled edge and nodal finite element spaces

A 1-form function space, containing curl free fields in certain regions WcC (DomainCC) of W, which are the complementary part of Wc (DomainC) in W, can be explicitly characterized by

h = Sum [ hk * sk, for all e in Ec ] + Sum [ phin * vn, for all n in NcC ], h in S1(W)

where Ec is the set of inner edges of W, NcC is the set of nodes inside WcC and on its boundary dWcC, sk is an edge basis function and vn is a vector nodal function. Such a space, coupling a vector field with a scalar potential, can be defined by

 
FunctionSpace {
  { Name Hcurl_hphi; Type Form1;
    BasisFunction {
      { Name sk; NameOfCoef hk; Function BF_Edge; 
        Support DomainC; Entity EdgesOf[All, Not SkinDomainC]; }
      { Name vn; NameOfCoef phin; Function BF_GradNode; 
        Support DomainCC; Entity NodesOf[All]; }
      { Name vn; NameOfCoef phic; Function BF_GroupOfEdges; 
        Support DomainC; Entity GroupsOfEdgesOnNodesOf[SkinDomainC];}
    }
    Constraint {
      { NameOfCoef hk;
        EntityType EdgesOf; NameOfConstraint MagneticField; }
      { NameOfCoef phin;
        EntityType NodesOf; NameOfConstraint MagneticScalarPotential; }
      { NameOfCoef phic;
        EntityType NodesOf; NameOfConstraint MagneticScalarPotential; }
    }
  }
}

This example points out the definition of a piecewise defined basis function, e.g., function vn being defined with BF_GradNode in DomainCC and BF_GroupOfEdges in DomainC. This leads to an easy coupling between these regions.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.5.7 Coupled edge and nodal finite element spaces for multiply connected domains

In case a multiply connected domain WcC is considered, basis functions associated with cuts (SurfaceCut) have to be added to the previous basis functions, which gives the function space below:

 
Group {
  _TransitionLayer_SkinDomainC_ =
    ElementsOf[SkinDomainC, OnOneSideOf SurfaceCut];
}

FunctionSpace {
  { Name Hcurl_hphi; Type Form1;
    BasisFunction {

      ... same as above ...

      { Name sc; NameOfCoef Ic; Function BF_GradGroupOfNodes;
        Support ElementsOf[DomainCC, OnOneSideOf SurfaceCut];
        Entity GroupsOfNodesOf[SurfaceCut]; }
      { Name sc; NameOfCoef Icc; Function BF_GroupOfEdges;
        Support DomainC;
        Entity GroupsOfEdgesOf
                 [SurfaceCut,
                  InSupport _TransitionLayer_SkinDomainC_]; }
    }
    GlobalQuantity {
      { Name I; Type AliasOf       ; NameOfCoef Ic; }
      { Name U; Type AssociatedWith; NameOfCoef Ic; }
    }
    Constraint {

      ... same as above ...

      { NameOfCoef Ic;
        EntityType GroupsOfNodesOf; NameOfConstraint Current; }
      { NameOfCoef Icc;
        EntityType GroupsOfNodesOf; NameOfConstraint Current; }
      { NameOfCoef U;
        EntityType GroupsOfNodesOf; NameOfConstraint Voltage; }
    }
  }
}

Global quantities associated with the cuts, i.e., currents and voltages if h is the magnetic field, have also been defined.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.6 Jacobian examples

A simple Jacobian method is for volume transformations (of n-D regions in n-D geometries; n = 1, 2 or 3), e.g., in region Domain,

 
Jacobian {
  { Name Vol;
    Case {
      { Region Domain; Jacobian Vol; }
    }
  }
}

Jacobian VolAxi would define a volume Jacobian for axisymmetrical problems.

A Jacobian method can also be piecewise defined, in DomainInf, where an infinite geometrical transformation has to be made using two constant parameters (inner and outer radius of a spherical shell), and in all the other regions (All, being the default); in each case, a volume Jacobian is used. This method is defined by:

 
Jacobian {
  { Name Vol;
    Case {
      { Region DomainInf; Jacobian VolSphShell {Val_Rint, Val_Rext}; }
      { Region All; Jacobian Vol; }
    }
  }
}


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.7 Integration examples

A commonly used numerical integration method is the Gauss integration, with a number of integration points (NumberOfPoints) depending on geometrical element types (GeoElement), i.e.

 
Integration {
  { Name Int_1;
    Case { {Type Gauss;
            Case { { GeoElement Triangle   ; NumberOfPoints 4; }
                   { GeoElement Quadrangle ; NumberOfPoints 4; }
                   { GeoElement Tetrahedron; NumberOfPoints 4; }
                   { GeoElement Hexahedron ; NumberOfPoints 6; }
                   { GeoElement Prism      ; NumberOfPoints 9; } }
           }
         }
  }
}

The method above is valid for both 2D and 3D problems, for different kinds of elements.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.8 Formulation examples

5.8.1 Electrostatic scalar potential formulation  
5.8.2 Electrostatic scalar potential formulation with floating potentials and electric charges  
5.8.3 Magnetostatic 3D vector potential formulation  
5.8.4 Magnetodynamic 3D or 2D magnetic field and magnetic scalar potential formulation  
5.8.5 Nonlinearities, Mixed formulations, ...  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.8.1 Electrostatic scalar potential formulation

An electrostatic formulation using an electric scalar potential v, i.e.

( epsr grad v, grad vp ) W = 0, for all vp in S0(W),

is expressed by

 
Formulation {
  { Name Electrostatics_v; Type FemEquation;
    Quantity {
      { Name v; Type Local; NameOfSpace Hgrad_v; }
    }
    Equation { 
      Galerkin { [ epsr[] * Dof{Grad v} , {Grad v} ];
                 In Domain; Jacobian Vol; Integration Int_1; }
    }
  }
}

The density of the Galerkin term is a copy of the symbolic form of the formulation, i.e., the product of a relative permittivity function epsr[] by a vector of degrees of freedom (Dof{.}); the scalar product of this with the gradient of test function v results in a symmetrical matrix.

Note that another Quantity could be defined for test functions, e.g., vp defined by { Name vp; Type Local; NameOfSpace Hgrad_v; }. However, its use would result in the computation of a full matrix and consequently in a loss of efficiency.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.8.2 Electrostatic scalar potential formulation with floating potentials and electric charges

An extension of the formulation above can be made to take floating potentials and electrical charges into account (the latter being defined in FunctionSpace Hgrad_v_floating), i.e.

 
Formulation {
  { Name Electrostatics_v_floating; Type FemEquation;
    Quantity {
      { Name v; Type Local; NameOfSpace Hgrad_v_floating; }
      { Name V; Type Global; 
        NameOfSpace Hgrad_v_floating [GlobalElectricPotential]; }
      { Name Q; Type Global; 
        NameOfSpace Hgrad_v_floating [GlobalElectricCharge]; }
    }
    Equation { 
      Galerkin { [ epsr[] * Dof{Grad v} , {Grad v} ];
                 In Domain; Jacobian Vol; Integration Int_1; }
      GlobalTerm { [ - Dof{Q}/eps0 , {V} ]; In SkinDomainC; }
    }
  }
}

with the predefinition Function { eps0 = 8.854187818e-12; }.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.8.3 Magnetostatic 3D vector potential formulation

A magnetostatic 3D vector potential formulation

( nu curl a , curl ap ) W - ( js , ap ) Ws = 0, for all ap in S1(W) with gauge condition,

with a source current density js in inductors Ws, is expressed by

 
Formulation {
  { Name Magnetostatics_a_3D; Type FemEquation;
    Quantity {
      { Name a; Type Local; NameOfSpace Hcurl_a_Gauge; }
    }
    Equation { 
      Galerkin { [ nu[] * Dof{Curl a} , {Curl a} ];
                 In Domain; Jacobian Vol; Integration Int_1; }
      Galerkin { [ - SourceCurrentDensity[] , {a} ];
                 In DomainWithSourceCurrentDensity;
                 Jacobian Vol; Integration Int_1; }
    }
  }
}

Note that js is here given by a function SourceCurrentDensity[], but could also be given by data computed from another problem, e.g., from an electrokinetic problem (coupling of formulations) or from a fully fixed function space (constraints fixing the density, which is usually more efficient in time domain analyses).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.8.4 Magnetodynamic 3D or 2D magnetic field and magnetic scalar potential formulation

A magnetodynamic 3D or 2D h-phi formulation, i.e., coupling the magnetic field h with a magnetic scalar potential phi,

Dt ( mu h , hp ) W + ( ro curl h , curl hp ) Wc = 0, for all hp in S1(W),

can be expressed by

 
Formulation {
  { Name Magnetodynamics_hphi; Type FemEquation;
    Quantity {
      { Name h; Type Local; NameOfSpace Hcurl_hphi; }
    }
    Equation { 
      Galerkin { Dt [ mu[] * Dof{h} , {h} ];
                 In Domain; Jacobian Vol; Integration Int_1; }
      Galerkin { [ rho[] * Dof{Curl h} , {Curl h} ];
                 In DomainC; Jacobian Vol; Integration Int_1; }
    }
  }
}


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.8.5 Nonlinearities, Mixed formulations, ...

In case nonlinear physical characteristics are considered, arguments are used for associated functions, e.g., mu[{h}]. Several test functions can be considered in an Equation field. Consequently, mixed formulations can be defined.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.9 Resolution examples

5.9.1 Static resolution (electrostatic problem)  
5.9.2 Frequency domain resolution (magnetodynamic problem)  
5.9.3 Time domain resolution (magnetodynamic problem)  
5.9.4 Nonlinear time domain resolution (magnetodynamic problem)  
5.9.5 Coupled formulations  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.9.1 Static resolution (electrostatic problem)

A static resolution, e.g., for the electrostatic formulation (see section 5.8 Formulation examples), can be defined by

 
Resolution {
  { Name Electrostatics_v;
    System {
      { Name Sys_Ele; NameOfFormulation Electrostatics_v; }
    }
    Operation { 
      Generate[Sys_Ele]; Solve[Sys_Ele]; SaveSolution[Sys_Ele];
    }
  }
}

The generation (Generate) of the matrix of the system Sys_Ele will be made with the formulation Electrostatics_v, followed by its solving (Solve) and the saving of the solution (SaveSolution).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.9.2 Frequency domain resolution (magnetodynamic problem)

A frequency domain resolution, e.g., for the magnetodynamic h-phi formulation (see section 5.8 Formulation examples), is given by

 
Resolution {
  { Name Magnetodynamics_hphi;
    System {
      { Name Sys_Mag; NameOfFormulation Magnetodynamics_hphi;
        Frequency Freq; }
    }
    Operation {
      Generate[Sys_Mag]; Solve[Sys_Mag]; SaveSolution[Sys_Mag];
    }
  }
}

preceded by the definition of constant Freq, e.g.,

 
Function {
  Freq = 50.;
}


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.9.3 Time domain resolution (magnetodynamic problem)

A time domain resolution, e.g., for the same magnetodynamic h-phi formulation (see section 5.8 Formulation examples), is given by

 
Resolution {
  { Name Magnetodynamics_hphi_Time;
    System {
      { Name Sys_Mag; NameOfFormulation Magnetodynamics_hphi; }
    }
    Operation {
      InitSolution[Sys_Mag]; SaveSolution[Sys_Mag];
      TimeLoopTheta[Mag_Time0, Mag_TimeMax, Mag_DTime[], Mag_Theta[]] {
        Generate[Sys_Mag]; Solve[Sys_Mag]; SaveSolution[Sys_Mag]; 
      }
    }
  }
}

If, e.g., the Resolution above is preceded by the constant and function definitions below

 
Function {
  Tc = 10.e-3;
  Mag_Time0 = 0.; Mag_TimeMax = 2.*Tc; Mag_DTime[] = Tc/20.;
  Mag_Theta[] = 1./2.;
}

the performed time analysis will be a Crank-Nicolson scheme (theta-scheme with Theta = 0.5) with initial time 0 ms, end time 20 ms and time step 1 ms.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.9.4 Nonlinear time domain resolution (magnetodynamic problem)

In case a nonlinear problem is solved, an iterative loop has to be defined in an appropriate level of the recursive resolution operations, e.g., for the magnetodynamic problem above,

 
...
    Operation {
      InitSolution[Sys_Mag]; SaveSolution[Sys_Mag];
      TimeLoopTheta[Mag_Time0, Mag_TimeMax, Mag_DTime[], Mag_Theta[]] {
        IterativeLoop[NL_NbrMax, NL_Eps, NL_Relax] {
          GenerateJac[Sys_Mag]; SolveJac[Sys_Mag]; 
        }
        SaveSolution[Sys_Mag];
      }
    }
...

preceded by constant definitions, e.g.,

 
Function {
  NL_Eps = 1.e-4; NL_Relax = 1.; NL_NbrMax = 80;
}


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.9.5 Coupled formulations

A coupled problem, e.g., magnetodynamic (in frequency domain; Frequency Freq) - thermal (in time domain) coupling, with temperature dependent characteristics (e.g., rho[{T}], ...), can be defined by:

 
Resolution {
  { Name MagnetoThermalCoupling_hphi_T;
    System {
      { Name Sys_Mag; NameOfFormulation Magnetodynamics_hphi; 
        Frequency Freq; }
      { Name Sys_The; NameOfFormulation Thermal_T; }
    }
    Operation {
      InitSolution[Sys_Mag]; InitSolution[Sys_The];
      IterativeLoop[NL_NbrMax, NL_Eps, NL_Relax] {
        GenerateJac[Sys_Mag]; SolveJac[Sys_Mag];
        GenerateJac[Sys_The]; SolveJac[Sys_The];
      }
      SaveSolution[Sys_Mag]; SaveSolution[Sys_The];
    }
  }
}

Two systems of equations, Sys_Mag and Sys_The, will be solved iteratively until convergence (Criterion), using a relaxation factor (RelaxationFactor).

It can be seen through these examples that many resolutions can be linked or nested directly by the user, which gives a great freedom for coupled problems.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.10 PostProcessing examples

The quantities to be post-computed based on a solution of a resolution are defined, e.g., for the electrostatic problem (see section 5.8 Formulation examples; see section 5.9 Resolution examples), for the solution associated with the formulation Electrostatics_v, by

 
PostProcessing {
  { Name EleSta_v; NameOfFormulation Electrostatics_v;
    Quantity {
      { Name v; Value { Local { [ {v} ]; In Domain; } } }
      { Name e; Value { Local { [ -{Grad v} ]; In Domain; } } }
      { Name d; Value { Local { [ -eps0*epsr[] *{Grad v} ];
                                  In Domain; } } }
    }
  }
}

The electric scalar potential v (v), the electric field e (e) and the electric flux density d (d) can all be computed from the solution. They are all defined in the region Domain.

The quantities for the solution associated with the formulation Electrostatics_v_floating are defined by

 
PostProcessing {
  { Name EleSta_vf; NameOfFormulation Electrostatics_v_floating;
    Quantity {

      ... same as above ...

      { Name Q; Value { Local { [ {Q} ]; In SkinDomainC; } } }
      { Name V; Value { Local { [ {V} ]; In SkinDomainC; } } }
    }
  }
}

which points out the way to define post-quantities based on global quantities.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

5.11 PostOperation examples

The simplest post-processing operation is the generation of maps of local quantities, i.e., the display of the computed fields on the mesh. For example, using the PostProcessing defined in 5.10 PostProcessing examples, the maps of the electric scalar potential and of the electric field on the elements of the region Domain are defined as:

 
PostOperation {
  { Name Map_v_e; NameOfPostProcessing EleSta_v ;
     Operation {
       Print [ v, OnElementsOf Domain, File "map_v.pos" ];
       Print [ e, OnElementsOf Domain, File "map_e.pos" ];
     }
  }
}

It is also possible to display local quantities on sections of the mesh, here for example on the plane containing the points (0,0,1), (1,0,1) and (0,1,1):

 
Print [ v, OnSection { {0,0,1} {1,0,1} {0,1,1} }, File "sec_v.pos" ];

Finally, local quantities can also be interpolated on another mesh than the one on which they have been computed. Six types of grids can be specified for this interpolation: a single point, a set of points evenly distributed on a line, a set of points evenly distributed on a plane, a set of points evenly distributed in a box, a set of points defined by a parametric equation, and a set of elements belonging to a different mesh than the original one:

 
Print [ e, OnPoint {0,0,1} ];
Print [ e, OnLine { {0,0,1} {1,0,1} } {125} ];
Print [ e, OnPlane { {0,0,1} {1,0,1} {0,1,1} } {125, 75} ];
Print [ e, OnBox { {0,0,1} {1,0,1} {0,1,1} {0,0,2} } {125, 75, 85} ];
Print [ e, OnGrid {$A, $B, 1} { 0:1:1/125, 0:1:1/75, 0 } ];
Print [ e, OnGrid Domain2 ];

Many options can be used to modify the aspect of all these maps, as well as the default behaviour of the Print commands. See 4.10 Types for PostOperation, to get the list of all these options. For example, to obtain a map of the scalar potential at the barycenters of the elements on the boundary of the region Domain, in a table oriented format appended to an already existing file out.txt, the operation would be:

 
Print [ v, OnElementsOf Domain, Depth 0, Skin, Format Table, 
        File >> "out.txt" ];

Global quantities, which are associated with regions (and not with the elements of the mesh of these regions), are displayed thanks to the OnRegion operation. For example, the global potential and charge on the region SkinDomainC can be displayed with:

 
PostOperation {
  { Name Val_V_Q; NameOfPostProcessing EleSta_vf ;
     Operation {
       Print [ V, OnRegion SkinDomainC ];
       Print [ Q, OnRegion SkinDomainC ];
     }
  }
}


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

6. Complete examples

This chapter presents complete examples that can be run "as is" with GetDP (see section 7. Running GetDP).

Many other ready-to-use examples are available in the GetDP wiki at the following address: http://www.geuz.org/getdp/wiki/FrontPage (username=getdp; password=wiki).

6.1 Electrostatic problem  
6.2 Magnetostatic problem  
6.3 Magnetodynamic problem  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

6.1 Electrostatic problem

Let us first consider a simple electrostatic problem. The formulation used is an electric scalar potential formulation (file `EleSta_v.pro', including files `Jacobian_Lib.pro' and `Integration_Lib.pro'). It is applied to a microstrip line (file `mStrip.pro'), whose geometry is defined in the file `mStrip.geo' (see section C. Gmsh examples). The geometry is two-dimensional and by symmetry only one half of the structure is modeled.

Strip

Note that the structure of the following files points out the separation of the data describing the particular problem and the method used to solve it (see section 1.1 Numerical tools as objects), and therefore how it is possible to build black boxes adapted to well defined categories of problems. The files are commented (see section 1.3 Comments) and can be run without any modification.

/* -------------------------------------------------------------------
   File "mStrip.pro"

   This file defines the problem dependent data structures for the
   microstrip problem.
   
   To compute the solution: 
       getdp mStrip -solve EleSta_v

   To compute post-results:
       getdp mStrip -pos Map
    or getdp mStrip -pos Cut
   ------------------------------------------------------------------- */

Group {
   
  /* Let's start by defining the interface (i.e. elementary groups)
     between the mesh file and GetDP (no mesh object is defined, so 
     the default mesh will be assumed to be in GMSH format and located 
     in "mStrip.msh") */

  Air = Region[101]; Diel1 = Region[111];
  Ground = Region[120]; Line = Region[121]; 
  SurfInf = Region[130];

  /* We can then define a global group (used in "EleSta_v.pro",
     the file containing the function spaces and formulations) */

  DomainCC_Ele = Region[{Air, Diel1}];

}

Function {

  /* The relative permittivity (needed in the formulation) is piecewise
     defined in elementary groups */

  epsr[Air] = 1.;
  epsr[Diel1] = 9.8;

}

Constraint {

  /* Now, some Dirichlet conditions are defined. The name 
     'ElectricScalarPotential' refers to the constraint name given in
     the function space */

  { Name ElectricScalarPotential; Type Assign;
    Case {
      { Region Region[{Ground, SurfInf}]; Value 0.; }
      { Region Line; Value 1.e-3; }
    }
  }

}

/* The formulation used and its tools, considered as being 
   in a black box, can now be included */

Include "Jacobian_Lib.pro"
Include "Integration_Lib.pro"
Include "EleSta_v.pro"

/* Finally, we can define some operations to output results */

e = 1.e-7;

PostOperation {
  { Name Map; NameOfPostProcessing EleSta_v;
     Operation {
       Print [ v, OnElementsOf DomainCC_Ele, File "mStrip_v.pos" ];
       Print [ e, OnElementsOf DomainCC_Ele, File "mStrip_e.pos" ];
     }
  }
  { Name Cut; NameOfPostProcessing EleSta_v;
     Operation {
       Print [ e, OnLine {{e,e,0}{10.e-3,e,0}} {500}, File "Cut_e" ];
     }
  }

}

/* -------------------------------------------------------------------
   File "EleSta_v.pro"

   Electrostatics - Electric scalar potential v formulation
   ------------------------------------------------------------------- 

   I N P U T
   ---------

   Global Groups :  (Extension '_Ele' is for Electric problem)
   -------------
   Domain_Ele               Whole electric domain (not used)
   DomainCC_Ele             Nonconducting regions
   DomainC_Ele              Conducting regions (not used)

   Function :
   --------
   epsr[]                   Relative permittivity

   Constraint :
   ----------
   ElectricScalarPotential  Fixed electric scalar potential
                            (classical boundary condition)

   Physical constants :
   ------------------                                               */

   eps0 = 8.854187818e-12;

Group {
  DefineGroup[ Domain_Ele, DomainCC_Ele, DomainC_Ele ];
}

Function {
  DefineFunction[ epsr ];
}

FunctionSpace {
  { Name Hgrad_v_Ele; Type Form0;
    BasisFunction {
      // v = v  s   ,  for all nodes
      //      n  n
      { Name sn; NameOfCoef vn; Function BF_Node;
        Support DomainCC_Ele; Entity NodesOf[ All ]; }
    }
    Constraint {
      { NameOfCoef vn; EntityType NodesOf; 
        NameOfConstraint ElectricScalarPotential; }
    }
  }
}


Formulation {
  { Name Electrostatics_v; Type FemEquation;
    Quantity {
      { Name v; Type Local; NameOfSpace Hgrad_v_Ele; }
    }
    Equation {
      Galerkin { [ epsr[] * Dof{d v} , {d v} ]; In DomainCC_Ele; 
                 Jacobian Vol; Integration GradGrad; }
    }
  }
}


Resolution {
  { Name EleSta_v;
    System {
      { Name Sys_Ele; NameOfFormulation Electrostatics_v; }
    }
    Operation { 
      Generate[Sys_Ele]; Solve[Sys_Ele]; SaveSolution[Sys_Ele];
    }
  }
}


PostProcessing {
  { Name EleSta_v; NameOfFormulation Electrostatics_v;
    Quantity {
      { Name v; 
        Value { 
          Local { [ {v} ]; In DomainCC_Ele; Jacobian Vol; } 
        }
      }
      { Name e; 
        Value { 
          Local { [ -{d v} ]; In DomainCC_Ele; Jacobian Vol; }
        }
      }
      { Name d; 
        Value { 
          Local { [ -eps0*epsr[] * {d v} ]; In DomainCC_Ele; 
                                             Jacobian Vol; } 
        } 
      }
    }
  }
}

/* -------------------------------------------------------------------
   File "Jacobian_Lib.pro"

   Definition of a jacobian method
   -------------------------------------------------------------------

   I N P U T
   ---------

   GlobalGroup :
   -----------
   DomainInf                Regions with Spherical Shell Transformation

   Parameters :
   ----------
   Val_Rint, Val_Rext       Inner and outer radius of the Spherical Shell
                            of DomainInf
*/

Group {
  DefineGroup[ DomainInf ] ;
  DefineVariable[ Val_Rint, Val_Rext ] ;
}

Jacobian {
  { Name Vol ;
    Case { { Region DomainInf ;
             Jacobian VolSphShell {Val_Rint, Val_Rext} ; }
           { Region All ; Jacobian Vol ; }
    }
  }
}

/* -------------------------------------------------------------------
   File "Integration_Lib.pro"

   Definition of integration methods
   ------------------------------------------------------------------- */

Integration {
  { Name GradGrad ;
    Case { {Type Gauss ;
            Case { { GeoElement Triangle    ; NumberOfPoints  4 ; }
                   { GeoElement Quadrangle  ; NumberOfPoints  4 ; }
                   { GeoElement Tetrahedron ; NumberOfPoints  4 ; }
                   { GeoElement Hexahedron  ; NumberOfPoints  6 ; }
                   { GeoElement Prism       ; NumberOfPoints  9 ; } }
           }
         }
  }
  { Name CurlCurl ;
    Case { {Type Gauss ;
            Case { { GeoElement Triangle    ; NumberOfPoints  4 ; }
                   { GeoElement Quadrangle  ; NumberOfPoints  4 ; }
                   { GeoElement Tetrahedron ; NumberOfPoints  4 ; }
                   { GeoElement Hexahedron  ; NumberOfPoints  6 ; }
                   { GeoElement Prism       ; NumberOfPoints  9 ; } }
           }
         }
  }
}


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

6.2 Magnetostatic problem

We now consider a magnetostatic problem. The formulation used is a 2D magnetic vector potential formulation (see file `MagSta_a_2D.pro'). It is applied to a core-inductor system (file `CoreSta.pro'), whose geometry is defined in theh file `Core.geo' (see section C. Gmsh examples). The geometry is two-dimensional and, by symmetry, one fourth of the structure is modeled.

Core

The jacobian and integration methods used are the same as for the electrostatic problem presented in 6.1 Electrostatic problem.

/* -------------------------------------------------------------------
   File "CoreSta.pro"

   This file defines the problem dependent data structures for the
   static core-inductor problem.
   
   To compute the solution: 
       getdp CoreSta -msh Core.msh -solve MagSta_a_2D

   To compute post-results: 
       getdp CoreSta -msh Core.msh -pos Map_a
   ------------------------------------------------------------------- */

Group {

  Air    = Region[ 101 ];   Core   = Region[ 102 ];
  Ind    = Region[ 103 ];   AirInf = Region[ 111 ];

  SurfaceGh0  = Region[ 1100 ];  SurfaceGe0 = Region[ 1101 ];
  SurfaceGInf = Region[ 1102 ];

  Val_Rint = 200.e-3;
  Val_Rext = 250.e-3;

  DomainCC_Mag = Region[ {Air, AirInf, Core, Ind} ];
  DomainC_Mag  = Region[ {} ];
  DomainS_Mag  = Region[ {Ind} ]; // Stranded inductor
  DomainInf    = Region[ {AirInf} ];
  Domain_Mag   = Region[ {DomainCC_Mag, DomainC_Mag} ];

}

Function {

  mu0 = 4.e-7 * Pi;
  murCore = 100.;

  nu [ Region[{Air, Ind, AirInf}] ]  = 1. / mu0;
  nu [ Core ]  = 1. / (murCore * mu0);

  Sc[ Ind ] = 2.5e-2 * 5.e-2;

}

Constraint {

  { Name MagneticVectorPotential_2D;
    Case {
      { Region SurfaceGe0 ; Value 0.; }
      { Region SurfaceGInf; Value 0.; }
    }
  }

  Val_I_1_ = 0.01 * 1000.;

  { Name SourceCurrentDensityZ;
    Case {
      { Region Ind; Value Val_I_1_/Sc[]; }
    }
  }

}

Include "Jacobian_Lib.pro"
Include "Integration_Lib.pro"
Include "MagSta_a_2D.pro"

e = 1.e-5;
p1 = {e,e,0};
p2 = {0.12,e,0};

PostOperation {

  { Name Map_a; NameOfPostProcessing MagSta_a_2D;
    Operation {
      Print[ az, OnElementsOf Domain_Mag, File "CoreSta_a.pos" ];
      Print[ b, OnLine{{List[p1]}{List[p2]}} {1000}, File "k_a" ];
    }
  }

}

/* -------------------------------------------------------------------
   File "MagSta_a_2D.pro"

    Magnetostatics - Magnetic vector potential a formulation (2D)
   ------------------------------------------------------------------- 

   I N P U T
   ---------

   GlobalGroup :  (Extension '_Mag' is for Magnetic problem)
   -----------
   Domain_Mag               Whole magnetic domain
   DomainS_Mag              Inductor regions (Source)

   Function :
   --------
   nu[]                     Magnetic reluctivity

   Constraint :
   ----------
   MagneticVectorPotential_2D
                            Fixed magnetic vector potential (2D)
                            (classical boundary condition)
   SourceCurrentDensityZ    Fixed source current density (in Z direction)
*/

Group {
  DefineGroup[ Domain_Mag, DomainS_Mag ];
}

Function {
  DefineFunction[ nu ];
}

FunctionSpace {

  // Magnetic vector potential a (b = curl a)
  { Name Hcurl_a_Mag_2D; Type Form1P;
    BasisFunction {
      // a = a  s
      //      e  e
      { Name se; NameOfCoef ae; Function BF_PerpendicularEdge;
        Support Domain_Mag; Entity NodesOf[ All ]; }
    }
    Constraint {
      { NameOfCoef ae; EntityType NodesOf;
        NameOfConstraint MagneticVectorPotential_2D; }
    }
  }

  // Source current density js (fully fixed space)
  { Name Hregion_j_Mag_2D; Type Vector;
    BasisFunction {
      { Name sr; NameOfCoef jsr; Function BF_RegionZ;
        Support DomainS_Mag; Entity DomainS_Mag; }
    }
    Constraint {
      { NameOfCoef jsr; EntityType Region;
        NameOfConstraint SourceCurrentDensityZ; }
    }
  }

}

Formulation {
  { Name Magnetostatics_a_2D; Type FemEquation;
    Quantity {
      { Name a ; Type Local; NameOfSpace Hcurl_a_Mag_2D; }
      { Name js; Type Local; NameOfSpace Hregion_j_Mag_2D; }
    }
    Equation {
      Galerkin { [ nu[] * Dof{d a} , {d a} ]; In Domain_Mag;
                 Jacobian Vol; Integration CurlCurl; }
      Galerkin { [ - Dof{js} , {a} ]; In DomainS_Mag;
                 Jacobian Vol; Integration CurlCurl; }
    }
  }
}

Resolution {
  { Name MagSta_a_2D;
    System {
      { Name Sys_Mag; NameOfFormulation Magnetostatics_a_2D; }
    }
    Operation {
      Generate[Sys_Mag]; Solve[Sys_Mag]; SaveSolution[Sys_Mag];
    }
  }
}

PostProcessing {
  { Name MagSta_a_2D; NameOfFormulation Magnetostatics_a_2D;
    Quantity {
      { Name a; 
        Value { 
          Local { [ {a} ]; In Domain_Mag; Jacobian Vol; } 
        }
      }
      { Name az; 
        Value { 
          Local { [ CompZ[{a}] ]; In Domain_Mag; Jacobian Vol; }
        }
      }
      { Name b; 
        Value { 
          Local { [ {d a} ]; In Domain_Mag; Jacobian Vol; }
        }
      }
      { Name h; 
        Value { 
          Local { [ nu[] * {d a} ]; In Domain_Mag; Jacobian Vol; }
        }
      }
    }
  }
}


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

6.3 Magnetodynamic problem

As a third example we consider a magnetodynamic problem. The formulation is a two-dimensional a-v formulation (see file `MagDyn_av_2D.pro', which includes the same jacobian and integration library files as in 6.1 Electrostatic problem). It is applied to a core-inductor system (defined in file `CoreMassive.pro'), whose geometry has already been defined in file `Core.geo'.

/* -------------------------------------------------------------------
   File "CoreMassive.pro"

   This file defines the problem dependent data structures for the
   dynamic core-inductor problem.
   
   To compute the solution: 
        getdp CoreMassive -msh Core.msh -solve MagDyn_av_2D

   To compute post-results: 
        getdp CoreMassive -msh Core.msh -pos Map_a
        getdp CoreMassive -msh Core.msh -pos U_av
   ------------------------------------------------------------------- */

Group {

  Air    = Region[ 101 ];   Core   = Region[ 102 ];
  Ind    = Region[ 103 ];   AirInf = Region[ 111 ];

  SurfaceGh0  = Region[ 1100 ];  SurfaceGe0 = Region[ 1101 ];
  SurfaceGInf = Region[ 1102 ];

  Val_Rint = 200.e-3;
  Val_Rext = 250.e-3;

  DomainCC_Mag = Region[ {Air, AirInf} ];
  DomainC_Mag  = Region[ {Ind, Core} ]; // Massive inductor + conducting core
  DomainB_Mag  = Region[ {} ];
  DomainS_Mag  = Region[ {} ];
  DomainInf    = Region[ {AirInf} ];
  Domain_Mag   = Region[ {DomainCC_Mag, DomainC_Mag} ];

}

Function {

  mu0 = 4.e-7 * Pi;

  murCore = 100.;

  nu [ #{Air, Ind, AirInf} ]  = 1. / mu0;
  nu [ Core ]  = 1. / (murCore * mu0);
  sigma [ Ind ] = 5.9e7;
  sigma [ Core ] = 2.5e7;

  Freq = 1.;

}

Constraint {

  { Name MagneticVectorPotential_2D;
    Case {
      { Region SurfaceGe0 ; Value 0.; }
      { Region SurfaceGInf; Value 0.; }
    }
  }

  { Name SourceCurrentDensityZ;
    Case {
    }
  }

  Val_I_ = 0.01 * 1000.;

  { Name Current_2D;
    Case {
      { Region Ind; Value Val_I_; }
    }
  }

  { Name Voltage_2D;
    Case {
      { Region Core; Value 0.; }
    }
  }

}

Include "Jacobian_Lib.pro"
Include "Integration_Lib.pro"
Include "MagDyn_av_2D.pro"

PostOperation {
  { Name Map_a; NameOfPostProcessing MagDyn_av_2D;
    Operation {
      Print[ az, OnElementsOf Domain_Mag, File "Core_m_a.pos" ];
      Print[ j, OnElementsOf Domain_Mag, File "Core_m_j.pos" ];
    }
  }
  { Name U_av;  NameOfPostProcessing MagDyn_av_2D;
    Operation {
      Print[ U, OnRegion Ind ];
      Print[ I, OnRegion Ind ];
    }
  }
}

/* -------------------------------------------------------------------
   File "MagDyn_av_2D.pro"

    Magnetodynamics - Magnetic vector potential and electric scalar 
                      potential a-v formulation (2D)
   ------------------------------------------------------------------- 

   I N P U T
   ---------

   GlobalGroup :  (Extension '_Mag' is for Magnetic problem)
   -----------
   Domain_Mag               Whole magnetic domain
   DomainCC_Mag             Nonconducting regions (not used)
   DomainC_Mag              Conducting regions
   DomainS_Mag              Inductor regions (Source)
   DomainV_Mag              All regions in movement (for speed term)

   Function :
   --------
   nu[]                     Magnetic reluctivity
   sigma[]                  Electric conductivity

   Velocity[]               Velocity of regions

   Constraint :
   ----------
   MagneticVectorPotential_2D
                            Fixed magnetic vector potential (2D)
                            (classical boundary condition)
   SourceCurrentDensityZ    Fixed source current density (in Z direction)

   Voltage_2D               Fixed voltage
   Current_2D               Fixed Current

   Parameters :
   ----------

   Freq                     Frequency (Hz)

   Parameters for time loop with theta scheme :
   Mag_Time0, Mag_TimeMax, Mag_DTime
                            Initial time, Maximum time, Time step  (s)
   Mag_Theta                Theta  (e.g. 1.  : Implicit Euler,
                                         0.5 : Cranck Nicholson)
*/

Group {
  DefineGroup[ Domain_Mag, DomainCC_Mag, DomainC_Mag,
               DomainS_Mag, DomainV_Mag ];
}

Function {
  DefineFunction[ nu, sigma ];
  DefineFunction[ Velocity ];
  DefineVariable[ Freq ];
  DefineVariable[ Mag_Time0, Mag_TimeMax, Mag_DTime, Mag_Theta ];
}

FunctionSpace {

  // Magnetic vector potential a (b = curl a)
  { Name Hcurl_a_Mag_2D; Type Form1P;
    BasisFunction {
      // a = a  s
      //      e  e
      { Name se; NameOfCoef ae; Function BF_PerpendicularEdge;
        Support Domain_Mag; Entity NodesOf[ All ]; }
    }
    Constraint {
      { NameOfCoef ae; EntityType NodesOf;
        NameOfConstraint MagneticVectorPotential_2D; }
    }
  }

  // Gradient of Electric scalar potential (2D)
  { Name Hregion_u_Mag_2D; Type Form1P;
    BasisFunction {
      { Name sr; NameOfCoef ur; Function BF_RegionZ;
        Support DomainC_Mag; Entity DomainC_Mag; }
    }
    GlobalQuantity {
      { Name U; Type AliasOf       ; NameOfCoef ur; }
      { Name I; Type AssociatedWith; NameOfCoef ur; }
    }
    Constraint {
      { NameOfCoef U; EntityType Region;
        NameOfConstraint Voltage_2D; }
      { NameOfCoef I; EntityType Region;
        NameOfConstraint Current_2D; }
    }
  }

  // Source current density js (fully fixed space)
  { Name Hregion_j_Mag_2D; Type Vector;
    BasisFunction {
      { Name sr; NameOfCoef jsr; Function BF_RegionZ;
        Support DomainS_Mag; Entity DomainS_Mag; }
    }
    Constraint {
      { NameOfCoef jsr; EntityType Region;
        NameOfConstraint SourceCurrentDensityZ; }
    }
  }

}


Formulation {
  { Name Magnetodynamics_av_2D; Type FemEquation;
    Quantity {
      { Name a ; Type Local ; NameOfSpace Hcurl_a_Mag_2D; }
      { Name ur; Type Local ; NameOfSpace Hregion_u_Mag_2D; }
      { Name I ; Type Global; NameOfSpace Hregion_u_Mag_2D [I]; }
      { Name U ; Type Global; NameOfSpace Hregion_u_Mag_2D [U]; }
      { Name js; Type Local ; NameOfSpace Hregion_j_Mag_2D; }
    }
    Equation {
      Galerkin { [ nu[] * Dof{d a} , {d a} ]; In Domain_Mag;
                 Jacobian Vol; Integration CurlCurl; }

      Galerkin { DtDof [ sigma[] * Dof{a} , {a} ]; In DomainC_Mag;
                 Jacobian Vol; Integration CurlCurl; }
      Galerkin { [ sigma[] * Dof{ur} , {a} ]; In DomainC_Mag;
                 Jacobian Vol; Integration CurlCurl; }

      Galerkin { [ - sigma[] * (Velocity[] *^ Dof{d a}) , {a} ];
                 In DomainV_Mag;
                 Jacobian Vol; Integration CurlCurl; }

      Galerkin { [ - Dof{js} , {a} ]; In DomainS_Mag;
                 Jacobian Vol;
                 Integration CurlCurl; }

      Galerkin { DtDof [ sigma[] * Dof{a} , {ur} ]; In DomainC_Mag;
                 Jacobian Vol; Integration CurlCurl; }
      Galerkin { [ sigma[] * Dof{ur} , {ur} ]; In DomainC_Mag;
                 Jacobian Vol; Integration CurlCurl; }
      GlobalTerm { [ Dof{I} , {U} ]; In DomainC_Mag; }
    }
  }
}


Resolution {
  { Name MagDyn_av_2D;
    System {
      { Name Sys_Mag; NameOfFormulation Magnetodynamics_av_2D;
        Type ComplexValue; Frequency Freq; }
    }
    Operation {
      Generate[Sys_Mag]; Solve[Sys_Mag]; SaveSolution[Sys_Mag];
    }
  }

  { Name MagDyn_t_av_2D;
    System {
      { Name Sys_Mag; NameOfFormulation Magnetodynamics_av_2D; }
    }
    Operation {
      InitSolution[Sys_Mag]; SaveSolution[Sys_Mag];
      TimeLoopTheta[Mag_Time0, Mag_TimeMax, Mag_DTime, Mag_Theta] {
        Generate[Sys_Mag]; Solve[Sys_Mag]; SaveSolution[Sys_Mag]; 
      }
    }
  }

}


PostProcessing {
  { Name MagDyn_av_2D; NameOfFormulation Magnetodynamics_av_2D;
    Quantity {
      { Name a; 
        Value { 
          Local { [ {a} ]; In Domain_Mag; Jacobian Vol; } 
        } 
      }
      { Name az; 
        Value { 
          Local { [ CompZ[{a}] ]; In Domain_Mag; Jacobian Vol; }
        }
      }
      { Name b; 
        Value { 
          Local { [ {d a} ]; In Domain_Mag; Jacobian Vol; }
        }
      }
      { Name h; 
        Value { 
          Local { [ nu[] * {d a} ]; In Domain_Mag; Jacobian Vol; }
        }
      }
      { Name j; 
        Value { 
          Local { [ - sigma[]*(Dt[{a}]+{ur}) ]; In DomainC_Mag; 
                  Jacobian Vol; } 
        } 
      }
      { Name jz; 
        Value { 
          Local { [ - sigma[]*CompZ[Dt[{a}]+{ur}] ]; In DomainC_Mag; 
                  Jacobian Vol; } 
        } 
      }
      { Name roj2;
        Value { 
          Local { [ sigma[]*SquNorm[Dt[{a}]+{ur}] ]; In DomainC_Mag; 
                  Jacobian Vol; } 
        } 
      }
      { Name U; Value { Local { [ {U} ]; In DomainC_Mag; } } }
      { Name I; Value { Local { [ {I} ]; In DomainC_Mag; } } }
    }
  }
}


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

7. Running GetDP

GetDP has no graphical interface(1). It is a command-line driven program that reads a problem definition file once at the beginning of the processing. This problem definition file is a regular ASCII text file (see section 1.1 Numerical tools as objects), hence created with whatever text editor you like.

If you just type the program name at your shell prompt (without any argument), you will get a short help on how to run GetDP. All GetDP calls look like

 
getdp filename options 

where filename is the ASCII file containing the problem definition, i.e., the structures this user's guide has taught you to create. This file can include other files (see section 1.4 Includes), so that only one problem definition file should always be given on the command line. The input files containing the problem definition structure are usually given the `.pro' extension (if so, there is no need to specify the extension on the command line). The name of this file (without the extension) is used as a basis for the creation of intermediate files during the pre-processing and the processing stages.

The options are a combination of the following commands (in any order):

-pre
resolution-id

Performs the pre-processing associated with the resolution resolution-id. In the pre-processing stage, GetDP creates the geometric database (from the mesh file), identifies the degrees of freedom (the unknowns) of the problem and sets up the constraints on these degrees of freedom. The pre-processing creates a file with a `.pre' extension. If resolution-id is omitted, the list of available choices is displayed.

-cal

Performs the processing. This requires that a pre-processing has been performed previously, or that a -pre option is given on the same command line. The performed resolution is the one given as an argument to the -pre option. In the processing stage, GetDP executes all the commands given in the Operation field of the selected Resolution object (such as matrix assemblies, system resolutions, ...).

-pos
post-operation-id ...

Performs the operations in the PostOperation(s) selected by the post-operation-id(s). This requires that a processing has been performed previously, or that a -cal option is given on the same command line. If post-operation-id is omitted, the list of available choices is displayed.

-ipos
post-processing-id ...

Enters an interactive post-processing mode, permetting to manually type PostOperation-like commands. These commands are based on the PostProcessing object(s) selected by the post-processing-id(s). If post-processing-id is omitted, the list of available choices is displayed.

-msh
filename

Reads the mesh (in .msh format) from filename (see section 8. File formats) rather than from the default problem file name (with the `.msh' extension appended).

-split

Saves processing results in separate files (one for each timestep).

-res
filename ...

Loads processing results from file(s).

-name
string

Uses string as the default generic file name for input or output of mesh, pre-processing and processing files.

-restart

Restarts processing of a time stepping resolution interrupted before being complete.

-solve
resolution-id

Same as -pre resolution-id -cal.

-adapt
file

Reads adaptation constraints from file.

-order
real

Specifies the maximum interpolation order.

-bin

Selects binary format for output files.

-log

Saves all processing history in a log file (the input file name with an appended .log extension).

-socket
string

Communicates through socket string.

-check

Lets you check the problem structure interactively.

-v
-verbose
integer

Sets the verbosity level. A value of 0 means that no information will be displayed during the processing.

-p
-progress
integer

Sets the progress update rate. This controls the refreshment rate of the counter indicating the progress of the current computation (in %).

-info

Displays the version information.

-version

Displays the version number.

-help

Displays a message listing basic usage and available options.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

8. File formats

This chapter describes the file formats that cannot be modified by the user. The format of the problem definition structure is explained in 3. Objects, and 4. Types for objects. The format of the post-processing files is explained in 4.10 Types for PostOperation.

8.1 Input file format  
8.2 Output file format  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

8.1 Input file format

The native mesh format read by GetDP is the mesh file format produced by Gmsh (http://www.geuz.org/gmsh/). The file is divided into two sections, defining the nodes and the elements in the mesh.

 
$NOD
number-of-nodes
node-number x-coord y-coord z-coord
...
$ENDNOD
$ELM
number-of-elements
elm-number elm-type elm-region unused number-of-nodes node-numbers
...
$ENDELM

All the syntactic variables stand for integers except x-coord, y-coord and z-coord which stand for floating point values. The elm-type value defines the geometrical type for the element:

elm-type:

1
Line (2 nodes, 1 edge).
2
Triangle (3 nodes, 3 edges).
3
Quadrangle (4 nodes, 4 edges).
4
Tetrahedron (4 nodes, 6 edges, 4 facets).
5
Hexahedron (8 nodes, 12 edges, 6 facets).
6
Prism (6 nodes, 9 edges, 5 facets).
7
Pyramid (5 nodes, 8 edges, 5 facets).
15
Point (1 node).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

8.2 Output file format

8.2.1 File `.pre'  
8.2.2 File `.res'  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

8.2.1 File `.pre'

The `.pre' file is generated by the pre-processing stage. It contains all the information about the degrees of freedom to be considered during the processing stage for a given resolution (i.e., unknowns, fixed values, initial values, etc.).

 
$Resolution /* 'resolution-id' */
main-resolution-number number-of-dofdata
$EndResolution
$DofData /* #dofdata-number */
resolution-number system-number
number-of-function-spaces function-space-number ...
number-of-time-functions time-function-number ...
number-of-partitions partition-index ...
number-of-any-dof number-of-dof
dof-basis-function-number dof-entity dof-harmonic dof-type dof-data
...
$EndDofData
...

with

 
dof-data:
  equation-number nnz
    (dof-type: 1; unknown) |
  dof-value dof-time-function-number
    (dof-type: 2; fixed value) |
  dof-associate-dof-number dof-value dof-time-function-number
    (dof-type: 3; associated degree of freedom) |
  equation-number dof-value
    (dof-type: 5; initial value for an unknown)

Notes:

  1. There is one $DofData field for each system of equations considered in the resolution (including those considered in pre-resolutions).
  2. The dofdata-number of a $DofData field is determined by the order of this field in the `.pre' file.
  3. number-of-dof is the dimension of the considered system of equations, while number-of-any-dof is the total number of degrees of freedom before the application of constraints.
  4. Each degree of freedom is coded with three integer values, which are the associated basis function, entity and harmonic numbers, i.e., dof-basis-function-number, dof-entity and dof-harmonic.
  5. nnz is not used at the moment.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

8.2.2 File `.res'

The `.res' file is generated by the processing stage. It contains the solution of the problem (or a part of it in case of program interruption).

 
$ResFormat /* GetDP vgetdp-version-number, string-for-format */
1.1 file-res-format
$EndResFormat
$Solution  /* DofData #dofdata-number */
dofdata-number time-value time-imag-value time-step-number
solution-value
...
$EndSolution
...

Notes:

  1. A $Solution field contains the solution associated with a $DofData field.
  2. There is one $Solution field for each time step, of which the time is time-value (0 for non time dependent or non modal analyses) and the imaginary time is time-imag-value (0 for non time dependent or non modal analyses).
  3. The order of the solution-values in a $Solution field follows the numbering of the equations given in the `.pre' file (one floating point value for each degree of freedom).


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

9. Bugs, versions and credits

9.1 Bugs  
9.2 Versions  
9.3 Credits  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

9.1 Bugs

If you think you have found a bug in GetDP, you can report it by electronic mail to the GetDP mailing list at getdp@geuz.org. Please send as precise a description of the problem as you can, including sample input files that produce the bug (problem definition and mesh files). Don't forget to mention both the version of GetDP and the version of your operation system (see section 7. Running GetDP to see how to get this information).

See the `TODO' file in the distribution to check the problems we already know about.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

9.2 Versions

$Id: VERSIONS,v 1.58 2006/03/18 05:05:59 geuzaine Exp $

New since 1.2: small bug fixes (binary .res read, Newmark -restart).

New in 1.2: Windows versions do not depend on Cygwin anymore; major
parser cleanup (loops & co).

New in 1.1: New eigensolver based on Arpack (EigenSolve); generalized
old Lanczos solver to work with GSL+lapack; reworked PETSc interface,
which now requires PETSc 2.3; documented many previously undocumented
features (loops, conditionals, strings, link constraints, etc.);
various improvements and bug fixes.

New in 1.0: New license (GNU GPL); added support for latest Gmsh mesh
file format; more code cleanups.

New in 0.91: Merged moving band and multi-harmonic code; new loops and
conditionals in the parser; removed old readline code (just use GNU
readline if available); upgraded to latest Gmsh post-processing
format; various small enhancements and bug fixes.

New in 0.89: Code cleanup.

New in 0.88: Integrated FMM code.

New in 0.87: Fixed major performance problem on Windows (matrix
assembly and post-processing can be up to 3-4 times faster with 0.87
compared to 0.86, bringing performance much closer to Unix versions);
fixed stack overflow on Mac OS X; Re-introduced face basis functions
mistakenly removed in 0.86; fixed post-processing bug with pyramidal
basis functions; new build system based on autoconf.

New in 0.86: Updated Gmsh output format; many small bug fixes.

New in 0.85: Upgraded communication interface with Gmsh; new
ChangeOfValues option in PostOperation; many internal changes.

New in 0.84: New ChangeOfCoordinate option in PostOperation; fixed
crash in InterpolationAkima; improved interactive postprocessing
(-ipos); changed syntax of parametric OnGrid ($S, $T -> $A, $B, $C);
corrected Skin for non simplicial meshes; fixed floating point
exception in diagonal matrix scaling; many other small fixes and
cleanups.

New in 0.83: Fixed bugs in SaveSolutions[] and InitSolution[]; fixed
corrupted binary post-processing files in the harmonic case for the
Gmsh format; output files are now created relatively to the input file
directory; made solver options available on the command line; added
optional matrix scaling and changed default parameter file name to
'solver.par' (Warning: please check the scaling definition in your old
SOLVER.PAR files); generalized syntax for lists (start:[incr]end ->
start:end:incr); updated reference guide; added a new short
presentation on the web site; OnCut -> OnSection; new functional
syntax for resolution operations (e.g. Generate X -> Generate[X]);
many other small fixes and cleanups.

New in 0.82: Added communication socket for interactive use with Gmsh;
corrected (again) memory problem (leak + seg. fault) in time stepping
schemes; corrected bug in Update[].

New in 0.81: Generalization of transformation jacobians (spherical and
rectangular, with optional parameters); changed handling of missing
command line arguments; enhanced Print OnCut; fixed memory leak for
time domain analysis of coupled problems; -name option; fixed
seg. fault in ILUK.

New in 0.80: Fixed computation of time derivatives on first time step
(in post-processing); added tolerance in transformation jacobians;
fixed parsing of DOS files (carriage return problems); automatic
memory reallocation in ILUD/ILUK.

New in 0.79: Various bug fixes (mainly for the post-processing of
intergal quantities); automatic treatment of degenerated cases in
axisymmetrical problems.

New in 0.78: Various bug fixes.

New in 0.77: Changed syntax for PostOperations (Plot suppressed in
favour of Print; Plot OnRegion becomes Print OnElementsOf); changed
table oriented post-processing formats; new binary formats; new error
diagnostics.

New in 0.76: Reorganized high order shape functions; optimization of
the post-processing (faster and less bloated); lots of internal
cleanups.

New in 0.74: High order shape functions; lots of small bug fixes.

New in 0.73: Eigen value problems (Lanczos); minor corrections.

New in 0.7: constraint syntax; fourier transform; unary minus
correction; complex integral quantity correction; separate iteration
matrix generation.

New in 0.6: Second order time derivatives; Newton nonlinear scheme;
Newmark time stepping scheme; global quantity syntax; interactive
post-processing; tensors; integral quantities; post-processing
facilities.

New in 0.3: First distributed version.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

9.3 Credits

$Id: CREDITS,v 1.17 2006/03/13 20:49:51 geuzaine Exp $

             GetDP is copyright (C) 1997-2006

                      Patrick Dular
              <patrick.dular at ulg.ac.be>

                           and 

                   Christophe Geuzaine 
             <christophe.geuzaine at case.edu> 

Major code contributions to GetDP have been provided by Johan
Gyselinck <johan.gyselinck at ulb.ac.be> and Ruth Sabariego
<r.sabariego at ulg.ac.be>.

Other code contributors include: David Colignon <david.colignon at
ulg.ac.be>, Tuan Ledinh <tuan.ledinh at student.ulg.ac.be>, Patrick
Lefevre <patrick.lefevre at umh.ac.be>, Andre Nicolet <andre.nicolet
at fresnel.fr>, <ohyeahq at yahoo.co.jp>, Jean-Francois Remacle
<remacle at gce.ucl.ac.be>, Timo Tarhasaari <timo.tarhasaari at
tut.fi>, Christophe Trophime <trophime at labs.polycnrs-gre.fr> and
Marc Ume <Marc.Ume at digitalgraphics.be>. See the source code for
more details.

The AVL tree code (DataStr/avl.*) is copyright (C) 1988-1993, 1995 The
Regents of the University of California. Permission to use, copy,
modify, and distribute this software and its documentation for any
purpose and without fee is hereby granted, provided that the above
copyright notice appear in all copies and that both that copyright
notice and this permission notice appear in supporting documentation,
and that the name of the University of California not be used in
advertising or publicity pertaining to distribution of the software
without specific, written prior permission.  The University of
California makes no representations about the suitability of this
software for any purpose.  It is provided "as is" without express or
implied warranty.

The Arpack code (in the Arpack subdirectory) was written by Danny
Sorensen, Richard Lehoucq, Chao Yang and Kristi Maschhoff from the
Dept. of Computational & Applied Mathematics at Rice University,
Houston, Texas, USA. See http://www.caam.rice.edu/software/ARPACK/ for
more info.

This version of GetDP may contain code (in the Sparskit subdirectory)
copyright (C) 1990 Yousef Saad: check the configuration options.

Thanks to the following folks who have contributed by providing fresh
ideas on theoretical or programming topics, who have sent patches,
requests for changes or improvements, or who gave us access to exotic
machines for testing GetDP: Olivier Adam <o.adam at ulg.ac.be>,
Alejandro Angulo <aacusm at yahoo.com>, Geoffrey Deliege <deliege at
mailserv.esat.kuleuven.ac.be>, Mark Evans <evans at gte.net>, Philippe
Geuzaine <geuzaine at gnat.colorado.edu>, Eric Godard <godard at
montefiore.ulg.ac.be>, Sebastien Guenneau <guenneau at
liverpool.ac.uk>, Francois Henrotte <fhenrott at esat.kuleuven.ac.be>,
Daniel Kedzierski <kedzierski at uol.com.br>, Samuel Kvasnica
<kvasnica at iaee.tuwien.ac.at>, Benoit Meys <bmeys at
techspace-aero.be>, Uwe Pahner <uwe.pahner at esat.kuleuven.ac.be>,
Georgia Psoni and Robert Struijs <gpsoni at free.fr>, Ahmed Rassili
<a.rassili at ulg.ac.be>, Thierry Scordilis <t.scordilis at atral.fr>,
Herve Tortel <tortel at loe.u-3mrs.fr>, Jose Geraldo A. Brito Neto
<jgabrito at terra.com.br>, Matthias Fenner <m.fenner at gmx.net>,
Daryl Van Vorst <dvpub@telus.net>.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

A. Tips and tricks


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

B. Frequently asked questions

$Id: FAQ,v 1.11 2006/02/25 08:02:58 geuzaine Exp $

This is the GetDP FAQ

********************************************************************

Section 1: The basics

* 1.1 What is GetDP?

GetDP is a scientific software environment for the numerical solution
of integro-differential equations, open to the coupling of physical
problems (electromagnetic, thermal, mechanical, etc) as well as of
numerical methods (finite element method, integral methods, etc). It
can deal with such problems of various dimensions (1D, 2D, 2D
axisymmetric or 3D) and time states (static, transient or
harmonic). The main feature of GetDP is the closeness between the
organization of data defining discrete problems (written by the user
in ASCII data files) and the symbolic mathematical expressions of
these problems.

* 1.2 What are the terms and conditions of use?

GetDP is distributed under the terms of the GNU General Public
License. See the file doc/LICENSE for more information, or go to the
GNU foundation's web site at http://www.gnu.org.

* 1.3 What does 'GetDP' mean?

It's an acronym for a "General environment for the treatment of
Discrete Problems".

* 1.4 Where can I find more information?

<http://www.geuz.org/getdp/> is the primary site to obtain information
about GetDP. You will find a short presentation, a complete reference
guide as well as a searchable archive of the GetDP mailing list
(<getdp@geuz.org>) on this site.

********************************************************************

Section 2: Installation

* 2.1 Which OSes does GetDP run on?

GetDP is known to run on Windows 95/98/NT/2000/XP, Linux, Mac OS X,
Compaq Tru64 Unix (aka OSF1, aka Digital Unix), Sun OS, IBM AIX, SGI
IRIX, FreeBSD and HP-UX. It should compile on any Unix-like operating
system, provided that you have access to a recent C and Fortran 77
compiler.

* 2.2 What do I need to compile GetDP from the sources?

You need a C and a Fortran 77 compiler (e.g. the GNU compilers gcc and
g77) as well as the GSL (version 1.2 or higher; freely available from
http://sources.redhat.com/gsl/).

* 2.3 How to I compile GetDP?

Just type 

./configure; make; make install

If you change some configuration options (type ./configure --help to
get the list of all available choices), don't forget to do 'make
clean' before rebuilding GetDP.

* 2.4 GetDP [from a binary distribution] complains about missing
libraries.

Try 'ldd getdp' (or 'otool -L getdp' on Mac OS X) to check if all the
required shared libraries are installed on your system. If not,
install them. If it still doesn't work, recompile GetDP from the
sources.

********************************************************************

Section 3: Usage

* 3.1 How can I provide a mesh to GetDP?

The only meshing format accepted by this version of GetDP is the 'msh'
format created by Gmsh (<http://www.geuz.org/gmsh>). This format being
very simple (see the reference GetDP manual for more details), it
should be straightforward to write a converter from your mesh format
to the 'msh' format.

* 3.2 How can I visualize the results produced by GetDP?  

You can specify a format in all post-processing operations. Available
formats include 'Table', 'TimeTable' and 'Gmsh'. 'Table' and
'TimeTable' output lists of numbers easily readable by
Excel/gnuplot/Caleida Graph/etc. 'Gmsh' outputs post-processing views
directly loadable by Gmsh.

* 3.3 How do I change the linear solver used by GetDP?

Edit the 'solver.par' file in the current working directory. You can
also remove the file: GetDP will give you the opportunity to create it
dynamically next time you perform a linear system solving.

If you don't like the default linear solvers (based on Yousef Saad's
Sparskit 2.0), you can configure and compile GetDP to use PETSc
(http://www.mcs.anl.gov/petsc) instead: run ./configure --help for
more info.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

C. Gmsh examples

Gmsh is a three-dimensional finite element mesh generator with simple CAD and post-processing capabilities that can be used as a graphical front-end for GetDP. Gmsh can be downloaded from http://www.geuz.org/gmsh/.

This appendix reproduces verbatim the input files needed by Gmsh to produce the mesh files `mStrip.msh' and `Core.msh' used in the examples of 6. Complete examples.

/* -------------------------------------------------------------------
   File "mStrip.geo"

   This file is the geometrical description used by GMSH to produce
   the file "mStrip.msh".
   ------------------------------------------------------------------- */

/* Definition of some parameters for geometrical dimensions, i.e.
   h (height of 'Diel1'), w (width of 'Line'), t (thickness of 'Line')
   xBox (width of the air box) and yBox (height of the air box) */

h = 1.e-3 ; w = 4.72e-3 ;  t = 0.035e-3 ;
xBox = w/2. * 6. ;  yBox = h * 12. ;

/* Definition of parameters for local mesh dimensions */

s = 1. ;
p0 = h / 10. * s ;
pLine0 = w/2. / 10. * s ;  pLine1 = w/2. / 50. * s ;
pxBox = xBox / 10. * s ;  pyBox = yBox / 8. * s ;

/* Definition of gemetrical points */

Point(1) = { 0   , 0, 0, p0} ;
Point(2) = { xBox, 0, 0, pxBox} ;
Point(3) = { xBox, h, 0, pxBox} ;
Point(4) = { 0   , h, 0, pLine0} ;
Point(5) = { w/2., h, 0, pLine1} ;
Point(6) = { 0   , h+t, 0, pLine0} ;
Point(7) = { w/2., h+t, 0, pLine1} ;
Point(8) = { 0   , yBox, 0, pyBox} ;
Point(9) = { xBox, yBox, 0, pyBox} ;

/* Definition of gemetrical lines */

Line(1) = {1,2};   Line(2) = {2,3};  Line(3) = {3,9};
Line(4) = {9,8};   Line(5) = {8,6};  Line(7) = {4,1};
Line(8) = {5,3};   Line(9) = {4,5};  Line(10) = {6,7};
Line(11) = {5,7};

/* Definition of geometrical surfaces */

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 */

Physical Surface (101) = {15} ;   /* Air */
Physical Surface (111) = {13} ;   /* Diel1 */

Physical Line (120) = {1} ;       /* Ground */
Physical Line (121) = {9,10,11} ; /* Line */
Physical Line (130) = {2,3,4} ;   /* SurfInf */

/* -------------------------------------------------------------------
   File "Core.geo"

   This file is the geometrical description used by GMSH to produce
   the file "Core.msh".
   ------------------------------------------------------------------- */

dxCore =  50.e-3; dyCore = 100.e-3;
xInd   =  75.e-3; dxInd  =  25.e-3; dyInd  =  50.e-3;
rInt   = 200.e-3; rExt   = 250.e-3;

s       =  1.;
p0      = 12.e-3 *s;
pCorex  =  4.e-3 *s; pCorey0 =  8.e-3 *s; pCorey  =  4.e-3 *s;
pIndx   =  5.e-3 *s; pIndy   =  5.e-3 *s;
pInt    = 12.5e-3*s; pExt    = 12.5e-3*s;

Point(1) = {0,0,0,p0};
Point(2) = {dxCore,0,0,pCorex};
Point(3) = {dxCore,dyCore,0,pCorey};
Point(4) = {0,dyCore,0,pCorey0};
Point(5) = {xInd,0,0,pIndx};
Point(6) = {xInd+dxInd,0,0,pIndx};
Point(7) = {xInd+dxInd,dyInd,0,pIndy};
Point(8) = {xInd,dyInd,0,pIndy};
Point(9) = {rInt,0,0,pInt};
Point(10) = {rExt,0,0,pExt};
Point(11) = {0,rInt,0,pInt};
Point(12) = {0,rExt,0,pExt};

Line(1) = {1,2};  Line(2) = {2,5};   Line(3) = {5,6};
Line(4) = {6,9};  Line(5) = {9,10};  Line(6) = {1,4};
Line(7) = {4,11}; Line(8) = {11,12}; Line(9) = {2,3};
Line(10) = {3,4}; Line(11) = {6,7};  Line(12) = {7,8};
Line(13) = {8,5};

Circle(14) = {9,1,11};  Circle(15) = {10,1,12};

Line Loop(16) = {-6,1,9,10};                 Plane Surface(17) = {16};
Line Loop(18) = {11,12,13,3};                Plane Surface(19) = {18};
Line Loop(20) = {7,-14,-4,11,12,13,-2,9,10}; Plane Surface(21) = {20};
Line Loop(22) = {8,-15,-5,14};               Plane Surface(23) = {22};

Physical Surface(101) = {21};  /* Air */
Physical Surface(102) = {17};  /* Core */
Physical Surface(103) = {19};  /* Ind */
Physical Surface(111) = {23};  /* AirInf */

Physical Line(1000) = {1,2};       /* Cut */
Physical Line(1001) = {2};         /* CutAir */
Physical Line(202)  = {9,10};      /* SkinCore */
Physical Line(203)  = {11,12,13};  /* SkinInd */
Physical Line(1100) = {1,2,3,4,5}; /* SurfaceGh0 */
Physical Line(1101) = {6,7,8};     /* SurfaceGe0 */
Physical Line(1102) = {15};        /* SurfaceGInf */


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

D. License

GNU General Public License

Version 2, June 1991

 
Copyright © 1989, 1991 Free Software Foundation, Inc.
59 Temple Place - Suite 330, Boston, MA  02111-1307, USA

Everyone is permitted to copy and distribute verbatim copies
of this license document, but changing it is not allowed.

Preamble

The licenses for most software are designed to take away your freedom to share and change it. By contrast, the GNU General Public License is intended to guarantee your freedom to share and change free software--to make sure the software is free for all its users. This General Public License applies to most of the Free Software Foundation's software and to any other program whose authors commit to using it. (Some other Free Software Foundation software is covered by the GNU Library General Public License instead.) You can apply it to your programs, too.

When we speak of free software, we are referring to freedom, not price. Our General Public Licenses are designed to make sure that you have the freedom to distribute copies of free software (and charge for this service if you wish), that you receive source code or can get it if you want it, that you can change the software or use pieces of it in new free programs; and that you know you can do these things.

To protect your rights, we need to make restrictions that forbid anyone to deny you these rights or to ask you to surrender the rights. These restrictions translate to certain responsibilities for you if you distribute copies of the software, or if you modify it.

For example, if you distribute copies of such a program, whether gratis or for a fee, you must give the recipients all the rights that you have. You must make sure that they, too, receive or can get the source code. And you must show them these terms so they know their rights.

We protect your rights with two steps: (1) copyright the software, and (2) offer you this license which gives you legal permission to copy, distribute and/or modify the software.

Also, for each author's protection and ours, we want to make certain that everyone understands that there is no warranty for this free software. If the software is modified by someone else and passed on, we want its recipients to know that what they have is not the original, so that any problems introduced by others will not reflect on the original authors' reputations.

Finally, any free program is threatened constantly by software patents. We wish to avoid the danger that redistributors of a free program will individually obtain patent licenses, in effect making the program proprietary. To prevent this, we have made it clear that any patent must be licensed for everyone's free use or not licensed at all.

The precise terms and conditions for copying, distribution and modification follow.

TERMS AND CONDITIONS FOR COPYING, DISTRIBUTION AND MODIFICATION

  1. This License applies to any program or other work which contains a notice placed by the copyright holder saying it may be distributed under the terms of this General Public License. The "Program", below, refers to any such program or work, and a "work based on the Program" means either the Program or any derivative work under copyright law: that is to say, a work containing the Program or a portion of it, either verbatim or with modifications and/or translated into another language. (Hereinafter, translation is included without limitation in the term "modification".) Each licensee is addressed as "you".

    Activities other than copying, distribution and modification are not covered by this License; they are outside its scope. The act of running the Program is not restricted, and the output from the Program is covered only if its contents constitute a work based on the Program (independent of having been made by running the Program). Whether that is true depends on what the Program does.

  2. You may copy and distribute verbatim copies of the Program's source code as you receive it, in any medium, provided that you conspicuously and appropriately publish on each copy an appropriate copyright notice and disclaimer of warranty; keep intact all the notices that refer to this License and to the absence of any warranty; and give any other recipients of the Program a copy of this License along with the Program.

    You may charge a fee for the physical act of transferring a copy, and you may at your option offer warranty protection in exchange for a fee.

  3. You may modify your copy or copies of the Program or any portion of it, thus forming a work based on the Program, and copy and distribute such modifications or work under the terms of Section 1 above, provided that you also meet all of these conditions:

    1. You must cause the modified files to carry prominent notices stating that you changed the files and the date of any change.

    2. You must cause any work that you distribute or publish, that in whole or in part contains or is derived from the Program or any part thereof, to be licensed as a whole at no charge to all third parties under the terms of this License.

    3. If the modified program normally reads commands interactively when run, you must cause it, when started running for such interactive use in the most ordinary way, to print or display an announcement including an appropriate copyright notice and a notice that there is no warranty (or else, saying that you provide a warranty) and that users may redistribute the program under these conditions, and telling the user how to view a copy of this License. (Exception: if the Program itself is interactive but does not normally print such an announcement, your work based on the Program is not required to print an announcement.)

    These requirements apply to the modified work as a whole. If identifiable sections of that work are not derived from the Program, and can be reasonably considered independent and separate works in themselves, then this License, and its terms, do not apply to those sections when you distribute them as separate works. But when you distribute the same sections as part of a whole which is a work based on the Program, the distribution of the whole must be on the terms of this License, whose permissions for other licensees extend to the entire whole, and thus to each and every part regardless of who wrote it.

    Thus, it is not the intent of this section to claim rights or contest your rights to work written entirely by you; rather, the intent is to exercise the right to control the distribution of derivative or collective works based on the Program.

    In addition, mere aggregation of another work not based on the Program with the Program (or with a work based on the Program) on a volume of a storage or distribution medium does not bring the other work under the scope of this License.

  4. You may copy and distribute the Program (or a work based on it, under Section 2) in object code or executable form under the terms of Sections 1 and 2 above provided that you also do one of the following:

    1. Accompany it with the complete corresponding machine-readable source code, which must be distributed under the terms of Sections 1 and 2 above on a medium customarily used for software interchange; or,

    2. Accompany it with a written offer, valid for at least three years, to give any third party, for a charge no more than your cost of physically performing source distribution, a complete machine-readable copy of the corresponding source code, to be distributed under the terms of Sections 1 and 2 above on a medium customarily used for software interchange; or,

    3. Accompany it with the information you received as to the offer to distribute corresponding source code. (This alternative is allowed only for noncommercial distribution and only if you received the program in object code or executable form with such an offer, in accord with Subsection b above.)

    The source code for a work means the preferred form of the work for making modifications to it. For an executable work, complete source code means all the source code for all modules it contains, plus any associated interface definition files, plus the scripts used to control compilation and installation of the executable. However, as a special exception, the source code distributed need not include anything that is normally distributed (in either source or binary form) with the major components (compiler, kernel, and so on) of the operating system on which the executable runs, unless that component itself accompanies the executable.

    If distribution of executable or object code is made by offering access to copy from a designated place, then offering equivalent access to copy the source code from the same place counts as distribution of the source code, even though third parties are not compelled to copy the source along with the object code.

  5. You may not copy, modify, sublicense, or distribute the Program except as expressly provided under this License. Any attempt otherwise to copy, modify, sublicense or distribute the Program is void, and will automatically terminate your rights under this License. However, parties who have received copies, or rights, from you under this License will not have their licenses terminated so long as such parties remain in full compliance.

  6. You are not required to accept this License, since you have not signed it. However, nothing else grants you permission to modify or distribute the Program or its derivative works. These actions are prohibited by law if you do not accept this License. Therefore, by modifying or distributing the Program (or any work based on the Program), you indicate your acceptance of this License to do so, and all its terms and conditions for copying, distributing or modifying the Program or works based on it.

  7. Each time you redistribute the Program (or any work based on the Program), the recipient automatically receives a license from the original licensor to copy, distribute or modify the Program subject to these terms and conditions. You may not impose any further restrictions on the recipients' exercise of the rights granted herein. You are not responsible for enforcing compliance by third parties to this License.

  8. If, as a consequence of a court judgment or allegation of patent infringement or for any other reason (not limited to patent issues), conditions are imposed on you (whether by court order, agreement or otherwise) that contradict the conditions of this License, they do not excuse you from the conditions of this License. If you cannot distribute so as to satisfy simultaneously your obligations under this License and any other pertinent obligations, then as a consequence you may not distribute the Program at all. For example, if a patent license would not permit royalty-free redistribution of the Program by all those who receive copies directly or indirectly through you, then the only way you could satisfy both it and this License would be to refrain entirely from distribution of the Program.

    If any portion of this section is held invalid or unenforceable under any particular circumstance, the balance of the section is intended to apply and the section as a whole is intended to apply in other circumstances.

    It is not the purpose of this section to induce you to infringe any patents or other property right claims or to contest validity of any such claims; this section has the sole purpose of protecting the integrity of the free software distribution system, which is implemented by public license practices. Many people have made generous contributions to the wide range of software distributed through that system in reliance on consistent application of that system; it is up to the author/donor to decide if he or she is willing to distribute software through any other system and a licensee cannot impose that choice.

    This section is intended to make thoroughly clear what is believed to be a consequence of the rest of this License.

  9. If the distribution and/or use of the Program is restricted in certain countries either by patents or by copyrighted interfaces, the original copyright holder who places the Program under this License may add an explicit geographical distribution limitation excluding those countries, so that distribution is permitted only in or among countries not thus excluded. In such case, this License incorporates the limitation as if written in the body of this License.

  10. The Free Software Foundation may publish revised and/or new versions of the General Public License from time to time. Such new versions will be similar in spirit to the present version, but may differ in detail to address new problems or concerns.

    Each version is given a distinguishing version number. If the Program specifies a version number of this License which applies to it and "any later version", you have the option of following the terms and conditions either of that version or of any later version published by the Free Software Foundation. If the Program does not specify a version number of this License, you may choose any version ever published by the Free Software Foundation.

  11. If you wish to incorporate parts of the Program into other free programs whose distribution conditions are different, write to the author to ask for permission. For software which is copyrighted by the Free Software Foundation, write to the Free Software Foundation; we sometimes make exceptions for this. Our decision will be guided by the two goals of preserving the free status of all derivatives of our free software and of promoting the sharing and reuse of software generally.

    NO WARRANTY

  12. BECAUSE THE PROGRAM IS LICENSED FREE OF CHARGE, THERE IS NO WARRANTY FOR THE PROGRAM, TO THE EXTENT PERMITTED BY APPLICABLE LAW. EXCEPT WHEN OTHERWISE STATED IN WRITING THE COPYRIGHT HOLDERS AND/OR OTHER PARTIES PROVIDE THE PROGRAM "AS IS" WITHOUT WARRANTY OF ANY KIND, EITHER EXPRESSED OR IMPLIED, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE. THE ENTIRE RISK AS TO THE QUALITY AND PERFORMANCE OF THE PROGRAM IS WITH YOU. SHOULD THE PROGRAM PROVE DEFECTIVE, YOU ASSUME THE COST OF ALL NECESSARY SERVICING, REPAIR OR CORRECTION.

  13. IN NO EVENT UNLESS REQUIRED BY APPLICABLE LAW OR AGREED TO IN WRITING WILL ANY COPYRIGHT HOLDER, OR ANY OTHER PARTY WHO MAY MODIFY AND/OR REDISTRIBUTE THE PROGRAM AS PERMITTED ABOVE, BE LIABLE TO YOU FOR DAMAGES, INCLUDING ANY GENERAL, SPECIAL, INCIDENTAL OR CONSEQUENTIAL DAMAGES ARISING OUT OF THE USE OR INABILITY TO USE THE PROGRAM (INCLUDING BUT NOT LIMITED TO LOSS OF DATA OR DATA BEING RENDERED INACCURATE OR LOSSES SUSTAINED BY YOU OR THIRD PARTIES OR A FAILURE OF THE PROGRAM TO OPERATE WITH ANY OTHER PROGRAMS), EVEN IF SUCH HOLDER OR OTHER PARTY HAS BEEN ADVISED OF THE POSSIBILITY OF SUCH DAMAGES.

END OF TERMS AND CONDITIONS

How to Apply These Terms to Your New Programs

If you develop a new program, and you want it to be of the greatest possible use to the public, the best way to achieve this is to make it free software which everyone can redistribute and change under these terms.

To do so, attach the following notices to the program. It is safest to attach them to the start of each source file to most effectively convey the exclusion of warranty; and each file should have at least the "copyright" line and a pointer to where the full notice is found.

 
one line to give the program's name and a brief idea of what it does.
Copyright (C) yyyy  name of author

This program is free software; you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation; either version 2 of the License, or
(at your option) any later version.

This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
GNU General Public License for more details.

You should have received a copy of the GNU General Public License
along with this program; if not, write to the Free Software
Foundation, Inc., 59 Temple Place - Suite 330, Boston,
MA  02111-1307, USA.

Also add information on how to contact you by electronic and paper mail.

If the program is interactive, make it output a short notice like this when it starts in an interactive mode:

 
Gnomovision version 69, Copyright (C) 19yy name of author
Gnomovision comes with ABSOLUTELY NO WARRANTY; for details type
`show w'.
This is free software, and you are welcome to redistribute it
under certain conditions; type `show c' for details.

The hypothetical commands `show w' and `show c' should show the appropriate parts of the General Public License. Of course, the commands you use may be called something other than `show w' and `show c'; they could even be mouse-clicks or menu items--whatever suits your program.

You should also get your employer (if you work as a programmer) or your school, if any, to sign a "copyright disclaimer" for the program, if necessary. Here is a sample; alter the names:

 
Yoyodyne, Inc., hereby disclaims all copyright interest in the program
`Gnomovision' (which makes passes at compilers) written by James Hacker.

signature of Ty Coon, 1 April 1989
Ty Coon, President of Vice

This General Public License does not permit incorporating your program into proprietary programs. If your program is a subroutine library, you may consider it more useful to permit linking proprietary applications with the library. If this is what you want to do, use the GNU Library General Public License instead of this License.


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

Concept index

Jump to:   .  
A   B   C   D   E   F   G   H   I   J   K   L   M   N   O   P   Q   R   S   T   U   V   W  

Index Entry Section

.
`.msh' file8.1 Input file format
`.pre' file8.2.1 File `.pre'
`.res' file8.2.2 File `.res'

A
Acknowledgments9.3 Credits
Analytical integration3.6 Integration: defining integration methods
Approximation spaces3.4 FunctionSpace: building function spaces
Arguments2.4 Functions
Arguments, definition2.6 Arguments
Authors, e-mail9.1 Bugs
Axisymmetric, transformation3.5 Jacobian: defining jacobian methods

B
Basis Functions3.4 FunctionSpace: building function spaces
Binary operators2.3.1 Operator types
Boundary conditions3.3 Constraint: specifying constraints on function spaces and formulations
Boundary Element Method1.5 Which problems can GetDP actually solve?
Bugs, known9.1 Bugs
Bugs, reporting9.1 Bugs
Built-in functions2.4 Functions

C
Change of coordinates3.5 Jacobian: defining jacobian methods
Changelog9.2 Versions
Circuit equations3.3 Constraint: specifying constraints on function spaces and formulations
Command line options7. Running GetDP
Comments1.3 Comments
Complete examples6. Complete examples
Complex-valued, system3.8 Resolution: solving systems of equations
Concepts, indexConcept index
Conditionals2.9 Loops and conditionals
Constant, definition2.2 Constants
Constant, evaluation2.2 Constants
Constraint, definition3.3 Constraint: specifying constraints on function spaces and formulations
Constraint, examples5.4 Constraint examples
Constraint, types4.3 Types for Constraint
Contributors, list9.3 Credits
Coordinate change3.5 Jacobian: defining jacobian methods
CopyrightCopying conditions
Credits9.3 Credits
Curl2.8 Fields
Current values2.5 Current values

D
Dependences, objects1.1 Numerical tools as objects
Derivative, exterior2.8 Fields
Derivative, time3.7 Formulation: building equations
Developments, future1.5 Which problems can GetDP actually solve?
Differential operators2.8 Fields
Discrete function spaces3.4 FunctionSpace: building function spaces
Discrete quantities2.8 Fields
Discretized Geometry3.1 Group: defining topological entities
Divergence2.8 Fields
Document syntax1.2 Syntactic rules used in this document
DownloadCopying conditions

E
E-mail, authors9.1 Bugs
Edge element space, example5.8.3 Magnetostatic 3D vector potential formulation
Efficiency, tipsA. Tips and tricks
Electromagnetism1.5 Which problems can GetDP actually solve?
Electrostatic formulation5.8.1 Electrostatic scalar potential formulation
Elementary matrices3.7 Formulation: building equations
Entities, topological3.1 Group: defining topological entities
Equations3.7 Formulation: building equations
Evaluation mechanism2.2 Constants
Evaluation, order2.3.2 Evaluation order
Examples, complete6. Complete examples
Examples, short5. Short examples
Exporting results3.10 PostOperation: exporting results
Expression, definition2.1 Definition
Exterior derivative2.8 Fields

F
FAQB. Frequently asked questions
Fields2.8 Fields
File, `.msh'8.1 Input file format
File, `.pre'8.2.1 File `.pre'
File, `.res'8.2.2 File `.res'
File, comment1.3 Comments
File, include1.4 Includes
File, mesh8.1 Input file format
File, pre-processing8.2.1 File `.pre'
File, result8.2.2 File `.res'
Finite Difference Method1.5 Which problems can GetDP actually solve?
Finite Element Method1.5 Which problems can GetDP actually solve?
Finite Volume Method1.5 Which problems can GetDP actually solve?
Floating point numbers2.2 Constants
Floating potential, example5.8.2 Electrostatic scalar potential formulation with floating potentials and electric charges
Format, output3.10 PostOperation: exporting results
Formulation, definition3.7 Formulation: building equations
Formulation, electrostatics5.8.1 Electrostatic scalar potential formulation
Formulation, examples5.8 Formulation examples
Formulation, types4.7 Types for Formulation
Frequency3.8 Resolution: solving systems of equations
Frequently asked questionsB. Frequently asked questions
Function groups3.1 Group: defining topological entities
Function space, definition3.4 FunctionSpace: building function spaces
Function space, examples5.5 FunctionSpace examples
Function space, types4.4 Types for FunctionSpace
Function, definition2.4 Functions
Function, definition3.2 Function: defining global and piecewise expressions
Function, examples5.3 Function examples
Future developments1.5 Which problems can GetDP actually solve?

G
Gauss, integration3.6 Integration: defining integration methods
Geometric transformations3.5 Jacobian: defining jacobian methods
Global quantity3.7 Formulation: building equations
Global quantity, example5.8.2 Electrostatic scalar potential formulation with floating potentials and electric charges
Gmsh, examplesC. Gmsh examples
Gmsh, file format8.1 Input file format
GNU General Public LicenseD. License
Gradient2.8 Fields
Grid3.1 Group: defining topological entities
Group, definition3.1 Group: defining topological entities
Group, examples5.2 Group examples
Group, types4.1 Types for Group

H
Hierarchical basis functions3.4 FunctionSpace: building function spaces
History, versions9.2 Versions

I
Includes1.4 Includes
Index, conceptsConcept index
Index, metasyntactic variablesMetasyntactic variable index
Index, syntaxSyntax index
Input file format8.1 Input file format
Integer numbers2.2 Constants
Integral Equation Method1.5 Which problems can GetDP actually solve?
Integral quantity3.7 Formulation: building equations
Integration, definition3.6 Integration: defining integration methods
Integration, examples5.7 Integration examples
Integration, types4.6 Types for Integration
Internet addressCopying conditions
Interpolation2.8 Fields
Interpolation3.4 FunctionSpace: building function spaces
IntroductionIntroduction
Iterative loop3.8 Resolution: solving systems of equations

J
Jacobian, definition3.5 Jacobian: defining jacobian methods
Jacobian, examples5.6 Jacobian examples
Jacobian, types4.5 Types for Jacobian

K
Keywords, indexSyntax index
Known bugs9.1 Bugs

L
LicenseCopying conditions
LicenseD. License
Linear system solving3.8 Resolution: solving systems of equations
Linking, objects1.1 Numerical tools as objects
Local quantity3.7 Formulation: building equations
Loops2.9 Loops and conditionals

M
Mailing listCopying conditions
Maps3.10 PostOperation: exporting results
Matrices, elementary3.7 Formulation: building equations
Mechanics1.5 Which problems can GetDP actually solve?
Mesh3.1 Group: defining topological entities
Mesh, examplesC. Gmsh examples
Mesh, file format8.1 Input file format
Metasyntactic variables, indexMetasyntactic variable index
Method of Moments1.5 Which problems can GetDP actually solve?

N
Networks3.3 Constraint: specifying constraints on function spaces and formulations
Newmark, time scheme3.8 Resolution: solving systems of equations
Newton, nonlinear scheme3.8 Resolution: solving systems of equations
Nodal function space, example5.8.1 Electrostatic scalar potential formulation
Nonlinear system solving3.8 Resolution: solving systems of equations
Numbers, integer2.2 Constants
Numbers, real2.2 Constants
Numerical integration3.6 Integration: defining integration methods
Numerical tools, overview1. Overview

O
Objects, definition3. Objects
Objects, dependences1.1 Numerical tools as objects
Objects, overview1. Overview
Objects, types4. Types for objects
Operating system7. Running GetDP
Operation, priorities2.3.2 Evaluation order
Operators, definition2.3.1 Operator types
Operators, differential2.8 Fields
Options, command line7. Running GetDP
Order of evaluation2.3.2 Evaluation order
Output file format8.2 Output file format
Overview1. Overview

P
Parameters2.4 Functions
Philosophy, general1.1 Numerical tools as objects
Physical problems1.5 Which problems can GetDP actually solve?
Picard, nonlinear scheme3.8 Resolution: solving systems of equations
Piecewise functions2.4 Functions
Piecewise functions3.2 Function: defining global and piecewise expressions
Platforms7. Running GetDP
Post-operation, definition3.10 PostOperation: exporting results
Post-operation, examples5.11 PostOperation examples
Post-operation, types4.10 Types for PostOperation
Post-processing, definition3.9 PostProcessing: exploiting computational results
Post-processing, examples5.10 PostProcessing examples
Post-processing, types4.9 Types for PostProcessing
Priorities, operations2.3.2 Evaluation order
Processing cycle1.1 Numerical tools as objects

Q
Quantities, discrete2.8 Fields
Quantity, global3.7 Formulation: building equations
Quantity, integral3.7 Formulation: building equations
Quantity, local3.7 Formulation: building equations
Quantity, post-processing3.9 PostProcessing: exploiting computational results
Questions, frequently askedB. Frequently asked questions

R
Reading, guidelinesHow to read this manual
Real numbers2.2 Constants
Region groups3.1 Group: defining topological entities
Registers, definition2.7 Registers
Relaxation factor3.8 Resolution: solving systems of equations
Reporting bugs9.1 Bugs
Resolution, definition3.8 Resolution: solving systems of equations
Resolution, examples5.9 Resolution examples
Resolution, types4.8 Types for Resolution
Results, exploitation3.9 PostProcessing: exploiting computational results
Results, export3.10 PostOperation: exporting results
Rules, syntactic1.2 Syntactic rules used in this document
Running GetDP7. Running GetDP

S
Scope of GetDP1.5 Which problems can GetDP actually solve?
Sections3.10 PostOperation: exporting results
Short examples5. Short examples
Solving, system3.8 Resolution: solving systems of equations
Spaces, discrete3.4 FunctionSpace: building function spaces
String2.2 Constants
Symmetry, integral kernel3.7 Formulation: building equations
Syntax, indexSyntax index
Syntax, rules1.2 Syntactic rules used in this document
System, definition3.8 Resolution: solving systems of equations

T
Ternary operators2.3.1 Operator types
Thermics1.5 Which problems can GetDP actually solve?
Theta, time scheme3.8 Resolution: solving systems of equations
Time derivative3.7 Formulation: building equations
Time stepping3.8 Resolution: solving systems of equations
Time, discretization3.8 Resolution: solving systems of equations
TipsA. Tips and tricks
Tools, order of definition1.1 Numerical tools as objects
Tools, overview1. Overview
Topology3.1 Group: defining topological entities
Transformations, geometric3.5 Jacobian: defining jacobian methods
Tree3.1 Group: defining topological entities
TricksA. Tips and tricks
Types, definition4. Types for objects

U
Unary operators2.3.1 Operator types
User-defined functions3.2 Function: defining global and piecewise expressions

V
Values, current2.5 Current values
Variables, indexMetasyntactic variable index
Versions9.2 Versions

W
Web siteCopying conditions
Wiki6. Complete examples

Jump to:   .  
A   B   C   D   E   F   G   H   I   J   K   L   M   N   O   P   Q   R   S   T   U   V   W  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

Metasyntactic variable index

Jump to:   :   <   |  
A   B   C   D   E   F   G   I   J   L   M   O   P   Q   R   S   T  

Index Entry Section

:
:1.2 Syntactic rules used in this document

<
<, >1.2 Syntactic rules used in this document

|
|1.2 Syntactic rules used in this document

A
affectation2.2 Constants
argument2.6 Arguments

B
basis-function-id3.4 FunctionSpace: building function spaces
basis-function-list3.4 FunctionSpace: building function spaces
basis-function-type3.4 FunctionSpace: building function spaces
basis-function-type4.4 Types for FunctionSpace
built-in-function-id2.4 Functions

C
coef-id3.4 FunctionSpace: building function spaces
constant-def2.2 Constants
constant-id2.2 Constants
constraint-case-id3.3 Constraint: specifying constraints on function spaces and formulations
constraint-case-val3.3 Constraint: specifying constraints on function spaces and formulations
constraint-id3.3 Constraint: specifying constraints on function spaces and formulations
constraint-type3.3 Constraint: specifying constraints on function spaces and formulations
constraint-type4.3 Types for Constraint
constraint-val3.3 Constraint: specifying constraints on function spaces and formulations
coord-function-id4.2.5 Coordinate functions

D
...1.2 Syntactic rules used in this document

E
element-type3.6 Integration: defining integration methods
element-type4.6 Types for Integration
etc1.2 Syntactic rules used in this document
expression2.1 Definition
expression-char2.2 Constants
expression-cst2.2 Constants
expression-cst-list2.2 Constants
expression-cst-list-item2.2 Constants
expression-list2.1 Definition
extended-math-function-id4.2.2 Extended math functions

F
formulation-id3.7 Formulation: building equations
formulation-list3.8 Resolution: solving systems of equations
formulation-type3.7 Formulation: building equations
formulation-type4.7 Types for Formulation
function-id3.2 Function: defining global and piecewise expressions
function-space-id3.4 FunctionSpace: building function spaces
function-space-type3.4 FunctionSpace: building function spaces
function-space-type4.4 Types for FunctionSpace

G
global-quantity-id3.4 FunctionSpace: building function spaces
global-quantity-type3.4 FunctionSpace: building function spaces
global-quantity-type4.4 Types for FunctionSpace
green-function-id4.2.3 Green functions
group-def3.1 Group: defining topological entities
group-id3.1 Group: defining topological entities
group-list3.1 Group: defining topological entities
group-list-item3.1 Group: defining topological entities
group-sub-type3.1 Group: defining topological entities
group-type3.1 Group: defining topological entities
group-type4.1 Types for Group

I
integer2.2 Constants
integral-value3.9 PostProcessing: exploiting computational results
integration-id3.6 Integration: defining integration methods
integration-type3.6 Integration: defining integration methods
integration-type4.6 Types for Integration

J
jacobian-id3.5 Jacobian: defining jacobian methods
jacobian-type3.5 Jacobian: defining jacobian methods
jacobian-type4.5 Types for Jacobian

L
local-term-type3.7 Formulation: building equations
local-term-type4.7 Types for Formulation
local-value3.9 PostProcessing: exploiting computational results
loop2.9 Loops and conditionals

M
math-function-id4.2.1 Math functions
misc-function-id4.2.6 Miscellaneous functions

O
operator-binary2.3.1 Operator types
operator-ternary-left2.3.1 Operator types
operator-ternary-right2.3.1 Operator types
operator-unary2.3.1 Operator types

P
post-operation-fmt3.10 PostOperation: exporting results
post-operation-fmt4.10 Types for PostOperation
post-operation-id3.10 PostOperation: exporting results
post-operation-op3.10 PostOperation: exporting results
post-processing-id3.9 PostProcessing: exploiting computational results
post-quantity-id3.9 PostProcessing: exploiting computational results
post-quantity-type3.9 PostProcessing: exploiting computational results
post-value3.9 PostProcessing: exploiting computational results
post-value4.9 Types for PostProcessing
print-option3.10 PostOperation: exporting results
print-option4.10 Types for PostOperation
print-support3.10 PostOperation: exporting results
print-support4.10 Types for PostOperation

Q
quantity2.8 Fields
quantity-dof2.8 Fields
quantity-id2.8 Fields
quantity-operator2.8 Fields
quantity-type3.7 Formulation: building equations
quantity-type4.7 Types for Formulation

R
real2.2 Constants
register-value-get2.7 Registers
register-value-set2.7 Registers
resolution-id3.8 Resolution: solving systems of equations
resolution-op3.8 Resolution: solving systems of equations
resolution-op4.8 Types for Resolution

S
string2.2 Constants
string-id2.2 Constants
sub-space-id3.4 FunctionSpace: building function spaces
system-id3.8 Resolution: solving systems of equations
system-type3.8 Resolution: solving systems of equations

T
term-op-type3.7 Formulation: building equations
term-op-type4.7 Types for Formulation
type-function-id4.2.4 Type manipulation functions

Jump to:   :   <   |  
A   B   C   D   E   F   G   I   J   L   M   O   P   Q   R   S   T  


[ < ] [ > ]   [ << ] [ Up ] [ >> ]         [Top] [Contents] [Index] [ ? ]

Syntax index

Jump to:   !   #   $   %   &   (   *   +   -   /   0   1   2   3   <   =   >   ?   ^   |   ~  
A   B   C   D   E   F   G   H   I   J   L   N   O   P   Q   R   S   T   U   V   X   Y   Z  

Index Entry Section

!
!2.3.1 Operator types
!=2.3.1 Operator types

#
#include1.4 Includes
#integer2.7 Registers

$
$A2.5 Current values
$B2.5 Current values
$C2.5 Current values
$DTime2.5 Current values
$EigenvalueImag2.5 Current values
$EigenvalueReal2.5 Current values
$integer2.6 Arguments
$Iteration2.5 Current values
$Theta2.5 Current values
$Time2.5 Current values
$TimeStep2.5 Current values
$X2.5 Current values
$XS2.5 Current values
$Y2.5 Current values
$YS2.5 Current values
$Z2.5 Current values
$ZS2.5 Current values

%
%2.3.1 Operator types

&
&&2.3.1 Operator types

(
()2.3.2 Evaluation order

*
*2.3.1 Operator types

+
+2.3.1 Operator types

-
-2.3.1 Operator types
-2.3.1 Operator types
-adapt7. Running GetDP
-adapt7. Running GetDP
-bin7. Running GetDP
-bin7. Running GetDP
-cal7. Running GetDP
-cal7. Running GetDP
-check7. Running GetDP
-check7. Running GetDP
-help7. Running GetDP
-help7. Running GetDP
-info7. Running GetDP
-info7. Running GetDP
-ipos7. Running GetDP
-ipos7. Running GetDP
-log7. Running GetDP
-log7. Running GetDP
-msh7. Running GetDP
-msh7. Running GetDP
-name7. Running GetDP
-name7. Running GetDP
-order7. Running GetDP
-order7. Running GetDP
-p7. Running GetDP
-p7. Running GetDP
-pos7. Running GetDP
-pos7. Running GetDP
-pre7. Running GetDP
-pre7. Running GetDP
-progress7. Running GetDP
-progress7. Running GetDP
-res7. Running GetDP
-res7. Running GetDP
-restart7. Running GetDP
-restart7. Running GetDP
-socket7. Running GetDP
-socket7. Running GetDP
-solve7. Running GetDP
-solve7. Running GetDP
-split7. Running GetDP
-split7. Running GetDP
-v7. Running GetDP
-v7. Running GetDP
-verbose7. Running GetDP
-verbose7. Running GetDP
-version7. Running GetDP
-version7. Running GetDP

/
/2.3.1 Operator types
/*, */1.3 Comments
//1.3 Comments
/\2.3.1 Operator types

0
0D2.2 Constants

1
1D2.2 Constants

2
2D2.2 Constants

3
3D2.2 Constants

<
<2.3.1 Operator types
<=2.3.1 Operator types

=
=2.2 Constants
=3.1 Group: defining topological entities
=3.2 Function: defining global and piecewise expressions
==2.3.1 Operator types

>
>2.3.1 Operator types
>=2.3.1 Operator types

?
?:2.3.1 Operator types

^
^2.3.1 Operator types

|
||2.3.1 Operator types

~
~2.2 Constants

A
Acos4.2.1 Math functions
Acos4.2.1 Math functions
Adapt4.10 Types for PostOperation
Adapt4.10 Types for PostOperation
Adaptation4.10 Types for PostOperation
Adaptation4.10 Types for PostOperation
AliasOf4.4 Types for FunctionSpace
AliasOf4.4 Types for FunctionSpace
All3.5 Jacobian: defining jacobian methods
Analytic3.6 Integration: defining integration methods
Asin4.2.1 Math functions
Asin4.2.1 Math functions
Assign4.3 Types for Constraint
Assign4.3 Types for Constraint
AssignFromResolution4.3 Types for Constraint
AssignFromResolution4.3 Types for Constraint
AssociatedWith4.4 Types for FunctionSpace
AssociatedWith4.4 Types for FunctionSpace
Atan4.2.1 Math functions
Atan4.2.1 Math functions
Atan24.2.1 Math functions
Atan24.2.1 Math functions

B
BasisFunction3.4 FunctionSpace: building function spaces
BF2.8 Fields
BF2.8 Fields
BF_CurlEdge4.4 Types for FunctionSpace
BF_CurlEdge4.4 Types for FunctionSpace
BF_CurlGroupOfEdges4.4 Types for FunctionSpace
BF_CurlGroupOfEdges4.4 Types for FunctionSpace
BF_CurlGroupOfPerpendicularEdge4.4 Types for FunctionSpace
BF_CurlGroupOfPerpendicularEdge4.4 Types for FunctionSpace
BF_CurlPerpendicularEdge4.4 Types for FunctionSpace
BF_CurlPerpendicularEdge4.4 Types for FunctionSpace
BF_dGlobal4.4 Types for FunctionSpace
BF_dGlobal4.4 Types for FunctionSpace
BF_DivFacet4.4 Types for FunctionSpace
BF_DivFacet4.4 Types for FunctionSpace
BF_DivPerpendicularFacet4.4 Types for FunctionSpace
BF_DivPerpendicularFacet4.4 Types for FunctionSpace
BF_Edge4.4 Types for FunctionSpace
BF_Edge4.4 Types for FunctionSpace
BF_Facet4.4 Types for FunctionSpace
BF_Facet4.4 Types for FunctionSpace
BF_Global4.4 Types for FunctionSpace
BF_Global4.4 Types for FunctionSpace
BF_GradGroupOfNodes4.4 Types for FunctionSpace
BF_GradGroupOfNodes4.4 Types for FunctionSpace
BF_GradNode4.4 Types for FunctionSpace
BF_GradNode4.4 Types for FunctionSpace
BF_GroupOfEdges4.4 Types for FunctionSpace
BF_GroupOfEdges4.4 Types for FunctionSpace
BF_GroupOfNodes4.4 Types for FunctionSpace
BF_GroupOfNodes4.4 Types for FunctionSpace
BF_GroupOfPerpendicularEdge4.4 Types for FunctionSpace
BF_GroupOfPerpendicularEdge4.4 Types for FunctionSpace
BF_Node4.4 Types for FunctionSpace
BF_Node4.4 Types for FunctionSpace
BF_NodeX4.4 Types for FunctionSpace
BF_NodeX4.4 Types for FunctionSpace
BF_NodeY4.4 Types for FunctionSpace
BF_NodeY4.4 Types for FunctionSpace
BF_NodeZ4.4 Types for FunctionSpace
BF_NodeZ4.4 Types for FunctionSpace
BF_One4.4 Types for FunctionSpace
BF_One4.4 Types for FunctionSpace
BF_PerpendicularEdge4.4 Types for FunctionSpace
BF_PerpendicularEdge4.4 Types for FunctionSpace
BF_PerpendicularFacet4.4 Types for FunctionSpace
BF_PerpendicularFacet4.4 Types for FunctionSpace
BF_Region4.4 Types for FunctionSpace
BF_Region4.4 Types for FunctionSpace
BF_RegionX4.4 Types for FunctionSpace
BF_RegionX4.4 Types for FunctionSpace
BF_RegionY4.4 Types for FunctionSpace
BF_RegionY4.4 Types for FunctionSpace
BF_RegionZ4.4 Types for FunctionSpace
BF_RegionZ4.4 Types for FunctionSpace
BF_Volume4.4 Types for FunctionSpace
BF_Volume4.4 Types for FunctionSpace
BF_Zero4.4 Types for FunctionSpace
BF_Zero4.4 Types for FunctionSpace
Break4.8 Types for Resolution
Break4.8 Types for Resolution

C
Case3.3 Constraint: specifying constraints on function spaces and formulations
Case3.5 Jacobian: defining jacobian methods
Case3.6 Integration: defining integration methods
ChangeOfCoordinates4.10 Types for PostOperation
ChangeOfCoordinates4.10 Types for PostOperation
ChangeOfValues4.10 Types for PostOperation
ChangeOfValues4.10 Types for PostOperation
Complex4.2.4 Type manipulation functions
Complex4.2.4 Type manipulation functions
CompX4.2.4 Type manipulation functions
CompX4.2.4 Type manipulation functions
CompXX4.2.4 Type manipulation functions
CompXX4.2.4 Type manipulation functions
CompXY4.2.4 Type manipulation functions
CompXY4.2.4 Type manipulation functions
CompXZ4.2.4 Type manipulation functions
CompXZ4.2.4 Type manipulation functions
CompY4.2.4 Type manipulation functions
CompY4.2.4 Type manipulation functions
CompYX4.2.4 Type manipulation functions
CompYX4.2.4 Type manipulation functions
CompYY4.2.4 Type manipulation functions
CompYY4.2.4 Type manipulation functions
CompYZ4.2.4 Type manipulation functions
CompYZ4.2.4 Type manipulation functions
CompZ4.2.4 Type manipulation functions
CompZ4.2.4 Type manipulation functions
CompZX4.2.4 Type manipulation functions
CompZX4.2.4 Type manipulation functions
CompZY4.2.4 Type manipulation functions
CompZY4.2.4 Type manipulation functions
CompZZ4.2.4 Type manipulation functions
CompZZ4.2.4 Type manipulation functions
Constraint3.3 Constraint: specifying constraints on function spaces and formulations
Constraint3.4 FunctionSpace: building function spaces
Cos4.2.1 Math functions
Cos4.2.1 Math functions
Cosh4.2.1 Math functions
Cosh4.2.1 Math functions
Criterion3.6 Integration: defining integration methods
Cross4.2.2 Extended math functions
Cross4.2.2 Extended math functions
Curl2.8 Fields
Curl2.8 Fields
CurlInv2.8 Fields
CurlInv2.8 Fields

D
d2.8 Fields
d2.8 Fields
DecomposeInSimplex4.10 Types for PostOperation
DecomposeInSimplex4.10 Types for PostOperation
DefineConstant2.2 Constants
DefineFunction3.2 Function: defining global and piecewise expressions
DefineGroup3.1 Group: defining topological entities
Depth4.10 Types for PostOperation
Depth4.10 Types for PostOperation
deRham4.7 Types for Formulation
deRham4.7 Types for Formulation
DestinationSystem3.8 Resolution: solving systems of equations
Dimension4.10 Types for PostOperation
Dimension4.10 Types for PostOperation
dInterpolationAkima4.2.6 Miscellaneous functions
dInterpolationAkima4.2.6 Miscellaneous functions
dInterpolationLinear4.2.6 Miscellaneous functions
dInterpolationLinear4.2.6 Miscellaneous functions
dInv2.8 Fields
dInv2.8 Fields
Div2.8 Fields
Div2.8 Fields
DivInv2.8 Fields
DivInv2.8 Fields
Dof2.8 Fields
Dof2.8 Fields
Dt4.7 Types for Formulation
Dt4.7 Types for Formulation
DtDof4.7 Types for Formulation
DtDof4.7 Types for Formulation
DtDt4.7 Types for Formulation
DtDt4.7 Types for Formulation
DtDtDof4.7 Types for Formulation
DtDtDof4.7 Types for Formulation
DualEdgesOf4.1 Types for Group
DualEdgesOf4.1 Types for Group
DualFacetsOf4.1 Types for Group
DualFacetsOf4.1 Types for Group
DualNodesOf4.1 Types for Group
DualNodesOf4.1 Types for Group
DualVolumesOf4.1 Types for Group
DualVolumesOf4.1 Types for Group

E
EdgesOf4.1 Types for Group
EdgesOf4.1 Types for Group
EdgesOfTreeIn4.1 Types for Group
EdgesOfTreeIn4.1 Types for Group
EigenSolve4.8 Types for Resolution
EigenSolve4.8 Types for Resolution
EigenvalueLegend4.10 Types for PostOperation
EigenvalueLegend4.10 Types for PostOperation
ElementsOf4.1 Types for Group
ElementsOf4.1 Types for Group
Else4.8 Types for Resolution
Else4.8 Types for Resolution
EndFor2.9 Loops and conditionals
EndFor2.9 Loops and conditionals
EndIf2.9 Loops and conditionals
EndIf2.9 Loops and conditionals
Entity3.4 FunctionSpace: building function spaces
EntitySubType3.4 FunctionSpace: building function spaces
EntityType3.4 FunctionSpace: building function spaces
Equation3.7 Formulation: building equations
Evaluate4.8 Types for Resolution
Evaluate4.8 Types for Resolution
Exp4.2.1 Math functions
Exp4.2.1 Math functions

F
F_CompElementNum4.2.6 Miscellaneous functions
F_CompElementNum4.2.6 Miscellaneous functions
F_Cos_wt_p4.2.2 Extended math functions
F_Cos_wt_p4.2.2 Extended math functions
F_Period4.2.2 Extended math functions
F_Period4.2.2 Extended math functions
F_Sin_wt_p4.2.2 Extended math functions
F_Sin_wt_p4.2.2 Extended math functions
Fabs4.2.1 Math functions
Fabs4.2.1 Math functions
FacetsOf4.1 Types for Group
FacetsOf4.1 Types for Group
FacetsOfTreeIn4.1 Types for Group
FacetsOfTreeIn4.1 Types for Group
FemEquation4.7 Types for Formulation
FemEquation4.7 Types for Formulation
File4.10 Types for PostOperation
File4.10 Types for PostOperation
File4.10 Types for PostOperation
File4.10 Types for PostOperation
File4.10 Types for PostOperation
File4.10 Types for PostOperation
Fmod4.2.1 Math functions
Fmod4.2.1 Math functions
For ( expression-cst : expression-cst )2.9 Loops and conditionals
For ( expression-cst : expression-cst )2.9 Loops and conditionals
For ( expression-cst : expression-cst : expression-cst )2.9 Loops and conditionals
For ( expression-cst : expression-cst : expression-cst )2.9 Loops and conditionals
For string In { expression-cst : expression-cst : expression-cst }2.9 Loops and conditionals
For string In { expression-cst : expression-cst : expression-cst }2.9 Loops and conditionals
For string In { expression-cst : expression-cst }2.9 Loops and conditionals
For string In { expression-cst : expression-cst }2.9 Loops and conditionals
Form04.4 Types for FunctionSpace
Form04.4 Types for FunctionSpace
Form14.4 Types for FunctionSpace
Form14.4 Types for FunctionSpace
Form1P4.4 Types for FunctionSpace
Form1P4.4 Types for FunctionSpace
Form24.4 Types for FunctionSpace
Form24.4 Types for FunctionSpace
Form2P4.4 Types for FunctionSpace
Form2P4.4 Types for FunctionSpace
Form34.4 Types for FunctionSpace
Form34.4 Types for FunctionSpace
Format3.10 PostOperation: exporting results
Format4.10 Types for PostOperation
Formulation3.4 FunctionSpace: building function spaces
Formulation3.7 Formulation: building equations
FourierTransform4.8 Types for Resolution
FourierTransform4.8 Types for Resolution
Frequency3.8 Resolution: solving systems of equations
Frequency4.10 Types for PostOperation
FrequencyLegend4.10 Types for PostOperation
FrequencyLegend4.10 Types for PostOperation
Function3.2 Function: defining global and piecewise expressions
Function3.4 FunctionSpace: building function spaces
FunctionSpace3.4 FunctionSpace: building function spaces

G
Galerkin4.7 Types for Formulation
Galerkin4.7 Types for Formulation
Gauss4.6 Types for Integration
Gauss4.6 Types for Integration
GaussLegendre4.6 Types for Integration
GaussLegendre4.6 Types for Integration
Generate4.8 Types for Resolution
Generate4.8 Types for Resolution
GenerateJac4.8 Types for Resolution
GenerateJac4.8 Types for Resolution
GenerateOnly4.8 Types for Resolution
GenerateOnly4.8 Types for Resolution
GenerateOnlyJac4.8 Types for Resolution
GenerateOnlyJac4.8 Types for Resolution
GenerateSeparate4.8 Types for Resolution
GenerateSeparate4.8 Types for Resolution
GeoElement3.6 Integration: defining integration methods
Global4.1 Types for Group
Global4.1 Types for Group
Global4.7 Types for Formulation
Global4.7 Types for Formulation
GlobalEquation3.7 Formulation: building equations
GlobalQuantity3.4 FunctionSpace: building function spaces
GlobalTerm3.7 Formulation: building equations
Gmsh4.10 Types for PostOperation
Gmsh4.10 Types for PostOperation
GmshParsed4.10 Types for PostOperation
GmshParsed4.10 Types for PostOperation
Gnuplot4.10 Types for PostOperation
Gnuplot4.10 Types for PostOperation
Grad2.8 Fields
Grad2.8 Fields
GradHelmholtz4.2.3 Green functions
GradHelmholtz4.2.3 Green functions
GradInv2.8 Fields
GradInv2.8 Fields
GradLaplace4.2.3 Green functions
GradLaplace4.2.3 Green functions
Group3.1 Group: defining topological entities
Group3.4 FunctionSpace: building function spaces
GroupsOfEdgesOf4.1 Types for Group
GroupsOfEdgesOf4.1 Types for Group
GroupsOfEdgesOnNodesOf4.1 Types for Group
GroupsOfEdgesOnNodesOf4.1 Types for Group
GroupsOfNodesOf4.1 Types for Group
GroupsOfNodesOf4.1 Types for Group

H
HarmonicToTime4.10 Types for PostOperation
HarmonicToTime4.10 Types for PostOperation
Helmholtz4.2.3 Green functions
Helmholtz4.2.3 Green functions
Hexahedron4.6 Types for Integration
Hexahedron4.6 Types for Integration
Hypot4.2.2 Extended math functions
Hypot4.2.2 Extended math functions

I
If4.8 Types for Resolution
If4.8 Types for Resolution
If4.8 Types for Resolution
If4.8 Types for Resolution
If ( expression-cst )2.9 Loops and conditionals
If ( expression-cst )2.9 Loops and conditionals
Im4.2.4 Type manipulation functions
Im4.2.4 Type manipulation functions
In3.7 Formulation: building equations
In3.9 PostProcessing: exploiting computational results
Include1.4 Includes
IndexOfSystem3.7 Formulation: building equations
Init4.3 Types for Constraint
Init4.3 Types for Constraint
InitFromResolution4.3 Types for Constraint
InitFromResolution4.3 Types for Constraint
InitSolution4.8 Types for Resolution
InitSolution4.8 Types for Resolution
Integral3.9 PostProcessing: exploiting computational results
Integral4.7 Types for Formulation
Integral4.9 Types for PostProcessing
Integral4.9 Types for PostProcessing
Integration3.6 Integration: defining integration methods
Integration3.7 Formulation: building equations
Integration3.9 PostProcessing: exploiting computational results
InterpolationAkima4.2.6 Miscellaneous functions
InterpolationAkima4.2.6 Miscellaneous functions
InterpolationLinear4.2.6 Miscellaneous functions
InterpolationLinear4.2.6 Miscellaneous functions
Iso4.10 Types for PostOperation
Iso4.10 Types for PostOperation
Iso4.10 Types for PostOperation
Iso4.10 Types for PostOperation
IterativeLoop4.8 Types for Resolution
IterativeLoop4.8 Types for Resolution

J
JacNL4.7 Types for Formulation
JacNL4.7 Types for Formulation
Jacobian3.5 Jacobian: defining jacobian methods
Jacobian3.7 Formulation: building equations
Jacobian3.9 PostProcessing: exploiting computational results

L
Lanczos4.8 Types for Resolution
Lanczos4.8 Types for Resolution
Laplace4.2.3 Green functions
Laplace4.2.3 Green functions
LastTimeStepOnly4.10 Types for PostOperation
LastTimeStepOnly4.10 Types for PostOperation
Lin4.5 Types for Jacobian
Lin4.5 Types for Jacobian
Line4.6 Types for Integration
Line4.6 Types for Integration
Link4.3 Types for Constraint
Link4.3 Types for Constraint
LinkCplx4.3 Types for Constraint
LinkCplx4.3 Types for Constraint
List2.2 Constants
ListAlt2.2 Constants
Local3.9 PostProcessing: exploiting computational results
Local4.7 Types for Formulation
Local4.9 Types for PostProcessing
Local4.9 Types for PostProcessing
Log4.2.1 Math functions
Log4.2.1 Math functions
Log104.2.1 Math functions
Log104.2.1 Math functions
Loop3.7 Formulation: building equations

N
Name3.3 Constraint: specifying constraints on function spaces and formulations
Name3.4 FunctionSpace: building function spaces
Name3.5 Jacobian: defining jacobian methods
Name3.6 Integration: defining integration methods
Name3.7 Formulation: building equations
Name3.8 Resolution: solving systems of equations
Name3.9 PostProcessing: exploiting computational results
Name3.10 PostOperation: exporting results
NameOfBasisFunction3.4 FunctionSpace: building function spaces
NameOfCoef3.4 FunctionSpace: building function spaces
NameOfConstraint3.4 FunctionSpace: building function spaces
NameOfConstraint3.7 Formulation: building equations
NameOfFormulation3.8 Resolution: solving systems of equations
NameOfFormulation3.9 PostProcessing: exploiting computational results
NameOfMesh3.8 Resolution: solving systems of equations
NameOfPostProcessing3.10 PostOperation: exporting results
NameOfSpace3.7 Formulation: building equations
NameOfSystem3.9 PostProcessing: exploiting computational results
Network3.7 Formulation: building equations
Network4.3 Types for Constraint
NeverDt4.7 Types for Formulation
NeverDt4.7 Types for Formulation
Node3.7 Formulation: building equations
NodesOf4.1 Types for Group
NodesOf4.1 Types for Group
NoNewLine4.10 Types for PostOperation
NoNewLine4.10 Types for PostOperation
Norm4.2.2 Extended math functions
Norm4.2.2 Extended math functions
Normal4.2.6 Miscellaneous functions
Normal4.2.6 Miscellaneous functions
NormalSource4.2.6 Miscellaneous functions
NormalSource4.2.6 Miscellaneous functions
NumberOfPoints3.6 Integration: defining integration methods

O
OnBox4.10 Types for PostOperation
OnBox4.10 Types for PostOperation
OnElementsOf4.10 Types for PostOperation
OnElementsOf4.10 Types for PostOperation
OnGlobal4.10 Types for PostOperation
OnGlobal4.10 Types for PostOperation
OnGrid4.10 Types for PostOperation
OnGrid4.10 Types for PostOperation
OnGrid4.10 Types for PostOperation
OnGrid4.10 Types for PostOperation
OnLine4.10 Types for PostOperation
OnLine4.10 Types for PostOperation
OnPlane4.10 Types for PostOperation
OnPlane4.10 Types for PostOperation
OnPoint4.10 Types for PostOperation
OnPoint4.10 Types for PostOperation
OnRegion4.10 Types for PostOperation
OnRegion4.10 Types for PostOperation
OnSection4.10 Types for PostOperation
OnSection4.10 Types for PostOperation
Operation3.8 Resolution: solving systems of equations
Operation3.10 PostOperation: exporting results
Order4.2.6 Miscellaneous functions
Order4.2.6 Miscellaneous functions
OriginSystem3.8 Resolution: solving systems of equations

P
Pi2.2 Constants
Point4.6 Types for Integration
Point4.6 Types for Integration
PostOperation3.10 PostOperation: exporting results
PostOperation4.8 Types for Resolution
PostProcessing3.9 PostProcessing: exploiting computational results
Print3.10 PostOperation: exporting results
Print4.8 Types for Resolution
Print4.8 Types for Resolution
Print4.8 Types for Resolution
Printf4.2.6 Miscellaneous functions
Printf4.2.6 Miscellaneous functions
Prism4.6 Types for Integration
Prism4.6 Types for Integration
Pyramid4.6 Types for Integration
Pyramid4.6 Types for Integration

Q
Quadrangle4.6 Types for Integration
Quadrangle4.6 Types for Integration
Quantity3.4 FunctionSpace: building function spaces
Quantity3.7 Formulation: building equations
Quantity3.9 PostProcessing: exploiting computational results

R
Re4.2.4 Type manipulation functions
Re4.2.4 Type manipulation functions
Region3.3 Constraint: specifying constraints on function spaces and formulations
Region3.5 Jacobian: defining jacobian methods
Resolution3.4 FunctionSpace: building function spaces
Resolution3.8 Resolution: solving systems of equations
Rot2.8 Fields
Rot2.8 Fields
RotInv2.8 Fields
RotInv2.8 Fields

S
SaveSolution4.8 Types for Resolution
SaveSolution4.8 Types for Resolution
SaveSolutions4.8 Types for Resolution
SaveSolutions4.8 Types for Resolution
Scalar4.4 Types for FunctionSpace
Scalar4.4 Types for FunctionSpace
SetFrequency4.8 Types for Resolution
SetFrequency4.8 Types for Resolution
SetTime4.8 Types for Resolution
SetTime4.8 Types for Resolution
Sin4.2.1 Math functions
Sin4.2.1 Math functions
Sinh4.2.1 Math functions
Sinh4.2.1 Math functions
Skin4.10 Types for PostOperation
Skin4.10 Types for PostOperation
Smoothing4.10 Types for PostOperation
Smoothing4.10 Types for PostOperation
Solve4.8 Types for Resolution
Solve4.8 Types for Resolution
SolveJac4.8 Types for Resolution
SolveJac4.8 Types for Resolution
Solver3.8 Resolution: solving systems of equations
Sort4.10 Types for PostOperation
Sort4.10 Types for PostOperation
Sqrt4.2.1 Math functions
Sqrt4.2.1 Math functions
SquNorm4.2.2 Extended math functions
SquNorm4.2.2 Extended math functions
Store4.10 Types for PostOperation
Store4.10 Types for PostOperation
SubRegion3.3 Constraint: specifying constraints on function spaces and formulations
SubSpace3.4 FunctionSpace: building function spaces
Support3.4 FunctionSpace: building function spaces
Sur4.5 Types for Jacobian
Sur4.5 Types for Jacobian
SurAxi4.5 Types for Jacobian
SurAxi4.5 Types for Jacobian
Symmetry3.7 Formulation: building equations
System3.8 Resolution: solving systems of equations
SystemCommand4.8 Types for Resolution
SystemCommand4.8 Types for Resolution

T
Table4.10 Types for PostOperation
Table4.10 Types for PostOperation
Tan4.2.1 Math functions
Tan4.2.1 Math functions
Tanh4.2.1 Math functions
Tanh4.2.1 Math functions
Target4.10 Types for PostOperation
Target4.10 Types for PostOperation
Tensor4.2.4 Type manipulation functions
Tensor4.2.4 Type manipulation functions
TensorDiag4.2.4 Type manipulation functions
TensorDiag4.2.4 Type manipulation functions
TensorSym4.2.4 Type manipulation functions
TensorSym4.2.4 Type manipulation functions
TensorV4.2.4 Type manipulation functions
TensorV4.2.4 Type manipulation functions
Tetrahedron4.6 Types for Integration
Tetrahedron4.6 Types for Integration
TimeFunction3.3 Constraint: specifying constraints on function spaces and formulations
TimeLegend4.10 Types for PostOperation
TimeLegend4.10 Types for PostOperation
TimeLoopNewmark4.8 Types for Resolution
TimeLoopNewmark4.8 Types for Resolution
TimeLoopTheta4.8 Types for Resolution
TimeLoopTheta4.8 Types for Resolution
TimeStep4.10 Types for PostOperation
TimeStep4.10 Types for PostOperation
TimeTable4.10 Types for PostOperation
TimeTable4.10 Types for PostOperation
TransferInitSolution4.8 Types for Resolution
TransferInitSolution4.8 Types for Resolution
TransferSolution4.8 Types for Resolution
TransferSolution4.8 Types for Resolution
Transpose4.2.2 Extended math functions
Transpose4.2.2 Extended math functions
Triangle4.6 Types for Integration
Triangle4.6 Types for Integration
TTrace4.2.2 Extended math functions
TTrace4.2.2 Extended math functions
Type3.3 Constraint: specifying constraints on function spaces and formulations
Type3.4 FunctionSpace: building function spaces
Type3.6 Integration: defining integration methods
Type3.7 Formulation: building equations
Type3.8 Resolution: solving systems of equations

U
Unit4.2.2 Extended math functions
Unit4.2.2 Extended math functions
Update4.8 Types for Resolution
Update4.8 Types for Resolution
UpdateConstraint4.8 Types for Resolution
UpdateConstraint4.8 Types for Resolution
UsingPost3.10 PostOperation: exporting results

V
Value3.9 PostProcessing: exploiting computational results
Value4.10 Types for PostOperation
Vector4.2.4 Type manipulation functions
Vector4.2.4 Type manipulation functions
Vector4.4 Types for FunctionSpace
Vector4.4 Types for FunctionSpace
Vol4.5 Types for Jacobian
Vol4.5 Types for Jacobian
VolAxi4.5 Types for Jacobian
VolAxi4.5 Types for Jacobian
VolAxiRectShell4.5 Types for Jacobian
VolAxiRectShell4.5 Types for Jacobian
VolAxiSphShell4.5 Types for Jacobian
VolAxiSphShell4.5 Types for Jacobian
VolAxiSqu4.5 Types for Jacobian
VolAxiSqu4.5 Types for Jacobian
VolAxiSquRectShell4.5 Types for Jacobian
VolAxiSquRectShell4.5 Types for Jacobian
VolAxiSquSphShell4.5 Types for Jacobian
VolAxiSquSphShell4.5 Types for Jacobian
VolRectShell4.5 Types for Jacobian
VolRectShell4.5 Types for Jacobian
VolSphShell4.5 Types for Jacobian
VolSphShell4.5 Types for Jacobian
VolumesOf4.1 Types for Group
VolumesOf4.1 Types for Group

X
X4.2.5 Coordinate functions
X4.2.5 Coordinate functions
XYZ4.2.5 Coordinate functions
XYZ4.2.5 Coordinate functions

Y
Y4.2.5 Coordinate functions
Y4.2.5 Coordinate functions

Z
Z4.2.5 Coordinate functions
Z4.2.5 Coordinate functions

Jump to:   !   #   $   %   &   (   *   +   -   /   0   1   2   3   <   =   >   ?   ^   |   ~  
A   B   C   D   E   F   G   H   I   J   L   N   O   P   Q   R   S   T   U   V   X   Y   Z  


[Top] [Contents] [Index] [ ? ]

Footnotes

(1)

If you are looking for a graphical front-end to GetDP, you may consider using Gmsh (available at http://www.geuz.org/gmsh/). Gmsh permits to construct geometries, generate meshes, launch computations and visualize results directly from within a user-friendly graphical interface. The file formats used by Gmsh for mesh generation and post-processing are the default file formats accepted by GetDP (see 8.1 Input file format, and 4.10 Types for PostOperation).


[Top] [Contents] [Index] [ ? ]

Table of Contents


[Top] [Contents] [Index] [ ? ]

Short Table of Contents

Copying conditions
Introduction
1. Overview
2. Expressions
3. Objects
4. Types for objects
5. Short examples
6. Complete examples
7. Running GetDP
8. File formats
9. Bugs, versions and credits
A. Tips and tricks
B. Frequently asked questions
C. Gmsh examples
D. License
Concept index
Metasyntactic variable index
Syntax index

[Top] [Contents] [Index] [ ? ]

About this document

This document was generated by Charlie & on November, 29 2010 using texi2html

The buttons in the navigation panels have the following meaning:

Button Name Go to From 1.2.3 go to
[ < ] Back previous section in reading order 1.2.2
[ > ] Forward next section in reading order 1.2.4
[ << ] FastBack previous or up-and-previous section 1.1
[ Up ] Up up section 1.2
[ >> ] FastForward next or up-and-next section 1.3
[Top] Top cover (top) of document  
[Contents] Contents table of contents  
[Index] Index concept index  
[ ? ] About this page  

where the Example assumes that the current position is at Subsubsection One-Two-Three of a document of the following structure:

Back to geuz.org/getdp