// ============================================================================ // == == // == AMDiS - Adaptive multidimensional simulations == // == == // ============================================================================ // == == // == TU Dresden == // == == // == Institut für Wissenschaftliches Rechnen == // == Zellescher Weg 12-14 == // == 01069 Dresden == // == germany == // == == // ============================================================================ // == == // == https://gforge.zih.tu-dresden.de/projects/amdis/ == // == == // ============================================================================ /** \file BasisFunction.h */ #ifndef AMDIS_BASISFUNCTION_H #define AMDIS_BASISFUNCTION_H #include #include "AMDiS_fwd.h" #include "Global.h" #include "Boundary.h" #include "MatrixVector.h" namespace AMDiS { /// Function interface for evaluating basis functions. class BasFctType { public: BasFctType() {} virtual ~BasFctType() {} virtual double operator()(const DimVec&) const = 0; }; /// Function interface for evaluating gradients of basis functions. class GrdBasFctType { public: GrdBasFctType() {} virtual ~GrdBasFctType() {} virtual void operator()(const DimVec&, DimVec&) const = 0; }; /// Function interface for evaluating second derivative of basis functions. class D2BasFctType { public: D2BasFctType() {} virtual ~D2BasFctType() {} virtual void operator()(const DimVec&, DimMat&) const = 0; }; typedef BasFctType *BFptr; typedef GrdBasFctType *GBFptr; typedef D2BasFctType *DBFptr; /** \ingroup FEMSpace * \brief * Base class for finite element basis functions. In order to build up a * finite element space, we have to specify a set of local basis functions. * Together with the correspondig DOF administration and the underlying mesh, * the finite element space is given. * This class holds the local basis functions and their derivatives of the * reference element. They are evaluated at barycentric coordinates, so they * can be used on every element of the mesh. */ class BasisFunction { protected: /// Creates a BasisFunction object of given dim and degree BasisFunction(std::string name, int dim, int degree); /// destructor virtual ~BasisFunction(); public: /// compares two BasisFunction objects. virtual bool operator==(const BasisFunction& a) const { return a.getName() == name; } /// Returns !(*this == b) inline bool operator!=(const BasisFunction& b) const { return !(operator == (b)); } /// Used by \ref getDOFIndices and \ref getVec virtual int* orderOfPositionIndices(const Element* el, GeoIndex position, int positionIndex) const = 0; /** \brief * Pointer to a function which connects the set of local basis functions * with its global DOFs. * getDOFIndices(el, admin, dof) returns a pointer to a const vector of * length \ref nBasFcts where the i-th entry is the index of the DOF * associated to the i-th basis function; arguments are the actual element * el and the DOF admin admin of the corresponding finite element space * (these indices depend on all defined DOF admins and thus on the * corresponding admin); if the last argument dof is NULL, getDOFndices * has to provide memory for storing this vector, which is overwritten on the * next call of getDOFIndices; if dof is not NULL, dof is a pointer to a * vector which has to be filled; */ virtual const DegreeOfFreedom* getDOFIndices(const Element*, const DOFAdmin&, DegreeOfFreedom*) const = 0; /** \brief * The second argument 'bound' has to be a pointer to a vector which has * to be filled. Its length is \ref nBasFcts (the number of basis functions * in the used finite element space). After calling this function, the i-th * entry of the array is the boundary type of the i-th basis function of this * element. * * This function needs boundary information within the ElInfo object; thus, * all routines using this function on the elements need the FILL_BOUND * flag during mesh traversal; */ virtual void getBound(const ElInfo*, BoundaryType *) const {} /// Returns \ref degree of BasisFunction inline const int getDegree() const { return degree; } /// Returns \ref dim of BasisFunction inline const int getDim() const { return dim; } /// Returns \ref nBasFcts which is the number of local basis functions inline const int getNumber() const { return nBasFcts; } /// Returns \ref name of BasisFunction inline std::string getName() const { return name; } /// Returns \ref nDOF[i] int getNumberOfDOFs(int i) const; /// Returns \ref nDOF inline DimVec* getNumberOfDOFs() const { return nDOF; } /// Initialisation of the \ref nDOF vector. Must be implemented by sub classes virtual void setNDOF() = 0; /// Returns the barycentric coordinates of the i-th basis function. virtual DimVec *getCoords(int i) const = 0; /** \brief * Returns a pointer to a const vector with interpolation coefficients of the * function f; if indices is a pointer to NULL, the coefficient for all * basis functions are calculated and the i-th entry in the vector is the * coefficient of the i-th basis function; if indices is non NULL, only the * coefficients for a subset of the local basis functions has to be * calculated; n is the number of those basis functions, indices[0], . . . * , indices[n-1] are the local indices of the basis functions where the * coefficients have to be calculated, and the i-th entry in the return * vector is then the coefficient of the indices[i]-th basis function; coeff * may be a pointer to a vector which has to be filled * (compare the dof argument of \ref getDOFIndices()); * such a function usually needs vertex coordinate information; thus, all * routines using this function on the elements need the FILL COORDS flag * during mesh traversal. * Must be implemented by sub classes. */ virtual const double* interpol(const ElInfo *el_info, int n, const int *indices, AbstractFunction > *f, double *coeff) = 0; /// WorldVector valued interpol function. virtual const WorldVector* interpol(const ElInfo *el_info, int no, const int *b_no, AbstractFunction, WorldVector > *f, WorldVector *vec) = 0; /// Returns the i-th local basis function inline BasFctType *getPhi(int i) const { return (*phi)[i]; } /// Returns the gradient of the i-th local basis function inline GrdBasFctType *getGrdPhi(int i) const { return (*grdPhi)[i]; } /// Returns the second derivative of the i-th local basis function inline D2BasFctType *getD2Phi(int i) const { return (*d2Phi)[i]; } /** \brief * Approximates the L2 scalar products of a given function with all basis * functions by numerical quadrature and adds the corresponding values to a * DOF vector; * f is a pointer for the evaluation of the given function in world * coordinates x and returns the value of that function at x; if f is a NULL * pointer, nothing is done; * fh is the DOF vector where at the i-th entry the approximation of the L2 * scalar product of the given function with the i-th global basis function * of fh->feSpace is stored; * quad is the quadrature for the approximation of the integral on each leaf * element of fh->feSpace->mesh; if quad is a NULL pointer, a default * quadrature which is exact of degree 2*fh->feSpace->basFcts->degree-2 is * used. * The integrals are approximated by looping over all leaf elements, * computing the approximations to the element contributions and adding * these values to the vector fh by add element vec(). * The vector fh is not initialized with 0.0; only the new contributions are * added */ virtual void l2ScpFctBas(Quadrature*, AbstractFunction >* /*f*/, DOFVector* /*fh*/) {} /// WorldVector valued l2ScpFctBas function virtual void l2ScpFctBas(Quadrature* , AbstractFunction, WorldVector >* /*f*/, DOFVector >* /*fh*/) {} /// Interpolates a DOFIndexed after refinement virtual void refineInter(DOFIndexed *, RCNeighbourList*, int) {} /// Interpolates a DOFIndexed after coarsening virtual void coarseInter(DOFIndexed *, RCNeighbourList*, int) {} /// Restricts a DOFIndexed after coarsening virtual void coarseRestr(DOFIndexed *, RCNeighbourList*, int) {} /// Interpolates a DOFVector > after refinement virtual void refineInter(DOFVector >*, RCNeighbourList*, int) {} /// Interpolates a DOFVector > after coarsening virtual void coarseInter(DOFVector >*, RCNeighbourList*, int) {} /// Restricts a DOFVector > after coarsening virtual void coarseRestr(DOFVector >*, RCNeighbourList*, int) {} /// Returns local dof indices of the element for the given DOF admin. virtual const DegreeOfFreedom *getLocalIndices(const Element *el, const DOFAdmin *admin, DegreeOfFreedom *dofPtr) const { return NULL; } /// Calculates local dof indices of the element for the given DOF admin. virtual void getLocalIndices(const Element *el, const DOFAdmin *admin, std::vector &indices) const {} virtual void getLocalDofPtrVec(const Element *el, const DOFAdmin *admin, std::vector& vec) const {} /** \brief * Evaluates elements value at barycentric coordinates lambda with local * coefficient vector uh. */ double evalUh(const DimVec& lambda, const double* uh) const; /** \brief * Evaluates elements value at barycentric coordinates lambda with local * coefficient vector uh. If val is not NULL the result will be stored * there, otherwise a pointer to a static local variable is returned which * will be overwritten after the next call. */ const WorldVector& evalUh(const DimVec& lambda, const WorldVector* uh, WorldVector* val) const; /** \brief * Evaluates the gradient at barycentric coordinates lambda. Lambda is the * Jacobian of the barycentric coordinates. uh is the local coefficient * vector. If val is not NULL the result will be stored * there, otherwise a pointer to a static local variable is returned which * will be overwritten after the next call. */ const WorldVector& evalGrdUh(const DimVec& lambda, const DimVec >& Lambda, const double* uh, WorldVector* val) const; /** \brief * Evaluates the second derivative at barycentric coordinates lambda. * Lambda is the Jacobian of the barycentric coordinates. uh is the local * coefficient vector. If val is not NULL the result will be stored * there, otherwise a pointer to a static local variable is returned which * will be overwritten after the next call. */ const WorldMatrix& evalD2Uh(const DimVec& lambda, const DimVec >& Lambda, const double* uh, WorldMatrix* val) const; protected: /// Textual description std::string name; /// Number of basisfunctions on one Element int nBasFcts; /// Maximal degree of the basis functions int degree; /// Dimension of the basis functions int dim; /// Dimension of the world. int dow; /// Number of DOFs at the different positions DimVec *nDOF; /// Vector of the local functions std::vector *phi; /// Vector of gradients std::vector *grdPhi; /// Vector of second derivatives std::vector *d2Phi; /** \brief * Is used by function evalGrdUh. To make it thread safe, we need a * temporary DimVec for every thread. */ std::vector* > grdTmpVec1; /** \brief * Is used by function evalGrdUh. To make it thread safe, we need a * temporary DimVec for every thread. */ std::vector* > grdTmpVec2; }; } #endif // AMDIS_BASISFUNCTION_H