DFT-EFE
 
Loading...
Searching...
No Matches
DealiiFEEvaluationWrapper.h
Go to the documentation of this file.
1/******************************************************************************
2 * Copyright (c) 2021. *
3 * The Regents of the University of Michigan and DFT-EFE developers. *
4 * *
5 * This file is part of the DFT-EFE code. *
6 * *
7 * DFT-EFE is free software: you can redistribute it and/or modify *
8 * it under the terms of the Lesser GNU General Public License as *
9 * published by the Free Software Foundation, either version 3 of *
10 * the License, or (at your option) any later version. *
11 * *
12 * DFT-EFE is distributed in the hope that it will be useful, but *
13 * WITHOUT ANY WARRANTY; without even the implied warranty *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. *
15 * See the Lesser GNU General Public License for more details. *
16 * *
17 * You should have received a copy of the GNU Lesser General Public *
18 * License at the top level of DFT-EFE distribution. If not, see *
19 * <https://www.gnu.org/licenses/>. *
20 ******************************************************************************/
21
22/*
23 * @author DFTFE, Avirup Sircar
24 */
25
26#ifndef FEEvaluationWrapper_h
27#define FEEvaluationWrapper_h
28
29#include <boost/preprocessor.hpp>
30#include <cmath>
31#include <memory>
32
33#include <deal.II/matrix_free/matrix_free.h>
34#include <deal.II/matrix_free/fe_evaluation.h>
35
36namespace dftefe
37{
38 namespace basis
39 {
41 {
42 public:
43 template <typename T>
45 dealii::LinearAlgebra::distributed::Vector<T,
46 dealii::MemorySpace::Host>;
47
48 virtual ~FEEvaluationWrapperBase() = 0;
49
53 virtual unsigned int
55
59 virtual void
60 reinit(const unsigned int macrocell) = 0;
61
65 virtual void
67
71 virtual void
73
77 virtual void
78 evaluate(dealii::EvaluationFlags::EvaluationFlags evaluateFlags) = 0;
79
80 virtual void
82 dealii::VectorizedArray<double> &alpha) = 0;
83
84 // virtual void submitInterpolatedGradientsAndMultiply(
85 // dealii::AlignedVector<
86 // dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &alpha) = 0;
87
88 virtual void
90 dealii::VectorizedArray<double> &alpha) = 0;
91
92 virtual void
94
95 virtual void
97 dealii::AlignedVector<dealii::VectorizedArray<double>> &alpha) = 0;
98
99 virtual void
101 const dealii::VectorizedArray<double> &scaling,
102 dealii::AlignedVector<
103 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &alpha) = 0;
104
105 virtual void
107 dealii::AlignedVector<
108 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &alpha) = 0;
109
110 virtual void
112 dealii::AlignedVector<dealii::VectorizedArray<double>> &alpha) = 0;
113
114
115 virtual dealii::VectorizedArray<double>
117
118 virtual void
119 submitValueAtQuadpoint(const unsigned int iQuadPoint,
120 const dealii::VectorizedArray<double> &value) = 0;
121
122 virtual void
123 alphaTimesQuadValuesPlusYFromSubCell(const unsigned int subCellIndex,
124 const double alpha,
125 double * outputVector) = 0;
126
127 virtual void
128 getQuadGradientsForSubCell(const unsigned int subCellIndex,
129 const double alpha,
130 double * outputVector) = 0;
131
132 virtual void
133 getQuadHessianForSubCell(const unsigned int subCellIndex,
134 const double alpha,
135 double * outputVector) = 0;
136
137
138 virtual void
140 const dealii::VectorizedArray<double> &scaleValues,
141 const bool scaleValuesFlag,
142 const dealii::VectorizedArray<double> &scaleGradients,
143 const bool scaleGradientsFlag) = 0;
144
145 virtual void
146 integrate(dealii::EvaluationFlags::EvaluationFlags evaluateFlags) = 0;
147
148 virtual dealii::Point<3, dealii::VectorizedArray<double>>
149 getQuadraturePoint(const unsigned int iQuadPoint) = 0;
150
151 virtual void
153 dealii::AlignedVector<dealii::VectorizedArray<double>> &tempVec) = 0;
154
155 virtual void
157 };
158
159 template <int FEOrder,
160 unsigned int num_1d_quadPoints,
161 unsigned int n_components>
163 {
164 public:
166 const dealii::MatrixFree<3, double> &matrixFreeData,
167 const unsigned int matrixFreeVectorComponent,
168 const unsigned int matrixFreeQuadratureComponent);
169
171
172 std::unique_ptr<
173 dealii::FEEvaluation<3, FEOrder, num_1d_quadPoints, n_components>>
175
176
177 unsigned int
179
180 void
181 reinit(const unsigned int macrocell) override;
182
183 void
184 readDoFValues(const distributedCPUVec<double> &tempvec) override;
185
186 void
187 readDoFValuesPlain(const distributedCPUVec<double> &tempvec) override;
188
189 void
190 evaluate(dealii::EvaluationFlags::EvaluationFlags evaluateFlags) override;
191
192 void
194 dealii::VectorizedArray<double> &alpha) override;
195
196 // void submitInterpolatedGradientsAndMultiply(
197 // dealii::AlignedVector<
198 // dealii::Tensor<1, 3, dealii::VectorizedArray<double>>> &alpha)
199 // override;
200 void
202 dealii::VectorizedArray<double> &alpha) override;
203
204 void
206
207 void
209 dealii::AlignedVector<dealii::VectorizedArray<double>> &alpha) override;
210
211 void
212 submitValues(const dealii::VectorizedArray<double> &scaling,
213 dealii::AlignedVector<
214 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>>
215 &alpha) override;
216
217 void
218 submitGradients(dealii::AlignedVector<
219 dealii::Tensor<1, 3, dealii::VectorizedArray<double>>>
220 &alpha) override;
221
222 void
224 dealii::AlignedVector<dealii::VectorizedArray<double>> &alpha) override;
225
226 dealii::VectorizedArray<double>
227 integrateValue() override;
228
229 void
231 const unsigned int iQuadPoint,
232 const dealii::VectorizedArray<double> &value) override;
233
234 dealii::Point<3, dealii::VectorizedArray<double>>
235 getQuadraturePoint(const unsigned int iQuadPoint) override;
236
237 void
238 alphaTimesQuadValuesPlusYFromSubCell(const unsigned int subCellIndex,
239 const double alpha,
240 double *outputVector) override;
241
242 void
243 getQuadGradientsForSubCell(const unsigned int subCellIndex,
244 const double alpha,
245 double * outputVector) override;
246
247 void
248 getQuadHessianForSubCell(const unsigned int subCellIndex,
249 const double alpha,
250 double * outputVector) override;
251
252 void
254 const dealii::VectorizedArray<double> &scaleValues,
255 const bool scaleValuesFlag,
256 const dealii::VectorizedArray<double> &scaleGradients,
257 const bool scaleGradientsFlag) override;
258
259 void
260 integrate(
261 dealii::EvaluationFlags::EvaluationFlags evaluateFlags) override;
262
263 void
264 getValues(dealii::AlignedVector<dealii::VectorizedArray<double>> &tempVec)
265 override;
266
267 void
269 };
270
271
272
273 template <unsigned int numberOfComponents>
275 {
276 public:
288 unsigned int fe_degree,
289 unsigned int num_1d_quad,
290 const dealii::MatrixFree<3, double> &matrixFreeData,
291 const unsigned int matrixFreeVectorComponent,
292 const unsigned int matrixFreeQuadratureComponent);
293
295
298
299 private:
300 unsigned int d_feDegree;
301 unsigned int d_num1dQuad;
304
305 std::unique_ptr<FEEvaluationWrapperBase> d_feEvaluationBase;
306 const dealii::MatrixFree<3, double> * d_matrix_free_data;
307
308
309
310 }; // end of DealiiFEEvaluationWrapper
311
312 } // end of namespace basis
313
314} // end of namespace dftefe
315
316
317#endif // FEEvaluationWrapper_h
Definition: DealiiFEEvaluationWrapper.h:275
~DealiiFEEvaluationWrapper()
Definition: DealiiFEEvaluationWrapper.cpp:526
const dealii::MatrixFree< 3, double > * d_matrix_free_data
Definition: DealiiFEEvaluationWrapper.h:306
unsigned int d_matrixFreeVectorComponent
Definition: DealiiFEEvaluationWrapper.h:302
FEEvaluationWrapperBase & getFEEvaluationWrapperBase() const
Definition: DealiiFEEvaluationWrapper.cpp:536
unsigned int d_matrixFreeQuadratureComponent
Definition: DealiiFEEvaluationWrapper.h:303
std::unique_ptr< FEEvaluationWrapperBase > d_feEvaluationBase
Definition: DealiiFEEvaluationWrapper.h:305
unsigned int d_num1dQuad
Definition: DealiiFEEvaluationWrapper.h:301
unsigned int d_feDegree
Definition: DealiiFEEvaluationWrapper.h:300
Definition: DealiiFEEvaluationWrapper.h:41
virtual dealii::Point< 3, dealii::VectorizedArray< double > > getQuadraturePoint(const unsigned int iQuadPoint)=0
virtual void submitValueAtQuadpoint(const unsigned int iQuadPoint, const dealii::VectorizedArray< double > &value)=0
virtual void readDoFValuesPlain(const distributedCPUVec< double > &tempvec)=0
Calls dealii::FEEvaluation::read_dofs_values_plain.
virtual void alphaTimesQuadValuesPlusYFromSubCell(const unsigned int subCellIndex, const double alpha, double *outputVector)=0
virtual void getQuadHessianForSubCell(const unsigned int subCellIndex, const double alpha, double *outputVector)=0
virtual void submitInterpolatedValuesAndMultiplySquared()=0
virtual unsigned int totalNumberofQuadraturePoints()=0
Returns the total number of quadrature points in all 3 directions.
dealii::LinearAlgebra::distributed::Vector< T, dealii::MemorySpace::Host > distributedCPUVec
Definition: DealiiFEEvaluationWrapper.h:46
virtual void submitValues(const dealii::VectorizedArray< double > &scaling, dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &alpha)=0
virtual void readDoFValues(const distributedCPUVec< double > &tempvec)=0
Calls dealii::FEEvaluation::read_dof_values.
virtual void integrate(dealii::EvaluationFlags::EvaluationFlags evaluateFlags)=0
virtual void getValues(dealii::AlignedVector< dealii::VectorizedArray< double > > &tempVec)=0
virtual dealii::VectorizedArray< double > integrateValue()=0
virtual void submitInterpolatedValuesSubmitInterpolatedGradients(const dealii::VectorizedArray< double > &scaleValues, const bool scaleValuesFlag, const dealii::VectorizedArray< double > &scaleGradients, const bool scaleGradientsFlag)=0
virtual void reinit(const unsigned int macrocell)=0
reinits the dealii::FEEvaluation object for the macrocellIndex
virtual void submitInterpolatedValuesAndMultiply(dealii::VectorizedArray< double > &alpha)=0
virtual void submitInterpolatedValuesAndMultiply(dealii::AlignedVector< dealii::VectorizedArray< double > > &alpha)=0
virtual void getQuadGradientsForSubCell(const unsigned int subCellIndex, const double alpha, double *outputVector)=0
virtual void distributeLocalToGlobal(distributedCPUVec< double > &tempvec)=0
virtual void submitValues(dealii::AlignedVector< dealii::VectorizedArray< double > > &alpha)=0
virtual void submitGradients(dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &alpha)=0
virtual ~FEEvaluationWrapperBase()=0
Definition: DealiiFEEvaluationWrapper.cpp:57
virtual void submitInterpolatedGradientsAndMultiply(dealii::VectorizedArray< double > &alpha)=0
virtual void evaluate(dealii::EvaluationFlags::EvaluationFlags evaluateFlags)=0
Calls the dealii::FEEvaluation::evaluate.
Definition: DealiiFEEvaluationWrapper.h:163
void submitGradients(dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &alpha) override
Definition: DealiiFEEvaluationWrapper.cpp:257
void submitInterpolatedValuesAndMultiply(dealii::VectorizedArray< double > &alpha) override
Definition: DealiiFEEvaluationWrapper.cpp:128
void distributeLocalToGlobal(distributedCPUVec< double > &tempvec) override
Definition: DealiiFEEvaluationWrapper.cpp:451
void submitInterpolatedGradientsAndMultiply(dealii::VectorizedArray< double > &alpha) override
Definition: DealiiFEEvaluationWrapper.cpp:113
void submitInterpolatedValuesSubmitInterpolatedGradients(const dealii::VectorizedArray< double > &scaleValues, const bool scaleValuesFlag, const dealii::VectorizedArray< double > &scaleGradients, const bool scaleGradientsFlag) override
Definition: DealiiFEEvaluationWrapper.cpp:193
void submitValues(const dealii::VectorizedArray< double > &scaling, dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &alpha) override
Definition: DealiiFEEvaluationWrapper.cpp:176
std::unique_ptr< dealii::FEEvaluation< 3, FEOrder, num_1d_quadPoints, n_components > > d_dealiiFEEvaluation
Definition: DealiiFEEvaluationWrapper.h:174
dealii::Point< 3, dealii::VectorizedArray< double > > getQuadraturePoint(const unsigned int iQuadPoint) override
Definition: DealiiFEEvaluationWrapper.cpp:441
~FEEvaluationWrapperDerived()
Definition: DealiiFEEvaluationWrapper.cpp:52
void getQuadHessianForSubCell(const unsigned int subCellIndex, const double alpha, double *outputVector) override
Definition: DealiiFEEvaluationWrapper.cpp:328
void submitInterpolatedValuesAndMultiplySquared() override
Definition: DealiiFEEvaluationWrapper.cpp:143
void readDoFValuesPlain(const distributedCPUVec< double > &tempvec) override
Calls dealii::FEEvaluation::read_dofs_values_plain.
Definition: DealiiFEEvaluationWrapper.cpp:93
void reinit(const unsigned int macrocell) override
reinits the dealii::FEEvaluation object for the macrocellIndex
Definition: DealiiFEEvaluationWrapper.cpp:75
unsigned int totalNumberofQuadraturePoints() override
Returns the total number of quadrature points in all 3 directions.
Definition: DealiiFEEvaluationWrapper.cpp:65
void getValues(dealii::AlignedVector< dealii::VectorizedArray< double > > &tempVec) override
Definition: DealiiFEEvaluationWrapper.cpp:400
dealii::VectorizedArray< double > integrateValue() override
Definition: DealiiFEEvaluationWrapper.cpp:272
void alphaTimesQuadValuesPlusYFromSubCell(const unsigned int subCellIndex, const double alpha, double *outputVector) override
Definition: DealiiFEEvaluationWrapper.cpp:282
void evaluate(dealii::EvaluationFlags::EvaluationFlags evaluateFlags) override
Calls the dealii::FEEvaluation::evaluate.
Definition: DealiiFEEvaluationWrapper.cpp:102
void submitValueAtQuadpoint(const unsigned int iQuadPoint, const dealii::VectorizedArray< double > &value) override
Definition: DealiiFEEvaluationWrapper.cpp:429
void integrate(dealii::EvaluationFlags::EvaluationFlags evaluateFlags) override
Definition: DealiiFEEvaluationWrapper.cpp:419
void getQuadGradientsForSubCell(const unsigned int subCellIndex, const double alpha, double *outputVector) override
Definition: DealiiFEEvaluationWrapper.cpp:301
void readDoFValues(const distributedCPUVec< double > &tempvec) override
Calls dealii::FEEvaluation::read_dof_values.
Definition: DealiiFEEvaluationWrapper.cpp:84
dealii includes
Definition: AtomFieldDataSpherical.cpp:31