source: util/src/NanoCreator.c@ d97af9

Last change on this file since d97af9 was 21b8fa, checked in by Frederik Heber <heber@…>, 16 years ago

BUGFIX: in main() an else was missing when creating basename from filename.

  • Property mode set to 100644
File size: 64.0 KB
Line 
1#include <stdlib.h>
2#include <stdio.h>
3#include <string.h>
4#include <math.h>
5#include <time.h>
6
7#include <gsl/gsl_matrix.h>
8#include <gsl/gsl_vector.h>
9#include <gsl/gsl_eigen.h>
10#include <gsl/gsl_blas.h>
11
12#include "NanoCreator.h"
13
14
15// ---------------------------------- F U N C T I O N S ----------------------------------------------
16
17// ================================= File functions ==============================
18
19#define LINE_SIZE 80
20#define NDIM 3
21
22/** Allocs memory and prints a message on fail.
23 * \param *size size of alloc in bytes
24 * \param *msg error msg
25 * \return pointer to allocated memory
26 */
27void * Malloc (size_t size, const char *msg)
28{
29 void *ptr = malloc(size);
30 if (ptr == NULL) {
31 if (msg == NULL)
32 fprintf(stdout, "ERROR: Malloc\n");
33 else
34 fprintf(stdout, "ERROR: Malloc - %s\n", msg);
35 return NULL;
36 } else {
37 return ptr;
38 }
39}
40
41/** Callocs memory and prints a message on fail.
42 * \param *size size of alloc in bytes
43 * \param *value pointer to initial value
44 * \param *msg error msg
45 * \return pointer to allocated memory
46 */
47void * Calloc (size_t size, double value, const char *msg)
48{
49 void *ptr = calloc(size, value);
50 if (ptr == NULL) {
51 if (msg == NULL)
52 fprintf(stdout, "ERROR: Calloc\n");
53 else
54 fprintf(stdout, "ERROR: Calloc - %s\n", msg);
55 return NULL;
56 } else {
57 return ptr;
58 }
59}
60
61/** Frees memory only if ptr != NULL.
62 * \param *ptr pointer to array
63 * \param *msg
64 */
65void Free(void * ptr, const char *msg)
66{
67 if (ptr != NULL)
68 free(ptr);
69 else {
70 if (msg == NULL)
71 fprintf(stdout, "ERROR: Free\n");
72 else
73 fprintf(stdout, "ERROR: Free - %s\n", msg);
74 }
75}
76
77/** Allocs memory and prints a message on fail.
78 * \param **old_ptr pointer to old memory range
79 * \param *newsize new size of alloc in bytes
80 * \param *msg error msg
81 * \return pointer to allocated memory
82 */
83void * Realloc (void *old_ptr, size_t newsize, const char *msg)
84{
85 if (old_ptr == NULL) {
86 fprintf(stdout, "ERROR: Realloc - old_ptr NULL\n");
87 exit(255);
88 }
89 void *ptr = realloc(old_ptr, newsize);
90 if (ptr == NULL) {
91 if (msg == NULL)
92 fprintf(stdout, "ERROR: Realloc\n");
93 else
94 fprintf(stdout, "ERROR: Realloc - %s\n", msg);
95 return NULL;
96 } else {
97 return ptr;
98 }
99}
100
101/** Reads a file's contents into a char buffer of appropiate size.
102 * \param *filename name of file
103 * \param pointer to integer holding read/allocated buffer length
104 * \return pointer to allocated segment with contents
105 */
106char * ReadBuffer (char *filename, int *bufferlength)
107{
108 if ((filename == NULL) || (bufferlength == NULL)) {
109 fprintf(stderr, "ERROR: ReadBuffer - ptr to filename or bufferlength NULL\n");
110 exit(255);
111 }
112 char *buffer = Malloc(sizeof(char)*LINE_SIZE, "ReadBuffer: buffer");
113 int i = 0, number;
114 FILE *file = fopen(filename, "r");
115 if (file == NULL) {
116 fprintf(stderr, "File open error: %s", filename);
117 exit(255);
118 }
119 while ((number = fread(&buffer[i], sizeof(char), LINE_SIZE, file)) == LINE_SIZE) {
120 //fprintf(stdout, "%s", &buffer[i]);
121 i+= LINE_SIZE;
122 buffer = (char *) Realloc(buffer, i+LINE_SIZE, "ReadBuffer: buffer");
123 }
124 fclose(file);
125 fprintf(stdout, "Buffer length is %i\n", i);
126 *bufferlength = i+(number);
127 return buffer;
128}
129
130/** Extracts next line out of a buffer.
131 * \param *buffer buffer to parse for newline
132 * \param *line contains complete line on return
133 * \return length of line, 0 if end of file
134 */
135int GetNextline(char *buffer, char *line)
136{
137 if ((buffer == NULL) || (line == NULL)) {
138 fprintf(stderr, "ERROR: GetNextline - ptr to buffer or line NULL\n");
139 exit(255);
140 }
141 int length;
142 char *ptr = strchr(buffer, '\n');
143 //fprintf(stdout, "Newline at %p from %p\n", ptr, buffer);
144 if (ptr == NULL) { // buffer ends
145 return 0;
146 } else {
147 //fprintf(stdout, "length of line is %d\n", length);
148 length = (int)(ptr - buffer)/sizeof(char);
149 strncpy(line, buffer, length);
150 line[length]='\0';
151 return length+1;
152 }
153}
154
155/** Adds commentary stuff (needed for further stages) to Cell xyz files.
156 * \param *filename file name
157 * \param atomicnumner number of atoms in xyz
158 * \param **Vector list of three unit cell vectors
159 * \param **Recivector list of three reciprocal unit cell vectors
160 * \param atomicnumber number of atoms in cell
161 */
162void AddAtomicNumber(char *filename, int atomicnumber, double **Vector, double **Recivector)
163{
164 if ((filename == NULL) || (Vector == NULL) || (Recivector == NULL)) {
165 fprintf(stdout, "ERROR: AddAtomicNumber - ptr to filename, Vector or Recivector NULL\n");
166 exit(255);
167 }
168 int bufferlength;
169 char *buffer = ReadBuffer(filename, &bufferlength);
170 FILE *file2 = fopen(filename, "w+");
171 if (file2 == NULL) {
172 fprintf(stdout, "ERROR: AddAtomicNumber: %s can't open for writing\n", filename);
173 exit(255);
174 }
175 double volume = Determinant(Vector);
176 time_t now;
177
178 now = time((time_t *)NULL); // Get the system time and put it into 'now' as 'calender time'
179 // open for writing and prepend
180 fprintf(file2,"%d\n", atomicnumber); // 2
181 fprintf(file2,"\tgenerated with Nanotube creator on %s", ctime(&now));
182 fwrite(buffer, sizeof(char), bufferlength, file2); // append buffer
183
184 // Add primitive vectors as comment
185 fprintf(file2,"\n****************************************\n\n");
186 fprintf(file2,"Primitive vectors\n");
187 fprintf(file2,"a(1) = %f\t%f\t%f\n", Vector[0][0], Vector[0][1], Vector[0][2]);
188 fprintf(file2,"a(2) = %f\t%f\t%f\n", Vector[1][0], Vector[1][1], Vector[1][2]);
189 fprintf(file2,"a(3) = %f\t%f\t%f\n", Vector[2][0], Vector[2][1], Vector[2][2]);
190 fprintf(file2,"\nVolume = %f", volume);
191 fprintf(file2,"\nReciprocal Vectors\n");
192 fprintf(file2,"b(1) = %f\t%f\t%f\n", Recivector[0][0], Recivector[0][1], Recivector[0][2]);
193 fprintf(file2,"b(2) = %f\t%f\t%f\n", Recivector[1][0], Recivector[1][1], Recivector[1][2]);
194 fprintf(file2,"b(3) = %f\t%f\t%f\n", Recivector[2][0], Recivector[2][1], Recivector[2][2]);
195
196 fclose(file2); // close file
197 Free(buffer, "AddAtomicNumber: buffer");
198}
199
200/** Adds commentary stuff (needed for further stages) to Sheet xyz files.
201 * \param *filename file name
202 * \param *axis array with major, minor and no axis
203 * \param *chiral pointer to array with both chiral values
204 * \param *factors pointer to array with length and radius factor
205 * \param seed random number seed
206 * \param numbercell number of atoms in unit cell, needed as length of \a *randomness
207 * \param *randomness for each atom in unit cell a factor between 0..1, giving its probability of appearance
208 */
209void AddSheetInfo(char *filename, int *axis, int *chiral, int *factors, int seed, int numbercell, double *randomness)
210{
211 int i;
212 if ((filename == NULL) || (axis == NULL) || (chiral == NULL) || (factors == NULL)) {
213 fprintf(stdout, "ERROR: AddSheetInfo - ptr to filename, axis, chiral or factors NULL\n");
214 exit(255);
215 }
216 // open for writing and append
217 FILE *file2 = fopen(filename,"a");
218 if (file2 == NULL) {
219 fprintf(stderr, "ERROR: AddSheetInfo - can't open %s for appending\n", filename);
220 exit(255);
221 }
222 // Add primitive vectors as comment
223 fprintf(file2,"\n****************************************\n\n");
224 fprintf(file2,"Axis %d\t%d\t%d\n", axis[0], axis[1], axis[2]);
225 fprintf(file2,"(n,m) %d\t%d\n", chiral[0], chiral[1]);
226 fprintf(file2,"factors %d\t%d\n", factors[0], factors[1]);
227 fprintf(file2,"seed %d\n", seed);
228 fprintf(file2,"\nRandomness\n");
229 for (i=0; i<numbercell; i++) {
230 fprintf(file2,"%d %g\n", i, randomness[i]);
231 }
232 fclose(file2);
233}
234
235
236// ================================= Vector functions ==============================
237
238/** Transforms a vector b with a matrix A: Ab = x.
239 * \param **matrixref reference to NDIMxNDIM matrix A
240 * \param *vectorref reference to NDIM vector b
241 * \return reference to resulting NDIM vector Ab = x
242 */
243double *MatrixTrafo(double **matrixref, double *vectorref)
244{
245 if ((matrixref == NULL) || (vectorref == NULL)) {
246 fprintf(stderr, "ERROR: MatrixTrafo: ptr to matrixref or vectorref NULL\n");
247 exit(255);
248 }
249 //double *returnvector = Calloc(sizeof(double)*NDIM, 0., "MatrixTrafo: returnvector");
250 double *returnvector = calloc(sizeof(double)*NDIM, 0.);
251 if (returnvector == NULL) {
252 fprintf(stderr, "ERROR: MatrixTrafo - returnvector\n");
253 exit(255);
254 }
255 int i,j;
256
257 for (i=0;i<NDIM;i++)
258 for (j=0;j<NDIM;j++)
259 returnvector[j] += matrixref[i][j] * vectorref[i];
260
261 return returnvector;
262}
263double *MatrixTrafoInverse(double *vectorref, double **matrixref)
264{
265 if ((matrixref == NULL) || (vectorref == NULL)) {
266 fprintf(stderr, "ERROR: MatrixTrafo: ptr to matrixref or vectorref NULL\n");
267 exit(255);
268 }
269 //double *returnvector = Calloc(sizeof(double)*NDIM, 0., "MatrixTrafo: returnvector");
270 double *returnvector = calloc(sizeof(double)*NDIM, 0.);
271 if (returnvector == NULL) {
272 fprintf(stderr, "ERROR: MatrixTrafo - returnvector\n");
273 exit(255);
274 }
275 int i,j;
276
277 for (i=0;i<NDIM;i++)
278 for (j=0;j<NDIM;j++)
279 returnvector[i] += matrixref[i][j] * vectorref[j];
280
281 return returnvector;
282}
283
284/** Inverts a NDIMxNDIM matrix.
285 * \param **matrix to be inverted
286 * \param **inverse allocated space for inverse of \a **matrix
287 */
288void MatrixInversion(double **matrix, double **inverse)
289{
290 if ((matrix == NULL) || (inverse == NULL)) {
291 fprintf(stderr, "ERROR: MatrixInversion: ptr to matrix or inverse NULL\n");
292 exit(255);
293 }
294 // determine inverse
295 double det = Determinant(matrix);
296 inverse[0][0] = (matrix[1][1]*matrix[2][2] - matrix[1][2]*matrix[2][1])/det;
297 inverse[1][0] = (matrix[0][2]*matrix[2][1] - matrix[0][1]*matrix[2][2])/det;
298 inverse[2][0] = (matrix[0][1]*matrix[1][2] - matrix[0][2]*matrix[1][1])/det;
299 inverse[0][1] = (matrix[1][2]*matrix[2][0] - matrix[1][0]*matrix[2][2])/det;
300 inverse[1][1] = (matrix[0][0]*matrix[2][2] - matrix[0][2]*matrix[2][0])/det;
301 inverse[2][1] = (matrix[0][2]*matrix[1][0] - matrix[0][0]*matrix[1][2])/det;
302 inverse[0][2] = (matrix[1][0]*matrix[2][1] - matrix[1][1]*matrix[2][0])/det;
303 inverse[1][2] = (matrix[0][1]*matrix[2][0] - matrix[0][0]*matrix[2][1])/det;
304 inverse[2][2] = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0])/det;
305
306 // check inverse
307 int flag = 0;
308 int i,j,k;
309 double tmp;
310 fprintf(stdout, "Checking inverse ... ");
311 for (i=0;i<NDIM;i++)
312 for (j=0;j<NDIM;j++) {
313 tmp = 0.;
314 for (k=0;k<NDIM;k++)
315 tmp += matrix[i][k]*inverse[j][k];
316 if (!flag) {
317 if (i == j) {
318 flag = (fabs(1.-tmp) > MYEPSILON) ? 1 : 0;
319 } else {
320 flag = (fabs(tmp) > MYEPSILON) ? 1 : 0;
321 }
322 }
323 }
324 if (!flag)
325 fprintf(stdout, "ok\n");
326 else
327 fprintf(stdout, "false\n");
328}
329
330/** Flips to double numbers in place.
331 * \param *number1 pointer to first double
332 * \param *number2 pointer to second double
333 */
334void flip(double *number1, double *number2)
335{
336 if ((number1 == NULL) || (number2 == NULL)) {
337 fprintf(stderr, "ERROR: flip: ptr to number1 or number2 NULL\n");
338 exit(255);
339 }
340 double tmp = *number1;
341 *number1 = *number2;
342 *number2 = tmp;
343}
344
345/** Transposes a matrix.
346 * \param **matrix pointer to NDIMxNDIM-matrix array
347 */
348void Transpose(double **matrix)
349{
350 if (matrix == NULL) {
351 fprintf(stderr, "ERROR: Transpose: ptr to matrix NULL\n");
352 exit(255);
353 }
354 int i,j;
355 for (i=0;i<NDIM;i++)
356 for (j=0;j<i;j++)
357 flip(&matrix[i][j],&matrix[j][i]);
358}
359
360
361/** Computes the determinant of a NDIMxNDIM matrix.
362 * \param **matrix pointer to matrix array
363 * \return det(matrix)
364 */
365double Determinant(double **matrix) {
366 if (matrix == NULL) {
367 fprintf(stderr, "ERROR: Determinant: ptr to Determinant NULL\n");
368 exit(255);
369 }
370 double det = matrix[0][0] * (matrix[1][1]*matrix[2][2] - matrix[1][2]*matrix[2][1])
371 - matrix[1][1] * (matrix[0][0]*matrix[2][2] - matrix[0][2]*matrix[2][0])
372 + matrix[2][2] * (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
373 return det;
374}
375
376/** Adds \a *vector1 onto \a *vector2 coefficient-wise.
377 * \param *vector1 first vector, on return contains sum
378 * \param *vector2 vector which is projected
379 * \return sum of the two vectors
380 */
381double * VectorAdd(double *vector1, double *vector2)
382{
383 if ((vector1 == NULL) || (vector2 == NULL)) {
384 fprintf(stderr, "ERROR: VectorAdd: ptr to vector1 or vector2 NULL\n");
385 exit(255);
386 }
387 //double *returnvector = Calloc(sizeof(double)*NDIM, 0., "VectorAdd: returnvector");
388 double *returnvector = calloc(sizeof(double)*NDIM, 0.);
389 if (returnvector == NULL) {
390 fprintf(stderr, "ERROR: VectorAdd - returnvector\n");
391 exit(255);
392 }
393 int i;
394
395 for (i=0;i<NDIM;i++)
396 returnvector[i] = vector1[i] + vector2[i];
397
398 return returnvector;
399}
400
401/** Fixed GramSchmidt-Orthogonalization for NDIM vectors
402 * \param @orthvector reference to NDIMxNDIM matrix
403 * \param @orthbetrag reference to NDIM vector with vector magnitudes
404 * \param @axis major-, minor- and noaxis for specific order for the GramSchmidt procedure
405 */
406void Orthogonalize(double **orthvector, int *axis)
407{
408 if ((orthvector == NULL) || (axis == NULL)) {
409 fprintf(stderr, "ERROR: Orthogonalize: ptr to orthvector or axis NULL\n");
410 exit(255);
411 }
412 double betrag;
413 int i;
414
415 // first vector is untouched
416 // second vector
417 betrag = Projection(orthvector[axis[1]], orthvector[axis[0]]);
418 fprintf(stdout,"%lg\t",betrag);
419 for (i=0;i<NDIM;i++)
420 orthvector[axis[1]][i] -= orthvector[axis[0]][i] * betrag;
421 // third vector
422 betrag = Projection(orthvector[axis[0]], orthvector[axis[2]]);
423 fprintf(stdout,"%lg\t",betrag);
424 for (i=0;i<NDIM;i++)
425 orthvector[axis[2]][i] -= orthvector[axis[0]][i] * betrag;
426 betrag = Projection(orthvector[axis[1]], orthvector[axis[2]]);
427 fprintf(stdout,"%lg\n",betrag);
428 for (i=0;i<NDIM;i++)
429 orthvector[axis[2]][i] -= orthvector[axis[1]][i] * betrag;
430}
431
432/** Computes projection of \a *vector2 onto \a *vector1.
433 * \param *vector1 reference vector
434 * \param *vector2 vector which is projected
435 * \return projection
436 */
437double Projection(double *vector1, double *vector2)
438{
439 if ((vector1 == NULL) || (vector2 == NULL)) {
440 fprintf(stderr, "ERROR: Projection: ptr to vector1 or vector2 NULL\n");
441 exit(255);
442 }
443 return (ScalarProduct(vector1, vector2)/Norm(vector1)/Norm(vector2));
444}
445
446/** Determine scalar product between two vectors.
447 * \param *vector1 first vector
448 * \param *vector2 second vector
449 * \return scalar product
450 */
451double ScalarProduct(double *vector1, double *vector2)
452{
453 if ((vector1 == NULL) || (vector2 == NULL)) {
454 fprintf(stderr, "ERROR: ScalarProduct: ptr to vector1 or vector2 NULL\n");
455 exit(255);
456 }
457 double scp = 0.;
458 int i;
459
460 for (i=0;i<NDIM;i++)
461 scp += vector1[i] * vector2[i];
462
463 return scp;
464}
465
466/** Computes norm of \a *vector.
467 * \param *vector pointer to NDIM vector
468 * \return norm of \a *vector
469 */
470double Norm(double *vector)
471{
472 if (vector == NULL) {
473 fprintf(stderr, "ERROR: Norm: ptr to vector NULL\n");
474 exit(255);
475 }
476 return sqrt(ScalarProduct(vector, vector));
477}
478
479/** Prints vector to \a *file.
480 * \param *file file or e.g. stdout
481 * \param *vector vector to be printed
482 */
483void PrintVector(FILE *file, double *vector)
484{
485 if ((file == NULL) || (vector == NULL)) {
486 fprintf(stderr, "ERROR: PrintVector: ptr to file or vector NULL\n");
487 exit(255);
488 }
489 int i;
490 for (i=0;i<NDIM;i++)
491 fprintf(file, "%5.5f\t", vector[i]);
492 fprintf(file, "\n");
493}
494
495/** Prints matrix to \a *file.
496 * \param *file file or e.g. stdout
497 * \param **matrix matrix to be printed
498 */
499void PrintMatrix(FILE *file, double **matrix)
500{
501 if ((file == NULL) || (matrix == NULL)) {
502 fprintf(stderr, "ERROR: PrintMatrix: ptr to file or matrix NULL\n");
503 exit(255);
504 }
505 int i,j;
506 for (i=0;i<NDIM;i++) {
507 for (j=0;j<NDIM;j++)
508 fprintf(file, "%5.5f\t", matrix[i][j]);
509 fprintf(file, "\n");
510 }
511}
512
513/** Returns greatest common denominator.
514 * \param a first integer
515 * \param b second integer
516 * \return GCD of a and b
517 */
518int GCD(int a, int b)
519{
520 int c;
521 do {
522 c = a % b; /* Rest of integer divison */
523 a = b; b = c; /* flip the two values */
524 } while( c != 0);
525 return a;
526}
527
528/** Determines the biggest diameter of a sheet.
529 * \param **matrix reference to NDIMxNDIM matrix with row vectors
530 * \param *axis reference to NDIM vector with permutation of axis indices [0,1,2]
531 * \param *factors factorsfor axis[0] and axis[1]
532 * \return biggest diameter of sheet
533*/
534double DetermineBiggestDiameter(double **matrix, int *axis, int *factors)
535{
536 if ((axis == NULL) || (factors == NULL) || (matrix == NULL)) {
537 fprintf(stderr, "ERROR: DetermineBiggestDiameter: ptr to factors, axis or matrix NULL\n");
538 exit(255);
539 }
540 double diameter[2] = {0.,0.};
541 int i, biggest;
542
543 for (i=0;i<NDIM;i++) {
544 diameter[0] += (matrix[axis[0]][i]*factors[0] - matrix[axis[1]][i]*factors[1]) * (matrix[axis[0]][i]*factors[0] - matrix[axis[1]][i]*factors[1]);
545 diameter[1] += (matrix[axis[0]][i]*factors[0] + matrix[axis[1]][i]*factors[1]) * (matrix[axis[0]][i]*factors[0] + matrix[axis[1]][i]*factors[1]);
546 }
547 if ((diameter[0] - diameter[1]) > MYEPSILON) {
548 biggest = 0;
549 } else {
550 biggest = 1;
551 }
552 diameter[0] = sqrt(diameter[0]);
553 diameter[1] = sqrt(diameter[1]);
554 fprintf(stdout, "\n\nMajor diameter of the sheet is %5.5f, minor diameter is %5.5f.\n",diameter[biggest],diameter[(biggest+1)%2]);
555
556 return diameter[biggest];
557}
558
559/** Determines the center of gravity of atoms in a buffer \a bufptr with given \a number
560 * \param *bufptr pointer to char buffer with atoms in (name x y z)-manner
561 * \param number number of atoms/lines to scan
562 * \return NDIM vector (allocated doubles) pointing back to center of gravity
563 */
564double * CenterOfGravity(char *bufptr, int number)
565{
566 if (bufptr == NULL) {
567 fprintf(stderr, "ERROR: CenterOfGravity - bufptr NULL\n");
568 exit(255);
569 }
570 double *cog = calloc(sizeof(double)*NDIM, 0.);
571 if (cog == NULL) {
572 fprintf(stderr, "ERROR: CenterOfGravity - cog\n");
573 exit(255);
574 }
575 double *atom = Malloc(sizeof(double)*NDIM, "CenterOfGravity: atom");
576 char name[255], line[255];
577 int i,j;
578
579 // Determine center of gravity
580 for (i=0;i<number;i++) {
581 bufptr += GetNextline(bufptr, line);
582 sscanf(line, "%s %lg %lg %lg", name, &atom[0], &atom[1], &atom[2]);
583 //fprintf(stdout, "Read Atom %s %lg %lg %lg\n", name, atom[0], atom[1], atom[2]);
584 for (j=0;j<NDIM;j++)
585 cog[j] += atom[j];
586 }
587 for (i=0;i<NDIM;i++)
588 cog[i] /= -number;
589
590 Free(atom, "CenterOfGravity: atom");
591 return cog;
592}
593
594/** Creates orthogonal vectors to directions axis[0] and axis[1], assuming that axis[2] is always orthogonal.
595 * \param **OrthoVector contains vectors set on return with axis[2] equal to Vector[axis[2]]
596 * \param **Vector vectors to orthogonalize against
597 * \param *axis lookup for which direction is which.
598 */
599void CreateOrthogonalAxisVectors(double **OrthoVector, double **Vector, int *axis) {
600 int i,j;
601 double factor;
602 // allocate memory
603 int *TempAxis = (int *) Malloc(sizeof(int)*NDIM, "Main: *TempAxis");
604 double **TempVectors = (double **) Malloc(sizeof(double *)*NDIM, "Main: *TempVectors");
605 for (i=0; i<NDIM; i++)
606 TempVectors[i] = (double *) Malloc(sizeof(double)*NDIM, "Main: TempVectors");
607
608 for (i=0; i<NDIM; i++)
609 for (j=0; j<NDIM; j++)
610 TempVectors[i][j] = Vector[i][j];
611 // GramSchmidt generates Vector[1] orthogonal to Vector[0] and Vector[2] ortho. to Vector[1] and Vector[0]!
612 TempAxis[0] = axis[2]; // (axis 2 is the orthogonal plane axis)
613 TempAxis[1] = axis[0];
614 TempAxis[2] = axis[1];
615 Orthogonalize(TempVectors, TempAxis);
616 factor = Norm(Vector[axis[0]])/Norm(TempVectors[TempAxis[2]]);
617 factor *= (Projection(TempVectors[TempAxis[2]], Vector[axis[0]]) > 0) ? 1. : -1.;
618 for (i=0; i<NDIM; i++)
619 OrthoVector[axis[0]][i] = TempVectors[TempAxis[2]][i]*factor;
620
621 TempAxis[1] = axis[1];
622 TempAxis[2] = axis[0];
623 for (i=0; i<NDIM; i++)
624 for (j=0; j<NDIM; j++)
625 TempVectors[i][j] = Vector[i][j];
626 Orthogonalize(TempVectors, TempAxis);
627 factor = Norm(Vector[axis[1]])/Norm(TempVectors[TempAxis[2]]);
628 factor *= (Projection(TempVectors[TempAxis[2]], Vector[axis[1]]) > 0) ? 1. : -1.;
629 for (i=0; i<NDIM; i++)
630 OrthoVector[axis[1]][i] = TempVectors[TempAxis[2]][i]*factor;
631
632 for (i=0; i<NDIM; i++)
633 OrthoVector[axis[2]][i] = Vector[axis[2]][i];
634 //print vectors
635 fprintf(stdout, "Orthogonal vectors are: \n");
636 for (i=0; i<NDIM; i++) {
637 for (j=0; j<NDIM; j++)
638 fprintf(stdout, "%lg\t", OrthoVector[axis[i]][j]);
639 fprintf(stdout, "\n");
640 }
641
642 // free memory
643 Free(TempAxis, "CreateOrthogonalAxisVectors: *TempAxis");
644 for (i=0; i<NDIM; i++ )
645 Free(TempVectors[i], "CreateOrthogonalAxisVectors: TempVectors");
646 Free(TempVectors, "CreateOrthogonalAxisVectors: *TempVectors");
647};
648
649// ================================= other functions ==============================
650
651/** Prints a debug message.
652 * \param *msg debug message
653 */
654void Debug(char *msg)
655{
656 if (msg == NULL) {
657 fprintf(stderr, "ERROR: Debug: ptr to msg NULL\n");
658 exit(255);
659 }
660 fprintf(stdout, "%s", msg);
661}
662
663
664// --------------------------------------- M A I N ---------------------------------------------------
665int main(int argc, char **argv) {
666 char *filename = NULL;
667 char *CellFilename = NULL, *SheetFilename = NULL, *TubeFilename = NULL, *TorusFilename = NULL;
668 char *SheetFilenameAligned = NULL, *TubeFilenameAligned = NULL;
669 double **Vector, **OrthoVector, **Recivector = NULL, **Tubevector = NULL, **TubevectorInverse = NULL;
670 double *atom = NULL, *atom_transformed = NULL;
671 double *x = NULL, *coord = NULL, *tempvector = NULL, *offset = NULL;
672 //double *cog = NULL, *cog_projected = NULL;
673 char stage[6];
674 int axis[NDIM] = {0,1,2}, chiral[2] = {1,1}, factors[NDIM] = {1,1,1};
675 char name[255], line[255], input = 'n';
676 char *CellBuffer = NULL, *SheetBuffer = NULL, *TubeBuffer = NULL, *bufptr = NULL;
677 double *randomness = NULL, percentage; // array with percentages for presence in sheet and beyond
678 int i,j, ggT;
679 int length;
680
681 // Read command line arguments
682 if (argc <= 2) {
683 fprintf(stdout, "Usage: %s <file> <stage>\n\tWhere <file> specifies a file to start from <stage> or a basename\n\t<stage> is either None, Cell, Sheet, Tube, Torus and specifies where to start the rolling up from.\n\tNote: The .Aligned. files can't be used (rotation is essential).\n", argv[0]);
684 exit(255);
685 } else {
686 // store variables
687 filename = argv[1];
688 strncpy(stage, argv[2], 5);
689 fprintf(stdout, "I will begin at stage %s.\n", stage);
690
691 // build filenames
692 char *ptr = strrchr(filename, '.');
693 if (ptr == NULL) {
694 ptr = filename;
695 } else {
696 length = strlen(filename);
697 if (ptr != NULL) {
698 length = ((int)(ptr-filename) >= length-5) ? (int)(ptr-filename) : length;
699 filename[length] = '\0'; // eventueller
700 }
701 }
702 fprintf(stdout, "I will use \'%s' as base for the filenames.\n\n", filename);
703 CellFilename = Malloc(sizeof(char)*(length+10), "Main: CellFilename");
704 SheetFilename = Malloc(sizeof(char)*(length+11), "Main: SheetFilename");
705 TubeFilename = Malloc(sizeof(char)*(length+10), "Main: TubeFilename");
706 TorusFilename = Malloc(sizeof(char)*(length+11), "Main: TorusFilename");
707 SheetFilenameAligned = Malloc(sizeof(char)*(length+20), "Main: SheetFilenameAligned");
708 TubeFilenameAligned = Malloc(sizeof(char)*(length+19), "Main: TubeFilenameAligned");
709 sprintf(CellFilename, "%s.Cell.xyz", filename);
710 sprintf(SheetFilename, "%s.Sheet.xyz", filename);
711 sprintf(TubeFilename, "%s.Tube.xyz", filename);
712 sprintf(TorusFilename, "%s.Torus.xyz", filename);
713 sprintf(SheetFilenameAligned, "%s.Sheet.Aligned.xyz", filename);
714 sprintf(TubeFilenameAligned, "%s.Tube.Aligned.xyz", filename);
715 }
716
717 // Allocating memory
718 Debug ("Allocating memory\n");
719 atom = (double *) Malloc(sizeof(double)*NDIM, "Main: atom");
720 Vector = (double **) Malloc(sizeof(double *)*NDIM, "Main: *Vector");
721 OrthoVector = (double **) Malloc(sizeof(double *)*NDIM, "Main: *OrthoVector");
722 Recivector = (double **) Malloc(sizeof(double *)*NDIM, "Main: *Recivector");
723 Tubevector = (double **) Malloc(sizeof(double *)*NDIM, "Main: *Tubevector");
724 TubevectorInverse = (double **) Malloc(sizeof(double *)*NDIM, "Main: *TubevectorInverse");
725 for (i=0; i<NDIM; i++ ) {
726 Vector[i] = (double *) Malloc(sizeof(double)*NDIM, "Main: Vector");
727 OrthoVector[i] = (double *) Malloc(sizeof(double)*NDIM, "Main: OrthoVector");
728 Recivector[i] = (double *) Malloc(sizeof(double)*NDIM, "Main: Recivector");
729 Tubevector[i] = (double *) Malloc(sizeof(double)*NDIM, "Main: Tubevector");
730 TubevectorInverse[i] = (double *) Malloc(sizeof(double)*NDIM, "Main: TubevectorInverse");
731 }
732
733 // ======================== STAGE: Cell ==============================
734 // The cell is simply created by transforming relative coordinates within the cell
735 // into cartesian ones using the unit cell vectors.
736
737 double volume;
738 int numbercell;
739 FILE *CellFile;
740
741 Debug ("STAGE: None\n");
742 // Read cell vectors from stdin or from file
743 if (!strncmp(stage, "Non", 3)) {
744 fprintf(stdout, "You will give the unit cell of the given substance.\nAfterwards, the programme will create a Sheet, a Tube and a Torus, each with their own xyz-file named accordingly.\n\n");
745 fprintf(stdout, "Enter 1st primitive vector: ");
746 fscanf(stdin, "%lg %lg %lg", &Vector[0][0], &Vector[0][1], &Vector[0][2]);
747 fprintf(stdout, "Enter 2nd primitive vector: ");
748 fscanf(stdin, "%lg %lg %lg", &Vector[1][0], &Vector[1][1], &Vector[1][2]);
749 fprintf(stdout, "Enter 3rd primitive vector: ");
750 fscanf(stdin, "%lg %lg %lg", &Vector[2][0], &Vector[2][1], &Vector[2][2]);
751 fprintf(stdout,"Unit vectors are\n");
752 PrintMatrix(stdout, Vector);
753 } else {
754 char *ptr = NULL;
755 char dummy[10];
756 CellBuffer = bufptr = ReadBuffer(CellFilename, &length);
757 for (i=0;i<NDIM;i++) {
758 sprintf(dummy, "a(%i) = ", i+1);
759 fprintf(stdout, "%s", dummy);
760 while ((length = GetNextline(bufptr, line)) != -1) {
761 bufptr += (length)*sizeof(char);
762 //fprintf(stdout, "LINE at %p: %s\n", bufptr, line);
763 if ((ptr = strstr(line, dummy)) != NULL)
764 break;
765 }
766 ptr += strlen(dummy);
767 sscanf(ptr, "%lg %lg %lg", &Vector[i][0], &Vector[i][1], &Vector[i][2]);
768 fprintf(stdout, "%5.5lg %5.5lg %5.5lg\n", Vector[i][0], Vector[i][1], Vector[i][2]);
769 }
770 }
771
772 volume = Determinant(Vector);
773 fprintf(stdout,"Volume is %lg\n", volume);
774 MatrixInversion(Vector, Recivector);
775 //Transpose(Recivector); // inverse's got row vectors if normal matrix' got column ones
776 fprintf(stdout, "Reciprocal vector is ");
777 PrintMatrix(stdout, Recivector);
778 fprintf(stdout, "Reciprocal volume is %lg\n", Determinant(Recivector));
779
780 fprintf(stdout, "Vector magnitudes: %5.5lg %5.5lg %5.5lg\n", Norm(Vector[0]), Norm(Vector[1]), Norm(Vector[2]));
781
782 Debug ("STAGE: Cell\n");
783 if (!strncmp(stage, "Non", 3)) {
784 fprintf(stdout, "\nHow many atoms are in the unit cell: ");
785 fscanf(stdin, "%i", &numbercell);
786 CellFile = fopen(CellFilename, "w");
787 if (CellFile == NULL) {
788 fprintf(stderr, "ERROR: main - can't open %s for writing\n", CellFilename);
789 exit(255);
790 }
791 fprintf(stdout, "\nNext, you have to enter each atom in the cell as follows, e.g.\n");
792 fprintf(stdout, "Enter \'ChemicalSymbol X Y Z\' (relative to primitive vectors): C 0.5 0.25 0.5\n\n");
793 for (i = 0; i < numbercell; i++) {
794 fprintf(stdout, "Enter for atom %i \'ChemicalSymbol X Y Z\': ", i+1);
795 fscanf(stdin, "%s %lg %lg %lg", name, &atom[0], &atom[1], &atom[2]);
796 tempvector = MatrixTrafo(Vector, atom);
797 fprintf(stdout, "Atom %i: %s %5.5lg %5.5lg %5.5lg\n", i, name, tempvector[0], tempvector[1], tempvector[2]);
798 fprintf(stdout, "Probe: %s %5.5lg %5.5lg %5.5lg\n", name,
799 atom[0]*Vector[0][0]+atom[1]*Vector[1][0]+atom[2]*Vector[2][0],
800 atom[0]*Vector[0][1]+atom[1]*Vector[1][1]+atom[2]*Vector[2][1],
801 atom[0]*Vector[0][2]+atom[1]*Vector[1][2]+atom[2]*Vector[2][2]
802 );
803 fprintf(CellFile, "%s %lg %lg %lg\n", name, tempvector[0], tempvector[1], tempvector[2]);
804 Free(tempvector, "Main: At stage Cell - tempvector");
805 }
806 fflush(CellFile);
807 fclose(CellFile);
808 AddAtomicNumber(CellFilename, numbercell, Vector, Recivector);
809
810 CellBuffer = ReadBuffer(CellFilename, &length);
811
812 sprintf(stage, "Cell");
813 } else {
814 bufptr = CellBuffer;
815 GetNextline(bufptr, line);
816 sscanf(line, "%i", &numbercell);
817 }
818
819 fprintf(stdout, "\nThere are %i atoms in the cell.\n", numbercell);
820
821 // ======================== STAGE: Sheet =============================
822 // The sheet is a bit more complex. We read the cell in cartesian coordinates
823 // from the file. Next, we have to rotate the unit cell vectors by the so called
824 // chiral angle. This gives a slanted and elongated section upon the sheet of
825 // periodically repeated original unit cells. It only matches up if the factors
826 // were all integer! (That's why the rotation is discrete and the chiral angle
827 // specified not as (cos alpha, sin alpha) but as (n,m)) Also, we want this
828 // section to be rectangular, thus we orthogonalize the original unit vectors
829 // to gain our (later) tube axis.
830 // By looking at the biggest possible diameter we know whose original cells to
831 // look at and check if their respective compounds (contained atoms) still reside
832 // in the rotated, elongated section we need for the later tube.
833 // Then in a for loop we go through every cell. By projecting the vector leading
834 // from the origin to the specific atom down onto the major and minor axis we
835 // know if it's still within the boundaries spanned by these rotated and elongated
836 // (radius-, length factor) unit vectors or not. If yes, its coordinates are
837 // written to file.
838
839 int numbersheet, biggestdiameter, sheetnr[NDIM], tmp, seed;
840 double x1,x2,x3, angle;
841 char flag = 'n';
842 FILE *SheetFile = NULL;
843 FILE *SheetFileAligned = NULL;
844
845 Debug ("STAGE: Sheet\n");
846 if (!strncmp(stage, "Cell", 4)) {
847 fprintf(stdout, "\nEnter seed unequal 0 if any of the bonds shall have a randomness in their being: ");
848 fscanf(stdin, "%d", &seed);
849 if (seed != 0)
850 input = 'y';
851 randomness = (double *) Malloc(sizeof(double)*numbercell, "Main: at sheet - randomness");
852 for(i=0;i<numbercell;i++)
853 randomness[i] = 0.;
854 i = 0;
855 fprintf(stdout, "\n");
856 while (input == 'y') {
857 fprintf(stdout, "Enter atom number (-1 0 to end) and percentage (0.0...1.0): ");
858 fscanf(stdin, "%d %lg", &i, &percentage);
859 if (i == -1) { input = 'n'; fprintf(stdout, "Breaking\n"); }
860 else { randomness[i] = 1.-percentage; }
861 };
862
863 fprintf(stdout, "\nSpecify the axis permutation that is going to be perpendicular to the sheet [tubeaxis, torusaxis, noaxis]: ");
864 fscanf(stdin, "%d %d %d", &axis[0], &axis[1], &axis[2]);
865 fprintf(stdout, "axis: %d %d %d\n", axis[0], axis[1], axis[2]);
866
867 // create orthogonal vectors individually for each unit cell vector
868 CreateOrthogonalAxisVectors(OrthoVector, Vector, axis);
869 fprintf(stdout, "Orthogonal vector axis[0] %lg\n", Projection(Vector[axis[0]], OrthoVector[axis[0]]));
870 fprintf(stdout, "Orthogonal vector axis[1] %lg\n", Projection(Vector[axis[1]], OrthoVector[axis[1]]));
871
872 do {
873 fprintf(stdout, "\nNow specify the two natural numbers (m n) defining the chiral angle, \nif the result is crap, try flipping to (m,n): ");
874 fscanf(stdin, "%d %d", &chiral[0], &chiral[1]);
875 ggT = GCD(2*chiral[1]+chiral[0],2*chiral[0]+chiral[1]);
876 fprintf(stdout, "Greatest Common Denominator of (2n+m, 2m+n) is %d\n", ggT);
877 fprintf(stdout, "chiral0: %d\tchiral1: %d\n", chiral[0], chiral[1]);
878 for (i=0;i<NDIM;i++) {
879 Tubevector[axis[0]][i] = (double)chiral[0] * Vector[axis[0]][i] + (double)chiral[1] * Vector[axis[1]][i];
880 //Tubevector[axis[0]][i] = chiral[0] * Vector[axis[0]][i] + chiral[1] * Vector[axis[1]][i];
881 //Tubevector[axis[0]][i] = (2.*chiral[0]+chiral[1])/(double)ggT * Vector[axis[0]][i] + (-chiral[0]-2.*chiral[1])/(double)ggT * Vector[axis[1]][i];
882 //Tubevector[axis[1]][i] = -chiral[1] * Vector[axis[0]][i] + chiral[0] * Vector[axis[1]][i];
883 Tubevector[axis[1]][i] = (double)chiral[0] * OrthoVector[axis[0]][i] - (double)chiral[1] * OrthoVector[axis[1]][i];
884 //Tubevector[axis[1]][i] = (-chiral[0]-2.*chiral[1])/(double)ggT * Vector[axis[0]][i] + (2.*chiral[0]+chiral[1])/(double)ggT * Vector[axis[1]][i];
885// fprintf(stderr, "Tubevector[axis[0]][i] = (double)chiral[0] * Vector[axis[0]][i] + (double)chiral[1] * Vector[axis[1]][i]\n = %lg * %lg + %lg * %lg = %lg + %lg = %lg\n\n",
886// (double)chiral[0], Vector[axis[0]][i], (double)chiral[1], Vector[axis[1]][i],
887// (double)chiral[0] * Vector[axis[0]][i], (double)chiral[1] * Vector[axis[1]][i],
888// Tubevector[axis[0]][i]);
889 Tubevector[axis[2]][i] = Vector[axis[2]][i];
890 }
891 // here we assume, that Vector[axis[2]] is along z direction!
892 gsl_matrix *M = gsl_matrix_alloc(2,2);
893 gsl_matrix *C = gsl_matrix_alloc(2,2);
894 gsl_matrix *evec = gsl_matrix_alloc(2,2);
895 gsl_vector *eval = gsl_vector_alloc(2);
896 gsl_vector *v = gsl_vector_alloc(2);
897 gsl_vector *u = gsl_vector_alloc(2);
898 gsl_eigen_symmv_workspace *w = gsl_eigen_symmv_alloc(2);
899 gsl_matrix_set(C, 0,0, Vector[axis[0]][0]);
900 gsl_matrix_set(C, 1,0, Vector[axis[0]][1]);
901 gsl_matrix_set(C, 0,1, Vector[axis[1]][0]);
902 gsl_matrix_set(C, 1,1, Vector[axis[1]][1]);
903 gsl_blas_dgemm(CblasTrans,CblasNoTrans, 1.0, C, C, 0.0, M);
904 fprintf(stdout, "M: \t%lg\t%lg\n\t%lg\t%lg\n", gsl_matrix_get(M,0,0), gsl_matrix_get(M,0,1), gsl_matrix_get(M,1,0), gsl_matrix_get(M,1,1));
905 gsl_eigen_symmv(M, eval, evec, w);
906 gsl_eigen_symmv_sort(eval,evec,GSL_EIGEN_SORT_ABS_DESC);
907 fprintf(stdout, "Eigenvalues: %lg\t%lg\n", gsl_vector_get(eval,0), gsl_vector_get(eval,1));
908 fprintf(stdout, "Eigenvectors: \t%lg\t%lg\n\t\t%lg\t%lg\n", gsl_matrix_get(evec,0,0), gsl_matrix_get(evec,0,1), gsl_matrix_get(evec,1,0), gsl_matrix_get(evec,1,1));
909 gsl_matrix_set(M, 0,0, 0.);
910 gsl_matrix_set(M, 1,0, 1.);
911 gsl_matrix_set(M, 0,1, -gsl_vector_get(eval,1)/gsl_vector_get(eval,0));
912 gsl_matrix_set(M, 1,1, 0.);
913 gsl_vector_set(v,0,(double)chiral[0]);
914 gsl_vector_set(v,1,(double)chiral[1]);
915 gsl_blas_dgemm(CblasNoTrans,CblasNoTrans, 1.0, evec, M, 0.0, C);
916 gsl_blas_dgemm(CblasNoTrans,CblasTrans, 1.0, C, evec, 0.0, M);
917 fprintf(stdout, "M: \t%lg\t%lg\n\t%lg\t%lg\n", gsl_matrix_get(M,0,0), gsl_matrix_get(M,0,1), gsl_matrix_get(M,1,0), gsl_matrix_get(M,1,1));
918 gsl_blas_dgemv(CblasNoTrans, 1.0, M, v, 0.0, u);
919 fprintf(stdout, "Looking for factor to integer...\n");
920 for(i=1;i<(chiral[0]+chiral[1])*(chiral[0]+chiral[1]);i++) {
921 x1 = gsl_vector_get(u,0)*(double)i;
922 x2 = gsl_vector_get(u,1)*(double)i;
923 x3 =
924 fprintf(stdout, "%d: %d\t%d vs. %lg\t%lg\n",i, ((int)(x1+x1/fabs(x1)*.5)), ((int)(x2+x2/fabs(x2)*.5)), (x1), (x2));
925 if (( fabs( ((int)(x1+x1/fabs(x1)*.5)) - (x1) ) < 1e-6) && ( fabs( ((int)(x2+x2/fabs(x2)*.5)) - (x2) ) < 1e-6 )) {
926 gsl_blas_dscal((double)i, u);
927 break;
928 }
929 }
930 fprintf(stdout, "(c,d) = (%lg,%lg)\n",gsl_vector_get(u,0), gsl_vector_get(u,1));
931
932 // get length
933 double x[NDIM];
934 for (i=0;i<NDIM;i++)
935 x[i] = gsl_vector_get(u,0) * Vector[axis[0]][i] + gsl_vector_get(u,1) * Vector[axis[1]][i];
936 angle = Norm(x)/Norm(Tubevector[axis[1]]) ;//ScalarProduct(x,Tubevector[axis[1]])/Norm(Tubevector[axis[1]]);
937 fprintf(stdout, "angle is %lg\n", angle);
938 for (i=0;i<NDIM;i++) {
939 Tubevector[axis[1]][i] = gsl_vector_get(u,0) * Vector[axis[0]][i] + gsl_vector_get(u,1) * Vector[axis[1]][i];
940 }
941
942 // Probe
943 gsl_matrix_set(M, 0,0, Vector[axis[0]][0]);
944 gsl_matrix_set(M, 1,0, Vector[axis[0]][1]);
945 gsl_matrix_set(M, 0,1, Vector[axis[1]][0]);
946 gsl_matrix_set(M, 1,1, Vector[axis[1]][1]);
947 gsl_vector_set(v,0,(double)chiral[0]);
948 gsl_vector_set(v,1,(double)chiral[1]);
949 gsl_blas_dgemv(CblasNoTrans, 1.0, M, u, 0.0, eval);
950 gsl_blas_dgemv(CblasNoTrans, 1.0, M, v, 0.0, u);
951 x1=1.;
952 gsl_blas_ddot(u,eval,&x1);
953 fprintf(stdout, "Testing (c,d): (a,b) M^t M (c,d)^t = 0 ? : %lg\n", x1);
954
955 gsl_matrix_free(M);
956 gsl_matrix_free(C);
957 gsl_matrix_free(evec);
958 gsl_vector_free(eval);
959 gsl_vector_free(v);
960 gsl_vector_free(u);
961 gsl_eigen_symmv_free(w);
962
963 if (fabs(x1) > 1e-6) {
964 fprintf(stderr,"Resulting TubeVectors of axis %d and %d and not orthogonal, aborting.\n", axis[0], axis[1]);
965 return(128);
966 }
967
968
969 angle = Projection(Tubevector[axis[1]], Vector[axis[0]]);
970 fprintf(stdout, "Projection Tubevector1 axis[0] %lg %lg\n", angle, 1./angle);
971 angle = Projection(Tubevector[axis[1]], Vector[axis[1]]);
972 fprintf(stdout, "Projection Tubevector1 axis[1] %lg %lg\n", angle, 1./angle);
973
974/* fprintf(stdout, "Vector\n");
975 PrintMatrix(stdout, Vector);
976 fprintf(stdout, "Tubevector\n");
977 PrintMatrix(stdout, Tubevector);
978 for (i=0;i<NDIM;i++) {
979 fprintf(stdout, "Tubevector %d in Unit cell vectors:\t", axis[i]);
980 tempvector = MatrixTrafoInverse(Tubevector[axis[i]], Recivector);
981 PrintVector(stdout, tempvector);
982 Free(tempvector, "Main:tempvector");
983 }*/
984
985 // Give info for length and radius factors
986 fprintf(stdout, "\nThe chiral angle then is %lg degrees with tube radius %5.5f A and length %5.5f A, i.e. torus radius of %5.5f A.\n",
987 acos(Projection(Vector[axis[0]], Tubevector[axis[0]]))/M_PI*180.,
988 Norm(Tubevector[axis[0]])/(2.*M_PI),
989 Norm(Tubevector[axis[1]]),
990 Norm(Tubevector[axis[1]])/(2.*M_PI)
991 );
992 fprintf(stdout, "\nGive integer factors for length and radius of tube (multiple of %d suggested) : ", ggT);
993 fscanf(stdin, "%d %d", &factors[1], &factors[0]);
994 fprintf(stdout, "\nThe chiral angle then is %5.5f degrees with tube radius %5.5f A and length %5.5f A, i.e. torus radius of %5.5f A.\n",
995 acos(Projection(Vector[axis[0]], Tubevector[axis[0]]))/M_PI*180.,
996 (double)factors[0]*Norm(Tubevector[axis[0]])/(2.*M_PI),
997 (double)factors[1]*Norm(Tubevector[axis[1]]),
998 (double)factors[1]*Norm(Tubevector[axis[1]])/(2.*M_PI)
999 );
1000 fprintf(stdout, "Satisfied? [yn] ");
1001 fscanf(stdin, "%c", &flag);
1002 fscanf(stdin, "%c", &flag);
1003 } while (flag != 'y');
1004 } else {
1005 char *ptr = NULL;
1006 char dummy[10];
1007 double dummydouble;
1008
1009 SheetBuffer = bufptr = ReadBuffer(SheetFilename, &length);
1010 bufptr += (GetNextline(bufptr, line))*sizeof(char);
1011 sscanf(line, "%d", &numbersheet);
1012
1013 // retrieve axis permutation
1014 sprintf(dummy, "Axis");
1015 fprintf(stdout, "%s ", dummy);
1016 while ((length = GetNextline(bufptr, line)) != 0) {
1017 bufptr += (length)*sizeof(char);
1018 if ((ptr = strstr(line, dummy)) != NULL)
1019 break;
1020 }
1021 if (length == 0) {
1022 fprintf(stderr, "ERROR: Main at stage Sheet - could not find %s in %s\n", dummy, SheetFilename);
1023 exit(255);
1024 }
1025 ptr += strlen(dummy);
1026 sscanf(ptr, "%d %d %d", &axis[0], &axis[1], &axis[2]);
1027 fprintf(stdout, "%d %d %d\n", axis[0], axis[1], axis[2]);
1028
1029 // retrieve chiral numbers
1030 sprintf(dummy, "(n,m)");
1031 fprintf(stdout, "%s ", dummy);
1032 while ((length = GetNextline(bufptr, line)) != 0) {
1033 bufptr += (length)*sizeof(char);
1034 if ((ptr = strstr(line, dummy)) != NULL)
1035 break;
1036 }
1037 if (length == 0) {
1038 fprintf(stderr, "ERROR: Main at stage Sheet - could not find %s in %s\n", dummy, SheetFilename);
1039 exit(255);
1040 }
1041 ptr += strlen(dummy);
1042 sscanf(ptr, "%d %d", &chiral[0], &chiral[1]);
1043 fprintf(stdout, "%d %d\n", chiral[0], chiral[1]);
1044 ggT = GCD(2*chiral[1]+chiral[0],2*chiral[0]+chiral[1]);
1045 fprintf(stdout, "Greatest Common Denominator of (2n+m, 2m+n) is %d\n", ggT);
1046
1047 // retrieve length and radius factors
1048 sprintf(dummy, "factors");
1049 fprintf(stdout, "%s ", dummy);
1050 while ((length = GetNextline(bufptr, line)) != 0) {
1051 bufptr += (length)*sizeof(char);
1052 if ((ptr = strstr(line, dummy)) != NULL)
1053 break;
1054 }
1055 if (length == 0) {
1056 fprintf(stderr, "ERROR: Main at stage Sheet - could not find %s in %s\n", dummy, SheetFilename);
1057 exit(255);
1058 }
1059 ptr += strlen(dummy);
1060 sscanf(ptr, "%d %d %d", &factors[0], &factors[1], &factors[2]);
1061 fprintf(stdout, "%d %d %d\n", factors[0], factors[1], factors[2]);
1062
1063 // create orthogonal vectors individually for each unit cell vector
1064 CreateOrthogonalAxisVectors(OrthoVector, Vector, axis);
1065 fprintf(stdout, "Orthogonal vector axis[0] %lg", Projection(Vector[axis[0]], OrthoVector[axis[0]]));
1066 fprintf(stdout, "Orthogonal vector axis[1] %lg", Projection(Vector[axis[1]], OrthoVector[axis[1]]));
1067 // create Tubevectors
1068 for (i=0;i<NDIM;i++) {
1069 Tubevector[axis[0]][i] = chiral[0] * Vector[axis[0]][i] + chiral[1] * Vector[axis[1]][i];
1070 //Tubevector[axis[0]][i] = (2.*chiral[0]+chiral[1])/(double)ggT * Vector[axis[0]][i] + (-chiral[0]-2.*chiral[1])/(double)ggT * Vector[axis[1]][i];
1071 //Tubevector[axis[1]][i] = -chiral[1] * Vector[axis[0]][i] + chiral[0] * Vector[axis[1]][i];
1072 Tubevector[axis[1]][i] = chiral[0] * OrthoVector[axis[0]][i] + chiral[1] * OrthoVector[axis[1]][i];
1073 //Tubevector[axis[1]][i] = (-chiral[0]-2.*chiral[1])/(double)ggT * Vector[axis[0]][i] + (2.*chiral[0]+chiral[1])/(double)ggT * Vector[axis[1]][i];
1074 Tubevector[axis[2]][i] = Vector[axis[2]][i];
1075 }
1076 // here we assume, that Vector[axis[2]] is along z direction!
1077 gsl_matrix *M = gsl_matrix_alloc(2,2);
1078 gsl_matrix *C = gsl_matrix_alloc(2,2);
1079 gsl_matrix *evec = gsl_matrix_alloc(2,2);
1080 gsl_vector *eval = gsl_vector_alloc(2);
1081 gsl_vector *v = gsl_vector_alloc(2);
1082 gsl_vector *u = gsl_vector_alloc(2);
1083 gsl_eigen_symmv_workspace *w = gsl_eigen_symmv_alloc(2);
1084 gsl_matrix_set(C, 0,0, Vector[axis[0]][0]);
1085 gsl_matrix_set(C, 1,0, Vector[axis[0]][1]);
1086 gsl_matrix_set(C, 0,1, Vector[axis[1]][0]);
1087 gsl_matrix_set(C, 1,1, Vector[axis[1]][1]);
1088 gsl_blas_dgemm(CblasTrans,CblasNoTrans, 1.0, C, C, 0.0, M);
1089 fprintf(stdout, "M: \t%lg\t%lg\n\t%lg\t%lg\n", gsl_matrix_get(M,0,0), gsl_matrix_get(M,0,1), gsl_matrix_get(M,1,0), gsl_matrix_get(M,1,1));
1090 gsl_eigen_symmv(M, eval, evec, w);
1091 gsl_eigen_symmv_sort(eval,evec,GSL_EIGEN_SORT_ABS_DESC);
1092 fprintf(stdout, "Eigenvalues: %lg\t%lg\n", gsl_vector_get(eval,0), gsl_vector_get(eval,1));
1093 fprintf(stdout, "Eigenvectors: \t%lg\t%lg\n\t\t%lg\t%lg\n", gsl_matrix_get(evec,0,0), gsl_matrix_get(evec,0,1), gsl_matrix_get(evec,1,0), gsl_matrix_get(evec,1,1));
1094 gsl_matrix_set(M, 0,0, 0.);
1095 gsl_matrix_set(M, 1,0, 1.);
1096 gsl_matrix_set(M, 0,1, -gsl_vector_get(eval,1)/gsl_vector_get(eval,0));
1097 gsl_matrix_set(M, 1,1, 0.);
1098 gsl_vector_set(v,0,(double)chiral[0]);
1099 gsl_vector_set(v,1,(double)chiral[1]);
1100 gsl_blas_dgemm(CblasNoTrans,CblasNoTrans, 1.0, evec, M, 0.0, C);
1101 gsl_blas_dgemm(CblasNoTrans,CblasTrans, 1.0, C, evec, 0.0, M);
1102 fprintf(stdout, "M: \t%lg\t%lg\n\t%lg\t%lg\n", gsl_matrix_get(M,0,0), gsl_matrix_get(M,0,1), gsl_matrix_get(M,1,0), gsl_matrix_get(M,1,1));
1103 gsl_blas_dgemv(CblasNoTrans, 1.0, M, v, 0.0, u);
1104 fprintf(stdout, "Looking for factor to integer...\n");
1105 for(i=1;i<(chiral[0]+chiral[1])*(chiral[0]+chiral[1]);i++) {
1106 x1 = gsl_vector_get(u,0)*(double)i;
1107 x2 = gsl_vector_get(u,1)*(double)i;
1108 x3 =
1109 fprintf(stdout, "%d: %d\t%d vs. %lg\t%lg\n",i, ((int)(x1+x1/fabs(x1)*.5)), ((int)(x2+x2/fabs(x2)*.5)), (x1), (x2));
1110 if (( fabs( ((int)(x1+x1/fabs(x1)*.5)) - (x1) ) < 1e-6) && ( fabs( ((int)(x2+x2/fabs(x2)*.5)) - (x2) ) < 1e-6 )) {
1111 gsl_blas_dscal((double)i, u);
1112 break;
1113 }
1114 }
1115 fprintf(stdout, "(c,d) = (%lg,%lg)\n",gsl_vector_get(u,0), gsl_vector_get(u,1));
1116
1117 // get length
1118 double x[NDIM];
1119 for (i=0;i<NDIM;i++)
1120 x[i] = gsl_vector_get(u,0) * Vector[axis[0]][i] + gsl_vector_get(u,1) * Vector[axis[1]][i];
1121 angle = Norm(x)/Norm(Tubevector[axis[1]]) ;//ScalarProduct(x,Tubevector[axis[1]])/Norm(Tubevector[axis[1]]);
1122 fprintf(stdout, "angle is %lg\n", angle);
1123 for (i=0;i<NDIM;i++) {
1124 Tubevector[axis[1]][i] = gsl_vector_get(u,0) * Vector[axis[0]][i] + gsl_vector_get(u,1) * Vector[axis[1]][i];
1125 }
1126
1127 // Probe
1128 gsl_matrix_set(M, 0,0, Vector[axis[0]][0]);
1129 gsl_matrix_set(M, 1,0, Vector[axis[0]][1]);
1130 gsl_matrix_set(M, 0,1, Vector[axis[1]][0]);
1131 gsl_matrix_set(M, 1,1, Vector[axis[1]][1]);
1132 gsl_vector_set(v,0,(double)chiral[0]);
1133 gsl_vector_set(v,1,(double)chiral[1]);
1134 gsl_blas_dgemv(CblasNoTrans, 1.0, M, u, 0.0, eval);
1135 gsl_blas_dgemv(CblasNoTrans, 1.0, M, v, 0.0, u);
1136 x1=1.;
1137 gsl_blas_ddot(u,eval,&x1);
1138 fprintf(stdout, "Testing (c,d): (a,b) M^t M (c,d)^t = 0 ? : %lg\n", x1);
1139
1140 gsl_matrix_free(M);
1141 gsl_matrix_free(C);
1142 gsl_matrix_free(evec);
1143 gsl_vector_free(eval);
1144 gsl_vector_free(v);
1145 gsl_vector_free(u);
1146 gsl_eigen_symmv_free(w);
1147
1148 if (fabs(x1) > 1e-6) {
1149 fprintf(stderr,"Resulting TubeVectors of axis %d and %d and not orthogonal, aborting.\n", axis[0], axis[1]);
1150 return(128);
1151 }
1152
1153 // retrieve seed ...
1154 randomness = (double *) Calloc(sizeof(double)*numbercell, 0., "Main: at sheet - randomness");
1155 sprintf(dummy, "seed");
1156 fprintf(stdout, "%s ", dummy);
1157 while ((length = GetNextline(bufptr, line)) != 0) {
1158 bufptr += (length)*sizeof(char);
1159 if ((ptr = strstr(line, dummy)) != NULL)
1160 break;
1161 }
1162 if (length == 0) {
1163 fprintf(stderr, "ERROR: Main at stage Sheet - could not find %s in %s\n", dummy, SheetFilename);
1164 exit(255);
1165 }
1166 ptr += strlen(dummy);
1167 sscanf(ptr, "%d", &seed);
1168 fprintf(stdout, "%d\n", seed);
1169
1170 // ... and randomness
1171 if (seed != 0) { // only parse for values if a seed, i.e. randomness wanted, was specified
1172 sprintf(dummy, "Randomness");
1173 fprintf(stdout, "%s\n", dummy);
1174 while ((length = GetNextline(bufptr, line)) != 0) {
1175 bufptr += (length)*sizeof(char);
1176 if ((ptr = strstr(line, dummy)) != NULL)
1177 break;
1178 }
1179 if (length == 0) {
1180 fprintf(stderr, "ERROR: Main at stage Sheet - could not find %s in %s\n", dummy, SheetFilename);
1181 exit(255);
1182 }
1183 sprintf(dummy, "probability values");
1184 for (i=0;i<numbercell;i++) {
1185 length = GetNextline(bufptr, line);
1186 if (length == 0) {
1187 fprintf(stderr, "ERROR: Main at stage Sheet - could not find %s in %s\n", dummy, SheetFilename);
1188 exit(255);
1189 }
1190 bufptr += (length)*sizeof(char);
1191 sscanf(line, "%d %lg", &j, &dummydouble);
1192 randomness[j] = dummydouble;
1193 fprintf(stdout, "%d %g\n", j, randomness[j]);
1194 }
1195 }
1196 }
1197
1198 //int OrthOrder[NDIM] = { axis[2], axis[0], axis[1] };
1199 //Orthogonalize(Tubevector,OrthOrder);
1200 angle = acos(Projection(Vector[axis[0]], Vector[axis[1]])); // calcs angle between shanks in unit cell
1201 fprintf(stdout, "The basic angle between the two shanks of the unit cell is %lg %lg\n", angle/M_PI*180., Projection(Vector[axis[0]], Vector[axis[1]]));
1202 if ( angle/M_PI*180. > 90 ) {
1203 fprintf(stderr, "There seems to be something wrong with the unit cell! for nanotube the angle should be 60 degrees for example!\n");
1204 return 1;
1205 }
1206 angle = acos(Projection(Tubevector[axis[0]], Tubevector[axis[1]])); // calcs angle between shanks in unit cell
1207 fprintf(stdout, "The basic angle between the two shanks of the tube unit cell is %lg %lg\n", angle/M_PI*180., Projection(Tubevector[axis[0]], Tubevector[axis[1]]));
1208 //angle -= acos(Projection(Vector[axis[0]], Tubevector[axis[0]]));
1209 //angle = 30./180.*M_PI - acos(Projection(Vector[axis[0]], Tubevector[axis[0]]));
1210 //angle = acos(Projection(Tubevector[axis[0]], Vector[axis[0]]));
1211 fprintf(stdout, "The relative alignment rotation angle then is %lg\n", angle/M_PI*180.);
1212 if (fabs(Tubevector[axis[0]][0]) > MYEPSILON)
1213 angle = -M_PI/2. + acos(Tubevector[axis[0]][0]/Norm(Tubevector[axis[0]]));
1214 else
1215 angle = 0.;
1216 fprintf(stdout, "The absolute alignment rotation angle then is %lg %lg\n", angle/M_PI*180., Tubevector[axis[0]][0]/Norm(Tubevector[axis[0]]));
1217 fprintf(stdout, "\nThe chiral angle then is %5.5f degrees with tube radius %5.5f A and length %5.5f A, i.e. final torus radius of %5.5f A.\n",
1218 acos(Projection(Vector[axis[0]], Tubevector[axis[0]]))/M_PI*180.,
1219 (double)factors[0]*Norm(Tubevector[axis[0]])/(2.*M_PI),
1220 (double)factors[1]*Norm(Tubevector[axis[1]]),
1221 (double)factors[1]*Norm(Tubevector[axis[1]])/(2.*M_PI)
1222 );
1223 Orthogonalize(Tubevector, axis); // with correct translational vector, not needed anymore (? what's been done here. Hence, re-inserted)
1224 fprintf(stdout, "Tubevector magnitudes: %5.5lg %5.5lg %5.5lg\n", Norm(Tubevector[0]), Norm(Tubevector[1]), Norm(Tubevector[2]));
1225 fprintf(stdout, "Tubevectors are \n");
1226 PrintMatrix(stdout, Tubevector);
1227 MatrixInversion(Tubevector, TubevectorInverse);
1228 //Transpose(TubevectorInverse);
1229 fprintf(stdout, "Vector\n");
1230 PrintMatrix(stdout, Vector);
1231 fprintf(stdout, "TubevectorInverse\n");
1232 PrintMatrix(stdout, TubevectorInverse);
1233 for (i=0;i<NDIM;i++) {
1234 fprintf(stdout, "Vector %d in TubeVectorInverse vectors:\t", axis[i]);
1235 tempvector = MatrixTrafoInverse(Vector[axis[i]], TubevectorInverse);
1236 PrintVector(stdout, tempvector);
1237 Free(tempvector, "Main:tempvector");
1238 }
1239 fprintf(stdout, "Reciprocal Tubebvectors are \n");
1240 PrintMatrix(stdout, TubevectorInverse);
1241 fprintf(stdout, "Tubevector magnitudes: %5.5lg %5.5lg %5.5lg\n", Norm(Tubevector[0]), Norm(Tubevector[1]), Norm(Tubevector[2]));
1242
1243 biggestdiameter = DetermineBiggestDiameter(Tubevector, axis, factors);
1244 for (i=0;i<NDIM;i++) {
1245 sheetnr[i] = 0;
1246 }
1247 for (i=0;i<NDIM;i++) {
1248 for (j=0;j<NDIM;j++) {
1249// sheetnr[j] = ceil(biggestdiameter/Norm(Vector[j]));
1250 if (fabs(Vector[i][j]) > MYEPSILON) {
1251 tmp = ceil(biggestdiameter/fabs(Vector[i][j]));
1252 } else {
1253 tmp = 0;
1254 }
1255 sheetnr[j] = sheetnr[j] > tmp ? sheetnr[j] : tmp;
1256 }
1257 }
1258 fprintf(stdout, "Maximum indices to regard: %d %d %d\n", sheetnr[0], sheetnr[1], sheetnr[2]);
1259 for (i=0;i<NDIM;i++) {
1260 fprintf(stdout, "For axis %d: (%5.5lg\t%5.5lg\t%5.5lg) with %5.5lg\n", i, (Vector[i][0]*sheetnr[i]), (Vector[i][1]*sheetnr[i]), (Vector[i][2]*sheetnr[i]), Norm(Vector[i]));
1261 }
1262
1263 //if (!strncmp(stage, "Cell", 4)) {
1264 // parse in atoms for quicker processing
1265 struct Atoms *atombuffer = malloc(sizeof(struct Atoms)*numbercell);
1266 bufptr = CellBuffer;
1267 bufptr += GetNextline(bufptr, line)*sizeof(char);
1268 bufptr += GetNextline(bufptr, line)*sizeof(char);
1269 for (i=0;i<numbercell;i++) {
1270 if ((length = GetNextline(bufptr, line)) != 0) {
1271 bufptr += length*sizeof(char);
1272 sscanf(line, "%s %lg %lg %lg", atombuffer[i].name, &(atombuffer[i].x[0]), &(atombuffer[i].x[1]), &(atombuffer[i].x[2]));
1273 fprintf(stdout, "Read Atombuffer Nr %i: %s %5.5lg %5.5lg %5.5lg\n", i+1, atombuffer[i].name, atombuffer[i].x[0], atombuffer[i].x[1], atombuffer[i].x[2]);
1274 } else {
1275 fprintf(stdout, "Error reading Atom Nr. %i\n", i+1);
1276 break;
1277 }
1278 }
1279 SheetFile = fopen(SheetFilename, "w");
1280 if (SheetFile == NULL) {
1281 fprintf(stderr, "ERROR: main - can't open %s for writing\n", SheetFilename);
1282 exit(255);
1283 }
1284 SheetFileAligned = fopen(SheetFilenameAligned, "w");
1285 if (SheetFile == NULL) {
1286 fprintf(stderr, "ERROR: main - can't open %s for writing\n", SheetFilenameAligned);
1287 exit(255);
1288 }
1289 // Now create the sheet
1290 double index[NDIM];
1291 int nr;//, nummer = 0;
1292 numbersheet = 0;
1293 index[axis[2]] = 0;
1294 // initialise pseudo random number generator with given seed
1295 fprintf(stdout, "Initialising pseudo random number generator with given seed %d.\n", seed);
1296 srand(seed);
1297 //for (index[axis[0]] = 0; index[axis[0]] <= sheetnr[axis[0]]; index[axis[0]]++) { // NOTE: minor axis may start from 0! Check on this later ...
1298 for (index[axis[0]] = -sheetnr[axis[0]]+1; index[axis[0]] < sheetnr[axis[0]]; index[axis[0]]++) { // NOTE: minor axis may start from 0! Check on this later ...
1299 //for (index[axis[1]] = 0; index[axis[1]] <= sheetnr[axis[1]]; index[axis[1]]++) { // These are all the cells that need be checked on
1300 for (index[axis[1]] = -sheetnr[axis[1]]+1; index[axis[1]] < sheetnr[axis[1]]; index[axis[1]]++) { // These are all the cells that need be checked on
1301 // Calculate offset in cartesian coordinates
1302 offset = MatrixTrafo(Vector, index);
1303
1304 //fprintf(stdout, "Now dealing with numbercell atoms in unit cell at R = (%lg,%lg,%lg)\n", offset[0], offset[1], offset[2]);
1305 for (nr = 0; nr < numbercell; nr++) {
1306 percentage = rand()/(RAND_MAX+1.0);
1307 //fprintf(stdout, "Lucky number for %d is %lg >? %lg\n", nr, percentage, randomness[nr]);
1308 if (percentage >= randomness[nr]) {
1309 // Create coordinates at atom site
1310 coord = VectorAdd(atombuffer[nr].x, offset);
1311 //fprintf(stdout, "Atom Nr. %i: ", (numbersheet+1));
1312 //PrintVector(stdout, coord);
1313 // project down on major and minor Tubevectors and check for length if out of sheet
1314 tempvector = MatrixTrafoInverse(coord, TubevectorInverse);
1315 if (((tempvector[axis[0]] + MYEPSILON) > 0) && ((factors[0] - tempvector[axis[0]]) > MYEPSILON) &&
1316 ((tempvector[axis[1]] + MYEPSILON) > 0) && ((factors[1] - tempvector[axis[1]]) > MYEPSILON) &&
1317 ((tempvector[axis[2]] + MYEPSILON) > 0) && ((factors[2] - tempvector[axis[2]]) > MYEPSILON)) { // check if within rotated cell numbersheet++;
1318 //if (nummer >= 2) strcpy(atombuffer[nr].name, "O");
1319 //nummer++;
1320 fprintf(SheetFile, "%s\t%5.5lg\t%5.5lg\t%5.5lg\n", atombuffer[nr].name, coord[0], coord[1], coord[2]);
1321 // rotate to align the sheet in xy plane
1322 x1 = coord[0]*cos(-angle) + coord[1] * sin(-angle);
1323 x2 = coord[0]*(-sin(-angle)) + coord[1] * cos(-angle);
1324 x3 = coord[2];
1325 fprintf(SheetFileAligned, "%s\t%5.5lg\t%5.5lg\t%5.5lg\n", atombuffer[nr].name, x1, x2, x3);
1326 //fprintf(SheetFile, "O\t%5.5lg\t%5.5lg\t%5.5lg\n", coord[0], coord[1], coord[2]);
1327 //fprintf(stdout, "%s/%d\t(%lg\t%lg\t%lg)\t", atombuffer[nr].name, numbersheet+1, coord[0], coord[1], coord[2]);
1328 //PrintVector(stdout, tempvector);
1329 numbersheet++;
1330 //fprintf(stdout, "%i,", nr);
1331 } //else {
1332 //numbersheet++;
1333 //fprintf(SheetFile, "B\t%lg\t%lg\t%lg\n", coord[0], coord[1], coord[2]);
1334 //fprintf(stdout, "O \t(%lg\t%lg\t%lg)\n", coord[0], coord[1], coord[2]);
1335 //fprintf(stdout, "!!%i!!, ", nr);
1336 //}
1337 Free(tempvector, "Main: At stage Sheet - tempvector");
1338 Free(coord, "Main: At stage Sheet - coord");
1339 }
1340 }
1341 Free(offset, "Main: At stage Sheet - offset");
1342 }
1343 //fprintf(stdout, "\n";
1344 }
1345
1346 fclose(SheetFile);
1347 fclose(SheetFileAligned);
1348 AddAtomicNumber(SheetFilename,numbersheet, Vector, Recivector); // prepend atomic number and comment
1349 AddAtomicNumber(SheetFilenameAligned,numbersheet, Vector, Recivector); // prepend atomic number and comment
1350 AddSheetInfo(SheetFilename,axis,chiral, factors, seed, numbercell, randomness);
1351 fprintf(stdout, "\nThere are %i atoms in the created sheet.\n", numbersheet);
1352
1353 strncpy(stage, "Sheet", 5);
1354 //}
1355 SheetBuffer = ReadBuffer(SheetFilename, &length);
1356
1357
1358 // ======================== STAGE: Tube ==============================
1359 // The tube starts with the rectangular (due to the orthogonalization) sheet
1360 // just created (or read). Along the minor axis it is rolled up, i.e. projected
1361 // from a 2d surface onto a cylindrical surface (x,y,z <-> r,alpha,z). The only
1362 // thing that's a bit complex is that the sheet it not aligned along the cartesian
1363 // axis but along major and minor. That's why we have to transform the atomic
1364 // cartesian coordinates into the orthogonal tubevector base, do the rolling up
1365 // there (and regard that minor and major axis must not necessarily be of equal
1366 // length) and afterwards transform back again (where we need the $halfaxis due to
1367 // the above possible inequality).
1368
1369 FILE *TubeFile = NULL;
1370 FILE *TubeFileAligned = NULL;
1371
1372 Debug ("STAGE: Tube\n");
1373 if (!strncmp(stage, "Sheet", 4)) {
1374 TubeFile = fopen(TubeFilename, "w");
1375 if (TubeFile == NULL) {
1376 fprintf(stderr, "ERROR: Main - can't open %s for writing\n", TubeFilename);
1377 exit(255);
1378 }
1379 TubeFileAligned = fopen(TubeFilenameAligned, "w");
1380 if (TubeFile == NULL) {
1381 fprintf(stderr, "ERROR: Main - can't open %s for writing\n", TubeFilenameAligned);
1382 exit(255);
1383 }
1384 bufptr = SheetBuffer;
1385 bufptr += GetNextline(bufptr, line); // write numbers to file
1386 bufptr += GetNextline(bufptr, line); // write comment to file
1387
1388 //cog = CenterOfGravity(bufptr, numbersheet);
1389 //cog_projected = MatrixTrafoInverse(cog, TubevectorInverse);
1390 //fprintf(stdout, "\nCenter of Gravity at (%5.5lg\t%5.5lg\t%5.5lg) and projected at (%5.5lg\t%5.5lg\t%5.5lg)\n", cog[0], cog[1], cog[2], cog_projected[0], cog_projected[1], cog_projected[2]);
1391
1392 // restart
1393 bufptr = SheetBuffer;
1394 bufptr += GetNextline(bufptr, line); // write numbers to file
1395 bufptr += GetNextline(bufptr, line); // write numbers to file
1396
1397 // determine half axis as tube vector not necessarily have the same length
1398 double halfaxis[NDIM];
1399 for (i=0;i<NDIM;i++)
1400 halfaxis[i] = factors[0]*Norm(Tubevector[axis[0]])/Norm(Tubevector[i]);
1401
1402 double arg, radius;
1403 for (i=0;i<numbersheet;i++) {
1404 // scan next atom
1405 bufptr += GetNextline(bufptr, line);
1406 sscanf(line, "%s %lg %lg %lg", name, &atom[0], &atom[1], &atom[2]);
1407
1408 // transform atom coordinates in cartesian system to the axis eigensystem
1409 x = MatrixTrafoInverse(atom, TubevectorInverse);
1410 //x = VectorAdd(y, cog_projected);
1411 //free(y);
1412
1413 // roll up (project (x,y,z) on cylindrical coordinates (radius,arg,z))
1414 arg = 2.*M_PI*x[axis[0]]/(factors[0]) - M_PI; // is angle
1415 radius = 1./(2.*M_PI); // is length of sheet in units of axis vector, divide by pi to get radius (from circumference)
1416 // fprintf(stdout, "arg: %5.2f (c%2.2f,s%2.2f)\t",$arg, cos($arg), sin($arg));
1417 x[axis[0]] = cos(arg)*halfaxis[axis[0]]*(radius+x[axis[2]]/halfaxis[axis[2]]); // as both vectors are not normalized additional betrag has to be taken into account!
1418 x[axis[2]] = sin(arg)*halfaxis[axis[2]]*(radius+x[axis[2]]/halfaxis[axis[2]]); // due to the back-transformation from eigensystem to cartesian one
1419 //fprintf(stdout, "rotated: (%5.2f,%5.2f,%5.2f)\n",x[0],x[1],x[2]);
1420 atom_transformed = MatrixTrafo(Tubevector, x);
1421 fprintf(TubeFile, "%s\t%lg\t%lg\t%lg\n", name, atom_transformed[0], atom_transformed[1], atom_transformed[2]);
1422 // rotate and flip to align tube in z-direction
1423 x1 = atom_transformed[0]*cos(-angle) + atom_transformed[1] * sin(-angle);
1424 x2 = atom_transformed[0]*(-sin(-angle)) + atom_transformed[1] * cos(-angle);
1425 x3 = atom_transformed[2];
1426 fprintf(TubeFileAligned, "%s\t%lg\t%lg\t%lg\n", name, x3, x2, x1); // order so that symmetry is along z axis
1427 //fprintf(stdout, "%s\t%5.5lg\t%5.5lg\t%5.5lg\n", name, atom_transformed[0], atom_transformed[1] ,atom_transformed[2]);
1428
1429 Free(x, "Main: at stage Tube - x");
1430 Free(atom_transformed, "Main: at stage Tube - atom_transformed");
1431 }
1432
1433
1434 fclose(TubeFile);
1435 fclose(TubeFileAligned);
1436 //free(cog);
1437 //free(cog_projected);
1438 AddAtomicNumber(TubeFilename,numbersheet, Vector, Recivector); // prepend atomic number and comment
1439 AddAtomicNumber(TubeFilenameAligned,numbersheet, Vector, Recivector); // prepend atomic number and comment
1440 AddSheetInfo(TubeFilename,axis,chiral, factors, seed, numbercell, randomness);
1441 fprintf(stdout, "\nThere are %i atoms in the created tube.\n", numbersheet);
1442
1443 strncpy(stage, "Tube", 4);
1444 } else {
1445 }
1446
1447 TubeBuffer = ReadBuffer(TubeFilename, &length);
1448
1449 // ======================== STAGE: Torus =============================
1450 // The procedure for the torus is very much alike to the one used to make the
1451 // tube. Only the projection is not from 2d surface onto a cylindrical one but
1452 // from a cylindrial onto a torus surface
1453 // (x,y,z) <-> (cos(s)*(R+r*cos(t)), sin(s)*(R+rcos(t)), r*sin(t)).
1454 // Here t is the angle within the tube with radius r, s is the torus angle with
1455 // radius R. We get R from the tubelength (that's why we need lengthfactor to
1456 // make it long enough). And due to fact that we have it already upon a cylindrical
1457 // surface, r*cos(t) and r*sin(t) already reside in $minoraxis and $noaxis.
1458
1459 FILE *TorusFile;
1460
1461 Debug ("STAGE: Torus\n");
1462 if (!strncmp(stage, "Tube", 4)) {
1463 TorusFile = fopen(TorusFilename, "w");
1464 if (TorusFile == NULL) {
1465 fprintf(stderr, "ERROR: main - can't open %s for writing\n", TorusFilename);
1466 exit(255);
1467 }
1468 bufptr = TubeBuffer;
1469 bufptr += GetNextline(bufptr, line); // write numbers to file
1470 bufptr += GetNextline(bufptr, line); // write comment to file
1471
1472 //cog = CenterOfGravity(bufptr, numbersheet);
1473 //cog_projected = MatrixTrafoInverse(cog, TubevectorInverse);
1474 //fprintf(stdout, "\nCenter of Gravity at (%5.5lg\t%5.5lg\t%5.5lg) and projected at (%5.5lg\t%5.5lg\t%5.5lg)\n", cog[0], cog[1], cog[2], cog_projected[0], cog_projected[1], cog_projected[2]);
1475
1476 // determine half axis as tube vectors not necessarily have same length
1477 double halfaxis[NDIM];
1478 for (i=0;i<NDIM;i++)
1479 halfaxis[i] = Norm(Tubevector[axis[1]])/Norm(Tubevector[i]);
1480
1481 double arg, radius;
1482 for (i=0;i<numbersheet;i++) {
1483 // scan next atom
1484 bufptr += GetNextline(bufptr, line);
1485 sscanf(line, "%s %lg %lg %lg", name, &atom[0], &atom[1], &atom[2]);
1486
1487 // transform atom coordinates in cartesian system to the axis eigensystem
1488 x = MatrixTrafoInverse(atom, TubevectorInverse);
1489 //x = VectorAdd(y, cog_projected);
1490 //free(y);
1491
1492 // roll up (project (x,y,z) on cylindrical coordinates (radius,arg,z))
1493 arg = 2.*M_PI*x[axis[1]]/(factors[1]) - M_PI; // is angle
1494 radius = (factors[1])/(2.*M_PI) + x[axis[0]]/halfaxis[axis[0]]; // is length of sheet in units of axis vector, divide by pi to get radius (from circumference)
1495 // fprintf(stdout, "arg: %5.2f (c%2.2f,s%2.2f)\t",$arg, cos($arg), sin($arg));
1496 x[axis[0]] = cos(arg)*halfaxis[axis[0]]*radius; // as both vectors are not normalized additional betrag has to be taken into account!
1497 x[axis[1]] = sin(arg)*halfaxis[axis[1]]*radius; // due to the back-transformation from eigensystem to cartesian one
1498 //fprintf(stdout, "rotated: (%5.2f,%5.2f,%5.2f)\n",x[0],x[1],x[2]);
1499 atom_transformed = MatrixTrafo(Tubevector, x);
1500 fprintf(TorusFile, "%s\t%lg\t%lg\t%lg\n", name, atom_transformed[0], atom_transformed[1] ,atom_transformed[2]);
1501 //fprintf(stdout, "%s\t%5.5lg\t%5.5lg\t%5.5lg\n", name, atom_transformed[0], atom_transformed[1] ,atom_transformed[2]);
1502
1503 Free(x, "Main: at stage Torus - x");
1504 Free(atom_transformed, "Main: at stage Torus - atom_transformed");
1505 }
1506
1507 fclose(TorusFile);
1508 //free(cog);
1509 //free(cog_projected);
1510 AddAtomicNumber(TorusFilename,numbersheet, Vector, Recivector); // prepend atomic number and comment
1511 AddSheetInfo(TorusFilename,axis,chiral, factors, seed, numbercell, randomness);
1512 fprintf(stdout, "\nThere are %i atoms in the created torus.\n", numbersheet);
1513
1514 strncpy(stage, "Torus", 5);
1515 } else {
1516 }
1517
1518 // Free memory
1519 for (i=0; i<NDIM; i++ ) {
1520 Free(Vector[i], "Main: end of stages - *Vector");
1521 Free(Recivector[i], "Main: end of stages - *Recivector");
1522 Free(Tubevector[i], "Main: end of stages - *Tubevector");
1523 Free(TubevectorInverse[i], "Main: end of stages - *TubevectorInverse");
1524 }
1525 Free(atom, "Main: end of stages - atom");
1526 Free(Vector, "Main: end of stages - Vector");
1527 Free(Recivector, "Main: end of stages - Recivector");
1528 Free(Tubevector, "Main: end of stages - Tubevector");
1529 Free(TubevectorInverse, "Main: end of stages - TubevectorInverse");
1530 Free(randomness, "Main: at stage Sheet - randomness");
1531
1532 if (CellBuffer != NULL) Free(CellBuffer, "Main: end of stages - CellBuffer");
1533 if (SheetBuffer != NULL) Free(SheetBuffer, "Main: end of stages - SheetBuffer");
1534 if (TubeBuffer != NULL) Free(TubeBuffer, "Main: end of stages - TubeBuffer");
1535
1536 Free(CellFilename, "Main: end of stafes - CellFilename");
1537 Free(SheetFilename, "Main: end of stafes - CellFilename");
1538 Free(TubeFilename, "Main: end of stafes - CellFilename");
1539 Free(TorusFilename, "Main: end of stafes - CellFilename");
1540 Free(SheetFilenameAligned, "Main: end of stafes - CellFilename");
1541 Free(TubeFilenameAligned, "Main: end of stafes - CellFilename");
1542
1543 // exit
1544 exit(0);
1545}
Note: See TracBrowser for help on using the repository browser.