From bfe8b97c7645d6b8078ac63b7570011225a88be4 Mon Sep 17 00:00:00 2001 From: Nicolas Essing <n.essing@fz-juelich.de> Date: Thu, 28 Oct 2021 11:17:37 +0200 Subject: [PATCH] Cleaned up input parameters for nonco bfields. Remove some inputs that were moved to atom-wise switches in nonco_angle.dat and some for functionality that was not implemented. Changed the rest to form a more consistent interface. Keeping a switch to turn everything off. If that switch is turned on, this implies calculating the torque and constraint magnetism in general. The fields are calculated for all atoms where that is specified in the nonco_angle.dat. It does not anymore imply external fields: There is a new switch for that. If that is turned on, a file with external fields must be provided. Introduced a verbosity parameter for bfields, select between no output, summary and detailed information. See comments and documentation for details. Changed output depending on the new input parameter. --- source/KKRnano/source/CalculationData_mod.F90 | 9 +- .../source/NonCollinearMagnetism_mod.F90 | 22 ++- .../source/ScatteringCalculation_mod.F90 | 2 +- source/KKRnano/source/bfield/bfield.f90 | 142 ++++++++++++------ source/KKRnano/source/bfield/torque.f90 | 5 +- .../source/datastructures/InputParamsNew.txt | 16 +- .../source/datastructures/InputParams_mod.F90 | 58 ++----- source/KKRnano/source/wrappers_mod.F90 | 2 +- 8 files changed, 134 insertions(+), 122 deletions(-) diff --git a/source/KKRnano/source/CalculationData_mod.F90 b/source/KKRnano/source/CalculationData_mod.F90 index 7c2e39df9..0ff9ee443 100644 --- a/source/KKRnano/source/CalculationData_mod.F90 +++ b/source/KKRnano/source/CalculationData_mod.F90 @@ -295,6 +295,7 @@ module CalculationData_mod integer, intent(in) :: voronano integer :: atom_id, ila, irmd + integer :: verbosity call create(self%lattice_vectors, arrays%bravais) ! createLatticeVectors @@ -385,8 +386,14 @@ module CalculationData_mod call setup_iguess(self, dims, arrays%nofks, kmesh) ! setup storage for iguess if (params%noncobfield) then + ! Output only as master, in that case copy the input parameter + verbosity = -1 + if (mp%isMasterRank) verbosity = params%bfield_verbosity + ! Initialize the noncolinear magnetic field. If present, read from disk - call load_bfields_from_disk(self%bfields, params%constr_field, dims%naez, self%atom_ids) + call load_bfields_from_disk(self%bfields, params%external_bfield, & + verbosity, dims%naez, self%atom_ids, & + self%noco_data%angle_fix_mode) ! Initialize the fields do ila = 1, self%num_local_atoms call init_bfield(self%bfields(ila), dims%lmaxd, & diff --git a/source/KKRnano/source/NonCollinearMagnetism_mod.F90 b/source/KKRnano/source/NonCollinearMagnetism_mod.F90 index cf586ae2e..6b5da416e 100644 --- a/source/KKRnano/source/NonCollinearMagnetism_mod.F90 +++ b/source/KKRnano/source/NonCollinearMagnetism_mod.F90 @@ -33,7 +33,7 @@ SUBROUTINE tmat_newsolver(ie,nspin,lmax,zat,socscale, & ! lly, & lmpotd,irmd_new,TmatN,soc,enable_quad_prec, & bfield, imt, iteration_number, itscf0, itscf1, & - lbfield, lbfield_constr, lbfield_trans, lbfield_mt) + lbfield, lbfield_trans, lbfield_mt) ! Code converted using TO_F90 by Alan Miller ! Date: 2016-04-18 Time: 14:58:02 @@ -68,7 +68,7 @@ type(bfield_data), intent(in) :: bfield !! Information on the noncoll integer, intent(in) :: imt !! MT radius (index in cheb mesh) integer, intent(in) :: iteration_number !! Current iteration number integer, intent(in) :: itscf0, itscf1 !! Window of iterations when to apply fields -logical, intent(in) :: lbfield, lbfield_constr !! Wheter to apply the fields +logical, intent(in) :: lbfield !! Wheter to apply nonco bfields logical, intent(in) :: lbfield_trans ! Apply only transveral logical, intent(in) :: lbfield_mt ! Apply only up do MT radius @@ -167,8 +167,7 @@ allocate(ull(nsra*lmmaxso,lmmaxso,irmdnew)) ! add bfield if (lbfield) then call add_bfield(bfield, vnspll1, theta, phi, imt, iteration_number, & - itscf0, itscf1, lbfield_constr, lbfield_trans, & - lbfield_mt, .false.) + itscf0, itscf1, lbfield_trans, lbfield_mt, .false.) end if !c extend matrix for the SRA treatment @@ -339,7 +338,7 @@ SUBROUTINE rhovalnew(ldorhoef,ielast,nsra,nspin,lmax,ez,wez,zat, & den_out,espv,rho2ns,r2nef,gmatn, muorb, & lpotd,lmaxd,irmd,irmd_new,iemxd,soc,enable_quad_prec, & bfield, imt, iteration_number, itscf0, itscf1, & - lbfield, lbfield_constr, ltorque, lbfield_trans, lbfield_mt) + lbfield, lbfield_trans, lbfield_mt) ! Code converted using TO_F90 by Alan Miller ! Date: 2016-04-21 Time: 11:39:57 @@ -395,8 +394,7 @@ type(bfield_data), intent(inout) :: bfield !! Information on the noncoll integer, intent(in) :: imt !! MT radius (index in cheb mesh) integer, intent(in) :: iteration_number !! Current iteration number integer, intent(in) :: itscf0, itscf1 !! Window of iterations when to apply fields -logical, intent(in) :: lbfield, lbfield_constr !! Wheter to apply the fields -logical, intent(in) :: ltorque !! Whether to calculate torques +logical, intent(in) :: lbfield !! Wheter to apply nonco bfields logical, intent(in) :: lbfield_trans ! Apply only transveral logical, intent(in) :: lbfield_mt ! Apply only up do MT radius @@ -559,8 +557,7 @@ DO ie=1,ielast ! Add bfield if (lbfield) then call add_bfield(bfield, vnspll1, theta, phi, imt, iteration_number, & - itscf0, itscf1, lbfield_constr, lbfield_trans, & - lbfield_mt, .false.) + itscf0, itscf1, lbfield_trans, lbfield_mt, .false.) end if !c extend matrix for the SRA treatment @@ -625,8 +622,7 @@ DO ie=1,ielast ! Add bfield if (lbfield) then call add_bfield(bfield, vnspll1, theta, phi, imt, iteration_number, & - itscf0, itscf1, lbfield_constr, lbfield_trans, & - lbfield_mt, .true.) + itscf0, itscf1, lbfield_trans, lbfield_mt, .true.) end if !c extend matrix for the SRA treatment @@ -823,10 +819,10 @@ END DO ! Calculate torque ! Compare this position to KKhost: calling after the rho2nsc_loop, before cheb2oldgrid -if (ltorque) then +if (lbfield) then call calc_torque(bfield, vins, rho2nsc, theta, phi, lmax, rpan_intervall, ipan_intervall, & npan_tot, ncheb, imt1, iend, icleb, cleb, ifunm, thetasnew, & - lbfield_constr, lbfield_mt, itscf0, itscf1, iteration_number, angle_fix_mode) + lbfield_mt, itscf0, itscf1, iteration_number, angle_fix_mode) end if allocate(rhotemp(irmdnew,lmpotd)) diff --git a/source/KKRnano/source/ScatteringCalculation_mod.F90 b/source/KKRnano/source/ScatteringCalculation_mod.F90 index 7e33e5c34..31a4b6c8e 100644 --- a/source/KKRnano/source/ScatteringCalculation_mod.F90 +++ b/source/KKRnano/source/ScatteringCalculation_mod.F90 @@ -274,7 +274,7 @@ implicit none kkr(ila)%TmatN(:,:,ispin),params%soc,params%enable_quad_prec, & calc%bfields(ila), get_muffin_tin_index(atomdata%chebmesh_ptr), & iter, params%itbfield0, params%itbfield1, & - params%noncobfield, params%constr_field, params%trans_bfield, params%mt_bfield) + params%noncobfield, params%trans_bfield, params%mt_bfield) call rotatematrix(kkr(ila)%TmatN(:,:,ispin),noco%theta_noco(i1),noco%phi_noco(i1),lmmaxd,0) else diff --git a/source/KKRnano/source/bfield/bfield.f90 b/source/KKRnano/source/bfield/bfield.f90 index 3dc12d5b5..4dd30e397 100644 --- a/source/KKRnano/source/bfield/bfield.f90 +++ b/source/KKRnano/source/bfield/bfield.f90 @@ -38,16 +38,18 @@ module mod_bfield contains - !> Load the external noncollinear magnetic field (if a file is present) and - !> the initial guess for the constraining fields (if constraint magnetism is - !> used and a file is present) from disk. - !> This subroutine loads information for the locally treated atoms into an - !> array of structs. - subroutine load_bfields_from_disk(bfields, lbfield_constr, naez, atom_ids) + !> Load the external noncollinear magnetic field (if used) and the initial + !> guess for the constraining fields from disk. + !> This subroutine loads information for the locally treated atoms into a + !> given array of structs of the right size. + subroutine load_bfields_from_disk(bfields, external_bfield, verbosity, naez, & + atom_ids, fix_angle_modes) type(bfield_data), intent(inout) :: bfields(:) !! Array for the read fields - logical, intent(in) :: lbfield_constr !! Wheter constraint magnetism is used + logical, intent(in) :: external_bfield !! Whether to apply external nonco fields + integer, intent(in) :: verbosity !! Output verbosity integer, intent(in) :: naez !! Number of atoms in the unit cell integer(kind=4), dimension(:), intent(in) :: atom_ids !! Locally treated atoms + integer(kind=1), dimension(:), intent(in) :: fix_angle_modes integer :: num_local_atoms ! number of atoms treated in this MPI rank integer :: ila ! atom index @@ -61,10 +63,8 @@ contains allocate(all_bfields(naez)) ! Read the files - call read_bfield(all_bfields) - if (lbfield_constr) then - call read_bconstr(all_bfields) - end if + call read_external_bfield(all_bfields, external_bfield, verbosity) + call read_constr_bfield_initial_guess(all_bfields, verbosity, fix_angle_modes) ! Copy the locally treated atoms to the output do ila = 1, num_local_atoms @@ -105,18 +105,19 @@ contains end subroutine init_bfield - !TODO save_bconst() - !------------------------------------------------------------------------------- - !> Summary: Reads the atom-wise constraining field from bconstr.dat + !> Summary: Reads the atom-wise initial guess for the constraining magnetic + !> field from bconstr.dat !> Author: MdSD, Nicolas Essing !> !> The file should contain three floating point values in each line and a line !> for each atom. The first line is treated as a header and skipped. !> Intepreted as initial constraining bfield in cartesian coordinates, in Ry. !------------------------------------------------------------------------------- - subroutine read_bconstr(bfields) + subroutine read_constr_bfield_initial_guess(bfields, verbosity, fix_angle_modes) type(bfield_data), dimension(:), intent(inout) :: bfields + integer, intent(in) :: verbosity + integer(kind=1), dimension(:), intent(in) :: fix_angle_modes integer :: number_of_atoms integer :: iatom @@ -136,6 +137,19 @@ contains write(*,*) "Error reading bconstr_in.dat" stop end if + + ! In case constraint magnetism is not used for this atom, set the initial + ! guess to zero. Will not be updated in that case and can be added to the + ! potential without further checking on input parameters. + ! If an initial guess was provided that is not zero, give a warning. + if (.not. (fix_angle_modes(iatom) == 2 .or. fix_angle_modes(iatom) == 3)) then + if (any(bfields(iatom)%bfield_constr /= 0) .and. verbosity >= 0) then + write(*,'(2A,I3,2A)') 'Warning: Initial guess for constraint magnetic field ', & + 'for atom ', iatom, ' was not zero, but no constraint magnetism is ', & + 'used for this atom. Will be set to zero.' + end if + bfields(iatom)%bfield_constr(:) = 0 + end if end do close(57493215) else @@ -145,43 +159,65 @@ contains end do end if - write(*,*) ' ############################################################' - write(*,*) ' input constraining fields' - write(*,*) ' ############################################################' - write(*,*) ' iatom Bx By Bz (in Ry)' - do iatom = 1, number_of_atoms - write(*,'(2X,I4,3(E16.8))') iatom, bfields(iatom)%bfield_constr(:) - end do - end subroutine read_bconstr + if (verbosity >= 3) then + ! Write detailed information + write(*,'(79("#"))') + write(*,'(2X,A)') 'Initial guess for constraining fields' + write(*,'(79("#"))') + write(*,'(2X,A4,3(7X,A7,3X))') 'atom','Bx [Ry]','By [Ry]','Bz [Ry]' + do iatom = 1, number_of_atoms + write(*,'(2X,I4,3(2X,E15.8))') iatom, bfields(iatom)%bfield_constr(:) + end do + write(*,'(79("#"))') + else if (verbosity >= 2) then + ! Give an overview + if (file_exists) then + write(*,'(A)') 'Initial constraining magnetic fields read from file.' + else + write(*,'(A)') 'Initial constraining magnetic fields set to zero.' + end if + else + ! No output + end if + end subroutine !------------------------------------------------------------------------------- - !> Summary: Reads the atom-wise magnetic field from bfield.dat + !> Summary: Reads the atom-wise external magnetic field from bfield.dat, + !> if it should be used according to the input parameters. !> Author: Nicolas Essing !> + !> If external fields shall not be used, the external field is set to zero and + !> an eventually present file is ignored. + !> If external fields shall be used and no file is present, the program stops + !> with an error. !> The file should contain three numbers per line and one line per atom. !> Read as 'theta phi bfield_strength' with the angles in degrees and the - !> strength in Ry ( ! might get changed to Tesla ! ) + !> strength in Ry. !> Lines can be commented out with a # as first character. !>------------------------------------------------------------------------------ - subroutine read_bfield(bfields) + subroutine read_external_bfield(bfields, external_bfield, verbosity) type(bfield_data), dimension(:), intent(inout) :: bfields + logical, intent(in) :: external_bfield + integer, intent(in) :: verbosity integer :: number_of_atoms integer :: iatom, iostat character(256) :: linebuffer - logical :: file_exists double precision, dimension(:), allocatable :: phi, theta, strength number_of_atoms = size(bfields) - inquire(file='bfield.dat', exist=file_exists) - allocate(phi(number_of_atoms), theta(number_of_atoms), strength(number_of_atoms)) - if (file_exists) then + if (external_bfield) then open(unit=57493215, file='bfield.dat', iostat=iostat) - + if (iostat /= 0) then + write(*,*) "I/O Error opening bfield.dat" + write(*,*) "If you do not want external nonco bfields, turn it of in the input.conf" + stop + end if + ! manual loop over iatom because comments are allowed iatom = 1 do while (iatom <= number_of_atoms) @@ -214,14 +250,28 @@ contains bfields(iatom)%bfield_ext(3) = strength(iatom) * cos(theta(iatom)) end do - write(*,*) ' ###############################################' - write(*,*) ' external non-collinear magnetic fields' - write(*,*) ' ###############################################' - write(*,*) ' iatom theta phi bfield (in Ry)' - do iatom = 1, number_of_atoms - write(*,'(2X,I4,3(E16.8))') iatom, theta(iatom), phi(iatom), strength(iatom) - end do - end subroutine read_bfield + if (verbosity >= 3) then + ! Write detailed information + write(*,'(79("#"))') + write(*,'(16X,A)') 'external non-collinear magnetic fields' + write(*,'(79("#"))') + write(*,'(2X,A4,3(5X,A11,1X))') 'atom', ' theta ', ' phi ', 'bfield [Ry]' + do iatom = 1, number_of_atoms + write(*,'(2X,I4,3(2X,E15.8))') iatom, theta(iatom), phi(iatom), strength(iatom) + end do + write(*,'(79("#"))') + else if (verbosity >= 2) then + ! Give an overview + if (external_bfield) then + write(*,'(A)') 'External non-collinear magnetic fields loaded from file.' + write(*,'(2X,A,E16.8)') 'Mean magnitude [Ry]:', sum(abs(strength)) / number_of_atoms + else + write(*,'(A)') 'Not applying external non-collinear magnetic fields.' + end if + else + ! No output + end if + end subroutine !>------------------------------------------------------------------------------ @@ -234,7 +284,7 @@ contains !> pauli matrices and B the combined bfield. !>------------------------------------------------------------------------------ subroutine add_bfield(bfield, vnspll, theta, phi, imt, iteration_number, & - itscf0, itscf1, lbfield_constr, lbfield_trans, & + itscf0, itscf1, lbfield_trans, & lbfield_mt, transpose_bfield) type(bfield_data), intent(in) :: bfield double complex, dimension(:,:,:), intent(inout) :: vnspll ! The potential to add to @@ -242,7 +292,6 @@ contains integer, intent(in) :: imt ! MT radius (index in cheb mesh) integer, intent(in) :: iteration_number !TODO this, or just a logical and do the check outside? integer, intent(in) :: itscf0, itscf1 !TODO ^ - logical, intent(in) :: lbfield_constr !! Wheter to use constraint fields logical, intent(in) :: lbfield_trans ! Apply only transveral logical, intent(in) :: lbfield_mt ! Apply only up do MT radius logical, intent(in) :: transpose_bfield ! Transpose the bfield (for left solutions) @@ -262,11 +311,10 @@ contains lmmax = size(bfield%thetallmat, 1) ! size(vnspll, 1) is 2*lmmax irmd = size(vnspll, 3) - combined_bfields(:) = bfield%bfield_ext(:) ! start with external, is zero if unused - if (lbfield_constr) then - ! Add constraint field - combined_bfields(:) = combined_bfields(:) + bfield%bfield_constr(:) - end if + ! Add external and constraint field. If one of them is tured off by input + ! parameters or mode, it is zero, so this distinction does not have to be + ! done here. + combined_bfields(:) = bfield%bfield_ext(:) + bfield%bfield_constr(:) if (lbfield_trans) then dir(1) = sin(theta) * cos(phi) @@ -312,8 +360,6 @@ contains !------------------------------------------------------------------------------------ !> Summary: Shape function LL' expansion !> Author: Sascha Brinker - !> Category: KKRhost, geometry, new-mesh, shapefun - !> Deprecated: False !> Calculates the LL' expansion of the shape function similarly to vllmat_new.f90 !> @note The input shapefunction (single L) uses pointer arrays for the lm index. !> The output does not need pointers! diff --git a/source/KKRnano/source/bfield/torque.f90 b/source/KKRnano/source/bfield/torque.f90 index 7a3a704f3..485515e17 100644 --- a/source/KKRnano/source/bfield/torque.f90 +++ b/source/KKRnano/source/bfield/torque.f90 @@ -23,7 +23,7 @@ contains !------------------------------------------------------------------------------- subroutine calc_torque(bfield, vpot, rho2nsc, theta, phi, lmax, rpan_intervall, ipan_intervall, & npan_tot, ncheb, imt, iend, icleb, cleb, ifunm, thetasnew, & - lbfield_constr, lbfield_mt, itscf0, itscf1, iteration, constr_mode) + lbfield_mt, itscf0, itscf1, iteration, constr_mode) type(bfield_data), intent(inout) :: bfield !! Information on the magnetic field double precision, dimension(:,:,:), intent(in) :: vpot !! The potential double complex, dimension(:,:,:), intent(in) :: rho2nsc !! complex density matrix @@ -39,7 +39,6 @@ contains double precision, dimension(:), intent(in) :: cleb !! Gaunt coefficients integer , dimension(:), intent(in) :: ifunm !! pointer array for shapefun double precision, dimension(:, :), intent(in) :: thetasnew !! shapefun on the Cheby mesh - logical, intent(in) :: lbfield_constr !! Use constraint fields logical, intent(in) :: lbfield_mt !! Use magnetic fields only inside the muffin tin integer, intent(in) :: itscf0, itscf1 !! Apply magnetic fields between these iterations integer, intent(in) :: iteration !! Current iteration @@ -159,7 +158,7 @@ contains end if ! Scf-cycle for constraint fields, based either on torque or on fields alone - if (lbfield_constr .and. itscf0 <= iteration .and. iteration <= itscf1) then + if (itscf0 <= iteration .and. iteration <= itscf1) then if (constr_mode == 3) then bfield%bfield_constr(:) = bfield%bfield_constr(:) - bfield%mag_torque(:) / mag_mom_len else if (constr_mode == 2) then diff --git a/source/KKRnano/source/datastructures/InputParamsNew.txt b/source/KKRnano/source/datastructures/InputParamsNew.txt index 20f704c84..14f919210 100644 --- a/source/KKRnano/source/datastructures/InputParamsNew.txt +++ b/source/KKRnano/source/datastructures/InputParamsNew.txt @@ -126,21 +126,15 @@ d a_log 0.025D0 l enable_quad_prec .FALSE. ### [bfields] turn on or off everything concerning non-collinear magnetic fields (called lbfield in KKRhost code) l noncobfield .FALSE. -### [bfields] constraining fields, non-collinear magnetic field (called lbfield_constr in KKRhost code) -l constr_field .FALSE. -### [bfields] apply same field to all atoms (True) or individual fields to each atom (called lbfield_all in KKRhost code) -l same_bfield .FALSE. +### [bfields] turn on or off external non-collinear magnetic fields +l external_bfield .FALSE. ### [bfields] apply only transversal bfield (called lbfield_trans in KKRhost code) l trans_bfield .FALSE. ### [bfields] apply magnetic field only in the muffin-tin (called lbfield_mt in KKRhost code) l mt_bfield .FALSE. -### [bfields] calculate magnetic torque (called ltorque in KKRhost code) -l torque .FALSE. -### [bfields] spin (0), orbital (1), spin+orbial (2) fields (called ibfield in KKRhost code) -i ibfield 0 -### [bfields] type of contraint (0 = torque, 1 = magnetic moment) (called ibfield_constr in KKRhost code) -i iconstr 0 -### [bfields] start magnetic field at iteration itscf0 (called ibfield_itscf0 in KKRhost code) +### [bfields] start applying magnetic field at iteration (called ibfield_itscf0 in KKRhost code) i itbfield0 0 ### [bfields] stop applying magnetic field after iteration itscf1 (called ibfield_itscf1 in KKRhost code) i itbfield1 10000 +### [bfields] control verbosity concerning non-collinear magnetic fields +i bfield_verbosity 2 diff --git a/source/KKRnano/source/datastructures/InputParams_mod.F90 b/source/KKRnano/source/datastructures/InputParams_mod.F90 index 73ebc527a..a20d14824 100644 --- a/source/KKRnano/source/datastructures/InputParams_mod.F90 +++ b/source/KKRnano/source/datastructures/InputParams_mod.F90 @@ -79,15 +79,12 @@ module InputParams_mod double precision :: a_log logical :: enable_quad_prec logical :: noncobfield - logical :: constr_field - logical :: same_bfield + logical :: external_bfield logical :: trans_bfield logical :: mt_bfield - logical :: torque - integer :: ibfield - integer :: iconstr integer :: itbfield0 integer :: itbfield1 + integer :: bfield_verbosity endtype ! InputParams @@ -657,21 +654,12 @@ integer function getValues(filename, self) result(ierror) destroy_and_return endif - ierror = getValue(cr, "constr_field", self%constr_field , def=.FALSE.) + ierror = getValue(cr, "external_bfield", self%external_bfield , def=.FALSE.) if (ierror == use_default) then - write(*,*) "WARNING: Bad/no value given for constr_field. Set constr_field to .FALSE." + write(*,*) "WARNING: Bad/no value given for external_bfield. Set external_bfield to .FALSE." ierror = 0 ! ok, no error elseif (ierror /= 0) then - write(*,*) "Bad/no value given for constr_field." - destroy_and_return - endif - - ierror = getValue(cr, "same_bfield", self%same_bfield , def=.FALSE.) - if (ierror == use_default) then - write(*,*) "WARNING: Bad/no value given for same_bfield. Set same_bfield to .FALSE." - ierror = 0 ! ok, no error - elseif (ierror /= 0) then - write(*,*) "Bad/no value given for same_bfield." + write(*,*) "Bad/no value given for external_bfield." destroy_and_return endif @@ -693,33 +681,6 @@ integer function getValues(filename, self) result(ierror) destroy_and_return endif - ierror = getValue(cr, "torque", self%torque , def=.FALSE.) - if (ierror == use_default) then - write(*,*) "WARNING: Bad/no value given for torque. Set torque to .FALSE." - ierror = 0 ! ok, no error - elseif (ierror /= 0) then - write(*,*) "Bad/no value given for torque." - destroy_and_return - endif - - ierror = getValue(cr, "ibfield", self%ibfield , def=0) - if (ierror == use_default) then - write(*,*) "WARNING: Bad/no value given for ibfield. Set ibfield to 0" - ierror = 0 ! ok, no error - elseif (ierror /= 0) then - write(*,*) "Bad/no value given for ibfield." - destroy_and_return - endif - - ierror = getValue(cr, "iconstr", self%iconstr , def=0) - if (ierror == use_default) then - write(*,*) "WARNING: Bad/no value given for iconstr. Set iconstr to 0" - ierror = 0 ! ok, no error - elseif (ierror /= 0) then - write(*,*) "Bad/no value given for iconstr." - destroy_and_return - endif - ierror = getValue(cr, "itbfield0", self%itbfield0 , def=0) if (ierror == use_default) then write(*,*) "WARNING: Bad/no value given for itbfield0. Set itbfield0 to 0" @@ -738,6 +699,15 @@ integer function getValues(filename, self) result(ierror) destroy_and_return endif + ierror = getValue(cr, "bfield_verbosity", self%bfield_verbosity , def=2) + if (ierror == use_default) then + write(*,*) "WARNING: Bad/no value given for bfield_verbosity. Set bfield_verbosity to 2" + ierror = 0 ! ok, no error + elseif (ierror /= 0) then + write(*,*) "Bad/no value given for bfield_verbosity." + destroy_and_return + endif + write(*,*) "Finished reading information from input.conf" destroy_and_return #undef destroy_and_return diff --git a/source/KKRnano/source/wrappers_mod.F90 b/source/KKRnano/source/wrappers_mod.F90 index 88f360b1b..261a31d99 100644 --- a/source/KKRnano/source/wrappers_mod.F90 +++ b/source/KKRnano/source/wrappers_mod.F90 @@ -99,7 +99,7 @@ module wrappers_mod den,espv,rho2ns,r2nef, gmatn(:,:,:,1), muorb, & ! just one spin component of gmatn needed atomdata%potential%lpot,lmaxd,mesh%irmd,chebmesh%irmd_new,iemxd, params%soc,params%enable_quad_prec, & bfield, imt, iteration_number, params%itbfield0, params%itbfield1, & - params%noncobfield, params%constr_field, params%torque, params%trans_bfield, params%mt_bfield) + params%noncobfield, params%trans_bfield, params%mt_bfield) ! calculate correct orbital moment do ispin=1,nspind -- GitLab