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 virtual bool
68 hasPeriodicNeighbor(const unsigned int i) const override;
69
70 double
71 diameter() const override;
72
73 void
74 center(dftefe::utils::Point &centerPoint) const override;
75 void
76 setRefineFlag() override;
77
78 void
79 clearRefineFlag() override;
80
81 double
82 minimumVertexDistance() const override;
83
84 double
85 distanceToUnitCell(dftefe::utils::Point &parametricPoint) const override;
86
87 void
88 setCoarsenFlag() override;
89
90 void
91 clearCoarsenFlag() override;
92
93 bool
94 isActive() const override;
95
96 bool
97 isLocallyOwned() const override;
98
99 bool
100 isGhost() const override;
101
102 bool
103 isArtificial() const override;
104
106 getDim() const override;
107
108 void
109 getParametricPoint(const utils::Point & realPoint,
110 const CellMappingBase &cellMapping,
111 utils::Point & parametricPoint) const override;
112
113 void
114 getRealPoint(const utils::Point & parametricPoint,
115 const CellMappingBase &cellMapping,
116 utils::Point & realPoint) const override;
117
118 void
120 std::vector<global_size_type> &vecId) const override;
121
123 getFaceBoundaryId(size_type faceId) const override;
124
125 void
127 size_type faceId,
128 std::vector<global_size_type> &vecNodeId) const override;
129
130
132 getFEOrder() const override;
133
136
137 private:
139
140
141 }; // end of class FECellDealii
142 } // namespace basis
143} // namespace dftefe
144#include "FECellDealii.t.cpp"
145#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:135
void center(dftefe::utils::Point &centerPoint) const override
Definition: FECellDealii.t.cpp:125
bool isArtificial() const override
Definition: FECellDealii.t.cpp:201
void getVertex(size_type i, utils::Point &point) const override
Definition: FECellDealii.t.cpp:62
double minimumVertexDistance() const override
Definition: FECellDealii.t.cpp:149
void getFaceDoFGlobalIndices(size_type faceId, std::vector< global_size_type > &vecNodeId) const override
Definition: FECellDealii.t.cpp:253
size_type getFEOrder() const override
Definition: FECellDealii.t.cpp:262
bool isLocallyOwned() const override
Definition: FECellDealii.t.cpp:187
void getParametricPoint(const utils::Point &realPoint, const CellMappingBase &cellMapping, utils::Point &parametricPoint) const override
Definition: FECellDealii.t.cpp:215
void setCoarsenFlag() override
Definition: FECellDealii.t.cpp:166
void clearRefineFlag() override
Definition: FECellDealii.t.cpp:142
DealiiFECellIterator d_dealiiFECellIter
Definition: FECellDealii.h:138
bool isGhost() const override
Definition: FECellDealii.t.cpp:194
size_type getFaceBoundaryId(size_type faceId) const override
Definition: FECellDealii.t.cpp:246
bool isActive() const override
Definition: FECellDealii.t.cpp:180
void clearCoarsenFlag() override
Definition: FECellDealii.t.cpp:173
void getRealPoint(const utils::Point &parametricPoint, const CellMappingBase &cellMapping, utils::Point &realPoint) const override
Definition: FECellDealii.t.cpp:228
DealiiFECellIterator & getDealiiFECellIter()
Definition: FECellDealii.t.cpp:269
virtual bool hasPeriodicNeighbor(const unsigned int i) const override
Definition: FECellDealii.t.cpp:105
double diameter() const override
Definition: FECellDealii.t.cpp:118
bool isAtBoundary() const override
Definition: FECellDealii.t.cpp:112
double distanceToUnitCell(dftefe::utils::Point &parametricPoint) const override
Definition: FECellDealii.t.cpp:156
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:238
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:208
Definition: PointImpl.h:13
dealii includes
Definition: AtomFieldDataSpherical.cpp:31
unsigned int size_type
Definition: TypeConfig.h:8