DFT-FE 1.3.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 // resize offset vectors
73 offsetVectors.resize(3);
74
75 for (dftfe::Int i = 0; i < 3; ++i)
76 {
77 for (dftfe::Int j = 0; j < 3; ++j)
78 {
79 offsetVectors[i][j] = -latticeVectors[i][j];
80 }
81 }
82 }
83
84
85 inline double
86 getCosineAngle(dealii::Tensor<1, 3> &Vector1, std::vector<double> &Vector2)
87 {
88 double dotProduct = Vector1[0] * Vector2[0] + Vector1[1] * Vector2[1] +
89 Vector1[2] * Vector2[2];
90 double lengthVector1Sq =
91 (Vector1[0] * Vector1[0] + Vector1[1] * Vector1[1] +
92 Vector1[2] * Vector1[2]);
93 double lengthVector2Sq =
94 (Vector2[0] * Vector2[0] + Vector2[1] * Vector2[1] +
95 Vector2[2] * Vector2[2]);
96
97 double angle = dotProduct / sqrt(lengthVector1Sq * lengthVector2Sq);
98
99 return angle;
100 }
101
102
103 inline double
104 getCosineAngle(std::vector<double> &Vector1, std::vector<double> &Vector2)
105 {
106 double dotProduct = Vector1[0] * Vector2[0] + Vector1[1] * Vector2[1] +
107 Vector1[2] * Vector2[2];
108 double lengthVector1Sq =
109 (Vector1[0] * Vector1[0] + Vector1[1] * Vector1[1] +
110 Vector1[2] * Vector1[2]);
111 double lengthVector2Sq =
112 (Vector2[0] * Vector2[0] + Vector2[1] * Vector2[1] +
113 Vector2[2] * Vector2[2]);
114
115 double angle = dotProduct / sqrt(lengthVector1Sq * lengthVector2Sq);
116
117 return angle;
118 }
119
120
121 inline void
123 dealii::Triangulation<3, 3> &triangulation,
124 std::vector<std::vector<double>> &latticeVectors,
125 const MPI_Comm &mpiCommParent,
126 const dftParameters &dftParams)
127 {
128 dealii::ConditionalOStream pcout(
129 std::cout,
130 (dealii::Utilities::MPI::this_mpi_process(mpiCommParent) == 0));
131 std::vector<std::vector<double>> periodicFaceNormals;
132 std::vector<dealii::Tensor<1, 3>> offsetVectors;
133 // bool periodicX = dftParameters::periodicX, periodicY =
134 // dftParameters::periodicY, periodicZ=dftParameters::periodicZ;
135
136 // compute periodic face normals from lattice vector information
137 computePeriodicFaceNormals(latticeVectors, periodicFaceNormals);
138
139 // compute offset vectors such that lattice vectors plus offset vectors
140 // are alligned along spatial x,y,z directions
141 computeOffsetVectors(latticeVectors, offsetVectors);
142
143 if (dftParams.verbosity >= 4)
144 {
145 pcout << "Periodic Face Normals 1: " << periodicFaceNormals[0][0]
146 << " " << periodicFaceNormals[0][1] << " "
147 << periodicFaceNormals[0][2] << std::endl;
148 pcout << "Periodic Face Normals 2: " << periodicFaceNormals[1][0]
149 << " " << periodicFaceNormals[1][1] << " "
150 << periodicFaceNormals[1][2] << std::endl;
151 pcout << "Periodic Face Normals 3: " << periodicFaceNormals[2][0]
152 << " " << periodicFaceNormals[2][1] << " "
153 << periodicFaceNormals[2][2] << std::endl;
154 }
155
156 dealii::QGauss<2> quadratureFace_formula(2);
157 dealii::FESystem<3> FE(dealii::FE_Q<3>(dealii::QGaussLobatto<1>(2)), 1);
158 dealii::FEFaceValues<3> feFace_values(FE,
159 quadratureFace_formula,
160 dealii::update_normal_vectors);
161
162 typename dealii::Triangulation<3, 3>::active_cell_iterator cell, endc;
163 //
164 // mark faces
165 //
166 // const dftfe::uInt px=dftParameters::periodicX,
167 // py=dftParameters::periodicX, pz=dftParameters::periodicX;
168 //
169 cell = triangulation.begin_active(), endc = triangulation.end();
170 const std::array<dftfe::Int, 3> periodic = {dftParams.periodicX,
171 dftParams.periodicY,
172 dftParams.periodicZ};
173 for (; cell != endc; ++cell)
174 {
175 for (dftfe::uInt f = 0; f < dealii::GeometryInfo<3>::faces_per_cell;
176 ++f)
177 {
178 const dealii::Point<3> face_center = cell->face(f)->center();
179 if (cell->face(f)->at_boundary())
180 {
181 feFace_values.reinit(cell, f);
182 dealii::Tensor<1, 3> faceNormalVector =
183 feFace_values.normal_vector(0);
184
185 // std::cout<<"Face normal vector: "<<faceNormalVector[0]<<"
186 // "<<faceNormalVector[1]<<"
187 // "<<faceNormalVector[2]<<std::endl; pcout<<"Angle :
188 // "<<getCosineAngle(faceNormalVector,periodicFaceNormals[0])<<"
189 // "<<getCosineAngle(faceNormalVector,periodicFaceNormals[1])<<"
190 // "<<getCosineAngle(faceNormalVector,periodicFaceNormals[2])<<std::endl;
191
192 dftfe::uInt i = 1;
193
194 for (dftfe::uInt d = 0; d < 3; ++d)
195 {
196 if (periodic[d] == 1)
197 {
198 if (std::abs(getCosineAngle(faceNormalVector,
199 periodicFaceNormals[d]) -
200 1.0) < 1.0e-05)
201 cell->face(f)->set_boundary_id(i);
202 else if (std::abs(
203 getCosineAngle(faceNormalVector,
204 periodicFaceNormals[d]) +
205 1.0) < 1.0e-05)
206 cell->face(f)->set_boundary_id(i + 1);
207 i = i + 2;
208 }
209 }
210 }
211 }
212 }
213
214 std::vector<dftfe::Int> periodicDirectionVector;
215
216 for (dftfe::uInt d = 0; d < 3; ++d)
217 {
218 if (periodic[d] == 1)
219 {
220 periodicDirectionVector.push_back(d);
221 }
222 }
223
224 // pcout << "Done with Boundary Flags\n";
225 std::vector<dealii::GridTools::PeriodicFacePair<
226 typename dealii::Triangulation<3, 3>::cell_iterator>>
227 periodicity_vector;
228 for (dftfe::Int i = 0;
229 i < std::accumulate(periodic.begin(), periodic.end(), 0);
230 ++i)
231 {
232 dealii::GridTools::collect_periodic_faces(
233 triangulation,
234 /*b_id1*/ 2 * i + 1,
235 /*b_id2*/ 2 * i + 2,
236 /*direction*/ periodicDirectionVector[i],
237 periodicity_vector,
238 offsetVectors[periodicDirectionVector[i]]);
239 // dealii::GridTools::collect_periodic_faces(triangulation, /*b_id1*/
240 // 2*i+1,
241 // /*b_id2*/ 2*i+2,/*direction*/ i, periodicity_vector);
242 }
243 triangulation.add_periodicity(periodicity_vector);
244
245 if (dftParams.verbosity >= 4)
246 pcout << "Periodic Facepairs size: " << periodicity_vector.size()
247 << std::endl;
248 /*
249 for(dftfe::uInt i=0; i< periodicity_vector.size(); ++i)
250 {
251 if (!periodicity_vector[i].cell[0]->active() ||
252 !periodicity_vector[i].cell[1]->active()) continue; if
253 (periodicity_vector[i].cell[0]->is_artificial() ||
254 periodicity_vector[i].cell[1]->is_artificial()) continue;
255
256 std::cout << "matched face pairs: "<<
257 periodicity_vector[i].cell[0]->face(periodicity_vector[i].face_idx[0])->boundary_id()
258 << " "<<
259 periodicity_vector[i].cell[1]->face(periodicity_vector[i].face_idx[1])->boundary_id()<<std::endl;
260 }
261 */
262 }
263
264 } // namespace meshGenUtils
265} // namespace dftfe
266#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:64
bool periodicY
Definition dftParameters.h:64
dftfe::Int verbosity
Definition dftParameters.h:105
bool periodicZ
Definition dftParameters.h:64
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:86
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:122
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10
std::int32_t Int
Definition TypeConfig.h:11