DFT-EFE
 
Loading...
Searching...
No Matches
TriangulationDealiiParallel.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 Vishal Subramanian
24 */
25
26#ifndef dftefeTriangulationDealiiParallel_h
27#define dftefeTriangulationDealiiParallel_h
28
29#include "TriangulationBase.h"
30#include <utils/TypeConfig.h>
32#include <deal.II/distributed/tria.h>
33namespace dftefe
34{
35 namespace basis
36 {
37 template <unsigned int dim>
39 {
40 public:
41 TriangulationDealiiParallel(const MPI_Comm &mpi_communicator);
43
44 void
46 void
48 void
50 const std::vector<unsigned int> &subdivisions,
51 const std::vector<utils::Point> &domainVectors,
52 const std::vector<bool> & isPeriodicFlags) override;
53 void
55 const std::vector<utils::Point> &vertices) override;
56 void
57 shiftTriangulation(const utils::Point &origin) override;
58 void
59 refineGlobal(const unsigned int times = 1) override;
60 void
61 coarsenGlobal(const unsigned int times = 1) override;
62 void
63 clearUserFlags() override;
64 void
66 unsigned int
67 nLocallyOwnedCells() const override;
68 double
69 maxElementLength() const override;
70 double
71 minElementLength() const override;
73 nGlobalCells() const override;
75 nLocalCells() const override;
82 std::vector<size_type>
83 getBoundaryIds() const override;
85 beginLocal() override;
87 endLocal() override;
89 beginLocal() const override;
91 endLocal() const override;
92 unsigned int
93 getDim() const override;
94 std::vector<bool>
95 getPeriodicFlags() const override;
96 void
97 saveRefineFlags(std::vector<bool> &v) const override;
98 void
99 writeToVtkFile(std::ostream &out) const override;
100 std::vector<utils::Point>
101 getDomainVectors() const override;
102
103 // Class specific member function
104
105
106 const dealii::parallel::distributed::Triangulation<dim> &
107 returnDealiiTria() const;
108
109 private:
117 void
118 markPeriodicFaces(const std::vector<bool> & isPeriodicFlags,
119 const std::vector<utils::Point> &domainVectors);
120
121 private:
124 dealii::parallel::distributed::Triangulation<dim> d_triangulationDealii;
125 std::vector<std::shared_ptr<TriangulationCellBase>> d_triaVectorCell;
126 std::vector<bool> d_isPeriodicFlags;
127 std::vector<utils::Point> d_domainVectors;
130 const MPI_Comm &d_mpiDomainCommunicator;
131
132 }; // end of class TriangulationDealiiParallel
133
134 } // end of namespace basis
135
136} // end of namespace dftefe
138#endif // # ifndef dftefeTriangulationDealiiParallel_h
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: TriangulationDealiiParallel.h:39
void shiftTriangulation(const utils::Point &origin) override
Definition: TriangulationDealiiParallel.t.cpp:159
const dealii::parallel::distributed::Triangulation< dim > & returnDealiiTria() const
Definition: TriangulationDealiiParallel.t.cpp:412
void writeToVtkFile(std::ostream &out) const override
Definition: TriangulationDealiiParallel.t.cpp:397
void executeCoarseningAndRefinement() override
Definition: TriangulationDealiiParallel.t.cpp:246
TriangulationBase::TriangulationCellIterator endLocal() override
Definition: TriangulationDealiiParallel.t.cpp:354
std::vector< bool > getPeriodicFlags() const override
Definition: TriangulationDealiiParallel.t.cpp:382
bool isFinalized
Definition: TriangulationDealiiParallel.h:123
void createSingleCellTriangulation(const std::vector< utils::Point > &vertices) override
Definition: TriangulationDealiiParallel.t.cpp:137
unsigned int nLocallyOwnedCells() const override
Definition: TriangulationDealiiParallel.t.cpp:257
unsigned int getDim() const override
Definition: TriangulationDealiiParallel.t.cpp:375
double d_minElemLength
Definition: TriangulationDealiiParallel.h:129
std::vector< size_type > getBoundaryIds() const override
Definition: TriangulationDealiiParallel.t.cpp:283
void clearUserFlags() override
Definition: TriangulationDealiiParallel.t.cpp:235
std::vector< utils::Point > getDomainVectors() const override
Definition: TriangulationDealiiParallel.t.cpp:405
~TriangulationDealiiParallel()
Definition: TriangulationDealiiParallel.t.cpp:26
bool isInitialized
Definition: TriangulationDealiiParallel.h:122
size_type nLocalCells() const override
Definition: TriangulationDealiiParallel.t.cpp:273
const MPI_Comm & d_mpiDomainCommunicator
Definition: TriangulationDealiiParallel.h:130
size_type nGlobalCells() const override
Definition: TriangulationDealiiParallel.t.cpp:265
double d_maxElemLength
Definition: TriangulationDealiiParallel.h:128
std::vector< utils::Point > d_domainVectors
Definition: TriangulationDealiiParallel.h:127
void finalizeTriangulationConstruction() override
Definition: TriangulationDealiiParallel.t.cpp:49
std::vector< bool > d_isPeriodicFlags
Definition: TriangulationDealiiParallel.h:126
double minElementLength() const override
Definition: TriangulationDealiiParallel.t.cpp:323
void saveRefineFlags(std::vector< bool > &v) const override
Definition: TriangulationDealiiParallel.t.cpp:389
TriangulationBase::TriangulationCellIterator beginLocal() override
Definition: TriangulationDealiiParallel.t.cpp:347
void initializeTriangulationConstruction() override
Definition: TriangulationDealiiParallel.t.cpp:33
void createUniformParallelepiped(const std::vector< unsigned int > &subdivisions, const std::vector< utils::Point > &domainVectors, const std::vector< bool > &isPeriodicFlags) override
Definition: TriangulationDealiiParallel.t.cpp:103
double maxElementLength() const override
Definition: TriangulationDealiiParallel.t.cpp:299
dealii::parallel::distributed::Triangulation< dim > d_triangulationDealii
Definition: TriangulationDealiiParallel.h:124
void markPeriodicFaces(const std::vector< bool > &isPeriodicFlags, const std::vector< utils::Point > &domainVectors)
Definition: TriangulationDealiiParallel.t.cpp:175
void refineGlobal(const unsigned int times=1) override
Definition: TriangulationDealiiParallel.t.cpp:211
void coarsenGlobal(const unsigned int times=1) override
Definition: TriangulationDealiiParallel.t.cpp:223
std::vector< std::shared_ptr< TriangulationCellBase > > d_triaVectorCell
Definition: TriangulationDealiiParallel.h:125
Definition: PointImpl.h:13
dealii includes
Definition: AtomFieldDataSpherical.cpp:31
unsigned int size_type
Definition: TypeConfig.h:8