hp_coarsentest.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_HP_COARSENTEST_H
21 #define LIBMESH_HP_COARSENTEST_H
22 
23 // Local Includes
24 #include "libmesh/auto_ptr.h" // deprecated
25 #include "libmesh/dense_matrix.h"
26 #include "libmesh/dense_vector.h"
27 #include "libmesh/hp_selector.h"
28 #include "libmesh/id_types.h"
29 #include "libmesh/libmesh_common.h"
30 #include "libmesh/fe.h" // MipsPro requires fe.h and quadrature.h
31 #include "libmesh/quadrature.h" // Required for inline deletion std::unique_ptrs<> in destructor
32 
33 // C++ includes
34 #include <vector>
35 #include <memory>
36 
37 #ifdef LIBMESH_ENABLE_AMR
38 
39 namespace libMesh
40 {
41 
42 // Forward Declarations
43 class Elem;
44 class Point;
45 class System;
46 template <typename T> class TensorValue;
47 template <typename T> class VectorValue;
52 
53 
68 class HPCoarsenTest : public HPSelector
69 {
70 public:
71 
76  {
77  libmesh_experimental();
78  }
79 
85  HPCoarsenTest (const HPCoarsenTest &) = delete;
86  HPCoarsenTest & operator= (const HPCoarsenTest &) = delete;
87 
91  HPCoarsenTest (HPCoarsenTest &&) = default;
92  HPCoarsenTest & operator= (HPCoarsenTest &&) = default;
93  virtual ~HPCoarsenTest() = default;
94 
101  virtual void select_refinement (System & system) override;
102 
108 
109 protected:
114  void add_projection(const System &, const Elem *, unsigned int var);
115 
119  const Elem * coarse;
120 
124  std::vector<dof_id_type> dof_indices;
125 
129  std::unique_ptr<FEBase> fe, fe_coarse;
130 
134  const std::vector<std::vector<Real>> * phi, * phi_coarse;
135  const std::vector<std::vector<RealGradient>> * dphi, * dphi_coarse;
136  const std::vector<std::vector<RealTensor>> * d2phi, * d2phi_coarse;
137 
141  const std::vector<Real> * JxW;
142 
146  const std::vector<Point> * xyz_values;
147  std::vector<Point> coarse_qpoints;
148 
152  std::unique_ptr<QBase> qrule;
153 
165 };
166 
167 } // namespace libMesh
168 
169 #endif // #ifdef LIBMESH_ENABLE_AMR
170 
171 #endif // LIBMESH_HP_COARSENTEST_H
RealVectorValue RealGradient
const std::vector< std::vector< Real > > * phi
void add_projection(const System &, const Elem *, unsigned int var)
DenseVector< Number > Fe
HPCoarsenTest & operator=(const HPCoarsenTest &)=delete
virtual void select_refinement(System &system) override
VectorValue< Real > RealVectorValue
std::unique_ptr< FEBase > fe_coarse
const std::vector< std::vector< Real > > * phi_coarse
The base class for all geometric element types.
Definition: elem.h:100
const std::vector< std::vector< RealGradient > > * dphi_coarse
RealTensorValue RealTensor
const std::vector< Point > * xyz_values
DenseVector< Number > Uc
const std::vector< Real > * JxW
TensorValue< Real > RealTensorValue
virtual ~HPCoarsenTest()=default
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:92
DenseVector< Number > Up
std::unique_ptr< FEBase > fe
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< std::vector< RealTensor > > * d2phi_coarse
const std::vector< std::vector< RealTensor > > * d2phi
const std::vector< std::vector< RealGradient > > * dphi
std::vector< dof_id_type > dof_indices
std::vector< Point > coarse_qpoints
std::unique_ptr< QBase > qrule
DenseMatrix< Number > Ke