DFT-EFE
 
Loading...
Searching...
No Matches
TriangulationDealiiSerial.h
Go to the documentation of this file.
1#ifndef dftefeTriangulationDealiiSerial_h
2#define dftefeTriangulationDealiiSerial_h
3
4#include "TriangulationBase.h"
5#include <utils/TypeConfig.h>
7#include <deal.II/grid/tria.h>
8namespace dftefe
9{
10 namespace basis
11 {
12 template <unsigned int dim>
14 {
15 public:
18
19 void
21 void
23 void
25 const std::vector<unsigned int> &subdivisions,
26 const std::vector<utils::Point> &domainVectors,
27 const std::vector<bool> & isPeriodicFlags) override;
28 void
30 const std::vector<utils::Point> &vertices) override;
31 void
32 shiftTriangulation(const utils::Point &origin) override;
33 void
34 refineGlobal(const unsigned int times = 1) override;
35 void
36 coarsenGlobal(const unsigned int times = 1) override;
37 void
38 clearUserFlags() override;
39 void
41 unsigned int
42 nLocallyOwnedCells() const override;
43 double
44 maxElementLength() const override;
45 double
46 minElementLength() const override;
48 nGlobalCells() const override;
50 nLocalCells() const override;
57 std::vector<size_type>
58 getBoundaryIds() const override;
60 beginLocal() override;
62 endLocal() override;
64 beginLocal() const override;
66 endLocal() const override;
67 unsigned int
68 getDim() const override;
69 std::vector<bool>
70 getPeriodicFlags() const override;
71 void
72 saveRefineFlags(std::vector<bool> &v) const override;
73 void
74 writeToVtkFile(std::ostream &out) const override;
75 std::vector<utils::Point>
76 getDomainVectors() const override;
77
78 // class specific member function
79 const dealii::Triangulation<dim> &
80 returnDealiiTria() const;
81
82 private:
90 void
91 markPeriodicFaces(const std::vector<bool> & isPeriodicFlags,
92 const std::vector<utils::Point> &domainVectors);
93
94 private:
97 dealii::Triangulation<dim> d_triangulationDealii;
98 std::vector<std::shared_ptr<TriangulationCellBase>> d_triaVectorCell;
99 std::vector<bool> d_isPeriodicFlags;
100 std::vector<utils::Point> d_domainVectors;
101
102 }; // end of class TriangulationDealiiSerial
103
104 } // end of namespace basis
105
106} // end of namespace dftefe
108#endif
An abstract class for the triangulation class. The derived class specialises this class to dealii and...
Definition: TriangulationBase.h:17
std::vector< std::shared_ptr< TriangulationCellBase > >::const_iterator const_TriangulationCellIterator
Definition: TriangulationBase.h:23
std::vector< std::shared_ptr< TriangulationCellBase > >::iterator TriangulationCellIterator
Definition: TriangulationBase.h:21
Definition: TriangulationDealiiSerial.h:14
double maxElementLength() const override
Definition: TriangulationDealiiSerial.t.cpp:264
size_type nGlobalCells() const override
Definition: TriangulationDealiiSerial.t.cpp:234
~TriangulationDealiiSerial()
Definition: TriangulationDealiiSerial.t.cpp:21
void shiftTriangulation(const utils::Point &origin) override
Definition: TriangulationDealiiSerial.t.cpp:129
dealii::Triangulation< dim > d_triangulationDealii
Definition: TriangulationDealiiSerial.h:97
unsigned int nLocallyOwnedCells() const override
Definition: TriangulationDealiiSerial.t.cpp:227
std::vector< utils::Point > getDomainVectors() const override
Definition: TriangulationDealiiSerial.t.cpp:353
void clearUserFlags() override
Definition: TriangulationDealiiSerial.t.cpp:204
TriangulationBase::TriangulationCellIterator endLocal() override
Definition: TriangulationDealiiSerial.t.cpp:303
void finalizeTriangulationConstruction() override
Definition: TriangulationDealiiSerial.t.cpp:44
void initializeTriangulationConstruction() override
Definition: TriangulationDealiiSerial.t.cpp:28
bool isInitialized
Definition: TriangulationDealiiSerial.h:95
unsigned int getDim() const override
Definition: TriangulationDealiiSerial.t.cpp:324
void createSingleCellTriangulation(const std::vector< utils::Point > &vertices) override
Definition: TriangulationDealiiSerial.t.cpp:107
void coarsenGlobal(const unsigned int times=1) override
Definition: TriangulationDealiiSerial.t.cpp:192
void refineGlobal(const unsigned int times=1) override
Definition: TriangulationDealiiSerial.t.cpp:180
TriangulationBase::TriangulationCellIterator beginLocal() override
Definition: TriangulationDealiiSerial.t.cpp:296
void saveRefineFlags(std::vector< bool > &v) const override
Definition: TriangulationDealiiSerial.t.cpp:338
void createUniformParallelepiped(const std::vector< unsigned int > &subdivisions, const std::vector< utils::Point > &domainVectors, const std::vector< bool > &isPeriodicFlags) override
Definition: TriangulationDealiiSerial.t.cpp:73
std::vector< bool > d_isPeriodicFlags
Definition: TriangulationDealiiSerial.h:99
void writeToVtkFile(std::ostream &out) const override
Definition: TriangulationDealiiSerial.t.cpp:345
std::vector< std::shared_ptr< TriangulationCellBase > > d_triaVectorCell
Definition: TriangulationDealiiSerial.h:98
std::vector< utils::Point > d_domainVectors
Definition: TriangulationDealiiSerial.h:100
double minElementLength() const override
Definition: TriangulationDealiiSerial.t.cpp:280
size_type nLocalCells() const override
Definition: TriangulationDealiiSerial.t.cpp:241
void executeCoarseningAndRefinement() override
Definition: TriangulationDealiiSerial.t.cpp:215
TriangulationDealiiSerial()
Definition: TriangulationDealiiSerial.t.cpp:14
void markPeriodicFaces(const std::vector< bool > &isPeriodicFlags, const std::vector< utils::Point > &domainVectors)
Definition: TriangulationDealiiSerial.t.cpp:145
const dealii::Triangulation< dim > & returnDealiiTria() const
Definition: TriangulationDealiiSerial.t.cpp:360
std::vector< bool > getPeriodicFlags() const override
Definition: TriangulationDealiiSerial.t.cpp:331
bool isFinalized
Definition: TriangulationDealiiSerial.h:96
std::vector< size_type > getBoundaryIds() const override
Definition: TriangulationDealiiSerial.t.cpp:248
Definition: PointImpl.h:13
dealii includes
Definition: AtomFieldDataSpherical.cpp:31
unsigned int size_type
Definition: TypeConfig.h:8