DFT-EFE
 
Loading...
Searching...
No Matches
CFEBasisDofHandlerDealii.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, Avirup Sircar
24 */
25
26#ifndef dftefeCFEBasisDofHandlerDealii_h
27#define dftefeCFEBasisDofHandlerDealii_h
28
29#include <utils/TypeConfig.h>
30#include <utils/Point.h>
32#include <memory>
33#include <deal.II/fe/fe_q.h>
34
37#include <memory>
38#include <map>
39
41#include <deal.II/dofs/dof_handler.h>
42#include <deal.II/matrix_free/matrix_free.h>
43
44namespace dftefe
45{
46 namespace basis
47 {
52 template <typename ValueTypeBasisCoeff,
53 utils::MemorySpace memorySpace,
54 size_type dim>
56 : public FEBasisDofHandler<ValueTypeBasisCoeff, memorySpace, dim>
57 {
58 public:
59 using FECellIterator = typename FEBasisDofHandler<ValueTypeBasisCoeff,
60 memorySpace,
65
67 std::shared_ptr<const TriangulationBase> triangulation,
68 const size_type feOrder,
69 const utils::mpi::MPIComm & mpiComm);
70
72 std::shared_ptr<const TriangulationBase> triangulation,
73 const size_type feOrder);
74
76
77 double
78 getBasisFunctionValue(const size_type basisId,
79 const utils::Point &point) const override;
80 std::vector<double>
82 const size_type basisId,
83 const utils::Point &point,
84 const size_type derivativeOrder = 1) const override;
85
87 void
88 reinit(std::shared_ptr<const TriangulationBase> triangulation,
89 const size_type feOrder,
90 const utils::mpi::MPIComm & mpiComm);
91
92 void
93 reinit(std::shared_ptr<const TriangulationBase> triangulation,
94 const size_type feOrder);
95
96 std::shared_ptr<const TriangulationBase>
97 getTriangulation() const override;
98
100 nLocalCells() const override;
102 nLocallyOwnedCells() const override;
103
105 nGlobalCells() const override;
107 getFEOrder(size_type cellId) const override;
108
110 nCellDofs(size_type cellId) const override;
111
112 bool
113 isVariableDofsPerCell() const override;
114
116 nLocalNodes() const override;
117
118 std::vector<std::pair<global_size_type, global_size_type>>
119 getLocallyOwnedRanges() const override;
120
121 std::vector<std::pair<global_size_type, global_size_type>>
122 getGlobalRanges() const override;
123
124 std::map<BasisIdAttribute, size_type>
125 getBasisAttributeToRangeIdMap() const override;
126
128 nGlobalNodes() const override;
129
130 std::vector<size_type>
131 getLocalNodeIds(size_type cellId) const override;
132
133 std::vector<size_type>
134 getGlobalNodeIds() const override;
135
136 void
138 size_type cellId,
139 std::vector<global_size_type> &vecGlobalNodeId) const override;
140
141 const std::vector<global_size_type> &
142 getBoundaryIds() const override;
143
145 beginLocallyOwnedCells() override;
146
148 endLocallyOwnedCells() override;
149
151 beginLocallyOwnedCells() const override;
152
154 endLocallyOwnedCells() const override;
155
157 beginLocalCells() override;
159 endLocalCells() override;
161 beginLocalCells() const override;
163 endLocalCells() const override;
164 unsigned int
165 getDim() const override;
166
168 nCumulativeLocallyOwnedCellDofs() const override;
169
171 nCumulativeLocalCellDofs() const override;
172
174 totalRanges() const override;
175
176 // This assumes a linear cell mapping
177 void
179 std::map<global_size_type, utils::Point> &dofCoords) const override;
180
181 // Additional functions for getting geometric constriants matrix
182 // Additional functions for getting the communication pattern object
183 // for MPI case
184
185 std::shared_ptr<const ConstraintsLocal<ValueTypeBasisCoeff, memorySpace>>
186 getIntrinsicConstraints() const override;
187
188 // use this to add extra constraints on top of geometric constraints
189 std::shared_ptr<ConstraintsLocal<ValueTypeBasisCoeff, memorySpace>>
190 createConstraintsStart() const override;
191
192 // call this after calling start
193 void
196 constraintsLocal) const override;
197
198 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
199 getMPIPatternP2P() const override;
200
201 bool
202 isDistributed() const override;
203
204
205 //
206 // dealii specific functions
207 //
208 std::shared_ptr<const dealii::DoFHandler<dim>>
209 getDoFHandler() const;
210
211 const dealii::FiniteElement<dim> &
212 getReferenceFE(const size_type cellId) const;
213
214 private:
215 std::shared_ptr<const TriangulationBase> d_triangulation;
216 std::shared_ptr<dealii::DoFHandler<dim>> d_dofHandler;
218 std::vector<std::shared_ptr<FECellBase>> d_localCells;
219 std::vector<std::shared_ptr<FECellBase>> d_locallyOwnedCells;
223
224 std::vector<global_size_type> d_boundaryIds;
226 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
228 std::shared_ptr<const ConstraintsLocal<ValueTypeBasisCoeff, memorySpace>>
230 std::vector<std::pair<global_size_type, global_size_type>>
232
233 void
235 const std::pair<global_size_type, global_size_type> &locallyOwnedRanges,
236 std::vector<std::pair<global_size_type, global_size_type>>
237 & allOwnedRanges,
238 const utils::mpi::MPIComm &mpiComm) const;
239
240 }; // end of CFEBasisDofHandlerDealii
241 } // end of namespace basis
242} // end of namespace dftefe
244#endif // dftefeCFEBasisDofHandlerDealii_h
Definition: CFEBasisDofHandlerDealii.h:57
std::shared_ptr< const dealii::DoFHandler< dim > > getDoFHandler() const
Definition: CFEBasisDofHandlerDealii.t.cpp:1061
FECellIterator beginLocallyOwnedCells() override
Definition: CFEBasisDofHandlerDealii.t.cpp:961
size_type d_numCumulativeLocalCellDofs
Definition: CFEBasisDofHandlerDealii.h:221
const dealii::FiniteElement< dim > & getReferenceFE(const size_type cellId) const
Definition: CFEBasisDofHandlerDealii.t.cpp:1071
std::vector< global_size_type > d_boundaryIds
Definition: CFEBasisDofHandlerDealii.h:224
size_type nCumulativeLocalCellDofs() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1122
const std::vector< global_size_type > & getBoundaryIds() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:950
global_size_type nGlobalNodes() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:811
typename FEBasisDofHandler< ValueTypeBasisCoeff, memorySpace, dim >::const_FECellIterator const_FECellIterator
Definition: CFEBasisDofHandlerDealii.h:64
FECellIterator endLocallyOwnedCells() override
Definition: CFEBasisDofHandlerDealii.t.cpp:972
std::vector< size_type > getGlobalNodeIds() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:921
std::shared_ptr< const ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > getIntrinsicConstraints() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1146
std::shared_ptr< const TriangulationBase > getTriangulation() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:730
std::vector< std::shared_ptr< FECellBase > > d_locallyOwnedCells
Definition: CFEBasisDofHandlerDealii.h:219
void getBasisCenters(std::map< global_size_type, utils::Point > &dofCoords) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1095
size_type totalRanges() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1132
unsigned int getDim() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1047
std::shared_ptr< const ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > d_constraintsLocal
Definition: CFEBasisDofHandlerDealii.h:229
void getCellDofsGlobalIds(size_type cellId, std::vector< global_size_type > &vecGlobalNodeId) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:937
std::vector< double > getBasisFunctionDerivative(const size_type basisId, const utils::Point &point, const size_type derivativeOrder=1) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:713
void reinit(std::shared_ptr< const TriangulationBase > triangulation, const size_type feOrder, const utils::mpi::MPIComm &mpiComm)
Definition: CFEBasisDofHandlerDealii.t.cpp:133
size_type d_totalRanges
Definition: CFEBasisDofHandlerDealii.h:222
bool d_isDistributed
Definition: CFEBasisDofHandlerDealii.h:225
size_type getFEOrder(size_type cellId) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:770
std::shared_ptr< ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > createConstraintsStart() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1156
std::vector< std::pair< global_size_type, global_size_type > > getGlobalRanges() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:879
std::shared_ptr< const TriangulationBase > d_triangulation
Definition: CFEBasisDofHandlerDealii.h:215
std::vector< size_type > getLocalNodeIds(size_type cellId) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:906
std::vector< std::pair< global_size_type, global_size_type > > getLocallyOwnedRanges() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:821
void createConstraintsEnd(std::shared_ptr< ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > constraintsLocal) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1178
size_type nCumulativeLocallyOwnedCellDofs() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1112
FECellIterator beginLocalCells() override
Definition: CFEBasisDofHandlerDealii.t.cpp:1005
size_type nCellDofs(size_type cellId) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:780
double getBasisFunctionValue(const size_type basisId, const utils::Point &point) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:699
std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > getMPIPatternP2P() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1190
std::vector< std::pair< global_size_type, global_size_type > > d_locallyOwnedRanges
Definition: CFEBasisDofHandlerDealii.h:231
size_type nLocalNodes() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:801
size_type nLocalCells() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:740
size_type d_numCumulativeLocallyOwnedCellDofs
Definition: CFEBasisDofHandlerDealii.h:220
FECellIterator endLocalCells() override
Definition: CFEBasisDofHandlerDealii.t.cpp:1016
std::map< BasisIdAttribute, size_type > getBasisAttributeToRangeIdMap() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:894
bool isDistributed() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1200
std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > d_mpiPatternP2P
Definition: CFEBasisDofHandlerDealii.h:227
bool d_isVariableDofsPerCell
Definition: CFEBasisDofHandlerDealii.h:217
size_type nLocallyOwnedCells() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:750
bool isVariableDofsPerCell() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:791
size_type nGlobalCells() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:760
std::vector< std::shared_ptr< FECellBase > > d_localCells
Definition: CFEBasisDofHandlerDealii.h:218
std::shared_ptr< dealii::DoFHandler< dim > > d_dofHandler
Definition: CFEBasisDofHandlerDealii.h:216
typename FEBasisDofHandler< ValueTypeBasisCoeff, memorySpace, dim >::FECellIterator FECellIterator
Definition: CFEBasisDofHandlerDealii.h:61
void getAllOwnedClassicalRanges(const std::pair< global_size_type, global_size_type > &locallyOwnedRanges, std::vector< std::pair< global_size_type, global_size_type > > &allOwnedRanges, const utils::mpi::MPIComm &mpiComm) const
Definition: CFEBasisDofHandlerDealii.t.cpp:831
Definition: ConstraintsLocal.h:43
Definition: FEBasisDofHandler.h:57
Definition: PointImpl.h:13
int MPIComm
Definition: MPITypes.h:83
MemorySpace
Definition: MemorySpaceType.h:37
dealii includes
Definition: AtomFieldDataSpherical.cpp:31
unsigned int size_type
Definition: TypeConfig.h:8
unsigned long int global_size_type
Definition: TypeConfig.h:9