DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
FECell.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 Vishal Subramanian, Bikash Kanungo
24 */
25
26#ifndef DFTFE_FECELL_H
27#define DFTFE_FECELL_H
28
29#include "Cell.h"
30
31namespace dftfe
32{
33 namespace utils
34 {
35 /**
36 * @brief This class provides the specialisation for Finite element cell.
37 * This class provides the function defintions for shape function
38 * calculation at an arbitrary point inside the cell, checks if a point lies
39 * within a cell by computing its parametric points
40 *
41 * @author Vishal Subramanian, Bikash Kanungo
42 */
43 template <size_type dim>
44 class FECell : public Cell<dim>
45 {
46 public:
48 typename dealii::DoFHandler<dim>::active_cell_iterator;
49
51 typename dealii::DoFHandler<dim>::active_cell_iterator dealiiFECellIter,
52 const dealii::FiniteElement<dim, dim> & fe);
53
54 void
55 reinit(DealiiFECellIterator dealiiFECellIter);
56
57 void
58 getVertices(std::vector<std::vector<double>> &points) const override;
59
60 void
61 getVertex(size_type i, std::vector<double> &point) const override;
62
63 std::pair<std::vector<double>, std::vector<double>>
64 getBoundingBox() const override;
65
66 bool
67 isPointInside(const std::vector<double> &point,
68 const double tol) const override;
69
70 std::vector<double>
71 getParametricPoint(const std::vector<double> &realPoint) const override;
72
73
74 std::vector<double>
76 unsigned int numPoints,
77 const std::vector<double> &realPoint) const;
78
81
82 void
84 unsigned int numPointsInCell,
85 const std::vector<double> & parametricPoints,
86 std::vector<dataTypes::number> &shapeFuncValues,
87 unsigned int cellShapeFuncStartIndex,
88 unsigned int numDofsPerElement) const;
89
90 void
91 getShapeFuncValues(unsigned int numPointsInCell,
92 const std::vector<double> &coordinatesOfPointsInCell,
93 std::vector<dataTypes::number> &shapeFuncValues,
94 unsigned int cellShapeFuncStartIndex,
95 unsigned int numDofsPerElement) const override;
96
97 private:
98 std::vector<double> d_lowerLeft;
99 std::vector<double> d_upperRight;
100
102
103 dealii::MappingQ1<dim, dim> d_mappingQ1;
104
105 const dealii::FiniteElement<dim, dim> &d_feCell;
106
107 }; // end of class FECell
108 } // end of namespace utils
109} // end of namespace dftfe
110
111#include "../utils/FECell.t.cc"
112
113//#include "../utils/FECell.t.cc"
114#endif // DFTFE_FECELL_H
DealiiFECellIterator & getDealiiFECellIter()
void getVertex(size_type i, std::vector< double > &point) const override
DealiiFECellIterator d_dealiiFECellIter
Definition FECell.h:101
void reinit(DealiiFECellIterator dealiiFECellIter)
std::pair< std::vector< double >, std::vector< double > > getBoundingBox() const override
dealii::MappingQ1< dim, dim > d_mappingQ1
Definition FECell.h:103
void getShapeFuncValuesFromParametricPoints(unsigned int numPointsInCell, const std::vector< double > &parametricPoints, std::vector< dataTypes::number > &shapeFuncValues, unsigned int cellShapeFuncStartIndex, unsigned int numDofsPerElement) const
std::vector< double > getParametricPointForAllPoints(unsigned int numPoints, const std::vector< double > &realPoint) const
void getVertices(std::vector< std::vector< double > > &points) const override
const dealii::FiniteElement< dim, dim > & d_feCell
Definition FECell.h:105
FECell(typename dealii::DoFHandler< dim >::active_cell_iterator dealiiFECellIter, const dealii::FiniteElement< dim, dim > &fe)
typename dealii::DoFHandler< dim >::active_cell_iterator DealiiFECellIterator
Definition FECell.h:47
std::vector< double > d_lowerLeft
Definition FECell.h:98
bool isPointInside(const std::vector< double > &point, const double tol) const override
std::vector< double > d_upperRight
Definition FECell.h:99
std::vector< double > getParametricPoint(const std::vector< double > &realPoint) const override
void getShapeFuncValues(unsigned int numPointsInCell, const std::vector< double > &coordinatesOfPointsInCell, std::vector< dataTypes::number > &shapeFuncValues, unsigned int cellShapeFuncStartIndex, unsigned int numDofsPerElement) const override
Definition Cell.h:36
Definition pseudoPotentialToDftfeConverter.cc:34
unsigned int size_type
Definition TypeConfig.h:6