~siesta-maint/siesta/trunk

« back to all changes in this revision

Viewing changes to Src/m_ts_full_scat.F90

  • Committer: Alberto Garcia
  • Date: 2018-12-07 21:49:47 UTC
  • mfrom: (746.1.31 trunk-clean)
  • Revision ID: albertog@icmab.es-20181207214947-cybuwjlg6fyuvz01
Compliance of code to F2003 standard plus further cleaning

* Removal of non-f2003 constructs

   -- non-standard intrinsics: dcmplx, dconjg, dreal, etc
   -- non-standard extensions for POSIX calls (chdir, system, etc)
      These are now implemented using the 'iso_c_binding' module.
   -- non-standard variable declarations (real*8, etc)

* Avoid use of 'mpif.h' in MPI compilation. Replace by 'use mpi'

* Avoid compiler checks for reallocation on assignment: whole array
  assignments now use array-section syntax in the left-hand side.

* Fix of a number of potential interface problems

* Some dead code has been removed

* The code should now compile with only 'benign' warnings with the
  gfortran flags:

      FFLAGS= -std=f2003 -Wall -Wextra -Wrealloc-lhs-all

A sizable number of warnings refer to unused arguments and variables,
which have not been removed in this pass.


Show diffs side-by-side

added added

removed removed

Lines of Context:
61
61
! *********************
62
62
! * LOCAL variables   *
63
63
! *********************
64
 
    complex(dp), parameter :: z0 = dcmplx(0._dp, 0._dp)
65
 
    complex(dp), parameter :: z1 = dcmplx(1._dp, 0._dp)
66
 
    complex(dp), parameter :: zi = dcmplx(0._dp, 1._dp)
 
64
    complex(dp), parameter :: z0 = cmplx(0._dp, 0._dp,dp)
 
65
    complex(dp), parameter :: z1 = cmplx(1._dp, 0._dp,dp)
 
66
    complex(dp), parameter :: zi = cmplx(0._dp, 1._dp,dp)
67
67
 
68
68
    integer :: i, NB, ind, iB
69
69
 
82
82
       ! Collect the top row of complex conjugated Gf
83
83
       ind = no_u_TS * no * iB + 1
84
84
       do i = 1 , no
85
 
          GGG(ind:ind-1+no) = dconjg(Gf(iB*no+1:(iB+1)*no,i))
 
85
          GGG(ind:ind-1+no) = conjg(Gf(iB*no+1:(iB+1)*no,i))
86
86
          ind = ind + no
87
87
       end do
88
88
       ind = no_u_TS * no * iB + 1
121
121
       ind = no_u_TS * no * NB + 1
122
122
       do i = 1 , no
123
123
          ! So this is the complex conjugated of the iB'th block
124
 
          GGG(ind:ind-1+iB) = dconjg(Gf(NB*no+1:NB*no+iB,i))
 
124
          GGG(ind:ind-1+iB) = conjg(Gf(NB*no+1:NB*no+iB,i))
125
125
          ind = ind + iB
126
126
       end do
127
127
       ind = no_u_TS * no * NB + 1
173
173
    do j = 1 , N
174
174
       do i = 1 , j
175
175
!          if(ionode)print *,M(j,i),M(i,j)
176
 
          M(j,i) = dimag(M(i,j))
 
176
          M(j,i) = aimag(M(i,j))
177
177
          M(i,j) = M(j,i)
178
178
 
179
179
       end do
287
287
         call die('GFB: Wrong size of Green function')
288
288
 
289
289
    ! Create the RHS for inversion...
290
 
    GF(:) = dcmplx(0._dp,0._dp)
 
290
    GF(:) = cmplx(0._dp,0._dp,dp)
291
291
 
292
292
    o = 0
293
293
    do iEl = 1 , N_Elec
294
294
       i = Elecs(iEl)%idx_o
295
295
       off_row = i - orb_offset(i) - 1
296
296
       do i = 1 , TotUsedOrbs(Elecs(iEl))
297
 
          GF(o*no_u_TS+off_row+i) = dcmplx(1._dp,0._dp)
 
297
          GF(o*no_u_TS+off_row+i) = cmplx(1._dp,0._dp,dp)
298
298
          o = o + 1
299
299
       end do
300
300
    end do
371
371
         call die('GFP: Wrong size of Green function')
372
372
 
373
373
    ! initialize
374
 
    GF(:) = dcmplx(0._dp,0._dp)
 
374
    GF(:) = cmplx(0._dp,0._dp,dp)
375
375
 
376
376
    i = 0
377
377
    do j = 0 , no - 1
380
380
          i = i + 1
381
381
          if ( .not. any(OrbInElec(Elecs,i) .and. Elecs(:)%DM_update == 0) ) then
382
382
             
383
 
             GF(j*no_u_TS+i-orb_offset(i)) = dcmplx(1._dp,0._dp)
 
383
             GF(j*no_u_TS+i-orb_offset(i)) = cmplx(1._dp,0._dp,dp)
384
384
             found = .true.
385
385
          end if
386
386
       end do