atomes 1.1.16
atomes: an atomic scale modeling tool box
Loading...
Searching...
No Matches
pdb.F90
Go to the documentation of this file.
1! This file is part of the 'atomes' software.
2!
3! 'atomes' is free software: you can redistribute it and/or modify it under the terms
4! of the GNU Affero General Public License as published by the Free Software Foundation,
5! either version 3 of the License, or (at your option) any later version.
6!
7! 'atomes' is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
8! without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
9! See the GNU General Public License for more details.
10!
11! You should have received a copy of the GNU Affero General Public License along with 'atomes'.
12! If not, see <https://www.gnu.org/licenses/>
13!
14! Copyright (C) 2022-2024 by CNRS and University of Strasbourg
15!
20
21INTEGER (KIND=c_int) FUNCTION read_pdb (pdb_f, lpdb) bind (C,NAME='read_pdb_')
22
23!
24! Lecture of atom types and coordinates - pdb files
25!
26
27USE parameters
28USE mendeleiev
29
30IMPLICIT NONE
31
32LOGICAL :: isok
33INTEGER (KIND=c_int), INTENT(IN) :: lpdb
34INTEGER, DIMENSION(:), ALLOCATABLE :: linea
35LOGICAL :: with_cell = .false.
36CHARACTER (KIND=c_char), DIMENSION(*), INTENT(IN) :: pdb_f
37CHARACTER (LEN=lpdb) :: pdb_file
38CHARACTER (LEN=4) :: aname
39CHARACTER (LEN=6) :: rid
40INTEGER :: asi
41CHARACTER (LEN=1) :: alb
42CHARACTER (LEN=3) :: acc
43CHARACTER (LEN=1) :: cic
44INTEGER :: rsn
45CHARACTER (LEN=1) :: ic
46DOUBLE PRECISION :: oca
47DOUBLE PRECISION thf
48CHARACTER (LEN=2) :: sym
49CHARACTER (LEN=2) :: cha
50REAL, DIMENSION(6) :: cell
51CHARACTER (LEN=11) :: spg
52INTEGER :: zval
53
54do i=1, lpdb
55 pdb_file(i:i) = pdb_f(i)
56enddo
57cell(:) = 0.0
58
59inquire (file=pdb_file, exist=isok)
60if (isok) then
61
62 open (unit=20, file=pdb_file, action='read', status='old')
63 ns=1
64 na=0
65 l = 0
66 do while (.true.)
67 read (20, '(A6)', end=999) rid
68 l = l+1
69 if (rid.eq.'HETATM' .or. rid.eq.'ATOM ') na=na+1
70 if (rid.eq. 'CRYST1') with_cell = .true.
71 enddo
72 999 continue
73 if (na .ne. 0) then
74 rewind(20)
75 if (allocated(linea)) deallocate(linea)
76 allocate (linea(l), stat=err)
77 if (err .ne. 0) then
78 call show_error ("Impossible to allocate memory"//char(0), &
79 "Function: read_pdb"//char(0), "Table: LINEA"//char(0))
80 endif
81 l=0
82 do while (.true.)
83 read (20, '(A6)', end=888) rid
84 l = l+1
85 linea(l)=0
86 if (rid.eq.'HETATM' .or. rid.eq.'ATOM ') linea(l) = 1
87 if (rid.eq.'CRYST1') linea(l) = 2
88 enddo
89 888 continue
90 rewind(20)
91 if (allocated(fullpos)) deallocate(fullpos)
92 allocate (fullpos(na,3,ns), stat=err)
93 if (err .ne. 0) then
94 call show_error ("Impossible to allocate memory"//char(0), &
95 "Function: read_pdb"//char(0), "Table: FULLPOS"//char(0))
96 endif
97 if (allocated(tab_of_type)) deallocate(tab_of_type)
98 allocate (tab_of_type(na), stat=err)
99 if (err .ne. 0) then
100 call show_error ("Impossible to allocate memory"//char(0), &
101 "Function: read_pdb"//char(0), "Table: TAB_OF_TYPE"//char(0))
102 endif
103 do i=1, ns
104 m = 0
105 do j=1, l-1
106 if (linea(j) .eq. 1) then
107 m = m+1
108 read (20, 001) rid, &
109 asi, &
110 aname, &
111 alb, &
112 acc, &
113 cic, &
114 rsn, &
115 ic, &
116 fullpos(m,1,i), &
117 fullpos(m,2,i), &
118 fullpos(m,3,i), &
119 oca, &
120 thf, &
121 sym, &
122 cha
123 ! write (6, '("m= ",i4," x= ",f15.10,", y= ",f15.10,", z= ",f15.10)') m, FULLPOS(m,1,i),FULLPOS(m,2,i),FULLPOS(m,3,i)
124 tab_of_type(m)=atompdb(j, sym)
125 else if (linea(j) .eq. 2) then
126 read (20, 002) rid, &
127 cell(1), &
128 cell(2), &
129 cell(3), &
130 cell(4), &
131 cell(5), &
132 cell(6), &
133 spg, &
134 zval
135 else
136 read (20, *)
137 endif
138 enddo
139 enddo
140 if (allocated(linea)) deallocate(linea)
141 if (with_cell) then
142 call cell_data_from_pdb (cell(1), cell(2), cell(3), cell(4), cell(5), cell(6))
143 endif
144 read_pdb=0
145 else
146 read_pdb=2
147 endif
148 close(20)
149else
150 read_pdb=1
151endif
152
153001 FORMAT (a6,i5,1x,a4,a1,a3,1x,a1,i4,a1,3x,f8.3,f8.3,f8.3,f6.2,f6.2,10x,a2,a2)
154002 FORMAT (a6,f9.3,f9.3,f9.3,f7.2,f7.2,f7.2,1x,a11,i4)
155
156!
157! PDB format v3.30
158!
159! HETATM / ATOM keyword
160!
161!---------------------------------------------------------------------------
162!Field | Column | FORTRAN |
163! No. | range | format | Description
164!---------------------------------------------------------------------------
165! 1. | 1 - 6 | A6 | "ATOM " or "HETATM"
166! 2. | 7 - 11 | I5 | Atom serial number
167! - | 12 - 12 | 1X | Blank
168! 3. | 13 - 16 | A4 | Atom name
169! 4. | 17 - 17 | A1 | Alternative location code (if any)
170! 5. | 18 - 20 | A3 | Standard 3-letter amino acid code for residue
171! - | 21 - 21 | 1X | Blank
172! 6. | 22 - 22 | A1 | Chain identifier code
173! 7. | 23 - 26 | I4 | Residue sequence number
174! 8. | 27 - 27 | A1 | Insertion code (if any)
175! - | 28 - 30 | 3X | Blank
176! 9. | 31 - 38 | F8.3 | Atom's x-coordinate
177! 10. | 39 - 46 | F8.3 | Atom's y-coordinate
178! 11. | 47 - 54 | F8.3 | Atom's z-coordinate
179! 12. | 55 - 60 | F6.2 | Occupancy value for atom
180! 13. | 61 - 66 | F6.2 | B-value (thermal factor)
181! - | 67 - 76 | 10X | Blank
182! 14. | 77 - 78 | A2 | Element symbol
183! 15. | 79 - 80 | A2 | Charge
184!---------------------------------------------------------------------------
185
186! CRYST1 keyword
187!
188!---------------------------------------------------------------------------
189!Field | Column | FORTRAN |
190! No. | range | format | Description
191!---------------------------------------------------------------------------
192! 1. | 1 - 6 | A6 | "CRYST1"
193! 2. | 7 - 15 | F9.3 | a
194! 3. | 16 - 24 | F9.3 | b
195! 4. | 25 - 33 | F9.3 | c
196! 5. | 34 - 40 | F7.2 | alpha
197! 6. | 41 - 47 | F7.2 | beta
198! 7. | 48 - 54 | F7.2 | gamma
199! - | 55 - 55 | 1X | Blank
200! 8. | 56 - 66 | A11 | Space group
201! 9. | 67 - 70 | I4 | z
202!---------------------------------------------------------------------------
203
204CONTAINS
205
206CHARACTER (LEN=2) FUNCTION atompdb(LIN, PDBSTR)
207
208INTEGER :: pa, pb
209CHARACTER (LEN=2), INTENT(IN) :: pdbstr
210INTEGER, INTENT(IN) :: lin
211
212atompdb=' '
213pb=0
214do pa=1, 2
215 if (pdbstr(pa:pa) .ne. ' ') then
216 pb= pb+1
217 atompdb(pb:pb)=pdbstr(pa:pa)
218 endif
219enddo
220
221if (pb .eq. 2) then
222 if (atompdb(2:2).eq.'A') atompdb(2:2)='a'
223 if (atompdb(2:2).eq.'B') atompdb(2:2)='b'
224 if (atompdb(2:2).eq.'C') atompdb(2:2)='c'
225 if (atompdb(2:2).eq.'D') atompdb(2:2)='d'
226 if (atompdb(2:2).eq.'E') atompdb(2:2)='e'
227 if (atompdb(2:2).eq.'F') atompdb(2:2)='f'
228 if (atompdb(2:2).eq.'G') atompdb(2:2)='g'
229 if (atompdb(2:2).eq.'H') atompdb(2:2)='h'
230 if (atompdb(2:2).eq.'I') atompdb(2:2)='i'
231 if (atompdb(2:2).eq.'J') atompdb(2:2)='j'
232 if (atompdb(2:2).eq.'K') atompdb(2:2)='k'
233 if (atompdb(2:2).eq.'L') atompdb(2:2)='l'
234 if (atompdb(2:2).eq.'M') atompdb(2:2)='m'
235 if (atompdb(2:2).eq.'N') atompdb(2:2)='n'
236 if (atompdb(2:2).eq.'O') atompdb(2:2)='o'
237 if (atompdb(2:2).eq.'P') atompdb(2:2)='p'
238 if (atompdb(2:2).eq.'Q') atompdb(2:2)='q'
239 if (atompdb(2:2).eq.'R') atompdb(2:2)='r'
240 if (atompdb(2:2).eq.'S') atompdb(2:2)='s'
241 if (atompdb(2:2).eq.'T') atompdb(2:2)='t'
242 if (atompdb(2:2).eq.'U') atompdb(2:2)='u'
243 if (atompdb(2:2).eq.'V') atompdb(2:2)='v'
244 if (atompdb(2:2).eq.'W') atompdb(2:2)='w'
245 if (atompdb(2:2).eq.'X') atompdb(2:2)='x'
246 if (atompdb(2:2).eq.'Y') atompdb(2:2)='y'
247 if (atompdb(2:2).eq.'Z') atompdb(2:2)='z'
248endif
249
250pb=0
251do pa=1, maxe
252 if (atsym(pa) .eq. atompdb) then
253 pb=1
254 endif
255enddo
256if (pb .eq. 0) then
257 write (6, '("Warning reading PDB file: atom symbol not found line: ",i10)') lin
258endif
259
260END FUNCTION
261
262END FUNCTION
action
Definition glview.h:189
void show_error(char *error, int val, GtkWidget *win)
show error message
Definition interface.c:293
character(len=2) function atompdb(lin, pdbstr)
Definition pdb.F90:207
integer(kind=c_int) function read_pdb(pdb_f, lpdb)
Definition pdb.F90:22