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 ~FEBasisDofHandler() = default;
67
68 virtual double
70 const utils::Point &point) const = 0;
71 virtual std::vector<double>
73 const utils::Point &point,
74 const size_type derivativeOrder = 1) const = 0;
75
77 // virtual void
78 // reinit(std::shared_ptr<const TriangulationBase> triangulation,
79 // const size_type feOrder) = 0;
80
81 virtual std::shared_ptr<const TriangulationBase>
82 getTriangulation() const = 0;
83
84 virtual size_type
85 nLocalCells() const = 0;
86 virtual size_type
87 nLocallyOwnedCells() const = 0;
88 virtual size_type
89 nGlobalCells() const = 0;
90 virtual size_type
91 getFEOrder(size_type cellId) const = 0;
92 virtual size_type
93 nCellDofs(size_type cellId) const = 0;
94 virtual bool
96
97 virtual std::vector<std::pair<global_size_type, global_size_type>>
99
100 virtual std::vector<std::pair<global_size_type, global_size_type>>
101 getGlobalRanges() const = 0;
102
103 virtual std::map<BasisIdAttribute, size_type>
105
106 virtual size_type
107 nLocalNodes() const = 0;
108 virtual global_size_type
109 nGlobalNodes() const = 0;
110 virtual std::vector<size_type>
111 getLocalNodeIds(size_type cellId) const = 0;
112 virtual std::vector<size_type>
113 getGlobalNodeIds() const = 0;
114 virtual void
116 size_type cellId,
117 std::vector<global_size_type> &vecGlobalNodeId) const = 0;
118 virtual const std::vector<global_size_type> &
119 getBoundaryIds() const = 0;
120 virtual FECellIterator
122 virtual FECellIterator
128 virtual FECellIterator
130 virtual FECellIterator
133 beginLocalCells() const = 0;
135 endLocalCells() const = 0;
136
137 virtual size_type
139
140 virtual size_type
142
143 // This assumes a linear cell mapping
144 virtual void
146 std::map<global_size_type, utils::Point> &dofCoords) const = 0;
147
148 virtual size_type
149 totalRanges() const = 0;
150
151 virtual unsigned int
152 getDim() const = 0;
153
154 // Additional functions for getting geometric constriants matrix
155 // Additional functions for getting the communication pattern object
156 // for MPI case
157
158 virtual std::shared_ptr<
161
162 // use this to add extra constraints on top of geometric constraints
163 virtual std::shared_ptr<
166
167 // call this after calling start
168 virtual void
171 constraintsLocal) const = 0;
172
173 virtual std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
174 getMPIPatternP2P() const = 0;
175
176 virtual bool
177 isDistributed() const = 0;
178
179 }; // end of FEBasisDofHandler
180 } // end of namespace basis
181} // end of namespace dftefe
182#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 ~FEBasisDofHandler()=default
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