NGSolve  4.9
Classes | Enumerations | Functions
ngcomp Namespace Reference

namespace for NGS-components. More...

Classes

class  BilinearForm
 A bilinear-form. More...
class  S_BilinearForm
 We specify the scalar (double or Complex) of the biform. More...
class  T_BilinearForm
class  T_BilinearFormSymmetric
class  T_BilinearFormDiagonal
class  BilinearFormApplication
 This objects provide the bilinear-form application as matrix vector product. More...
class  LinearizedBilinearFormApplication
 Applies the matrix-vector product of linearized matrix. More...
class  ElementByElement_BilinearForm
 This bilinearform stores the element-matrices. More...
class  FacetFESpace
class  FESpace
 Base class for finite element space. More...
class  NodalFESpace
 A space of continuous finite elements. More...
class  NonconformingFESpace
class  ElementFESpace
class  SurfaceElementFESpace
 Non-continous fe space on boundary. More...
class  CompoundFESpace
 A combination of fe-spaces. More...
class  FESpaceClasses
 Registered FESpace classes. More...
class  RegisterFESpace
 template for registration of finite element spaces. More...
class  ParallelMeshDofs
class  GridFunction
 Grid-functions. More...
class  S_GridFunction
class  T_GridFunction
class  S_ComponentGridFunction
class  GridFunctionCoefficientFunction
class  VisualizeGridFunction
class  VisualizeCoefficientFunction
class  H1HighOrderFESpace
 High Order Finite Element Space. More...
class  NedelecFESpace
 Lowest order Nedelec space (edge elements) More...
class  NedelecFESpace2
class  HCurlHighOrderFESpace
 HCurl High Order Finite Element Space. More...
class  RaviartThomasFESpace
 Finite Element Space for H(div) More...
class  HDivHighOrderFESpace
 HDiv High Order Finite Element Space. More...
class  L2HighOrderFESpace
 High Order Finite Element Space for L2 (element by element) More...
class  L2SurfaceHighOrderFESpace
class  LinearForm
 Linearform. More...
class  S_LinearForm
class  T_LinearForm
 Template argument specifies vector type. More...
class  MeshAccess
 Access to mesh topology and geometry. More...
class  ProgressOutput
 Controls the progress - output. More...
class  NGS_Object
 NGSolve base class. More...
class  TCreateVecObjectS
class  TCreateVecObjectS< Object, Base, SCAL, ARG, 1 >
class  TCreateVecObject3S
class  TCreateVecObject3S< Object, Base, SCAL, ARG, ARG2, ARG3, 1 >
class  TCreateSymMatObject3S
class  TCreateSymMatObject3S< Object, Base, SCAL, ARG, ARG2, ARG3, 1 >
class  Preconditioner
 Base class for preconditioners. More...
class  MGPreconditioner
 Multigrid preconditioner. More...
class  LocalPreconditioner
 Local (Block-Jacobi or Block-Gauss-Seidel) preconditioner. More...
class  TwoLevelPreconditioner
class  ComplexPreconditioner
class  ChebychevPreconditioner
class  CommutingAMGPreconditioner
class  NonsymmetricPreconditioner
class  PreconditionerClasses
 Registered Preconditioner classes. More...
class  RegisterPreconditioner
class  VectorFacetFESpace

Enumerations

enum  TRANSFORM_TYPE {
  TRANSFORM_MAT_LEFT = 1, TRANSFORM_MAT_RIGHT = 2, TRANSFORM_MAT_LEFT_RIGHT = 3, TRANSFORM_RHS = 4,
  TRANSFORM_SOL = 8
}
 transformation from local to global orientation used for low order Nedelec elements
enum  COUPLING_TYPE {
  UNUSED_DOF = 0, LOCAL_DOF = 1, INTERFACE_DOF = 2, NONWIREBASKET_DOF = 3,
  WIREBASKET_DOF = 4, EXTERNAL_DOF = 6, ANY_DOF = 7
}
 coupling types: Each degree of freedom is either More...

Functions

NGS_DLL_HEADER BilinearFormCreateBilinearForm (const FESpace *space, const string &name, const Flags &flags)
 Allocates a bilinearform.
NGS_DLL_HEADER FESpaceClassesGetFESpaceClasses ()
 returns createion object
NGS_DLL_HEADER FESpaceCreateFESpace (const string &type, const MeshAccess &ma, const Flags &flags)
 creates a fespace of that type
NGS_DLL_HEADER GridFunctionCreateGridFunction (const FESpace *space, const string &name, const Flags &flags)
NGS_DLL_HEADER LinearFormCreateLinearForm (const FESpace *space, const string &name, const Flags &flags)
ELEMENT_TYPE ConvertElementType (NG_ELEMENT_TYPE type)
 Converts element-type from Netgen to element-types of NGSolve.
template<typename T >
void AllReduceNodalData (NODE_TYPE nt, Array< T > &data, MPI_Op op, const MeshAccess &ma)
 old style, compatibility for a while
