source: util/NanoCreator.c@ bd7b85

Last change on this file since bd7b85 was a0bcf1, checked in by Frederik Heber <heber@…>, 17 years ago

-initial commit
-Minimum set of files needed from ESPACK SVN repository
-Switch to three tantamount package parts instead of all relating to pcp (as at some time Ralf's might find inclusion as well)

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