TABLE OF CONTENTS


ABINIT/testfi [ Functions ]

[ Top ] [ Functions ]

NAME

 testfi

FUNCTION

 Routine "Final test" for generation of the test report in the status file:
 if it appears that the run was a "Build-in Test", then
 compare the final values of different quantities to the reference
 values, here hard-coded.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA,XG,GMR)
 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 .

INPUTS

  builtintest=number of the builtintest, from the input file.
  etotal=total energy (sum of 7 contributions) (hartree)
  filstat=name of the status file
  fred(3,natom)=symmetrized gradient of etotal with respect to tn
  natom=number of atoms in cell.
  strten(6)=components of the stress tensor (hartree/bohr^3)
  xred(3,natom)=reduced atomic coordinates

OUTPUT

  (only writing)

PARENTS

      abinit

CHILDREN

SOURCE

 39 #if defined HAVE_CONFIG_H
 40 #include "config.h"
 41 #endif
 42 
 43 #include "abi_common.h"
 44 
 45 
 46 subroutine testfi(builtintest,etotal,filstat,fred,natom,strten,xred)
 47 
 48  use defs_basis
 49  use m_profiling_abi
 50  use m_errors
 51 
 52  use m_io_tools,     only : open_file
 53 
 54 !This section has been created automatically by the script Abilint (TD).
 55 !Do not modify the following lines by hand.
 56 #undef ABI_FUNC
 57 #define ABI_FUNC 'testfi'
 58 !End of the abilint section
 59 
 60  implicit none
 61 
 62 !Arguments ------------------------------------
 63 !scalars
 64  integer,intent(in) :: builtintest,natom
 65  real(dp),intent(in) :: etotal
 66  character(len=fnlen),intent(in) :: filstat
 67 !arrays
 68  real(dp),intent(in) :: fred(3,natom),strten(6),xred(3,natom)
 69 
 70 !Local variables-------------------------------
 71  character(len=fnlen) :: testname(7)='         '
 72  character(len=*), parameter :: format01000="(a,d22.14,a,d12.4)"
 73  character(len=500) :: msg
 74 !scalars
 75  integer,parameter :: mtest=7
 76  integer :: iatom,ii,problem,tok,temp_unit
 77  real(dp) :: etot_mxdev,etot_ref,fred_mxdev,strten_mxdev,xred_mxdev
 78 !arrays
 79  integer,parameter :: natom_test(mtest)=(/2,1,2,1,1,1,2/)
 80  real(dp) :: fred_ref(3,2),strten_ref(6),xred_ref(3,2)
 81 
 82 ! ***********************************************************************
 83 
 84  testname(1)='fast'
 85  testname(2)='v1'
 86  testname(3)='v5'
 87  testname(4)='bigdft'
 88  testname(5)='etsf_io'
 89  testname(6)='libxc'
 90  testname(7)='wannier90'
 91 
 92 !---------------------------------------------------------
 93 
 94 !Now, open the status file, and either delete it, or produce a report
 95  if (open_file(filstat,msg,newunit=temp_unit,form='formatted',status='unknown') /= 0) then
 96    MSG_ERROR(msg)
 97  end if
 98 
 99  if(builtintest==0)then
