The Elk Code
modtest.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2009 J. K. Dewhurst, S. Sharma and E. K. U. Gross
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 module modtest
7 
8 use modmpi
9 
10 ! if test is .true. then the test variables are written to file
11 logical test
12 
13 contains
14 
15 subroutine writetest(id,descr,nv,iv,iva,tol,rv,rva,zv,zva)
16 implicit none
17 ! arguments
18 integer, intent(in) :: id
19 character(*), intent(in) :: descr
20 integer, optional, intent(in) :: nv
21 integer, optional, intent(in) :: iv
22 integer, optional, intent(in) :: iva(*)
23 real(8), optional, intent(in) :: tol
24 real(8), optional, intent(in) :: rv
25 real(8), optional, intent(in) :: rva(*)
26 complex(8), optional, intent(in) :: zv
27 complex(8), optional, intent(in) :: zva(*)
28 ! local variables
29 integer j
30 character(256) fname
31 if (.not.test) return
32 if (.not.mp_mpi) return
33 if ((id < 0).or.(id > 999)) then
34  write(*,*)
35  write(*,'("Error(writetest): id out of range : ",I8)') id
36  write(*,*)
37  stop
38 end if
39 if (present(iva).or.present(rva).or.present(zva)) then
40  if (.not.present(nv)) then
41  write(*,*)
42  write(*,'("Error(writetest): missing argument nv")')
43  write(*,*)
44  stop
45  else
46  if (nv < 1) then
47  write(*,*)
48  write(*,'("Error(writetest): nv < 1 : ",I8)') nv
49  write(*,*)
50  stop
51  end if
52  end if
53 end if
54 if (present(rv).or.present(rva).or.present(zv).or.present(zva)) then
55  if (.not.present(tol)) then
56  write(*,*)
57  write(*,'("Error(writetest): missing argument tol")')
58  write(*,*)
59  stop
60  end if
61 end if
62 write(fname,'("TEST_",I3.3,".OUT")') id
63 !$OMP CRITICAL(writetest_)
64 open(90,file=trim(fname),form='FORMATTED',action='WRITE')
65 write(90,'("''",A,"''")') trim(descr)
66 if (present(iv)) then
67  write(90,'(2I8)') 1,1
68  write(90,'(2I8)') 1,iv
69 else if (present(rv)) then
70  write(90,'(2I8)') 2,1
71  write(90,'(G24.14)') tol
72  write(90,'(I8,G24.14)') 1,rv
73 else if (present(zv)) then
74  write(90,'(2I8)') 3,1
75  write(90,'(G24.14)') tol
76  write(90,'(I8,2G24.14)') 1,dble(zv),aimag(zv)
77 else if (present(iva)) then
78  write(90,'(2I8)') 1,nv
79  do j=1,nv
80  write(90,'(2I8)') j,iva(j)
81  end do
82 else if (present(rva)) then
83  write(90,'(2I8)') 2,nv
84  write(90,'(G24.14)') tol
85  do j=1,nv
86  write(90,'(I8,G24.14)') j,rva(j)
87  end do
88 else if (present(zva)) then
89  write(90,'(2I8)') 3,nv
90  write(90,'(G24.14)') tol
91  do j=1,nv
92  write(90,'(I8,2G24.14)') j,dble(zva(j)),aimag(zva(j))
93  end do
94 end if
95 close(90)
96 !$OMP END CRITICAL(writetest_)
97 end subroutine
98 
99 end module
100 
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
logical mp_mpi
Definition: modmpi.f90:17
logical test
Definition: modtest.f90:11
Definition: modmpi.f90:6