MFEM v2.0
gridfunc.hpp
Go to the documentation of this file.
00001 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
00002 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
00003 // reserved. See file COPYRIGHT for details.
00004 //
00005 // This file is part of the MFEM library. For more information and source code
00006 // availability see http://mfem.googlecode.com.
00007 //
00008 // MFEM is free software; you can redistribute it and/or modify it under the
00009 // terms of the GNU Lesser General Public License (as published by the Free
00010 // Software Foundation) version 2.1 dated February 1999.
00011 
00012 #ifndef MFEM_GRIDFUNC
00013 #define MFEM_GRIDFUNC
00014 
00016 class GridFunction : public Vector
00017 {
00018 protected:
00020    FiniteElementSpace *fes;
00021 
00023    FiniteElementCollection *fec;
00024 
00025    void SaveSTLTri(ostream &out, double p1[], double p2[], double p3[]);
00026 
00027    void GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh);
00028 
00029 public:
00030 
00031    GridFunction() { fes = NULL; fec = NULL; }
00032 
00034    GridFunction(FiniteElementSpace *f) : Vector(f->GetVSize())
00035    { fes = f; fec = NULL; }
00036 
00037    GridFunction(Mesh *m, istream &input);
00038 
00039    GridFunction(Mesh *m, GridFunction *gf_array[], int num_pieces);
00040 
00042    void MakeOwner(FiniteElementCollection *_fec)
00043    { fec = _fec; }
00044 
00045    FiniteElementCollection *OwnFEC() { return fec; }
00046 
00047    int VectorDim() const;
00048 
00050    void GetNodalValues(int i, Array<double> &nval, int vdim = 1) const;
00051 
00052    double GetValue(int i, const IntegrationPoint &ip, int vdim = 1) const;
00053 
00054    void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const;
00055 
00056    void GetValues(int i, const IntegrationRule &ir, Vector &vals,
00057                   DenseMatrix &tr, int vdim = 1) const;
00058 
00059    int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals,
00060                      DenseMatrix &tr, int vdim = 1) const;
00061 
00062    void GetVectorValues(int i, const IntegrationRule &ir,
00063                         DenseMatrix &vals, DenseMatrix &tr) const;
00064 
00065    int GetFaceVectorValues(int i, int side, const IntegrationRule &ir,
00066                            DenseMatrix &vals, DenseMatrix &tr) const;
00067 
00068    void GetValuesFrom(GridFunction &);
00069 
00070    void GetBdrValuesFrom(GridFunction &);
00071 
00072    void GetVectorFieldValues(int i, const IntegrationRule &ir,
00073                              DenseMatrix &vals,
00074                              DenseMatrix &tr, int comp = 0) const;
00075 
00077    void ReorderByNodes();
00078 
00080    void GetNodalValues(Vector &nval, int vdim = 1) const;
00081 
00082    void GetVectorFieldNodalValues(Vector &val, int comp) const;
00083 
00084    void ProjectVectorFieldOn(GridFunction &vec_field, int comp = 0);
00085 
00086    void GetDerivative(int comp, int der_comp, GridFunction &der);
00087 
00088    double GetDivergence(ElementTransformation &tr);
00089 
00090    void GetGradient(ElementTransformation &tr, Vector &grad);
00091 
00092    void GetGradients(const int elem, const IntegrationRule &ir,
00093                      DenseMatrix &grad);
00094 
00095    void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad);
00096 
00100    void GetElementAverages(GridFunction &avgs);
00101 
00102    void ProjectCoefficient(Coefficient &coeff);
00103 
00104    // call fes -> BuildDofToArrays() before using this projection
00105    void ProjectCoefficient(Coefficient &coeff, Array<int> &dofs, int vd = 0);
00106 
00107    void ProjectCoefficient(VectorCoefficient &vcoeff);
00108 
00109    void ProjectCoefficient(Coefficient *coeff[]);
00110 
00111    void ProjectBdrCoefficient(Coefficient &coeff, Array<int> &attr)
00112    {
00113       Coefficient *coeff_p = &coeff;
00114       ProjectBdrCoefficient(&coeff_p, attr);
00115    }
00116 
00117    void ProjectBdrCoefficient(Coefficient *coeff[], Array<int> &attr);
00118 
00122    void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff,
00123                                     Array<int> &bdr_attr);
00124 
00128    void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff,
00129                                      Array<int> &bdr_attr);
00130 
00131    double ComputeL2Error(Coefficient &exsol,
00132                          const IntegrationRule *irs[] = NULL) const
00133    {
00134       Coefficient *exsol_p = &exsol;
00135       return ComputeL2Error(&exsol_p, irs);
00136    }
00137 
00138    double ComputeL2Error(Coefficient *exsol[],
00139                          const IntegrationRule *irs[] = NULL) const;
00140 
00141    double ComputeL2Error(VectorCoefficient &exsol,
00142                          const IntegrationRule *irs[] = NULL,
00143                          Array<int> *elems = NULL) const;
00144 
00145    double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
00146                          Coefficient *ell_coef, double Nu,
00147                          int norm_type) const;
00148 
00149    double ComputeMaxError(Coefficient *exsol[],
00150                           const IntegrationRule *irs[] = NULL) const;
00151 
00152    double ComputeMaxError(VectorCoefficient &exsol,
00153                           const IntegrationRule *irs[] = NULL) const;
00154 
00155    double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad,
00156                           int norm_type, Array<int> *elems = NULL,
00157                           const IntegrationRule *irs[] = NULL) const;
00158 
00159    double ComputeL1Error(VectorCoefficient &exsol,
00160                          const IntegrationRule *irs[] = NULL) const;
00161 
00163    GridFunction &operator=(double value);
00164 
00165    GridFunction &operator=(const Vector &v);
00166 
00167    GridFunction &operator=(const GridFunction &v);
00168 
00169    FiniteElementSpace *FESpace() { return fes; }
00170 
00171    void Update() { SetSize(fes->GetVSize()); }
00172 
00173    void Update(FiniteElementSpace *f);
00174 
00175    void Update(FiniteElementSpace *f, Vector &v, int v_offset);
00176 
00178    virtual void Save(ostream &out);
00179 
00183    void SaveVTK(ostream &out, const string &field_name, int ref);
00184 
00185    void SaveSTL(ostream &out, int TimesToRefine = 1);
00186 
00188    virtual ~GridFunction();
00189 };
00190 
00191 void ComputeFlux(BilinearFormIntegrator &blfi,
00192                  GridFunction &u,
00193                  GridFunction &flux, int wcoef = 1, int sd = -1);
00194 
00195 void ZZErrorEstimator(BilinearFormIntegrator &blfi,
00196                       GridFunction &u,
00197                       GridFunction &flux, Vector &ErrorEstimates,
00198                       int wsd = 1);
00199 
00200 #endif
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines