TABLE OF CONTENTS


ABINIT/m_slc_primitive_potential [ Modules ]

[ Top ] [ Modules ]

NAME

 m_slc_primitive_potential

FUNCTION

 This module contains the atomic structures and the spin-lattice coupling terms inside the primitive cell
 which can be directly mapped to the netcdf file. It is not for the calculation, but for constructing
 the hamiltonian in supercell.

 Datatypes:
  slc_primitive_potential_t

 Subroutines:

COPYRIGHT

 Copyright (C) 2001-2024 ABINIT group (hexu,nehelbig)
 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

 24 #if defined HAVE_CONFIG_H
 25 #include "config.h"
 26 #endif
 27 #include "abi_common.h"
 28 
 29 module m_slc_primitive_potential
 30 
 31   use m_nctk
 32 #if defined HAVE_NETCDF
 33   use netcdf
 34 #endif
 35   use m_abicore
 36   use defs_basis
 37   use m_errors
 38 
 39   use m_slc_potential
 40   use m_spmat_ndcoo
 41 
 42   use m_abstract_potential, only: abstract_potential_t
 43   use m_dynamic_array, only: int2d_array_type
 44   use m_mpi_scheduler, only: init_mpi_info
 45   use m_multibinit_cell, only: mbcell_t, mbsupercell_t
 46   use m_multibinit_dataset, only: multibinit_dtset_type
 47   use m_primitive_potential, only: primitive_potential_t
 48 
 49   use m_supercell_maker, only: supercell_maker_t
 50   use m_xmpi
 51 
 52   implicit none
 53 
 54   type, public, extends(primitive_potential_t) :: slc_primitive_potential_t
 55      integer :: natom, nspin        ! on every mpi node
 56      logical :: has_bilin=.False.   ! bilinear coupling term, i.e. liu        
 57      logical :: has_linquad=.False. ! spin first then lattice, i.e. niuv
 58      logical :: has_quadlin=.False. ! spin first then lattice, i.e. oiju
 59      logical :: has_biquad=.False.  ! biquadratic coupling term, i.e. tijuv
 60 
 61      !four possible coupling terms, bilinear, quadratic-linear, linear-quadratic, and biquadratic
 62      !all the following are only on master node 
 63      type(ndcoo_mat_t) :: liu           ! parameter values
 64      type(ndcoo_mat_t) :: niuv          ! parameter values
 65      type(ndcoo_mat_t) :: oiju          ! parameter values
 66      type(ndcoo_mat_t) :: tijuv         ! parameter values
 67      !R vector of i index is always zero
 68      type(int2d_array_type) :: lRulist  ! R vector for u index of liu 
 69      type(int2d_array_type) :: nRulist  ! R vector for u index of niuv 
 70      type(int2d_array_type) :: nRvlist  ! R vector for v index of niuv 
 71      type(int2d_array_type) :: oRjlist  ! R vector for j index of oiju 
 72      type(int2d_array_type) :: oRulist  ! R vector for u index of oiju
 73      type(int2d_array_type) :: tRjlist  ! R vector for j index of tijuv
 74      type(int2d_array_type) :: tRulist  ! R vector for u index of tijuv
 75      type(int2d_array_type) :: tRvlist  ! R vector for v index of tijuv
 76 
 77    contains
 78      procedure:: initialize
 79      procedure:: finalize
 80      procedure :: set_liu
 81      procedure :: set_liu_1term
 82      procedure :: set_niuv
 83      procedure :: set_niuv_1term
 84      procedure :: set_oiju
 85      procedure :: set_oiju_1term
 86      procedure :: set_tijuv
 87      procedure :: set_tijuv_1term
 88      procedure :: load_from_files
 89      procedure :: read_netcdf
 90      procedure :: read_liu
 91      procedure :: read_niuv
 92      procedure :: read_oiju
 93      procedure :: read_tijuv
 94      procedure :: fill_supercell
 95      procedure :: fill_liu
 96      procedure :: fill_niuv
 97      procedure :: fill_oiju
 98      procedure :: fill_tijuv
 99      procedure :: set_liu_sc
