source: util/src/NanoCreator.c@ b53a7e

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

NanoCreator - length and radius parameters are now parsed as doubles and rounded to integers.

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