DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
meshGenUtils.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16//
17// @author Sambit Das, Phani Motamarri
18//
19
20#ifndef meshGenUtils_H_
21#define meshGenUtils_H_
22
23#include <headers.h>
24#include "dftParameters.h"
25
26namespace dftfe
27{
28 namespace meshGenUtils
29 {
30 inline void
31 cross_product(std::vector<double> &a,
32 std::vector<double> &b,
33 std::vector<double> &crossProduct)
34 {
35 crossProduct.resize(a.size(), 0.0);
36
37 crossProduct[0] = a[1] * b[2] - a[2] * b[1];
38 crossProduct[1] = a[2] * b[0] - a[0] * b[2];
39 crossProduct[2] = a[0] * b[1] - a[1] * b[0];
40 }
41
42 inline void
44 std::vector<std::vector<double>> &latticeVectors,
45 std::vector<std::vector<double>> &periodicFaceNormals)
46 {
47 // resize periodic face normal
48 periodicFaceNormals.resize(3);
49
50
51 // evaluate cross product between first two lattice vectors
52 cross_product(latticeVectors[0],
53 latticeVectors[1],
54 periodicFaceNormals[2]);
55
56 // evaluate cross product between second two lattice vectors
57 cross_product(latticeVectors[1],
58 latticeVectors[2],
59 periodicFaceNormals[0]);
60
61
62 // evaluate cross product between third and first two lattice vectors
63 cross_product(latticeVectors[2],
64 latticeVectors[0],
65 periodicFaceNormals[1]);
66 }
67
68 inline void
69 computeOffsetVectors(std::vector<std::vector<double>> &latticeVectors,
70 std::vector<dealii::Tensor<1, 3>> &offsetVectors)
71 {
72 // create unitVectorsXYZ
73 std::vector<std::vector<double>> unitVectorsXYZ;
74 unitVectorsXYZ.resize(3);
75
76 for (dftfe::Int i = 0; i < 3; ++i)
77 {
78 unitVectorsXYZ[i].resize(3, 0.0);
79 unitVectorsXYZ[i][i] = 0.0;
80 }
81
82
83 // resize offset vectors
84 offsetVectors.resize(3);
85
86 for (dftfe::Int i = 0; i < 3; ++i)
87 {
88 for (dftfe::Int j = 0; j < 3; ++j)
89 {
90 offsetVectors[i][j] = unitVectorsXYZ[i][j] - latticeVectors[i][j];
91 }
92 }
93 }
94
95
96 inline double
97 getCosineAngle(dealii::Tensor<1, 3> &Vector1, std::vector<double> &Vector2)
98 {
99 double dotProduct = Vector1[0] * Vector2[0] + Vector1[1] * Vector2[1] +
100 Vector1[2] * Vector2[2];
101 double lengthVector1Sq =
102 (Vector1[0] * Vector1[0] + Vector1[1] * Vector1[1] +
103 Vector1[2] * Vector1[2]);
104 double lengthVector2Sq =
105 (Vector2[0] * Vector2[0] + Vector2[1] * Vector2[1] +
106 Vector2[2] * Vector2[2]);
107
108 double angle = dotProduct / sqrt(lengthVector1Sq * lengthVector2Sq);
109
110 return angle;
111 }
112
113
114 inline double
115 getCosineAngle(std::vector<double> &Vector1, std::vector<double> &Vector2)
116 {
117 double dotProduct = Vector1[0] * Vector2[0] + Vector1[1] * Vector2[1] +
118 Vector1[2] * Vector2[2];
119 double lengthVector1Sq =
120 (Vector1[0] * Vector1[0] + Vector1[1] * Vector1[1] +
121 Vector1[2] * Vector1[2]);
122 double lengthVector2Sq =
123 (Vector2[0] * Vector2[0] + Vector2[1] * Vector2[1] +
124 Vector2[2] * Vector2[2]);
125
126 double angle = dotProduct / sqrt(lengthVector1Sq * lengthVector2Sq);
127
128 return angle;
129 }
130
131
132 inline void
134 dealii::Triangulation<3, 3> &triangulation,
135 std::vector<std::vector<double>> &latticeVectors,
136 const MPI_Comm &mpiCommParent,
137 const dftParameters &dftParams)
138 {
139 dealii::ConditionalOStream pcout(
140 std::cout,
141 (dealii::Utilities::MPI::this_mpi_process(mpiCommParent) == 0));
142 std::vector<std::vector<double>> periodicFaceNormals;
143 std::vector<dealii::Tensor<1, 3>> offsetVectors;
144 // bool periodicX = dftParameters::periodicX, periodicY =
145 // dftParameters::periodicY, periodicZ=dftParameters::periodicZ;
146
147 // compute periodic face normals from lattice vector information
148 computePeriodicFaceNormals(latticeVectors, periodicFaceNormals);
149
150 // compute offset vectors such that lattice vectors plus offset vectors
151 // are alligned along spatial x,y,z directions
152 computeOffsetVectors(latticeVectors, offsetVectors);
153
154 if (dftParams.verbosity >= 4)
155 {
156 pcout << "Periodic Face Normals 1: " << periodicFaceNormals[0][0]
157 << " " << periodicFaceNormals[0][1] << " "
158 << periodicFaceNormals[0][2] << std::endl;
159 pcout << "Periodic Face Normals 2: " << periodicFaceNormals[1][0]
160 << " " << periodicFaceNormals[1][1] << " "
161 << periodicFaceNormals[1][2] << std::endl;
162 pcout << "Periodic Face Normals 3: " << periodicFaceNormals[2][0]
163 << " " << periodicFaceNormals[2][1] << " "
164 << periodicFaceNormals[2][2] << std::endl;
165 }
166
167 dealii::QGauss<2> quadratureFace_formula(2);
168 dealii::FESystem<3> FE(dealii::FE_Q<3>(dealii::QGaussLobatto<1>(2)), 1);
169 dealii::FEFaceValues<3> feFace_values(FE,
170 quadratureFace_formula,
171 dealii::update_normal_vectors);
172
173 typename dealii::Triangulation<3, 3>::active_cell_iterator cell, endc;
174 //
175 // mark faces
176 //
177 // const dftfe::uInt px=dftParameters::periodicX,
178 // py=dftParameters::periodicX, pz=dftParameters::periodicX;
179 //
180 cell = triangulation.begin_active(), endc = triangulation.end();
181 const std::array<dftfe::Int, 3> periodic = {dftParams.periodicX,
182 dftParams.periodicY,
183 dftParams.periodicZ};
184 for (; cell != endc; ++cell)
185 {
186 for (dftfe::uInt f = 0; f < dealii::GeometryInfo<3>::faces_per_cell;
187 ++f)
188 {
189 const dealii::Point<3> face_center = cell->face(f)->center();
190 if (cell->face(f)->at_boundary())
191 {
192 feFace_values.reinit(cell, f);
193 dealii::Tensor<1, 3> faceNormalVector =
194 feFace_values.normal_vector(0);
195
196 // std::cout<<"Face normal vector: "<<faceNormalVector[0]<<"
197 // "<<faceNormalVector[1]<<"
198 // "<<faceNormalVector[2]<<std::endl; pcout<<"Angle :
199 // "<<getCosineAngle(faceNormalVector,periodicFaceNormals[0])<<"
200 // "<<getCosineAngle(faceNormalVector,periodicFaceNormals[1])<<"
201 // "<<getCosineAngle(faceNormalVector,periodicFaceNormals[2])<<std::endl;
202
203 dftfe::uInt i = 1;
204
205 for (dftfe::uInt d = 0; d < 3; ++d)
206 {
207 if (periodic[d] == 1)
208 {
209 if (std::abs(getCosineAngle(faceNormalVector,
210 periodicFaceNormals[d]) -
211 1.0) < 1.0e-05)
212 cell->face(f)->set_boundary_id(i);
213 else if (std::abs(
214 getCosineAngle(faceNormalVector,
215 periodicFaceNormals[d]) +
216 1.0) < 1.0e-05)
217 cell->face(f)->set_boundary_id(i + 1);
218 i = i + 2;
219 }
220 }
221 }
222 }
223 }
224
225 std::vector<dftfe::Int> periodicDirectionVector;
226
227 for (dftfe::uInt d = 0; d < 3; ++d)
228 {
229 if (periodic[d] == 1)
230 {
231 periodicDirectionVector.push_back(d);
232 }
233 }
234
235 // pcout << "Done with Boundary Flags\n";
236 std::vector<dealii::GridTools::PeriodicFacePair<
237 typename dealii::Triangulation<3, 3>::cell_iterator>>
238 periodicity_vector;
239 for (dftfe::Int i = 0;
240 i < std::accumulate(periodic.begin(), periodic.end(), 0);
241 ++i)
242 {
243 dealii::GridTools::collect_periodic_faces(
244 triangulation,
245 /*b_id1*/ 2 * i + 1,
246 /*b_id2*/ 2 * i + 2,
247 /*direction*/ periodicDirectionVector[i],
248 periodicity_vector,
249 offsetVectors[periodicDirectionVector[i]]);
250 // dealii::GridTools::collect_periodic_faces(triangulation, /*b_id1*/
251 // 2*i+1,
252 // /*b_id2*/ 2*i+2,/*direction*/ i, periodicity_vector);
253 }
254 triangulation.add_periodicity(periodicity_vector);
255
256 if (dftParams.verbosity >= 4)
257 pcout << "Periodic Facepairs size: " << periodicity_vector.size()
258 << std::endl;
259 /*
260 for(dftfe::uInt i=0; i< periodicity_vector.size(); ++i)
261 {
262 if (!periodicity_vector[i].cell[0]->active() ||
263 !periodicity_vector[i].cell[1]->active()) continue; if
264 (periodicity_vector[i].cell[0]->is_artificial() ||
265 periodicity_vector[i].cell[1]->is_artificial()) continue;
266
267 std::cout << "matched face pairs: "<<
268 periodicity_vector[i].cell[0]->face(periodicity_vector[i].face_idx[0])->boundary_id()
269 << " "<<
270 periodicity_vector[i].cell[1]->face(periodicity_vector[i].face_idx[1])->boundary_id()<<std::endl;
271 }
272 */
273 }
274
275 } // namespace meshGenUtils
276} // namespace dftfe
277#endif
Namespace which declares the input parameters and the functions to parse them from the input paramete...
Definition dftParameters.h:36
bool periodicX
Definition dftParameters.h:62
bool periodicY
Definition dftParameters.h:62
dftfe::Int verbosity
Definition dftParameters.h:103
bool periodicZ
Definition dftParameters.h:62
Definition meshGenUtils.h:29
void cross_product(std::vector< double > &a, std::vector< double > &b, std::vector< double > &crossProduct)
Definition meshGenUtils.h:31
void computeOffsetVectors(std::vector< std::vector< double > > &latticeVectors, std::vector< dealii::Tensor< 1, 3 > > &offsetVectors)
Definition meshGenUtils.h:69
double getCosineAngle(dealii::Tensor< 1, 3 > &Vector1, std::vector< double > &Vector2)
Definition meshGenUtils.h:97
void computePeriodicFaceNormals(std::vector< std::vector< double > > &latticeVectors, std::vector< std::vector< double > > &periodicFaceNormals)
Definition meshGenUtils.h:43
void markPeriodicFacesNonOrthogonal(dealii::Triangulation< 3, 3 > &triangulation, std::vector< std::vector< double > > &latticeVectors, const MPI_Comm &mpiCommParent, const dftParameters &dftParams)
Definition meshGenUtils.h:133
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10
std::int32_t Int
Definition TypeConfig.h:11