100 
101    close (temp_unit,status='delete')
102 
103  else
104 
105 !  Note : all processors have their own file, so no special
106 !  attention must be paid to the parallel case.
107 
108    write(temp_unit,*)
109    write(temp_unit,*)'Status file, reporting on built-in test ',trim(testname(builtintest))
110    write(temp_unit,*)
111 
112 !  Define reference values, as well as maximum tolerable deviation
113    if(builtintest==1)then
114 
115      etot_ref=-1.05814441948188d+00
116      etot_mxdev=1.0d-9
117      xred_ref(1:3,1)=(/ -0.65048430042634D-01 , 0.0_dp , 0.0_dp /)
118      xred_ref(1:3,2)=(/  0.65048430042634D-01 , 0.0_dp , 0.0_dp /)
119 !    xred(*,*) are reduced coordinates
120      xred_mxdev=1.0d-6
121      fred_ref(1:3,1:2)= 0.0_dp
122 !    fred(*,*) are gradients with respect to reduced coordinates
123      fred_mxdev=5.0d-4
124      strten_ref(1:6)=(/ 0.149D-04  , 0.560D-04 , 0.560D-04 ,&
125 &     0.0_dp , 0.0_dp , 0.0_dp /)
126      strten_mxdev=1.0d-5
127 
128    else if(builtintest==2)then
129 
130 !    This value of etot is accurate to about 1.0d-12
131      etot_ref=-.70811958266295D+02
132 !    Initialisation conditions might give fluctuations
133      etot_mxdev=3.0d-7
134      xred_ref(1:3,1)=(/ 0.0_dp , 0.0_dp , 0.0_dp /)
135      xred_mxdev=1.0d-12
136      fred_ref(1:3,1)=(/ 0.0_dp , 0.0_dp , 0.0_dp /)
137      fred_mxdev=1.0d-12
138 !    This value of strten is accurate to at least 1.0d-8
139      strten_ref(1:3)= 5.09324870E-03
140      strten_ref(4:6)= 0.0_dp
141      strten_mxdev=1.0d-6
142 
143    else if(builtintest==3)then
144 
145      etot_ref=-.892746696311772D+01
146      etot_mxdev=1.0d-8
147      xred_ref(1:3,1)=(/ -0.125_dp , 0.0_dp , 0.0_dp /)
148      xred_ref(1:3,2)=(/  0.125_dp , 0.0_dp , 0.0_dp /)
149      xred_mxdev=1.0d-12
150      fred_ref(1:3,1)=(/ -.140263620278D+00 , 0.0_dp , 0.0_dp /)
151      fred_ref(1:3,2)=(/  .140013483725D+00 , 0.0_dp , 0.0_dp /)
152      fred_mxdev=1.0d-3
153      strten_ref(1:6)=(/  1.3949d-3 ,  1.3643d-3 , 1.3643d-3 ,&
154 &     .0_dp ,  .0_dp  ,  .0_dp    /)
155      strten_mxdev=1.0d-5
156 
157 !    Bigdft
158    else if(builtintest==4)then
159 
160 !    This value of etot is accurate to about 1.0d-12
161      etot_ref=-0.56231990141D+00
162 !    Initialisation conditions might give fluctuations
163      etot_mxdev=3.0d-7
164      xred_ref(1:3,1)=(/ 0.5_dp , 0.5_dp , 0.5_dp /)
165      xred_mxdev=1.0d-12
166      fred_ref(1:3,1)=(/ 0.0_dp , 0.0_dp , 0.0_dp /)
167      fred_mxdev=1.0d-12
168      strten_ref(1)=-0.22537238382D-01
169      strten_ref(2)=-0.22536232141D-01
170      strten_ref(3)=-0.22529038043D-01
171      strten_ref(4)= 0.39104928264D-06
172      strten_ref(5)= 0.62823846408D-06
173      strten_ref(6)= 0.38414125548E-09
174      strten_mxdev=1.0d-7
175 
176 !    ETSF_IO
177    else if(builtintest==5)then
178 
179 !    This value of etot is accurate to about 1.0d-12
180      etot_ref=-0.33307825915795D+02
181 !    Initialisation conditions might give fluctuations
182      etot_mxdev=3.0d-7
183      xred_ref(1:3,1)=(/ 0.0_dp , 0.0_dp , 0.0_dp /)
184      xred_mxdev=1.0d-12
185      fred_ref(1:3,1)=(/ 0.0_dp , -0.000000120877_dp , -0.000000164287_dp /)
186      fred_mxdev=1.0d-7
187 !    This value of strten is accurate to at least 1.0d-8
188      strten_ref(1)= 0.13569940015175D-01
189      strten_ref(2)= 0.33822108610352D-01
190      strten_ref(3)= 0.37991262607028D-01
191      strten_ref(4:6)= 0.0_dp
192      strten_mxdev=1.0d-6
193 
194 !    LibXC
195    else if(builtintest==6)then
196 
197      etot_ref=-0.55196688523897D+01
198      etot_mxdev=5.0d-7
199      xred_ref(1:3,1)=(/ 0.0_dp , 0.0_dp , 0.0_dp /)
200      xred_mxdev=1.0d-12
201      fred_ref(1:3,1)=(/ 0.0_dp , 0.0_dp , 0.0_dp /)
202      fred_mxdev=1.0d-12
203      strten_ref(1:3)= 0.13246699375127D-04
204      strten_ref(4:6)= 0.0_dp
205      strten_mxdev=1.0d-8
206 
207 !    Wannier90
208    else if(builtintest==7)then
209 
210 !    This value of etot is accurate to about 1.0d-12
211      etot_ref=-0.10620085133544D+02
212 !    Initialisation conditions might give fluctuations
213      etot_mxdev=3.0d-7
214      xred_ref(1:3,1)=(/ 0.0_dp , 0.0_dp , 0.0_dp /)
215      xred_ref(1:3,2)=(/ 0.25_dp , 0.25_dp , 0.25_dp /)
216      xred_mxdev=1.0d-12
217      fred_ref(1:3,1)=(/ 0.0_dp , 0.0_dp , 0.0_dp /)
218      fred_ref(1:3,2)=(/ 0.0_dp , 0.0_dp , 0.0_dp /)
219      fred_mxdev=1.0d-12
220 !    This value of strten is accurate to at least 1.0d-8
221      strten_ref(1:3)=0.31413922197317D-03
222      strten_ref(4:6)= 0.0_dp
223      strten_mxdev=1.0d-6
224 
225 !    End of the reference value set up, for different tests
226    end if
227 
228 !  Compare reference and actual values
229 
230    problem=0
231 
232 !  Take care of total energy
233    if(abs(etot_ref-etotal)>etot_mxdev)then
234      problem=1
235      write(temp_unit,'(a)')' Error for total energy : '
236      write(temp_unit,format01000)'        expected ',etot_ref,&
237 &     '  with maximum   deviation',etot_mxdev
238      write(temp_unit,format01000)'        computed ',etotal,&
239 &     '  with effective deviation',abs(etotal-etot_ref)
240    end if
241 
242 !  Take care of nuclei positions
243    tok=1
244    do iatom=1,natom
245      do ii=1,3
246        if(abs(xred_ref(ii,iatom)-xred(ii,iatom))>xred_mxdev)then
247          tok=0
248          write(temp_unit, '(a,i1,a,i1,a)' )&
249 &         ' Error for nuclei position xred(',ii,',',iatom,')'
250          write(temp_unit,format01000)'        expected ',xred_ref(ii,iatom),&
251 &         '  with maximum   deviation',xred_mxdev
252          write(temp_unit,format01000)'        computed ',xred(ii,iatom),&
253 &         '  with effective deviation',&
254 &         abs( xred(ii,iatom)-xred_ref(ii,iatom) )
255        end if
256      end do
257    end do
258    if(tok==0)problem=1
259 
260 !  Take care of forces
261    tok=1
262    do iatom=1,natom
263      do ii=1,3
264        if(abs(fred_ref(ii,iatom)-fred(ii,iatom))&
265 &       >fred_mxdev)then
266          tok=0
267          write(temp_unit, '(a,i1,a,i1,a)' )&
268 &         ' Error for force fred(',ii,',',iatom,')'
269          write(temp_unit,format01000)'        expected ',fred_ref(ii,iatom),&
270 &         '  with maximum   deviation',fred_mxdev
271          write(temp_unit,format01000)'        computed ',fred(ii,iatom),&
272 &         '  with effective deviation',&
273 &         abs( fred(ii,iatom)-fred_ref(ii,iatom) )
274        end if
275      end do
276    end do
277    if(tok==0)problem=1
278 
279 !  Take care of stress
280    tok=1
281    do ii=1,6
282      if(abs(strten_ref(ii)-strten(ii))>strten_mxdev)then
283        tok=0
284        write(temp_unit,'(a,i1,a)')' Error for stress strten(',ii,')'
285        write(temp_unit,format01000)'        expected ',strten_ref(ii),&
286 &       '  with maximum   deviation',strten_mxdev
287        write(temp_unit,format01000)'        computed ',strten(ii),&
288 &       '  with effective deviation',&
289 &       abs( strten(ii)-strten_ref(ii) )
290      end if
291    end do
292    if(tok==0)problem=1
293 
294    if(problem==0)then
295      write(temp_unit,'(a)')' ==> The run finished cleanly.'
296      write(temp_unit,'(a,a)')'     Moreover, comparison of the total energy, and other (few) ',&
297 &     'relevant quantities with reference values has been successful.'
298      write(temp_unit,'(a)')'     This does not mean that no problem is present, however. '
299      write(temp_unit,'(a)')'     Please run the complete set of ABINIT tests to gain a better confidence in your installation.'
300    end if
301 
302    write(temp_unit,*)
303 
304    close(temp_unit)
305 
306  end if !  End of the choice between produce a report, and produce no report
307 
308 end subroutine testfi