libMesh
hp_coarsentest.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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/dense_matrix.h"
25 #include "libmesh/dense_vector.h"
26 #include "libmesh/hp_selector.h"
27 #include "libmesh/id_types.h"
28 #include "libmesh/libmesh_common.h"
29 #include "libmesh/fe.h" // MipsPro requires fe.h and quadrature.h
30 #include "libmesh/quadrature.h" // Required for inline deletion std::unique_ptrs<> in destructor
31 
32 // C++ includes
33 #include <vector>
34 #include <memory>
35 
36 #ifdef LIBMESH_ENABLE_AMR
37 
38 namespace libMesh
39 {
40 
41 // Forward Declarations
42 class Elem;
43 class Point;
44 class System;
45 template <typename T> class TensorValue;
46 template <typename T> class VectorValue;
51 
52 
67 class HPCoarsenTest : public HPSelector
68 {
69 public:
70 
75  {
76  libmesh_experimental();
77  }
78 
84  HPCoarsenTest (const HPCoarsenTest &) = delete;
85  HPCoarsenTest & operator= (const HPCoarsenTest &) = delete;
86 
90  HPCoarsenTest (HPCoarsenTest &&) = default;
91  HPCoarsenTest & operator= (HPCoarsenTest &&) = default;
92  virtual ~HPCoarsenTest() = default;
93 
100  virtual void select_refinement (System & system) override;
101 
107 
108 protected:
113  void add_projection(const System &, const Elem *, unsigned int var);
114 
119 
123  std::vector<dof_id_type> dof_indices;
124 
128  std::unique_ptr<FEBase> fe, fe_coarse;
129 
133  const std::vector<std::vector<Real>> * phi, * phi_coarse;
134  const std::vector<std::vector<RealGradient>> * dphi, * dphi_coarse;
135  const std::vector<std::vector<RealTensor>> * d2phi, * d2phi_coarse;
136 
140  const std::vector<Real> * JxW;
141 
145  const std::vector<Point> * xyz_values;
146  std::vector<Point> coarse_qpoints;
147 
151  std::unique_ptr<QBase> qrule;
152 
164 };
165 
166 } // namespace libMesh
167 
168 #endif // #ifdef LIBMESH_ENABLE_AMR
169 
170 #endif // LIBMESH_HP_COARSENTEST_H
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::HPCoarsenTest::Fe
DenseVector< Number > Fe
Definition: hp_coarsentest.h:157
libMesh::HPCoarsenTest::Up
DenseVector< Number > Up
Definition: hp_coarsentest.h:163
libMesh::HPCoarsenTest::fe_coarse
std::unique_ptr< FEBase > fe_coarse
Definition: hp_coarsentest.h:128
libMesh::RealGradient
RealVectorValue RealGradient
Definition: hp_coarsentest.h:49
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::RealTensor
RealTensorValue RealTensor
Definition: hp_coarsentest.h:50
libMesh::HPCoarsenTest::JxW
const std::vector< Real > * JxW
Mapping jacobians.
Definition: hp_coarsentest.h:140
libMesh::DenseMatrix< Number >
libMesh::HPCoarsenTest::xyz_values
const std::vector< Point > * xyz_values
Quadrature locations.
Definition: hp_coarsentest.h:145
libMesh::RealVectorValue
VectorValue< Real > RealVectorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
Definition: hp_coarsentest.h:46
libMesh::VectorValue
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
Definition: exact_solution.h:54
libMesh::HPCoarsenTest::qrule
std::unique_ptr< QBase > qrule
The quadrature rule for the fine element.
Definition: hp_coarsentest.h:151
libMesh::HPCoarsenTest::Ke
DenseMatrix< Number > Ke
Linear system for projections.
Definition: hp_coarsentest.h:156
libMesh::HPCoarsenTest::coarse
Elem * coarse
The coarse element on which a solution projection is cached.
Definition: hp_coarsentest.h:118
libMesh::RealTensorValue
TensorValue< Real > RealTensorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
Definition: hp_coarsentest.h:48
libMesh::HPSelector
Subclasses of this abstract base class choose between h refining and p elevation.
Definition: hp_selector.h:48
libMesh::TensorValue
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
Definition: exact_solution.h:53
libMesh::HPCoarsenTest::~HPCoarsenTest
virtual ~HPCoarsenTest()=default
libMesh::HPCoarsenTest::HPCoarsenTest
HPCoarsenTest()
Constructor.
Definition: hp_coarsentest.h:74
libMesh::HPCoarsenTest::Uc
DenseVector< Number > Uc
Coefficients for projected coarse and projected p-derefined solutions.
Definition: hp_coarsentest.h:162
libMesh::HPCoarsenTest::add_projection
void add_projection(const System &, const Elem *, unsigned int var)
The helper function which adds individual fine element data to the coarse element projection.
Definition: hp_coarsentest.C:47
libMesh::HPCoarsenTest::fe
std::unique_ptr< FEBase > fe
The finite element objects for fine and coarse elements.
Definition: hp_coarsentest.h:128
libMesh::HPCoarsenTest::dof_indices
std::vector< dof_id_type > dof_indices
Global DOF indices for fine elements.
Definition: hp_coarsentest.h:123
libMesh::HPCoarsenTest::d2phi
const std::vector< std::vector< RealTensor > > * d2phi
Definition: hp_coarsentest.h:135
libMesh::HPCoarsenTest::p_weight
Real p_weight
Because the coarsening test seems to always choose p refinement, we're providing an option to make h ...
Definition: hp_coarsentest.h:106
libMesh::HPCoarsenTest::phi_coarse
const std::vector< std::vector< Real > > * phi_coarse
Definition: hp_coarsentest.h:133
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::HPCoarsenTest::dphi
const std::vector< std::vector< RealGradient > > * dphi
Definition: hp_coarsentest.h:134
libMesh::HPCoarsenTest::phi
const std::vector< std::vector< Real > > * phi
The shape functions and their derivatives.
Definition: hp_coarsentest.h:133
libMesh::HPCoarsenTest::dphi_coarse
const std::vector< std::vector< RealGradient > > * dphi_coarse
Definition: hp_coarsentest.h:134
libMesh::HPCoarsenTest
This class uses the error estimate given by different types of derefinement (h coarsening or p reduct...
Definition: hp_coarsentest.h:67
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::HPCoarsenTest::select_refinement
virtual void select_refinement(System &system) override
This pure virtual function must be redefined in derived classes to take a mesh flagged for h refineme...
Definition: hp_coarsentest.C:137
libMesh::HPCoarsenTest::operator=
HPCoarsenTest & operator=(const HPCoarsenTest &)=delete
libMesh::HPCoarsenTest::coarse_qpoints
std::vector< Point > coarse_qpoints
Definition: hp_coarsentest.h:146
libMesh::DenseVector< Number >
libMesh::HPCoarsenTest::d2phi_coarse
const std::vector< std::vector< RealTensor > > * d2phi_coarse
Definition: hp_coarsentest.h:135