TABLE OF CONTENTS
- ABINIT/m_spin_model
- m_spin_model/spin_model_t
- m_spin_model/spin_model_t_finalize
- m_spin_model/spin_model_t_initialize
- m_spin_model/spin_model_t_make_supercell
- m_spin_model/spin_model_t_read_xml
- m_spin_model/spin_model_t_run
- m_spin_model/spin_model_t_run_MvH
- m_spin_model/spin_model_t_run_MvT
- m_spin_model/spin_model_t_run_one_step
- m_spin_model/spin_model_t_run_time
- m_spin_model/spin_model_t_set_initial_spin
- m_spin_model/spin_model_t_set_params
ABINIT/m_spin_model [ Modules ]
NAME
m_spin_model
FUNCTION
This module contains the manager of spin dynamics. It call other related modules to run spin dynamics. Itself does nothing in detail. Datatypes: * spin_model_t Subroutines: * spin_model_t * spin_model_t_run * spin_model_t_initialize * spin_model_t_set_initial_spin * spin_model_t_finalize * spin_model_t_read_xml * spin_model_t_make_supercell * spin_model_t_run_one_step * spin_model_t_run_time
COPYRIGHT
Copyright (C) 2001-2017 ABINIT group (hexu) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
SOURCE
38 #if defined HAVE_CONFIG_H 39 #include "config.h" 40 #endif 41 42 #include "abi_common.h" 43 44 module m_spin_model 45 46 use defs_basis 47 use m_abicore 48 use m_errors 49 use m_xmpi 50 51 use m_multibinit_dataset, only: multibinit_dtset_type 52 use m_spin_terms, only: spin_terms_t, spin_terms_t_finalize, spin_terms_t_set_external_hfield 53 use m_spin_model_primitive, only: spin_model_primitive_t, & 54 & spin_model_primitive_t_initialize, & 55 & spin_model_primitive_t_print_terms, & 56 & spin_model_primitive_t_finalize, & 57 & spin_model_primitive_t_read_xml, & 58 & spin_model_primitive_t_make_supercell 59 use m_spin_hist, only: spin_hist_t, spin_hist_t_set_vars, spin_hist_t_init, spin_hist_t_get_s, spin_hist_t_free, & 60 & spin_hist_t_set_params 61 use m_spin_mover, only: spin_mover_t, spin_mover_t_initialize, spin_mover_t_finalize, & 62 & spin_mover_t_run_time, spin_mover_t_run_one_step 63 use m_spin_ncfile, only: spin_ncfile_t, spin_ncfile_t_init, spin_ncfile_t_close, spin_ncfile_t_def_sd, & 64 & spin_ncfile_t_write_primitive_cell, spin_ncfile_t_write_supercell, spin_ncfile_t_write_parameters, & 65 & spin_ncfile_t_write_one_step 66 implicit none
m_spin_model/spin_model_t [ Types ]
[ Top ] [ m_spin_model ] [ Types ]
NAME
spin_model_t
FUNCTION
This type contains all the data for spin dynamics. It contains: spin_primitive = information of the primitive cell, which is a map to the xml file. spin_calculator = the calculator for spin dynamics, include information of supercell. spin_mover = include information for how to run spin dynamics params = parameters from input file spin_ncfile = wrapper for spin hist netcdf file. nmatom= number of magnetic atoms
SOURCE
86 type spin_model_t 87 type(spin_model_primitive_t) :: spin_primitive 88 type(spin_terms_t) :: spin_calculator 89 type(spin_hist_t):: spin_hist 90 type(spin_mover_t):: spin_mover 91 type(multibinit_dtset_type) :: params 92 type(spin_ncfile_t) :: spin_ncfile 93 integer :: nmatoms 94 ! CONTAINS 95 ! procedure :: run => spin_model_t_run 96 ! procedure :: initialize=>spin_model_t_initialize 97 ! procedure :: set_initial_spin => spin_model_t_set_initial_spin 98 ! procedure :: finalize => spin_model_t_finalize 99 ! procedure :: read_xml => spin_model_t_read_xml 100 ! procedure :: make_supercell => spin_model_t_make_supercell 101 ! procedure :: run_one_step => spin_model_t_run_one_step 102 ! procedure :: run_time => spin_model_t_run_time 103 ! procedure :: run_MvT => spin_model_t_run_MvT 104 ! procedure :: run_MvH => spin_model_t_run_MvH 105 end type spin_model_t
m_spin_model/spin_model_t_finalize [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_finalize
FUNCTION
clean memory for spin dynamics
INPUTS
OUTPUT
NOTES
PARENTS
CHILDREN
SOURCE
245 subroutine spin_model_t_finalize(self) 246 247 248 !This section has been created automatically by the script Abilint (TD). 249 !Do not modify the following lines by hand. 250 #undef ABI_FUNC 251 #define ABI_FUNC 'spin_model_t_finalize' 252 !End of the abilint section 253 254 class(spin_model_t), intent(inout) :: self 255 !call self%spin_primitive%finalize() 256 call spin_model_primitive_t_finalize(self%spin_primitive) 257 call spin_terms_t_finalize(self%spin_calculator) 258 call spin_hist_t_free(self%spin_hist) 259 call spin_mover_t_finalize(self%spin_mover) 260 call spin_ncfile_t_close(self%spin_ncfile) 261 ! TODO: finalize others 262 end subroutine spin_model_t_finalize
m_spin_model/spin_model_t_initialize [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_initialize
FUNCTION
initialize spin dynamics from input
INPUTS
xml_fname= file of xml file params = multibinit dataset from input file
OUTPUT
self <spin_model_t>
NOTES
PARENTS
CHILDREN
SOURCE
165 subroutine spin_model_t_initialize(self, xml_fname, params) 166 167 168 !This section has been created automatically by the script Abilint (TD). 169 !Do not modify the following lines by hand. 170 #undef ABI_FUNC 171 #define ABI_FUNC 'spin_model_t_initialize' 172 !End of the abilint section 173 174 class(spin_model_t), intent(inout) :: self 175 character(len=*), intent(in) :: xml_fname 176 integer :: sc_matrix(3,3) 177 type(multibinit_dtset_type), intent(in) :: params 178 179 ! read input 180 !call self%spin_primitive%initialize() 181 call spin_model_primitive_t_initialize(self%spin_primitive) 182 self%params=params 183 !call self%read_xml(xml_fname) 184 call spin_model_t_read_xml(self, xml_fname) 185 !call self%spin_primitive%print_terms() 186 call spin_model_primitive_t_print_terms(self%spin_primitive) 187 188 ! make supercell 189 sc_matrix(:,:)=0.0_dp 190 sc_matrix(1,1)=params%ncell(1) 191 sc_matrix(2,2)=params%ncell(2) 192 sc_matrix(3,3)=params%ncell(3) 193 !call self%make_supercell(sc_matrix) 194 call spin_model_t_make_supercell(self, sc_matrix) 195 196 ! set parameters to hamiltonian and mover 197 self%nmatoms= self%spin_calculator%nmatoms 198 199 call spin_model_t_set_params(self) 200 201 202 !TODO hexu: mxhist, has_latt, natom should be input with their true values when lattice part also added 203 call spin_hist_t_init(hist=self%spin_hist, nmatom=self%nmatoms, mxhist=3, has_latt=.False.) 204 call spin_hist_t_set_params(self%spin_hist, spin_nctime=self%params%spin_nctime, & 205 & spin_temperature=self%params%spin_temperature) 206 !TODO 207 !call spin_hist_t_set_atomic_structure(self%spin_hist, acell, rprimd, xred, spin_index, ntypat, typat, znucl) 208 209 !call self%set_initial_spin(mode=1) 210 call spin_model_t_set_initial_spin(self, mode=0) 211 212 !call self%spin_mover%initialize(self%nmatoms, dt=params%dtspin, total_time=params%dtspin*params%ntime_spin, temperature=self%params%self) 213 call spin_mover_t_initialize(self%spin_mover, self%nmatoms, dt=params%spin_dt, & 214 & total_time=params%spin_dt*params%spin_ntime, temperature=self%params%spin_temperature) 215 216 call spin_ncfile_t_init(self%spin_ncfile, 'spinhist.nc') 217 call spin_ncfile_t_def_sd(self%spin_ncfile, self%spin_hist ) 218 !call spin_ncfile_t_write_primitive_cell(self%spin_ncfile, self%spin_primitive) 219 call spin_ncfile_t_write_supercell(self%spin_ncfile, self%spin_calculator) 220 call spin_ncfile_t_write_parameters(self%spin_ncfile, self%params) 221 222 call spin_ncfile_t_write_one_step(self%spin_ncfile, self%spin_hist) 223 end subroutine spin_model_t_initialize
m_spin_model/spin_model_t_make_supercell [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_make_supercell
FUNCTION
make supercell and prepare calculator
INPUTS
sc_mat(3,3) = supercell matrix
OUTPUT
NOTES
PARENTS
CHILDREN
SOURCE
362 subroutine spin_model_t_make_supercell(self, sc_mat) 363 364 365 !This section has been created automatically by the script Abilint (TD). 366 !Do not modify the following lines by hand. 367 #undef ABI_FUNC 368 #define ABI_FUNC 'spin_model_t_make_supercell' 369 !End of the abilint section 370 371 class(spin_model_t), intent(inout) :: self 372 integer , intent(in):: sc_mat(3, 3) 373 !call self%spin_primitive%make_supercell(sc_mat, self%spin_calculator) 374 call spin_model_primitive_t_make_supercell(self%spin_primitive, sc_mat, self%spin_calculator) 375 end subroutine spin_model_t_make_supercell
m_spin_model/spin_model_t_read_xml [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_read_xml
FUNCTION
read xml file and set primitive cell parameters.
INPUTS
OUTPUT
NOTES
PARENTS
CHILDREN
SOURCE
326 subroutine spin_model_t_read_xml(self, xml_fname) 327 328 329 !This section has been created automatically by the script Abilint (TD). 330 !Do not modify the following lines by hand. 331 #undef ABI_FUNC 332 #define ABI_FUNC 'spin_model_t_read_xml' 333 !End of the abilint section 334 335 class(spin_model_t), intent(inout) :: self 336 character(len=*), intent(in) :: xml_fname 337 !call self%spin_primitive%read_xml(xml_fname) 338 call spin_model_primitive_t_read_xml(self%spin_primitive, xml_fname) 339 end subroutine spin_model_t_read_xml
m_spin_model/spin_model_t_run [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_run
FUNCTION
run the whole spin dynamics.
INPUTS
OUTPUT
NOTES
currently it just call run_time. It will be able to call e.g. run_MvT to do other things.
PARENTS
CHILDREN
SOURCE
129 subroutine spin_model_t_run(self) 130 131 132 !This section has been created automatically by the script Abilint (TD). 133 !Do not modify the following lines by hand. 134 #undef ABI_FUNC 135 #define ABI_FUNC 'spin_model_t_run' 136 !End of the abilint section 137 138 class(spin_model_t), intent(inout) :: self 139 integer :: i 140 !call spin_mover_t_run_time(self%spin_mover, self%spin_calculator, self%spin_hist) 141 call spin_model_t_run_time(self) 142 end subroutine spin_model_t_run
m_spin_model/spin_model_t_run_MvH [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_run_MvT
FUNCTION
run spin vs external magnetic field
INPUTS
OUTPUT
NOTES
Not yet implemented.
PARENTS
CHILDREN
SOURCE
564 subroutine spin_model_t_run_MvH(self) 565 566 567 !This section has been created automatically by the script Abilint (TD). 568 !Do not modify the following lines by hand. 569 #undef ABI_FUNC 570 #define ABI_FUNC 'spin_model_t_run_MvH' 571 !End of the abilint section 572 573 class(spin_model_t), intent(inout) :: self 574 write(std_err, *) "MvH is not yet implemented. " 575 end subroutine spin_model_t_run_MvH
m_spin_model/spin_model_t_run_MvT [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_run_MvT
FUNCTION
run spin vs temperature
INPUTS
OUTPUT
NOTES
Not yet implemented.
PARENTS
CHILDREN
SOURCE
530 subroutine spin_model_t_run_MvT(self) 531 532 533 !This section has been created automatically by the script Abilint (TD). 534 !Do not modify the following lines by hand. 535 #undef ABI_FUNC 536 #define ABI_FUNC 'spin_model_t_run_MvT' 537 !End of the abilint section 538 539 class(spin_model_t), intent(inout) :: self 540 !TODO 541 write(std_err, *) "MvH is not yet implemented. " 542 end subroutine spin_model_t_run_MvT
m_spin_model/spin_model_t_run_one_step [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_run_one_step
FUNCTION
run one spin dynamics step
INPUTS
OUTPUT
NOTES
This is not called by run_time. insted run_time call mover%run_time (in the spirit that this module should do nothing in detail!)
PARENTS
CHILDREN
SOURCE
460 subroutine spin_model_t_run_one_step(self) 461 462 463 !This section has been created automatically by the script Abilint (TD). 464 !Do not modify the following lines by hand. 465 #undef ABI_FUNC 466 #define ABI_FUNC 'spin_model_t_run_one_step' 467 !End of the abilint section 468 469 class(spin_model_t), intent(inout) :: self 470 real(dp) :: S_tmp(3,self%nmatoms), etot 471 call spin_mover_t_run_one_step(self%spin_mover, self%spin_calculator, & 472 spin_hist_t_get_S(self%spin_hist),S_tmp, etot) 473 call spin_hist_t_set_vars(self%spin_hist, S=S_tmp, etot=etot, inc=.False.) 474 end subroutine spin_model_t_run_one_step
m_spin_model/spin_model_t_run_time [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_run_time
FUNCTION
run all spin time step
INPUTS
OUTPUT
NOTES
PARENTS
CHILDREN
SOURCE
496 subroutine spin_model_t_run_time(self) 497 498 499 !This section has been created automatically by the script Abilint (TD). 500 !Do not modify the following lines by hand. 501 #undef ABI_FUNC 502 #define ABI_FUNC 'spin_model_t_run_time' 503 !End of the abilint section 504 505 class(spin_model_t), intent(inout) :: self 506 call spin_mover_t_run_time(self%spin_mover,self%spin_calculator,self%spin_hist, self%spin_ncfile) 507 end subroutine spin_model_t_run_time
m_spin_model/spin_model_t_set_initial_spin [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_set_initial_spin
FUNCTION
set initial spin state
INPUTS
mode= 0 : all along z. 1: random
OUTPUT
NOTES
PARENTS
CHILDREN
SOURCE
397 subroutine spin_model_t_set_initial_spin(self, mode) 398 399 400 !This section has been created automatically by the script Abilint (TD). 401 !Do not modify the following lines by hand. 402 #undef ABI_FUNC 403 #define ABI_FUNC 'spin_model_t_set_initial_spin' 404 !End of the abilint section 405 406 class(spin_model_t), intent(inout) :: self 407 integer, intent(in) :: mode 408 integer :: i 409 real(dp) :: S(3, self%nmatoms) 410 character(len=500) :: msg 411 if(mode==0) then 412 ! set all spin to z direction. 413 S(1,:)=0.0d0 414 S(2,:)=0.0d0 415 S(3,:)=1.0d0 416 else if (mode==1) then 417 ! randomize S using uniform random number 418 ! print *, "Initial spin set to random value" 419 write(msg,*) "Initial spin set to random value." 420 call wrtout(ab_out,msg,'COLL') 421 call wrtout(std_out,msg,'COLL') 422 call random_number(S) 423 S=S-0.5 424 do i=1, self%nmatoms 425 S(:,i)=S(:,i)/sqrt(sum(S(:, i)**2)) 426 end do 427 else 428 write(msg,*) "Error: Set initial spin: mode should be 0 (FM) or 1 (random)" 429 call wrtout(ab_out,msg,'COLL') 430 call wrtout(std_out,msg,'COLL') 431 432 end if 433 434 call spin_hist_t_set_vars(self%spin_hist, S=S, time=0.0_dp, ihist_latt=0, inc=.True.) 435 !print *, "initial spin", self%spin_hist%current_S 436 end subroutine spin_model_t_set_initial_spin
m_spin_model/spin_model_t_set_params [ Functions ]
[ Top ] [ m_spin_model ] [ Functions ]
NAME
spin_model_t_set_params
FUNCTION
set parameters for spin dynamics from input.
INPUTS
OUTPUT
NOTES
Most are already done in initialize. (should be moved to here or remove this function)
PARENTS
CHILDREN
SOURCE
284 subroutine spin_model_t_set_params(self) 285 286 287 !This section has been created automatically by the script Abilint (TD). 288 !Do not modify the following lines by hand. 289 #undef ABI_FUNC 290 #define ABI_FUNC 'spin_model_t_set_params' 291 !End of the abilint section 292 293 class(spin_model_t), intent(inout) :: self 294 real(dp):: mfield(3, self%nmatoms) 295 integer :: i 296 ! params -> mover 297 ! params -> calculator 298 do i=1, self%nmatoms 299 mfield(:, i)=self%params%spin_mag_field(:) 300 end do 301 call spin_terms_t_set_external_hfield(self%spin_calculator, mfield) 302 ! params -> hist 303 304 end subroutine spin_model_t_set_params