100      procedure :: set_niuv_sc
101      procedure :: set_oiju_sc
102      procedure :: set_tijuv_sc
103   end type slc_primitive_potential_t
104 
105 contains
106 
107   subroutine initialize(self, primcell)
108     class(slc_primitive_potential_t), intent(inout) :: self
109     type(mbcell_t), target, intent(inout) :: primcell
110     !integer, intent(in) :: nspin
111     self%primcell=>primcell
112     self%label="SLC_primitive_potential"
113     self%has_spin=.True.
114     self%has_displacement=.True.
115     self%has_strain=.False.
116     self%has_lwf=.False.
117   end subroutine initialize
118 
119 
120   subroutine finalize(self)
121     class(slc_primitive_potential_t), intent(inout) :: self
122     call self%oiju%finalize()
123     call self%oRjlist%finalize()
124     call self%oRulist%finalize()
125     call self%tijuv%finalize()
126     call self%tRjlist%finalize()
127     call self%tRulist%finalize()
128     call self%tRvlist%finalize()
129     nullify(self%primcell)
130     self%nspin=0
131     self%natom=0
132     self%label="Destroyed SLC_primitive_potential"
133     call self%primitive_potential_t%finalize()
134   end subroutine finalize
135 
136   subroutine load_from_files(self, params, fnames)
137     class(slc_primitive_potential_t), intent(inout) :: self
138     type(multibinit_dtset_type), intent(in) :: params
139     character(len=fnlen), intent(in) :: fnames(:)
140     character(len=fnlen) :: ncdf_fname
141     character(len=500) :: message
142     integer :: ii
143     ABI_UNUSED(fnames)
144     if (xmpi_comm_rank(xmpi_world)==0) then
145        ncdf_fname=params%slc_pot_fname
146        write(message,'(a,(81a, 80a),3a)') ch10,('=',ii=1,80),ch10,ch10,&
147             &     '- Reading spin-lattice coupling terms from ', trim(ncdf_fname)
148        call wrtout(ab_out,message,'COLL')
149        call wrtout(std_out,message,'COLL')
150     endif
151     call self%read_netcdf(trim(ncdf_fname)//char(0))
152     
153     ABI_UNUSED_A(params)
154   end subroutine load_from_files
155 
156   !-----------------------------------
157   ! reading parameters from ncdf file
158   !-----------------------------------
159   subroutine read_netcdf(self, ncdf_fname)
160     class(slc_primitive_potential_t), intent(inout) :: self
161     character (len = *),              intent(in)    :: ncdf_fname
162 
163     integer:: ncid, ncerr, comm
164 
165 #if defined HAVE_NETCDF
166 
167     comm = xmpi_world
168 
169     ncerr = nctk_open_read(ncid, ncdf_fname, comm)
170 
171     ! read primcell info
172     ncerr=nctk_get_dim(ncid, "natom", self%natom)
173     ncerr=nctk_get_dim(ncid, "nspin", self%nspin)
174 
175     ! read different coupling terms if they are present in the file
176     call self%read_liu(ncid)
177     call self%read_niuv(ncid)
178     call self%read_oiju(ncid)
179     call self%read_tijuv(ncid)
180     
181     ncerr = nf90_close(ncid)
182     if(ncerr /= NF90_NOERR) then
183       write(std_out,'(A25)') 'Could not close netcdf file'
184     endif
185 #else
186     ABI_ERROR("Multibint should be installed with netcdf enabled to run this.")
187 
188 #endif
189 
190   end subroutine read_netcdf
191 
192   !---------------------------------------
193   ! reading liu parameters from ncdf file
194   ! TODO: test
195   !---------------------------------------
196   subroutine read_liu(self, ncid)
197     class(slc_primitive_potential_t), intent(inout) :: self
198     integer,                          intent(in)    :: ncid    
199 
200     integer :: ncerr, dimid, ndata, varid
201     integer, allocatable :: ilist(:), ulist(:,:)
202     real(dp), allocatable :: vallist(:)
203 #if defined HAVE_NETCDF
204 
205     ncerr = nf90_inq_dimid(ncid, "spin_lattice_Liu_number_of_entries", dimid)
206     if(ncerr /= NF90_NOERR) then
207       ndata = 0
208     else
209       self%has_bilin=.True.
210       ncerr = nctk_get_dim(ncid, "spin_lattice_Liu_number_of_entries", ndata)
211       ABI_MALLOC(ilist, (ndata))
212       ABI_MALLOC(ulist, (4,ndata))
213       ABI_MALLOC(vallist, (ndata))
214 
215       varid = nctk_idname(ncid, "spin_lattice_Liu_ilist")
216       ncerr = nf90_get_var(ncid, varid, ilist)
217       call netcdf_check(ncerr, "when reading Liu_ilist") 
218 
219       varid = nctk_idname(ncid, "spin_lattice_Liu_ulist")
220       ncerr = nf90_get_var(ncid, varid, ulist)
221       call netcdf_check(ncerr, "when reading Liu_ulist")
222 
223       varid = nctk_idname(ncid, "spin_lattice_Liu_valuelist")
224       ncerr = nf90_get_var(ncid, varid, vallist)
225       call netcdf_check(ncerr, "when reading Liu_valuelist")
226 
227       !change units from eV to Ha and Ang to Bohr
228       vallist(:) = vallist(:)*eV_Ha*Bohr_Ang
229   
230       !fill the sparse matrix for liu parameters
231       call self%set_liu(ndata, ilist, ulist, vallist)
232 
233       ABI_SFREE(ilist)
234       ABI_SFREE(ulist)
235       ABI_SFREE(vallist)
236     endif
237 
238     write(std_out,'(A8,I10,A11)') 'L_iu:  ', ndata, 'terms read'  
239 #else
240     ABI_ERROR("Multibinit should be installed with netcdf.")
241 #endif
242 
243   end subroutine read_liu
244 
245   !----------------------------------------
246   ! reading niuv parameters from ncdf file
247   ! TODO: test
248   !----------------------------------------
249   subroutine read_niuv(self, ncid)
250     class(slc_primitive_potential_t), intent(inout) :: self
251     integer,                          intent(in)    :: ncid    
252 
253     integer :: ncerr, dimid, ndata, varid
254     integer, allocatable :: ilist(:), ulist(:,:), vlist(:,:)
255     real(dp), allocatable :: vallist(:)
256 #if defined HAVE_NETCDF
257     ncerr = nf90_inq_dimid(ncid, "spin_lattice_Niuv_number_of_entries", dimid)
258     if(ncerr /= NF90_NOERR) then
259       ndata = 0
260     else
261       self%has_linquad=.True.
262       ncerr = nctk_get_dim(ncid, "spin_lattice_Niuv_number_of_entries", ndata)
263       ABI_MALLOC(ilist, (ndata))
264       ABI_MALLOC(ulist, (4,ndata))
265       ABI_MALLOC(vlist, (4,ndata))
266       ABI_MALLOC(vallist, (ndata))
267 
268       varid = nctk_idname(ncid, "spin_lattice_Niuv_ilist")
269       ncerr = nf90_get_var(ncid, varid, ilist)
270       call netcdf_check(ncerr, "when reading Niuv_ilist") 
271 
272       varid = nctk_idname(ncid, "spin_lattice_Niuv_ulist")
273       ncerr = nf90_get_var(ncid, varid, ulist)
274       call netcdf_check(ncerr, "when reading Niuv_ulist")
275 
276       varid = nctk_idname(ncid, "spin_lattice_Niuv_vlist")
277       ncerr = nf90_get_var(ncid, varid, ulist)
278       call netcdf_check(ncerr, "when reading Niuv_vlist")
279 
280       varid = nctk_idname(ncid, "spin_lattice_Niuv_valuelist")
281       ncerr = nf90_get_var(ncid, varid, vallist)
282       call netcdf_check(ncerr, "when reading Niuv_valuelist")
283 
284       !change units from eV to Ha and Ang to Bohr
285       vallist(:) = vallist(:)*eV_Ha*(Bohr_Ang*Bohr_Ang)
286   
287       !fill the sparse matrix for liu parameters
288       call self%set_niuv(ndata, ilist, ulist, vlist, vallist)
289 
290       ABI_SFREE(ilist)
291       ABI_SFREE(ulist)
292       ABI_SFREE(vlist)
293       ABI_SFREE(vallist)
294     endif
295 
296     write(std_out,'(A8,I10,A11)') 'N_iuv: ', ndata, 'terms read'  
297 #else
298     ABI_ERROR("Multibinit should be installed with netcdf") 
299 #endif
300 
301   end subroutine read_niuv
302 
303   !---------------------------------------
304   ! reading oiju parameters from ncdf file
305   !---------------------------------------
306   subroutine read_oiju(self, ncid)  
307     class(slc_primitive_potential_t), intent(inout) :: self
308     integer,                          intent(in)    :: ncid    
309 
310     integer :: ncerr, dimid, ndata, varid
311     integer, allocatable :: ilist(:), jlist(:,:), ulist(:,:)
312     real(dp), allocatable :: vallist(:)
313 
314 #if defined HAVE_NETCDF
315     ncerr = nf90_inq_dimid(ncid, "spin_lattice_Oiju_number_of_entries", dimid)
316     if(ncerr /= NF90_NOERR) then
317       write(std_out,'(A20)') 'No O_iju term found'
318     else
319       self%has_quadlin=.True.
320       ncerr = nctk_get_dim(ncid, "spin_lattice_Oiju_number_of_entries", ndata)
321       ABI_MALLOC(ilist, (ndata))
322       ABI_MALLOC(jlist, (4,ndata))
323       ABI_MALLOC(ulist, (4,ndata))
324       ABI_MALLOC(vallist, (ndata))
325 
326       varid = nctk_idname(ncid, "spin_lattice_Oiju_ilist")
327       ncerr = nf90_get_var(ncid, varid, ilist)
328       call netcdf_check(ncerr, "when reading Oiju_ilist") 
329 
330       varid = nctk_idname(ncid, "spin_lattice_Oiju_jlist")
331       ncerr = nf90_get_var(ncid, varid, jlist)
332       call netcdf_check(ncerr, "when reading Oiju_jlist")
333 
334       varid = nctk_idname(ncid, "spin_lattice_Oiju_ulist")
335       ncerr = nf90_get_var(ncid, varid, ulist)
336       call netcdf_check(ncerr, "when reading Oiju_ulist")
337 
338       varid = nctk_idname(ncid, "spin_lattice_Oiju_valuelist")
339       ncerr = nf90_get_var(ncid, varid, vallist)
340       call netcdf_check(ncerr, "when reading spin_lattice_Oiju_valuelist")
341 
342       write(std_out,'(A8,I10,A11)') 'O_iju: ', ndata, 'terms read'  
343 
344       !change units from eV to Ha and Ang to Bohr
345       vallist(:) = vallist(:)*eV_Ha*Bohr_Ang
346 
347       !fill the sparse matrix for oiju parameters
348       call self%set_oiju(ndata, ilist, jlist, ulist, vallist)
349 
350       ABI_SFREE(ilist)
351       ABI_SFREE(jlist)
352       ABI_SFREE(ulist)
353       ABI_SFREE(vallist)
354     endif
355 #else
356     ABI_ERROR('Multibinit should be install with netcdf to run this.')
357 #endif
358 
359   end subroutine read_oiju
360 
361   !-----------------------------------------
362   ! reading tijuv parameters from ncdf file
363   !-----------------------------------------
364   subroutine read_tijuv(self, ncid)
365     class(slc_primitive_potential_t), intent(inout) :: self
366     integer,                          intent(in)    :: ncid    
367 
368     integer :: ncerr, dimid, ndata, varid
369     integer, allocatable :: ilist(:), jlist(:,:), ulist(:,:), vlist(:,:)
370     real(dp), allocatable :: vallist(:)
371 
372 #if defined HAVE_NETCDF
373     ncerr = nf90_inq_dimid(ncid, "spin_lattice_Tijuv_number_of_entries", dimid)
374     if(ncerr /= NF90_NOERR) then
375       ndata = 0
376     else
377       self%has_biquad=.True.
378       ncerr = nctk_get_dim(ncid, "spin_lattice_Tijuv_number_of_entries", ndata)
379       ABI_MALLOC(ilist, (ndata))
380       ABI_MALLOC(jlist, (4, ndata))
381       ABI_MALLOC(ulist, (4, ndata))
382       ABI_MALLOC(vlist, (4, ndata))
383       ABI_MALLOC(vallist, (ndata))
384 
385       varid = nctk_idname(ncid, "spin_lattice_Tijuv_ilist")
386       ncerr = nf90_get_var(ncid, varid, ilist)
387       call netcdf_check(ncerr, "when reading Tijuv_ilist") 
388 
389       varid = nctk_idname(ncid, "spin_lattice_Tijuv_jlist")
390       ncerr = nf90_get_var(ncid, varid, jlist)
391       call netcdf_check(ncerr, "when reading Tijuv_jlist")
392 
393       varid = nctk_idname(ncid, "spin_lattice_Tijuv_ulist")
394       ncerr = nf90_get_var(ncid, varid, ulist)
395       call netcdf_check(ncerr, "when reading Tijuv_ulist")
396 
397       varid = nctk_idname(ncid, "spin_lattice_Tijuv_vlist")
398       ncerr = nf90_get_var(ncid, varid, vlist)
399       call netcdf_check(ncerr, "when reading Tijuv_vlist")
400 
401       varid = nctk_idname(ncid, "spin_lattice_Tijuv_valuelist")
402       ncerr = nf90_get_var(ncid, varid, vallist)
403       call netcdf_check(ncerr, "when reading Tijuv_valuelist")
404 
405       write(std_out,'(A8,I10,A11)') 'T_ijuv:', ndata, 'terms read'  
406 
407       !change units from eV to Ha and Ang to Bohr
408       vallist(:) = vallist(:)*eV_Ha*(Bohr_Ang*Bohr_Ang)
409   
410       !fill the sparse matrix for tijuv parameters
411       call self%set_tijuv(ndata, ilist, jlist, ulist, vlist, vallist)
412 
413       ABI_SFREE(ilist)
414       ABI_SFREE(jlist)
415       ABI_SFREE(ulist)
416       ABI_SFREE(vlist)
417       ABI_SFREE(vallist)
418     endif
419 
420 #else
421     ABI_ERROR('Multibinit should be install with netcdf to run this.')
422 #endif
423 
424 
425   end subroutine read_tijuv
426 
427   !---------------------------------------
428   ! store liu parameters in sparse matrix
429   ! TODO: test
430   !---------------------------------------
431   subroutine set_liu(self, nn, ilist, ulist, vallist)
432 
433     class(slc_primitive_potential_t), intent(inout) :: self
434     integer,                          intent(inout) :: nn
435     integer,                          intent(in)    :: ilist(nn)
436     integer,                          intent(in)    :: ulist(4,nn)
437     real(dp),                         intent(in)    :: vallist(nn)
438 
439     integer :: idx
440 
441     call self%liu%initialize(mshape=[-1, self%nspin*3, self%natom*3])
442     
443     if (xmpi_comm_rank(xmpi_world)==0) then
444       do idx=1, nn
445         call self%set_liu_1term(ilist(idx), ulist(1,idx), ulist(2:4,idx), vallist(idx))
446       end do
447     endif
448   end subroutine set_liu
449 
450   !------------------------------------
451   ! add one entry to sparse liu matrix
452   ! TODO: test
453   !------------------------------------
454   subroutine set_liu_1term(self, ii, uu, Ru, val)
455     class(slc_primitive_potential_t), intent(inout) :: self
456     integer,                          intent(in)    :: ii    
457     integer,                          intent(in)    :: uu
458     integer,                          intent(in)    :: Ru(3)
459     real(dp),                         intent(in)    :: val
460 
461     integer :: indRu
462     
463     call self%lRulist%push_unique(Ru, position=indRu)
464     
465     call self%liu%add_entry(ind=[indRu, ii, uu], val=val)
466 
467   end subroutine set_liu_1term
468 
469   !----------------------------------------
470   ! store niuv parameters in sparse matrix
471   ! TODO: test
472   !----------------------------------------
473   subroutine set_niuv(self, nn, ilist, ulist, vlist, vallist)
474 
475     class(slc_primitive_potential_t), intent(inout) :: self
476     integer,                          intent(inout) :: nn
477     integer,                          intent(in)    :: ilist(nn)
478     integer,                          intent(in)    :: ulist(4,nn)
479     integer,                          intent(in)    :: vlist(4,nn)
480     real(dp),                         intent(in)    :: vallist(nn)
481 
482     integer :: idx
483    
484     call self%niuv%initialize(mshape=[-1, -1, self%nspin*3, self%natom*3, self%natom*3])
485     
486     if (xmpi_comm_rank(xmpi_world)==0) then
487       do idx=1, nn
488         call self%set_niuv_1term(ilist(idx), ulist(1,idx), vlist(1,idx), ulist(2:4,idx), vlist(2:4,idx), vallist(idx))
489       end do
490     endif
491 
492   end subroutine set_niuv
493 
494   !-------------------------------------
495   ! add one entry to sparse niuv matrix
496   !-------------------------------------
497   subroutine set_niuv_1term(self, ii, uu, vv, Ru, Rv, val)
498     class(slc_primitive_potential_t), intent(inout) :: self
499     integer,                          intent(in)    :: ii    
500     integer,                          intent(in)    :: uu
501     integer,                          intent(in)    :: vv
502     integer,                          intent(in)    :: Ru(3)
503     integer,                          intent(in)    :: Rv(3)
504     real(dp),                         intent(in)    :: val
505 
506     integer :: indRu, indRv
507     
508     call self%nRulist%push_unique(Ru, position=indRu)
509     call self%nRvlist%push_unique(Rv, position=indRv)
510     
511     call self%niuv%add_entry(ind=[indRu, indRv, ii, uu, vv], val=val)
512 
513   end subroutine set_niuv_1term
514 
515 
516   !----------------------------------------
517   ! store oiju parameters in sparse matrix
518   ! TODO: test
519   !----------------------------------------
520   subroutine set_oiju(self, nn, ilist, jlist, ulist, vallist)
521 
522     class(slc_primitive_potential_t), intent(inout) :: self
523     integer,                          intent(inout) :: nn
524     integer,                          intent(in)    :: ilist(nn)
525     integer,                          intent(in)    :: jlist(4,nn)
526     integer,                          intent(in)    :: ulist(4,nn)
527     real(dp),                         intent(in)    :: vallist(nn)
528 
529     integer :: idx
530    
531     call self%oiju%initialize(mshape=[-1, -1, self%nspin*3, self%nspin*3, self%natom*3])
532     
533     if (xmpi_comm_rank(xmpi_world)==0) then
534       do idx=1, nn
535         call self%set_oiju_1term(ilist(idx), jlist(1,idx), ulist(1,idx), jlist(2:4,idx), ulist(2:4,idx), vallist(idx))
536       end do
537     endif
538 
539     call self%oiju%sum_duplicates()
540   end subroutine set_oiju
541 
542   !-------------------------------------
543   ! add one entry to sparse oiju matrix
544   ! TODO: test
545   !-------------------------------------
546   subroutine set_oiju_1term(self, ii, jj, uu, Rj, Ru, val)
547     class(slc_primitive_potential_t), intent(inout) :: self
548     integer,                          intent(in)    :: ii    
549     integer,                          intent(in)    :: jj
550     integer,                          intent(in)    :: uu
551     integer,                          intent(in)    :: Rj(3)
552     integer,                          intent(in)    :: Ru(3)
553     real(dp),                         intent(in)    :: val
554 
555     integer :: indRj, indRu
556     
557     call self%oRjlist%push_unique(Rj, position=indRj)
558     call self%oRulist%push_unique(Ru, position=indRu)
559     
560     call self%oiju%add_entry(ind=[indRj, indRu, ii, jj, uu], val=val)
561 
562   end subroutine set_oiju_1term
563 
564   !----------------------------------------
565   ! store tijuv parameters in sparse matrix
566   !----------------------------------------
567   subroutine set_tijuv(self, nn, ilist, jlist, ulist, vlist, vallist)
568 
569     class(slc_primitive_potential_t), intent(inout) :: self
570     integer,                          intent(inout) :: nn
571     integer,                          intent(in)    :: ilist(nn)
572     integer,                          intent(in)    :: jlist(4, nn)
573     integer,                          intent(in)    :: ulist(4, nn)
574     integer,                          intent(in)    :: vlist(4, nn)
575     real(dp),                         intent(in)    :: vallist(nn)
576 
577     integer :: idx
578    
579     call self%tijuv%initialize(mshape=[-1, -1, -1, self%nspin*3, self%nspin*3, self%natom*3, self%natom*3])
580     
581     if (xmpi_comm_rank(xmpi_world)==0) then
582       do idx=1, nn
583         call self%set_tijuv_1term(ilist(idx), jlist(1, idx), ulist(1, idx), vlist(1, idx), & 
584                                 & jlist(2:4,idx), ulist(2:4,idx), vlist(2:4,idx), vallist(idx))
585       end do
586     endif
587 
588     call self%tijuv%sum_duplicates()
589 
590   end subroutine set_tijuv
591 
592   !-------------------------------------
593   ! add one entry to sparse tijuv matrix
594   !-------------------------------------
595   subroutine set_tijuv_1term(self, ii, jj, uu, vv, Rj, Ru, Rv, val)
596     class(slc_primitive_potential_t), intent(inout) :: self
597     integer,                          intent(in)    :: ii    
598     integer,                          intent(in)    :: jj
599     integer,                          intent(in)    :: uu
600     integer,                          intent(in)    :: vv
601     integer,                          intent(in)    :: Rj(3)
602     integer,                          intent(in)    :: Ru(3)
603     integer,                          intent(in)    :: Rv(3)
604     real(dp),                         intent(in)    :: val
605 
606     integer :: indRj, indRu, indRv
607     
608     call self%tRjlist%push_unique(Rj, position=indRj)
609     call self%tRulist%push_unique(Ru, position=indRu)
610     call self%tRvlist%push_unique(Rv, position=indRv)
611     
612     call self%tijuv%add_entry(ind=[indRj, indRu, indRv, ii, jj, uu, vv], val=val)
613 
614   end subroutine set_tijuv_1term
615 
616 
617   !-----------------------------------------------------------------
618   ! transfer parameter information from primitive cell to supercell
619   ! TODO: test
620   !-----------------------------------------------------------------
621   subroutine fill_supercell(self, scmaker, params, scpot, supercell)
622     class(slc_primitive_potential_t) , intent(inout) :: self
623     type(supercell_maker_t),           intent(inout) :: scmaker
624     type(multibinit_dtset_type),       intent(inout) :: params
625     class(abstract_potential_t), pointer, intent(inout) :: scpot
626     type(mbsupercell_t), target :: supercell
627 
628     integer :: nspin, sc_nspin, natom, sc_natom
629     integer :: master, my_rank, comm, nproc, ierr
630     logical :: iam_master
631 
632     call init_mpi_info(master, iam_master, my_rank, comm, nproc) 
633 
634     nspin=self%nspin
635     natom=self%natom
636     sc_nspin=nspin * scmaker%ncells
637     sc_natom=natom * scmaker%ncells
638     
639     call xmpi_bcast(sc_nspin, master, comm, ierr)
640     call xmpi_bcast(sc_natom, master, comm, ierr)
641     ABI_MALLOC_TYPE_SCALAR(slc_potential_t, scpot)
642 
643     select type(scpot) ! use select type because properties only defined for slc_potential are used
644     type is (slc_potential_t) 
645       call scpot%initialize(sc_nspin, sc_natom)
646       call scpot%set_supercell(supercell)
647       call scpot%set_params(params)
648       ! fill different coupling terms
649       if (iam_master) then
650         call self%fill_liu(scpot, scmaker)
651         call self%fill_oiju(scpot, scmaker)
652         call self%fill_niuv(scpot, scmaker)
653         call self%fill_tijuv(scpot, scmaker)
654 
655         !Write information which terms are used
656         write(std_out,'(A55)') 'Using the following terms for the spin-lattice coupling'
657         if(scpot%has_bilin)   write(std_out,'(A19)') 'Bilinear term: Liu'
658         if(scpot%has_quadlin) write(std_out,'(A28)') 'Quadratic-linear term: Oiju'
659         if(scpot%has_linquad) write(std_out,'(A28)') 'Linear-quadratic term: Niuv'
660         if(scpot%has_biquad)  write(std_out,'(A24)') 'Biquadratic term: Tijuv'
661       endif
662     end select
663   end subroutine fill_supercell
664 
665   !--------------------------------------------------------
666   ! Check for each term if it is needed in the calculation
667   ! and present in the ncdf input then put into supercell
668   !--------------------------------------------------------
669   subroutine fill_liu(self, scpot, scmaker)
670     class(slc_primitive_potential_t) , intent(inout) :: self
671     type(supercell_maker_t),           intent(inout) :: scmaker
672     class(slc_potential_t),            intent(inout) :: scpot
673 
674     if(scpot%has_bilin) then
675       if(self%has_bilin) then !did we actually find this in the netcdf file
676         call self%liu%sum_duplicates()
677         call self%set_liu_sc(scpot, scmaker)
678       else
679         ABI_ERROR("No parameters for bilinear coupling available. Check your input and parameter files.")
680         scpot%has_bilin = .False.
681       endif
682     endif
683   end subroutine fill_liu
684   
685   subroutine fill_niuv(self, scpot, scmaker)
686     class(slc_primitive_potential_t) , intent(inout) :: self
687     type(supercell_maker_t),           intent(inout) :: scmaker
688     class(slc_potential_t),            intent(inout) :: scpot
689 
690     if(scpot%has_linquad) then
691       if(self%has_linquad) then
692         call self%niuv%sum_duplicates()
693         call self%set_niuv_sc(scpot, scmaker)
694       else
695         ABI_ERROR("No parameters for linear-quadratic coupling available. Check your input and parameter files.")
696         scpot%has_linquad = .False.
697       endif
698     endif
699   end subroutine fill_niuv
700 
701   subroutine fill_oiju(self, scpot, scmaker)
702     class(slc_primitive_potential_t) , intent(inout) :: self
703     type(supercell_maker_t),           intent(inout) :: scmaker
704     class(slc_potential_t),            intent(inout) :: scpot
705 
706     if(scpot%has_quadlin) then
707       if(self%has_quadlin) then
708         call self%oiju%sum_duplicates()
709         call self%set_oiju_sc(scpot, scmaker)
710        else
711          ABI_ERROR("No parameters for quadratic-linear coupling available. Check your input and parameter files.")
712          scpot%has_quadlin = .False.
713        endif
714      endif
715   end subroutine fill_oiju
716 
717   subroutine fill_tijuv(self, scpot, scmaker)
718     class(slc_primitive_potential_t) , intent(inout) :: self
719     type(supercell_maker_t),           intent(inout) :: scmaker
720     class(slc_potential_t),            intent(inout) :: scpot
721 
722     if(scpot%has_biquad) then
723       if(self%has_biquad) then
724         call self%tijuv%sum_duplicates()
725         call self%set_tijuv_sc(scpot, scmaker)
726       else
727         ABI_ERROR("No parameters for biquadratic coupling available. Check your input and parameter files.")
728         scpot%has_biquad = .False.
729       endif
730     endif
731   end subroutine fill_tijuv
732 
733   !------------------------------
734   ! fill liu terms in supercell
735   ! TODO: test
736   !------------------------------
737   subroutine set_liu_sc(self, scpot, scmaker)
738     class(slc_primitive_potential_t), intent(inout) :: self
739     type(slc_potential_t),            intent(inout) :: scpot
740     type(supercell_maker_t),          intent(inout) :: scmaker
741         
742     integer :: icell, Ru_prim(3), liu_ind(3), iRu, i_prim, u_prim, inz
743     integer :: ngroup
744     integer, allocatable :: i_sc(:), u_sc(:), Ru_sc(:,:)
745     integer, allocatable :: i1list(:), ise(:)
746     real(dp) :: val_sc(scmaker%ncells)
747 
748     integer :: master, my_rank, comm, nproc
749     logical :: iam_master
750 
751     call init_mpi_info(master, iam_master, my_rank, comm, nproc) 
752 
753     if(iam_master) then
754       call scpot%liu_sc%initialize(mshape=[self%nspin*3, self%natom*3])
755     endif
756 
757     do inz=1, self%liu%nnz
758       liu_ind=self%liu%get_ind_inz(inz)
759       iRu=liu_ind(1)
760       i_prim=liu_ind(2)
761       u_prim=liu_ind(3)
762       Ru_prim=self%lRulist%data(:,iRu)
763       call scmaker%trans_i(nbasis=self%nspin*3, i=i_prim, i_sc=i_sc) 
764       call scmaker%trans_j_and_Rj(nbasis=self%natom*3, j=u_prim, Rj=Ru_prim, j_sc=u_sc, Rj_sc=Ru_sc)
765       val_sc(:)= self%liu%val%data(inz)
766       do icell=1, scmaker%ncells
767         call scpot%add_liu_term(i_sc(icell), u_sc(icell), val_sc(icell))
768       end do
769       ABI_SFREE(i_sc)
770       ABI_SFREE(u_sc)
771       ABI_SFREE(Ru_sc)
772     end do
773     
774     call scpot%liu_sc%group_by_1dim(ngroup, i1list, ise)
775     ABI_SFREE(i1list)
776     ABI_SFREE(ise)
777 
778   end subroutine set_liu_sc
779 
780   !------------------------------
781   ! fill niuv terms in supercell
782   ! TODO: test
783   !------------------------------
784   subroutine set_niuv_sc(self, scpot, scmaker)
785     class(slc_primitive_potential_t), intent(inout) :: self
786     type(slc_potential_t),            intent(inout) :: scpot
787     type(supercell_maker_t),          intent(inout) :: scmaker
788         
789     integer :: icell, Ru_prim(3), Rv_prim(3), niuv_ind(5), iRu, iRv, i_prim, u_prim, v_prim, inz
790     integer :: ngroup
791     integer, allocatable :: i_sc(:), u_sc(:), v_sc(:), Ru_sc(:,:), Rv_sc(:,:)
792     integer, allocatable :: i1list(:), ise(:)
793     real(dp) :: val_sc(scmaker%ncells)
794 
795     integer :: master, my_rank, comm, nproc
796     logical :: iam_master
797 
798     call init_mpi_info(master, iam_master, my_rank, comm, nproc) 
799 
800     if(iam_master) then
801       call scpot%niuv_sc%initialize(mshape=[self%nspin*3, self%natom*3, self%natom*3])
802     endif
803 
804     do inz=1, self%niuv%nnz
805       niuv_ind=self%niuv%get_ind_inz(inz)
806       iRu=niuv_ind(1)
807       iRv=niuv_ind(2)
808       i_prim=niuv_ind(3)
809       u_prim=niuv_ind(4)
810       v_prim=niuv_ind(5)
811       Ru_prim=self%nRulist%data(:,iRu)
812       Rv_prim=self%nRvlist%data(:,iRv)
813       call scmaker%trans_i(nbasis=self%nspin*3, i=i_prim, i_sc=i_sc) 
814       call scmaker%trans_j_and_Rj(nbasis=self%natom*3, j=u_prim, Rj=Ru_prim, j_sc=u_sc, Rj_sc=Ru_sc)
815       call scmaker%trans_j_and_Rj(nbasis=self%natom*3, j=v_prim, Rj=Rv_prim, j_sc=v_sc, Rj_sc=Rv_sc)
816       val_sc(:)= self%niuv%val%data(inz)
817       do icell=1, scmaker%ncells
818         !call scpot%add_niuv_term(i_sc(icell), u_sc(icell), v_sc(icell), val_sc(icell))
819       end do
820       ABI_SFREE(i_sc)
821       ABI_SFREE(u_sc)
822       ABI_SFREE(v_sc)
823       ABI_SFREE(Ru_sc)
824       ABI_SFREE(Rv_sc)
825     end do
826     
827     call scpot%niuv_sc%group_by_1dim(ngroup, i1list, ise)
828     ABI_SFREE(i1list)
829     ABI_SFREE(ise)
830 
831   end subroutine set_niuv_sc
832 
833 
834   !------------------------------
835   ! fill oiju terms in supercell
836   ! TODO: test
837   !------------------------------
838   subroutine set_oiju_sc(self, scpot, scmaker)
839     class(slc_primitive_potential_t), intent(inout) :: self
840     class(slc_potential_t),            intent(inout) :: scpot
841     type(supercell_maker_t),          intent(inout) :: scmaker
842         
843     integer :: icell, Rj_prim(3), Ru_prim(3), oiju_ind(5), iRj, iRu, i_prim, j_prim, u_prim, inz
844     integer :: ngroup
845     integer, allocatable :: i_sc(:), j_sc(:), u_sc(:), Rj_sc(:, :), Ru_sc(:,:)
846     integer, allocatable :: i1list(:), ise(:)
847     real(dp) :: val_sc(scmaker%ncells)
848 
849     integer :: master, my_rank, comm, nproc
850     logical :: iam_master
851 
852     call init_mpi_info(master, iam_master, my_rank, comm, nproc) 
853 
854     if(iam_master) then
855       call scpot%oiju_sc%initialize(mshape=[self%nspin*3, self%nspin*3, self%natom*3])
856     endif
857 
858     do inz=1, self%oiju%nnz
859       oiju_ind=self%oiju%get_ind_inz(inz)
860       iRj=oiju_ind(1)
861       iRu=oiju_ind(2)
862       i_prim=oiju_ind(3)
863       j_prim=oiju_ind(4)
864       u_prim=oiju_ind(5)
865       Rj_prim=self%oRjlist%data(:,iRj)
866       Ru_prim=self%oRulist%data(:,iRu)
867       call scmaker%trans_i(nbasis=self%nspin*3, i=i_prim, i_sc=i_sc) 
868       call scmaker%trans_j_and_Rj(nbasis=self%nspin*3, j=j_prim, Rj=Rj_prim, j_sc=j_sc, Rj_sc=Rj_sc)
869       call scmaker%trans_j_and_Rj(nbasis=self%natom*3, j=u_prim, Rj=Ru_prim, j_sc=u_sc, Rj_sc=Ru_sc)
870       val_sc(:)= self%oiju%val%data(inz)
871       do icell=1, scmaker%ncells
872         call scpot%add_oiju_term(i_sc(icell), j_sc(icell), u_sc(icell), val_sc(icell))
873       end do
874       ABI_SFREE(i_sc)
875       ABI_SFREE(j_sc)
876       ABI_SFREE(u_sc)
877       ABI_SFREE(Rj_sc)
878       ABI_SFREE(Ru_sc)
879     end do
880     
881     call scpot%oiju_sc%group_by_1dim(ngroup, i1list, ise)
882     ABI_SFREE(i1list)
883     ABI_SFREE(ise)
884 
885   end subroutine set_oiju_sc
886 
887   subroutine set_tijuv_sc(self, scpot, scmaker)
888     class(slc_primitive_potential_t), intent(inout) :: self
889     type(slc_potential_t),            intent(inout) :: scpot
890     type(supercell_maker_t),          intent(inout) :: scmaker
891         
892     integer :: icell, Rj(3), Ru(3), Rv(3), tijuv_ind(7), iRj, iRu, iRv, ii, ij, iu, iv, inz
893     integer, allocatable :: i_sc(:), j_sc(:), u_sc(:), v_sc(:), Rj_sc(:, :), Ru_sc(:,:), Rv_sc(:,:)
894     real(dp) :: val_sc(scmaker%ncells)
895 
896     integer :: master, my_rank, comm, nproc
897     logical :: iam_master
898 
899     call init_mpi_info(master, iam_master, my_rank, comm, nproc) 
900 
901     if(iam_master) then
902       call scpot%tijuv_sc%initialize(mshape=[self%nspin*3, self%nspin*3, self%natom*3, self%natom*3])
903       call scpot%tuvij_sc%initialize(mshape=[self%natom*3, self%natom*3, self%nspin*3, self%nspin*3])
904     endif
905 
906     do inz=1, self%tijuv%nnz
907       tijuv_ind=self%tijuv%get_ind_inz(inz)
908       iRj=tijuv_ind(1)
909       iRu=tijuv_ind(2)
910       iRv=tijuv_ind(3)
911       ii=tijuv_ind(4)
912       ij=tijuv_ind(5)
913       iu=tijuv_ind(6)
914       iv=tijuv_ind(7)
915       Rj=self%tRjlist%data(:,iRj)
916       Ru=self%tRulist%data(:,iRu)
917       Rv=self%tRvlist%data(:,iRv)
918       call scmaker%trans_i(nbasis=self%nspin*3, i=ii, i_sc=i_sc)
919       call scmaker%trans_j_and_Rj(nbasis=self%nspin*3, j=ij, Rj=Rj, j_sc=j_sc, Rj_sc=Rj_sc)
920       call scmaker%trans_j_and_Rj(nbasis=self%natom*3, j=iu, Rj=Ru, j_sc=u_sc, Rj_sc=Ru_sc)
921       call scmaker%trans_j_and_Rj(nbasis=self%natom*3, j=iv, Rj=Rv, j_sc=v_sc, Rj_sc=Rv_sc)
922       val_sc(:)= self%tijuv%val%data(inz)
923       do icell=1, scmaker%ncells
924         call scpot%add_tijuv_term(i_sc(icell), j_sc(icell), u_sc(icell), v_sc(icell), val_sc(icell))
925       end do
926       ABI_SFREE(i_sc)
927       ABI_SFREE(j_sc)
928       ABI_SFREE(u_sc)
929       ABI_SFREE(v_sc)
930       ABI_SFREE(Rj_sc)
931       ABI_SFREE(Ru_sc)
932       ABI_SFREE(Rv_sc)
933     end do
934     
935     call scpot%tijuv_sc%group_by_pair()
936     call scpot%tuvij_sc%group_by_pair()
937 
938   end subroutine set_tijuv_sc