The Elk Code
 
Loading...
Searching...
No Matches
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
6module modtest
7
8use modmpi
9
10! if test is .true. then the test variables are written to file
11logical test
12
13contains
14
15subroutine writetest(id,descr,nv,iv,iva,tol,rv,rva,zv,zva)
16implicit none
17! arguments
18integer, intent(in) :: id
19character(*), intent(in) :: descr
20integer, optional, intent(in) :: nv
21integer, optional, intent(in) :: iv
22integer, optional, intent(in) :: iva(*)
23real(8), optional, intent(in) :: tol
24real(8), optional, intent(in) :: rv
25real(8), optional, intent(in) :: rva(*)
26complex(8), optional, intent(in) :: zv
27complex(8), optional, intent(in) :: zva(*)
28! local variables
29integer j
30character(256) fname
31if (.not.test) return
32if (.not.mp_mpi) return
33if ((id < 0).or.(id > 999)) then
34 write(*,*)
35 write(*,'("Error(writetest): id out of range : ",I8)') id
36 write(*,*)
37 stop
38end if
39if (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
53end if
54if (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
61end if
62write(fname,'("TEST_",I3.3,".OUT")') id
63!$OMP CRITICAL(writetest_)
64open(90,file=trim(fname),form='FORMATTED',action='WRITE')
65write(90,'("''",A,"''")') trim(descr)
66if (present(iv)) then
67 write(90,'(2I8)') 1,1
68 write(90,'(2I8)') 1,iv
69else if (present(rv)) then
70 write(90,'(2I8)') 2,1
71 write(90,'(G24.14)') tol
72 write(90,'(I8,G24.14)') 1,rv
73else 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)
77else 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
82else 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
88else 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
94end if
95close(90)
96!$OMP END CRITICAL(writetest_)
97end subroutine
98
99end module
100
logical mp_mpi
Definition modmpi.f90:17
logical test
Definition modtest.f90:11
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16