valib
Vortex Analysis LIBrary
Loading...
Searching...
No Matches
linalg3.h
Go to the documentation of this file.
1
5#ifndef VA_LINALG3_H
6#define VA_LINALG3_H
7
8#include "common_defs.h"
9
10/***************************************************************************/
15/***************************************************************************/
20VA_DEVICE_FUN VA_REAL valib_sign(VA_REAL a, VA_REAL b)
21{
22 return (b < (VA_REAL) 0.) ? -fabs(a) : fabs(a);
23}
24
25/***************************************************************************/
29VA_DEVICE_FUN VA_REAL valib_cabs(VA_DOUBLE_COMPLEX a)
30{
31#if defined(CUDA)
32 return abs(a);
33#elif defined(HIP)
34 VA_REAL ar = hipCreal(a);
35 VA_REAL ai = hipCimag(a);
36 return sqrt(ar*ar + ai*ai);
37#else
38 return cabs(a);
39#endif
40}
41
42/***************************************************************************/
46VA_DEVICE_FUN VA_DOUBLE_COMPLEX valib_cmultiply(VA_DOUBLE_COMPLEX a, VA_DOUBLE_COMPLEX b)
47{
48#if defined(CUDA)
49 return a*b;
50#elif defined(HIP)
51 double ar = hipCreal(a);
52 double ai = hipCimag(a);
53 double br = hipCreal(b);
54 double bi = hipCimag(b);
55
56 double cr = ar*br - ai*bi;
57 double ci = ar*bi + ai*br;
58 VA_DOUBLE_COMPLEX c = {cr,ci};
59 return c;
60#else
61 return a*b;
62#endif
63}
64
65/***************************************************************************/
69VA_DEVICE_FUN VA_REAL valib_cubic_root(VA_REAL x)
70{
71 return valib_sign(pow(fabs(x), 1.0/3.0), x);
72}
73
74/***************************************************************************/
79VA_DEVICE_FUN VA_REAL valib_max(VA_REAL a, VA_REAL b)
80{
81 return (a > b) ? a : b;
82}
83
84/***************************************************************************/
89VA_DEVICE_FUN VA_REAL valib_min(VA_REAL a, VA_REAL b)
90{
91 return (a < b) ? a : b;
92}
93
94/***************************************************************************/
98VA_DEVICE_FUN VA_REAL valib_square(VA_REAL a)
99{
100 return (a * a);
101}
102
103/***************************************************************************/
107VA_DEVICE_FUN void valib_swap_values(VA_REAL *a, VA_REAL *b)
108{
109 VA_REAL buf;
110 buf = *a;
111 *a = *b;
112 *b = buf;
113}
114
115/***************************************************************************/
123VA_DEVICE_FUN void valib_solve_cubic_equation(VA_REAL a, VA_REAL b,
124 VA_REAL c, VA_REAL d,
125 VA_COMPLEX *x1,
126 VA_COMPLEX *x2,
127 VA_COMPLEX *x3,
128 int *info)
129{
130 int debug = 0;
131
132 // This computation only works in double precision.
133 double ad = a;
134 double bd = b;
135 double cd = c;
136 double dd = d;
137
138 double q = pow(bd/(3.*ad),3) - (bd*cd)/(6.*ad*ad) + dd/(2.*ad);
139 double p = (3*ad*cd - bd*bd) / (9*ad*ad);
140
141 double delta = - (p*p*p + q*q);
142
143 if (debug) {
144 printf("q: %lf \n", q);
145 printf("p: %lf \n", p);
146 printf("delta: %lf \n", delta);
147 }
148
149 if (delta <= 0) {
150 if (debug) {
151 printf("This is a good case with two complex roots.\n");
152 }
153 }
154 else {
155 if (debug) {
156 printf("Discriminant is positive, this is not a good way to get roots.\n");
157 }
158 *info = -1;
159 return;
160 }
161
162 double sr = sqrt(-delta);
163
164 if (debug) {
165 printf("sr: %lf \n", sr);
166 }
167
168 double mqpsr = -q + sr;
169 double mqmsr = -q - sr;
170
171 if (debug) {
172 printf("mqpsr: %lf \n", mqpsr);
173 printf("mqmsr: %lf \n", mqmsr);
174 }
175
176 double u = valib_cubic_root(mqpsr);
177 double v = valib_cubic_root(mqmsr);
178
179 if (debug) {
180 printf("fabs(u*v + p): %lf \n", fabs(u*v + p));
181 }
182 if (fabs(u*v + p) > 1.e-4) {
183 if (debug) {
184 printf("Check of uv = -p failed, check = %lf.", fabs(u*v + p));
185 }
186 *info = -2;
187 return;
188 }
189
190// define imaginary unit for CUDA or HIP
191#if defined(CUDA)
192 VA_DOUBLE_COMPLEX e1(-0.5, sqrt(3.)/2.);
193 VA_DOUBLE_COMPLEX e2(-0.5, -sqrt(3.)/2.);
194
195 VA_DOUBLE_COMPLEX uc = u;
196 VA_DOUBLE_COMPLEX vc = v;
197#elif defined(HIP)
198 VA_DOUBLE_COMPLEX e1 = {-0.5, sqrt(3.)/2.};
199 VA_DOUBLE_COMPLEX e2 = {-0.5, -sqrt(3.)/2.};
200
201 VA_DOUBLE_COMPLEX uc = {u, 0};
202 VA_DOUBLE_COMPLEX vc = {v, 0};
203#else
204 VA_DOUBLE_COMPLEX e1 = -0.5 + sqrt(3.)/2.*I;
205 VA_DOUBLE_COMPLEX e2 = -0.5 - sqrt(3.)/2.*I;
206
207 VA_DOUBLE_COMPLEX uc = u + 0.*I;
208 VA_DOUBLE_COMPLEX vc = v + 0.*I;
209#endif
210
211 VA_DOUBLE_COMPLEX y1 = uc + vc;
212 VA_DOUBLE_COMPLEX y2 = valib_cmultiply(e1,uc) + valib_cmultiply(e2,vc);
213 VA_DOUBLE_COMPLEX y3 = valib_cmultiply(e2,uc) + valib_cmultiply(e1,vc);
214
215 double shift = -bd/(3.*ad);
216#if defined(CUDA)
217 VA_DOUBLE_COMPLEX shiftc = shift;
218#elif defined(HIP)
219 VA_DOUBLE_COMPLEX shiftc = {shift,0};
220#else
221 VA_DOUBLE_COMPLEX shiftc = shift + 0.*I;
222#endif
223
224 *x1 = y1 + shiftc;
225 *x2 = y2 + shiftc;
226 *x3 = y3 + shiftc;
227
228#if defined(CUDA)
229 VA_DOUBLE_COMPLEX ac = a;
230 VA_DOUBLE_COMPLEX bc = b;
231 VA_DOUBLE_COMPLEX cc = c;
232 VA_DOUBLE_COMPLEX dc = d;
233#elif defined(HIP)
234 VA_DOUBLE_COMPLEX ac = {a,0};
235 VA_DOUBLE_COMPLEX bc = {b,0};
236 VA_DOUBLE_COMPLEX cc = {c,0};
237 VA_DOUBLE_COMPLEX dc = {d,0};
238#else
239 VA_DOUBLE_COMPLEX ac = a + 0.*I;
240 VA_DOUBLE_COMPLEX bc = b + 0.*I;
241 VA_DOUBLE_COMPLEX cc = c + 0.*I;
242 VA_DOUBLE_COMPLEX dc = d + 0.*I;
243#endif
244
245 // check that ax^3 + bx^2 + cx + d is close to zero for the roots
246 VA_DOUBLE_COMPLEX check1 = valib_cmultiply(ac,valib_cmultiply(valib_cmultiply((*x1),(*x1)),(*x1))) + valib_cmultiply(bc,valib_cmultiply((*x1),(*x1))) + valib_cmultiply(cc,(*x1)) + dc;
247 VA_DOUBLE_COMPLEX check2 = valib_cmultiply(ac,valib_cmultiply(valib_cmultiply((*x2),(*x2)),(*x2))) + valib_cmultiply(bc,valib_cmultiply((*x2),(*x2))) + valib_cmultiply(cc,(*x2)) + dc;
248 VA_DOUBLE_COMPLEX check3 = valib_cmultiply(ac,valib_cmultiply(valib_cmultiply((*x3),(*x3)),(*x3))) + valib_cmultiply(bc,valib_cmultiply((*x3),(*x3))) + valib_cmultiply(cc,(*x3)) + dc;
249
250 if (valib_cabs(check1) > 1.e-5) {
251 if (debug) {
252 printf("Check of x1 failed.");
253 }
254 }
255 if (valib_cabs(check2) > 1.e-5) {
256 if (debug) {
257 printf("Check of x2 failed.");
258 }
259 }
260 if (valib_cabs(check3) > 1.e-5) {
261 if (debug) {
262 printf("Check of x3 failed.");
263 }
264 }
265
266 // finished correctly
267 *info = 0;
268}
269
270/***************************************************************************/
277/***************************************************************************/
281VA_DEVICE_FUN void valib_vec_zero3(VA_REAL *v)
282{
283 int i;
284 for (i = 0; i < 3; i++) {
285 v[i] = 0.;
286 }
287}
288
289/***************************************************************************/
293VA_DEVICE_FUN VA_REAL valib_scalar_prod3(VA_REAL *v1, VA_REAL *v2)
294{
295 VA_REAL result = (VA_REAL) 0.;
296 int i;
297 for (i = 0; i < 3; i++)
298 {
299 result = result + v1[i]*v2[i];
300 }
301 return result;
302}
303
304/***************************************************************************/
308VA_DEVICE_FUN void valib_cross_prod3(VA_REAL *v1, VA_REAL *v2, VA_REAL *v3)
309{
310 v3[0] = v1[1] * v2[2] - v1[2] * v2[1];
311 v3[1] = v1[2] * v2[0] - v1[0] * v2[2];
312 v3[2] = v1[0] * v2[1] - v1[1] * v2[0];
313}
314
315/***************************************************************************/
322/***************************************************************************/
326VA_DEVICE_FUN void valib_mat_zero3(VA_REAL *A)
327{
328 int i, j;
329 for (j = 0; j < 3; j++) {
330 for (i = 0; i < 3; i++) {
331 A[j*3+i] = 0.;
332 }
333 }
334}
335
336/***************************************************************************/
340VA_DEVICE_FUN void valib_mat_copy3(VA_REAL *A, VA_REAL *B)
341{
342 int i, j;
343 for (j = 0; j < 3; j++) {
344 for (i = 0; i < 3; i++) {
345 B[j*3+i] = A[j*3+i];
346 }
347 }
348}
349
350/***************************************************************************/
354VA_DEVICE_FUN void valib_norm_frobenius3(VA_REAL *A, VA_REAL *norm)
355{
356 *norm = 0.0;
357 int i, j;
358 for (j = 0; j < 3; j++) {
359 for (i = 0; i < 3; i++) {
360 *norm = *norm + A[j*3+i]*A[j*3+i];
361 }
362 }
363 *norm = sqrt(*norm);
364}
365
366/***************************************************************************/
370VA_DEVICE_FUN void valib_trace3(VA_REAL *A, VA_REAL *tr3)
371{
372 *tr3 = A[0] + A[4] + A[8];
373}
374
375/***************************************************************************/
379VA_DEVICE_FUN void valib_second_invariant3(VA_REAL *A,
380 VA_REAL *second_invariant)
381{
382 *second_invariant = ( A[0]*A[4] + A[4]*A[8] + A[0]*A[8] )
383 - ( A[1]*A[3] + A[2]*A[6] + A[5]*A[7] );
384}
385
386/***************************************************************************/
390VA_DEVICE_FUN void valib_determinant3(VA_REAL *A, VA_REAL *det)
391{
392 *det = A[0] * A[4] * A[8]
393 + A[3] * A[7] * A[2]
394 + A[6] * A[1] * A[5]
395 - A[0] * A[7] * A[5]
396 - A[3] * A[1] * A[8]
397 - A[6] * A[4] * A[2];
398}
399
400/***************************************************************************/
404VA_DEVICE_FUN void valib_transpose3(VA_REAL *A)
405{
406 valib_swap_values(&A[1], &A[3]);
407 valib_swap_values(&A[2], &A[6]);
408 valib_swap_values(&A[5], &A[7]);
409}
410
411/***************************************************************************/
417VA_DEVICE_FUN void valib_sym_antisym3(VA_REAL *A)
418{
419 VA_REAL ud, ld;
420
421 int i, j;
422 for (i = 0; i < 3; i++) {
423 for (j = i+1; j < 3; j++) {
424 ud = A[i + j*3];
425 ld = A[j + i*3];
426
427 A[i + j*3] = (VA_REAL) 0.5 * (ud + ld);
428 A[j + i*3] = (VA_REAL) 0.5 * (ld - ud);
429 }
430 }
431}
432
433/***************************************************************************/
441VA_DEVICE_FUN void valib_deviatorise3(VA_REAL *A)
442{
443 // get trace(A)
444 VA_REAL trace;
445 valib_trace3(A, &trace);
446
447 // subtract a multiple of the third of the trace from the diagonal
448 A[0] = A[0] - (VA_REAL) (1./3.) * trace;
449 A[4] = A[4] - (VA_REAL) (1./3.) * trace;
450 A[8] = A[8] - (VA_REAL) (1./3.) * trace;
451}
452
453/***************************************************************************/
458VA_DEVICE_FUN void valib_matvec3_prod(int trans_a, VA_REAL *A,
459 VA_REAL *x, VA_REAL *y)
460{
461 int i, j;
462 for (i = 0; i < 3; i++) {
463 y[i] = 0.0;
464 for (j = 0; j < 3; j++) {
465 if (trans_a == 0) {
466 y[i] = y[i] + A[j*3+i]*x[j];
467 }
468 else {
469 y[i] = y[i] + A[i*3+j]*x[j];
470 }
471 }
472 }
473}
474
475/***************************************************************************/
480VA_DEVICE_FUN void valib_rank1_update3(VA_REAL alpha,
481 VA_REAL *x, VA_REAL *y,
482 VA_REAL *A)
483{
484 int ind;
485 VA_REAL tmp;
486
487 ind = 0;
488 int i, j;
489 for (j = 0; j < 3; j++) {
490 tmp = alpha * y[j];
491 for (i = 0; i < 3; i++) {
492 A[ind] = A[ind] + x[i]*tmp;
493 ind = ind + 1;
494 }
495 }
496}
497
498/***************************************************************************/
503VA_DEVICE_FUN void valib_tatt3(VA_REAL *A, VA_REAL *T, VA_REAL *TATT)
504{
505 unsigned i, j, ind;
506
507 // zero TATT
508 valib_mat_zero3(TATT);
509
510 // perform series of rank one updates
511 ind = 0;
512 for (i = 0; i < 3; i++) {
513 for (j = 0; j < 3; j++) {
514 // call rank one update on each element of A
515 valib_rank1_update3(A[ind], &T[j*3], &T[i*3], TATT);
516 ind = ind + 1;
517 }
518 }
519}
520
521/***************************************************************************/
526VA_DEVICE_FUN void valib_matmat3_prod(int trans_a, int trans_b,
527 VA_REAL *A, VA_REAL *B, VA_REAL *C)
528{
529 int i, j, k;
530 for (i = 0; i < 3; i++) {
531 for (j = 0; j < 3; j++) {
532 for (k = 0; k < 3; k++) {
533 if (trans_a == 0 && trans_b == 0) {
534 C[j*3+i] = C[j*3+i] + A[k*3+i] * B[j*3+k];
535 }
536 else if (trans_a != 0 && trans_b == 0) {
537 C[j*3+i] = C[j*3+i] + A[i*3+k] * B[j*3+k];
538 }
539 else if (trans_a == 0 && trans_b != 0) {
540 C[j*3+i] = C[j*3+i] + A[k*3+i] * B[k*3+j];
541 }
542 else {//(trans_a != 0 && trans_b != 0)
543 C[j*3+i] = C[j*3+i] + A[i*3+k] * B[k*3+j];
544 }
545 }
546 }
547 }
548}
549
550/***************************************************************************/
555VA_DEVICE_FUN void valib_ttat3(VA_REAL *A, VA_REAL *T, VA_REAL *TTAT)
556{
557 // zero TATT
558 valib_mat_zero3(TTAT);
559
560 // first multiply A*T
561 VA_REAL AT[3*3];
562 valib_mat_zero3(AT);
563 valib_matmat3_prod(0, 0, A, T, AT);
564
565 // now premultiply T^T * (A*T)
566 valib_matmat3_prod(1, 0, T, AT, TTAT);
567}
568
569/***************************************************************************/
573VA_DEVICE_FUN void valib_gram_schmidt3(VA_REAL *Q)
574{
575 VA_REAL factor;
576
577 // column-wise process
578 int i, j, k;
579 for (j = 0; j < 3; j++) {
580 for (i = 0; i < j; i++) {
581 // compute the scalar product with previous columns
582 factor = valib_scalar_prod3(&Q[j*3], &Q[i*3]);
583
584 // subtract the projected vector
585 for (k = 0; k < 3; k++) {
586 Q[j*3+k] = Q[j*3+k] - factor * Q[i*3+k];
587 }
588 }
589
590 // normalize the column
591 factor = valib_scalar_prod3(&Q[j*3], &Q[j*3]);
592 for (k = 0; k < 3; k++) {
593 Q[j*3+k] = Q[j*3+k] / sqrt(factor);
594 }
595 }
596}
597
598/***************************************************************************/
602VA_DEVICE_FUN void valib_eigenvalues_sym3(VA_REAL *A, VA_REAL *eigval)
603{
604 // compute first invariant
605 VA_REAL Idev;
606 valib_trace3(A, &Idev);
607 VA_REAL IdevMod = Idev / 3.;
608
609 // compute deviator
610 VA_REAL dev[9];
611 valib_mat_copy3(A, dev);
612 dev[0] = dev[0] - IdevMod;
613 dev[4] = dev[4] - IdevMod;
614 dev[8] = dev[8] - IdevMod;
615
616 // compute second invariant of deviator
617 VA_REAL IIdev;
618 valib_second_invariant3(dev, &IIdev);
619 VA_REAL IIdevMod = -IIdev / 3.;
620
621 // compute half the determinant (third invariant) of deviator
622 VA_REAL IIIdev;
623 valib_determinant3(dev, &IIIdev);
624 VA_REAL IIIdevMod = 0.5 * IIIdev;
625
626 // compute the phi angle
627 VA_REAL aux = IIdevMod*IIdevMod*IIdevMod - IIIdevMod*IIIdevMod;
628 // assure nonnegativity
629 if (aux < 0.) {
630 aux = 0.;
631 }
632 VA_REAL phi;
633 phi = atan2( sqrt(aux), IIIdevMod ) / 3.;
634 // readjust to the right quadrant
635 if (phi < 0.) phi += VA_PI / 3.;
636
637 // eigenvalues
638 eigval[0] = IdevMod + 2.*sqrt(IIdevMod) * cos(phi);
639 eigval[1] = IdevMod - sqrt(IIdevMod) * (cos(phi) + sqrt(3.) * sin(phi));
640 eigval[2] = IdevMod - sqrt(IIdevMod) * (cos(phi) - sqrt(3.) * sin(phi));
641
642 // sort eigenvalues
643 if (eigval[0] > eigval[1] ) valib_swap_values( &eigval[0], &eigval[1] );
644 if (eigval[1] > eigval[2] ) valib_swap_values( &eigval[1], &eigval[2] );
645 if (eigval[0] > eigval[1] ) valib_swap_values( &eigval[0], &eigval[1] );
646
647 return;
648}
649
650/***************************************************************************/
656VA_DEVICE_FUN int valib_find_pivot3(VA_REAL *A, int jcol, int start)
657{
658 int debug = 0;
659 int ipiv = start;
660 for (int i = start+1; i<3; i++) {
661 if (fabs(A[i+jcol*3]) > fabs(A[ipiv+jcol*3]))
662 ipiv = i;
663 }
664
665 VA_REAL numerical_zero = 1.0e-5;
666 if (fabs(A[ipiv+jcol*3]) < numerical_zero) {
667 if (debug)
668 printf("Warning: no nonzero pivot in column %d.\n", jcol);
669
670 ipiv = -1;
671 }
672
673 return ipiv;
674}
675
676/***************************************************************************/
683VA_DEVICE_FUN void valib_find_nullspace3(VA_REAL *A, VA_REAL *x, int *info)
684{
685 int debug = 0;
686
687 // zero the vector
688 x[0] = 0.; x[1] = 0.; x[2] = 0.;
689
690 // handle special cases with zero column
691 int num_zero_cols = 0;
692 int i_zero_col = -1;
693 int ip;
694 for (int ji = 0; ji < 3; ji++) {
695 int j = 2-ji;
696 ip = valib_find_pivot3(A,j,0);
697 if (ip == -1) {
698 if (debug) {
699 printf("There is no meaningful pivot in the %d-th column.\n",j);
700 }
701 num_zero_cols++;
702 i_zero_col = j;
703 }
704 }
705
706 if (num_zero_cols > 0) {
707 if (num_zero_cols > 1) {
708 if (debug) {
709 printf("Oops! The matrix has more than one zero columns. I search only for 1-dimensional nullspace.\n");
710 }
711 *info = -1;
712 return;
713 }
714 else {
715 x[i_zero_col] = 1.;
716 *info = 0;
717 return;
718 }
719 }
720
721 // in the usual case, permute the rows and eliminate first column
722 if (ip != 0) {
723 for (int j = 0; j < 3; j++) {
724 valib_swap_values(&A[j*3], &A[ip+j*3]);
725 }
726 }
727
728 if (debug) {
729 printf("Matrix in GEM after 1st swap: \n");
730 for (int i = 0; i < 3; i++) {
731 for (int j = 0; j < 3; j++) {
732 printf(" %lf ", A[i+j*3]);
733 }
734 printf("\n");
735 }
736 }
737
738 VA_REAL coef;
739
740 // eliminate the first column using the pivot
741 for (int i = 1; i < 3; i++) {
742 coef = A[i] / A[0];
743 for (int j = 0; j < 3; j++) {
744 A[i+j*3] -= coef*A[j*3];
745 }
746 }
747
748 if (debug) {
749 printf("Matrix in GEM without 1st col: \n");
750 for (int i = 0; i < 3; i++) {
751 for (int j = 0; j < 3; j++) {
752 printf(" %lf ", A[i+j*3]);
753 }
754 printf("\n");
755 }
756 }
757
758 // now work on the second column
759 int jpiv, jfree;
760 ip = valib_find_pivot3(A,1,1);
761 if (ip == -1) {
762 // the second column is not a pivot one, we can search for the nullspace
763 jpiv = 2;
764 jfree = 1;
765 }
766 else {
767 jpiv = 1;
768 jfree = 2;
769 }
770
771 ip = valib_find_pivot3(A,jpiv,1);
772 if (ip == -1) {
773 if (debug) {
774 printf("Something is wrong with the matrix, the defect seems larger than one.");
775 }
776 *info = -2;
777 return;
778 }
779 if (ip == 2) {
780 // swap the rows
781 for (int j = 0; j < 3; j++) {
782 valib_swap_values(&A[1+j*3], &A[2+j*3]);
783 }
784 }
785
786 // eliminate the last row
787 coef = A[2+jpiv*3] / A[1+jpiv*3];
788 for (int j = 0; j < 3; j++) {
789 A[2+j*3] -= coef*A[1+j*3];
790 }
791 // eliminate the first row
792 coef = A[0+jpiv*3] / A[1+jpiv*3];
793 for (int j = 0; j < 3; j++) {
794 A[j*3] -= coef*A[1+j*3];
795 }
796
797 // check that the last row, in particular the last entry, is zero
798 if (valib_find_pivot3(A,2,2) != -1) {
799 if (debug) {
800 printf("Something is wrong, the matrix seems to be full rank.\n");
801 }
802 *info = -2;
803 return;
804 }
805
806 if (debug) {
807 printf("Matrix in GEM before dividing by diagonal: \n");
808 for (int i = 0; i < 3; i++) {
809 for (int j = 0; j < 3; j++) {
810 printf(" %lf ", A[i+j*3]);
811 }
812 printf("\n");
813 }
814 }
815
816 // finalize the Gauss-Jordan elimination by dividing with the diagonal values
817 coef = A[0];
818 for (int j = 0; j < 3; j++) {
819 A[j*3] /= coef;
820 }
821 coef = A[1+jpiv*3];
822 for (int j = 0; j < 3; j++) {
823 A[1+j*3] /= coef;
824 }
825
826 if (debug) {
827 printf("Matrix in GEM in the upper right echelon form:\n");
828 for (int i = 0; i < 3; i++) {
829 for (int j = 0; j < 3; j++) {
830 printf(" %d %lf ", i+j*3, A[i+j*3]);
831 }
832 printf("\n");
833 }
834 }
835
836 // fill the nullspace basis vector
837 x[jfree] = 1.;
838 x[jpiv] = -A[jpiv+jfree*3];
839 x[0] = -A[jfree*3];
840
841 // normalize the basis vector
842 VA_REAL norm = sqrt(valib_scalar_prod3(&x[0], &x[0]));
843 for (int i = 0; i < 3; i++) {
844 x[i] /= norm;
845 }
846
847 *info = 0;
848}
849
850
851/***************************************************************************/
857VA_DEVICE_FUN void valib_constructQ3(VA_REAL sina, VA_REAL cosa,
858 VA_REAL sinb, VA_REAL cosb,
859 VA_REAL sing, VA_REAL cosg,
860 VA_REAL *Q )
861{
862 Q[0] = cosa*cosb*cosg - sina*sing;
863 Q[1] = -cosa*cosb*sing - sina*cosg;
864 Q[2] = cosa*sinb;
865 Q[3] = sina*cosb*cosg + cosa*sing;
866 Q[4] = -sina*cosb*sing + cosa*cosg;
867 Q[5] = sina*sinb;
868 Q[6] = -sinb*cosg;
869 Q[7] = sinb*sing;
870 Q[8] = cosb;
871}
872
873/***************************************************************************/
879VA_DEVICE_FUN void valib_constructQfromN3(VA_REAL *n, VA_REAL *Q)
880{
881 VA_REAL vec[3];
882 VA_REAL factor;
883
884 // initialize the first row of Q with n
885 Q[0] = n[0]; Q[1] = n[1]; Q[2] = n[2];
886
887 // fill some entries to initialize the other entries in matrix Q
888 Q[3] = n[1]; Q[4] = n[2]; Q[5] = n[0];
889 Q[6] = n[2]; Q[7] = n[0]; Q[8] = n[1];
890
891 // perform Gram-Schmidt orthogonalization
893
894 // swap the first and the third column - at the moment, the matrix contains
895 // | z y x |
896 int i;
897 for (i = 0; i < 3; i++) {
898 valib_swap_values(&Q[i], &Q[6+i]);
899 }
900
901 // verify that the matrix has a right-hand oriented coordinate system,
902 // Q_1 x Q_2 = Q_3
903 valib_cross_prod3(&Q[0], &Q[3], vec);
904
905 // revert sign of the new x-axis to keep right-hand oriented coordinate
906 // system
907 factor = valib_scalar_prod3(vec, &Q[6]);
908 if (factor < 0.) {
909 for (i = 0; i < 3; i++) {
910 Q[i] = -Q[i];
911 }
912 }
913
914 // transpose the matrix
916}
917
918/***************************************************************************/
922#endif // VA_LINALG3_H
VA_DEVICE_FUN VA_REAL valib_cubic_root(VA_REAL x)
Definition linalg3.h:69
VA_DEVICE_FUN void valib_cross_prod3(VA_REAL *v1, VA_REAL *v2, VA_REAL *v3)
Definition linalg3.h:308
VA_DEVICE_FUN void valib_tatt3(VA_REAL *A, VA_REAL *T, VA_REAL *TATT)
Definition linalg3.h:503
VA_DEVICE_FUN void valib_solve_cubic_equation(VA_REAL a, VA_REAL b, VA_REAL c, VA_REAL d, VA_COMPLEX *x1, VA_COMPLEX *x2, VA_COMPLEX *x3, int *info)
Definition linalg3.h:123
VA_DEVICE_FUN void valib_rank1_update3(VA_REAL alpha, VA_REAL *x, VA_REAL *y, VA_REAL *A)
Definition linalg3.h:480
VA_DEVICE_FUN void valib_ttat3(VA_REAL *A, VA_REAL *T, VA_REAL *TTAT)
Definition linalg3.h:555
VA_DEVICE_FUN void valib_deviatorise3(VA_REAL *A)
Definition linalg3.h:441
VA_DEVICE_FUN VA_REAL valib_scalar_prod3(VA_REAL *v1, VA_REAL *v2)
Definition linalg3.h:293
VA_DEVICE_FUN void valib_vec_zero3(VA_REAL *v)
Definition linalg3.h:281
VA_DEVICE_FUN void valib_matvec3_prod(int trans_a, VA_REAL *A, VA_REAL *x, VA_REAL *y)
Definition linalg3.h:458
VA_DEVICE_FUN VA_REAL valib_min(VA_REAL a, VA_REAL b)
Definition linalg3.h:89
VA_DEVICE_FUN void valib_sym_antisym3(VA_REAL *A)
Definition linalg3.h:417
VA_DEVICE_FUN void valib_transpose3(VA_REAL *A)
Definition linalg3.h:404
VA_DEVICE_FUN void valib_mat_copy3(VA_REAL *A, VA_REAL *B)
Definition linalg3.h:340
VA_DEVICE_FUN VA_DOUBLE_COMPLEX valib_cmultiply(VA_DOUBLE_COMPLEX a, VA_DOUBLE_COMPLEX b)
Definition linalg3.h:46
VA_DEVICE_FUN void valib_swap_values(VA_REAL *a, VA_REAL *b)
Definition linalg3.h:107
VA_DEVICE_FUN VA_REAL valib_square(VA_REAL a)
Definition linalg3.h:98
VA_DEVICE_FUN void valib_second_invariant3(VA_REAL *A, VA_REAL *second_invariant)
Definition linalg3.h:379
VA_DEVICE_FUN void valib_find_nullspace3(VA_REAL *A, VA_REAL *x, int *info)
Definition linalg3.h:683
VA_DEVICE_FUN void valib_gram_schmidt3(VA_REAL *Q)
Definition linalg3.h:573
VA_DEVICE_FUN void valib_constructQfromN3(VA_REAL *n, VA_REAL *Q)
Definition linalg3.h:879
VA_DEVICE_FUN int valib_find_pivot3(VA_REAL *A, int jcol, int start)
Definition linalg3.h:656
VA_DEVICE_FUN void valib_mat_zero3(VA_REAL *A)
Definition linalg3.h:326
VA_DEVICE_FUN VA_REAL valib_sign(VA_REAL a, VA_REAL b)
Definition linalg3.h:20
VA_DEVICE_FUN VA_REAL valib_cabs(VA_DOUBLE_COMPLEX a)
Definition linalg3.h:29
VA_DEVICE_FUN void valib_eigenvalues_sym3(VA_REAL *A, VA_REAL *eigval)
Definition linalg3.h:602
VA_DEVICE_FUN void valib_constructQ3(VA_REAL sina, VA_REAL cosa, VA_REAL sinb, VA_REAL cosb, VA_REAL sing, VA_REAL cosg, VA_REAL *Q)
Definition linalg3.h:857
VA_DEVICE_FUN void valib_matmat3_prod(int trans_a, int trans_b, VA_REAL *A, VA_REAL *B, VA_REAL *C)
Definition linalg3.h:526
VA_DEVICE_FUN void valib_trace3(VA_REAL *A, VA_REAL *tr3)
Definition linalg3.h:370
VA_DEVICE_FUN void valib_norm_frobenius3(VA_REAL *A, VA_REAL *norm)
Definition linalg3.h:354
VA_DEVICE_FUN void valib_determinant3(VA_REAL *A, VA_REAL *det)
Definition linalg3.h:390
VA_DEVICE_FUN VA_REAL valib_max(VA_REAL a, VA_REAL b)
Definition linalg3.h:79