Link between finite elements classes and EllipticProblem

We describe here how finite element classes are handled in the class VarFiniteElement. The class EllipticProblem derives from this class. The class VarFiniteElement stores a set a finite element classes. Different finite elements may be needed in the simulation for the following reasons

We are using this class to construct and store finite element objects.

Methods of class VarFiniteElement (and VarFiniteElement_Base)

AddFiniteElement1D adds one-dimensional finite element
GetTypeIntegrationBoundary retrieves the type of quadrature rules for surface integrals
ClearFiniteElement clears the finite elements stored in the class
AddSurfaceFiniteElement adds surface finite elements
AddFiniteElement adds volume finite elements
UpdateInterpolationElement updates finite elements
CopyFiniteElement copies finite elements from another class
SetFiniteElement changes the finite elements to use for a part of elements
GetReferenceElement returns the reference finite element associated with a given element
GetSurfaceFiniteElement returns the surface finite element associated with a given surface
SetAxisymGeometry informs that the 2-D mesh is used for axisymmetric computation
GetReferenceElementH1 returns the reference finite element associated with a given element (for H1 space)
GetSurfaceFiniteElementH1 returns the surface finite element associated with a given face (for H1 space)
GetReferenceElementHcurl returns the reference finite element associated with a given element (for H(curl) space)
GetSurfaceFiniteElementHcurl returns the surface finite element associated with a given face (for H(curl) space)
GetReferenceElementHdiv returns the reference finite element associated with a given element (for H(div) space)
GetNeighborReferenceElement returns the reference finite element associated with a given neighbor element (parallel)

The class VarFiniteElementEnum lists the finite elements in Montjoie and makes the link between the name (e.g. QUADRANGLE_LOBATTO) and its identification number (an integer). In regular use, this class does not need to be used, only VarFiniteElement should be used.

Methods of class VarFiniteElementEnum

GetIdentityNumber returns an identification number of the finite element from its name
InitStaticData initializes the arrays containing identification numbers
GetNewReferenceElement creates a new instance of the required finite element
GetNewSurfaceElement creates a new instance of the required surface finite element

AddFiniteElement1D

Syntax

void AddFiniteElement1D(const string& name_elt, int order, VectBool& change_elt, Mesh<Dimension1>& mesh)

This method constructs a 1-D finite element following the name given as argument. Some elements can keep the old finite element. The name of the finite element has to be chosen as detailed in the description of TypeElement keyword.

Parameters

name_elt (in)
name of the finite element
order (in)
order of approximation
change_elt (in)
if change_elt(i) is true, the element i is affected with the created finite element
mesh (in)
1-D mesh

Example :

