1 | /** \file mymath.c
|
---|
2 | * Linear algebra mathematical routines.
|
---|
3 | * Small library of often needed mathematical routines such as hard-coded
|
---|
4 | * vector VP3(), scalar SP(), matrix products RMat33Vec3(), RMatMat33(), RVec3Mat33(),
|
---|
5 | * multiplication with scalar SM(), euclidian distance Dist(),inverse RMatReci3(),
|
---|
6 | * transposed RTranspose3(), modulo Rest(), nullifying NV(), SetArrayToDouble0(),
|
---|
7 | * gamma function gammln(), gaussian error function derf(), integration via
|
---|
8 | * Simpsons Rule Simps().\n
|
---|
9 | * Also for printing matrixes PrintCMat330(), PrintRMat330() and vectors
|
---|
10 | * PrintCVec30(), PrintRVec30() to screen.\n
|
---|
11 | * All specialized for 3x3 real or complex ones.\n
|
---|
12 | * Rather specialized is RotateToAlign() which is needed in transforming the whole coordinate
|
---|
13 | * system in order to align a certain vector.
|
---|
14 | *
|
---|
15 | Project: ParallelCarParrinello
|
---|
16 | \author Jan Hamaekers
|
---|
17 | \date 2000
|
---|
18 |
|
---|
19 | File: helpers.c
|
---|
20 | $Id: mymath.c,v 1.25 2007-03-29 13:38:30 foo Exp $
|
---|
21 | */
|
---|
22 |
|
---|
23 | #include<stdlib.h>
|
---|
24 | #include<stdio.h>
|
---|
25 | #include<stddef.h>
|
---|
26 | #include<math.h>
|
---|
27 | #include<string.h>
|
---|
28 | #include"mymath.h"
|
---|
29 |
|
---|
30 | // use double precision fft when we have it
|
---|
31 | #ifdef HAVE_CONFIG_H
|
---|
32 | #include <config.h>
|
---|
33 | #endif
|
---|
34 |
|
---|
35 | #ifdef HAVE_DFFTW_H
|
---|
36 | #include "dfftw.h"
|
---|
37 | #else
|
---|
38 | #include "fftw.h"
|
---|
39 | #endif
|
---|
40 |
|
---|
41 | #ifdef HAVE_GSL_GSL_SF_ERF_H
|
---|
42 | #include "gsl/gsl_sf_erf.h"
|
---|
43 | #endif
|
---|
44 |
|
---|
45 |
|
---|
46 | /** efficiently compute x^n
|
---|
47 | * \param x argument
|
---|
48 | * \param n potency
|
---|
49 | * \return \f$x^n\f$
|
---|
50 | */
|
---|
51 | #ifdef HAVE_INLINE
|
---|
52 | inline double tpow(double x, int n)
|
---|
53 | #else
|
---|
54 | double tpow(double x, int n)
|
---|
55 | #endif
|
---|
56 | {
|
---|
57 | double y = 1;
|
---|
58 | int neg = (n < 0);
|
---|
59 |
|
---|
60 | if (neg) n = -n;
|
---|
61 |
|
---|
62 | while (n) {
|
---|
63 | if (n & 1) y *= x;
|
---|
64 | x *= x;
|
---|
65 | n >>= 1;
|
---|
66 | }
|
---|
67 | return neg ? 1.0/y : y;
|
---|
68 | }
|
---|
69 |
|
---|
70 |
|
---|
71 | /** Modulo function.
|
---|
72 | * Normal modulo operation, yet return value is >=0
|
---|
73 | * \param n denominator
|
---|
74 | * \param m divisor
|
---|
75 | * \return modulo >=0
|
---|
76 | */
|
---|
77 | #ifdef HAVE_INLINE
|
---|
78 | inline int Rest(int n, int m) /* normale modulo-Funktion, Ausgabe>=0 */
|
---|
79 | #else
|
---|
80 | int Rest(int n, int m) /* normale modulo-Funktion, Ausgabe>=0 */
|
---|
81 | #endif
|
---|
82 | {
|
---|
83 | int q = n%m;
|
---|
84 | if (q >= 0) return (q);
|
---|
85 | return ((q) + m);
|
---|
86 | }
|
---|
87 |
|
---|
88 | /* Rechnungen */
|
---|
89 |
|
---|
90 | /** Real 3x3 inverse of matrix.
|
---|
91 | * Calculates the inverse of a matrix by b_ij = A_ij/det(A), where
|
---|
92 | * is A_ij is the matrix with row j and column i removed.
|
---|
93 | * \param B inverse matrix array (set by function)
|
---|
94 | * \param A matrix array to be inverted
|
---|
95 | * \return 0 - error: det A == 0, 1 - success
|
---|
96 | */
|
---|
97 | #ifdef HAVE_INLINE
|
---|
98 | inline int RMatReci3(double B[NDIM_NDIM], const double A[NDIM_NDIM])
|
---|
99 | #else
|
---|
100 | int RMatReci3(double B[NDIM_NDIM], const double A[NDIM_NDIM])
|
---|
101 | #endif
|
---|
102 | {
|
---|
103 | double detA = RDET3(A);
|
---|
104 | double detAReci;
|
---|
105 | if (detA == 0.0) return 1; // RDET3(A) yields precisely zero if A irregular
|
---|
106 | detAReci = 1./detA;
|
---|
107 | B[0] = detAReci*RDET2(A[4],A[5],A[7],A[8]); // A_11
|
---|
108 | B[1] = -detAReci*RDET2(A[1],A[2],A[7],A[8]); // A_12
|
---|
109 | B[2] = detAReci*RDET2(A[1],A[2],A[4],A[5]); // A_13
|
---|
110 | B[3] = -detAReci*RDET2(A[3],A[5],A[6],A[8]); // A_21
|
---|
111 | B[4] = detAReci*RDET2(A[0],A[2],A[6],A[8]); // A_22
|
---|
112 | B[5] = -detAReci*RDET2(A[0],A[2],A[3],A[5]); // A_23
|
---|
113 | B[6] = detAReci*RDET2(A[3],A[4],A[6],A[7]); // A_31
|
---|
114 | B[7] = -detAReci*RDET2(A[0],A[1],A[6],A[7]); // A_32
|
---|
115 | B[8] = detAReci*RDET2(A[0],A[1],A[3],A[4]); // A_33
|
---|
116 | return 0;
|
---|
117 | }
|
---|
118 |
|
---|
119 | /** Real 3x3 Matrix multiplication.
|
---|
120 | * Hard-coded falk scheme for multiplication of matrix1 * matrix2
|
---|
121 | * \param C product matrix
|
---|
122 | * \param A matrix1 array
|
---|
123 | * \param B matrix2 array
|
---|
124 | */
|
---|
125 | #ifdef HAVE_INLINE
|
---|
126 | inline void RMatMat33(double C[NDIM*NDIM], const double A[NDIM*NDIM], const double B[NDIM*NDIM])
|
---|
127 | #else
|
---|
128 | void RMatMat33(double C[NDIM*NDIM], const double A[NDIM*NDIM], const double B[NDIM*NDIM])
|
---|
129 | #endif
|
---|
130 | {
|
---|
131 | C[0] = A[0]*B[0]+A[3]*B[1]+A[6]*B[2];
|
---|
132 | C[1] = A[1]*B[0]+A[4]*B[1]+A[7]*B[2];
|
---|
133 | C[2] = A[2]*B[0]+A[5]*B[1]+A[8]*B[2];
|
---|
134 | C[3] = A[0]*B[3]+A[3]*B[4]+A[6]*B[5];
|
---|
135 | C[4] = A[1]*B[3]+A[4]*B[4]+A[7]*B[5];
|
---|
136 | C[5] = A[2]*B[3]+A[5]*B[4]+A[8]*B[5];
|
---|
137 | C[6] = A[0]*B[6]+A[3]*B[7]+A[6]*B[8];
|
---|
138 | C[7] = A[1]*B[6]+A[4]*B[7]+A[7]*B[8];
|
---|
139 | C[8] = A[2]*B[6]+A[5]*B[7]+A[8]*B[8];
|
---|
140 | }
|
---|
141 |
|
---|
142 | /** Real 3x3 Matrix vector multiplication.
|
---|
143 | * hard-coded falk scheme for multiplication of matrix * vector
|
---|
144 | * \param C resulting vector
|
---|
145 | * \param M matrix array
|
---|
146 | * \param V vector array
|
---|
147 | */
|
---|
148 | #ifdef HAVE_INLINE
|
---|
149 | inline void RMat33Vec3(double C[NDIM], const double M[NDIM*NDIM], const double V[NDIM])
|
---|
150 | #else
|
---|
151 | void RMat33Vec3(double C[NDIM], const double M[NDIM*NDIM], const double V[NDIM])
|
---|
152 | #endif
|
---|
153 | {
|
---|
154 | C[0] = M[0]*V[0]+M[3]*V[1]+M[6]*V[2];
|
---|
155 | C[1] = M[1]*V[0]+M[4]*V[1]+M[7]*V[2];
|
---|
156 | C[2] = M[2]*V[0]+M[5]*V[1]+M[8]*V[2];
|
---|
157 | }
|
---|
158 |
|
---|
159 | /** Real 3x3 vector Matrix multiplication.
|
---|
160 | * hard-coded falk scheme for multiplication of vector * matrix
|
---|
161 | * \param C resulting vector
|
---|
162 | * \param V vector array
|
---|
163 | * \param M matrix array
|
---|
164 | */
|
---|
165 | #ifdef HAVE_INLINE
|
---|
166 | inline void RVec3Mat33(double C[NDIM], const double V[NDIM], const double M[NDIM*NDIM])
|
---|
167 | #else
|
---|
168 | void RVec3Mat33(double C[NDIM], const double V[NDIM], const double M[NDIM*NDIM])
|
---|
169 | #endif
|
---|
170 | {
|
---|
171 | C[0] = V[0]*M[0]+V[1]*M[1]+V[2]*M[2];
|
---|
172 | C[1] = V[0]*M[3]+V[1]*M[4]+V[2]*M[5];
|
---|
173 | C[2] = V[0]*M[6]+V[1]*M[7]+V[2]*M[8];
|
---|
174 | }
|
---|
175 |
|
---|
176 | /** Real 3x3 vector product.
|
---|
177 | * vector product of vector1 x vector 2
|
---|
178 | * \param V resulting orthogonal vector
|
---|
179 | * \param A vector1 array
|
---|
180 | * \param B vector2 array
|
---|
181 | */
|
---|
182 | #ifdef HAVE_INLINE
|
---|
183 | inline void VP3(double V[NDIM], double A[NDIM], double B[NDIM])
|
---|
184 | #else
|
---|
185 | void VP3(double V[NDIM], double A[NDIM], double B[NDIM])
|
---|
186 | #endif
|
---|
187 | {
|
---|
188 | V[0] = A[1]*B[2]-A[2]*B[1];
|
---|
189 | V[1] = A[2]*B[0]-A[0]*B[2];
|
---|
190 | V[2] = A[0]*B[1]-A[1]*B[0];
|
---|
191 | }
|
---|
192 |
|
---|
193 | /** Real transposition of 3x3 Matrix.
|
---|
194 | * \param *A Matrix
|
---|
195 | */
|
---|
196 | #ifdef HAVE_INLINE
|
---|
197 | inline void RTranspose3(double *A)
|
---|
198 | #else
|
---|
199 | void RTranspose3(double *A)
|
---|
200 | #endif
|
---|
201 | {
|
---|
202 | double dummy = A[1];
|
---|
203 | A[1] = A[3];
|
---|
204 | A[3] = dummy;
|
---|
205 | dummy = A[2];
|
---|
206 | A[2] = A[6];
|
---|
207 | A[6] = dummy;
|
---|
208 | dummy = A[5];
|
---|
209 | A[5] = A[7];
|
---|
210 | A[7] = dummy;
|
---|
211 | }
|
---|
212 |
|
---|
213 | /** Scalar product.
|
---|
214 | * \param *a first vector
|
---|
215 | * \param *b second vector
|
---|
216 | * \param n dimension
|
---|
217 | * \return scalar product of a with b
|
---|
218 | */
|
---|
219 | #ifdef HAVE_INLINE
|
---|
220 | inline double SP(const double *a, const double *b, const int n)
|
---|
221 | #else
|
---|
222 | double SP(const double *a, const double *b, const int n)
|
---|
223 | #endif
|
---|
224 | {
|
---|
225 | int i;
|
---|
226 | double dummySP;
|
---|
227 | dummySP = 0;
|
---|
228 | for (i = 0; i < n; i++) {
|
---|
229 | dummySP += ((a[i]) * (b[i]));
|
---|
230 | }
|
---|
231 | return dummySP;
|
---|
232 | }
|
---|
233 |
|
---|
234 | /** Euclidian distance.
|
---|
235 | * \param *a first vector
|
---|
236 | * \param *b second vector
|
---|
237 | * \param n dimension
|
---|
238 | * \return sqrt(a-b)
|
---|
239 | */
|
---|
240 | #ifdef HAVE_INLINE
|
---|
241 | inline double Dist(const double *a, const double *b, const int n)
|
---|
242 | #else
|
---|
243 | double Dist(const double *a, const double *b, const int n)
|
---|
244 | #endif
|
---|
245 | {
|
---|
246 | int i;
|
---|
247 | double dummyDist = 0;
|
---|
248 | for (i = 0; i < n; i++) {
|
---|
249 | dummyDist += (a[i]-b[i])*(a[i]-b[i]);
|
---|
250 | }
|
---|
251 | return (sqrt(dummyDist));
|
---|
252 | }
|
---|
253 |
|
---|
254 |
|
---|
255 | /** Multiplication with real scalar.
|
---|
256 | * \param *a vector (changed)
|
---|
257 | * \param c scalar
|
---|
258 | * \param n dimension
|
---|
259 | */
|
---|
260 | #ifdef HAVE_INLINE
|
---|
261 | inline void SM(double *a, const double c, const int n)
|
---|
262 | #else
|
---|
263 | void SM(double *a, const double c, const int n)
|
---|
264 | #endif
|
---|
265 | {
|
---|
266 | int i;
|
---|
267 | for (i = 0; i < n; i++) a[i] *= c;
|
---|
268 | }
|
---|
269 |
|
---|
270 | /** nullify vector.
|
---|
271 | * sets all components of vector /a a to zero.
|
---|
272 | * \param *a vector (changed)
|
---|
273 | * \param n dimension
|
---|
274 | */
|
---|
275 | #ifdef HAVE_INLINE
|
---|
276 | inline void NV(double *a, const int n)
|
---|
277 | #else
|
---|
278 | void NV(double *a, const int n)
|
---|
279 | #endif
|
---|
280 | {
|
---|
281 | int i;
|
---|
282 | for (i = 0; i < n; i++) a[i] = 0;
|
---|
283 | }
|
---|
284 |
|
---|
285 | /** Differential step sum.
|
---|
286 | * Sums up entries from array *dx, taking each \a incx of it, \a n times.
|
---|
287 | * \param n number of steps
|
---|
288 | * \param *dx incremental value array
|
---|
289 | * \param incx step width
|
---|
290 | * \return sum_i+=incx dx[i]
|
---|
291 | * \sa Simps
|
---|
292 | */
|
---|
293 | #ifdef HAVE_INLINE
|
---|
294 | inline double dSum(int n, double *dx, int incx)
|
---|
295 | #else
|
---|
296 | double dSum(int n, double *dx, int incx)
|
---|
297 | #endif
|
---|
298 | {
|
---|
299 | int i;
|
---|
300 | double res;
|
---|
301 | if (n <= 0) return(0.0);
|
---|
302 | res = dx[0];
|
---|
303 | for(i = incx+1; i <= n*incx; i +=incx)
|
---|
304 | res += dx[i-1];
|
---|
305 | return (res);
|
---|
306 | }
|
---|
307 |
|
---|
308 | /** Simpson formula for integration.
|
---|
309 | * \a f is replaced by a polynomial of 2nd degree in order
|
---|
310 | * to approximate the integral
|
---|
311 | * \param n number of sampling points
|
---|
312 | * \param *f function value array
|
---|
313 | * \param h half the width of the integration interval
|
---|
314 | * \return \f$\int_a^b f(x) dx = \frac{h}{3} (y_0 + 4 y_1 + 2 y_2 + 4 y_3 + ... + 2 y_{n-2} + 4 y_{n-1} + y_n)\f$
|
---|
315 | * \sa dSum() - used by this function.
|
---|
316 | */
|
---|
317 | #ifdef HAVE_INLINE
|
---|
318 | inline double Simps(int n, double *f, double h)
|
---|
319 | #else
|
---|
320 | double Simps(int n, double *f, double h)
|
---|
321 | #endif
|
---|
322 | {
|
---|
323 | double res;
|
---|
324 | int nm12=(n-1)/2;
|
---|
325 | if (nm12*2 != n-1) {
|
---|
326 | fprintf(stderr,"Simps: wrong n in Simps");
|
---|
327 | }
|
---|
328 | res = 4.*dSum(nm12,&f[1],2)+2.*dSum(nm12-1,&f[2],2)+f[0]+f[n-1];
|
---|
329 | return(res*h/3.);
|
---|
330 | }
|
---|
331 |
|
---|
332 | /* derf */
|
---|
333 |
|
---|
334 | #ifndef HAVE_GSL_GSL_SF_ERF_H
|
---|
335 | /** Logarithm of Gamma function.
|
---|
336 | * \param xx x-value for function
|
---|
337 | * \return ln(gamma(xx))
|
---|
338 | * \note formula and coefficients are taken from "Numerical Receipes in C"
|
---|
339 | */
|
---|
340 | static double gammln(double xx) {
|
---|
341 | int j;
|
---|
342 | double x,tmp,ser;
|
---|
343 | double stp = 2.50662827465;
|
---|
344 | double cof[6] = { 76.18009173,-86.50532033,24.01409822,-1.231739516,.120858003e-2,-.536382e-5 };
|
---|
345 | x = xx -1.;
|
---|
346 | tmp = x+5.5;
|
---|
347 | tmp = (x+0.5)*log(tmp)-tmp;
|
---|
348 | ser = 1.;
|
---|
349 | for(j=0;j<6;j++) {
|
---|
350 | x+=1.0;
|
---|
351 | ser+=cof[j]/x;
|
---|
352 | }
|
---|
353 | return(tmp+log(stp*ser));
|
---|
354 | }
|
---|
355 |
|
---|
356 | /** Series used by gammp().
|
---|
357 | * \param a
|
---|
358 | * \param x
|
---|
359 | * \bug when x equals 0 is 0 returned?
|
---|
360 | * \note formula and coefficients are taken from "Numerical Receipes in C"
|
---|
361 | * \warning maximum precision 1e-7
|
---|
362 | */
|
---|
363 | static double gser(double a, double x) {
|
---|
364 | double gln = gammln(a);
|
---|
365 | double ap,sum,del;
|
---|
366 | int n;
|
---|
367 | if (x <= 0.) {
|
---|
368 | if (x < 0.) {
|
---|
369 | return(0.0);
|
---|
370 | }
|
---|
371 | }
|
---|
372 | ap=a;
|
---|
373 | sum=1./a;
|
---|
374 | del=sum;
|
---|
375 | for (n=1;n<=100;n++) {
|
---|
376 | ap += 1.;
|
---|
377 | del *=x/ap;
|
---|
378 | sum += del;
|
---|
379 | if(fabs(del) < fabs(sum)*1.e-7) {
|
---|
380 | return(sum*exp(-x+a*log(x)-gln));
|
---|
381 | }
|
---|
382 | }
|
---|
383 | return(sum*exp(-x+a*log(x)-gln));
|
---|
384 | }
|
---|
385 |
|
---|
386 | /** Continued fraction used by gammp().
|
---|
387 | * \param a
|
---|
388 | * \param x
|
---|
389 | * \note formula and coefficients are taken from "Numerical Receipes in C"
|
---|
390 | */
|
---|
391 | static double gcf(double a, double x) {
|
---|
392 | double gln = gammln(a);
|
---|
393 | double gold = 0.0;
|
---|
394 | double a0 = 1.;
|
---|
395 | double a1 = x;
|
---|
396 | double b0 = 0.;
|
---|
397 | double b1 = 1.;
|
---|
398 | double fac = 1.;
|
---|
399 | double an,ana,anf,g=0.0;
|
---|
400 | int n;
|
---|
401 | for (n=1; n <= 100; n++) {
|
---|
402 | an = n;
|
---|
403 | ana = an-a;
|
---|
404 | a0=(a1+a0*ana)*fac;
|
---|
405 | b0=(b1+b0*ana)*fac;
|
---|
406 | anf=an*fac;
|
---|
407 | a1=x*a0+anf*a1;
|
---|
408 | b1=x*b0+anf*b1;
|
---|
409 | if(a1 != 0.) {
|
---|
410 | fac=1./a1;
|
---|
411 | g=b1*fac;
|
---|
412 | if (fabs((g-gold)/g)<1.e-7) {
|
---|
413 | return(exp(-x+a*log(x)-gln)*g);
|
---|
414 | }
|
---|
415 | }
|
---|
416 | }
|
---|
417 | return(exp(-x+a*log(x)-gln)*g);
|
---|
418 | }
|
---|
419 |
|
---|
420 | /** Incomplete gamma function.
|
---|
421 | * Either calculated via series gser() or via continued fraction gcf()
|
---|
422 | * Needed by derf()
|
---|
423 | * \f[
|
---|
424 | * gammp(a,x) = \frac{1}{\gamma(a)} \int_x^\infty t^{a-1} \exp(-t) dt
|
---|
425 | * \f]
|
---|
426 | * \param a
|
---|
427 | * \param x
|
---|
428 | * \return f(a,x) = (x < 1+a) ? gser(a,x) : 1-gcf(a,x)
|
---|
429 | * \note formula and coefficients are taken from "Numerical Receipes in C"
|
---|
430 | */
|
---|
431 | static double gammp(double a, double x) {
|
---|
432 | double res;
|
---|
433 | if (x < a+1.) {
|
---|
434 | res = gser(a,x);
|
---|
435 | } else {
|
---|
436 | res = 1.-gcf(a,x);
|
---|
437 | }
|
---|
438 | return(res);
|
---|
439 | }
|
---|
440 | #endif
|
---|
441 |
|
---|
442 | /** Error function of integrated normal distribution.
|
---|
443 | * Either realized via GSL function gsl_sf_erf or via gammp()
|
---|
444 | * \f[
|
---|
445 | erf(x) = \frac{2}{\sqrt{\pi}} \int^x_0 \exp(-t^2) dt
|
---|
446 | = \pi^{-1/2} \gamma(\frac{1}{2},x^2)
|
---|
447 | * \f]
|
---|
448 | * \param x
|
---|
449 | * \return f(x) = sign(x) * gammp(0.5,x^2)
|
---|
450 | * \sa gammp
|
---|
451 | */
|
---|
452 | #ifdef HAVE_INLINE
|
---|
453 | inline double derf(double x)
|
---|
454 | #else
|
---|
455 | double derf(double x)
|
---|
456 | #endif
|
---|
457 | {
|
---|
458 | double res;
|
---|
459 | #ifdef HAVE_GSL_GSL_SF_ERF_H
|
---|
460 | // call gsl instead of numerical recipes routines
|
---|
461 | res = gsl_sf_erf(x);
|
---|
462 | #else
|
---|
463 | if (x < 0) {
|
---|
464 | res = -gammp(0.5,x*x);
|
---|
465 | } else {
|
---|
466 | res = gammp(0.5,x*x);
|
---|
467 | }
|
---|
468 | #endif
|
---|
469 | return(res);
|
---|
470 | }
|
---|
471 |
|
---|
472 | /** Sets array to zero.
|
---|
473 | * \param *a pointer to the double array
|
---|
474 | * \param n number of array elements
|
---|
475 | */
|
---|
476 | #ifdef HAVE_INLINE
|
---|
477 | inline void SetArrayToDouble0(double *a, int n)
|
---|
478 | #else
|
---|
479 | void SetArrayToDouble0(double *a, int n)
|
---|
480 | #endif
|
---|
481 | {
|
---|
482 | int i;
|
---|
483 | for(i=0;i<n;i++) a[i] = 0.0;
|
---|
484 | }
|
---|
485 |
|
---|
486 | /** Print complex 3x3 matrix.
|
---|
487 | * Checks if matrix has only zero entries, if not print each to screen: (re, im) ...
|
---|
488 | * \param M matrix array
|
---|
489 | */
|
---|
490 | void PrintCMat330(fftw_complex M[NDIM_NDIM])
|
---|
491 | {
|
---|
492 | int i,p=0;
|
---|
493 | for (i=0;i<NDIM_NDIM;i++)
|
---|
494 | if (M[i].re != 0.0 || M[i].im != 0.0) p++;
|
---|
495 | if (p) {
|
---|
496 | for (i=0;i<NDIM_NDIM;i++) fprintf(stderr," (%f %f)", M[i].re, M[i].im);
|
---|
497 | fprintf(stderr,"\n");
|
---|
498 | }
|
---|
499 | }
|
---|
500 |
|
---|
501 | /** Print real 3x3 matrix.
|
---|
502 | * Checks if matrix has only zero entries, if not print each to screen: re ...
|
---|
503 | * \param M matrix array
|
---|
504 | */
|
---|
505 | void PrintRMat330(fftw_real M[NDIM_NDIM])
|
---|
506 | {
|
---|
507 | int i,p=0;
|
---|
508 | for (i=0;i<NDIM_NDIM;i++)
|
---|
509 | if (M[i] != 0.0) p++;
|
---|
510 | if (p) {
|
---|
511 | for (i=0;i<NDIM_NDIM;i++) fprintf(stderr," %f", M[i]);
|
---|
512 | fprintf(stderr,"\n");
|
---|
513 | }
|
---|
514 | }
|
---|
515 |
|
---|
516 | /** Print complex 3-dim vector.
|
---|
517 | * Checks if vector has only zero entries, if not print each to screen: (re, im) ...
|
---|
518 | * \param M vector array
|
---|
519 | */
|
---|
520 | void PrintCVec30(fftw_complex M[NDIM])
|
---|
521 | {
|
---|
522 | int i,p=0;
|
---|
523 | for (i=0;i<NDIM;i++)
|
---|
524 | if (M[i].re != 0.0 || M[i].im != 0.0) p++;
|
---|
525 | if (p) {
|
---|
526 | for (i=0;i<NDIM;i++) fprintf(stderr," (%f %f)", M[i].re, M[i].im);
|
---|
527 | fprintf(stderr,"\n");
|
---|
528 | }
|
---|
529 | }
|
---|
530 |
|
---|
531 | /** Print real 3-dim vector.
|
---|
532 | * Checks if vector has only zero entries, if not print each to screen: re ...
|
---|
533 | * \param M matrix array
|
---|
534 | */
|
---|
535 | void PrintRVec30(fftw_real M[NDIM])
|
---|
536 | {
|
---|
537 | int i,p=0;
|
---|
538 | for (i=0;i<NDIM;i++)
|
---|
539 | if (M[i] != 0.0) p++;
|
---|
540 | if (p) {
|
---|
541 | for (i=0;i<NDIM;i++) fprintf(stderr," %f", M[i]);
|
---|
542 | fprintf(stderr,"\n");
|
---|
543 | }
|
---|
544 | }
|
---|
545 |
|
---|
546 | /** Rotates \a matrix, such that simultaneously given \a vector is aligned with z axis.
|
---|
547 | * Is used to rotate the unit cell in case of an external magnetic field. This field
|
---|
548 | * is rotated so that it aligns with z axis in order to simplify necessary perturbation
|
---|
549 | * calculations (only one component of each perturbed wave function necessary then).
|
---|
550 | * \param vector which is aligned with z axis by rotation \a Q
|
---|
551 | * \param Q return rotation matrix
|
---|
552 | * \param matrix which is transformed under the above rotation \a Q
|
---|
553 | */
|
---|
554 | void RotateToAlign(fftw_real Q[NDIM_NDIM], fftw_real matrix[NDIM_NDIM], fftw_real vector[NDIM]) {
|
---|
555 | double tmp[NDIM_NDIM], Q1[NDIM_NDIM], Qtmp[NDIM_NDIM];
|
---|
556 | double alpha, beta, new_y;
|
---|
557 | int i,j ;
|
---|
558 |
|
---|
559 | // calculate rotation angles
|
---|
560 | if (vector[0] < MYEPSILON) {
|
---|
561 | alpha = 0;
|
---|
562 | } else if (vector[1] > MYEPSILON) {
|
---|
563 | alpha = atan(-vector[0]/vector[1]);
|
---|
564 | } else alpha = PI/2;
|
---|
565 | new_y = -sin(alpha)*vector[0]+cos(alpha)*vector[1];
|
---|
566 | if (new_y < MYEPSILON) {
|
---|
567 | beta = 0;
|
---|
568 | } else if (vector[2] > MYEPSILON) {
|
---|
569 | beta = atan(-new_y/vector[2]);//asin(-vector[1]/vector[2]);
|
---|
570 | } else beta = PI/2;
|
---|
571 |
|
---|
572 | // create temporary matrix copy
|
---|
573 | // set Q to identity
|
---|
574 | for (i=0;i<NDIM;i++)
|
---|
575 | for (j=0;j<NDIM;j++) {
|
---|
576 | Q[i*NDIM+j] = (i == j) ? 1 : 0;
|
---|
577 | tmp[i*NDIM+j] = matrix[i*NDIM+j];
|
---|
578 | }
|
---|
579 |
|
---|
580 | // construct rotation matrices
|
---|
581 | Q1[0] = cos(alpha);
|
---|
582 | Q1[1] = sin(alpha);
|
---|
583 | Q1[2] = 0;
|
---|
584 | Q1[3] = -sin(alpha);
|
---|
585 | Q1[4] = cos(alpha);
|
---|
586 | Q1[5] = 0;
|
---|
587 | Q1[6] = 0;
|
---|
588 | Q1[7] = 0;
|
---|
589 | Q1[8] = 1;
|
---|
590 | // apply rotation and store
|
---|
591 | RMatMat33(tmp,Q1,matrix);
|
---|
592 | RMatMat33(Qtmp,Q1,Q);
|
---|
593 |
|
---|
594 | Q1[0] = 1;
|
---|
595 | Q1[1] = 0;
|
---|
596 | Q1[2] = 0;
|
---|
597 | Q1[3] = 0;
|
---|
598 | Q1[4] = cos(beta);
|
---|
599 | Q1[5] = sin(beta);
|
---|
600 | Q1[6] = 0;
|
---|
601 | Q1[7] = -sin(beta);
|
---|
602 | Q1[8] = cos(beta);
|
---|
603 | // apply rotation and store
|
---|
604 | RMatMat33(matrix,Q1,tmp);
|
---|
605 | RMatMat33(Q,Q1,Qtmp);
|
---|
606 |
|
---|
607 | // in order to avoid unncessary calculations, set everything below epsilon to zero
|
---|
608 | for (i=0;i<NDIM_NDIM;i++) {
|
---|
609 | matrix[i] = (fabs(matrix[i]) > MYEPSILON) ? matrix[i] : 0;
|
---|
610 | Q[i] = (fabs(Q[i]) > MYEPSILON) ? Q[i] : 0;
|
---|
611 | }
|
---|
612 | }
|
---|