AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
analytic1d.hpp
Go to the documentation of this file.
1 // Copyright 2015-2018 The ALMA Project Developers
2 //
3 // Licensed under the Apache License, Version 2.0 (the "License");
4 // you may not use this file except in compliance with the License.
5 // You may obtain a copy of the License at
6 //
7 // http://www.apache.org/licenses/LICENSE-2.0
8 //
9 // Unless required by applicable law or agreed to in writing, software
10 // distributed under the License is distributed on an "AS IS" BASIS,
11 // WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or
12 // implied. See the License for the specific language governing
13 // permissions and limitations under the License.
14 
15 #pragma once
16 
20 
21 #include <constants.hpp>
22 #include <structures.hpp>
23 #include <qpoint_grid.hpp>
24 #include <Eigen/Dense>
25 
26 namespace alma {
27 
28 namespace analytic1D {
29 
35 
37 public:
40  const alma::Gamma_grid* grid,
41  const Eigen::ArrayXXd* w,
42  double T);
43 
45  void setDirection(const Eigen::Vector3d unitvector);
46 
49  void setLogMFPbins(double MFPmin, double MFPmax, int Nbins);
50 
54  void setAutoMFPbins(int Nbins);
55 
59  void setAutoProjMFPbins(int Nbins);
60 
63  void setLogRTbins(double taumin, double taumax, int Nbins);
64 
69  void setAutoRTbins(int Nbins);
70 
73  void setLinOmegabins(double omegamin, double omegamax, int Nbins);
74 
78  void setAutoOmegabins(int Nbins);
79 
81  Eigen::VectorXd getMFPbins();
82 
84  Eigen::VectorXd getRTbins();
85 
87  Eigen::VectorXd getOmegabins();
88 
91  void setBulk();
92 
97  void setInPlaneFilm(double filmthickness,
98  const Eigen::Vector3d normal,
99  double specularity = 0.0);
100  void setCrossPlaneFilm(double filmthickness);
101 
103  double getConductivity();
104 
106  double getSpectralConductivity();
107 
109  double getCapacity();
110 
112  double getDiffusivity();
113 
115  double getDominantProjMFP();
116 
118  double getDominantRT();
119 
121  double getAnisotropyIndex();
122 
124  void resolveByMFP();
125 
127  void resolveByProjMFP();
128 
130  void resolveByRT();
131 
133  void resolveByOmega();
134 
137  Eigen::VectorXd getCumulativeConductivity();
138 
141  Eigen::VectorXd getCumulativeCapacity();
142 
144  void setDOSgridsize(int nbins);
145 
147  Eigen::VectorXd getDOSgrid();
148 
150  Eigen::VectorXd getDOS();
151 
152 private:
154  const alma::Crystal_structure* poscar;
156  const alma::Gamma_grid* grid;
158  const Eigen::ArrayXXd* w;
160  double T;
161 
163  Eigen::Vector3d unitvector;
165  Eigen::Vector3d normalvector;
167  Eigen::VectorXd MFPbins;
169  Eigen::VectorXd taubins;
171  Eigen::VectorXd omegabins;
173  double internal_MFPmax;
175  double internal_MFPprojmax;
177  double internal_taumax;
179  double internal_omegamax;
181  bool thinfilm;
182  bool crossplane;
183  double filmthickness;
184  double specularity; // for in-plane films
186  double kappa;
187  double Cv;
190  double dominantProjMFP;
193  double dominantRT;
197  Eigen::Matrix3d FuchsRotation;
199  void initFuchsRotation();
201  void updateMe();
203  Eigen::VectorXd kappacumul;
205  Eigen::VectorXd Cvcumul;
206 
208  enum kappacumulID {
209  resolve_by_MFP,
210  resolve_by_ProjMFP,
211  resolve_by_RT,
212  resolve_by_omega
213  };
214  int kappacumulIdentifier;
216  int DOS_Nbins;
217 };
218 
225 
227 public:
230  const alma::Gamma_grid* grid,
231  const Eigen::ArrayXXd* w,
232  double T);
233 
235  void setDirection(const Eigen::Vector3d unitvector);
236 
239  void setLinGrid(double ximin, double ximax, int Nxi);
240 
243  void setLogGrid(double ximin, double ximax, int Nxi);
244 
246  void setXiGrid(const Eigen::Ref<const Eigen::VectorXd> xigrid);
247 
249  Eigen::VectorXd getSpatialFrequencies();
250 
253  void normaliseOutput(bool norm);
254 
256  double getDiffusivity();
257 
259  Eigen::VectorXd getPsi();
260 
261 private:
263  const alma::Crystal_structure* poscar;
265  const alma::Gamma_grid* grid;
267  const Eigen::ArrayXXd* w;
269  double T;
270 
272  Eigen::Vector3d unitvector;
274  Eigen::VectorXd xi;
276  bool scaleoutput;
278  double Dbulk;
279  void updateDiffusivity();
281  Eigen::VectorXd psi;
282 };
283 
285 
289 
291 public:
294  const alma::Gamma_grid* grid,
295  const Eigen::ArrayXXd* w,
296  double T);
297 
299  void setDirection(const Eigen::Vector3d unitvector);
300 
303  void setLinSpatialGrid(double ximin, double ximax, int Nxi);
304 
307  void setLogSpatialGrid(double ximin, double ximax, int Nxi);
308 
312  void setLinTemporalGrid(double fmin, double fmax, int Nf);
313 
317  void setLogTemporalGrid(double fmin, double fmax, int Nf);
318 
320  Eigen::VectorXd getSpatialFrequencies();
321 
323  Eigen::VectorXd getTemporalFrequencies();
324 
326  Eigen::MatrixXcd getSPR();
327 
328 private:
330  const alma::Crystal_structure* poscar;
332  const alma::Gamma_grid* grid;
334  const Eigen::ArrayXXd* w;
336  double T;
337 
339  Eigen::Vector3d unitvector;
341  Eigen::VectorXd xi;
343  Eigen::VectorXd f;
344  Eigen::VectorXcd s;
346  Eigen::MatrixXcd Pxis;
347 };
348 
350 
361 
363 public:
366  const alma::Gamma_grid* grid,
367  const Eigen::ArrayXXd* w,
368  double T);
369 
371  void setDirection(const Eigen::Vector3d unitvector);
372 
375  void setLinGrid(double xmin, double xmax, int Nx);
376 
379  void setLogGrid(double xmin, double xmax, int Nx);
380 
385  void declareGridNormalised(bool norm);
386 
387  // Optional normalisation of calculation output
388  // by sqrt(4*pi*Dbulk*t)
389  void normaliseOutput(bool norm);
390 
392  void setTime(double t);
393 
395  void setLogMFPbins(double MFPmin, double MFPmax, int Nbins);
396 
398  Eigen::VectorXd getGrid();
399  Eigen::VectorXd getNormalisedGrid();
400 
402  double getTime();
403 
405  Eigen::VectorXd getMFPbins();
406 
408  double getDiffusivity();
409 
411  Eigen::VectorXd getSourceTransient(
412  const Eigen::Ref<Eigen::VectorXd> timegrid);
413 
415  Eigen::VectorXd getSPR();
416 
418  Eigen::MatrixXd resolveSPRbyMFP();
419 
420 private:
422  const alma::Crystal_structure* poscar;
424  const alma::Gamma_grid* grid;
426  const Eigen::ArrayXXd* w;
428  double T;
429 
431  Eigen::Vector3d unitvector;
433  Eigen::Vector3d filmnormal;
435  Eigen::VectorXd x;
436  bool gridIsNormalised;
438  bool normaliseoutput;
440  Eigen::VectorXd MFPbins;
442  double t;
444  double Dbulk;
445  void updateDiffusivity();
447  Eigen::VectorXd Pxt;
449  Eigen::MatrixXd Pmodes;
450 };
451 
453 
464 
466 public:
469  const alma::Gamma_grid* grid,
470  const Eigen::ArrayXXd* w,
471  double T);
472 
474  void setDirection(const Eigen::Vector3d unitvector);
475 
478  void setLinGrid(double fmin, double fmax, int Nf);
479 
482  void setLogGrid(double fmin, double fmax, int Nf);
483 
485  void setLaplaceGrid(const Eigen::Ref<const Eigen::VectorXcd> sgrid);
486 
488  Eigen::VectorXd getGrid();
489  Eigen::VectorXcd getLaplaceGrid();
490 
492  Eigen::VectorXcd getMSD();
493 
494 private:
496  const alma::Crystal_structure* poscar;
498  const alma::Gamma_grid* grid;
500  const Eigen::ArrayXXd* w;
502  double T;
503 
505  Eigen::Vector3d unitvector;
507  Eigen::VectorXd f;
509  Eigen::VectorXcd s;
511  Eigen::VectorXcd MSD;
512 };
513 
515 
522 
524 public:
527  const alma::Gamma_grid* grid,
528  const Eigen::ArrayXXd* w,
529  double T);
530 
532  void setDirection(const Eigen::Vector3d unitvector);
533 
536  void setLinGrid(double tmin, double tmax, int Nt);
537 
540  void setLogGrid(double tmin, double tmax, int Nt);
541 
543  void setTimeGrid(const Eigen::Ref<const Eigen::VectorXd> tgrid);
544 
546  Eigen::VectorXd getGrid();
547 
550  void normaliseOutput(bool norm);
551 
553  double getDiffusivity();
554 
556  Eigen::VectorXd getMSD();
557 
558 private:
560  const alma::Crystal_structure* poscar;
562  const alma::Gamma_grid* grid;
564  const Eigen::ArrayXXd* w;
566  double T;
567 
569  Eigen::Vector3d unitvector;
571  Eigen::Vector3d filmnormal;
573  Eigen::VectorXd t;
575  double Dbulk;
576  void updateDiffusivity();
578  bool normaliseoutput;
580  Eigen::VectorXd MSD;
581 
583  int GS_depth;
584  Eigen::VectorXd GS_coeffs;
586  double factorial(int n);
587 };
588 } // end namespace analytic1D
589 } // end namespace alma
Eigen::VectorXd getDOSgrid()
Obtain the DOS energy grid (in meV)
Definition: analytic1d.cpp:434
void setDOSgridsize(int nbins)
Set the number of bins to be used in DOS evaluation.
Definition: analytic1d.cpp:430
Class for computing the RTA propagator function psi(xi) of a medium.
Definition: analytic1d.hpp:226
Definition: analytic1d.hpp:26
void setAutoProjMFPbins(int Nbins)
Construct a logarithmic grid of MFP bins witn Nbins elements.
Definition: analytic1d.cpp:72
Class for computing basic thermal properties (kappa, Cv, diffusivity) and cumulative functions (resol...
Definition: analytic1d.hpp:36
Class for computing the exact analytical RTA single pulse energy density response of the infinite bul...
Definition: analytic1d.hpp:290
Class for computing the exact analytical RTA solution for mean square thermal energy displacement...
Definition: analytic1d.hpp:465
double getConductivity()
Retrieve thermal conductivity.
Definition: analytic1d.cpp:492
void resolveByRT()
Choose resolving cumulative curves by relaxation time.
Definition: analytic1d.cpp:564
void setLogRTbins(double taumin, double taumax, int Nbins)
Construct a logarithmic grid of tau bins from taumin to taumax with Nbins elements.
Definition: analytic1d.cpp:88
double getDominantRT()
Retrieve the dominant phonon relaxation time.
Definition: analytic1d.cpp:541
void setAutoMFPbins(int Nbins)
Construct a logarithmic grid of MFP bins witn Nbins elements.
Definition: analytic1d.cpp:78
double getDominantProjMFP()
Retrieve the dominant projected MFP as defined below.
Definition: analytic1d.cpp:526
Eigen::VectorXd getCumulativeCapacity()
Compute the cumulative heat capacity curve with respect to the MFP bins and retrieve it...
Definition: analytic1d.cpp:744
Eigen::VectorXd getDOS()
Calculate and obtain the DOS.
Definition: analytic1d.cpp:442
Classes and functions used to manipulate grids in reciprocal space.
void setAutoOmegabins(int Nbins)
Construct a linear grid of omega bins witn Nbins elements.
Definition: analytic1d.cpp:100
Eigen::VectorXd getRTbins()
Retrieve the relaxation time bins.
Definition: analytic1d.cpp:109
Eigen::VectorXd getMFPbins()
Retrieve the MFP bins.
Definition: analytic1d.cpp:105
void setDirection(const Eigen::Vector3d unitvector)
Function setting the thermal transport axis.
Definition: analytic1d.cpp:60
void setLogMFPbins(double MFPmin, double MFPmax, int Nbins)
Construct a logarithmic grid of MFP bins from MFPmin to MFPmax with Nbins elements.
Definition: analytic1d.cpp:66
double getSpectralConductivity()
Obtain spectrally computed thermal conductivity.
Definition: analytic1d.cpp:312
Eigen::VectorXd getOmegabins()
Retrieve the relaxation time bins.
Definition: analytic1d.cpp:113
void resolveByProjMFP()
Choose resolving cumulative curves by projected MFP.
Definition: analytic1d.cpp:560
Class for computing the exact analytical RTA solution for mean square thermal energy displacement in ...
Definition: analytic1d.hpp:523
void setAutoRTbins(int Nbins)
Construct a logarithmic grid of tau bins witn Nbins elements.
Definition: analytic1d.cpp:83
Eigen::VectorXd getCumulativeConductivity()
Compute the cumulative thermal conductivity curve with respect to the MFP bins and retrieve it...
Definition: analytic1d.cpp:572
BasicProperties_calculator(const alma::Crystal_structure *poscar, const alma::Gamma_grid *grid, const Eigen::ArrayXXd *w, double T)
Constructor: initialise internal variables.
Definition: analytic1d.cpp:36
void resolveByOmega()
Choose resolving cumulative curves by angular frequency.
Definition: analytic1d.cpp:568
double getDiffusivity()
Compute the Fourier diffusivity and retrieve it.
Definition: analytic1d.cpp:511
Definitions of the basic data-handling classes in ALMA.
Objects of this class represent a regular grid with the Gamma point in one corner.
Definition: qpoint_grid.hpp:32
Hold information about a crystal structure.
Definition: structures.hpp:51
void setBulk()
Specify that the medium should be treated as infinite bulk (default option)
Definition: analytic1d.cpp:118
Physical, mathematical and miscellaneous constants used in alma.
Class for computing the approximate analytical RTA single pulse energy density response in real space...
Definition: analytic1d.hpp:362
void resolveByMFP()
Choose resolving cumulative curves by MFP.
Definition: analytic1d.cpp:556
void setInPlaneFilm(double filmthickness, const Eigen::Vector3d normal, double specularity=0.0)
Specify that the medium should be treated as a thin film with provided thickness. ...
Definition: analytic1d.cpp:130
double getAnisotropyIndex()
Compute the anisotropy index kappa_max/kappa_min and retrieve it.
Definition: analytic1d.cpp:483
double getCapacity()
Compute the heat capacity and retrieve it.
Definition: analytic1d.cpp:507
void setLinOmegabins(double omegamin, double omegamax, int Nbins)
Construct a linear grid of omega bins from omegamin to omegamax with Nbins elements.
Definition: analytic1d.cpp:94