<!DOCTYPE html>
< html class = "writer-html5" lang = "en" data-content_root = "../../" >
< head >
< meta charset = "utf-8" / >
< meta name = "viewport" content = "width=device-width, initial-scale=1.0" / >
< title > grogupy.grogu — grogupy 1.0.0 documentation< / title >
< link rel = "stylesheet" type = "text/css" href = "../../_static/pygments.css?v=80d5e7a1" / >
< link rel = "stylesheet" type = "text/css" href = "../../_static/css/theme.css?v=e59714d7" / >
< link rel = "stylesheet" type = "text/css" href = "../../_static/sphinx-design.min.css?v=95c83b7e" / >
< script src = "../../_static/jquery.js?v=5d32c60e" > < / script >
< script src = "../../_static/_sphinx_javascript_frameworks_compat.js?v=2cd50e6c" > < / script >
< script src = "../../_static/documentation_options.js?v=8d563738" > < / script >
< script src = "../../_static/doctools.js?v=9a2dae69" > < / script >
< script src = "../../_static/sphinx_highlight.js?v=dc90522c" > < / script >
< script src = "../../_static/design-tabs.js?v=f930bc37" > < / script >
< script src = "../../_static/js/theme.js" > < / script >
< link rel = "index" title = "Index" href = "../../genindex.html" / >
< link rel = "search" title = "Search" href = "../../search.html" / >
< / head >
< body class = "wy-body-for-nav" >
< div class = "wy-grid-for-nav" >
< nav data-toggle = "wy-nav-shift" class = "wy-nav-side" >
< div class = "wy-side-scroll" >
< div class = "wy-side-nav-search" >
< a href = "../../index.html" class = "icon icon-home" >
grogupy
< / a >
< div role = "search" >
< form id = "rtd-search-form" class = "wy-form" action = "../../search.html" method = "get" >
< input type = "text" name = "q" placeholder = "Search docs" aria-label = "Search docs" / >
< input type = "hidden" name = "check_keywords" value = "yes" / >
< input type = "hidden" name = "area" value = "default" / >
< / form >
< / div >
< / div > < div class = "wy-menu wy-menu-vertical" data-spy = "affix" role = "navigation" aria-label = "Navigation menu" >
< p class = "caption" role = "heading" > < span class = "caption-text" > Getting started< / span > < / p >
< ul >
< li class = "toctree-l1" > < a class = "reference internal" href = "../../introduction.html" > Introduction< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "../../quickstart/index.html" > Quickstart< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "../../cite.html" > Citing grogupy< / a > < / li >
< / ul >
< p class = "caption" role = "heading" > < span class = "caption-text" > User Guide< / span > < / p >
< ul >
< li class = "toctree-l1" > < a class = "reference internal" href = "../../tutorials/index.html" > Tutorials< / a > < / li >
< / ul >
< p class = "caption" role = "heading" > < span class = "caption-text" > Advanced usage< / span > < / p >
< ul >
< li class = "toctree-l1" > < a class = "reference internal" href = "../../implementation/modules.html" > src< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "../../implementation/grogupy.html" > grogupy package< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "../../implementation/environment.html" > Environment variables< / a > < / li >
< / ul >
< p class = "caption" role = "heading" > < span class = "caption-text" > Development< / span > < / p >
< ul >
< li class = "toctree-l1" > < a class = "reference internal" href = "../../dev/index.html" > Contributing to grogupy< / a > < / li >
< / ul >
< p class = "caption" role = "heading" > < span class = "caption-text" > Extras< / span > < / p >
< ul >
< li class = "toctree-l1" > < a class = "reference internal" href = "../../changelog/index.html" > Changelog< / a > < / li >
< li class = "toctree-l1" > < a class = "reference internal" href = "../../bibliography.html" > Bibliography< / a > < / li >
< / ul >
< / div >
< / div >
< / nav >
< section data-toggle = "wy-nav-shift" class = "wy-nav-content-wrap" > < nav class = "wy-nav-top" aria-label = "Mobile navigation menu" >
< i data-toggle = "wy-nav-top" class = "fa fa-bars" > < / i >
< a href = "../../index.html" > grogupy< / a >
< / nav >
< div class = "wy-nav-content" >
< div class = "rst-content" >
< div role = "navigation" aria-label = "Page navigation" >
< ul class = "wy-breadcrumbs" >
< li > < a href = "../../index.html" class = "icon icon-home" aria-label = "Home" > < / a > < / li >
< li class = "breadcrumb-item" > < a href = "../index.html" > Module code< / a > < / li >
< li class = "breadcrumb-item active" > grogupy.grogu< / li >
< li class = "wy-breadcrumbs-aside" >
< / li >
< / ul >
< hr / >
< / div >
< div role = "main" class = "document" itemscope = "itemscope" itemtype = "http://schema.org/Article" >
< div itemprop = "articleBody" >
< h1 > Source code for grogupy.grogu< / h1 > < div class = "highlight" > < pre >
< span > < / span > < span class = "c1" > # Copyright (c) [2024] []< / span >
< span class = "c1" > #< / span >
< span class = "c1" > # Permission is hereby granted, free of charge, to any person obtaining a copy< / span >
< span class = "c1" > # of this software and associated documentation files (the " Software" ), to deal< / span >
< span class = "c1" > # in the Software without restriction, including without limitation the rights< / span >
< span class = "c1" > # to use, copy, modify, merge, publish, distribute, sublicense, and/or sell< / span >
< span class = "c1" > # copies of the Software, and to permit persons to whom the Software is< / span >
< span class = "c1" > # furnished to do so, subject to the following conditions:< / span >
< span class = "c1" > #< / span >
< span class = "c1" > # The above copyright notice and this permission notice shall be included in all< / span >
< span class = "c1" > # copies or substantial portions of the Software.< / span >
< span class = "c1" > #< / span >
< span class = "c1" > # THE SOFTWARE IS PROVIDED " AS IS" , WITHOUT WARRANTY OF ANY KIND, EXPRESS OR< / span >
< span class = "c1" > # IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,< / span >
< span class = "c1" > # FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE< / span >
< span class = "c1" > # AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER< / span >
< span class = "c1" > # LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,< / span >
< span class = "c1" > # OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE< / span >
< span class = "c1" > # SOFTWARE.< / span >
< span class = "sd" > " " " Docstring in grogupy.< / span >
< span class = "sd" > " " " < / span >
< span class = "kn" > import< / span > < span class = "nn" > warnings< / span >
< span class = "kn" > from< / span > < span class = "nn" > argparse< / span > < span class = "kn" > import< / span > < span class = "n" > ArgumentParser< / span >
< span class = "kn" > from< / span > < span class = "nn" > sys< / span > < span class = "kn" > import< / span > < span class = "n" > getsizeof< / span >
< span class = "kn" > from< / span > < span class = "nn" > timeit< / span > < span class = "kn" > import< / span > < span class = "n" > default_timer< / span > < span class = "k" > as< / span > < span class = "n" > timer< / span >
< span class = "c1" > # use numpy number of threads one< / span >
< span class = "k" > try< / span > < span class = "p" > :< / span >
< span class = "kn" > from< / span > < span class = "nn" > threadpoolctl< / span > < span class = "kn" > import< / span > < span class = "n" > threadpool_info< / span > < span class = "p" > ,< / span > < span class = "n" > threadpool_limits< / span >
< span class = "n" > user_api< / span > < span class = "o" > =< / span > < span class = "n" > threadpool_info< / span > < span class = "p" > ()[< / span > < span class = "s2" > " user_api" < / span > < span class = "p" > ]< / span >
< span class = "n" > threadpool_limits< / span > < span class = "p" > (< / span > < span class = "n" > limits< / span > < span class = "o" > =< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "n" > user_api< / span > < span class = "o" > =< / span > < span class = "n" > user_api< / span > < span class = "p" > )< / span >
< span class = "k" > except< / span > < span class = "p" > :< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " Warning: threadpoolctl could not make numpy use single thread!" < / span > < span class = "p" > )< / span >
< span class = "kn" > import< / span > < span class = "nn" > numpy< / span > < span class = "k" > as< / span > < span class = "nn" > np< / span >
< span class = "kn" > import< / span > < span class = "nn" > sisl< / span >
< span class = "kn" > from< / span > < span class = "nn" > mpi4py< / span > < span class = "kn" > import< / span > < span class = "n" > MPI< / span >
< span class = "k" > try< / span > < span class = "p" > :< / span >
< span class = "kn" > from< / span > < span class = "nn" > tqdm< / span > < span class = "kn" > import< / span > < span class = "n" > tqdm< / span >
< span class = "n" > tqdm_imported< / span > < span class = "o" > =< / span > < span class = "kc" > True< / span >
< span class = "k" > except< / span > < span class = "p" > :< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " Please install tqdm for nice progress bar." < / span > < span class = "p" > )< / span >
< span class = "n" > tqdm_imported< / span > < span class = "o" > =< / span > < span class = "kc" > False< / span >
< span class = "kn" > from< / span > < span class = "nn" > grogupy< / span > < span class = "kn" > import< / span > < span class = "o" > *< / span >
< div class = "viewcode-block" id = "parse_command_line" >
< a class = "viewcode-back" href = "../../implementation/grogupy.html#grogupy.grogu.parse_command_line" > [docs]< / a >
< span class = "k" > def< / span > < span class = "nf" > parse_command_line< / span > < span class = "p" > ():< / span >
< span class = "w" > < / span > < span class = "sd" > " " " This function can read input from the command line." " " < / span >
< span class = "n" > parser< / span > < span class = "o" > =< / span > < span class = "n" > ArgumentParser< / span > < span class = "p" > ()< / span >
< span class = "n" > parser< / span > < span class = "o" > .< / span > < span class = "n" > add_argument< / span > < span class = "p" > (< / span >
< span class = "s2" > " -i" < / span > < span class = "p" > ,< / span >
< span class = "s2" > " --input" < / span > < span class = "p" > ,< / span >
< span class = "n" > dest< / span > < span class = "o" > =< / span > < span class = "s2" > " infile" < / span > < span class = "p" > ,< / span >
< span class = "n" > default< / span > < span class = "o" > =< / span > < span class = "kc" > None< / span > < span class = "p" > ,< / span >
< span class = "nb" > type< / span > < span class = "o" > =< / span > < span class = "nb" > str< / span > < span class = "p" > ,< / span >
< span class = "n" > help< / span > < span class = "o" > =< / span > < span class = "s2" > " Input file name" < / span > < span class = "p" > ,< / span >
< span class = "n" > required< / span > < span class = "o" > =< / span > < span class = "kc" > True< / span > < span class = "p" > ,< / span >
< span class = "p" > )< / span >
< span class = "n" > parser< / span > < span class = "o" > .< / span > < span class = "n" > add_argument< / span > < span class = "p" > (< / span >
< span class = "s2" > " -o" < / span > < span class = "p" > ,< / span >
< span class = "s2" > " --output" < / span > < span class = "p" > ,< / span >
< span class = "n" > dest< / span > < span class = "o" > =< / span > < span class = "s2" > " outfile" < / span > < span class = "p" > ,< / span >
< span class = "n" > default< / span > < span class = "o" > =< / span > < span class = "kc" > None< / span > < span class = "p" > ,< / span >
< span class = "nb" > type< / span > < span class = "o" > =< / span > < span class = "nb" > str< / span > < span class = "p" > ,< / span >
< span class = "n" > help< / span > < span class = "o" > =< / span > < span class = "s2" > " Output file name" < / span > < span class = "p" > ,< / span >
< span class = "p" > )< / span >
< span class = "c1" > # parser.add_argument(' --scf-orientation' , dest = ' scf_xcf_orientation' , default = None, help = ' Output file name' )< / span >
< span class = "c1" > # parser.add_argument(' --ref-orientation' , dest = ' ref_xcf_orientations' , default = None, help = ' Output file name' )< / span >
< span class = "n" > parser< / span > < span class = "o" > .< / span > < span class = "n" > add_argument< / span > < span class = "p" > (< / span >
< span class = "s2" > " --kset" < / span > < span class = "p" > ,< / span >
< span class = "n" > dest< / span > < span class = "o" > =< / span > < span class = "s2" > " kset" < / span > < span class = "p" > ,< / span >
< span class = "n" > default< / span > < span class = "o" > =< / span > < span class = "kc" > None< / span > < span class = "p" > ,< / span >
< span class = "nb" > type< / span > < span class = "o" > =< / span > < span class = "nb" > int< / span > < span class = "p" > ,< / span >
< span class = "n" > help< / span > < span class = "o" > =< / span > < span class = "s2" > " k-space resolution of calculation" < / span > < span class = "p" > ,< / span >
< span class = "p" > )< / span >
< span class = "n" > parser< / span > < span class = "o" > .< / span > < span class = "n" > add_argument< / span > < span class = "p" > (< / span >
< span class = "s2" > " --kdirs" < / span > < span class = "p" > ,< / span >
< span class = "n" > dest< / span > < span class = "o" > =< / span > < span class = "s2" > " kdirs" < / span > < span class = "p" > ,< / span >
< span class = "n" > default< / span > < span class = "o" > =< / span > < span class = "kc" > None< / span > < span class = "p" > ,< / span >
< span class = "nb" > type< / span > < span class = "o" > =< / span > < span class = "nb" > str< / span > < span class = "p" > ,< / span >
< span class = "n" > help< / span > < span class = "o" > =< / span > < span class = "s2" > " Definition of k-space dimensionality" < / span > < span class = "p" > ,< / span >
< span class = "p" > )< / span >
< span class = "n" > parser< / span > < span class = "o" > .< / span > < span class = "n" > add_argument< / span > < span class = "p" > (< / span >
< span class = "s2" > " --ebot" < / span > < span class = "p" > ,< / span >
< span class = "n" > dest< / span > < span class = "o" > =< / span > < span class = "s2" > " ebot" < / span > < span class = "p" > ,< / span >
< span class = "n" > default< / span > < span class = "o" > =< / span > < span class = "kc" > None< / span > < span class = "p" > ,< / span >
< span class = "nb" > type< / span > < span class = "o" > =< / span > < span class = "nb" > float< / span > < span class = "p" > ,< / span >
< span class = "n" > help< / span > < span class = "o" > =< / span > < span class = "s2" > " Bottom energy of the contour" < / span > < span class = "p" > ,< / span >
< span class = "p" > )< / span >
< span class = "n" > parser< / span > < span class = "o" > .< / span > < span class = "n" > add_argument< / span > < span class = "p" > (< / span >
< span class = "s2" > " --eset" < / span > < span class = "p" > ,< / span >
< span class = "n" > dest< / span > < span class = "o" > =< / span > < span class = "s2" > " eset" < / span > < span class = "p" > ,< / span >
< span class = "n" > default< / span > < span class = "o" > =< / span > < span class = "kc" > None< / span > < span class = "p" > ,< / span >
< span class = "nb" > type< / span > < span class = "o" > =< / span > < span class = "nb" > int< / span > < span class = "p" > ,< / span >
< span class = "n" > help< / span > < span class = "o" > =< / span > < span class = "s2" > " Number of energy points on the contour" < / span > < span class = "p" > ,< / span >
< span class = "p" > )< / span >
< span class = "n" > parser< / span > < span class = "o" > .< / span > < span class = "n" > add_argument< / span > < span class = "p" > (< / span >
< span class = "s2" > " --eset-p" < / span > < span class = "p" > ,< / span >
< span class = "n" > dest< / span > < span class = "o" > =< / span > < span class = "s2" > " esetp" < / span > < span class = "p" > ,< / span >
< span class = "n" > default< / span > < span class = "o" > =< / span > < span class = "kc" > None< / span > < span class = "p" > ,< / span >
< span class = "nb" > type< / span > < span class = "o" > =< / span > < span class = "nb" > int< / span > < span class = "p" > ,< / span >
< span class = "n" > help< / span > < span class = "o" > =< / span > < span class = "s2" > " Parameter tuning the distribution on the contour" < / span > < span class = "p" > ,< / span >
< span class = "p" > )< / span >
< span class = "n" > parser< / span > < span class = "o" > .< / span > < span class = "n" > add_argument< / span > < span class = "p" > (< / span >
< span class = "s2" > " --parallel-green" < / span > < span class = "p" > ,< / span >
< span class = "n" > dest< / span > < span class = "o" > =< / span > < span class = "s2" > " parallel_solver_for_Gk" < / span > < span class = "p" > ,< / span >
< span class = "n" > default< / span > < span class = "o" > =< / span > < span class = "kc" > None< / span > < span class = "p" > ,< / span >
< span class = "nb" > type< / span > < span class = "o" > =< / span > < span class = "nb" > bool< / span > < span class = "p" > ,< / span >
< span class = "n" > help< / span > < span class = "o" > =< / span > < span class = "s2" > " Whether to use the parallel or sequential solver for Greens function" < / span > < span class = "p" > ,< / span >
< span class = "p" > )< / span >
< span class = "n" > parser< / span > < span class = "o" > .< / span > < span class = "n" > add_argument< / span > < span class = "p" > (< / span >
< span class = "s2" > " --padawan-mode" < / span > < span class = "p" > ,< / span >
< span class = "n" > dest< / span > < span class = "o" > =< / span > < span class = "s2" > " padawan_mode" < / span > < span class = "p" > ,< / span >
< span class = "n" > default< / span > < span class = "o" > =< / span > < span class = "kc" > None< / span > < span class = "p" > ,< / span >
< span class = "nb" > type< / span > < span class = "o" > =< / span > < span class = "nb" > bool< / span > < span class = "p" > ,< / span >
< span class = "n" > help< / span > < span class = "o" > =< / span > < span class = "s2" > " If it is on it turns on extra helpful information for new users" < / span > < span class = "p" > ,< / span >
< span class = "p" > )< / span >
< span class = "n" > cmd_line_args< / span > < span class = "o" > =< / span > < span class = "n" > parser< / span > < span class = "o" > .< / span > < span class = "n" > parse_args< / span > < span class = "p" > ()< / span >
< span class = "k" > return< / span > < span class = "n" > cmd_line_args< / span > < / div >
< div class = "viewcode-block" id = "main" >
< a class = "viewcode-back" href = "../../implementation/grogupy.html#grogupy.grogu.main" > [docs]< / a >
< span class = "k" > def< / span > < span class = "nf" > main< / span > < span class = "p" > (< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > ,< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > ,< / span > < span class = "n" > pairs< / span > < span class = "p" > ):< / span >
< span class = "c1" > # runtime information< / span >
< span class = "n" > times< / span > < span class = "o" > =< / span > < span class = "nb" > dict< / span > < span class = "p" > ()< / span >
< span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s2" > " start_time" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > timer< / span > < span class = "p" > ()< / span >
< span class = "c1" > # MPI parameters< / span >
< span class = "n" > comm< / span > < span class = "o" > =< / span > < span class = "n" > MPI< / span > < span class = "o" > .< / span > < span class = "n" > COMM_WORLD< / span >
< span class = "n" > size< / span > < span class = "o" > =< / span > < span class = "n" > comm< / span > < span class = "o" > .< / span > < span class = "n" > Get_size< / span > < span class = "p" > ()< / span >
< span class = "n" > rank< / span > < span class = "o" > =< / span > < span class = "n" > comm< / span > < span class = "o" > .< / span > < span class = "n" > Get_rank< / span > < span class = "p" > ()< / span >
< span class = "n" > root_node< / span > < span class = "o" > =< / span > < span class = "mi" > 0< / span >
< span class = "k" > if< / span > < span class = "n" > rank< / span > < span class = "o" > ==< / span > < span class = "n" > root_node< / span > < span class = "p" > :< / span >
< span class = "c1" > # include parallel size in simulation parameters< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " parallel_size" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > size< / span >
< span class = "c1" > # check versions for debugging< / span >
< span class = "k" > try< / span > < span class = "p" > :< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " sisl version: " < / span > < span class = "p" > ,< / span > < span class = "n" > sisl< / span > < span class = "o" > .< / span > < span class = "n" > __version__< / span > < span class = "p" > )< / span >
< span class = "k" > except< / span > < span class = "p" > :< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " sisl version unknown." < / span > < span class = "p" > )< / span >
< span class = "k" > try< / span > < span class = "p" > :< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " numpy version: " < / span > < span class = "p" > ,< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > __version__< / span > < span class = "p" > )< / span >
< span class = "k" > except< / span > < span class = "p" > :< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " numpy version unknown." < / span > < span class = "p" > )< / span >
< span class = "c1" > # rename outfile< / span >
< span class = "k" > if< / span > < span class = "ow" > not< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " outfile" < / span > < span class = "p" > ]< / span > < span class = "o" > .< / span > < span class = "n" > endswith< / span > < span class = "p" > (< / span > < span class = "s2" > " .pickle" < / span > < span class = "p" > ):< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " outfile" < / span > < span class = "p" > ]< / span > < span class = "o" > +=< / span > < span class = "s2" > " .pickle" < / span >
< span class = "c1" > # if ebot is not given put it 0.1 eV under the smallest energy< / span >
< span class = "k" > if< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " ebot" < / span > < span class = "p" > ]< / span > < span class = "ow" > is< / span > < span class = "kc" > None< / span > < span class = "p" > :< / span >
< span class = "k" > try< / span > < span class = "p" > :< / span >
< span class = "n" > eigfile< / span > < span class = "o" > =< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " infile" < / span > < span class = "p" > ][:< / span > < span class = "o" > -< / span > < span class = "mi" > 3< / span > < span class = "p" > ]< / span > < span class = "o" > +< / span > < span class = "s2" > " EIG" < / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " ebot" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > read_siesta_emin< / span > < span class = "p" > (< / span > < span class = "n" > eigfile< / span > < span class = "p" > )< / span > < span class = "o" > -< / span > < span class = "mf" > 0.1< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " automatic_ebot" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "kc" > True< / span >
< span class = "k" > except< / span > < span class = "p" > :< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " Could not determine ebot." < / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " Parameter was not given and .EIG file was not found." < / span > < span class = "p" > )< / span >
< span class = "k" > else< / span > < span class = "p" > :< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " automatic_ebot" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "kc" > False< / span >
< span class = "c1" > # read sile< / span >
< span class = "n" > fdf< / span > < span class = "o" > =< / span > < span class = "n" > sisl< / span > < span class = "o" > .< / span > < span class = "n" > get_sile< / span > < span class = "p" > (< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " infile" < / span > < span class = "p" > ])< / span >
< span class = "c1" > # read in hamiltonian< / span >
< span class = "n" > dh< / span > < span class = "o" > =< / span > < span class = "n" > fdf< / span > < span class = "o" > .< / span > < span class = "n" > read_hamiltonian< / span > < span class = "p" > ()< / span >
< span class = "c1" > # read unit cell vectors< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " cell" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > fdf< / span > < span class = "o" > .< / span > < span class = "n" > read_geometry< / span > < span class = "p" > ()< / span > < span class = "o" > .< / span > < span class = "n" > cell< / span >
< span class = "c1" > # unit cell index< / span >
< span class = "n" > uc_in_sc_idx< / span > < span class = "o" > =< / span > < span class = "n" > dh< / span > < span class = "o" > .< / span > < span class = "n" > lattice< / span > < span class = "o" > .< / span > < span class = "n" > sc_index< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])< / span >
< span class = "k" > if< / span > < span class = "n" > rank< / span > < span class = "o" > ==< / span > < span class = "n" > root_node< / span > < span class = "p" > :< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " < / span > < span class = "se" > \n\n\n\n\n< / span > < span class = "s2" > " < / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " #################################################################### JOB INFORMATION ###########################################################################" < / span >
< span class = "p" > )< / span >
< span class = "n" > print_job_description< / span > < span class = "p" > (< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " ################################################################################################################################################################" < / span >
< span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " < / span > < span class = "se" > \n\n\n\n\n< / span > < span class = "s2" > " < / span > < span class = "p" > )< / span >
< span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s2" > " setup_time" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > timer< / span > < span class = "p" > ()< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "sa" > f< / span > < span class = "s2" > " Setup done. Elapsed time: < / span > < span class = "si" > {< / span > < span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s1" > ' setup_time' < / span > < span class = "p" > ]< / span > < span class = "si" > }< / span > < span class = "s2" > s" < / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " ================================================================================================================================================================" < / span >
< span class = "p" > )< / span >
< span class = "n" > NO< / span > < span class = "o" > =< / span > < span class = "n" > dh< / span > < span class = "o" > .< / span > < span class = "n" > no< / span > < span class = "c1" > # shorthand for number of orbitals in the unit cell< / span >
< span class = "c1" > # reformat Hamltonian and Overlap matrix for manipulations< / span >
< span class = "n" > hh< / span > < span class = "p" > ,< / span > < span class = "n" > ss< / span > < span class = "o" > =< / span > < span class = "n" > build_hh_ss< / span > < span class = "p" > (< / span > < span class = "n" > dh< / span > < span class = "p" > )< / span >
< span class = "c1" > # symmetrizing Hamiltonian and Overlap matrix to make them hermitian< / span >
< span class = "k" > for< / span > < span class = "n" > i< / span > < span class = "ow" > in< / span > < span class = "nb" > range< / span > < span class = "p" > (< / span > < span class = "n" > dh< / span > < span class = "o" > .< / span > < span class = "n" > lattice< / span > < span class = "o" > .< / span > < span class = "n" > sc_off< / span > < span class = "o" > .< / span > < span class = "n" > shape< / span > < span class = "p" > [< / span > < span class = "mi" > 0< / span > < span class = "p" > ]):< / span >
< span class = "n" > j< / span > < span class = "o" > =< / span > < span class = "n" > dh< / span > < span class = "o" > .< / span > < span class = "n" > lattice< / span > < span class = "o" > .< / span > < span class = "n" > sc_index< / span > < span class = "p" > (< / span > < span class = "o" > -< / span > < span class = "n" > dh< / span > < span class = "o" > .< / span > < span class = "n" > lattice< / span > < span class = "o" > .< / span > < span class = "n" > sc_off< / span > < span class = "p" > [< / span > < span class = "n" > i< / span > < span class = "p" > ])< / span >
< span class = "n" > h1< / span > < span class = "p" > ,< / span > < span class = "n" > h1d< / span > < span class = "o" > =< / span > < span class = "n" > hh< / span > < span class = "p" > [< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > hh< / span > < span class = "p" > [< / span > < span class = "n" > j< / span > < span class = "p" > ]< / span >
< span class = "n" > hh< / span > < span class = "p" > [< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > hh< / span > < span class = "p" > [< / span > < span class = "n" > j< / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "p" > (< / span > < span class = "n" > h1< / span > < span class = "o" > +< / span > < span class = "n" > h1d< / span > < span class = "o" > .< / span > < span class = "n" > T< / span > < span class = "o" > .< / span > < span class = "n" > conj< / span > < span class = "p" > ())< / span > < span class = "o" > /< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "p" > (< / span > < span class = "n" > h1d< / span > < span class = "o" > +< / span > < span class = "n" > h1< / span > < span class = "o" > .< / span > < span class = "n" > T< / span > < span class = "o" > .< / span > < span class = "n" > conj< / span > < span class = "p" > ())< / span > < span class = "o" > /< / span > < span class = "mi" > 2< / span >
< span class = "n" > s1< / span > < span class = "p" > ,< / span > < span class = "n" > s1d< / span > < span class = "o" > =< / span > < span class = "n" > ss< / span > < span class = "p" > [< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > ss< / span > < span class = "p" > [< / span > < span class = "n" > j< / span > < span class = "p" > ]< / span >
< span class = "n" > ss< / span > < span class = "p" > [< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > ss< / span > < span class = "p" > [< / span > < span class = "n" > j< / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "p" > (< / span > < span class = "n" > s1< / span > < span class = "o" > +< / span > < span class = "n" > s1d< / span > < span class = "o" > .< / span > < span class = "n" > T< / span > < span class = "o" > .< / span > < span class = "n" > conj< / span > < span class = "p" > ())< / span > < span class = "o" > /< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "p" > (< / span > < span class = "n" > s1d< / span > < span class = "o" > +< / span > < span class = "n" > s1< / span > < span class = "o" > .< / span > < span class = "n" > T< / span > < span class = "o" > .< / span > < span class = "n" > conj< / span > < span class = "p" > ())< / span > < span class = "o" > /< / span > < span class = "mi" > 2< / span >
< span class = "c1" > # identifying TRS and TRB parts of the Hamiltonian< / span >
< span class = "n" > TAUY< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > kron< / span > < span class = "p" > (< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > eye< / span > < span class = "p" > (< / span > < span class = "n" > NO< / span > < span class = "p" > ),< / span > < span class = "n" > tau_y< / span > < span class = "p" > )< / span >
< span class = "n" > hTR< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "n" > TAUY< / span > < span class = "o" > @< / span > < span class = "n" > hh< / span > < span class = "p" > [< / span > < span class = "n" > i< / span > < span class = "p" > ]< / span > < span class = "o" > .< / span > < span class = "n" > conj< / span > < span class = "p" > ()< / span > < span class = "o" > @< / span > < span class = "n" > TAUY< / span > < span class = "k" > for< / span > < span class = "n" > i< / span > < span class = "ow" > in< / span > < span class = "nb" > range< / span > < span class = "p" > (< / span > < span class = "n" > dh< / span > < span class = "o" > .< / span > < span class = "n" > lattice< / span > < span class = "o" > .< / span > < span class = "n" > nsc< / span > < span class = "o" > .< / span > < span class = "n" > prod< / span > < span class = "p" > ())])< / span >
< span class = "n" > hTRS< / span > < span class = "o" > =< / span > < span class = "p" > (< / span > < span class = "n" > hh< / span > < span class = "o" > +< / span > < span class = "n" > hTR< / span > < span class = "p" > )< / span > < span class = "o" > /< / span > < span class = "mi" > 2< / span >
< span class = "n" > hTRB< / span > < span class = "o" > =< / span > < span class = "p" > (< / span > < span class = "n" > hh< / span > < span class = "o" > -< / span > < span class = "n" > hTR< / span > < span class = "p" > )< / span > < span class = "o" > /< / span > < span class = "mi" > 2< / span >
< span class = "c1" > # extracting the exchange field< / span >
< span class = "n" > traced< / span > < span class = "o" > =< / span > < span class = "p" > [< / span > < span class = "n" > spin_tracer< / span > < span class = "p" > (< / span > < span class = "n" > hTRB< / span > < span class = "p" > [< / span > < span class = "n" > i< / span > < span class = "p" > ])< / span > < span class = "k" > for< / span > < span class = "n" > i< / span > < span class = "ow" > in< / span > < span class = "nb" > range< / span > < span class = "p" > (< / span > < span class = "n" > dh< / span > < span class = "o" > .< / span > < span class = "n" > lattice< / span > < span class = "o" > .< / span > < span class = "n" > nsc< / span > < span class = "o" > .< / span > < span class = "n" > prod< / span > < span class = "p" > ())]< / span > < span class = "c1" > # equation 77< / span >
< span class = "n" > XCF< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > (< / span >
< span class = "p" > [< / span >
< span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "n" > f< / span > < span class = "p" > [< / span > < span class = "s2" > " x" < / span > < span class = "p" > ]< / span > < span class = "o" > /< / span > < span class = "mi" > 2< / span > < span class = "k" > for< / span > < span class = "n" > f< / span > < span class = "ow" > in< / span > < span class = "n" > traced< / span > < span class = "p" > ]),< / span >
< span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "n" > f< / span > < span class = "p" > [< / span > < span class = "s2" > " y" < / span > < span class = "p" > ]< / span > < span class = "o" > /< / span > < span class = "mi" > 2< / span > < span class = "k" > for< / span > < span class = "n" > f< / span > < span class = "ow" > in< / span > < span class = "n" > traced< / span > < span class = "p" > ]),< / span >
< span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "n" > f< / span > < span class = "p" > [< / span > < span class = "s2" > " z" < / span > < span class = "p" > ]< / span > < span class = "o" > /< / span > < span class = "mi" > 2< / span > < span class = "k" > for< / span > < span class = "n" > f< / span > < span class = "ow" > in< / span > < span class = "n" > traced< / span > < span class = "p" > ]),< / span >
< span class = "p" > ]< / span >
< span class = "p" > )< / span >
< span class = "c1" > # check if exchange field has scalar part< / span >
< span class = "n" > max_xcfs< / span > < span class = "o" > =< / span > < span class = "nb" > abs< / span > < span class = "p" > (< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > (< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "n" > f< / span > < span class = "p" > [< / span > < span class = "s2" > " c" < / span > < span class = "p" > ]< / span > < span class = "o" > /< / span > < span class = "mi" > 2< / span > < span class = "k" > for< / span > < span class = "n" > f< / span > < span class = "ow" > in< / span > < span class = "n" > traced< / span > < span class = "p" > ])))< / span > < span class = "o" > .< / span > < span class = "n" > max< / span > < span class = "p" > ()< / span >
< span class = "k" > if< / span > < span class = "n" > max_xcfs< / span > < span class = "o" > > < / span > < span class = "mf" > 1e-12< / span > < span class = "p" > :< / span >
< span class = "n" > warnings< / span > < span class = "o" > .< / span > < span class = "n" > warn< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " Exchange field has non negligible scalar part. Largest value is < / span > < span class = "si" > {< / span > < span class = "n" > max_xcfs< / span > < span class = "si" > }< / span > < span class = "s2" > " < / span >
< span class = "p" > )< / span >
< span class = "k" > if< / span > < span class = "n" > rank< / span > < span class = "o" > ==< / span > < span class = "n" > root_node< / span > < span class = "p" > :< / span >
< span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s2" > " H_and_XCF_time" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > timer< / span > < span class = "p" > ()< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " Hamiltonian and exchange field rotated. Elapsed time: < / span > < span class = "si" > {< / span > < span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s1" > ' H_and_XCF_time' < / span > < span class = "p" > ]< / span > < span class = "si" > }< / span > < span class = "s2" > s" < / span >
< span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " ================================================================================================================================================================" < / span >
< span class = "p" > )< / span >
< span class = "c1" > # initialize pairs and magnetic entities based on input information< / span >
< span class = "n" > pairs< / span > < span class = "p" > ,< / span > < span class = "n" > magnetic_entities< / span > < span class = "o" > =< / span > < span class = "n" > setup_pairs_and_magnetic_entities< / span > < span class = "p" > (< / span >
< span class = "n" > magnetic_entities< / span > < span class = "p" > ,< / span > < span class = "n" > pairs< / span > < span class = "p" > ,< / span > < span class = "n" > dh< / span > < span class = "p" > ,< / span > < span class = "n" > simulation_parameters< / span >
< span class = "p" > )< / span >
< span class = "k" > if< / span > < span class = "n" > rank< / span > < span class = "o" > ==< / span > < span class = "n" > root_node< / span > < span class = "p" > :< / span >
< span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s2" > " site_and_pair_dictionaries_time" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > timer< / span > < span class = "p" > ()< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " Site and pair dictionaries created. Elapsed time: < / span > < span class = "si" > {< / span > < span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s1" > ' site_and_pair_dictionaries_time' < / span > < span class = "p" > ]< / span > < span class = "si" > }< / span > < span class = "s2" > s" < / span >
< span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " ================================================================================================================================================================" < / span >
< span class = "p" > )< / span >
< span class = "c1" > # generate k space sampling< / span >
< span class = "n" > kset< / span > < span class = "o" > =< / span > < span class = "n" > make_kset< / span > < span class = "p" > (< / span >
< span class = "n" > dirs< / span > < span class = "o" > =< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " kdirs" < / span > < span class = "p" > ],< / span > < span class = "n" > NUMK< / span > < span class = "o" > =< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " kset" < / span > < span class = "p" > ]< / span >
< span class = "p" > )< / span >
< span class = "c1" > # generate weights for k points< / span >
< span class = "n" > wkset< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > ones< / span > < span class = "p" > (< / span > < span class = "nb" > len< / span > < span class = "p" > (< / span > < span class = "n" > kset< / span > < span class = "p" > ))< / span > < span class = "o" > /< / span > < span class = "nb" > len< / span > < span class = "p" > (< / span > < span class = "n" > kset< / span > < span class = "p" > )< / span >
< span class = "c1" > # split the k points based on MPI size< / span >
< span class = "n" > kpcs< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array_split< / span > < span class = "p" > (< / span > < span class = "n" > kset< / span > < span class = "p" > ,< / span > < span class = "n" > size< / span > < span class = "p" > )< / span >
< span class = "c1" > # use progress bar if available< / span >
< span class = "k" > if< / span > < span class = "n" > rank< / span > < span class = "o" > ==< / span > < span class = "n" > root_node< / span > < span class = "ow" > and< / span > < span class = "n" > tqdm_imported< / span > < span class = "p" > :< / span >
< span class = "n" > kpcs< / span > < span class = "p" > [< / span > < span class = "n" > root_node< / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > tqdm< / span > < span class = "p" > (< / span > < span class = "n" > kpcs< / span > < span class = "p" > [< / span > < span class = "n" > root_node< / span > < span class = "p" > ],< / span > < span class = "n" > desc< / span > < span class = "o" > =< / span > < span class = "s2" > " k loop" < / span > < span class = "p" > )< / span >
< span class = "k" > if< / span > < span class = "n" > rank< / span > < span class = "o" > ==< / span > < span class = "n" > root_node< / span > < span class = "p" > :< / span >
< span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s2" > " k_set_time" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > timer< / span > < span class = "p" > ()< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "sa" > f< / span > < span class = "s2" > " k set created. Elapsed time: < / span > < span class = "si" > {< / span > < span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s1" > ' k_set_time' < / span > < span class = "p" > ]< / span > < span class = "si" > }< / span > < span class = "s2" > s" < / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " ================================================================================================================================================================" < / span >
< span class = "p" > )< / span >
< span class = "c1" > # this will contain the three Hamiltonian in the< / span >
< span class = "c1" > # reference directions needed to calculate the energy< / span >
< span class = "c1" > # variations upon rotation< / span >
< span class = "n" > hamiltonians< / span > < span class = "o" > =< / span > < span class = "p" > []< / span >
< span class = "c1" > # iterate over the reference directions (quantization axes)< / span >
< span class = "k" > for< / span > < span class = "n" > i< / span > < span class = "p" > ,< / span > < span class = "n" > orient< / span > < span class = "ow" > in< / span > < span class = "nb" > enumerate< / span > < span class = "p" > (< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " ref_xcf_orientations" < / span > < span class = "p" > ]):< / span >
< span class = "c1" > # obtain rotated exchange field and Hamiltonian< / span >
< span class = "n" > R< / span > < span class = "o" > =< / span > < span class = "n" > RotMa2b< / span > < span class = "p" > (< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " scf_xcf_orientation" < / span > < span class = "p" > ],< / span > < span class = "n" > orient< / span > < span class = "p" > [< / span > < span class = "s2" > " o" < / span > < span class = "p" > ])< / span >
< span class = "n" > rot_XCF< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > einsum< / span > < span class = "p" > (< / span > < span class = "s2" > " ij,jklm-> iklm" < / span > < span class = "p" > ,< / span > < span class = "n" > R< / span > < span class = "p" > ,< / span > < span class = "n" > XCF< / span > < span class = "p" > )< / span >
< span class = "n" > rot_H_XCF< / span > < span class = "o" > =< / span > < span class = "nb" > sum< / span > < span class = "p" > (< / span >
< span class = "p" > [< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > kron< / span > < span class = "p" > (< / span > < span class = "n" > rot_XCF< / span > < span class = "p" > [< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > tau< / span > < span class = "p" > )< / span > < span class = "k" > for< / span > < span class = "n" > i< / span > < span class = "p" > ,< / span > < span class = "n" > tau< / span > < span class = "ow" > in< / span > < span class = "nb" > enumerate< / span > < span class = "p" > ([< / span > < span class = "n" > tau_x< / span > < span class = "p" > ,< / span > < span class = "n" > tau_y< / span > < span class = "p" > ,< / span > < span class = "n" > tau_z< / span > < span class = "p" > ])]< / span >
< span class = "p" > )< / span >
< span class = "n" > rot_H_XCF_uc< / span > < span class = "o" > =< / span > < span class = "n" > rot_H_XCF< / span > < span class = "p" > [< / span > < span class = "n" > uc_in_sc_idx< / span > < span class = "p" > ]< / span >
< span class = "c1" > # obtain total Hamiltonian with the rotated exchange field< / span >
< span class = "n" > rot_H< / span > < span class = "o" > =< / span > < span class = "n" > hTRS< / span > < span class = "o" > +< / span > < span class = "n" > rot_H_XCF< / span > < span class = "c1" > # equation 76< / span >
< span class = "c1" > # store the relevant information of the Hamiltonian< / span >
< span class = "n" > hamiltonians< / span > < span class = "o" > .< / span > < span class = "n" > append< / span > < span class = "p" > (< / span > < span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > orient< / span > < span class = "o" > =< / span > < span class = "n" > orient< / span > < span class = "p" > [< / span > < span class = "s2" > " o" < / span > < span class = "p" > ],< / span > < span class = "n" > H< / span > < span class = "o" > =< / span > < span class = "n" > rot_H< / span > < span class = "p" > ))< / span >
< span class = "c1" > # these are the rotations (for now) perpendicular to the quantization axis< / span >
< span class = "k" > for< / span > < span class = "n" > u< / span > < span class = "ow" > in< / span > < span class = "n" > orient< / span > < span class = "p" > [< / span > < span class = "s2" > " vw" < / span > < span class = "p" > ]:< / span >
< span class = "c1" > # section 2.H< / span >
< span class = "n" > Tu< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > kron< / span > < span class = "p" > (< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > eye< / span > < span class = "p" > (< / span > < span class = "n" > NO< / span > < span class = "p" > ,< / span > < span class = "n" > dtype< / span > < span class = "o" > =< / span > < span class = "nb" > int< / span > < span class = "p" > ),< / span > < span class = "n" > tau_u< / span > < span class = "p" > (< / span > < span class = "n" > u< / span > < span class = "p" > ))< / span >
< span class = "n" > Vu1< / span > < span class = "p" > ,< / span > < span class = "n" > Vu2< / span > < span class = "o" > =< / span > < span class = "n" > calc_Vu< / span > < span class = "p" > (< / span > < span class = "n" > rot_H_XCF_uc< / span > < span class = "p" > ,< / span > < span class = "n" > Tu< / span > < span class = "p" > )< / span >
< span class = "k" > for< / span > < span class = "n" > mag_ent< / span > < span class = "ow" > in< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > :< / span >
< span class = "n" > idx< / span > < span class = "o" > =< / span > < span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " spin_box_indices" < / span > < span class = "p" > ]< / span >
< span class = "c1" > # fill up the perturbed potentials (for now) based on the on-site projections< / span >
< span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " Vu1" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ]< / span > < span class = "o" > .< / span > < span class = "n" > append< / span > < span class = "p" > (< / span > < span class = "n" > onsite_projection< / span > < span class = "p" > (< / span > < span class = "n" > Vu1< / span > < span class = "p" > ,< / span > < span class = "n" > idx< / span > < span class = "p" > ,< / span > < span class = "n" > idx< / span > < span class = "p" > ))< / span >
< span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " Vu2" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ]< / span > < span class = "o" > .< / span > < span class = "n" > append< / span > < span class = "p" > (< / span > < span class = "n" > onsite_projection< / span > < span class = "p" > (< / span > < span class = "n" > Vu2< / span > < span class = "p" > ,< / span > < span class = "n" > idx< / span > < span class = "p" > ,< / span > < span class = "n" > idx< / span > < span class = "p" > ))< / span >
< span class = "k" > if< / span > < span class = "n" > rank< / span > < span class = "o" > ==< / span > < span class = "n" > root_node< / span > < span class = "p" > :< / span >
< span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s2" > " reference_rotations_time" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > timer< / span > < span class = "p" > ()< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " Rotations done perpendicular to quantization axis. Elapsed time: < / span > < span class = "si" > {< / span > < span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s1" > ' reference_rotations_time' < / span > < span class = "p" > ]< / span > < span class = "si" > }< / span > < span class = "s2" > s" < / span >
< span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " ================================================================================================================================================================" < / span >
< span class = "p" > )< / span >
< span class = "c1" > # provide helpful information to estimate the runtime and memory< / span >
< span class = "c1" > # requirements of the Greens function calculations< / span >
< span class = "k" > if< / span > < span class = "n" > rank< / span > < span class = "o" > ==< / span > < span class = "n" > root_node< / span > < span class = "p" > :< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " Starting matrix inversions." < / span > < span class = "p" > )< / span >
< span class = "k" > if< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " padawan_mode" < / span > < span class = "p" > ]:< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " Padawan mode: " < / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "sa" > f< / span > < span class = "s2" > " Total number of k points: < / span > < span class = "si" > {< / span > < span class = "n" > kset< / span > < span class = "o" > .< / span > < span class = "n" > shape< / span > < span class = "p" > [< / span > < span class = "mi" > 0< / span > < span class = "p" > ]< / span > < span class = "si" > }< / span > < span class = "s2" > " < / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " Number of energy samples per k point: < / span > < span class = "si" > {< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s1" > ' eset' < / span > < span class = "p" > ]< / span > < span class = "si" > }< / span > < span class = "s2" > " < / span >
< span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "sa" > f< / span > < span class = "s2" > " Total number of directions: < / span > < span class = "si" > {< / span > < span class = "nb" > len< / span > < span class = "p" > (< / span > < span class = "n" > hamiltonians< / span > < span class = "p" > )< / span > < span class = "si" > }< / span > < span class = "s2" > " < / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " Total number of matrix inversions: < / span > < span class = "si" > {< / span > < span class = "n" > kset< / span > < span class = "o" > .< / span > < span class = "n" > shape< / span > < span class = "p" > [< / span > < span class = "mi" > 0< / span > < span class = "p" > ]< / span > < span class = "w" > < / span > < span class = "o" > *< / span > < span class = "w" > < / span > < span class = "nb" > len< / span > < span class = "p" > (< / span > < span class = "n" > hamiltonians< / span > < span class = "p" > )< / span > < span class = "w" > < / span > < span class = "o" > *< / span > < span class = "w" > < / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s1" > ' eset' < / span > < span class = "p" > ]< / span > < span class = "si" > }< / span > < span class = "s2" > " < / span >
< span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " The shape of the Hamiltonian and the Greens function is < / span > < span class = "si" > {< / span > < span class = "n" > NO< / span > < span class = "si" > }< / span > < span class = "s2" > x< / span > < span class = "si" > {< / span > < span class = "n" > NO< / span > < span class = "si" > }< / span > < span class = "s2" > =< / span > < span class = "si" > {< / span > < span class = "n" > NO< / span > < span class = "o" > *< / span > < span class = "n" > NO< / span > < span class = "si" > }< / span > < span class = "s2" > " < / span >
< span class = "p" > )< / span >
< span class = "c1" > # https://stackoverflow.com/questions/70746660/how-to-predict-memory-requirement-for-np-linalg-inv< / span >
< span class = "c1" > # memory is O(64 n**2) for complex matrices< / span >
< span class = "n" > memory_size< / span > < span class = "o" > =< / span > < span class = "n" > getsizeof< / span > < span class = "p" > (< / span > < span class = "n" > hamiltonians< / span > < span class = "p" > [< / span > < span class = "mi" > 0< / span > < span class = "p" > ][< / span > < span class = "s2" > " H" < / span > < span class = "p" > ]< / span > < span class = "o" > .< / span > < span class = "n" > base< / span > < span class = "p" > )< / span > < span class = "o" > /< / span > < span class = "mi" > 1024< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " Memory taken by a single Hamiltonian is: < / span > < span class = "si" > {< / span > < span class = "n" > getsizeof< / span > < span class = "p" > (< / span > < span class = "n" > hamiltonians< / span > < span class = "p" > [< / span > < span class = "mi" > 0< / span > < span class = "p" > ][< / span > < span class = "s1" > ' H' < / span > < span class = "p" > ]< / span > < span class = "o" > .< / span > < span class = "n" > base< / span > < span class = "p" > )< / span > < span class = "w" > < / span > < span class = "o" > /< / span > < span class = "w" > < / span > < span class = "mi" > 1024< / span > < span class = "si" > }< / span > < span class = "s2" > KB" < / span >
< span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "sa" > f< / span > < span class = "s2" > " Expected memory usage per matrix inversion: < / span > < span class = "si" > {< / span > < span class = "n" > memory_size< / span > < span class = "w" > < / span > < span class = "o" > *< / span > < span class = "w" > < / span > < span class = "mi" > 32< / span > < span class = "si" > }< / span > < span class = "s2" > KB" < / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " Expected memory usage per k point for parallel inversion: < / span > < span class = "si" > {< / span > < span class = "n" > memory_size< / span > < span class = "w" > < / span > < span class = "o" > *< / span > < span class = "w" > < / span > < span class = "nb" > len< / span > < span class = "p" > (< / span > < span class = "n" > hamiltonians< / span > < span class = "p" > )< / span > < span class = "w" > < / span > < span class = "o" > *< / span > < span class = "w" > < / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s1" > ' eset' < / span > < span class = "p" > ]< / span > < span class = "w" > < / span > < span class = "o" > *< / span > < span class = "w" > < / span > < span class = "mi" > 32< / span > < span class = "si" > }< / span > < span class = "s2" > KB" < / span >
< span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " Expected memory usage on root node: < / span > < span class = "si" > {< / span > < span class = "nb" > len< / span > < span class = "p" > (< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array_split< / span > < span class = "p" > (< / span > < span class = "n" > kset< / span > < span class = "p" > ,< / span > < span class = "w" > < / span > < span class = "n" > size< / span > < span class = "p" > )[< / span > < span class = "mi" > 0< / span > < span class = "p" > ])< / span > < span class = "w" > < / span > < span class = "o" > *< / span > < span class = "w" > < / span > < span class = "n" > memory_size< / span > < span class = "w" > < / span > < span class = "o" > *< / span > < span class = "w" > < / span > < span class = "nb" > len< / span > < span class = "p" > (< / span > < span class = "n" > hamiltonians< / span > < span class = "p" > )< / span > < span class = "w" > < / span > < span class = "o" > *< / span > < span class = "w" > < / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s1" > ' eset' < / span > < span class = "p" > ]< / span > < span class = "w" > < / span > < span class = "o" > *< / span > < span class = "w" > < / span > < span class = "mi" > 32< / span > < span class = "w" > < / span > < span class = "o" > /< / span > < span class = "w" > < / span > < span class = "mi" > 1024< / span > < span class = "si" > }< / span > < span class = "s2" > MB" < / span >
< span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " ================================================================================================================================================================" < / span >
< span class = "p" > )< / span >
< span class = "c1" > # MPI barrier< / span >
< span class = "n" > comm< / span > < span class = "o" > .< / span > < span class = "n" > Barrier< / span > < span class = "p" > ()< / span >
< span class = "c1" > # make energy contour< / span >
< span class = "n" > cont< / span > < span class = "o" > =< / span > < span class = "n" > make_contour< / span > < span class = "p" > (< / span >
< span class = "n" > emin< / span > < span class = "o" > =< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " ebot" < / span > < span class = "p" > ],< / span >
< span class = "n" > enum< / span > < span class = "o" > =< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " eset" < / span > < span class = "p" > ],< / span >
< span class = "n" > p< / span > < span class = "o" > =< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " esetp" < / span > < span class = "p" > ],< / span >
< span class = "p" > )< / span >
< span class = "n" > eran< / span > < span class = "o" > =< / span > < span class = "n" > cont< / span > < span class = "o" > .< / span > < span class = "n" > ze< / span >
< span class = "c1" > # sampling the integrand on the contour and the BZ< / span >
< span class = "k" > for< / span > < span class = "n" > k< / span > < span class = "ow" > in< / span > < span class = "n" > kpcs< / span > < span class = "p" > [< / span > < span class = "n" > rank< / span > < span class = "p" > ]:< / span >
< span class = "c1" > # weight of k point in BZ integral< / span >
< span class = "n" > wk< / span > < span class = "o" > =< / span > < span class = "n" > wkset< / span > < span class = "p" > [< / span > < span class = "n" > rank< / span > < span class = "p" > ]< / span >
< span class = "c1" > # iterate over reference directions< / span >
< span class = "k" > for< / span > < span class = "n" > i< / span > < span class = "p" > ,< / span > < span class = "n" > hamiltonian_orientation< / span > < span class = "ow" > in< / span > < span class = "nb" > enumerate< / span > < span class = "p" > (< / span > < span class = "n" > hamiltonians< / span > < span class = "p" > ):< / span >
< span class = "c1" > # calculate Hamiltonian and Overlap matrix in a given k point< / span >
< span class = "n" > H< / span > < span class = "o" > =< / span > < span class = "n" > hamiltonian_orientation< / span > < span class = "p" > [< / span > < span class = "s2" > " H" < / span > < span class = "p" > ]< / span >
< span class = "n" > HK< / span > < span class = "p" > ,< / span > < span class = "n" > SK< / span > < span class = "o" > =< / span > < span class = "n" > hsk< / span > < span class = "p" > (< / span > < span class = "n" > H< / span > < span class = "p" > ,< / span > < span class = "n" > ss< / span > < span class = "p" > ,< / span > < span class = "n" > dh< / span > < span class = "o" > .< / span > < span class = "n" > sc_off< / span > < span class = "p" > ,< / span > < span class = "n" > k< / span > < span class = "p" > )< / span >
< span class = "k" > if< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " parallel_solver_for_Gk" < / span > < span class = "p" > ]:< / span >
< span class = "n" > Gk< / span > < span class = "o" > =< / span > < span class = "n" > parallel_Gk< / span > < span class = "p" > (< / span > < span class = "n" > HK< / span > < span class = "p" > ,< / span > < span class = "n" > SK< / span > < span class = "p" > ,< / span > < span class = "n" > eran< / span > < span class = "p" > ,< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " eset" < / span > < span class = "p" > ])< / span >
< span class = "k" > else< / span > < span class = "p" > :< / span >
< span class = "c1" > # solve Greens function sequentially for the energies, because of memory bound< / span >
< span class = "n" > Gk< / span > < span class = "o" > =< / span > < span class = "n" > sequential_GK< / span > < span class = "p" > (< / span > < span class = "n" > HK< / span > < span class = "p" > ,< / span > < span class = "n" > SK< / span > < span class = "p" > ,< / span > < span class = "n" > eran< / span > < span class = "p" > ,< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " eset" < / span > < span class = "p" > ])< / span >
< span class = "c1" > # store the Greens function slice of the magnetic entities< / span >
< span class = "k" > for< / span > < span class = "n" > mag_ent< / span > < span class = "ow" > in< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > :< / span >
< span class = "n" > idx< / span > < span class = "o" > =< / span > < span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " spin_box_indices" < / span > < span class = "p" > ]< / span >
< span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " Gii_tmp" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ]< / span > < span class = "o" > +=< / span > < span class = "n" > onsite_projection< / span > < span class = "p" > (< / span > < span class = "n" > Gk< / span > < span class = "p" > ,< / span > < span class = "n" > idx< / span > < span class = "p" > ,< / span > < span class = "n" > idx< / span > < span class = "p" > )< / span > < span class = "o" > *< / span > < span class = "n" > wk< / span >
< span class = "k" > for< / span > < span class = "n" > pair< / span > < span class = "ow" > in< / span > < span class = "n" > pairs< / span > < span class = "p" > :< / span >
< span class = "c1" > # add phase shift based on the cell difference< / span >
< span class = "n" > phase< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > exp< / span > < span class = "p" > (< / span > < span class = "mi" > 1< / span > < span class = "n" > j< / span > < span class = "o" > *< / span > < span class = "mi" > 2< / span > < span class = "o" > *< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > pi< / span > < span class = "o" > *< / span > < span class = "n" > k< / span > < span class = "o" > @< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " Ruc" < / span > < span class = "p" > ]< / span > < span class = "o" > .< / span > < span class = "n" > T< / span > < span class = "p" > )< / span >
< span class = "c1" > # get the pair orbital sizes from the magnetic entities< / span >
< span class = "n" > ai< / span > < span class = "o" > =< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > [< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " ai" < / span > < span class = "p" > ]][< / span > < span class = "s2" > " spin_box_indices" < / span > < span class = "p" > ]< / span >
< span class = "n" > aj< / span > < span class = "o" > =< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > [< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " aj" < / span > < span class = "p" > ]][< / span > < span class = "s2" > " spin_box_indices" < / span > < span class = "p" > ]< / span >
< span class = "c1" > # store the Greens function slice of the magnetic entities< / span >
< span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " Gij_tmp" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ]< / span > < span class = "o" > +=< / span > < span class = "n" > onsite_projection< / span > < span class = "p" > (< / span > < span class = "n" > Gk< / span > < span class = "p" > ,< / span > < span class = "n" > ai< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "p" > )< / span > < span class = "o" > *< / span > < span class = "n" > phase< / span > < span class = "o" > *< / span > < span class = "n" > wk< / span >
< span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " Gji_tmp" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ]< / span > < span class = "o" > +=< / span > < span class = "n" > onsite_projection< / span > < span class = "p" > (< / span > < span class = "n" > Gk< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "p" > ,< / span > < span class = "n" > ai< / span > < span class = "p" > )< / span > < span class = "o" > /< / span > < span class = "n" > phase< / span > < span class = "o" > *< / span > < span class = "n" > wk< / span >
< span class = "c1" > # summ reduce partial results of mpi nodes< / span >
< span class = "k" > for< / span > < span class = "n" > i< / span > < span class = "ow" > in< / span > < span class = "nb" > range< / span > < span class = "p" > (< / span > < span class = "nb" > len< / span > < span class = "p" > (< / span > < span class = "n" > hamiltonians< / span > < span class = "p" > )):< / span >
< span class = "k" > for< / span > < span class = "n" > mag_ent< / span > < span class = "ow" > in< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > :< / span >
< span class = "n" > comm< / span > < span class = "o" > .< / span > < span class = "n" > Reduce< / span > < span class = "p" > (< / span > < span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " Gii_tmp" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " Gii" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > root< / span > < span class = "o" > =< / span > < span class = "n" > root_node< / span > < span class = "p" > )< / span >
< span class = "k" > for< / span > < span class = "n" > pair< / span > < span class = "ow" > in< / span > < span class = "n" > pairs< / span > < span class = "p" > :< / span >
< span class = "n" > comm< / span > < span class = "o" > .< / span > < span class = "n" > Reduce< / span > < span class = "p" > (< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " Gij_tmp" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " Gij" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > root< / span > < span class = "o" > =< / span > < span class = "n" > root_node< / span > < span class = "p" > )< / span >
< span class = "n" > comm< / span > < span class = "o" > .< / span > < span class = "n" > Reduce< / span > < span class = "p" > (< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " Gji_tmp" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " Gji" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > root< / span > < span class = "o" > =< / span > < span class = "n" > root_node< / span > < span class = "p" > )< / span >
< span class = "k" > if< / span > < span class = "n" > rank< / span > < span class = "o" > ==< / span > < span class = "n" > root_node< / span > < span class = "p" > :< / span >
< span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s2" > " green_function_inversion_time" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > timer< / span > < span class = "p" > ()< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "sa" > f< / span > < span class = "s2" > " Calculated Greens functions. Elapsed time: < / span > < span class = "si" > {< / span > < span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s1" > ' green_function_inversion_time' < / span > < span class = "p" > ]< / span > < span class = "si" > }< / span > < span class = "s2" > s" < / span >
< span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " ================================================================================================================================================================" < / span >
< span class = "p" > )< / span >
< span class = "k" > if< / span > < span class = "n" > rank< / span > < span class = "o" > ==< / span > < span class = "n" > root_node< / span > < span class = "p" > :< / span >
< span class = "c1" > # iterate over the magnetic entities< / span >
< span class = "k" > for< / span > < span class = "n" > tracker< / span > < span class = "p" > ,< / span > < span class = "n" > mag_ent< / span > < span class = "ow" > in< / span > < span class = "nb" > enumerate< / span > < span class = "p" > (< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > ):< / span >
< span class = "c1" > # iterate over the quantization axes< / span >
< span class = "k" > for< / span > < span class = "n" > i< / span > < span class = "p" > ,< / span > < span class = "n" > Gii< / span > < span class = "ow" > in< / span > < span class = "nb" > enumerate< / span > < span class = "p" > (< / span > < span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " Gii" < / span > < span class = "p" > ]):< / span >
< span class = "n" > storage< / span > < span class = "o" > =< / span > < span class = "p" > []< / span >
< span class = "c1" > # iterate over the first and second order local perturbations< / span >
< span class = "k" > for< / span > < span class = "n" > Vu1< / span > < span class = "p" > ,< / span > < span class = "n" > Vu2< / span > < span class = "ow" > in< / span > < span class = "nb" > zip< / span > < span class = "p" > (< / span > < span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " Vu1" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ],< / span > < span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " Vu2" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ]):< / span >
< span class = "c1" > # The Szunyogh-Lichtenstein formula< / span >
< span class = "n" > traced< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > trace< / span > < span class = "p" > (< / span >
< span class = "p" > (< / span > < span class = "n" > Vu2< / span > < span class = "o" > @< / span > < span class = "n" > Gii< / span > < span class = "o" > +< / span > < span class = "mf" > 0.5< / span > < span class = "o" > *< / span > < span class = "n" > Gii< / span > < span class = "o" > @< / span > < span class = "n" > Vu1< / span > < span class = "o" > @< / span > < span class = "n" > Gii< / span > < span class = "p" > ),< / span > < span class = "n" > axis1< / span > < span class = "o" > =< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "n" > axis2< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span >
< span class = "p" > )< / span > < span class = "c1" > # this is the on site projection< / span >
< span class = "c1" > # evaluation of the contour integral< / span >
< span class = "n" > storage< / span > < span class = "o" > .< / span > < span class = "n" > append< / span > < span class = "p" > (< / span > < span class = "n" > int_de_ke< / span > < span class = "p" > (< / span > < span class = "n" > traced< / span > < span class = "p" > ,< / span > < span class = "n" > cont< / span > < span class = "o" > .< / span > < span class = "n" > we< / span > < span class = "p" > ))< / span >
< span class = "c1" > # fill up the magnetic entities dictionary with the energies< / span >
< span class = "n" > magnetic_entities< / span > < span class = "p" > [< / span > < span class = "n" > tracker< / span > < span class = "p" > ][< / span > < span class = "s2" > " energies" < / span > < span class = "p" > ]< / span > < span class = "o" > .< / span > < span class = "n" > append< / span > < span class = "p" > (< / span > < span class = "n" > storage< / span > < span class = "p" > )< / span >
< span class = "c1" > # convert to np array< / span >
< span class = "n" > magnetic_entities< / span > < span class = "p" > [< / span > < span class = "n" > tracker< / span > < span class = "p" > ][< / span > < span class = "s2" > " energies" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > (< / span >
< span class = "n" > magnetic_entities< / span > < span class = "p" > [< / span > < span class = "n" > tracker< / span > < span class = "p" > ][< / span > < span class = "s2" > " energies" < / span > < span class = "p" > ]< / span >
< span class = "p" > )< / span >
< span class = "c1" > # iterate over the pairs< / span >
< span class = "k" > for< / span > < span class = "n" > tracker< / span > < span class = "p" > ,< / span > < span class = "n" > pair< / span > < span class = "ow" > in< / span > < span class = "nb" > enumerate< / span > < span class = "p" > (< / span > < span class = "n" > pairs< / span > < span class = "p" > ):< / span >
< span class = "c1" > # iterate over the quantization axes< / span >
< span class = "k" > for< / span > < span class = "n" > i< / span > < span class = "p" > ,< / span > < span class = "p" > (< / span > < span class = "n" > Gij< / span > < span class = "p" > ,< / span > < span class = "n" > Gji< / span > < span class = "p" > )< / span > < span class = "ow" > in< / span > < span class = "nb" > enumerate< / span > < span class = "p" > (< / span > < span class = "nb" > zip< / span > < span class = "p" > (< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " Gij" < / span > < span class = "p" > ],< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " Gji" < / span > < span class = "p" > ])):< / span >
< span class = "n" > site_i< / span > < span class = "o" > =< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > [< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " ai" < / span > < span class = "p" > ]]< / span >
< span class = "n" > site_j< / span > < span class = "o" > =< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > [< / span > < span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " aj" < / span > < span class = "p" > ]]< / span >
< span class = "n" > storage< / span > < span class = "o" > =< / span > < span class = "p" > []< / span >
< span class = "c1" > # iterate over the first order local perturbations in all possible orientations for the two sites< / span >
< span class = "k" > for< / span > < span class = "n" > Vui< / span > < span class = "ow" > in< / span > < span class = "n" > site_i< / span > < span class = "p" > [< / span > < span class = "s2" > " Vu1" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ]:< / span >
< span class = "k" > for< / span > < span class = "n" > Vuj< / span > < span class = "ow" > in< / span > < span class = "n" > site_j< / span > < span class = "p" > [< / span > < span class = "s2" > " Vu1" < / span > < span class = "p" > ][< / span > < span class = "n" > i< / span > < span class = "p" > ]:< / span >
< span class = "c1" > # The Szunyogh-Lichtenstein formula< / span >
< span class = "n" > traced< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > trace< / span > < span class = "p" > (< / span >
< span class = "p" > (< / span > < span class = "n" > Vui< / span > < span class = "o" > @< / span > < span class = "n" > Gij< / span > < span class = "o" > @< / span > < span class = "n" > Vuj< / span > < span class = "o" > @< / span > < span class = "n" > Gji< / span > < span class = "p" > ),< / span > < span class = "n" > axis1< / span > < span class = "o" > =< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "n" > axis2< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span >
< span class = "p" > )< / span > < span class = "c1" > # this is the on site projection< / span >
< span class = "c1" > # evaluation of the contour integral< / span >
< span class = "n" > storage< / span > < span class = "o" > .< / span > < span class = "n" > append< / span > < span class = "p" > (< / span > < span class = "n" > int_de_ke< / span > < span class = "p" > (< / span > < span class = "n" > traced< / span > < span class = "p" > ,< / span > < span class = "n" > cont< / span > < span class = "o" > .< / span > < span class = "n" > we< / span > < span class = "p" > ))< / span >
< span class = "c1" > # fill up the pairs dictionary with the energies< / span >
< span class = "n" > pairs< / span > < span class = "p" > [< / span > < span class = "n" > tracker< / span > < span class = "p" > ][< / span > < span class = "s2" > " energies" < / span > < span class = "p" > ]< / span > < span class = "o" > .< / span > < span class = "n" > append< / span > < span class = "p" > (< / span > < span class = "n" > storage< / span > < span class = "p" > )< / span >
< span class = "c1" > # convert to np array< / span >
< span class = "n" > pairs< / span > < span class = "p" > [< / span > < span class = "n" > tracker< / span > < span class = "p" > ][< / span > < span class = "s2" > " energies" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > (< / span > < span class = "n" > pairs< / span > < span class = "p" > [< / span > < span class = "n" > tracker< / span > < span class = "p" > ][< / span > < span class = "s2" > " energies" < / span > < span class = "p" > ])< / span >
< span class = "c1" > # calculate magnetic parameters< / span >
< span class = "k" > for< / span > < span class = "n" > mag_ent< / span > < span class = "ow" > in< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > :< / span >
< span class = "n" > Kxx< / span > < span class = "p" > ,< / span > < span class = "n" > Kyy< / span > < span class = "p" > ,< / span > < span class = "n" > Kzz< / span > < span class = "p" > ,< / span > < span class = "n" > consistency< / span > < span class = "o" > =< / span > < span class = "n" > calculate_anisotropy_tensor< / span > < span class = "p" > (< / span > < span class = "n" > mag_ent< / span > < span class = "p" > )< / span >
< span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " K" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "n" > Kxx< / span > < span class = "p" > ,< / span > < span class = "n" > Kyy< / span > < span class = "p" > ,< / span > < span class = "n" > Kzz< / span > < span class = "p" > ])< / span > < span class = "o" > *< / span > < span class = "n" > sisl< / span > < span class = "o" > .< / span > < span class = "n" > unit_convert< / span > < span class = "p" > (< / span > < span class = "s2" > " eV" < / span > < span class = "p" > ,< / span > < span class = "s2" > " meV" < / span > < span class = "p" > )< / span >
< span class = "n" > mag_ent< / span > < span class = "p" > [< / span > < span class = "s2" > " K_consistency" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > consistency< / span >
< span class = "k" > for< / span > < span class = "n" > pair< / span > < span class = "ow" > in< / span > < span class = "n" > pairs< / span > < span class = "p" > :< / span >
< span class = "n" > J_iso< / span > < span class = "p" > ,< / span > < span class = "n" > J_S< / span > < span class = "p" > ,< / span > < span class = "n" > D< / span > < span class = "p" > ,< / span > < span class = "n" > J< / span > < span class = "o" > =< / span > < span class = "n" > calculate_exchange_tensor< / span > < span class = "p" > (< / span > < span class = "n" > pair< / span > < span class = "p" > )< / span >
< span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " J_iso" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > J_iso< / span > < span class = "o" > *< / span > < span class = "n" > sisl< / span > < span class = "o" > .< / span > < span class = "n" > unit_convert< / span > < span class = "p" > (< / span > < span class = "s2" > " eV" < / span > < span class = "p" > ,< / span > < span class = "s2" > " meV" < / span > < span class = "p" > )< / span >
< span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " J_S" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > J_S< / span > < span class = "o" > *< / span > < span class = "n" > sisl< / span > < span class = "o" > .< / span > < span class = "n" > unit_convert< / span > < span class = "p" > (< / span > < span class = "s2" > " eV" < / span > < span class = "p" > ,< / span > < span class = "s2" > " meV" < / span > < span class = "p" > )< / span >
< span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " D" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > D< / span > < span class = "o" > *< / span > < span class = "n" > sisl< / span > < span class = "o" > .< / span > < span class = "n" > unit_convert< / span > < span class = "p" > (< / span > < span class = "s2" > " eV" < / span > < span class = "p" > ,< / span > < span class = "s2" > " meV" < / span > < span class = "p" > )< / span >
< span class = "n" > pair< / span > < span class = "p" > [< / span > < span class = "s2" > " J" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > J< / span > < span class = "o" > *< / span > < span class = "n" > sisl< / span > < span class = "o" > .< / span > < span class = "n" > unit_convert< / span > < span class = "p" > (< / span > < span class = "s2" > " eV" < / span > < span class = "p" > ,< / span > < span class = "s2" > " meV" < / span > < span class = "p" > )< / span >
< span class = "n" > times< / span > < span class = "p" > [< / span > < span class = "s2" > " end_time" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > timer< / span > < span class = "p" > ()< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " < / span > < span class = "se" > \n\n\n\n\n< / span > < span class = "s2" > " < / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " ##################################################################### GROGU OUTPUT #############################################################################" < / span >
< span class = "p" > )< / span >
< span class = "n" > print_parameters< / span > < span class = "p" > (< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > )< / span >
< span class = "n" > print_atoms_and_pairs< / span > < span class = "p" > (< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > ,< / span > < span class = "n" > pairs< / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span >
< span class = "s2" > " ################################################################################################################################################################" < / span >
< span class = "p" > )< / span >
< span class = "n" > print_runtime_information< / span > < span class = "p" > (< / span > < span class = "n" > times< / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " " < / span > < span class = "p" > )< / span >
< span class = "c1" > # remove unwanted stuff before saving< / span >
< span class = "n" > pairs< / span > < span class = "p" > ,< / span > < span class = "n" > magnetic_entities< / span > < span class = "o" > =< / span > < span class = "n" > remove_clutter_for_save< / span > < span class = "p" > (< / span > < span class = "n" > pairs< / span > < span class = "p" > ,< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > )< / span >
< span class = "c1" > # create output dictionary with all the relevant data< / span >
< span class = "n" > results< / span > < span class = "o" > =< / span > < span class = "nb" > dict< / span > < span class = "p" > (< / span >
< span class = "n" > parameters< / span > < span class = "o" > =< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > ,< / span >
< span class = "n" > magnetic_entities< / span > < span class = "o" > =< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > ,< / span >
< span class = "n" > pairs< / span > < span class = "o" > =< / span > < span class = "n" > pairs< / span > < span class = "p" > ,< / span >
< span class = "n" > runtime< / span > < span class = "o" > =< / span > < span class = "n" > times< / span > < span class = "p" > ,< / span >
< span class = "p" > )< / span >
< span class = "c1" > # save results< / span >
< span class = "n" > save_pickle< / span > < span class = "p" > (< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " outfile" < / span > < span class = "p" > ],< / span > < span class = "n" > results< / span > < span class = "p" > )< / span >
< span class = "nb" > print< / span > < span class = "p" > (< / span > < span class = "s2" > " < / span > < span class = "se" > \n\n\n\n\n< / span > < span class = "s2" > " < / span > < span class = "p" > )< / span > < / div >
< span class = "k" > if< / span > < span class = "vm" > __name__< / span > < span class = "o" > ==< / span > < span class = "s2" > " __main__" < / span > < span class = "p" > :< / span >
< span class = "c1" > # loading parameters< / span >
< span class = "c1" > # it is not clear how to give grogu.fdf path...< / span >
< span class = "n" > command_line_arguments< / span > < span class = "o" > =< / span > < span class = "n" > parse_command_line< / span > < span class = "p" > ()< / span >
< span class = "n" > fdf_arguments< / span > < span class = "p" > ,< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > ,< / span > < span class = "n" > pairs< / span > < span class = "o" > =< / span > < span class = "n" > read_fdf< / span > < span class = "p" > (< / span > < span class = "n" > command_line_arguments< / span > < span class = "p" > [< / span > < span class = "s2" > " infile" < / span > < span class = "p" > ])< / span >
< span class = "c1" > # right now we do not use any of these input, but it shows< / span >
< span class = "c1" > # the order of priority when processing arguments< / span >
< span class = "n" > default_arguments< / span > < span class = "o" > =< / span > < span class = "kc" > False< / span >
< span class = "n" > fdf_arguments< / span > < span class = "o" > =< / span > < span class = "kc" > False< / span >
< span class = "n" > command_line_arguments< / span > < span class = "o" > =< / span > < span class = "kc" > False< / span >
< span class = "n" > simulation_parameters< / span > < span class = "o" > =< / span > < span class = "n" > process_input_args< / span > < span class = "p" > (< / span >
< span class = "n" > default_arguments< / span > < span class = "p" > ,< / span > < span class = "n" > fdf_arguments< / span > < span class = "p" > ,< / span > < span class = "n" > command_line_arguments< / span > < span class = "p" > ,< / span > < span class = "n" > ACCEPTED_INPUTS< / span >
< span class = "p" > )< / span >
< span class = "c1" > ####################################################################################################< / span >
< span class = "c1" > # This is the input file for now #< / span >
< span class = "c1" > ####################################################################################################< / span >
< span class = "n" > magnetic_entities< / span > < span class = "o" > =< / span > < span class = "p" > [< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > atom< / span > < span class = "o" > =< / span > < span class = "mi" > 3< / span > < span class = "p" > ,< / span > < span class = "n" > l< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > atom< / span > < span class = "o" > =< / span > < span class = "mi" > 4< / span > < span class = "p" > ,< / span > < span class = "n" > l< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > atom< / span > < span class = "o" > =< / span > < span class = "mi" > 5< / span > < span class = "p" > ,< / span > < span class = "n" > l< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ),< / span >
< span class = "p" > ]< / span >
< span class = "n" > pairs< / span > < span class = "o" > =< / span > < span class = "p" > [< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > ai< / span > < span class = "o" > =< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "o" > =< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "n" > Ruc< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > ai< / span > < span class = "o" > =< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "n" > Ruc< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > ai< / span > < span class = "o" > =< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "n" > Ruc< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > ai< / span > < span class = "o" > =< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "n" > Ruc< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "o" > -< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "o" > -< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > ai< / span > < span class = "o" > =< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "n" > Ruc< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "o" > -< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "o" > -< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > ai< / span > < span class = "o" > =< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "n" > Ruc< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "o" > -< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > ai< / span > < span class = "o" > =< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "n" > Ruc< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "o" > -< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > ai< / span > < span class = "o" > =< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "n" > Ruc< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "o" > -< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > ai< / span > < span class = "o" > =< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "n" > aj< / span > < span class = "o" > =< / span > < span class = "mi" > 2< / span > < span class = "p" > ,< / span > < span class = "n" > Ruc< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "o" > -< / span > < span class = "mi" > 3< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])),< / span >
< span class = "p" > ]< / span >
< span class = "n" > simulation_parameters< / span > < span class = "o" > =< / span > < span class = "nb" > dict< / span > < span class = "p" > ()< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " infile" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "p" > (< / span >
< span class = "s2" > " /Users/danielpozsar/Downloads/nojij/Fe3GeTe2/monolayer/soc/lat3_791/Fe3GeTe2.fdf" < / span >
< span class = "p" > )< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " outfile" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "s2" > " ./Fe3GeTe2_notebook" < / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " scf_xcf_orientation" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 1< / span > < span class = "p" > ])< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " ref_xcf_orientations" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "p" > [< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > o< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ]),< / span > < span class = "n" > vw< / span > < span class = "o" > =< / span > < span class = "p" > [< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ]),< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 1< / span > < span class = "p" > ])]),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > o< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ]),< / span > < span class = "n" > vw< / span > < span class = "o" > =< / span > < span class = "p" > [< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ]),< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 1< / span > < span class = "p" > ])]),< / span >
< span class = "nb" > dict< / span > < span class = "p" > (< / span > < span class = "n" > o< / span > < span class = "o" > =< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 1< / span > < span class = "p" > ]),< / span > < span class = "n" > vw< / span > < span class = "o" > =< / span > < span class = "p" > [< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ]),< / span > < span class = "n" > np< / span > < span class = "o" > .< / span > < span class = "n" > array< / span > < span class = "p" > ([< / span > < span class = "mi" > 0< / span > < span class = "p" > ,< / span > < span class = "mi" > 1< / span > < span class = "p" > ,< / span > < span class = "mi" > 0< / span > < span class = "p" > ])]),< / span >
< span class = "p" > ]< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " kset" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "mi" > 20< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " kdirs" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "s2" > " xy" < / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " ebot" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "kc" > None< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " eset" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "mi" > 600< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " esetp" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "mi" > 10000< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " parallel_solver_for_Gk" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "kc" > False< / span >
< span class = "n" > simulation_parameters< / span > < span class = "p" > [< / span > < span class = "s2" > " padawan_mode" < / span > < span class = "p" > ]< / span > < span class = "o" > =< / span > < span class = "kc" > True< / span >
< span class = "c1" > ####################################################################################################< / span >
< span class = "c1" > # This is the input file for now #< / span >
< span class = "c1" > ####################################################################################################< / span >
< span class = "n" > main< / span > < span class = "p" > (< / span > < span class = "n" > simulation_parameters< / span > < span class = "p" > ,< / span > < span class = "n" > magnetic_entities< / span > < span class = "p" > ,< / span > < span class = "n" > pairs< / span > < span class = "p" > )< / span >
< / pre > < / div >
< / div >
< / div >
< footer >
< hr / >
< div role = "contentinfo" >
< p > © Copyright 2024, grogupy.< / p >
< / div >
Built with < a href = "https://www.sphinx-doc.org/" > Sphinx< / a > using a
< a href = "https://github.com/readthedocs/sphinx_rtd_theme" > theme< / a >
provided by < a href = "https://readthedocs.org" > Read the Docs< / a > .
< / footer >
< / div >
< / div >
< / section >
< / div >
< script >
jQuery(function () {
SphinxRtdTheme.Navigation.enable(true);
});
< / script >
< / body >
< / html >