DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
meshMovementGaussian.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
18#ifndef meshMovementGaussian_H_
19#define meshMovementGaussian_H_
20#include "meshMovement.h"
21
22namespace dftfe
23{
24 /**
25 * @brief Class to move triangulation nodes using Gaussian functions attached to control points
26 *
27 * @author Sambit Das
28 */
30 {
31 public:
32 /** @brief Constructor
33 *
34 * @param mpi_comm_parent parent mpi communicator
35 * @param mpi_comm_domain mpi communicator for domain decomposition
36 */
37 meshMovementGaussianClass(const MPI_Comm & mpi_comm_parent,
38 const MPI_Comm & mpi_comm_domaim,
39 const dftParameters &dftParams);
40
41 /** @brief Moves the triangulation corresponding to Gaussians attached to control points
42 *
43 * This functions takes into account the hanging node and periodic
44 * constraints when computing the nodal increment field.
45 *
46 * @param controlPointLocations vector of coordinates of control points
47 * @param controlPointDisplacements vector of displacements of control points
48 * @ controllingParameter constant in the Gaussian function:
49 * exp(-(r/controllingParameter)^pow) where pow is controlled via the input
50 * file parameter (.prm)
51 * @return std::pair<bool,double> mesh quality metrics
52 * pair(bool for is negative jacobian, maximum jacobian ratio)
53 */
54 std::pair<bool, double>
55 moveMesh(const std::vector<dealii::Point<3>> &controlPointLocations,
56 const std::vector<dealii::Tensor<1, 3, double>>
57 & controlPointDisplacements,
58 const std::vector<double> &gaussianWidthParameter,
59 const std::vector<double> &flatTopWidthParameter,
60 const bool moveSubdivided = false);
61
62
63
64 std::pair<bool, double>
65 moveMeshTwoStep(const std::vector<dealii::Point<3>> &controlPointLocations1,
66 const std::vector<dealii::Point<3>> &controlPointLocations2,
67 const std::vector<dealii::Tensor<1, 3, double>>
68 &controlPointDisplacements1,
69 const std::vector<dealii::Tensor<1, 3, double>>
70 & controlPointDisplacements2,
71 const std::vector<double> &controllingParameter1,
72 const std::vector<double> &controllingParameter2,
73 const std::vector<double> &flatTopWidthParameter,
74 const bool moveSubdivided = false);
75
76
77 void
79
80
81 private:
82 /** @brief internal function which computes the nodal increment field in the local processor
83 *
84 */
85 void
86 computeIncrement(const std::vector<dealii::Point<3>> &controlPointLocations,
87 const std::vector<dealii::Tensor<1, 3, double>>
88 & controlPointDisplacements,
89 const std::vector<double> &gaussianWidthParameter,
90 const std::vector<double> &flatTopWidthParameter);
91
92 void
94 const std::vector<dealii::Point<3>> &controlPointLocations1,
95 const std::vector<dealii::Point<3>> &controlPointLocations2,
96 const std::vector<dealii::Tensor<1, 3, double>>
97 &controlPointDisplacements1,
98 const std::vector<dealii::Tensor<1, 3, double>>
99 & controlPointDisplacements2,
100 const std::vector<double> &gaussianWidthParameter1,
101 const std::vector<double> &gaussianWidthParameter2,
102 const std::vector<double> &flatTopWidthParameter);
103 };
104
105} // namespace dftfe
106#endif
Namespace which declares the input parameters and the functions to parse them from the input paramete...
Definition dftParameters.h:35
meshMovementClass(const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const dftParameters &dftParams)
Constructor.
std::pair< bool, double > moveMesh(const std::vector< dealii::Point< 3 > > &controlPointLocations, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements, const std::vector< double > &gaussianWidthParameter, const std::vector< double > &flatTopWidthParameter, const bool moveSubdivided=false)
Moves the triangulation corresponding to Gaussians attached to control points.
std::pair< bool, double > moveMeshTwoStep(const std::vector< dealii::Point< 3 > > &controlPointLocations1, const std::vector< dealii::Point< 3 > > &controlPointLocations2, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements1, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements2, const std::vector< double > &controllingParameter1, const std::vector< double > &controllingParameter2, const std::vector< double > &flatTopWidthParameter, const bool moveSubdivided=false)
void computeIncrementTwoStep(const std::vector< dealii::Point< 3 > > &controlPointLocations1, const std::vector< dealii::Point< 3 > > &controlPointLocations2, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements1, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements2, const std::vector< double > &gaussianWidthParameter1, const std::vector< double > &gaussianWidthParameter2, const std::vector< double > &flatTopWidthParameter)
meshMovementGaussianClass(const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domaim, const dftParameters &dftParams)
Constructor.
void computeIncrement(const std::vector< dealii::Point< 3 > > &controlPointLocations, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements, const std::vector< double > &gaussianWidthParameter, const std::vector< double > &flatTopWidthParameter)
internal function which computes the nodal increment field in the local processor
Definition pseudoPotentialToDftfeConverter.cc:34