GRASS Programmer's Manual  6.5.svn(2014)-r66266
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
la.c
Go to the documentation of this file.
1 
2 /******************************************************************************
3  * la.c
4  * wrapper modules for linear algebra problems
5  * linking to BLAS / LAPACK (and others?)
6 
7  * @Copyright David D.Gray <ddgray@armadce.demon.co.uk>
8  * 26th. Sep. 2000
9  * Last updated:
10  * 2006-11-23
11 
12  * This file is part of GRASS GIS. It is free software. You can
13  * redistribute it and/or modify it under the terms of
14  * the GNU General Public License as published by the Free Software
15  * Foundation; either version 2 of the License, or (at your option)
16  * any later version.
17 
18  * This program is distributed in the hope that it will be useful,
19  * but WITHOUT ANY WARRANTY; without even the implied warranty of
20  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21  * GNU General Public License for more details.
22 
23  ******************************************************************************/
24 
25 #include <stdio.h> /* needed here for ifdef/else */
26 #include <stdlib.h>
27 #include <string.h>
28 #include <math.h>
29 
30 #include <grass/config.h>
31 
32 #if defined(HAVE_LIBLAPACK) && defined(HAVE_LIBBLAS)
33 
34 #include <grass/gis.h>
35 #include <grass/glocale.h>
36 #include <grass/la.h>
37 
38 
39 static int egcmp(const void *pa, const void *pb);
40 
41 
55 mat_struct *G_matrix_init(int rows, int cols, int ldim)
56 {
57  mat_struct *tmp_arry;
58 
59  if (rows < 1 || cols < 1 || ldim < rows) {
60  G_warning(_("Matrix dimensions out of range"));
61  return NULL;
62  }
63 
64  tmp_arry = (mat_struct *) G_malloc(sizeof(mat_struct));
65  tmp_arry->rows = rows;
66  tmp_arry->cols = cols;
67  tmp_arry->ldim = ldim;
68  tmp_arry->type = MATRIX_;
69  tmp_arry->v_indx = -1;
70 
71  tmp_arry->vals = (doublereal *) G_calloc(ldim * cols, sizeof(doublereal));
72  tmp_arry->is_init = 1;
73 
74  return tmp_arry;
75 }
76 
77 
87 int G_matrix_zero(mat_struct * A)
88 {
89  if (!A->vals)
90  return 0;
91 
92  memset(A->vals, 0, sizeof(A->vals));
93 
94  return 1;
95 }
96 
97 
113 int G_matrix_set(mat_struct * A, int rows, int cols, int ldim)
114 {
115  if (rows < 1 || cols < 1 || ldim < 0) {
116  G_warning(_("Matrix dimensions out of range"));
117  return -1;
118  }
119 
120  A->rows = rows;
121  A->cols = cols;
122  A->ldim = ldim;
123  A->type = MATRIX_;
124  A->v_indx = -1;
125 
126  A->vals = (doublereal *) G_calloc(ldim * cols, sizeof(doublereal));
127  A->is_init = 1;
128 
129  return 0;
130 }
131 
132 
144 mat_struct *G_matrix_copy(const mat_struct * A)
145 {
146  mat_struct *B;
147 
148  if (!A->is_init) {
149  G_warning(_("Matrix is not initialised fully."));
150  return NULL;
151  }
152 
153  if ((B = G_matrix_init(A->rows, A->cols, A->ldim)) == NULL) {
154  G_warning(_("Unable to allocate space for matrix copy"));
155  return NULL;
156  }
157 
158  memcpy(&B->vals[0], &A->vals[0], A->cols * A->ldim * sizeof(doublereal));
159 
160  return B;
161 }
162 
163 
177 mat_struct *G_matrix_add(mat_struct * mt1, mat_struct * mt2)
178 {
179  return G__matrix_add(mt1, mt2, 1, 1);
180 }
181 
182 
196 mat_struct *G_matrix_subtract(mat_struct * mt1, mat_struct * mt2)
197 {
198  return G__matrix_add(mt1, mt2, 1, -1);
199 }
200 
201 
215 mat_struct *G_matrix_scale(mat_struct * mt1, const double c)
216 {
217  return G__matrix_add(mt1, NULL, c, 0);
218 }
219 
220 
236 mat_struct *G__matrix_add(mat_struct * mt1, mat_struct * mt2, const double c1,
237  const double c2)
238 {
239  mat_struct *mt3;
240  int i, j; /* loop variables */
241 
242  if (c1 == 0) {
243  G_warning(_("First scalar multiplier must be non-zero"));
244  return NULL;
245  }
246 
247  if (c2 == 0) {
248  if (!mt1->is_init) {
249  G_warning(_("One or both input matrices uninitialised"));
250  return NULL;
251  }
252  }
253 
254  else {
255 
256  if (!((mt1->is_init) && (mt2->is_init))) {
257  G_warning(_("One or both input matrices uninitialised"));
258  return NULL;
259  }
260 
261  if (mt1->rows != mt2->rows || mt1->cols != mt2->cols) {
262  G_warning(_("Matrix order does not match"));
263  return NULL;
264  }
265  }
266 
267  if ((mt3 = G_matrix_init(mt1->rows, mt1->cols, mt1->ldim)) == NULL) {
268  G_warning(_("Unable to allocate space for matrix sum"));
269  return NULL;
270  }
271 
272  if (c2 == 0) {
273 
274  for (i = 0; i < mt3->rows; i++) {
275  for (j = 0; j < mt3->cols; j++) {
276  mt3->vals[i + mt3->ldim * j] =
277  c1 * mt1->vals[i + mt1->ldim * j];
278  }
279  }
280  }
281 
282  else {
283 
284  for (i = 0; i < mt3->rows; i++) {
285  for (j = 0; j < mt3->cols; j++) {
286  mt3->vals[i + mt3->ldim * j] =
287  c1 * mt1->vals[i + mt1->ldim * j] + c2 * mt2->vals[i +
288  mt2->
289  ldim *
290  j];
291  }
292  }
293  }
294 
295  return mt3;
296 }
297 
298 
299 #if defined(HAVE_LIBBLAS)
300 
314 mat_struct *G_matrix_product(mat_struct * mt1, mat_struct * mt2)
315 {
316  mat_struct *mt3;
317  doublereal unity = 1, zero = 0;
318  integer rows, cols, interdim, lda, ldb;
319  integer1 no_trans = 'n';
320 
321  if (!((mt1->is_init) || (mt2->is_init))) {
322  G_warning(_("One or both input matrices uninitialised"));
323  return NULL;
324  }
325 
326  if (mt1->cols != mt2->rows) {
327  G_warning(_("Matrix order does not match"));
328  return NULL;
329  }
330 
331  if ((mt3 = G_matrix_init(mt1->rows, mt2->cols, mt1->ldim)) == NULL) {
332  G_warning(_("Unable to allocate space for matrix product"));
333  return NULL;
334  }
335 
336  /* Call the driver */
337 
338  rows = (integer) mt1->rows;
339  interdim = (integer) mt1->cols;
340  cols = (integer) mt2->cols;
341 
342  lda = (integer) mt1->ldim;
343  ldb = (integer) mt2->ldim;
344 
345  f77_dgemm(&no_trans, &no_trans, &rows, &cols, &interdim, &unity,
346  mt1->vals, &lda, mt2->vals, &ldb, &zero, mt3->vals, &lda);
347 
348  return mt3;
349 }
350 
351 #else /* defined(HAVE_LIBBLAS) */
352 #warning G_matrix_product() not compiled; requires BLAS library
353 #endif /* defined(HAVE_LIBBLAS) */
354 
368 mat_struct *G_matrix_transpose(mat_struct * mt)
369 {
370  mat_struct *mt1;
371  int ldim, ldo;
372  doublereal *dbo, *dbt, *dbx, *dby;
373  int cnt, cnt2;
374 
375  /* Word align the workspace blocks */
376  if (mt->cols % 2 == 0)
377  ldim = mt->cols;
378  else
379  ldim = mt->cols + 1;
380 
381  mt1 = G_matrix_init(mt->cols, mt->rows, ldim);
382 
383  /* Set initial values for reading arrays */
384  dbo = &mt->vals[0];
385  dbt = &mt1->vals[0];
386  ldo = mt->ldim;
387 
388  for (cnt = 0; cnt < mt->cols; cnt++) {
389  dbx = dbo;
390  dby = dbt;
391 
392  for (cnt2 = 0; cnt2 < ldo - 1; cnt2++) {
393  *dby = *dbx;
394  dby += ldim;
395  dbx++;
396  }
397 
398  *dby = *dbx;
399 
400  if (cnt < mt->cols - 1) {
401  dbo += ldo;
402  dbt++;
403  }
404  }
405 
406  return mt1;
407 }
408 
409 
410 #if defined(HAVE_LIBBLAS) && defined(HAVE_LIBLAPACK)
411 
435 /*** NOT YET COMPLETE: only some solutions' options available ***/
436 
437 int
438 G_matrix_LU_solve(const mat_struct * mt1, mat_struct ** xmat0,
439  const mat_struct * bmat, mat_type mtype)
440 {
441  mat_struct *wmat, *xmat, *mtx;
442 
443  if (mt1->is_init == 0 || bmat->is_init == 0) {
444  G_warning(_("Input: one or both data matrices uninitialised"));
445  return -1;
446  }
447 
448  if (mt1->rows != mt1->cols || mt1->rows < 1) {
449  G_warning(_("Principal matrix is not properly dimensioned"));
450  return -1;
451  }
452 
453  if (bmat->cols < 1) {
454  G_warning(_("Input: you must have at least one array to solve"));
455  return -1;
456  }
457 
458  /* Now create solution matrix by copying the original coefficient matrix */
459  if ((xmat = G_matrix_copy(bmat)) == NULL) {
460  G_warning(_("Could not allocate space for solution matrix"));
461  return -1;
462  }
463 
464  /* Create working matrix for the coefficient array */
465  if ((mtx = G_matrix_copy(mt1)) == NULL) {
466  G_warning(_("Could not allocate space for working matrix"));
467  return -1;
468  }
469 
470  /* Copy the contents of the data matrix, to preserve the
471  original information
472  */
473  if ((wmat = G_matrix_copy(bmat)) == NULL) {
474  G_warning(_("Could not allocate space for working matrix"));
475  return -1;
476  }
477 
478  /* Now call appropriate LA driver to solve equations */
479  switch (mtype) {
480 
481  case NONSYM:
482  {
483  integer *perm, res_info;
484  integer num_eqns, nrhs, lda, ldb;
485 
486  perm = (integer *) G_malloc(wmat->rows);
487 
488  /* Set fields to pass to fortran routine */
489  num_eqns = (integer) mt1->rows;
490  nrhs = (integer) wmat->cols;
491  lda = (integer) mt1->ldim;
492  ldb = (integer) wmat->ldim;
493 
494  /* Call LA driver */
495  f77_dgesv(&num_eqns, &nrhs, mtx->vals, &lda, perm, wmat->vals,
496  &ldb, &res_info);
497 
498  /* Copy the results from the modified data matrix, taking account
499  of pivot permutations ???
500  */
501 
502  /*
503  for(indx1 = 0; indx1 < num_eqns; indx1++) {
504  iperm = perm[indx1];
505  ptin = &wmat->vals[0] + indx1;
506  ptout = &xmat->vals[0] + iperm;
507 
508  for(indx2 = 0; indx2 < nrhs - 1; indx2++) {
509  *ptout = *ptin;
510  ptin += wmat->ldim;
511  ptout += xmat->ldim;
512  }
513 
514  *ptout = *ptin;
515  }
516  */
517 
518  memcpy(xmat->vals, wmat->vals,
519  wmat->cols * wmat->ldim * sizeof(doublereal));
520 
521  /* Free temp arrays */
522  G_free(perm);
523  G_matrix_free(wmat);
524  G_matrix_free(mtx);
525 
526  if (res_info > 0) {
527  G_warning(_("Matrix (or submatrix is singular). Solution undetermined"));
528  return 1;
529  }
530  else if (res_info < 0) {
531  G_warning(_("Problem in LA routine."));
532  return -1;
533  }
534  break;
535  }
536  default:
537  {
538  G_warning(_("Procedure not yet available for selected matrix type"));
539  return -1;
540  }
541  } /* end switch */
542 
543  *xmat0 = xmat;
544 
545  return 0;
546 }
547 
548 #else /* defined(HAVE_LIBBLAS) && defined(HAVE_LIBLAPACK) */
549 #warning G_matrix_LU_solve() not compiled; requires BLAS and LAPACK libraries
550 #endif /* defined(HAVE_LIBBLAS) && defined(HAVE_LIBLAPACK) */
551 
552 #if defined(HAVE_LIBBLAS) && defined(HAVE_LIBLAPACK)
553 
566 mat_struct *G_matrix_inverse(mat_struct * mt)
567 {
568  mat_struct *mt0, *res;
569  int i, j, k; /* loop */
570 
571  if (mt->rows != mt->cols) {
572  G_warning(_("Matrix is not square. Cannot determine inverse"));
573  return NULL;
574  }
575 
576  if ((mt0 = G_matrix_init(mt->rows, mt->rows, mt->ldim)) == NULL) {
577  G_warning(_("Unable to allocate space for matrix"));
578  return NULL;
579  }
580 
581  /* Set `B' matrix to unit matrix */
582  for (i = 0; i < mt0->rows - 1; i++) {
583  mt0->vals[i + i * mt0->ldim] = 1.0;
584 
585  for (j = i + 1; j < mt0->cols; j++) {
586  mt0->vals[i + j * mt0->ldim] = mt0->vals[j + i * mt0->ldim] = 0.0;
587  }
588  }
589 
590  mt0->vals[mt0->rows - 1 + (mt0->rows - 1) * mt0->ldim] = 1.0;
591 
592  /* Solve system */
593  if ((k = G_matrix_LU_solve(mt, &res, mt0, NONSYM)) == 1) {
594  G_warning(_("Matrix is singular"));
595  G_matrix_free(mt0);
596  return NULL;
597  }
598  else if (k < 0) {
599  G_warning(_("Problem in LA procedure."));
600  G_matrix_free(mt0);
601  return NULL;
602  }
603  else {
604  G_matrix_free(mt0);
605  return res;
606  }
607 }
608 
609 #else /* defined(HAVE_LIBBLAS) && defined(HAVE_LIBLAPACK) */
610 #warning G_matrix_inverse() not compiled; requires BLAS and LAPACK libraries
611 #endif /* defined(HAVE_LIBBLAS) && defined(HAVE_LIBLAPACK) */
612 
613 
625 void G_matrix_free(mat_struct * mt)
626 {
627  if (mt->is_init)
628  G_free(mt->vals);
629 
630  G_free(mt);
631 }
632 
633 
645 void G_matrix_print(mat_struct * mt)
646 {
647  int i, j;
648  char buf[64], numbuf[64];
649 
650  for (i = 0; i < mt->rows; i++) {
651  strcpy(buf, "");
652 
653  for (j = 0; j < mt->cols; j++) {
654 
655  sprintf(numbuf, "%14.6f", G_matrix_get_element(mt, i, j));
656  strcat(buf, numbuf);
657  if (j < mt->cols - 1)
658  strcat(buf, ", ");
659  }
660 
661  G_message("%s", buf);
662  }
663 
664  fprintf(stderr, "\n");
665 }
666 
667 
684 int G_matrix_set_element(mat_struct * mt, int rowval, int colval, double val)
685 {
686  if (!mt->is_init) {
687  G_warning(_("Element array has not been allocated"));
688  return -1;
689  }
690 
691  if (rowval >= mt->rows || colval >= mt->cols || rowval < 0 || colval < 0) {
692  G_warning(_("Specified element is outside array bounds"));
693  return -1;
694  }
695 
696  mt->vals[rowval + colval * mt->ldim] = (doublereal) val;
697 
698  return 0;
699 }
700 
701 
717 double G_matrix_get_element(mat_struct * mt, int rowval, int colval)
718 {
719  double val;
720 
721  /* Should do some checks, but this would require an error control
722  system: later? */
723 
724  return (val = (double)mt->vals[rowval + colval * mt->ldim]);
725 }
726 
727 
740 vec_struct *G_matvect_get_column(mat_struct * mt, int col)
741 {
742  int i; /* loop */
743  vec_struct *vc1;
744 
745  if (col < 0 || col >= mt->cols) {
746  G_warning(_("Specified matrix column index is outside range"));
747  return NULL;
748  }
749 
750  if (!mt->is_init) {
751  G_warning(_("Matrix is not initialised"));
752  return NULL;
753  }
754 
755  if ((vc1 = G_vector_init(mt->rows, mt->ldim, CVEC)) == NULL) {
756  G_warning(_("Could not allocate space for vector structure"));
757  return NULL;
758  }
759 
760  for (i = 0; i < mt->rows; i++)
761  G_matrix_set_element((mat_struct *) vc1, i, 0,
762  G_matrix_get_element(mt, i, col));
763 
764  return vc1;
765 }
766 
767 
781 vec_struct *G_matvect_get_row(mat_struct * mt, int row)
782 {
783  int i; /* loop */
784  vec_struct *vc1;
785 
786  if (row < 0 || row >= mt->cols) {
787  G_warning(_("Specified matrix row index is outside range"));
788  return NULL;
789  }
790 
791  if (!mt->is_init) {
792  G_warning(_("Matrix is not initialised"));
793  return NULL;
794  }
795 
796  if ((vc1 = G_vector_init(mt->cols, mt->ldim, RVEC)) == NULL) {
797  G_warning(_("Could not allocate space for vector structure"));
798  return NULL;
799  }
800 
801  for (i = 0; i < mt->cols; i++)
802  G_matrix_set_element((mat_struct *) vc1, 0, i,
803  G_matrix_get_element(mt, row, i));
804 
805  return vc1;
806 }
807 
808 
824 int G_matvect_extract_vector(mat_struct * mt, vtype vt, int indx)
825 {
826  if (vt == RVEC && indx >= mt->rows) {
827  G_warning(_("Specified row index is outside range"));
828  return -1;
829  }
830 
831  else if (vt == CVEC && indx >= mt->cols) {
832  G_warning(_("Specified column index is outside range"));
833  return -1;
834  }
835 
836  switch (vt) {
837 
838  case RVEC:
839  {
840  mt->type = ROWVEC_;
841  mt->v_indx = indx;
842  }
843 
844  case CVEC:
845  {
846  mt->type = COLVEC_;
847  mt->v_indx = indx;
848  }
849 
850  default:
851  {
852  G_warning(_("Unknown vector type."));
853  return -1;
854  }
855 
856  }
857 
858  return 0;
859 
860 }
861 
862 
874 int G_matvect_retrieve_matrix(vec_struct * vc)
875 {
876  /* We have to take the integrity of the vector structure
877  largely on trust
878  */
879 
880  vc->type = MATRIX_;
881  vc->v_indx = -1;
882 
883  return 0;
884 }
885 
886 
901 vec_struct *G_vector_init(int cells, int ldim, vtype vt)
902 {
903  vec_struct *tmp_arry;
904 
905  if ((cells < 1) || (vt == RVEC && ldim < 1)
906  || (vt == CVEC && ldim < cells) || ldim < 0) {
907  G_warning("Vector dimensions out of range.");
908  return NULL;
909  }
910 
911  tmp_arry = (vec_struct *) G_malloc(sizeof(vec_struct));
912 
913  if (vt == RVEC) {
914  tmp_arry->rows = 1;
915  tmp_arry->cols = cells;
916  tmp_arry->ldim = ldim;
917  tmp_arry->type = ROWVEC_;
918  }
919 
920  else if (vt == CVEC) {
921  tmp_arry->rows = cells;
922  tmp_arry->cols = 1;
923  tmp_arry->ldim = ldim;
924  tmp_arry->type = COLVEC_;
925  }
926 
927  tmp_arry->v_indx = 0;
928 
929  tmp_arry->vals = (doublereal *) G_calloc(ldim * tmp_arry->cols,
930  sizeof(doublereal));
931  tmp_arry->is_init = 1;
932 
933  return tmp_arry;
934 }
935 
936 
948 void G_vector_free(vec_struct * v)
949 {
950  if (v->is_init)
951  G_free(v->vals);
952 
953  G_free(v);
954 }
955 
956 
971 vec_struct *G_vector_sub(vec_struct * v1, vec_struct * v2, vec_struct * out)
972 {
973  int idx1, idx2, idx0;
974  int i;
975 
976  if (!out->is_init) {
977  G_warning(_("Output vector is uninitialized"));
978  return NULL;
979  }
980 
981  if (v1->type != v2->type) {
982  G_warning(_("Vectors are not of the same type"));
983  return NULL;
984  }
985 
986  if (v1->type != out->type) {
987  G_warning(_("Output vector is of incorrect type"));
988  return NULL;
989  }
990 
991  if (v1->type == MATRIX_) {
992  G_warning(_("Matrices not allowed"));
993  return NULL;
994  }
995 
996  if ((v1->type == ROWVEC_ && v1->cols != v2->cols) ||
997  (v1->type == COLVEC_ && v1->rows != v2->rows)) {
998  G_warning(_("Vectors have differing dimensions"));
999  return NULL;
1000  }
1001 
1002  if ((v1->type == ROWVEC_ && v1->cols != out->cols) ||
1003  (v1->type == COLVEC_ && v1->rows != out->rows)) {
1004  G_warning(_("Output vector has incorrect dimension"));
1005  return NULL;
1006  }
1007 
1008  idx1 = (v1->v_indx > 0) ? v1->v_indx : 0;
1009  idx2 = (v2->v_indx > 0) ? v2->v_indx : 0;
1010  idx0 = (out->v_indx > 0) ? out->v_indx : 0;
1011 
1012  if (v1->type == ROWVEC_) {
1013  for (i = 0; i < v1->cols; i++)
1014  G_matrix_set_element(out, idx0, i,
1015  G_matrix_get_element(v1, idx1, i) -
1016  G_matrix_get_element(v2, idx2, i));
1017  }
1018  else {
1019  for (i = 0; i < v1->rows; i++)
1020  G_matrix_set_element(out, i, idx0,
1021  G_matrix_get_element(v1, i, idx1) -
1022  G_matrix_get_element(v2, i, idx2));
1023  }
1024 
1025  return out;
1026 }
1027 
1028 
1046 int G_vector_set(vec_struct * A, int cells, int ldim, vtype vt, int vindx)
1047 {
1048  if ((cells < 1) || (vt == RVEC && ldim < 1)
1049  || (vt == CVEC && ldim < cells) || ldim < 0) {
1050  G_warning(_("Vector dimensions out of range"));
1051  return -1;
1052  }
1053 
1054  if ((vt == RVEC && vindx >= A->cols) || (vt == CVEC && vindx >= A->rows)) {
1055  G_warning(_("Row/column out of range"));
1056  return -1;
1057  }
1058 
1059  if (vt == RVEC) {
1060  A->rows = 1;
1061  A->cols = cells;
1062  A->ldim = ldim;
1063  A->type = ROWVEC_;
1064  }
1065  else {
1066  A->rows = cells;
1067  A->cols = 1;
1068  A->ldim = ldim;
1069  A->type = COLVEC_;
1070  }
1071 
1072  if (vindx < 0)
1073  A->v_indx = 0;
1074  else
1075  A->v_indx = vindx;
1076 
1077  A->vals = (doublereal *) G_calloc(ldim * A->cols, sizeof(doublereal));
1078  A->is_init = 1;
1079 
1080  return 0;
1081 
1082 }
1083 
1084 
1085 #if defined(HAVE_LIBBLAS)
1086 
1099 double G_vector_norm_euclid(vec_struct * vc)
1100 {
1101  integer incr, Nval;
1102  doublereal *startpt;
1103 
1104  if (!vc->is_init)
1105  G_fatal_error(_("Matrix is not initialised"));
1106 
1107  if (vc->type == ROWVEC_) {
1108  Nval = (integer) vc->cols;
1109  incr = (integer) vc->ldim;
1110  if (vc->v_indx < 0)
1111  startpt = vc->vals;
1112  else
1113  startpt = vc->vals + vc->v_indx;
1114  }
1115  else {
1116  Nval = (integer) vc->rows;
1117  incr = 1;
1118  if (vc->v_indx < 0)
1119  startpt = vc->vals;
1120  else
1121  startpt = vc->vals + vc->v_indx * vc->ldim;
1122  }
1123 
1124  /* Call the BLAS routine dnrm2_() */
1125  return (double)f77_dnrm2(&Nval, startpt, &incr);
1126 }
1127 
1128 #else /* defined(HAVE_LIBBLAS) */
1129 #warning G_vector_norm_euclid() not compiled; requires BLAS library
1130 #endif /* defined(HAVE_LIBBLAS) */
1131 
1132 
1150 double G_vector_norm_maxval(vec_struct * vc, int vflag)
1151 {
1152  doublereal xval, *startpt, *curpt;
1153  double cellval;
1154  int ncells, incr;
1155 
1156  if (!vc->is_init)
1157  G_fatal_error(_("Matrix is not initialised"));
1158 
1159  if (vc->type == ROWVEC_) {
1160  ncells = (integer) vc->cols;
1161  incr = (integer) vc->ldim;
1162  if (vc->v_indx < 0)
1163  startpt = vc->vals;
1164  else
1165  startpt = vc->vals + vc->v_indx;
1166  }
1167  else {
1168  ncells = (integer) vc->rows;
1169  incr = 1;
1170  if (vc->v_indx < 0)
1171  startpt = vc->vals;
1172  else
1173  startpt = vc->vals + vc->v_indx * vc->ldim;
1174  }
1175 
1176  xval = *startpt;
1177  curpt = startpt;
1178 
1179  while (ncells > 0) {
1180  if (curpt != startpt) {
1181  switch (vflag) {
1182 
1183  case MAX_POS:
1184  {
1185  if (*curpt > xval)
1186  xval = *curpt;
1187  break;
1188  }
1189 
1190  case MAX_NEG:
1191  {
1192  if (*curpt < xval)
1193  xval = *curpt;
1194  break;
1195  }
1196 
1197  case MAX_ABS:
1198  {
1199  cellval = (double)(*curpt);
1200  if (hypot(cellval, cellval) > (double)xval)
1201  xval = *curpt;
1202  }
1203  } /* switch */
1204  } /* if(curpt != startpt) */
1205 
1206  curpt += incr;
1207  ncells--;
1208  }
1209 
1210  return (double)xval;
1211 }
1212 
1213 
1225 double G_vector_norm1(vec_struct * vc)
1226 {
1227  double result = 0.0;
1228  int idx;
1229  int i;
1230 
1231  if (!vc->is_init) {
1232  G_warning(_("Matrix is not initialised"));
1233  return 0.0 / 0.0; /* NaN */
1234  }
1235 
1236  idx = (vc->v_indx > 0) ? vc->v_indx : 0;
1237 
1238  if (vc->type == ROWVEC_) {
1239  for (i = 0; i < vc->cols; i++)
1240  result += fabs(G_matrix_get_element(vc, idx, i));
1241  }
1242  else {
1243  for (i = 0; i < vc->rows; i++)
1244  result += fabs(G_matrix_get_element(vc, i, idx));
1245  }
1246 
1247  return result;
1248 }
1249 
1250 
1262 vec_struct *G_vector_copy(const vec_struct * vc1, int comp_flag)
1263 {
1264  vec_struct *tmp_arry;
1265  int incr1, incr2;
1266  doublereal *startpt1, *startpt2, *curpt1, *curpt2;
1267  int cnt;
1268 
1269  if (!vc1->is_init) {
1270  G_warning(_("Vector structure is not initialised"));
1271  return NULL;
1272  }
1273 
1274  tmp_arry = (vec_struct *) G_malloc(sizeof(vec_struct));
1275 
1276  if (comp_flag == DO_COMPACT) {
1277  if (vc1->type == ROWVEC_) {
1278  tmp_arry->rows = 1;
1279  tmp_arry->cols = vc1->cols;
1280  tmp_arry->ldim = 1;
1281  tmp_arry->type = ROWVEC_;
1282  tmp_arry->v_indx = 0;
1283  }
1284  else if (vc1->type == COLVEC_) {
1285  tmp_arry->rows = vc1->rows;
1286  tmp_arry->cols = 1;
1287  tmp_arry->ldim = vc1->ldim;
1288  tmp_arry->type = COLVEC_;
1289  tmp_arry->v_indx = 0;
1290  }
1291  else {
1292  G_warning("Type is not vector.");
1293  return NULL;
1294  }
1295  }
1296  else if (comp_flag == NO_COMPACT) {
1297  tmp_arry->v_indx = vc1->v_indx;
1298  tmp_arry->rows = vc1->rows;
1299  tmp_arry->cols = vc1->cols;
1300  tmp_arry->ldim = vc1->ldim;
1301  tmp_arry->type = vc1->type;
1302  }
1303  else {
1304  G_warning("Copy method must be specified: [DO,NO]_COMPACT.\n");
1305  return NULL;
1306  }
1307 
1308  tmp_arry->vals = (doublereal *) G_calloc(tmp_arry->ldim * tmp_arry->cols,
1309  sizeof(doublereal));
1310  if (comp_flag == DO_COMPACT) {
1311  if (tmp_arry->type == ROWVEC_) {
1312  startpt1 = tmp_arry->vals;
1313  startpt2 = vc1->vals + vc1->v_indx;
1314  curpt1 = startpt1;
1315  curpt2 = startpt2;
1316  incr1 = 1;
1317  incr2 = vc1->ldim;
1318  cnt = vc1->cols;
1319  }
1320  else if (tmp_arry->type == COLVEC_) {
1321  startpt1 = tmp_arry->vals;
1322  startpt2 = vc1->vals + vc1->v_indx * vc1->ldim;
1323  curpt1 = startpt1;
1324  curpt2 = startpt2;
1325  incr1 = 1;
1326  incr2 = 1;
1327  cnt = vc1->rows;
1328  }
1329  else {
1330  G_warning("Structure type is not vector.");
1331  return NULL;
1332  }
1333  }
1334  else if (comp_flag == NO_COMPACT) {
1335  startpt1 = tmp_arry->vals;
1336  startpt2 = vc1->vals;
1337  curpt1 = startpt1;
1338  curpt2 = startpt2;
1339  incr1 = 1;
1340  incr2 = 1;
1341  cnt = vc1->ldim * vc1->cols;
1342  }
1343  else {
1344  G_warning("Copy method must be specified: [DO,NO]_COMPACT.\n");
1345  return NULL;
1346  }
1347 
1348  while (cnt > 0) {
1349  memcpy(curpt1, curpt2, sizeof(doublereal));
1350  curpt1 += incr1;
1351  curpt2 += incr2;
1352  cnt--;
1353  }
1354 
1355  tmp_arry->is_init = 1;
1356 
1357  return tmp_arry;
1358 }
1359 
1360 
1375 int G_matrix_read(FILE * fp, mat_struct * out)
1376 {
1377  char buff[100];
1378  int rows, cols;
1379  int i, j, row;
1380  double val;
1381 
1382  /* skip comments */
1383  for (;;) {
1384  if (!G_getl(buff, sizeof(buff), fp))
1385  return -1;
1386  if (buff[0] != '#')
1387  break;
1388  }
1389 
1390  if (sscanf(buff, "Matrix: %d by %d", &rows, &cols) != 2) {
1391  G_warning(_("Input format error"));
1392  return -1;
1393  }
1394 
1395  G_matrix_set(out, rows, cols, rows);
1396 
1397  for (i = 0; i < rows; i++) {
1398  if (fscanf(fp, "row%d:", &row) != 1 || row != i) {
1399  G_warning(_("Input format error"));
1400  return -1;
1401  }
1402  for (j = 0; j < cols; j++) {
1403  if (fscanf(fp, "%lf:", &val) != 1) {
1404  G_warning(_("Input format error"));
1405  return -1;
1406  }
1407 
1408  G_matrix_set_element(out, i, j, val);
1409  }
1410  }
1411 
1412  return 0;
1413 }
1414 
1415 
1429 int G_matrix_stdin(mat_struct * out)
1430 {
1431  return G_matrix_read(stdin, out);
1432 }
1433 
1434 
1447 int G_matrix_eigen_sort(vec_struct * d, mat_struct * m)
1448 {
1449  mat_struct tmp;
1450  int i, j;
1451  int idx;
1452 
1453  G_matrix_set(&tmp, m->rows + 1, m->cols, m->ldim + 1);
1454 
1455  idx = (d->v_indx > 0) ? d->v_indx : 0;
1456 
1457  /* concatenate (vertically) m and d into tmp */
1458  for (i = 0; i < m->cols; i++) {
1459  for (j = 0; j < m->rows; j++)
1460  G_matrix_set_element(&tmp, j + 1, i,
1461  G_matrix_get_element(m, j, i));
1462  if (d->type == ROWVEC_)
1463  G_matrix_set_element(&tmp, 0, i, G_matrix_get_element(d, idx, i));
1464  else
1465  G_matrix_set_element(&tmp, 0, i, G_matrix_get_element(d, i, idx));
1466  }
1467 
1468  /* sort the combined matrix */
1469  qsort(tmp.vals, tmp.cols, tmp.ldim * sizeof(doublereal), egcmp);
1470 
1471  /* split tmp into m and d */
1472  for (i = 0; i < m->cols; i++) {
1473  for (j = 0; j < m->rows; j++)
1474  G_matrix_set_element(m, j, i,
1475  G_matrix_get_element(&tmp, j + 1, i));
1476  if (d->type == ROWVEC_)
1477  G_matrix_set_element(d, idx, i, G_matrix_get_element(&tmp, 0, i));
1478  else
1479  G_matrix_set_element(d, i, idx, G_matrix_get_element(&tmp, 0, i));
1480  }
1481 
1482  G_free(tmp.vals);
1483 
1484  return 0;
1485 }
1486 
1487 
1488 static int egcmp(const void *pa, const void *pb)
1489 {
1490  double a = *(doublereal * const)pa;
1491  double b = *(doublereal * const)pb;
1492 
1493  if (a > b)
1494  return 1;
1495  if (a < b)
1496  return -1;
1497 
1498  return 0;
1499 }
1500 
1501 
1502 #endif /* HAVE_BLAS && HAVE_LAPACK && HAVE_G2C */
sprintf(buf2,"%s", G3D_CATS_ELEMENT)
void G_free(void *buf)
Free allocated memory.
Definition: gis/alloc.c:142
float b
Definition: named_colr.c:8
char buff[1024]
Definition: g3dcats.c:89
void G_message(const char *msg,...)
Print a message to stderr.
Definition: lib/gis/error.c:74
if(!YY_CURRENT_BUFFER)
Definition: lex.yy.c:799
char buf[GNAME_MAX+sizeof(G3D_DIRECTORY)+2]
Definition: g3drange.c:62
int G_getl(char *buf, int n, FILE *fd)
gets a line of text from a file
Definition: getl.c:17
return NULL
Definition: dbfopen.c:1394
G_warning("category support for [%s] in mapset [%s] %s", name, mapset, type)
tuple cols
for(cat=0;;cat++)
Definition: g3dcats.c:140
int G_fatal_error(const char *msg,...)
Print a fatal error message to stderr.
string vtype
Definition: tools.py:3428