template<template< class T > class Object, class Base , class ARG >
Base * CreateVecObject (int dim, bool iscomplex, ARG &arg)
template<template< class T > class Object, class Base , class ARG , class ARG2 , class ARG3 >
Base * CreateVecObject (int dim, bool iscomplex, ARG &arg, ARG2 &arg2, ARG3 &arg3)
template<template< class T, class TV > class Object, class Base , class ARG , class ARG2 , class ARG3 >
Base * CreateSymMatObject (int dim, bool iscomplex, ARG &arg, ARG2 &arg2, ARG3 &arg3)
template<class SCAL >
NGS_DLL_HEADER void CalcFlux (const MeshAccess &ma, const S_GridFunction< SCAL > &u, S_GridFunction< SCAL > &flux, const BilinearFormIntegrator &bli, bool applyd, bool add, int domain)
NGS_DLL_HEADER void CalcFluxProject (const MeshAccess &ma, const GridFunction &u, GridFunction &flux, const BilinearFormIntegrator &bli, bool applyd, int domain, LocalHeap &lh)
template<class SCAL >
NGS_DLL_HEADER void CalcFluxProject (const MeshAccess &ma, const S_GridFunction< SCAL > &u, S_GridFunction< SCAL > &flux, const BilinearFormIntegrator &bli, bool applyd, const BitArray &domains, LocalHeap &lh)
NGS_DLL_HEADER void SetValues (const MeshAccess &ma, const CoefficientFunction &coef, GridFunction &u, bool bound, DifferentialOperator *diffop, LocalHeap &clh)
template<class SCAL >
NGS_DLL_HEADER int CalcPointFlux (const MeshAccess &ma, const GridFunction &u, const FlatVector< double > &point, FlatVector< SCAL > &flux, const BilinearFormIntegrator &bli, bool applyd, LocalHeap &lh, int component=0)
template<class SCAL >
NGS_DLL_HEADER int CalcPointFlux (const MeshAccess &ma, const GridFunction &u, const FlatVector< double > &point, const Array< int > &domains, FlatVector< SCAL > &flux, const BilinearFormIntegrator &bli, bool applyd, LocalHeap &lh, int component=0)
NGS_DLL_HEADER void CalcError (const MeshAccess &ma, const GridFunction &bu, const GridFunction &bflux, const BilinearFormIntegrator &bli, FlatVector< double > &err, int domain, LocalHeap &lh)
template<class SCAL >
NGS_DLL_HEADER void CalcError (const MeshAccess &ma, const S_GridFunction< SCAL > &u, const S_GridFunction< SCAL > &flux, const BilinearFormIntegrator &bli, FlatVector< double > &err, const BitArray &domains, LocalHeap &lh)
template<class SCAL >
NGS_DLL_HEADER void CalcDifference (const MeshAccess &ma, const S_GridFunction< SCAL > &u1, const S_GridFunction< SCAL > &u2, const BilinearFormIntegrator &bli1, const BilinearFormIntegrator &bli2, FlatVector< double > &diff, int domain, LocalHeap &lh)
NGS_DLL_HEADER void CalcDifference (const MeshAccess &ma, const GridFunction &u1, const BilinearFormIntegrator &bli1, const CoefficientFunction *coef, FlatVector< double > &diff, int domain, LocalHeap &lh)
template<class SCAL >
NGS_DLL_HEADER void CalcGradient (const MeshAccess &ma, const FESpace &fesh1, const S_BaseVector< SCAL > &vech1, const FESpace &feshcurl, S_BaseVector< SCAL > &vechcurl)
template<class SCAL >
NGS_DLL_HEADER void CalcGradientT (const MeshAccess &ma, const FESpace &feshcurl, const S_BaseVector< SCAL > &vechcurl, const FESpace &fesh1, S_BaseVector< SCAL > &vech1)
template<class SCAL >
NGS_DLL_HEADER void CalcErrorHierarchical (const MeshAccess &ma, const S_BilinearForm< SCAL > &bfa, const S_BilinearForm< SCAL > &bfa2, const S_LinearForm< SCAL > &lff, S_GridFunction< SCAL > &gfu, const FESpace &festest, FlatVector< double > &err, LocalHeap &lh)
NGS_DLL_HEADER
PreconditionerClasses
GetPreconditionerClasses ()

Detailed Description

namespace for NGS-components.

High Order Finite Element Space for Facets.

It contains the access to the mesh topology MeshAccess, finite element spaces derived from FESpace, and global objects such as Gridfunction, LinearForm, and BilinearForm, or a Preconditioner. Computational functions such as SetValues and CalcFlux.


Enumeration Type Documentation

coupling types: Each degree of freedom is either

  • a local degree of freedom
  • an interface degree of freedom or
  • a wirebasket degree of freedom

Function Documentation

ELEMENT_TYPE ngcomp::ConvertElementType ( NG_ELEMENT_TYPE  type) [inline]

Converts element-type from Netgen to element-types of NGSolve.

E.g. the Netgen-types NG_TRIG and NG_TRIG6 are merged to NGSolve type ET_TRIG.

NGS_DLL_HEADER BilinearForm* ngcomp::CreateBilinearForm ( const FESpace *  space,
const string &  name,
const Flags &  flags 
)

Allocates a bilinearform.

Some flags are: -symmetric ... assembles a symmetric matrix