DFT-EFE
 
Loading...
Searching...
No Matches
FECellDealii.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 Bikash Kanungo, Vishal Subramanian
24 */
25
26#ifndef dftefeFECellDealii_h
27#define dftefeFECellDealii_h
28
29#include <utils/TypeConfig.h>
30#include <deal.II/dofs/dof_accessor.h>
31#include <utils/Point.h>
32#include "FECellBase.h"
33namespace dftefe
34{
35 namespace basis
36 {
37 template <size_type dim>
38 class FECellDealii : public FECellBase
39 {
41 typename dealii::DoFHandler<dim>::active_cell_iterator;
42
43 public:
44 FECellDealii(DealiiFECellIterator dealiiFECellIter);
45
46 void
47 getVertices(std::vector<utils::Point> &points) const override;
48
49 void
50 getVertex(size_type i, utils::Point &point) const override;
51
52 std::vector<std::shared_ptr<utils::Point>>
53 getNodalPoints() const override;
54
56 getId() const override;
57
58 bool
59 isPointInside(const utils::Point &point) const override;
60
61 bool
62 isAtBoundary(const unsigned int i) const override;
63
64 bool
65 isAtBoundary() const override;
66
67
68 double
69 diameter() const override;
70
71 void
72 center(dftefe::utils::Point &centerPoint) const override;
73 void
74 setRefineFlag() override;
75
76 void
77 clearRefineFlag() override;
78
79 double
80 minimumVertexDistance() const override;
81
82 double
83 distanceToUnitCell(dftefe::utils::Point &parametricPoint) const override;
84
85 void
86 setCoarsenFlag() override;
87
88 void
89 clearCoarsenFlag() override;
90
91 bool
92 isActive() const override;
93
94 bool
95 isLocallyOwned() const override;
96
97 bool
98 isGhost() const override;
99
100 bool
101 isArtificial() const override;
102
104 getDim() const override;
105
106 void
107 getParametricPoint(const utils::Point & realPoint,
108 const CellMappingBase &cellMapping,
109 utils::Point & parametricPoint) const override;
110
111 void
112 getRealPoint(const utils::Point & parametricPoint,
113 const CellMappingBase &cellMapping,
114 utils::Point & realPoint) const override;
115
116 void
118 std::vector<global_size_type> &vecId) const override;
119
121 getFaceBoundaryId(size_type faceId) const override;
122
123 void
125 size_type faceId,
126 std::vector<global_size_type> &vecNodeId) const override;
127
128
130 getFEOrder() const override;
131
134
135 private:
137
138
139 }; // end of class FECellDealii
140 } // namespace basis
141} // namespace dftefe
142#include "FECellDealii.t.cpp"
143#endif // dftefeFECellDealii_h
An abstract class to map a real point to parametric point and vice-versa.
Definition: CellMappingBase.h:27
An abstract class for a finite element cell (can be of any dimension) This is created primarily to be...
Definition: FECellBase.h:22
Definition: FECellDealii.h:39
size_type getId() const override
Definition: FECellDealii.t.cpp:80
void setRefineFlag() override
Definition: FECellDealii.t.cpp:128
void center(dftefe::utils::Point &centerPoint) const override
Definition: FECellDealii.t.cpp:118
bool isArtificial() const override
Definition: FECellDealii.t.cpp:194
void getVertex(size_type i, utils::Point &point) const override
Definition: FECellDealii.t.cpp:62
double minimumVertexDistance() const override
Definition: FECellDealii.t.cpp:142
void getFaceDoFGlobalIndices(size_type faceId, std::vector< global_size_type > &vecNodeId) const override
Definition: FECellDealii.t.cpp:246
size_type getFEOrder() const override
Definition: FECellDealii.t.cpp:255
bool isLocallyOwned() const override
Definition: FECellDealii.t.cpp:180
void getParametricPoint(const utils::Point &realPoint, const CellMappingBase &cellMapping, utils::Point &parametricPoint) const override
Definition: FECellDealii.t.cpp:208
void setCoarsenFlag() override
Definition: FECellDealii.t.cpp:159
void clearRefineFlag() override
Definition: FECellDealii.t.cpp:135
DealiiFECellIterator d_dealiiFECellIter
Definition: FECellDealii.h:136
bool isGhost() const override
Definition: FECellDealii.t.cpp:187
size_type getFaceBoundaryId(size_type faceId) const override
Definition: FECellDealii.t.cpp:239
bool isActive() const override
Definition: FECellDealii.t.cpp:173
void clearCoarsenFlag() override
Definition: FECellDealii.t.cpp:166
void getRealPoint(const utils::Point &parametricPoint, const CellMappingBase &cellMapping, utils::Point &realPoint) const override
Definition: FECellDealii.t.cpp:221
DealiiFECellIterator & getDealiiFECellIter()
Definition: FECellDealii.t.cpp:262
double diameter() const override
Definition: FECellDealii.t.cpp:111
bool isAtBoundary() const override
Definition: FECellDealii.t.cpp:105
double distanceToUnitCell(dftefe::utils::Point &parametricPoint) const override
Definition: FECellDealii.t.cpp:149
std::vector< std::shared_ptr< utils::Point > > getNodalPoints() const override
Definition: FECellDealii.t.cpp:69
void cellNodeIdtoGlobalNodeId(std::vector< global_size_type > &vecId) const override
Definition: FECellDealii.t.cpp:231
void getVertices(std::vector< utils::Point > &points) const override
Definition: FECellDealii.t.cpp:45
bool isPointInside(const utils::Point &point) const override
Definition: FECellDealii.t.cpp:89
typename dealii::DoFHandler< dim >::active_cell_iterator DealiiFECellIterator
Definition: FECellDealii.h:41
size_type getDim() const override
Definition: FECellDealii.t.cpp:201
Definition: PointImpl.h:13
dealii includes
Definition: AtomFieldDataSpherical.cpp:31
unsigned int size_type
Definition: TypeConfig.h:8