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
75 double
76 getBasisFunctionValue(const size_type basisId,
77 const utils::Point &point) const override;
78 std::vector<double>
80 const size_type basisId,
81 const utils::Point &point,
82 const size_type derivativeOrder = 1) const override;
83
85 void
86 reinit(std::shared_ptr<const TriangulationBase> triangulation,
87 const size_type feOrder,
88 const utils::mpi::MPIComm & mpiComm);
89
90 void
91 reinit(std::shared_ptr<const TriangulationBase> triangulation,
92 const size_type feOrder);
93
94 std::shared_ptr<const TriangulationBase>
95 getTriangulation() const override;
96
98 nLocalCells() const override;
100 nLocallyOwnedCells() const override;
101
103 nGlobalCells() const override;
105 getFEOrder(size_type cellId) const override;
106
108 nCellDofs(size_type cellId) const override;
109
110 bool
111 isVariableDofsPerCell() const override;
112
114 nLocalNodes() const override;
115
116 std::vector<std::pair<global_size_type, global_size_type>>
117 getLocallyOwnedRanges() const override;
118
119 std::vector<std::pair<global_size_type, global_size_type>>
120 getGlobalRanges() const override;
121
122 std::map<BasisIdAttribute, size_type>
123 getBasisAttributeToRangeIdMap() const override;
124
126 nGlobalNodes() const override;
127
128 std::vector<size_type>
129 getLocalNodeIds(size_type cellId) const override;
130
131 std::vector<size_type>
132 getGlobalNodeIds() const override;
133
134 void
136 size_type cellId,
137 std::vector<global_size_type> &vecGlobalNodeId) const override;
138
139 const std::vector<global_size_type> &
140 getBoundaryIds() const override;
141
143 beginLocallyOwnedCells() override;
144
146 endLocallyOwnedCells() override;
147
149 beginLocallyOwnedCells() const override;
150
152 endLocallyOwnedCells() const override;
153
155 beginLocalCells() override;
157 endLocalCells() override;
159 beginLocalCells() const override;
161 endLocalCells() const override;
162 unsigned int
163 getDim() const override;
164
166 nCumulativeLocallyOwnedCellDofs() const override;
167
169 nCumulativeLocalCellDofs() const override;
170
172 totalRanges() const override;
173
174 // This assumes a linear cell mapping
175 void
177 std::map<global_size_type, utils::Point> &dofCoords) const override;
178
179 // Additional functions for getting geometric constriants matrix
180 // Additional functions for getting the communication pattern object
181 // for MPI case
182
183 std::shared_ptr<const ConstraintsLocal<ValueTypeBasisCoeff, memorySpace>>
184 getIntrinsicConstraints() const override;
185
186 // use this to add extra constraints on top of geometric constraints
187 std::shared_ptr<ConstraintsLocal<ValueTypeBasisCoeff, memorySpace>>
188 createConstraintsStart() const override;
189
190 // call this after calling start
191 void
194 constraintsLocal) const override;
195
196 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
197 getMPIPatternP2P() const override;
198
199 bool
200 isDistributed() const override;
201
202
203 //
204 // dealii specific functions
205 //
206 std::shared_ptr<const dealii::DoFHandler<dim>>
207 getDoFHandler() const;
208
209 const dealii::FiniteElement<dim> &
210 getReferenceFE(const size_type cellId) const;
211
212 private:
213 std::shared_ptr<const TriangulationBase> d_triangulation;
214 std::shared_ptr<dealii::DoFHandler<dim>> d_dofHandler;
216 std::vector<std::shared_ptr<FECellBase>> d_localCells;
217 std::vector<std::shared_ptr<FECellBase>> d_locallyOwnedCells;
221
222 std::vector<global_size_type> d_boundaryIds;
224 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
226 std::shared_ptr<const ConstraintsLocal<ValueTypeBasisCoeff, memorySpace>>
228 std::vector<std::pair<global_size_type, global_size_type>>
230
231 void
233 const std::pair<global_size_type, global_size_type> &locallyOwnedRanges,
234 std::vector<std::pair<global_size_type, global_size_type>>
235 & allOwnedRanges,
236 const utils::mpi::MPIComm &mpiComm) const;
237
238 }; // end of CFEBasisDofHandlerDealii
239 } // end of namespace basis
240} // end of namespace dftefe
242#endif // dftefeCFEBasisDofHandlerDealii_h
Definition: CFEBasisDofHandlerDealii.h:57
std::shared_ptr< const dealii::DoFHandler< dim > > getDoFHandler() const
Definition: CFEBasisDofHandlerDealii.t.cpp:1059
FECellIterator beginLocallyOwnedCells() override
Definition: CFEBasisDofHandlerDealii.t.cpp:959
size_type d_numCumulativeLocalCellDofs
Definition: CFEBasisDofHandlerDealii.h:219
const dealii::FiniteElement< dim > & getReferenceFE(const size_type cellId) const
Definition: CFEBasisDofHandlerDealii.t.cpp:1069
std::vector< global_size_type > d_boundaryIds
Definition: CFEBasisDofHandlerDealii.h:222
size_type nCumulativeLocalCellDofs() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1120
const std::vector< global_size_type > & getBoundaryIds() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:948
global_size_type nGlobalNodes() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:809
typename FEBasisDofHandler< ValueTypeBasisCoeff, memorySpace, dim >::const_FECellIterator const_FECellIterator
Definition: CFEBasisDofHandlerDealii.h:64
FECellIterator endLocallyOwnedCells() override
Definition: CFEBasisDofHandlerDealii.t.cpp:970
std::vector< size_type > getGlobalNodeIds() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:919
std::shared_ptr< const ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > getIntrinsicConstraints() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1144
std::shared_ptr< const TriangulationBase > getTriangulation() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:728
std::vector< std::shared_ptr< FECellBase > > d_locallyOwnedCells
Definition: CFEBasisDofHandlerDealii.h:217
void getBasisCenters(std::map< global_size_type, utils::Point > &dofCoords) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1093
size_type totalRanges() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1130
unsigned int getDim() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1045
std::shared_ptr< const ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > d_constraintsLocal
Definition: CFEBasisDofHandlerDealii.h:227
void getCellDofsGlobalIds(size_type cellId, std::vector< global_size_type > &vecGlobalNodeId) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:935
std::vector< double > getBasisFunctionDerivative(const size_type basisId, const utils::Point &point, const size_type derivativeOrder=1) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:711
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:220
bool d_isDistributed
Definition: CFEBasisDofHandlerDealii.h:223
size_type getFEOrder(size_type cellId) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:768
std::shared_ptr< ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > createConstraintsStart() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1154
std::vector< std::pair< global_size_type, global_size_type > > getGlobalRanges() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:877
std::shared_ptr< const TriangulationBase > d_triangulation
Definition: CFEBasisDofHandlerDealii.h:213
std::vector< size_type > getLocalNodeIds(size_type cellId) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:904
std::vector< std::pair< global_size_type, global_size_type > > getLocallyOwnedRanges() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:819
void createConstraintsEnd(std::shared_ptr< ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > constraintsLocal) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1176
size_type nCumulativeLocallyOwnedCellDofs() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1110
FECellIterator beginLocalCells() override
Definition: CFEBasisDofHandlerDealii.t.cpp:1003
size_type nCellDofs(size_type cellId) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:778
double getBasisFunctionValue(const size_type basisId, const utils::Point &point) const override
Definition: CFEBasisDofHandlerDealii.t.cpp:697
std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > getMPIPatternP2P() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1188
std::vector< std::pair< global_size_type, global_size_type > > d_locallyOwnedRanges
Definition: CFEBasisDofHandlerDealii.h:229
size_type nLocalNodes() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:799
size_type nLocalCells() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:738
size_type d_numCumulativeLocallyOwnedCellDofs
Definition: CFEBasisDofHandlerDealii.h:218
FECellIterator endLocalCells() override
Definition: CFEBasisDofHandlerDealii.t.cpp:1014
std::map< BasisIdAttribute, size_type > getBasisAttributeToRangeIdMap() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:892
bool isDistributed() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:1198
std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > d_mpiPatternP2P
Definition: CFEBasisDofHandlerDealii.h:225
bool d_isVariableDofsPerCell
Definition: CFEBasisDofHandlerDealii.h:215
size_type nLocallyOwnedCells() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:748
bool isVariableDofsPerCell() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:789
size_type nGlobalCells() const override
Definition: CFEBasisDofHandlerDealii.t.cpp:758
std::vector< std::shared_ptr< FECellBase > > d_localCells
Definition: CFEBasisDofHandlerDealii.h:216
std::shared_ptr< dealii::DoFHandler< dim > > d_dofHandler
Definition: CFEBasisDofHandlerDealii.h:214
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:829
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