atomes 1.1.16
atomes: an atomic scale modeling tool box
Loading...
Searching...
No Matches
read_trj.c
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
4of the GNU Affero General Public License as published by the Free Software Foundation,
5either 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;
8without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
9See the GNU General Public License for more details.
10
11You should have received a copy of the GNU Affero General Public License along with 'atomes'.
12If not, see <https://www.gnu.org/licenses/>
13
14Copyright (C) 2022-2024 by CNRS and University of Strasbourg */
15
22/*
23* This file: 'read_trj.c'
24*
25* Contains:
26*
27
28 - The functions to read CPMD atomic coordinates
29
30*
31* List of functions:
32
33 int trj_get_atom_coordinates ();
34 int open_trj_file (int linec);
35
36*/
37
38#include "global.h"
39#include "glview.h"
40#include "callbacks.h"
41#include "interface.h"
42#include "project.h"
43#include "bind.h"
44#include "readers.h"
45#ifdef OPENMP
46# include <omp.h>
47#endif
48
55{
56 int i, j, k, l;
57 gchar * lia[4] = {"a", "b", "c", "d"};
59#ifdef OPENMP
60 int res;
61 int numth = omp_get_max_threads ();
62 gboolean doatoms = FALSE;
63 gchar * saved_line;
64 if (active_project -> steps < numth)
65 {
66 if (numth >= 2*(active_project -> steps-1))
67 {
68 doatoms = TRUE;
69 }
70 else
71 {
72 numth = active_project -> steps;
73 }
74 }
75 if (doatoms)
76 {
77 // OpenMP on atoms
78 res = 0;
79 for (i=0; i<active_project -> steps; i++)
80 {
81 k = i*active_project -> natomes;
82 #pragma omp parallel for num_threads(numth) private(j,this_line,saved_line,this_word) shared(i,k,lia,coord_line,active_project,res)
83 for (j=0; j<active_project -> natomes; j++)
84 {
85 if (res == 2) goto enda;
86 this_line = g_strdup_printf ("%s", coord_line[k+j]);
87 saved_line = g_strdup_printf ("%s", this_line);
88 this_word = strtok_r (this_line, " ", & saved_line);
89 if (! this_word)
90 {
91 format_error (i+1, j+1, lia[0], k+j);
92 res = 2;
93 goto enda;
94 }
95 this_word = strtok_r (NULL, " ", & saved_line);
96 if (! this_word)
97 {
98 format_error (i+1, j+1, lia[1], k+j);
99 res = 2;
100 goto enda;
101 }
102 active_project -> atoms[i][j].x = string_to_double ((gpointer)this_word);
103 this_word = strtok_r (NULL, " ", & saved_line);
104 if (! this_word)
105 {
106 format_error (i+1, j+1, lia[2], k+j);
107 res = 2;
108 goto enda;
109 }
110 active_project -> atoms[i][j].y = string_to_double ((gpointer)this_word);
111 this_word = strtok_r (NULL, " ", & saved_line);
112 if (! this_word)
113 {
114 format_error (i+1, j+1, lia[3], k+j);
115 res = 2;
116 goto enda;
117 }
118 active_project -> atoms[i][j].z = string_to_double ((gpointer)this_word);
119 g_free (this_line);
120 enda:;
121 }
122 if (res == 2) break;
123 }
124 }
125 else
126 {
127 res = 0;
128 #pragma omp parallel for num_threads(numth) private(i,j,k,this_line,saved_line,this_word) shared(lia,coord_line,active_project,res)
129 for (i=0; i<active_project -> steps; i++)
130 {
131 if (res == 2) goto ends;
132 k = i*active_project -> natomes;
133 for (j=0; j<active_project -> natomes; j++)
134 {
135 this_line = g_strdup_printf ("%s", coord_line[k+j]);
136 saved_line = g_strdup_printf ("%s", this_line);
137 this_word = strtok_r (this_line, " ", & saved_line);
138 if (! this_word)
139 {
140 format_error (i+1, j+1, lia[0], k+j);
141 res = 2;
142 goto ends;
143 }
144 this_word = strtok_r (NULL, " ", & saved_line);
145 if (! this_word)
146 {
147 format_error (i+1, j+1, lia[1], k+j);
148 res = 2;
149 goto ends;
150 }
151 active_project -> atoms[i][j].x = string_to_double ((gpointer)this_word) * 0.52917721;
152 this_word = strtok_r (NULL, " ", & saved_line);
153 if (! this_word)
154 {
155 format_error (i+1, j+1, lia[2], k+j);
156 res = 2;
157 goto ends;
158 }
159 active_project -> atoms[i][j].y = string_to_double ((gpointer)this_word) * 0.52917721;
160 this_word = strtok_r (NULL, " ", & saved_line);
161 if (! this_word)
162 {
163 format_error (i+1, j+1, lia[3], k+j);
164 res = 2;
165 goto ends;
166 }
167 active_project -> atoms[i][j].z = string_to_double ((gpointer)this_word) * 0.52917721;
168 g_free (this_line);
169 }
170 ends:;
171 }
172 }
173 g_free (coord_line);
174 if (res == 2) return 2;
175 for (i=1; i<active_project -> steps; i++)
176 {
177 for (j=0; j<active_project -> natomes; j++)
178 {
179 if (active_project -> atoms[i-1][j].sp != active_project -> atoms[i][j].sp)
180 {
181 add_reader_info (g_strdup_printf ("Error - chemical species changes between steps %d and %d, for atom %d !", i, i+1, j+1), 0);
182 return 2;
183 }
184 }
185 }
186#else
187 line_node * tmp_line;
188 tail = head;
189 k = 0;
190 for (i=0; i<active_project -> steps; i++)
191 {
192 for (j=0; j<active_project -> natomes; j++)
193 {
194 this_line = g_strdup_printf ("%s", tail -> line);
195 this_word = strtok (this_line, " ");
196 if (! this_word)
197 {
198 format_error (i+1, j+1, lia[0], k);
199 return 2;
200 }
201 this_word = strtok (NULL, " ");
202 if (! this_word)
203 {
204 format_error (i+1, j+1, lia[1], k);
205 return 2;
206 }
207 active_project -> atoms[i][j].x = string_to_double ((gpointer)this_word) * 0.52917721;
208 this_word = strtok (NULL, " ");
209 if (! this_word)
210 {
211 format_error (i+1, j+1, lia[2], k);
212 return 2;
213 }
214 active_project -> atoms[i][j].y = string_to_double ((gpointer)this_word) * 0.52917721;
215 this_word = strtok (NULL, " ");
216 if (! this_word)
217 {
218 format_error (i+1, j+1, lia[3], k);
219 return 2;
220 }
221 active_project -> atoms[i][j].z = string_to_double ((gpointer)this_word) * 0.52917721;
222 tmp_line = tail;
223 tail = tail -> next;
224 g_free (tmp_line);
225 k ++;
226 }
227 }
228#endif
229 i = 0;
230 for (j=0; j<this_reader -> nspec; j++)
231 {
232 for (k=0; k<this_reader -> nsps[j]; k++)
233 {
234 for (l=0; l<active_project -> steps; l++)
235 {
236 active_project -> atoms[l][i].sp = j;
237 }
238 i ++;
239 }
240 }
241 return 0;
242}
243
251int open_trj_file (int linec)
252{
253 if (linec%(this_reader -> natomes) != 0) return 2;
254 reader_info ("trj", "Number of atoms", this_reader -> natomes);
255 active_project -> steps = linec / this_reader -> natomes;
256 reader_info ("trj", "Number of steps", active_project -> steps);
257 active_project -> natomes = this_reader -> natomes;
258 return trj_get_atom_coordinates ();
259}
Binding to the Fortran90 subroutines.
Callback declarations for main window.
void allocatoms(project *this_proj)
allocate project data
Definition open_p.c:160
int atoms[NUM_STYLES][2]
double string_to_double(gpointer string)
convert string to double
Definition global.c:624
Global variable declarations Global convenience function declarations Global data structure defin...
project * active_project
Definition project.c:47
Variable declarations related to the OpenGL window Function declarations related to the OpenGL wind...
Messaging function declarations.
Function declarations for reading atomes project file Function declarations for saving atomes proje...
line_node * head
Definition read_coord.c:75
gchar ** coord_line
Definition read_coord.c:72
coord_file * this_reader
Definition read_coord.c:71
void reader_info(gchar *type, gchar *sinf, int val)
display reader information
Definition read_coord.c:101
line_node * tail
Definition read_coord.c:76
void add_reader_info(gchar *info, int mid)
append information message to the reader information
Definition read_coord.c:86
char * this_word
Definition read_coord.c:74
void format_error(int stp, int ato, gchar *mot, int line)
Message to display an error message.
Definition read_coord.c:116
gchar * this_line
Definition read_coord.c:73
int open_trj_file(int linec)
open CPMD file
Definition read_trj.c:251
int trj_get_atom_coordinates()
get the atomic coordinates from the CPMD file
Definition read_trj.c:54
Functions declaration to read atomic coordinates.
GtkWidget * res[2]
Definition w_encode.c:212