int main()
{
  VarFiniteElement<Dimension1> var; // it can also be a derived class

  // assuming that a 1-D mesh has been created
  Mesh<Dimension1> mesh; // to complete


  Vector<bool> change_elt(mesh.GetNbElt());
  change_elt.Fill(true); // if you want to affect the same finite element for all elements
  var.AddFiniteElement1D("EDGE_LOBATTO", 4, change_elt, mesh);

  // if one element has a different order :
  int num_elem = 3;
  change_elt.Fill(false); change_elt(num_elem) = true;
  var.AddFiniteElement1D("EDGE_LOBATTO", 6, change_elt, mesh);

  // to retrieve the finite element associated with an element :
  int i = 32;
  const ElementReference<Dimension1, 1>& Fb = var.GetReferenceElementH1(i);
  
   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

GetTypeIntegrationBoundary

Syntax

void GetTypeIntegrationBoundary(const Mesh& mesh, const MeshNumbering& mesh_num, const Vector<MeshNumbering*>& other_mesh_num,
int& type_integration_edge, int& type_integration_tri, int& type_integration_quad)

This method retrieves the type of quadrature rules that are used by finite elements of the mesh numberings given as argument. If these rules are not compatible (for example an element use Gauss rules on the edge, whereas another element use Gauss-Lobatto rules), the program is aborted.

Parameters

mesh (in)
considered mesh
mesh_num (in)
main mesh numbering
other_mesh_num (in)
additional mesh numberingsx
type_integration_edge (out)
quadrature rule used for edges (2-D only)
type_integration_tri (out)
quadrature rule used for triangles (3-D only)
type_integration_quad (out)
quadrature rule used for quadrilaterals (3-D only)

Example :

int main()
{
  
    EllipticProblem<HelmholtzEquationDG<Dimension3> > var;;

    // allocating the arrays for physical indexes
    var.InitIndices(100); // we choose N = 100
    var.SetTypeEquation("HELMHOLTZ");
    
    // then you can read the data file
    ReadInputFile("test.ini", var);

    // you can continue the computations
    var.ComputeMeshAndFiniteElement("QUADRANGLE_LOBATTO");

   // you can retrieve the type of quadrature rules
   Vector<MeshNumbering<Dimension3>* > other_mesh_num;
   int type_int_edge, type_int_tri, type_int_quad;
   var.GetTypeIntegrationBoundary(var.mesh, var.GetMeshNumbering(0), other_mesh_num,
                                  type_int_edge, type_int_tri, type_int_quad);

   DISP(type_int_tri); DISP(type_int_quad);
   
   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

AddFiniteElement

Syntax

void AddFiniteElement( const string& name_elt, TinyVector<IVect, 4>& order, VectBool& change_elt,
const Mesh& mesh, const MeshNumbering& mesh_num,
int r_over_quad, int dg_form, int type = 1)
void AddFiniteElement( const string& name_elt, TinyVector<IVect, 4>& order, int type,
const VectString& other_name, Vector<TinyVector<IVect, 4>>& other_order, IVect& other_type,
VectBool& change_elt, const Mesh& mesh, const MeshNumbering& mesh_num,
Vector<MeshNumbering*>& other_mesh_num, IVect r_over_quad, int dg_form)

This method creates a set of new finite elements (with orders specified in arguments) and affects them to the required elements. The name of the finite element has to be chosen as detailed in the description of TypeElement. order(0) lists the orders of approximation to create for triangles (in 2-D) and for tetrahedra (in 3-D). order(1) lists the orders of approximation for quadrilaterals in 2-D, pyramids in 3-D. order(2) and order(3) are used for wedges and hexahedra (in 3-D only). For the class EllipticProblem, the user should not call AddFiniteElement since it is alreay called by ComputeMeshAndFiniteElement. The second syntax is more complex and handles the case with several mesh numberings (with different finite elements for each numbering).

Parameters

name_elt (in)
name of the finite element
order (in)
orders for each type of element
change_elt (in)
if change_elt(i) is true, the element i is affected with the created finite element
mesh (in)
considered mesh
mesh_num (in)
main mesh numbering
r_over_quad (in)
additional accuracy for quadrature rules
dg_form (in)
type of formulation (continuous, discontinuous or hdg)
type (in)
type of finite element (H1, H(div) or H(curl))

Example :

int main()
{
  
  // starting from a mesh
  Mesh<Dimension3> mesh;
  mesh.Read("test.mesh");

  // and a mesh numbering
  MeshNumbering<Dimension3> mesh_num(mesh);
   // you have to fill orders for each element in mesh_num (uniform or variable)

   // here we retrieve the orders present in the mesh numbering
   TinyVector<IVect, 4> order;    
   this->mesh_num.GetOrder(order);

   // no over_integration (put a positive number to over-integrate)
   int r_over_quad = 0;

   // for a continuous formulation
   int dg_form = ElementReference_Base::CONTINUOUS;

   // we want to change all the elements
   VectBool change_elt(mesh.GetNbElt()); change_elt.Fill(true);
   
   // for H^1 elements, type = 1
   VarFiniteElement<Dimension3> var;
   var.AddFiniteElement("HEXAHEDRON_LOBATTO", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

   // to retrieve the finite element for a given element
   int num_elem = 20;
   const ElementReference<Dimension3, 1>& Fb = var.GetReferenceElementH1(num_elem);
   
   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

ClearFiniteElement

Syntax

void ClearFiniteElement()

This method clear previous finite elements created.

Example :

int main()
{
  
  // starting from a mesh
  Mesh<Dimension3> mesh;
  mesh.Read("test.mesh");

  // and a mesh numbering
  MeshNumbering<Dimension3> mesh_num(mesh);
   // you have to fill orders for each element in mesh_num (uniform or variable)

   // here we retrieve the orders present in the mesh numbering
   TinyVector<IVect, 4> order;    
   this->mesh_num.GetOrder(order);

   // no over_integration (put a positive number to over-integrate)
   int r_over_quad = 0;

   // for a continuous formulation
   int dg_form = ElementReference_Base::CONTINUOUS;

   // we want to change all the elements
   VectBool change_elt(mesh.GetNbElt()); change_elt.Fill(true);
   
   // for H^1 elements, type = 1
   VarFiniteElement<Dimension3> var;
   var.AddFiniteElement("HEXAHEDRON_LOBATTO", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

   // to retrieve the finite element for a given element
   int num_elem = 20;
   const ElementReference<Dimension3, 1>& Fb = var.GetReferenceElementH1(num_elem);

   // if you want to clear finite elements
   var.ClearFiniteElement();

   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

AddSurfaceFiniteElement

Syntax

void AddSurfaceFiniteElement( const string& name_elt, TinyVector<IVect, 4>& order, VectBool& change_elt,
const Mesh& mesh, const MeshNumbering& mesh_num,
int r_over_quad, int dg_form, int type = 1)
void AddSurfaceFiniteElement( const string& name_elt, TinyVector<IVect, 4>& order, int type,
const VectString& other_name, Vector<TinyVector<IVect, 4>>& other_order, IVect& other_type,
VectBool& change_elt, const Mesh& mesh, const MeshNumbering& mesh_num,
Vector<MeshNumbering*>& other_mesh_num, IVect r_over_quad, int dg_form)

This method creates a set of new surface finite elements (with orders specified in arguments) and affects them to the required surfaces. The name of the finite element has to be chosen as detailed in the description of TypeElement. The methods creates the trace of the finite element on its boundaries. For the class EllipticProblem, the user should not call AddSurfaceFiniteElement since it is alreay called by ComputeMeshAndFiniteElement (for HDG formulation). The second syntax is more complex and handles the case with several mesh numberings (with different finite elements for each numbering).

Parameters

name_elt (in)
name of the finite element (volume finite element)
order (in)
orders for each type of surface
change_elt (in)
if change_elt(i) is true, the surface i is affected with the created finite element
mesh (in)
considered mesh
mesh_num (in)
main mesh numbering
r_over_quad (in)
additional accuracy for quadrature rules
dg_form (in)
type of formulation (continuous, discontinuous or hdg)
type (in)
type of finite element (H1, H(div) or H(curl))

Example :

int main()
{
  
  // starting from a mesh
  Mesh<Dimension3> mesh;
  mesh.Read("test.mesh");

  // and a mesh numbering
  MeshNumbering<Dimension3> mesh_num(mesh);
   // you have to fill orders for each element in mesh_num (uniform or variable)

   // here we retrieve the orders present in the mesh numbering
   TinyVector<IVect, 4> order, order_quad;    
   this->mesh_num.GetOrder(order);

   // no over_integration (put a positive number to over-integrate)
   int r_over_quad = 0;

   // for a continuous formulation
   int dg_form = ElementReference_Base::CONTINUOUS;

   // we want to change all the elements
   VectBool change_elt(mesh.GetNbElt()); change_elt.Fill(true);
   
   // for H^1 elements, type = 1
   VarFiniteElement<Dimension3> var;
   var.AddFiniteElement("HEXAHEDRON_LOBATTO", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

  // if you need also trace of finite elements (for unknown lambda in HDG)
  mesh_num.GetOrderQuadrature(order_quad);
  change_elt.Reallocate(mesh.GetNbBoundary()); change_elt.Fill(true);
  var.AddSurfaceFiniteElement("HEXAHEDRON_LOBATTO", order_quad, change_elt, mesh,
                              mesh_num, r_over_quad, dg_form, 1);

   // to retrieve the finite element for a given element
   int num_elem = 20;
   const ElementReference<Dimension3, 1>& Fb = var.GetReferenceElementH1(num_elem);

   // and surface finite element for a given surface
   int num_surf = 11;
   const ElementReference<Dimension2, 1>& Fb_s = var.GetSurfaceFiniteElementH1(num_surf);
   
   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

UpdateInterpolationElement

Syntax

void UpdateInterpolationElement(const Mesh& mesh, const MeshNumbering& mesh_num)

This method updates the interpolation operators (used for variable orders) without reconstructing finite element classes. It can be useful if the same orders exist but with a different distribution (faces are affected with different orders of quadrature).

Parameters

mesh (in)
considered mesh
mesh_num (in)
main mesh numbering

Example :

int main()
{
  
  // starting from a mesh
  Mesh<Dimension3> mesh;
  mesh.Read("test.mesh");

  // and a mesh numbering
  MeshNumbering<Dimension3> mesh_num(mesh);
   // you have to fill orders for each element in mesh_num (uniform or variable)

   // here we retrieve the orders present in the mesh numbering
   TinyVector<IVect, 4> order, order_quad;    
   this->mesh_num.GetOrder(order);

   // no over_integration (put a positive number to over-integrate)
   int r_over_quad = 0;

   // for a continuous formulation
   int dg_form = ElementReference_Base::CONTINUOUS;

   // we want to change all the elements
   VectBool change_elt(mesh.GetNbElt()); change_elt.Fill(true);
   
   // for H^1 elements, type = 1
   VarFiniteElement<Dimension3> var;
   var.AddFiniteElement("HEXAHEDRON_LOBATTO", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

  // updating finite element (for a different distribution of orders)
  var.UpdateInterpolationElement();
   
   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

CopyFiniteElement

Syntax

void CopyFiniteElement(const VarFiniteElement& var, const Mesh& mesh, const MeshNumbering& mesh_num,
const Vector<MeshNumbering*>& other_mesh_num, int type, const IVect& other_type)

This method copies the finite elements stored in another class. This copy is a shallow copy such that only pointers are copied.

Parameters

var (in)
object containing finite elements to be copied
mesh (in)
considered mesh
mesh_num (in)
main mesh numbering
other_mesh_num (in)
other mesh numberings
type (in)
type of finite element (H1, H(div) or H(curl))
other_type (in)
type of finite element for other numberings

Example :

int main()
{
  
  // starting from a mesh
  Mesh<Dimension3> mesh;
  mesh.Read("test.mesh");

  // and a mesh numbering
  MeshNumbering<Dimension3> mesh_num(mesh);
   // you have to fill orders for each element in mesh_num (uniform or variable)

   // here we retrieve the orders present in the mesh numbering
   TinyVector<IVect, 4> order, order_quad;    
   this->mesh_num.GetOrder(order);

   // no over_integration (put a positive number to over-integrate)
   int r_over_quad = 0;

   // for a continuous formulation
   int dg_form = ElementReference_Base::CONTINUOUS;

   // we want to change all the elements
   VectBool change_elt(mesh.GetNbElt()); change_elt.Fill(true);
   
   // for H^1 elements, type = 1
   VarFiniteElement<Dimension3> var;
   var.AddFiniteElement("HEXAHEDRON_LOBATTO", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

   // if you want to copy finite elements to another class var2
   VarFiniteElement<Dimension3> var2;
   Vector<MeshNumbering<Dimension3>*> other_mesh_num; IVect other_type;
   var2.CopyFiniteElement(var, mesh, mesh_num, other_mesh_num, 1, other_type);
   
   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

SetFiniteElement

Syntax

void SetFiniteElement(const string& name_elt, const VectBool& change_elt, int type, int n = 0)

This method sets the finite elements to use for elements such that change_elt(i) is true. In this method, we assume that finite elements have been previously created (with AddFiniteElement), and we change only the pointers.

Parameters

name_elt (in)
name of the finite element to affect to selected elements
change_elt (in)
if change_elt(i) is true, the element i is affected with finite element name_elt
type (in)
type of finite element (H1, H(div) or H(curl))
n (optional)
number of the mesh numbering

Example :

int main()
{
  
  // starting from a mesh
  Mesh<Dimension3> mesh;
  mesh.Read("test.mesh");

  // and a mesh numbering
  MeshNumbering<Dimension3> mesh_num(mesh);
   // you have to fill orders for each element in mesh_num (uniform or variable)

   // here we retrieve the orders present in the mesh numbering
   TinyVector<IVect, 4> order, order_quad;    
   this->mesh_num.GetOrder(order);

   // no over_integration (put a positive number to over-integrate)
   int r_over_quad = 0;

   // for a continuous formulation
   int dg_form = ElementReference_Base::DISCONTINUOUS;

   // we want to change all the elements
   VectBool change_elt(mesh.GetNbElt()); change_elt.Fill(true);
   
   // for H^1 elements, type = 1
   VarFiniteElement<Dimension3> var;
   var.AddFiniteElement("HEXAHEDRON_LOBATTO", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

   // if you want to add other finite elements
   VarFiniteElement<Dimension3> var;
   var.AddFiniteElement("TETRAHEDRON_HIERARCHIC", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

   // you can change some elements to get back to hexahedron_lobatto
   // no finite elements are created, only pointers are modified
   change_elt(3) = false; change_elt(5) = false;
   var.SetFiniteElement("HEXAHEDRON_LOBATTO", change_elt, 1);
   
   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

GetReferenceElement, GetReferenceElementH1, GetReferenceElementHcurl, GetReferenceElementHdiv

Syntax

const ElementReference_Dim& GetReferenceElement(int i, int n = 0)
const ElementReference<Dimension, 1>& GetReferenceElementH1(int i, int n = 0)
const ElementReference<Dimension, 2>& GetReferenceElementHcurl(int i, int n = 0)
const ElementReference<Dimension, 3>& GetReferenceElementHdiv(int i, int n = 0)
void GetReferenceElement(Vector<const ElementReference_Dim*>& stored_elt)

This method returns the finite element associated with the element i. The second argument is optional and is used if there are several mesh numberings. The methods GetReferenceElementH1, GetReferenceElementHcurl, GetReferenceElementHdiv returns the finite element casted with the good type. The last syntax of GetReferenceElement allows to retrieve all the finite elements stored in the class.

Parameters

i (in)
element number
n (optional)
number of the mesh numbering

Example :

int main()
{
  
  // starting from a mesh
  Mesh<Dimension3> mesh;
  mesh.Read("test.mesh");

  // and a mesh numbering
  MeshNumbering<Dimension3> mesh_num(mesh);
   // you have to fill orders for each element in mesh_num (uniform or variable)

   // here we retrieve the orders present in the mesh numbering
   TinyVector<IVect, 4> order, order_quad;    
   this->mesh_num.GetOrder(order);

   // no over_integration (put a positive number to over-integrate)
   int r_over_quad = 0;

   // for a continuous formulation
   int dg_form = ElementReference_Base::DISCONTINUOUS;

   // we want to change all the elements
   VectBool change_elt(mesh.GetNbElt()); change_elt.Fill(true);
   
   // for H^1 elements, type = 1
   VarFiniteElement<Dimension3> var;
   var.AddFiniteElement("HEXAHEDRON_LOBATTO", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

   // loop over elements of the mesh
   for (int i = 0; i < mesh.GetNbElt(); i++)
     {
       // we retrieve the finite element associated with element i
       const ElementReference_Dim<Dimension3>& Fb = var.GetReferenceElement(i);

       // if you know that the element is scalar, you can retrieve the leaf class
       const ElementReference<Dimension3, 1>& Fb_s = var.GetReferenceElementH1(i);

       // GetReferenceElementHcurl for edge elements, GetReferenceElementHdiv for facet elements
     }

   // second syntax is less common, and serves to retrieve all the finite elements stored in the class
   // if there is only one finite element (e.g. the mesh is purely hexahedral and the order is uniform), the array will be of length 1
   Vector<const ElementReferece_Dim<Dimension3>* > all_Fb;
   var.GetReferenceElement(all_Fb);

   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

GetSurfaceFiniteElement, GetSurfaceFiniteElementH1, GetSurfaceFiniteElementHcurl

Syntax

const ElementReference_Dim& GetSurfaceFiniteElement(int i, int n = 0)
const ElementReference<Dimension, 1>& GetSurfaceFiniteElementH1(int i, int n = 0)
const ElementReference<Dimension, 2>& GetSurfaceFiniteElementHcurl(int i, int n = 0)
void GetSurfaceFiniteElement(Vector<const ElementReference_Dim*>& stored_elt)

This method returns the surface finite element associated with the face i. The second argument is optional and is used if there are several mesh numberings. The last syntax allows to retrieve all the surface finite elements stored in the class.

Parameters

i (in)
face number
n (optional)
number of the mesh numbering

Example :

int main()
{
  
  // starting from a mesh
  Mesh<Dimension3> mesh;
  mesh.Read("test.mesh");

  // and a mesh numbering
  MeshNumbering<Dimension3> mesh_num(mesh);
   // you have to fill orders for each element in mesh_num (uniform or variable)

   // here we retrieve the orders present in the mesh numbering
   TinyVector<IVect, 4> order, order_quad;    
   this->mesh_num.GetOrder(order);

   // no over_integration (put a positive number to over-integrate)
   int r_over_quad = 0;

   // for hdg formulation
   int dg_form = ElementReference_Base::HDG;

   // we want to change all the elements
   VectBool change_elt(mesh.GetNbElt()); change_elt.Fill(true);
   
   // for H^1 elements, type = 1
   VarFiniteElement<Dimension3> var;
   var.AddFiniteElement("HEXAHEDRON_LOBATTO", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

   // if you need also trace of finite elements (for unknown lambda in HDG)
   mesh_num.GetOrderQuadrature(order_quad);
   change_elt.Reallocate(mesh.GetNbBoundary()); change_elt.Fill(true);
   var.AddSurfaceFiniteElement("HEXAHEDRON_LOBATTO", order_quad, change_elt, mesh,
                               mesh_num, r_over_quad, dg_form, 1);

   // loop over faces of the mesh
   for (int i = 0; i < mesh.GetNbBoundary(); i++)
     {
       // we retrieve the finite element associated with face i
       const ElementReference_Dim<Dimension2>& Fb = var.GetSurfaceFiniteElement(i);

       // if you know that the element is scalar, you can retrieve the leaf class
       const ElementReference<Dimension2, 1>& Fb_s = var.GetSurfaceFiniteElementH1(i);

       // GetSurfaceFiniteElementHcurl for edge elements
     }

   // second syntax is less common, and serves to retrieve all the surface finite elements stored in the class
   // if there is only one finite element (e.g. the mesh is purely hexahedral and the order is uniform), the array will be of length 1
   Vector<const ElementReferece_Dim<Dimension2>* > all_Fb;
   var.GetSurfaceFiniteElement(all_Fb);

   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

SetAxisymGeometry

Syntax

void SetAxisymGeometry(bool axisym)

This method is called when the solved problem is axisymmetric. This is useful because the curvatures are different from a 2-D planar configuration.

Example :

int main()
{
  
  // starting from a mesh
  Mesh<Dimension2> mesh;
  mesh.Read("test.mesh");

  // and a mesh numbering
  MeshNumbering<Dimension2> mesh_num(mesh);
   // you have to fill orders for each element in mesh_num (uniform or variable)

   // here we retrieve the orders present in the mesh numbering
   TinyVector<IVect, 4> order, order_quad;    
   this->mesh_num.GetOrder(order);

   // no over_integration (put a positive number to over-integrate)
   int r_over_quad = 0;

   // for a continuous formulation
   int dg_form = ElementReference_Base::CONTINUOUS;

   // we want to change all the elements
   VectBool change_elt(mesh.GetNbElt()); change_elt.Fill(true);

   // the computation will be used for axisymmetric geometry
   VarFiniteElement<Dimension2> var;
   var.SetAxisymGeometry(true); // not needed for the common case

   // for H^1 elements, type = 1
   var.AddFiniteElement("QUADRANGLE_LOBATTO", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

GetNeighborReferenceElement

Syntax

const ElementReference_Dim& GetNeighborReferenceElement(int i, int n = 0)

This method returns the finite element associated with the neighboring element of face i. It is useful for parallel computations if you need to access to the finite element that is associated with the processor owing the neighbor element.

Parameters

i (in)
local face number
n (optional)
number of the mesh numbering

Example :

int main()
{
  
  // starting from a mesh
  Mesh<Dimension3> mesh;
  mesh.Read("test.mesh");

  // and a mesh numbering
  MeshNumbering<Dimension3> mesh_num(mesh);
   // you have to fill orders for each element in mesh_num (uniform or variable)

   // here we retrieve the orders present in the mesh numbering
   TinyVector<IVect, 4> order, order_quad;    
   this->mesh_num.GetOrder(order);

   // no over_integration (put a positive number to over-integrate)
   int r_over_quad = 0;

   // for a continuous formulation
   int dg_form = ElementReference_Base::DISCONTINUOUS;

   // we want to change all the elements
   VectBool change_elt(mesh.GetNbElt()); change_elt.Fill(true);
   
   // for H^1 elements, type = 1
   VarFiniteElement<Dimension3> var;
   var.AddFiniteElement("HEXAHEDRON_LOBATTO", order, change_elt, mesh, mesh_num,
                        r_over_quad, dg_form, 1);

   // loop over referenced faces of the mesh
   for (int i = 0; i < mesh.GetNbBoundaryRef(); i++)
     {
       // do we have a face on the interface between two processors ?
       int ref = mesh.BoundaryRef(i).GetReference();
       int cond = mesh.GetBoundaryCondition(ref);
       bool neighbor_face = (cond == BoundaryConditionEnum::LINE_NEIGHBOR);

       if (neighbor_face)
         {
           // local number of the face (between neighboring faces)
           int pos_loc = mesh_num.GetLocalEdgeNumberNeighborElement(i);
           
           // we retrieve the finite element associated with the neighbor element
           const ElementReference_Dim<Dimension3>& Fb = var.GetNeighborReferenceElement(pos_loc);
          }
     }

   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

GetIdentityNumber

Syntax

int GetIdentityNumber(string name, int type)

This method returns the integer associated with the name of a finite element. For a given dimension and type of finite element, the identification number is unique.

Example :

int main(int argc, char** argv)
{
  InitMontjoie(argc, argv);

  int id = VarFiniteElementEnum<Dimension2, 1>::GetIdentityNumber("QUADRANGLE_LOBATTO", 1);

  // then you can create a finite element with this id
  ElementReference<Dimension2, 1>* Fb = VarFiniteElementEnum<Dimension2, 1>::GetNewReferenceElement(id);

  // and constructs the finite element for a given order
   Fb->ConstructFiniteElement(4);

   return FinalizeMontjoie();
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

InitStaticData

Syntax

int InitStaticData()

This method initializes the mapping between name of finite elements and identification numbers. This method is automatically called by InitMontjoie, and should not be called in regular use.

Example :

int main(int argc, char** argv)
{
  // you do not call  InitMontjoie(argc, argv);
  // then you need to call InitStaticData
  VarFiniteElementEnum<Dimension2, 1>::InitStaticData();

  int id = VarFiniteElementEnum<Dimension2, 1>::GetIdentityNumber("QUADRANGLE_LOBATTO", 1);

  // then you can create a finite element with this id
  ElementReference<Dimension2, 1>* Fb = VarFiniteElementEnum<Dimension2, 1>::GetNewReferenceElement(id);

  // and constructs the finite element for a given order
   Fb->ConstructFiniteElement(4);

   return 0;
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

GetNewReferenceElement

Syntax

ElementReference* GetNewReferenceElement(int id)

This method creates a finite element object by giving its identification number. It is up to the user to delete the created object.

Example :

int main(int argc, char** argv)
{
  InitMontjoie(argc, argv);

  int id = VarFiniteElementEnum<Dimension2, 1>::GetIdentityNumber("QUADRANGLE_LOBATTO", 1);

  // then you can create a finite element with this id
  ElementReference<Dimension2, 1>* Fb = VarFiniteElementEnum<Dimension2, 1>::GetNewReferenceElement(id);

  // and constructs the finite element for a given order
   Fb->ConstructFiniteElement(4);


  // you have to release memory when you no longer need the object
  delete Fb;

   return FinalizeMontjoie();
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx

GetNewSurfaceElement

Syntax

ElementReference* GetNewSurfaceElement(int r, int r_geom, int r_quad, int id)

This method creates a surface finite element object by giving its identification number. The created surface finite element is the trace of the finite element on a face. The surface finite element is constructed. It is up to the user to delete the created object.

Parameters

r (in)
order of approximation
r_geom (in)
order of geometry approximation
r_quad (in)
order of quadrature rules
id (in)
order of approximation

Example :

int main(int argc, char** argv)
{
  InitMontjoie(argc, argv);

  int id = VarFiniteElementEnum<Dimension2, 1>::GetIdentityNumber("QUADRANGLE_LOBATTO", 1);

  // then you can create (and construct) a surface finite element with this id
  ElementReference<Dimension1, 1>* Fb = VarFiniteElementEnum<Dimension2, 1>::GetNewSurfaceElement(4, 4, 4, id);

  // you have to release memory when you no longer need the object
  delete Fb;

   return FinalizeMontjoie();
}

  

Location :

Harmonic/VarFiniteElement.hxx
Harmonic/VarFiniteElement.cxx