DFT-EFE
 
Loading...
Searching...
No Matches
FEBasisDofHandler.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 dftefeFEBasisDofHandler_h
27#define dftefeFEBasisDofHandler_h
28
29#include <utils/TypeConfig.h>
30#include <utils/Point.h>
33#include <basis/FECellBase.h>
35#include <utils/MPIPatternP2P.h>
36#include <map>
37namespace dftefe
38{
39 namespace basis
40 {
41 // add attribute to the classical and enriched ids for locallyOwnedRanges()
43 {
46 };
53 template <typename ValueTypeBasisCoeff,
55 size_type dim>
57 {
58 //
59 // Typedefs
60 //
61 public:
62 typedef std::vector<std::shared_ptr<FECellBase>>::iterator FECellIterator;
63 typedef std::vector<std::shared_ptr<FECellBase>>::const_iterator
65
66 virtual double
68 const utils::Point &point) const = 0;
69 virtual std::vector<double>
71 const utils::Point &point,
72 const size_type derivativeOrder = 1) const = 0;
73
75 // virtual void
76 // reinit(std::shared_ptr<const TriangulationBase> triangulation,
77 // const size_type feOrder) = 0;
78
79 virtual std::shared_ptr<const TriangulationBase>
80 getTriangulation() const = 0;
81
82 virtual size_type
83 nLocalCells() const = 0;
84 virtual size_type
85 nLocallyOwnedCells() const = 0;
86 virtual size_type
87 nGlobalCells() const = 0;
88 virtual size_type
89 getFEOrder(size_type cellId) const = 0;
90 virtual size_type
91 nCellDofs(size_type cellId) const = 0;
92 virtual bool
94
95 virtual std::vector<std::pair<global_size_type, global_size_type>>
97
98 virtual std::vector<std::pair<global_size_type, global_size_type>>
99 getGlobalRanges() const = 0;
100
101 virtual std::map<BasisIdAttribute, size_type>
103
104 virtual size_type
105 nLocalNodes() const = 0;
106 virtual global_size_type
107 nGlobalNodes() const = 0;
108 virtual std::vector<size_type>
109 getLocalNodeIds(size_type cellId) const = 0;
110 virtual std::vector<size_type>
111 getGlobalNodeIds() const = 0;
112 virtual void
114 size_type cellId,
115 std::vector<global_size_type> &vecGlobalNodeId) const = 0;
116 virtual const std::vector<global_size_type> &
117 getBoundaryIds() const = 0;
118 virtual FECellIterator
120 virtual FECellIterator
126 virtual FECellIterator
128 virtual FECellIterator
131 beginLocalCells() const = 0;
133 endLocalCells() const = 0;
134
135 virtual size_type
137
138 virtual size_type
140
141 // This assumes a linear cell mapping
142 virtual void
144 std::map<global_size_type, utils::Point> &dofCoords) const = 0;
145
146 virtual size_type
147 totalRanges() const = 0;
148
149 virtual unsigned int
150 getDim() const = 0;
151
152 // Additional functions for getting geometric constriants matrix
153 // Additional functions for getting the communication pattern object
154 // for MPI case
155
156 virtual std::shared_ptr<
159
160 // use this to add extra constraints on top of geometric constraints
161 virtual std::shared_ptr<
164
165 // call this after calling start
166 virtual void
169 constraintsLocal) const = 0;
170
171 virtual std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
172 getMPIPatternP2P() const = 0;
173
174 virtual bool
175 isDistributed() const = 0;
176
177 }; // end of FEBasisDofHandler
178 } // end of namespace basis
179} // end of namespace dftefe
180#endif // dftefeFEBasisDofHandler_h
Definition: BasisDofHandler.h:41
Definition: ConstraintsLocal.h:43
Definition: FEBasisDofHandler.h:57
virtual unsigned int getDim() const =0
virtual std::vector< std::pair< global_size_type, global_size_type > > getGlobalRanges() const =0
virtual FECellIterator endLocalCells()=0
virtual global_size_type nGlobalNodes() const =0
virtual FECellIterator endLocallyOwnedCells()=0
virtual std::vector< std::pair< global_size_type, global_size_type > > getLocallyOwnedRanges() const =0
virtual FECellIterator beginLocalCells()=0
virtual std::map< BasisIdAttribute, size_type > getBasisAttributeToRangeIdMap() const =0
virtual const_FECellIterator beginLocallyOwnedCells() const =0
virtual const std::vector< global_size_type > & getBoundaryIds() const =0
virtual size_type getFEOrder(size_type cellId) const =0
virtual size_type totalRanges() const =0
virtual size_type nGlobalCells() const =0
virtual size_type nCumulativeLocalCellDofs() const =0
virtual const_FECellIterator endLocalCells() const =0
std::vector< std::shared_ptr< FECellBase > >::iterator FECellIterator
Definition: FEBasisDofHandler.h:62
virtual std::shared_ptr< ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > createConstraintsStart() const =0
virtual FECellIterator beginLocallyOwnedCells()=0
virtual size_type nCellDofs(size_type cellId) const =0
virtual const_FECellIterator beginLocalCells() const =0
virtual size_type nLocalCells() const =0
virtual double getBasisFunctionValue(const size_type basisId, const utils::Point &point) const =0
std::vector< std::shared_ptr< FECellBase > >::const_iterator const_FECellIterator
Definition: FEBasisDofHandler.h:64
virtual void getCellDofsGlobalIds(size_type cellId, std::vector< global_size_type > &vecGlobalNodeId) const =0
virtual std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > getMPIPatternP2P() const =0
virtual void createConstraintsEnd(std::shared_ptr< ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > constraintsLocal) const =0
virtual std::vector< size_type > getGlobalNodeIds() const =0
virtual bool isDistributed() const =0
virtual std::vector< size_type > getLocalNodeIds(size_type cellId) const =0
virtual std::shared_ptr< const ConstraintsLocal< ValueTypeBasisCoeff, memorySpace > > getIntrinsicConstraints() const =0
virtual std::vector< double > getBasisFunctionDerivative(const size_type basisId, const utils::Point &point, const size_type derivativeOrder=1) const =0
virtual bool isVariableDofsPerCell() const =0
virtual size_type nLocalNodes() const =0
virtual size_type nLocallyOwnedCells() const =0
virtual const_FECellIterator endLocallyOwnedCells() const =0
virtual size_type nCumulativeLocallyOwnedCellDofs() const =0
virtual void getBasisCenters(std::map< global_size_type, utils::Point > &dofCoords) const =0
virtual std::shared_ptr< const TriangulationBase > getTriangulation() const =0
Definition: PointImpl.h:13
BasisIdAttribute
Definition: FEBasisDofHandler.h:43
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