Home | Doxygen Documentation | Tutorials | Developer Tools (restricted)

hp1D/formula.hh
Go to the documentation of this file.
00001 /* helper classes for bilinear forms and linear forms in hp1d
00002  */
00003 
00004 #ifndef formulahp1d_hh
00005 #define formulahp1d_hh
00006 
00007 #include "basics/operations.hh"
00008 #include "formula/formula.hh"
00009 #include "geometry/formula.hh"
00010 #include "geometry/arrays.hh"
00011 #include "function/vector.hh"
00012 #include "element.hh"
00013 
00014 namespace hp1D {
00015 
00016   // *********************************************** ElementFormulaInterpGrad **
00017 
00018   template<typename F>
00019   class ElementFormulaInterpGrad : public concepts::ElementFormula<F> {
00020   public:
00021     ElementFormulaInterpGrad(const concepts::ElementFormulaContainer<F> frm);
00022 
00023     virtual F operator() (const concepts::ElementWithCell<Real>& elm,
00024                           const Real p, const Real t = 0.0) const;
00025     virtual F operator() (const concepts::ElementWithCell<Real>& elm,
00026                           const concepts::Real2d& p, const Real t = 0.0) const
00027     {
00028       return 0.0;
00029     }
00030     virtual F operator() (const concepts::ElementWithCell<Real>& elm,
00031                           const concepts::Real3d& p, const Real t = 0.0) const
00032     {
00033       return 0.0;
00034     }
00035 
00037     virtual ElementFormulaInterpGrad<F>* clone() const {
00038       return new ElementFormulaInterpGrad<F>(frm_);
00039     }
00040   protected:
00041     virtual std::ostream& info(std::ostream& os) const;
00042   private:
00044     const concepts::ElementFormulaContainer<F> frm_;
00046     mutable concepts::Array<Real> x_;
00048     mutable const concepts::ElementWithCell<Real>* edge_;
00050     mutable concepts::Array<F> values_;
00051   };
00052 
00053   // *************************************************** ArrayElementFormula **
00054 
00056 
00057   template<class F = Real>
00058   class ArrayElementFormula : public concepts::Array<F> {
00059   public:
00064     ArrayElementFormula();
00066     ArrayElementFormula(const Element<Real>& elm,
00067                         const concepts::ElementFormulaContainer<F> frm);
00071     void compute(const Element<Real>& elm,
00072                  const concepts::ElementFormulaContainer<F> frm);
00073   };
00074 
00075   template<class F>
00076   ArrayElementFormula<F>::ArrayElementFormula() : concepts::Array<F>() {}
00077 
00078   template<class F>
00079   ArrayElementFormula<F>::ArrayElementFormula
00080   (const Element<Real>& elm, const concepts::ElementFormulaContainer<F> frm)
00081     : concepts::Array<F>(0) 
00082   {
00083     compute(elm, frm);
00084   }
00085 
00086   // ****************************************************** LinearFormHelper **
00087 
00092   template<uint i, class F = Real>
00093   class LinearFormHelper {
00094   };
00095 
00096   // *************************************************** LinearFormHelper<0> **
00097 
00112   template<class F>
00113   class LinearFormHelper<0,F> {
00114   public:
00115     LinearFormHelper(const concepts::ElementFormulaContainer<F> frm);
00116     virtual ~LinearFormHelper();
00117 
00118   protected:
00123     void computeIntermediate_(const Element<Real>& elm);
00124 
00129     concepts::Array<F> intermediateValue_;
00130 
00132     const concepts::ElementFormulaContainer<F> frm_;
00133   };
00134 
00135   // *************************************************** LinearFormHelper<1> **
00136 
00153   template<class F>
00154   class LinearFormHelper<1,F> {
00155   public:
00156     LinearFormHelper(const concepts::ElementFormulaContainer<F> frm);
00157     virtual ~LinearFormHelper();
00158 
00159   protected:
00164     void computeIntermediate_(const Element<Real>& elm);
00165 
00170     concepts::Array<F> intermediateValue_;
00171 
00173     const concepts::ElementFormulaContainer<F> frm_;
00174   };
00175 
00176 
00177   // **************************************************** BilinearFormHelper **
00178 
00183   template<uint i, uint j, class F = Real>
00184   class BilinearFormHelper {
00185   };
00186 
00187   // *********************************************** BilinearFormHelper<0,0> **
00188 
00189   template<class F>
00190   class BilinearFormHelper<0,0,F> : public LinearFormHelper<0,F> {
00191   public:
00192     BilinearFormHelper(const concepts::ElementFormulaContainer<F> frm)
00193       : LinearFormHelper<0,F>(frm) {}
00194 
00195     const concepts::ElementFormula<F>* formula() const { 
00196       return &this->frm_.frm();
00197     }
00198   };
00199 
00200   // *********************************************** BilinearFormHelper<1,1> **
00201 
00210   template<class F>
00211   class BilinearFormHelper<1,1,F> {
00212   public:
00213     BilinearFormHelper(const concepts::ElementFormulaContainer<F> 
00214                        frm = concepts::ElementFormulaContainer<F>());
00215     virtual ~BilinearFormHelper();
00216   protected:
00218     void computeIntermediate_(const Element<Real>& elm);
00223     concepts::Array<F> intermediateValue_;
00225     const concepts::ElementFormulaContainer<F> frm_;
00226   };
00227 
00228 } // namespace hp1D
00229 
00230 #endif //formulahp1d_hh

Home | Doxygen Documentation | Tutorials | Developer Tools (restricted)